CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DQMExample_Step1.cc
Go to the documentation of this file.
3 
4 
5 // Geometry
9 #include "TLorentzVector.h"
10 
11 
12 #include <iostream>
13 #include <iomanip>
14 #include <stdio.h>
15 #include <string>
16 #include <sstream>
17 #include <math.h>
18 
19 //
20 // -------------------------------------- Constructor --------------------------------------------
21 //
23 {
24  edm::LogInfo("DQMExample_Step1") << "Constructor DQMExample_Step1::DQMExample_Step1 " << std::endl;
25 
26  // Get parameters from configuration file
27  theElectronCollection_ = consumes<reco::GsfElectronCollection>(ps.getParameter<edm::InputTag>("electronCollection"));
28  theCaloJetCollection_ = consumes<reco::CaloJetCollection>(ps.getParameter<edm::InputTag>("caloJetCollection"));
29  thePfMETCollection_ = consumes<reco::PFMETCollection>(ps.getParameter<edm::InputTag>("pfMETCollection"));
30  theConversionCollection_ = consumes<reco::ConversionCollection>(ps.getParameter<edm::InputTag>("conversionsCollection"));
31  thePVCollection_ = consumes<reco::VertexCollection>(ps.getParameter<edm::InputTag>("PVCollection"));
32  theBSCollection_ = consumes<reco::BeamSpot>(ps.getParameter<edm::InputTag>("beamSpotCollection"));
33  triggerEvent_ = consumes<trigger::TriggerEvent>(ps.getParameter<edm::InputTag>("TriggerEvent"));
34  triggerResults_ = consumes<edm::TriggerResults>(ps.getParameter<edm::InputTag>("TriggerResults"));
35  triggerFilter_ = ps.getParameter<edm::InputTag>("TriggerFilter");
36  triggerPath_ = ps.getParameter<std::string>("TriggerPath");
37 
38 
39  // cuts:
40  ptThrL1_ = ps.getUntrackedParameter<double>("PtThrL1");
41  ptThrL2_ = ps.getUntrackedParameter<double>("PtThrL2");
42  ptThrJet_ = ps.getUntrackedParameter<double>("PtThrJet");
43  ptThrMet_ = ps.getUntrackedParameter<double>("PtThrMet");
44 
45  //DQMStore
47 }
48 
49 //
50 // -- Destructor
51 //
53 {
54  edm::LogInfo("DQMExample_Step1") << "Destructor DQMExample_Step1::~DQMExample_Step1 " << std::endl;
55 }
56 
57 //
58 // -------------------------------------- beginJob --------------------------------------------
59 //
61 {
62  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::beginJob " << std::endl;
63 }
64 //
65 // -------------------------------------- beginRun --------------------------------------------
66 //
68 {
69  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::beginRun" << std::endl;
70 
71  //book at beginRun
72  bookHistos(dbe_);
73 }
74 //
75 // -------------------------------------- beginLuminosityBlock --------------------------------------------
76 //
78  edm::EventSetup const& context)
79 {
80  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::beginLuminosityBlock" << std::endl;
81 }
82 
83 
84 //
85 // -------------------------------------- Analyze --------------------------------------------
86 //
88 {
89  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::analyze" << std::endl;
90 
91 
92  //-------------------------------
93  //--- Vertex Info
94  //-------------------------------
96  e.getByToken(thePVCollection_, vertexHandle);
97  if ( !vertexHandle.isValid() )
98  {
99  edm::LogError ("DQMClientExample") << "invalid collection: vertex" << "\n";
100  return;
101  }
102 
103  int vertex_number = vertexHandle->size();
104  reco::VertexCollection::const_iterator v = vertexHandle->begin();
105 
106  math::XYZPoint PVPoint(-999, -999, -999);
107  if(vertex_number != 0)
108  PVPoint = math::XYZPoint(v->position().x(), v->position().y(), v->position().z());
109 
110  PVPoint_=PVPoint;
111 
112  //-------------------------------
113  //--- MET
114  //-------------------------------
115  edm::Handle<reco::PFMETCollection> pfMETCollection;
116  e.getByToken(thePfMETCollection_, pfMETCollection);
117  if ( !pfMETCollection.isValid() )
118  {
119  edm::LogError ("DQMClientExample") << "invalid collection: MET" << "\n";
120  return;
121  }
122  //-------------------------------
123  //--- Electrons
124  //-------------------------------
125  edm::Handle<reco::GsfElectronCollection> electronCollection;
126  e.getByToken(theElectronCollection_, electronCollection);
127  if ( !electronCollection.isValid() )
128  {
129  edm::LogError ("DQMClientExample") << "invalid collection: electrons" << "\n";
130  return;
131  }
132 
133  float nEle=0;
134  int posEle=0, negEle=0;
135  const reco::GsfElectron* ele1 = NULL;
136  const reco::GsfElectron* ele2 = NULL;
137  for (reco::GsfElectronCollection::const_iterator recoElectron=electronCollection->begin(); recoElectron!=electronCollection->end(); ++recoElectron)
138  {
139  //decreasing pT
140  if( MediumEle(e,eSetup,*recoElectron) )
141  {
142  if(!ele1 && recoElectron->pt() > ptThrL1_)
143  ele1 = &(*recoElectron);
144 
145  else if(!ele2 && recoElectron->pt() > ptThrL2_)
146  ele2 = &(*recoElectron);
147 
148  }
149 
150  if(recoElectron->charge()==1)
151  posEle++;
152  else if(recoElectron->charge()==-1)
153  negEle++;
154 
155  } // end of loop over electrons
156 
157  nEle = posEle+negEle;
158 
159  //-------------------------------
160  //--- Jets
161  //-------------------------------
162  edm::Handle<reco::CaloJetCollection> caloJetCollection;
163  e.getByToken (theCaloJetCollection_,caloJetCollection);
164  if ( !caloJetCollection.isValid() )
165  {
166  edm::LogError ("DQMClientExample") << "invalid collection: jets" << "\n";
167  return;
168  }
169 
170  int nJet = 0;
171  const reco::CaloJet* jet1 = NULL;
172  const reco::CaloJet* jet2 = NULL;
173 
174  for (reco::CaloJetCollection::const_iterator i_calojet = caloJetCollection->begin(); i_calojet != caloJetCollection->end(); ++i_calojet)
175  {
176  //remove jet-ele matching
177  if(ele1)
178  if (Distance(*i_calojet,*ele1) < 0.3) continue;
179 
180  if(ele2)
181  if (Distance(*i_calojet,*ele2) < 0.3) continue;
182 
183  if (i_calojet->pt() < ptThrJet_) continue;
184 
185  nJet++;
186 
187  if (!jet1)
188  jet1 = &(*i_calojet);
189 
190  else if (!jet2)
191  jet2 = &(*i_calojet);
192  }
193 
194  // ---------------------------
195  // ---- Analyze Trigger Event
196  // ---------------------------
197 
198  //check what is in the menu
200  e.getByToken(triggerResults_,hltresults);
201 
202  if(!hltresults.isValid())
203  {
204  edm::LogError ("DQMClientExample") << "invalid collection: TriggerResults" << "\n";
205  return;
206  }
207 
208  bool hasFired = false;
209  const edm::TriggerNames& trigNames = e.triggerNames(*hltresults);
210  unsigned int numTriggers = trigNames.size();
211 
212  for( unsigned int hltIndex=0; hltIndex<numTriggers; ++hltIndex )
213  {
214  if (trigNames.triggerName(hltIndex)==triggerPath_ && hltresults->wasrun(hltIndex) && hltresults->accept(hltIndex))
215  hasFired = true;
216  }
217 
218 
219 
220  //access the trigger event
222  e.getByToken(triggerEvent_, triggerEvent);
223  if( triggerEvent.failedToGet() )
224  {
225  edm::LogError ("DQMClientExample") << "invalid collection: TriggerEvent" << "\n";
226  return;
227  }
228 
229 
230  reco::Particle* ele1_HLT = NULL;
231  int nEle_HLT = 0;
232 
233  size_t filterIndex = triggerEvent->filterIndex( triggerFilter_ );
234  trigger::TriggerObjectCollection triggerObjects = triggerEvent->getObjects();
235  if( !(filterIndex >= triggerEvent->sizeFilters()) )
236  {
237  const trigger::Keys& keys = triggerEvent->filterKeys( filterIndex );
238  std::vector<reco::Particle> triggeredEle;
239 
240  for( size_t j = 0; j < keys.size(); ++j )
241  {
242  trigger::TriggerObject foundObject = triggerObjects[keys[j]];
243  if( abs( foundObject.particle().pdgId() ) != 11 ) continue; //make sure that it is an electron
244 
245  triggeredEle.push_back( foundObject.particle() );
246  ++nEle_HLT;
247  }
248 
249  if( triggeredEle.size() >= 1 )
250  ele1_HLT = &(triggeredEle.at(0));
251  }
252 
253  //-------------------------------
254  //--- Fill the histos
255  //-------------------------------
256 
257  //vertex
258  h_vertex_number -> Fill( vertex_number );
259 
260  //met
261  h_pfMet -> Fill( pfMETCollection->begin()->et() );
262 
263  //multiplicities
264  h_eMultiplicity->Fill(nEle);
265  h_jMultiplicity->Fill(nJet);
266  h_eMultiplicity_HLT->Fill(nEle_HLT);
267 
268  //leading not matched
269  if(ele1)
270  {
271  h_ePt_leading->Fill(ele1->pt());
272  h_eEta_leading->Fill(ele1->eta());
273  h_ePhi_leading->Fill(ele1->phi());
274  }
275  if(ele1_HLT)
276  {
277  h_ePt_leading_HLT->Fill(ele1_HLT->pt());
278  h_eEta_leading_HLT->Fill(ele1_HLT->eta());
279  h_ePhi_leading_HLT->Fill(ele1_HLT->phi());
280  }
281  //leading Jet
282  if(jet1)
283  {
284  h_jPt_leading->Fill(jet1->pt());
285  h_jEta_leading->Fill(jet1->eta());
286  h_jPhi_leading->Fill(jet1->phi());
287  }
288 
289 
290  //fill only when the trigger candidate mathes with the reco one
291  if( ele1 && ele1_HLT && deltaR(*ele1_HLT,*ele1) < 0.3 && hasFired==true )
292  {
293  h_ePt_leading_matched->Fill(ele1->pt());
294  h_eEta_leading_matched->Fill(ele1->eta());
295  h_ePhi_leading_matched->Fill(ele1->phi());
296 
297  h_ePt_leading_HLT_matched->Fill(ele1_HLT->pt());
298  h_eEta_leading_HLT_matched->Fill(ele1_HLT->eta());
299  h_ePhi_leading_HLT_matched->Fill(ele1_HLT->phi());
300 
301  h_ePt_diff->Fill(ele1->pt()-ele1_HLT->pt());
302  }
303 }
304 //
305 // -------------------------------------- endLuminosityBlock --------------------------------------------
306 //
308 {
309  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::endLuminosityBlock" << std::endl;
310 }
311 
312 
313 //
314 // -------------------------------------- endRun --------------------------------------------
315 //
317 {
318  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::endRun" << std::endl;
319 }
320 
321 //
322 // -------------------------------------- endJob --------------------------------------------
323 //
325 {
326  edm::LogInfo("DQMExample_Step1") << "DQMExample_Step1::endJob" << std::endl;
327 }
328 
329 
330 //
331 // -------------------------------------- book histograms --------------------------------------------
332 //
334 {
335  dbe->cd();
336  dbe->setCurrentFolder("Physics/TopTest");
337 
338  h_vertex_number = dbe->book1D("Vertex_number", "Number of event vertices in collection", 40,-0.5, 39.5 );
339 
340  h_pfMet = dbe->book1D("pfMet", "Pf Missing E_{T}; GeV" , 20, 0.0 , 100);
341 
342  h_eMultiplicity = dbe_->book1D("NElectrons","# of electrons per event",10,0.,10.);
343  h_ePt_leading_matched = dbe_->book1D("ElePt_leading_matched","Pt of leading electron",50,0.,100.);
344  h_eEta_leading_matched = dbe_->book1D("EleEta_leading_matched","Eta of leading electron",50,-5.,5.);
345  h_ePhi_leading_matched = dbe_->book1D("ElePhi_leading_matched","Phi of leading electron",50,-3.5,3.5);
346 
347  h_ePt_leading = dbe_->book1D("ElePt_leading","Pt of leading electron",50,0.,100.);
348  h_eEta_leading = dbe_->book1D("EleEta_leading","Eta of leading electron",50,-5.,5.);
349  h_ePhi_leading = dbe_->book1D("ElePhi_leading","Phi of leading electron",50,-3.5,3.5);
350 
351  h_jMultiplicity = dbe_->book1D("NJets","# of electrons per event",10,0.,10.);
352  h_jPt_leading = dbe_->book1D("JetPt_leading","Pt of leading Jet",150,0.,300.);
353  h_jEta_leading = dbe_->book1D("JetEta_leading","Eta of leading Jet",50,-5.,5.);
354  h_jPhi_leading = dbe_->book1D("JetPhi_leading","Phi of leading Jet",50,-3.5,3.5);
355 
356  h_eMultiplicity_HLT = dbe_->book1D("NElectrons_HLT","# of electrons per event @HLT",10,0.,10.);
357  h_ePt_leading_HLT = dbe_->book1D("ElePt_leading_HLT","Pt of leading electron @HLT",50,0.,100.);
358  h_eEta_leading_HLT = dbe_->book1D("EleEta_leading_HLT","Eta of leading electron @HLT",50,-5.,5.);
359  h_ePhi_leading_HLT = dbe_->book1D("ElePhi_leading_HLT","Phi of leading electron @HLT",50,-3.5,3.5);
360 
361  h_ePt_leading_HLT_matched = dbe_->book1D("ElePt_leading_HLT_matched","Pt of leading electron @HLT",50,0.,100.);
362  h_eEta_leading_HLT_matched = dbe_->book1D("EleEta_leading_HLT_matched","Eta of leading electron @HLT",50,-5.,5.);
363  h_ePhi_leading_HLT_matched = dbe_->book1D("ElePhi_leading_HLT_matched","Phi of leading electron @HLT",50,-3.5,3.5);
364 
365  h_ePt_diff = dbe_->book1D("ElePt_diff_matched","pT(RECO) - pT(HLT) for mathed candidates",100,-10,10.);
366 
367  dbe->cd();
368 
369 }
370 
371 
372 //
373 // -------------------------------------- functions --------------------------------------------
374 //
376  return deltaR(c1,c2);
377 }
378 
380  return deltaPhi(c1.p4().phi(),c2.p4().phi());
381 }
382 
383 // This always returns only a positive deltaPhi
384 double DQMExample_Step1::calcDeltaPhi(double phi1, double phi2) {
385  double deltaPhi = phi1 - phi2;
386  if (deltaPhi < 0) deltaPhi = -deltaPhi;
387  if (deltaPhi > 3.1415926) {
388  deltaPhi = 2 * 3.1415926 - deltaPhi;
389  }
390  return deltaPhi;
391 }
392 
393 //
394 // -------------------------------------- electronID --------------------------------------------
395 //
397 {
398 
399  //********* CONVERSION TOOLS
401  iEvent.getByToken(theConversionCollection_, conversions_h);
402 
403  bool isMediumEle = false;
404 
405  float pt = electron.pt();
406  float eta = electron.eta();
407 
408  int isEB = electron.isEB();
409  float sigmaIetaIeta = electron.sigmaIetaIeta();
410  float DetaIn = electron.deltaEtaSuperClusterTrackAtVtx();
411  float DphiIn = electron.deltaPhiSuperClusterTrackAtVtx();
412  float HOverE = electron.hadronicOverEm();
413  float ooemoop = (1.0/electron.ecalEnergy() - electron.eSuperClusterOverP()/electron.ecalEnergy());
414 
415  int mishits = electron.gsfTrack()->trackerExpectedHitsInner().numberOfHits();
416  int nAmbiguousGsfTracks = electron.ambiguousGsfTracksSize();
417 
418  reco::GsfTrackRef eleTrack = electron.gsfTrack() ;
419  float dxy = eleTrack->dxy(PVPoint_);
420  float dz = eleTrack->dz (PVPoint_);
421 
423  iEvent.getByToken(theBSCollection_, BSHandle);
424  const reco::BeamSpot BS = *BSHandle;
425 
426  bool isConverted = ConversionTools::hasMatchedConversion(electron, conversions_h, BS.position());
427 
428  // default
429  if( (pt > 12.) && (fabs(eta) < 2.5) &&
430  ( ( (isEB == 1) && (fabs(DetaIn) < 0.004) ) || ( (isEB == 0) && (fabs(DetaIn) < 0.007) ) ) &&
431  ( ( (isEB == 1) && (fabs(DphiIn) < 0.060) ) || ( (isEB == 0) && (fabs(DphiIn) < 0.030) ) ) &&
432  ( ( (isEB == 1) && (sigmaIetaIeta < 0.010) ) || ( (isEB == 0) && (sigmaIetaIeta < 0.030) ) ) &&
433  ( ( (isEB == 1) && (HOverE < 0.120) ) || ( (isEB == 0) && (HOverE < 0.100) ) ) &&
434  ( ( (isEB == 1) && (fabs(ooemoop) < 0.050) ) || ( (isEB == 0) && (fabs(ooemoop) < 0.050) ) ) &&
435  ( ( (isEB == 1) && (fabs(dxy) < 0.020) ) || ( (isEB == 0) && (fabs(dxy) < 0.020) ) ) &&
436  ( ( (isEB == 1) && (fabs(dz) < 0.100) ) || ( (isEB == 0) && (fabs(dz) < 0.100) ) ) &&
437  ( ( (isEB == 1) && (!isConverted) ) || ( (isEB == 0) && (!isConverted) ) ) &&
438  ( mishits == 0 ) &&
439  ( nAmbiguousGsfTracks == 0 )
440  )
441  isMediumEle=true;
442 
443  return isMediumEle;
444 }
double calcDeltaPhi(double phi1, double phi2)
MonitorElement * h_ePhi_leading_matched
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:204
MonitorElement * h_jMultiplicity
math::XYZPoint PVPoint_
bool MediumEle(const edm::Event &iEvent, const edm::EventSetup &iESetup, const reco::GsfElectron &electron)
edm::EDGetTokenT< reco::GsfElectronCollection > theElectronCollection_
Jets made from CaloTowers.
Definition: CaloJet.h:29
float eSuperClusterOverP() const
Definition: GsfElectron.h:230
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:872
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
void cd(void)
go to top directory (ie. root)
Definition: DQMStore.cc:561
void endLuminosityBlock(edm::LuminosityBlock const &lumi, edm::EventSetup const &eSetup)
reco::Particle particle(reco::Particle::Charge q=0, const reco::Particle::Point &vertex=reco::Particle::Point(0, 0, 0), int status=0, bool integerCharge=true) const
Definition: TriggerObject.h:69
MonitorElement * h_ePhi_leading_HLT_matched
#define NULL
Definition: scimark2.h:8
MonitorElement * h_ePhi_leading
MonitorElement * h_eEta_leading_HLT
MonitorElement * h_jPhi_leading
edm::EDGetTokenT< reco::ConversionCollection > theConversionCollection_
Strings::size_type size() const
Definition: TriggerNames.cc:39
edm::EDGetTokenT< trigger::TriggerEvent > triggerEvent_
T eta() const
void bookHistos(DQMStore *dbe_)
double Distance(const reco::Candidate &c1, const reco::Candidate &c2)
bool isEB() const
Definition: GsfElectron.h:334
double eta() const
momentum pseudorapidity
Definition: Particle.cc:168
DQMExample_Step1(const edm::ParameterSet &ps)
void Fill(long long x)
MonitorElement * h_eEta_leading
edm::EDGetTokenT< reco::BeamSpot > theBSCollection_
MonitorElement * h_ePt_leading_HLT
Single trigger physics object (e.g., an isolated muon)
Definition: TriggerObject.h:22
double phi() const
momentum azimuthal angle
Definition: Particle.cc:159
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:234
int iEvent
Definition: GenABIO.cc:243
MonitorElement * h_ePt_diff
MonitorElement * h_ePt_leading_HLT_matched
float sigmaIetaIeta() const
Definition: GsfElectron.h:386
float hadronicOverEm() const
Definition: GsfElectron.h:410
MonitorElement * h_jEta_leading
void beginRun(edm::Run const &run, edm::EventSetup const &eSetup)
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:237
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
MonitorElement * h_eEta_leading_HLT_matched
MonitorElement * h_vertex_number
MonitorElement * h_ePt_leading
bool isValid() const
Definition: HandleBase.h:76
double pt() const
transverse momentum
Definition: Particle.cc:155
static bool hasMatchedConversion(const reco::GsfElectron &ele, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
MonitorElement * h_eMultiplicity
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
void analyze(edm::Event const &e, edm::EventSetup const &eSetup)
MonitorElement * h_pfMet
std::vector< TriggerObject > TriggerObjectCollection
collection of trigger physics objects (e.g., all isolated muons)
Definition: TriggerObject.h:81
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:74
bool failedToGet() const
Definition: HandleBase.h:80
edm::EDGetTokenT< edm::TriggerResults > triggerResults_
void endRun(edm::Run const &run, edm::EventSetup const &eSetup)
MonitorElement * h_jPt_leading
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:27
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< size_type > Keys
MonitorElement * h_ePt_leading_matched
virtual ~DQMExample_Step1()
GsfTrackRefVector::size_type ambiguousGsfTracksSize() const
Definition: GsfElectron.h:602
int pdgId() const
PDG identifier.
Definition: Particle.cc:246
edm::EDGetTokenT< reco::PFMETCollection > thePfMETCollection_
edm::InputTag triggerFilter_
float ecalEnergy() const
Definition: GsfElectron.h:741
MonitorElement * h_ePhi_leading_HLT
double DistancePhi(const reco::Candidate &c1, const reco::Candidate &c2)
std::string triggerPath_
const Point & position() const
position
Definition: BeamSpot.h:62
void beginLuminosityBlock(edm::LuminosityBlock const &lumi, edm::EventSetup const &eSetup)
MonitorElement * h_eMultiplicity_HLT
virtual float pt() const GCC11_FINAL
transverse momentum
MonitorElement * h_eEta_leading_matched
edm::EDGetTokenT< reco::VertexCollection > thePVCollection_
edm::EDGetTokenT< reco::CaloJetCollection > theCaloJetCollection_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:584
Definition: Run.h:41
virtual GsfTrackRef gsfTrack() const
reference to a GsfTrack
Definition: GsfElectron.h:170
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector