CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GEDGsfElectronFinalizer.cc
Go to the documentation of this file.
2 
8 
11 
12 #include <iostream>
13 #include <string>
14 
15 using namespace reco;
16 
18  {
19  previousGsfElectrons_ = consumes<reco::GsfElectronCollection>(cfg.getParameter<edm::InputTag>("previousGsfElectronsTag"));
20  pfCandidates_ = consumes<reco::PFCandidateCollection>(cfg.getParameter<edm::InputTag>("pfCandidatesTag"));
21  outputCollectionLabel_ = cfg.getParameter<std::string>("outputCollectionLabel");
22  edm::ParameterSet pfIsoVals(cfg.getParameter<edm::ParameterSet> ("pfIsolationValues"));
23 
24  tokenElectronIsoVals_.push_back(consumes<edm::ValueMap<double> >(pfIsoVals.getParameter<edm::InputTag>("pfSumChargedHadronPt")));
25  tokenElectronIsoVals_.push_back(consumes<edm::ValueMap<double> >(pfIsoVals.getParameter<edm::InputTag>("pfSumPhotonEt")));
26  tokenElectronIsoVals_.push_back(consumes<edm::ValueMap<double> >(pfIsoVals.getParameter<edm::InputTag>("pfSumNeutralHadronEt")));
27  tokenElectronIsoVals_.push_back(consumes<edm::ValueMap<double> >(pfIsoVals.getParameter<edm::InputTag>("pfSumPUPt")));
28 // std::vector<std::string> isoNames = pfIsoVals.getParameterNamesForType<edm::InputTag>();
29 // for(const std::string& name : isoNames) {
30 // edm::InputTag tag =
31 // pfIsoVals.getParameter<edm::InputTag>(name);
32 // tokenElectronIsoVals_.push_back(consumes<edm::ValueMap<double> >(tag));
33 // }
34 
35  nDeps_ = tokenElectronIsoVals_.size();
36 
37  produces<reco::GsfElectronCollection> (outputCollectionLabel_);
38 }
39 
41  {}
42 
43 // ------------ method called to produce the data ------------
45  {
46 
47  // Output collection
48  std::auto_ptr<reco::GsfElectronCollection> outputElectrons_p(new reco::GsfElectronCollection);
49 
50  // read input collections
51  // electrons
53  event.getByToken(previousGsfElectrons_,gedElectronHandle);
54 
55  // PFCandidates
57  event.getByToken(pfCandidates_,pfCandidateHandle);
58  // value maps
59  std::vector< edm::Handle< edm::ValueMap<double> > > isolationValueMaps(nDeps_);
60 
61  for(unsigned i=0; i < nDeps_ ; ++i) {
62  event.getByToken(tokenElectronIsoVals_[i],isolationValueMaps[i]);
63  }
64 
65  // prepare a map of PFCandidates having a valid GsfTrackRef to save time
66  std::map<reco::GsfTrackRef, const reco::PFCandidate* > gsfPFMap;
67  reco::PFCandidateCollection::const_iterator it = pfCandidateHandle->begin();
68  reco::PFCandidateCollection::const_iterator itend = pfCandidateHandle->end() ;
69  for(;it!=itend;++it) {
70  // First check that the GsfTrack is non null
71  if( it->gsfTrackRef().isNonnull()) {
72  if(abs(it->pdgId())==11) // consider only the electrons
73  gsfPFMap[it->gsfTrackRef()]=&(*it);
74  }
75  }
76 
77 
78  // Now loop on the electrons
79  unsigned nele=gedElectronHandle->size();
80  for(unsigned iele=0; iele<nele;++iele) {
81  reco::GsfElectronRef myElectronRef(gedElectronHandle,iele);
82 
83  reco::GsfElectron newElectron(*myElectronRef);
85  isoVariables.sumChargedHadronPt = (*(isolationValueMaps)[0])[myElectronRef];
86  isoVariables.sumPhotonEt = (*(isolationValueMaps)[1])[myElectronRef];
87  isoVariables.sumNeutralHadronEt = (*(isolationValueMaps)[2])[myElectronRef];
88  isoVariables.sumPUPt = (*(isolationValueMaps)[3])[myElectronRef];
89  newElectron.setPfIsolationVariables(isoVariables);
90 
91  // now set a status if not already done (in GEDGsfElectronProducer.cc)
92  // std::cout << " previous status " << newElectron.mvaOutput().status << std::endl;
93  if(newElectron.mvaOutput().status<=0) {
94  std::map<reco::GsfTrackRef, const reco::PFCandidate * >::const_iterator itcheck=gsfPFMap.find(newElectron.gsfTrack());
95  reco::GsfElectron::MvaOutput myMvaOutput(newElectron.mvaOutput());
96  if(itcheck!=gsfPFMap.end()) {
97  // it means that there is a PFCandidate with the same GsfTrack
98  myMvaOutput.status = 3; //as defined in PFCandidateEGammaExtra.h
99  }
100  else
101  myMvaOutput.status = 4 ; //
102 
103  newElectron.setMvaOutput(myMvaOutput);
104  }
105  outputElectrons_p->push_back(newElectron);
106  }
107 
108  event.put(outputElectrons_p,outputCollectionLabel_);
109  }
110 
111 
virtual void produce(edm::Event &, const edm::EventSetup &)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
float sumPUPt
sum pt of charged Particles not from PV (for Pu corrections)
Definition: GsfElectron.h:563
void setPfIsolationVariables(const PflowIsolationVariables &iso)
Definition: GsfElectron.h:607
void setMvaOutput(const MvaOutput &mo)
Definition: GsfElectron.h:609
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
float sumPhotonEt
sum pt of PF photons // old float photonIso ;
Definition: GsfElectron.h:558
GEDGsfElectronFinalizer(const edm::ParameterSet &)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float sumNeutralHadronEt
sum pt of neutral hadrons // old float neutralHadronIso ;
Definition: GsfElectron.h:557
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
const MvaOutput & mvaOutput() const
Definition: GsfElectron.h:603
float sumChargedHadronPt
sum-pt of charged Hadron // old float chargedHadronIso ;
Definition: GsfElectron.h:556
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
virtual GsfTrackRef gsfTrack() const
reference to a GsfTrack
Definition: GsfElectron.h:183