CMS 3D CMS Logo

pat::PATElectronProducer Class Reference

Produces pat::Electron's. More...

#include <PhysicsTools/PatAlgos/interface/PATElectronProducer.h>

Inheritance diagram for pat::PATElectronProducer:

edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Types

typedef edm::RefToBase
< ElectronType
ElectronBaseRef

Public Member Functions

 PATElectronProducer (const edm::ParameterSet &iConfig)
virtual void produce (edm::Event &iEvent, const edm::EventSetup &iSetup)
 ~PATElectronProducer ()

Private Types

typedef std::vector
< edm::Handle
< edm::Association
< reco::GenParticleCollection > > > 
GenAssociations
typedef std::pair< std::string,
edm::InputTag
NameTag
typedef std::vector
< edm::Handle
< edm::Association
< TriggerPrimitiveCollection > > > 
TrigAssociations

Private Member Functions

void FillElectron (Electron &aEl, const ElectronBaseRef &elecRef, const reco::CandidateBaseRef &baseRef, const GenAssociations &genMatches, const TrigAssociations &trigMatches) const

Private Attributes

bool addEfficiencies_
bool addElecID_
bool addElecShapes_
bool addGenMatch_
bool addResolutions_
bool addTrigMatch_
pat::helper::EfficiencyLoader efficiencyLoader_
std::vector< NameTagelecIDSrcs_
edm::InputTag electronSrc_
bool embedGenMatch_
bool embedGsfTrack_
bool embedPFCandidate_
bool embedSuperCluster_
bool embedTrack_
std::vector< edm::InputTaggenMatchSrc_
std::vector< std::pair
< pat::IsolationKeys,
edm::InputTag > > 
isoDepositLabels_
pat::helper::MultiIsolator isolator_
pat::helper::MultiIsolator::IsolationValuePairs isolatorTmpStorage_
std::auto_ptr
< EcalClusterLazyTools
lazyTools_
edm::InputTag pfElecSrc_
GreaterByPt< ElectronpTComparator_
edm::InputTag reducedBarrelRecHitCollection_
edm::InputTag reducedEndcapRecHitCollection_
std::vector< edm::InputTagtrigMatchSrc_
bool useParticleFlow_
 pflow specific
pat::PATUserDataHelper
< pat::Electron
userDataHelper_
bool useUserData_


Detailed Description

Produces pat::Electron's.

The PATElectronProducer produces analysis-level pat::Electron's starting from a collection of objects of ElectronType.

Author:
Steven Lowette, James Lamb
Version:
Id
PATElectronProducer.h,v 1.12.2.4 2009/02/06 23:54:55 pioppi Exp

Definition at line 51 of file PATElectronProducer.h.


Member Typedef Documentation

typedef edm::RefToBase<ElectronType> pat::PATElectronProducer::ElectronBaseRef

Definition at line 59 of file PATElectronProducer.h.

typedef std::vector<edm::Handle<edm::Association<reco::GenParticleCollection> > > pat::PATElectronProducer::GenAssociations [private]

Definition at line 81 of file PATElectronProducer.h.

typedef std::pair<std::string, edm::InputTag> pat::PATElectronProducer::NameTag [private]

Definition at line 92 of file PATElectronProducer.h.

typedef std::vector<edm::Handle<edm::Association<TriggerPrimitiveCollection> > > pat::PATElectronProducer::TrigAssociations [private]

Definition at line 83 of file PATElectronProducer.h.


Constructor & Destructor Documentation

PATElectronProducer::PATElectronProducer ( const edm::ParameterSet iConfig  )  [explicit]

Definition at line 28 of file PATElectronProducer.cc.

References addEfficiencies_, addElecID_, addElecShapes_, addGenMatch_, addResolutions_, addTrigMatch_, pat::ChargedParticleIso, pat::ECalIso, efficiencyLoader_, elecIDSrcs_, electronSrc_, embedGenMatch_, embedGsfTrack_, embedPFCandidate_, embedSuperCluster_, embedTrack_, edm::ParameterSet::exists(), edm::ParameterSet::existsAs(), pat::GammaParticleIso, genMatchSrc_, edm::ParameterSet::getParameter(), edm::ParameterSet::getParameterNamesForType(), pat::HCalIso, isoDepositLabels_, it, getDQMSummary::key, names, pat::NeutralParticleIso, pat::ParticleIso, pfElecSrc_, reducedBarrelRecHitCollection_, reducedEndcapRecHitCollection_, pat::TrackerIso, trigMatchSrc_, useParticleFlow_, pat::UserBaseIso, and useUserData_.

