Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:04:24

0001 //
0002 // template for implementation of jet associations
0003 // F Ratnikov, UMd, Sep. 5, 2007
0004 //
0005 
0006 #include "FWCore/Utilities/interface/Exception.h"
0007 #include "DataFormats/JetReco/interface/JetCollection.h"
0008 
0009 namespace JetAssociationTemplate {
0010 
0011   template <typename Container>
0012   inline unsigned findJet(const Container& fContainer, const reco::Jet& fJet) {
0013     for (unsigned i = 0; i < fContainer.size(); ++i) {
0014       if (&*(fContainer[i].first) == &fJet)
0015         return i;
0016     }
0017     throw cms::Exception("Invalid jet") << "JetAssociationTemplate-> inquire association with jet which is not "
0018                                            "available in the original jet collection";
0019   }
0020 
0021   template <typename Container, typename Value>
0022   inline bool setValue(Container* fContainer, const reco::JetBaseRef& fJet, const Value& fValue) {
0023     if (!fContainer)
0024       return false;
0025     fContainer->setValue(fJet.key(), fValue);
0026     return true;
0027   }
0028 
0029   template <typename Container, typename Value>
0030   inline bool setValue(Container& fContainer, const reco::JetBaseRef& fJet, const Value& fValue) {
0031     return setValue(&fContainer, fJet, fValue);
0032   }
0033 
0034   template <typename Container, typename Value>
0035   inline const Value& getValue(const Container& fContainer, const reco::JetBaseRef& fJet) {
0036     return fContainer[fJet];
0037   }
0038 
0039   template <typename Container, typename Value>
0040   inline const Value& getValue(const Container& fContainer, const reco::Jet& fJet) {
0041     unsigned i = findJet(fContainer, fJet);
0042     return fContainer[i].second;
0043   }
0044 
0045   template <typename Container>
0046   inline std::vector<reco::JetBaseRef> allJets(const Container& fContainer) {
0047     std::vector<reco::JetBaseRef> result;
0048     for (unsigned i = 0; i < fContainer.size(); ++i) {
0049       result.push_back(fContainer[i].first);
0050     }
0051     return result;
0052   }
0053 
0054   template <typename Container>
0055   inline bool hasJet(const Container& fContainer, const reco::JetBaseRef& fJet) {
0056     for (unsigned i = 0; i < fContainer.size(); ++i) {
0057       if (fContainer[i].first == fJet)
0058         return true;
0059     }
0060     return false;
0061   }
0062 
0063   template <typename Container>
0064   inline bool hasJet(const Container& fContainer, const reco::Jet& fJet) {
0065     for (unsigned i = 0; i < fContainer.size(); ++i) {
0066       if (&*(fContainer[i].first) == &fJet)
0067         return true;
0068     }
0069     return false;
0070   }
0071 }  // namespace JetAssociationTemplate