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