00028                                                                         :
00029   isolator_(iConfig.exists("isolation") ? iConfig.getParameter<edm::ParameterSet>("isolation") : edm::ParameterSet(), false) ,
00030   userDataHelper_ ( iConfig.getParameter<edm::ParameterSet>("userData") )
00031 {
00032 
00033   // general configurables
00034   electronSrc_      = iConfig.getParameter<edm::InputTag>( "electronSource" );
00035   embedGsfTrack_    = iConfig.getParameter<bool>         ( "embedGsfTrack" );
00036   embedSuperCluster_= iConfig.getParameter<bool>         ( "embedSuperCluster" );
00037   embedTrack_       = iConfig.getParameter<bool>         ( "embedTrack" );
00038 
00039   // pflow specific
00040   pfElecSrc_           = iConfig.getParameter<edm::InputTag>( "pfElectronSource" );
00041   useParticleFlow_        = iConfig.getParameter<bool>( "useParticleFlow" );
00042   embedPFCandidate_   = iConfig.getParameter<bool>( "embedPFCandidate" );
00043 
00044   
00045   // MC matching configurables
00046   addGenMatch_      = iConfig.getParameter<bool>          ( "addGenMatch" );
00047   if (addGenMatch_) {
00048       embedGenMatch_ = iConfig.getParameter<bool>         ( "embedGenMatch" );
00049       if (iConfig.existsAs<edm::InputTag>("genParticleMatch")) {
00050           genMatchSrc_.push_back(iConfig.getParameter<edm::InputTag>( "genParticleMatch" ));
00051       } else {
00052           genMatchSrc_ = iConfig.getParameter<std::vector<edm::InputTag> >( "genParticleMatch" );
00053       }
00054   }
00055   
00056   // trigger matching configurables
00057   addTrigMatch_     = iConfig.getParameter<bool>         ( "addTrigMatch" );
00058   trigMatchSrc_     = iConfig.getParameter<std::vector<edm::InputTag> >( "trigPrimMatch" );
00059 
00060   // resolution configurables
00061   addResolutions_   = iConfig.getParameter<bool>         ( "addResolutions" );
00062 
00063   // electron ID configurables
00064   addElecID_        = iConfig.getParameter<bool>         ( "addElectronID" );
00065   if (addElecID_) {
00066       // it might be a single electron ID
00067       if (iConfig.existsAs<edm::InputTag>("electronIDSource")) {
00068           elecIDSrcs_.push_back(NameTag("", iConfig.getParameter<edm::InputTag>("electronIDSource")));
00069       }
00070       // or there might be many of them
00071       if (iConfig.existsAs<edm::ParameterSet>("electronIDSources")) {
00072           // please don't configure me twice
00073           if (!elecIDSrcs_.empty()) throw cms::Exception("Configuration") << 
00074                 "PATElectronProducer: you can't specify both 'electronIDSource' and 'electronIDSources'\n";
00075           // read the different electron ID names
00076           edm::ParameterSet idps = iConfig.getParameter<edm::ParameterSet>("electronIDSources");
00077           std::vector<std::string> names = idps.getParameterNamesForType<edm::InputTag>();
00078           for (std::vector<std::string>::const_iterator it = names.begin(), ed = names.end(); it != ed; ++it) {
00079               elecIDSrcs_.push_back(NameTag(*it, idps.getParameter<edm::InputTag>(*it)));
00080           }
00081       }
00082       // but in any case at least once
00083       if (elecIDSrcs_.empty()) throw cms::Exception("Configuration") <<
00084             "PATElectronProducer: id addElectronID is true, you must specify either:\n" <<
00085             "\tInputTag electronIDSource = <someTag>\n" << "or\n" <<
00086             "\tPSet electronIDSources = { \n" <<
00087             "\t\tInputTag <someName> = <someTag>   // as many as you want \n " <<
00088             "\t}\n";
00089   }
00090   
00091   // construct resolution calculator
00092 
00093   // IsoDeposit configurables
00094   if (iConfig.exists("isoDeposits")) {
00095      edm::ParameterSet depconf = iConfig.getParameter<edm::ParameterSet>("isoDeposits");
00096      if (depconf.exists("tracker")) isoDepositLabels_.push_back(std::make_pair(TrackerIso, depconf.getParameter<edm::InputTag>("tracker")));
00097      if (depconf.exists("ecal"))    isoDepositLabels_.push_back(std::make_pair(ECalIso, depconf.getParameter<edm::InputTag>("ecal")));
00098      if (depconf.exists("hcal"))    isoDepositLabels_.push_back(std::make_pair(HCalIso, depconf.getParameter<edm::InputTag>("hcal")));
00099      if (depconf.exists("particle"))           isoDepositLabels_.push_back(std::make_pair(ParticleIso, depconf.getParameter<edm::InputTag>("particle")));
00100      if (depconf.exists("chargedparticle"))    isoDepositLabels_.push_back(std::make_pair(ChargedParticleIso, depconf.getParameter<edm::InputTag>("chargedparticle")));
00101      if (depconf.exists("neutralparticle")) isoDepositLabels_.push_back(std::make_pair(NeutralParticleIso,depconf.getParameter<edm::InputTag>("neutralparticle")));
00102      if (depconf.exists("gammaparticle"))    isoDepositLabels_.push_back(std::make_pair(GammaParticleIso, depconf.getParameter<edm::InputTag>("gammaparticle")));
00103 
00104 
00105      if (depconf.exists("user")) {
00106         std::vector<edm::InputTag> userdeps = depconf.getParameter<std::vector<edm::InputTag> >("user");
00107         std::vector<edm::InputTag>::const_iterator it = userdeps.begin(), ed = userdeps.end();
00108         int key = UserBaseIso;
00109         for ( ; it != ed; ++it, ++key) {
00110             isoDepositLabels_.push_back(std::make_pair(IsolationKeys(key), *it));
00111         }
00112      }
00113   }
00114 
00115   // Efficiency configurables
00116   addEfficiencies_ = iConfig.getParameter<bool>("addEfficiencies");
00117   if (addEfficiencies_) {
00118      efficiencyLoader_ = pat::helper::EfficiencyLoader(iConfig.getParameter<edm::ParameterSet>("efficiencies"));
00119   }
00120 
00121   // Check to see if the user wants to add user data
00122   useUserData_ = false;
00123   if ( iConfig.exists("userData") ) {
00124     useUserData_ = true;
00125   }
00126 
00127   // electron ID configurables
00128   addElecShapes_        = iConfig.getParameter<bool>("addElectronShapes" );
00129   reducedBarrelRecHitCollection_ = iConfig.getParameter<edm::InputTag>("reducedBarrelRecHitCollection") ;
00130   reducedEndcapRecHitCollection_ = iConfig.getParameter<edm::InputTag>("reducedEndcapRecHitCollection") ;
00131    
00132   // produces vector of muons
00133   produces<std::vector<Electron> >();
00134 
00135 }

