CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PatElectronAnalyzer.cc
Go to the documentation of this file.
1 #include "TH1.h"
2 #include "TH2.h"
3 #include "TMath.h"
4 
11 
14 
16 
17  public:
18 
19  explicit PatElectronAnalyzer(const edm::ParameterSet&);
21 
22  private:
23 
24  virtual void beginJob() override ;
25  virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
26  virtual void endJob() override ;
27 
28  // restrictions for the electron to be
29  // considered
30  double minPt_;
31  double maxEta_;
32 
33  // decide in which mode to run the analyzer
34  // 0 : genMatch, 1 : tagAndProbe are
35  // supported depending on the line comments
36  unsigned int mode_;
37 
38  // choose a given electronID for the electron
39  // in consideration; the following types are
40  // available:
41  // * eidRobustLoose
42  // * eidRobustTight
43  // * eidLoose
44  // * eidTight
45  // * eidRobustHighEnergy
47 
48  // source of electrons
50  // source of generator particles
52 
55 
56  // internal variables for genMatchMode and
57  // tagAndProbeMode
58  double maxDeltaR_;
59  double maxDeltaM_;
60  double maxTagIso_;
61 
62  // book histograms of interest
63  TH1I *nr_;
64  TH1F *pt_;
65  TH1F *eta_;
66  TH1F *phi_;
67  TH1F *genPt_;
68  TH1F *genEta_;
69  TH1F *genPhi_;
70  TH1F *deltaR_;
71  TH1F *isoTag_;
72  TH1F *invMass_;
73  TH1F *deltaPhi_;
74 };
75 
76 #include "Math/VectorUtil.h"
77 
79  minPt_ (cfg.getParameter<double>("minPt")),
80  maxEta_ (cfg.getParameter<double>("maxEta")),
81  mode_ (cfg.getParameter<unsigned int>("mode")),
82  electronID_ (cfg.getParameter<std::string>("electronID")),
83  electronSrcToken_ (consumes<std::vector<pat::Electron> >(cfg.getParameter<edm::InputTag>("electronSrc"))),
84  particleSrcToken_ (consumes<reco::GenParticleCollection>(cfg.getParameter<edm::InputTag>("particleSrc"))),
85  genMatchMode_(cfg.getParameter<edm::ParameterSet>("genMatchMode")),
86  tagAndProbeMode_(cfg.getParameter<edm::ParameterSet>("tagAndProbeMode"))
87 {
88  // complete the configuration of the analyzer
89  maxDeltaR_ = genMatchMode_ .getParameter<double>("maxDeltaR");
90  maxDeltaM_ = tagAndProbeMode_.getParameter<double>("maxDeltaM");
91  maxTagIso_ = tagAndProbeMode_.getParameter<double>("maxTagIso");
92 
93 
94  // register histograms to the TFileService
96  nr_ = fs->make<TH1I>("nr", "nr", 10, 0 , 10 );
97  pt_ = fs->make<TH1F>("pt", "pt", 20, 0., 100.);
98  eta_ = fs->make<TH1F>("eta", "eta", 30, -3., 3.);
99  phi_ = fs->make<TH1F>("phi", "phi", 35, -3.5, 3.5);
100  genPt_ = fs->make<TH1F>("genPt", "pt", 20, 0., 100.);
101  genEta_ = fs->make<TH1F>("genEta", "eta", 30, -3., 3.);
102  genPhi_ = fs->make<TH1F>("genPhi", "phi", 35, -3.5, 3.5);
103  deltaR_ = fs->make<TH1F>("deltaR", "log(dR)", 50, -5., 0.);
104  isoTag_ = fs->make<TH1F>("isoTag", "iso", 50, 0., 10.);
105  invMass_ = fs->make<TH1F>("invMass", "m", 100, 50., 150.);
106  deltaPhi_= fs->make<TH1F>("deltaPhi", "deltaPhi", 100, -3.5, 3.5);
107 }
108 
110 {
111 }
112 
113 void
115 {
116  // get electron collection
118  evt.getByToken(electronSrcToken_, electrons);
119  // get generator particle collection
121  evt.getByToken(particleSrcToken_, particles);
122 
123  nr_->Fill( electrons->size() );
124 
125  // ----------------------------------------------------------------------
126  //
127  // First Part Mode 0: genMatch
128  //
129  // ----------------------------------------------------------------------
130  if( mode_==0 ){
131  // loop generator particles
132  for(reco::GenParticleCollection::const_iterator part=particles->begin();
133  part!=particles->end(); ++part){
134  // only loop stable electrons
135  if( part->status()==1 && abs(part->pdgId())==11 ){
136  if( part->pt()>minPt_ && fabs(part->eta())<maxEta_ ){
137  genPt_ ->Fill( part->pt() );
138  genEta_->Fill( part->eta() );
139  genPhi_->Fill( part->phi() );
140  }
141  }
142  }
143 
144  // loop electrons
145  for( std::vector<pat::Electron>::const_iterator elec=electrons->begin(); elec!=electrons->end(); ++elec ){
146  if( elec->genLepton() ){
147  float deltaR = ROOT::Math::VectorUtil::DeltaR(elec->genLepton()->p4(), elec->p4());
148  deltaR_->Fill(TMath::Log10(deltaR));
149  if( deltaR<maxDeltaR_ ){
150  if( electronID_.compare("none")!=0 ){
151  if( elec->electronID(electronID_)<0.5 )
152  continue;
153  }
154  if( elec->pt()>minPt_ && fabs(elec->eta())<maxEta_ ){
155  pt_ ->Fill( elec->pt() );
156  eta_->Fill( elec->eta() );
157  phi_->Fill( elec->phi() );
158  }
159  }
160  }
161  }
162  }
163 
164  // ----------------------------------------------------------------------
165  //
166  // Second Part Mode 1: tagAndProbe
167  //
168  // ----------------------------------------------------------------------
169  if( mode_==1 ){
170  // loop tag electron
171  for( std::vector<pat::Electron>::const_iterator elec=electrons->begin(); elec!=electrons->end(); ++elec ){
172  isoTag_->Fill(elec->trackIso());
173  if( elec->trackIso()<maxTagIso_ && elec->electronID("eidTight")>0.5 ){
174  // loop probe electron
175  for( std::vector<pat::Electron>::const_iterator probe=electrons->begin(); probe!=electrons->end(); ++probe ){
176  // skip the tag electron itself
177  if( probe==elec ) continue;
178 
179  float zMass = (probe->p4()+elec->p4()).mass();
180  invMass_ ->Fill(zMass);
181  float deltaPhi = ROOT::Math::VectorUtil::DeltaPhi(elec->p4(), probe->p4());
182  deltaPhi_->Fill(deltaPhi);
183 
184  // check for the Z mass
185  if( fabs( zMass-90. )<maxDeltaM_ ){
186  if( electronID_.compare("none")!=0 ){
187  if( probe->electronID(electronID_)<0.5 )
188  continue;
189  }
190  if( probe->pt()>minPt_ && fabs(probe->eta())<maxEta_ ){
191  pt_ ->Fill( probe->pt() );
192  eta_->Fill( probe->eta() );
193  phi_->Fill( probe->phi() );
194  }
195  }
196  }
197  }
198  }
199  }
200 }
201 
203 {
204 }
205 
207 {
208 }
209 
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
T getParameter(std::string const &) const
tuple cfg
Definition: looper.py:293
edm::EDGetTokenT< std::vector< pat::Electron > > electronSrcToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void beginJob() override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< reco::GenParticleCollection > particleSrcToken_
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
edm::ParameterSet tagAndProbeMode_
part
Definition: HCALResponse.h:20
edm::ParameterSet genMatchMode_
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
PatElectronAnalyzer(const edm::ParameterSet &)
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
virtual void endJob() override