Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:01:05

0001 #ifndef CommonTools_ParticleFlow_ElectronIDPFCandidateSelectorDefinition
0002 #define CommonTools_ParticleFlow_ElectronIDPFCandidateSelectorDefinition
0003 
0004 /**
0005    \class    pf2pat::ElectronIDPFCandidateSelectorDefinition ElectronIDPFCandidateSelectorDefinition.h "CommonTools/ParticleFlow/interface/ElectronIDPFCandidateSelectorDefinition.h"
0006    \brief    Selects PFCandidates basing on cuts provided with string cut parser
0007 
0008    \author   Giovanni Petrucciani
0009    \version  $Id: ElectronIDPFCandidateSelectorDefinition.h,v 1.1 2011/01/28 20:56:44 srappocc Exp $
0010 */
0011 
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/Framework/interface/ConsumesCollector.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
0016 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0017 #include "DataFormats/Common/interface/ValueMap.h"
0018 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0019 #include "CommonTools/ParticleFlow/interface/PFCandidateSelectorDefinition.h"
0020 #include <algorithm>
0021 
0022 namespace pf2pat {
0023 
0024   struct ElectronIDPFCandidateSelectorDefinition : public PFCandidateSelectorDefinition {
0025     ElectronIDPFCandidateSelectorDefinition(const edm::ParameterSet& cfg, edm::ConsumesCollector&& iC)
0026         : electronsToken_(
0027               iC.consumes<reco::GsfElectronCollection>(cfg.getParameter<edm::InputTag>("recoGsfElectrons"))),
0028           electronIdToken_(iC.consumes<edm::ValueMap<float> >(cfg.getParameter<edm::InputTag>("electronIdMap"))) {
0029       if (cfg.exists("bitsToCheck")) {
0030         isBitMap_ = true;
0031         mask_ = 0;
0032         if (cfg.existsAs<std::vector<std::string> >("bitsToCheck")) {
0033           std::vector<std::string> strbits = cfg.getParameter<std::vector<std::string> >("bitsToCheck");
0034           for (std::vector<std::string>::const_iterator istrbit = strbits.begin(), estrbit = strbits.end();
0035                istrbit != estrbit;
0036                ++istrbit) {
0037             if (*istrbit == "id") {
0038               mask_ |= 1;
0039             } else if (*istrbit == "iso") {
0040               mask_ |= 2;
0041             } else if (*istrbit == "conv") {
0042               mask_ |= 4;
0043             } else if (*istrbit == "ip") {
0044               mask_ |= 8;
0045             } else
0046               throw cms::Exception("Configuration")
0047                   << "ElectronIDPFCandidateSelector: "
0048                   << "bitsToCheck allowed string values are only id(0), iso(1), conv(2), ip(3).\n"
0049                   << "Otherwise, use uint32_t bitmask).\n";
0050           }
0051         } else if (cfg.existsAs<uint32_t>("bitsToCheck")) {
0052           mask_ = cfg.getParameter<uint32_t>("bitsToCheck");
0053         } else {
0054           throw cms::Exception("Configuration")
0055               << "ElectronIDPFCandidateSelector: "
0056               << "bitsToCheck must be either a vector of strings, or a uint32 bitmask.\n";
0057         }
0058       } else {
0059         isBitMap_ = false;
0060         value_ = cfg.getParameter<double>("electronIdCut");
0061       }
0062     }
0063 
0064     void select(const HandleToCollection& hc, const edm::Event& e, const edm::EventSetup& s) {
0065       selected_.clear();
0066 
0067       edm::Handle<reco::GsfElectronCollection> electrons;
0068       e.getByToken(electronsToken_, electrons);
0069 
0070       edm::Handle<edm::ValueMap<float> > electronId;
0071       e.getByToken(electronIdToken_, electronId);
0072 
0073       unsigned key = 0;
0074       for (collection::const_iterator pfc = hc->begin(); pfc != hc->end(); ++pfc, ++key) {
0075         // Get GsfTrack for matching with reco::GsfElectron objects
0076         reco::GsfTrackRef PfTk = pfc->gsfTrackRef();
0077 
0078         // skip ones without GsfTrack: they won't be matched anyway
0079         if (PfTk.isNull())
0080           continue;
0081 
0082         int match = -1;
0083         // try first the non-ambiguous tracks
0084         for (auto it = electrons->begin(), ed = electrons->end(); it != ed; ++it) {
0085           if (it->gsfTrack() == PfTk) {
0086             match = it - electrons->begin();
0087             break;
0088           }
0089         }
0090         // then the ambiguous ones
0091         if (match == -1) {
0092           for (auto it = electrons->begin(), ed = electrons->end(); it != ed; ++it) {
0093             if (std::count(it->ambiguousGsfTracks().begin(), it->ambiguousGsfTracks().end(), PfTk) > 0) {
0094               match = it - electrons->begin();
0095               break;
0096             }
0097           }
0098         }
0099         // if found, make a GsfElectronRef and read electron id
0100         if (match != -1) {
0101           reco::GsfElectronRef ref(electrons, match);
0102           float eleId = (*electronId)[ref];
0103           bool pass = false;
0104           if (isBitMap_) {
0105             uint32_t thisval = eleId;
0106             pass = ((thisval & mask_) == mask_);
0107           } else {
0108             pass = (eleId > value_);
0109           }
0110           if (pass) {
0111             selected_.push_back(reco::PFCandidate(*pfc));
0112             reco::PFCandidatePtr ptrToMother(hc, key);
0113             selected_.back().setSourceCandidatePtr(ptrToMother);
0114           }
0115         }
0116       }
0117     }
0118 
0119   private:
0120     edm::EDGetTokenT<reco::GsfElectronCollection> electronsToken_;
0121     edm::EDGetTokenT<edm::ValueMap<float> > electronIdToken_;
0122     bool isBitMap_;
0123     uint32_t mask_;
0124     double value_;
0125   };
0126 }  // namespace pf2pat
0127 
0128 #endif