PATElectronProducer::~PATElectronProducer (  ) 

Definition at line 138 of file PATElectronProducer.cc.

00138                                           {
00139 }


Member Function Documentation

void PATElectronProducer::FillElectron ( Electron aEl,
const ElectronBaseRef elecRef,
const reco::CandidateBaseRef baseRef,
const GenAssociations genMatches,
const TrigAssociations trigMatches 
) const [private]

Definition at line 317 of file PATElectronProducer.cc.

References addElecShapes_, addGenMatch_, pat::PATObject< ObjectType >::addGenParticleRef(), pat::PATObject< ObjectType >::addTriggerMatch(), addTrigMatch_, efficiencyLoader_, embedGenMatch_, pat::PATObject< ObjectType >::embedGenParticle(), pat::Electron::embedGsfTrack(), embedGsfTrack_, pat::Electron::embedSuperCluster(), embedSuperCluster_, pat::Electron::embedTrack(), embedTrack_, pat::helper::EfficiencyLoader::enabled(), i, edm::Ref< C, T, F >::isAvailable(), edm::Ref< C, T, F >::isNonnull(), lazyTools_, n, pat::Electron::setClusterShapes(), pat::helper::EfficiencyLoader::setEfficiencies(), funct::sqrt(), and trigMatchSrc_.

