CMS 3D CMS Logo

Public Member Functions | Private Attributes

pf2pat::ElectronIDPFCandidateSelectorDefinition Class Reference

Selects PFCandidates basing on cuts provided with string cut parser. More...

#include <CommonTools/ParticleFlow/interface/ElectronIDPFCandidateSelectorDefinition.h>

Inheritance diagram for pf2pat::ElectronIDPFCandidateSelectorDefinition:
pf2pat::PFCandidateSelectorDefinition

List of all members.

Public Member Functions

 ElectronIDPFCandidateSelectorDefinition (const edm::ParameterSet &cfg)
void select (const HandleToCollection &hc, const edm::Event &e, const edm::EventSetup &s)

Private Attributes

edm::InputTag electronId_
edm::InputTag electrons_
bool isBitMap_
uint32_t mask_
double value_

Detailed Description

Selects PFCandidates basing on cuts provided with string cut parser.

Author:
Giovanni Petrucciani
Version:
Id:
ElectronIDPFCandidateSelectorDefinition.h,v 1.2 2011/04/06 12:12:38 rwolf Exp

Definition at line 21 of file ElectronIDPFCandidateSelectorDefinition.h.


Constructor & Destructor Documentation

pf2pat::ElectronIDPFCandidateSelectorDefinition::ElectronIDPFCandidateSelectorDefinition ( const edm::ParameterSet cfg) [inline]

Definition at line 23 of file ElectronIDPFCandidateSelectorDefinition.h.

References Exception, edm::ParameterSet::exists(), edm::ParameterSet::existsAs(), edm::ParameterSet::getParameter(), isBitMap_, mask_, and value_.

                                                                            :
      electrons_(  cfg.getParameter< edm::InputTag >( "recoGsfElectrons" ) ), 
      electronId_( cfg.getParameter< edm::InputTag >( "electronIdMap" ) )
    { 
        if (cfg.exists("bitsToCheck")) {
            isBitMap_ = true; 
            mask_ = 0;
            if (cfg.existsAs<std::vector<std::string> >("bitsToCheck")) {
                std::vector<std::string> strbits = cfg.getParameter<std::vector<std::string> >("bitsToCheck");
                for (std::vector<std::string>::const_iterator istrbit = strbits.begin(), estrbit = strbits.end(); 
                        istrbit != estrbit; ++istrbit) {
                    if      (*istrbit == "id" )  { mask_ |= 1; }
                    else if (*istrbit == "iso")  { mask_ |= 2; }
                    else if (*istrbit == "conv") { mask_ |= 4; }
                    else if (*istrbit == "ip")   { mask_ |= 8; }
                    else throw cms::Exception("Configuration") << "ElectronIDPFCandidateSelector: " <<
                        "bitsToCheck allowed string values are only id(0), iso(1), conv(2), ip(3).\n" << 
                            "Otherwise, use uint32_t bitmask).\n";
                }
            } else if (cfg.existsAs<uint32_t>("bitsToCheck")) {
                mask_ = cfg.getParameter<uint32_t>("bitsToCheck");
            } else {
                throw cms::Exception("Configuration") << "ElectronIDPFCandidateSelector: " <<
                        "bitsToCheck must be either a vector of strings, or a uint32 bitmask.\n";   
            }
        } else {
            isBitMap_ = false;
            value_ = cfg.getParameter<double>("electronIdCut");
        }
    }

Member Function Documentation

void pf2pat::ElectronIDPFCandidateSelectorDefinition::select ( const HandleToCollection hc,
const edm::Event e,
const edm::EventSetup s 
) [inline]

Definition at line 54 of file ElectronIDPFCandidateSelectorDefinition.h.

References prof2calltree::count, electronId_, HI_PhotonSkim_cff::electrons, electrons_, edm::Event::getByLabel(), isBitMap_, edm::Ref< C, T, F >::isNull(), combine::key, mask_, match(), pf2pat::PFCandidateSelectorDefinition::selected_, and value_.

                                         {
      selected_.clear();

      edm::Handle<reco::GsfElectronCollection> electrons;
      e.getByLabel(electrons_, electrons);

      edm::Handle<edm::ValueMap<float> > electronId;
      e.getByLabel(electronId_, electronId);

      unsigned key=0;
      for( collection::const_iterator pfc = hc->begin(); 
           pfc != hc->end(); ++pfc, ++key) {

        // Get GsfTrack for matching with reco::GsfElectron objects
        reco::GsfTrackRef PfTk = pfc->gsfTrackRef();

        // skip ones without GsfTrack: they won't be matched anyway     
        if (PfTk.isNull()) continue;

        int match = -1;
        // try first the non-ambiguous tracks
        for (reco::GsfElectronCollection::const_iterator it = electrons->begin(), ed = electrons->end(); it != ed; ++it) {
            if (it->gsfTrack() == PfTk) { match = it - electrons->begin(); break; }
        }
        // then the ambiguous ones
        if (match == -1) {
            for (reco::GsfElectronCollection::const_iterator it = electrons->begin(), ed = electrons->end(); it != ed; ++it) {
                if (std::count(it->ambiguousGsfTracksBegin(), it->ambiguousGsfTracksEnd(), PfTk) > 0) {
                    match = it - electrons->begin(); break;
                }
            }
        }
        // if found, make a GsfElectronRef and read electron id
        if (match != -1) {
            reco::GsfElectronRef ref(electrons,match);
            float eleId = (*electronId)[ref];
            bool pass = false;
            if (isBitMap_) {
                uint32_t thisval = eleId; 
                pass = ((thisval & mask_) == mask_);
            } else {
                pass = (eleId > value_);
            }
            if (pass) {
                selected_.push_back( reco::PFCandidate(*pfc) );
                reco::PFCandidatePtr ptrToMother( hc, key );
                selected_.back().setSourceCandidatePtr( ptrToMother );
            }
        }
      }
    }

Member Data Documentation

Definition at line 110 of file ElectronIDPFCandidateSelectorDefinition.h.

Referenced by select().

Definition at line 109 of file ElectronIDPFCandidateSelectorDefinition.h.

Referenced by select().