Referenced by produce().

00321                                                                                  {
00322   if (embedGsfTrack_) anElectron.embedGsfTrack();
00323   if (embedSuperCluster_) anElectron.embedSuperCluster();
00324   if (embedTrack_) anElectron.embedTrack();
00325 
00326   // store the match to the generated final state muons
00327   if (addGenMatch_) {
00328     for(size_t i = 0, n = genMatches.size(); i < n; ++i) {
00329       reco::GenParticleRef genElectron = (*genMatches[i])[elecsRef];
00330       anElectron.addGenParticleRef(genElectron);
00331     }
00332     if (embedGenMatch_) anElectron.embedGenParticle();
00333   }
00334     // matches to trigger primitives
00335     if ( addTrigMatch_ ) {
00336       for ( size_t i = 0; i < trigMatchSrc_.size(); ++i ) {
00337         TriggerPrimitiveRef trigPrim = (*trigMatches[i])[elecsRef];
00338         if ( trigPrim.isNonnull() && trigPrim.isAvailable() ) {
00339           anElectron.addTriggerMatch(*trigPrim);
00340         }
00341       }
00342     }
00343     
00344     if (efficiencyLoader_.enabled()) {
00345       efficiencyLoader_.setEfficiencies( anElectron, elecsRef );
00346     }
00347     //  add electron shapes info
00348     if (addElecShapes_) {
00349         std::vector<float> covariances = lazyTools_->covariances(*(elecsRef->superCluster()->seed())) ;
00350         std::vector<float> localCovariances = lazyTools_->localCovariances(*(elecsRef->superCluster()->seed())) ;
00351         float scSigmaEtaEta = sqrt(covariances[0]) ;
00352         float scSigmaIEtaIEta = sqrt(localCovariances[0]) ;
00353         float scE1x5 = lazyTools_->e1x5(*(elecsRef->superCluster()->seed()))  ;
00354         float scE2x5Max = lazyTools_->e2x5Max(*(elecsRef->superCluster()->seed()))  ;
00355         float scE5x5 = lazyTools_->e5x5(*(elecsRef->superCluster()->seed())) ;
00356         anElectron.setClusterShapes(scSigmaEtaEta,scSigmaIEtaIEta,scE1x5,scE2x5Max,scE5x5) ;
00357     }
00358 }

void PATElectronProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements edm::EDProducer.

Definition at line 142 of file PATElectronProducer.cc.

References pat::PATUserDataHelper< ObjectType >::add(), addElecID_, addElecShapes_, addGenMatch_, addTrigMatch_, pat::helper::MultiIsolator::beginEvent(), efficiencyLoader_, elecIDSrcs_, electronSrc_, pat::Electron::embedPFCandidate(), embedPFCandidate_, pat::helper::MultiIsolator::enabled(), pat::helper::EfficiencyLoader::enabled(), pat::helper::MultiIsolator::endEvent(), pat::helper::MultiIsolator::fill(), FillElectron(), first, genMatchSrc_, edm::Event::getByLabel(), i, index, isoDepositLabels_, isolator_, isolatorTmpStorage_, it, j, lazyTools_, reco::Matched, pat::helper::EfficiencyLoader::newEvent(), pfElecSrc_, python::pfElectrons_cfi::pfElectrons, pTComparator_, ptr, edm::Event::put(), reducedBarrelRecHitCollection_, reducedEndcapRecHitCollection_, edm::second(), pat::Electron::setElectronIDs(), pat::Lepton< LeptonType >::setIsoDeposit(), pat::Lepton< LeptonType >::setIsolation(), pat::Electron::setPFCandidateRef(), python::multivaluedict::sort(), trigMatchSrc_, useParticleFlow_, userDataHelper_, and useUserData_.

00142                                                                                  {
00143 
00144   // Get the collection of electrons from the event
00145   edm::Handle<edm::View<ElectronType> > electrons;
00146   iEvent.getByLabel(electronSrc_, electrons);
00147 
00148   if (isolator_.enabled()) isolator_.beginEvent(iEvent,iSetup);
00149 
00150   if (efficiencyLoader_.enabled()) efficiencyLoader_.newEvent(iEvent);
00151 
00152   std::vector<edm::Handle<edm::ValueMap<IsoDeposit> > > deposits(isoDepositLabels_.size());
00153   for (size_t j = 0, nd = deposits.size(); j < nd; ++j) {
00154     iEvent.getByLabel(isoDepositLabels_[j].second, deposits[j]);
00155   }
00156 
00157   // prepare the MC matching
00158   GenAssociations  genMatches(genMatchSrc_.size());
00159   if (addGenMatch_) {
00160         for (size_t j = 0, nd = genMatchSrc_.size(); j < nd; ++j) {
00161             iEvent.getByLabel(genMatchSrc_[j], genMatches[j]);
00162         }
00163   }
00164 
00165 
00166   // prepare the trigger matching
00167   TrigAssociations  trigMatches(trigMatchSrc_.size());
00168   if ( addTrigMatch_ ) {
00169     for ( size_t i = 0; i < trigMatchSrc_.size(); ++i ) {
00170       iEvent.getByLabel(trigMatchSrc_[i], trigMatches[i]);
00171     }
00172   }
00173 
00174   // prepare ID extraction 
00175   std::vector<edm::Handle<edm::ValueMap<float> > > idhandles;
00176   std::vector<pat::Electron::IdPair>               ids;
00177   if (addElecID_) {
00178      idhandles.resize(elecIDSrcs_.size());
00179      ids.resize(elecIDSrcs_.size());
00180      for (size_t i = 0; i < elecIDSrcs_.size(); ++i) {
00181         iEvent.getByLabel(elecIDSrcs_[i].second, idhandles[i]);
00182         ids[i].first = elecIDSrcs_[i].first;
00183      }
00184   }
00185 
00186 
00187   if (addElecShapes_) {
00188     lazyTools_ .reset(new EcalClusterLazyTools( iEvent , iSetup , reducedBarrelRecHitCollection_ , reducedEndcapRecHitCollection_ ));  
00189   }
00190 
00191   std::vector<Electron> * patElectrons = new std::vector<Electron>();
00192 
00193 
00194   if( useParticleFlow_ ) {
00195     edm::Handle< reco::PFCandidateCollection >  pfElectrons;
00196     iEvent.getByLabel(pfElecSrc_, pfElectrons);
00197     unsigned index=0;
00198   
00199     for( reco::PFCandidateConstIterator i = pfElectrons->begin(); 
00200          i != pfElectrons->end(); ++i, ++index) {
00201 
00202       reco::PFCandidateRef pfRef(pfElectrons,index);
00203       reco::PFCandidatePtr ptrToMother(pfElectrons,index);
00204       reco::CandidateBaseRef pfBaseRef( pfRef ); 
00205       
00206       reco::GsfTrackRef PfTk= i->gsfTrackRef(); 
00207 
00208       bool Matched=false;
00209       for (edm::View<ElectronType>::const_iterator itElectron = electrons->begin(); itElectron != electrons->end(); ++itElectron) {
00210         unsigned int idx = itElectron - electrons->begin();
00211         if (Matched) continue;
00212         reco::GsfTrackRef EgTk= itElectron->gsfTrack();
00213         if (itElectron->gsfTrack()==i->gsfTrackRef()){
00214           const ElectronBaseRef elecsRef = electrons->refAt(idx);
00215           Electron anElectron(elecsRef);
00216           FillElectron(anElectron,elecsRef,pfBaseRef, genMatches, trigMatches);
00217           Matched=true;
00218           anElectron.setPFCandidateRef( pfRef  );
00219           if( embedPFCandidate_ ) anElectron.embedPFCandidate();
00220 
00221           if (isolator_.enabled()){
00222             reco::CandidatePtr mother =  ptrToMother->sourceCandidatePtr(0);
00223             isolator_.fill(mother, isolatorTmpStorage_);
00224             typedef pat::helper::MultiIsolator::IsolationValuePairs IsolationValuePairs;
00225             for (IsolationValuePairs::const_reverse_iterator it = isolatorTmpStorage_.rbegin(), 
00226                    ed = isolatorTmpStorage_.rend(); it != ed; ++it) {
00227               anElectron.setIsolation(it->first, it->second);
00228             }
00229             for (size_t j = 0, nd = deposits.size(); j < nd; ++j) {
00230               anElectron.setIsoDeposit(isoDepositLabels_[j].first, (*deposits[j])[mother]);
00231             }
00232           }
00233           //Electron Id
00234           if (addElecID_) {
00235             //STANDARD EL ID 
00236             for (size_t i = 0; i < elecIDSrcs_.size(); ++i) {
00237               ids[i].second = (*idhandles[i])[elecsRef];    
00238             }
00239             //SPECIFIC PF ID
00240             ids.push_back(std::make_pair("pf_evspi",pfRef->mva_e_pi()));
00241             ids.push_back(std::make_pair("pf_evsmu",pfRef->mva_e_mu()));
00242             anElectron.setElectronIDs(ids);
00243           }
00244 
00245           patElectrons->push_back(anElectron);
00246 
00247         }
00248 
00249         
00250       }
00251  
00252       
00253     }
00254   }
00255 
00256   else{
00257   for (edm::View<ElectronType>::const_iterator itElectron = electrons->begin(); itElectron != electrons->end(); ++itElectron) {
00258     // construct the Electron from the ref -> save ref to original object
00259     unsigned int idx = itElectron - electrons->begin();
00260 
00261     const ElectronBaseRef elecsRef = electrons->refAt(idx);
00262     reco::CandidateBaseRef elecBaseRef(elecsRef);
00263 
00264     //    const edm::Ptr<ElectronType> electronPtr = electrons->ptrAt(idx);
00265     Electron anElectron(elecsRef);
00266 
00267     FillElectron(anElectron,elecsRef,elecBaseRef, genMatches, trigMatches);
00268     
00269 
00270     // add resolution info
00271     
00272     // Isolation
00273     if (isolator_.enabled()) {
00274         isolator_.fill(*electrons, idx, isolatorTmpStorage_);
00275         typedef pat::helper::MultiIsolator::IsolationValuePairs IsolationValuePairs;
00276         // better to loop backwards, so the vector is resized less times
00277         for (IsolationValuePairs::const_reverse_iterator it = isolatorTmpStorage_.rbegin(), ed = isolatorTmpStorage_.rend(); it != ed; ++it) {
00278             anElectron.setIsolation(it->first, it->second);
00279         }
00280     }
00281 
00282     for (size_t j = 0, nd = deposits.size(); j < nd; ++j) {
00283         anElectron.setIsoDeposit(isoDepositLabels_[j].first, (*deposits[j])[elecsRef]);
00284     }
00285 
00286     // add electron ID info
00287     if (addElecID_) {
00288         for (size_t i = 0; i < elecIDSrcs_.size(); ++i) {
00289             ids[i].second = (*idhandles[i])[elecsRef];    
00290         }
00291         anElectron.setElectronIDs(ids);
00292     }
00293     
00294 
00295     if ( useUserData_ ) {
00296       userDataHelper_.add( anElectron, iEvent, iSetup );
00297     }
00298     
00299     
00300     // add sel to selected
00301     patElectrons->push_back(anElectron);
00302   }
00303   }
00304   
00305   // sort electrons in pt
00306   std::sort(patElectrons->begin(), patElectrons->end(), pTComparator_);
00307 
00308   // add the electrons to the event output
00309   std::auto_ptr<std::vector<Electron> > ptr(patElectrons);
00310   iEvent.put(ptr);
00311 
00312   // clean up
00313   if (isolator_.enabled()) isolator_.endEvent();
00314 
00315 }


Member Data Documentation

bool pat::PATElectronProducer::addEfficiencies_ [private]

Definition at line 102 of file PATElectronProducer.h.

Referenced by PATElectronProducer().

bool pat::PATElectronProducer::addElecID_ [private]

Definition at line 74 of file PATElectronProducer.h.

Referenced by PATElectronProducer(), and produce().

bool pat::PATElectronProducer::addElecShapes_ [private]

Definition at line 109 of file PATElectronProducer.h.

Referenced by FillElectron(), PATElectronProducer(), and produce().

bool pat::PATElectronProducer::addGenMatch_ [private]

Definition at line 68 of file PATElectronProducer.h.

Referenced by FillElectron(), PATElectronProducer(), and produce().

bool pat::PATElectronProducer::addResolutions_ [private]

Definition at line 73 of file PATElectronProducer.h.

Referenced by PATElectronProducer().

bool pat::PATElectronProducer::addTrigMatch_ [private]

Definition at line 71 of file PATElectronProducer.h.

Referenced by FillElectron(), PATElectronProducer(), and produce().

pat::helper::EfficiencyLoader pat::PATElectronProducer::efficiencyLoader_ [private]

Definition at line 103 of file PATElectronProducer.h.

Referenced by FillElectron(), PATElectronProducer(), and produce().

std::vector<NameTag> pat::PATElectronProducer::elecIDSrcs_ [private]

Definition at line 93 of file PATElectronProducer.h.

Referenced by PATElectronProducer(), and produce().

edm::InputTag pat::PATElectronProducer::electronSrc_ [private]

Definition at line 64 of file PATElectronProducer.h.

Referenced by PATElectronProducer(), and produce().

bool pat::PATElectronProducer::embedGenMatch_ [private]

Definition at line 69 of file PATElectronProducer.h.

Referenced by FillElectron(), and PATElectronProducer().

bool pat::PATElectronProducer::embedGsfTrack_ [private]

Definition at line 65 of file PATElectronProducer.h.

Referenced by FillElectron(), and PATElectronProducer().

bool pat::PATElectronProducer::embedPFCandidate_ [private]

Definition at line 79 of file PATElectronProducer.h.

Referenced by PATElectronProducer(), and produce().

bool pat::PATElectronProducer::embedSuperCluster_ [private]

Definition at line 66 of file PATElectronProducer.h.

Referenced by FillElectron(), and PATElectronProducer().

bool pat::PATElectronProducer::embedTrack_ [private]

Definition at line 67 of file PATElectronProducer.h.

Referenced by FillElectron(), and PATElectronProducer().

std::vector<edm::InputTag> pat::PATElectronProducer::genMatchSrc_ [private]

Definition at line 70 of file PATElectronProducer.h.

Referenced by PATElectronProducer(), and produce().

std::vector<std::pair<pat::IsolationKeys,edm::InputTag> > pat::PATElectronProducer::isoDepositLabels_ [private]

Definition at line 100 of file PATElectronProducer.h.

Referenced by PATElectronProducer(), and produce().

pat::helper::MultiIsolator pat::PATElectronProducer::isolator_ [private]

Definition at line 98 of file PATElectronProducer.h.

Referenced by produce().

pat::helper::MultiIsolator::IsolationValuePairs pat::PATElectronProducer::isolatorTmpStorage_ [private]

Definition at line 99 of file PATElectronProducer.h.

Referenced by produce().

std::auto_ptr<EcalClusterLazyTools> pat::PATElectronProducer::lazyTools_ [private]

Definition at line 111 of file PATElectronProducer.h.

Referenced by FillElectron(), and produce().

edm::InputTag pat::PATElectronProducer::pfElecSrc_ [private]

Definition at line 78 of file PATElectronProducer.h.

Referenced by PATElectronProducer(), and produce().

GreaterByPt<Electron> pat::PATElectronProducer::pTComparator_ [private]

Definition at line 96 of file PATElectronProducer.h.

Referenced by produce().

edm::InputTag pat::PATElectronProducer::reducedBarrelRecHitCollection_ [private]

Definition at line 114 of file PATElectronProducer.h.

Referenced by PATElectronProducer(), and produce().

edm::InputTag pat::PATElectronProducer::reducedEndcapRecHitCollection_ [private]

Definition at line 115 of file PATElectronProducer.h.

Referenced by PATElectronProducer(), and produce().

std::vector<edm::InputTag> pat::PATElectronProducer::trigMatchSrc_ [private]

Definition at line 72 of file PATElectronProducer.h.

Referenced by FillElectron(), PATElectronProducer(), and produce().

bool pat::PATElectronProducer::useParticleFlow_ [private]

pflow specific

Definition at line 77 of file PATElectronProducer.h.

Referenced by PATElectronProducer(), and produce().

pat::PATUserDataHelper<pat::Electron> pat::PATElectronProducer::userDataHelper_ [private]

Definition at line 106 of file PATElectronProducer.h.

Referenced by produce().

bool pat::PATElectronProducer::useUserData_ [private]

Definition at line 105 of file PATElectronProducer.h.

Referenced by PATElectronProducer(), and produce().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:49:42 2009 for CMSSW by  doxygen 1.5.4