test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
QcdPhotonsDQM.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author Michael B. Anderson, University of Wisconsin Madison
5  */
6 
8 
12 
17 
19 
21 
22 // Physics Objects
25 
26 // Vertex
28 
29 // For removing ECAL Spikes
34 
37 
38 // geometry
39 //#include "Geometry/CaloGeometry/interface/CaloGeometry.h"
40 //#include "Geometry/Records/interface/CaloGeometryRecord.h"
41 //#include "Geometry/CaloEventSetup/interface/CaloTopologyRecord.h"
42 //#include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
43 //#include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
44 //#include "Geometry/EcalAlgo/interface/EcalPreshowerGeometry.h"
45 
46 // Math stuff
49 
50 #include <vector>
51 
52 #include <string>
53 #include <cmath>
54 using namespace std;
55 using namespace edm;
56 using namespace reco;
57 
59  // Get parameters from configuration file
60  theTriggerPathToPass_ = parameters.getParameter<string>("triggerPathToPass");
61  thePlotTheseTriggersToo_ =
62  parameters.getParameter<vector<string> >("plotTheseTriggersToo");
63  theJetCollectionLabel_ = parameters.getParameter<InputTag>("jetCollection");
64  trigTagToken_ = consumes<edm::TriggerResults>(
65  parameters.getUntrackedParameter<edm::InputTag>("trigTag"));
66  thePhotonCollectionToken_ = consumes<reco::PhotonCollection>(
67  parameters.getParameter<InputTag>("photonCollection"));
68  theJetCollectionToken_ = consumes<edm::View<reco::Jet> >(
69  parameters.getParameter<InputTag>("jetCollection"));
70  theVertexCollectionToken_ = consumes<reco::VertexCollection>(
71  parameters.getParameter<InputTag>("vertexCollection"));
72  theMinJetPt_ = parameters.getParameter<double>("minJetPt");
73  theMinPhotonEt_ = parameters.getParameter<double>("minPhotonEt");
74  theRequirePhotonFound_ = parameters.getParameter<bool>("requirePhotonFound");
75  thePlotPhotonMaxEt_ = parameters.getParameter<double>("plotPhotonMaxEt");
76  thePlotPhotonMaxEta_ = parameters.getParameter<double>("plotPhotonMaxEta");
77  thePlotJetMaxEta_ = parameters.getParameter<double>("plotJetMaxEta");
78  theBarrelRecHitTag_ = parameters.getParameter<InputTag>("barrelRecHitTag");
79  theEndcapRecHitTag_ = parameters.getParameter<InputTag>("endcapRecHitTag");
80  theBarrelRecHitToken_ = consumes<EcalRecHitCollection>(
81  parameters.getParameter<InputTag>("barrelRecHitTag"));
82  theEndcapRecHitToken_ = consumes<EcalRecHitCollection>(
83  parameters.getParameter<InputTag>("endcapRecHitTag"));
84 
85  // coverity says...
86  h_deltaEt_photon_jet = 0;
87  h_deltaPhi_jet_jet2 = 0;
88  h_deltaPhi_photon_jet = 0;
89  h_deltaPhi_photon_jet2 = 0;
90  h_deltaR_jet_jet2 = 0;
91  h_deltaR_photon_jet2 = 0;
92  h_jet2_eta = 0;
93  h_jet2_pt = 0;
94  h_jet2_ptOverPhotonEt = 0;
95  h_jet_count = 0;
96  h_jet_eta = 0;
97  h_jet_pt = 0;
98  h_photon_count_bar = 0;
99  h_photon_count_end = 0;
100  h_photon_et = 0;
101  h_photon_et_beforeCuts = 0;
102  h_photon_et_jetco = 0;
103  h_photon_et_jetcs = 0;
104  h_photon_et_jetfo = 0;
105  h_photon_et_jetfs = 0;
106  h_photon_eta = 0;
107  h_triggers_passed = 0;
108 
109 }
110 
112 
113 
115  edm::Run const &, edm::EventSetup const & ){
116 
117  logTraceName = "QcdPhotonAnalyzer";
118 
119  LogTrace(logTraceName) << "Parameters initialization";
120 
121  ibooker.setCurrentFolder("Physics/QcdPhotons"); // Use folder with name of PAG
122 
123  std::stringstream aStringStream;
124  std::string aString;
125  aStringStream << theMinJetPt_;
126  aString = aStringStream.str();
127 
128  // Monitor of triggers passed
129  int numOfTriggersToMonitor = thePlotTheseTriggersToo_.size();
130  h_triggers_passed = ibooker.book1D("triggers_passed",
131  "Events passing these trigger paths", numOfTriggersToMonitor,
132  0, numOfTriggersToMonitor);
133  for (int i = 0; i < numOfTriggersToMonitor; i++) {
134  h_triggers_passed->setBinLabel(i + 1, thePlotTheseTriggersToo_[i]);
135  }
136 
137  // Keep the number of plots and number of bins to a minimum!
138  h_photon_et_beforeCuts = ibooker.book1D("photon_et_beforeCuts",
139  "#gamma with highest E_{T};E_{T}(#gamma) (GeV)", 20, 0., thePlotPhotonMaxEt_);
140  h_photon_et = ibooker.book1D("photon_et",
141  "#gamma with highest E_{T};E_{T}(#gamma) (GeV)", 20, 0., thePlotPhotonMaxEt_);
142  h_photon_eta = ibooker.book1D("photon_eta",
143  "#gamma with highest E_{T};#eta(#gamma)", 40,
144  -thePlotPhotonMaxEta_, thePlotPhotonMaxEta_);
145  h_photon_count_bar = ibooker.book1D("photon_count_bar",
146  "Number of #gamma's passing selection (Barrel);Number of #gamma's", 8,
147  -0.5, 7.5);
148  h_photon_count_end = ibooker.book1D("photon_count_end",
149  "Number of #gamma's passing selection (Endcap);Number of #gamma's", 8,
150  -0.5, 7.5);
151  h_jet_pt = ibooker.book1D("jet_pt",
152  "Jet with highest p_{T} (from " + theJetCollectionLabel_.label() +
153  ");p_{T}(1^{st} jet) (GeV)", 20, 0., thePlotPhotonMaxEt_);
154  h_jet_eta = ibooker.book1D("jet_eta",
155  "Jet with highest p_{T} (from " + theJetCollectionLabel_.label() +
156  ");#eta(1^{st} jet)", 20, -thePlotJetMaxEta_, thePlotJetMaxEta_);
157  h_deltaPhi_photon_jet = ibooker.book1D("deltaPhi_photon_jet",
158  "#Delta#phi between Highest E_{T} #gamma and jet;#Delta#phi(#gamma,1^{st} jet)",
159  20, 0, 3.1415926);
160  h_deltaPhi_jet_jet2 = ibooker.book1D("deltaPhi_jet_jet2",
161  "#Delta#phi between Highest E_{T} jet and 2^{nd} "
162  "jet;#Delta#phi(1^{st} jet,2^{nd} jet)",
163  20, 0, 3.1415926);
164  h_deltaEt_photon_jet = ibooker.book1D("deltaEt_photon_jet",
165  "(E_{T}(#gamma)-p_{T}(jet))/E_{T}(#gamma) when #Delta#phi(#gamma,1^{st} "
166  "jet) > 2.8;#DeltaE_{T}(#gamma,1^{st} jet)/E_{T}(#gamma)",
167  20, -1.0, 1.0);
168  h_jet_count = ibooker.book1D("jet_count",
169  "Number of " + theJetCollectionLabel_.label() + " (p_{T} > " + aString +
170  " GeV);Number of Jets", 8, -0.5, 7.5);
171  h_jet2_pt = ibooker.book1D("jet2_pt",
172  "Jet with 2^{nd} highest p_{T} (from " + theJetCollectionLabel_.label() +
173  ");p_{T}(2^{nd} jet) (GeV)", 20, 0., thePlotPhotonMaxEt_);
174  h_jet2_eta = ibooker.book1D("jet2_eta",
175  "Jet with 2^{nd} highest p_{T} (from " + theJetCollectionLabel_.label() +
176  ");#eta(2^{nd} jet)", 20, -thePlotJetMaxEta_, thePlotJetMaxEta_);
177  h_jet2_ptOverPhotonEt = ibooker.book1D("jet2_ptOverPhotonEt",
178  "p_{T}(2^{nd} highest jet) / E_{T}(#gamma);p_{T}(2^{nd} Jet)/E_{T}(#gamma)",
179  20, 0.0, 4.0);
180  h_deltaPhi_photon_jet2 = ibooker.book1D("deltaPhi_photon_jet2",
181  "#Delta#phi between Highest E_{T} #gamma and 2^{nd} "
182  "highest jet;#Delta#phi(#gamma,2^{nd} jet)", 20, 0, 3.1415926);
183  h_deltaR_jet_jet2 = ibooker.book1D("deltaR_jet_jet2",
184  "#DeltaR between Highest Jet and 2^{nd} Highest;#DeltaR(1^{st} jet,2^{nd} jet)",
185  30, 0, 6.0);
186  h_deltaR_photon_jet2 = ibooker.book1D("deltaR_photon_jet2",
187  "#DeltaR between Highest E_{T} #gamma and 2^{nd} "
188  "jet;#DeltaR(#gamma, 2^{nd} jet)", 30, 0, 6.0);
189 
190  // Photon Et for different jet configurations
191  Float_t bins_et[] = {15, 20, 30, 50, 80};
192  int num_bins_et = 4;
193  h_photon_et_jetcs = ibooker.book1D("photon_et_jetcs",
194  "#gamma with highest E_{T} (#eta(jet)<1.45, "
195  "#eta(#gamma)#eta(jet)>0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
196  h_photon_et_jetco = ibooker.book1D("photon_et_jetco",
197  "#gamma with highest E_{T} (#eta(jet)<1.45, "
198  "#eta(#gamma)#eta(jet)<0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
199  h_photon_et_jetfs = ibooker.book1D("photon_et_jetfs",
200  "#gamma with highest E_{T} (1.55<#eta(jet)<2.5, "
201  "#eta(#gamma)#eta(jet)>0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
202  h_photon_et_jetfo = ibooker.book1D("photon_et_jetfo",
203  "#gamma with highest E_{T} (1.55<#eta(jet)<2.5, "
204  "#eta(#gamma)#eta(jet)<0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
205 
206  auto setSumw2 = [](MonitorElement* me) {
207  if (me->getTH1F()->GetSumw2N() == 0) {
208  me->getTH1F()->Sumw2();
209  }
210  };
211 
212  setSumw2(h_photon_et_jetcs);
213  setSumw2(h_photon_et_jetco);
214  setSumw2(h_photon_et_jetfs);
215  setSumw2(h_photon_et_jetfo);
216 }
217 
218 void QcdPhotonsDQM::analyze(const Event& iEvent, const EventSetup& iSetup) {
219  LogTrace(logTraceName) << "Analysis of event # ";
220 
222  // Did event pass HLT paths?
223  Handle<TriggerResults> HLTresults;
224  iEvent.getByToken(trigTagToken_, HLTresults);
225  if (!HLTresults.isValid()) {
226  // LogWarning("") << ">>> TRIGGER collection does not exist !!!";
227  return;
228  }
229  const edm::TriggerNames& trigNames = iEvent.triggerNames(*HLTresults);
230 
231  bool passed_HLT = false;
232 
233  // See if event passed trigger paths
234  // increment that bin in the trigger plot
235  for (unsigned int i = 0; i < thePlotTheseTriggersToo_.size(); i++) {
236  passed_HLT = false;
237  for (unsigned int ti = 0; (ti < trigNames.size()) && !passed_HLT; ++ti) {
238  size_t pos = trigNames.triggerName(ti).find(thePlotTheseTriggersToo_[i]);
239  if (pos == 0) passed_HLT = HLTresults->accept(ti);
240  }
241  if (passed_HLT) h_triggers_passed->Fill(i);
242  }
243 
244  // grab photons
245  Handle<PhotonCollection> photonCollection;
246  iEvent.getByToken(thePhotonCollectionToken_, photonCollection);
247 
248  // If photon collection is empty, exit
249  if (!photonCollection.isValid()) return;
250 
251  // Quit if the event did not pass the HLT path we care about
252  passed_HLT = false;
253  {
254  // bool found=false;
255  for (unsigned int ti = 0; ti < trigNames.size(); ++ti) {
256  size_t pos = trigNames.triggerName(ti).find(theTriggerPathToPass_);
257  if (pos == 0) {
258  passed_HLT = HLTresults->accept(ti);
259  // found=true;
260  break;
261  }
262  }
263 
264  // Assumption: reco photons are ordered by Et
265  for (PhotonCollection::const_iterator recoPhoton =
266  photonCollection->begin();
267  recoPhoton != photonCollection->end(); recoPhoton++) {
268  // stop looping over photons once we get to too low Et
269  if (recoPhoton->et() < theMinPhotonEt_) break;
270 
271  h_photon_et_beforeCuts->Fill(recoPhoton->et());
272  break; // leading photon only
273  }
274 
275  if (!passed_HLT) {
276  return;
277  }
278  }
279 
281 
282  // std::cout << "\tpassed main trigger (" << theTriggerPathToPass_ << ")" <<
283  // std::endl;
284 
286  // Does event have valid vertex?
287  // Get the primary event vertex
288  Handle<VertexCollection> vertexHandle;
289  iEvent.getByToken(theVertexCollectionToken_, vertexHandle);
290  VertexCollection vertexCollection = *(vertexHandle.product());
291  // double vtx_ndof = -1.0;
292  // double vtx_z = 0.0;
293  // bool vtx_isFake = true;
294  // if (vertexCollection.size()>0) {
295  // vtx_ndof = vertexCollection.begin()->ndof();
296  // vtx_z = vertexCollection.begin()->z();
297  // vtx_isFake = false;
298  //}
299  // if (vtx_isFake || fabs(vtx_z)>15 || vtx_ndof<4) return;
300 
301  int nvvertex = 0;
302  for (unsigned int i = 0; i < vertexCollection.size(); ++i) {
303  if (vertexCollection[i].isValid()) nvvertex++;
304  }
305  if (nvvertex == 0) return;
306 
308 
309  // std::cout << "\tpassed vertex selection" << std::endl;
310 
312  // Did the event pass certain L1 Technical Trigger bits?
313  // It's probably beam halo
314  // TODO: ADD code
316 
317  // For finding spikes
318  Handle<EcalRecHitCollection> EBReducedRecHits;
319  iEvent.getByToken(theBarrelRecHitToken_, EBReducedRecHits);
320  Handle<EcalRecHitCollection> EEReducedRecHits;
321  iEvent.getByToken(theEndcapRecHitToken_, EEReducedRecHits);
322  EcalClusterLazyTools lazyTool(iEvent, iSetup, theBarrelRecHitToken_,
323  theEndcapRecHitToken_);
324 
325  // Find the highest et "decent" photon
326  float photon_et = -9.0;
327  float photon_eta = -9.0;
328  float photon_phi = -9.0;
329  bool photon_passPhotonID = false;
330  bool found_lead_pho = false;
331  int photon_count_bar = 0;
332  int photon_count_end = 0;
333  // False Assumption: reco photons are ordered by Et
334  // find the photon with highest et
335  auto pho_maxet = std::max_element(
336  photonCollection->begin(), photonCollection->end(),
338  const PhotonCollection::value_type& b) { return a.et() < b.et(); });
339  if (pho_maxet != photonCollection->end() &&
340  pho_maxet->et() >= theMinPhotonEt_) {
341  /*
342  // Ignore ECAL Spikes
343  const reco::CaloClusterPtr seed = pho_maxet->superCluster()->seed();
344  DetId id = lazyTool.getMaximum(*seed).first; // Cluster shape variables
345  // float time = -999., outOfTimeChi2 = -999., chi2 = -999.; // UNUSED
346  int flags=-1, severity = -1;
347  const EcalRecHitCollection & rechits = ( pho_maxet->isEB() ?
348  *EBReducedRecHits : *EEReducedRecHits);
349  EcalRecHitCollection::const_iterator it = rechits.find( id );
350  if( it != rechits.end() ) {
351  // time = it->time(); // UNUSED
352  // outOfTimeChi2 = it->outOfTimeChi2(); // UNUSED
353  // chi2 = it->chi2(); // UNUSED
354  flags = it->recoFlag();
355 
356  edm::ESHandle<EcalSeverityLevelAlgo> sevlv;
357  iSetup.get<EcalSeverityLevelAlgoRcd>().get(sevlv);
358  severity = sevlv->severityLevel( id, rechits);
359  }
360  bool isNotSpike = ((pho_maxet->isEB() && (severity!=3 && severity!=4 ) &&
361  (flags != 2) ) || pho_maxet->isEE());
362  if (!isNotSpike) continue; // move on to next photon
363  // END of determining ECAL Spikes
364  */
365 
366  bool pho_current_passPhotonID = false;
367  bool pho_current_isEB = pho_maxet->isEB();
368  bool pho_current_isEE = pho_maxet->isEE();
369 
370  if (pho_current_isEB && (pho_maxet->sigmaIetaIeta() < 0.01 ||
371  pho_maxet->hadronicOverEm() < 0.05)) {
372  // Photon object in barrel passes photon ID
373  pho_current_passPhotonID = true;
374  photon_count_bar++;
375  } else if (pho_current_isEE && (pho_maxet->hadronicOverEm() < 0.05)) {
376  // Photon object in endcap passes photon ID
377  pho_current_passPhotonID = true;
378  photon_count_end++;
379  }
380 
381  if (!found_lead_pho) {
382  found_lead_pho = true;
383  photon_passPhotonID = pho_current_passPhotonID;
384  photon_et = pho_maxet->et();
385  photon_eta = pho_maxet->eta();
386  photon_phi = pho_maxet->phi();
387  }
388  }
389 
390  // If user requires a photon to be found, but none is, return.
391  // theRequirePhotonFound should pretty much always be set to 'True'
392  // except when running on qcd monte carlo just to see the jets.
393  if (theRequirePhotonFound_ &&
394  (!photon_passPhotonID || photon_et < theMinPhotonEt_))
395  return;
396 
398  // Find the highest et jet
399  Handle<View<Jet> > jetCollection;
400  iEvent.getByToken(theJetCollectionToken_, jetCollection);
401  if (!jetCollection.isValid()) return;
402 
403  float jet_pt = -8.0;
404  float jet_eta = -8.0;
405  float jet_phi = -8.0;
406  int jet_count = 0;
407  float jet2_pt = -9.0;
408  float jet2_eta = -9.0;
409  float jet2_phi = -9.0;
410  // Assumption: jets are ordered by Et
411  for (unsigned int i_jet = 0; i_jet < jetCollection->size(); i_jet++) {
412  const Jet* jet = &jetCollection->at(i_jet);
413 
414  float jet_current_pt = jet->pt();
415 
416  // don't care about jets that overlap with the lead photon
417  if (deltaR(jet->eta(), jet->phi(), photon_eta, photon_phi) < 0.5) continue;
418  // stop looping over jets once we get to too low Et
419  if (jet_current_pt < theMinJetPt_) break;
420 
421  jet_count++;
422  if (jet_current_pt > jet_pt) {
423  jet2_pt = jet_pt; // 2nd highest jet get's et from current highest
424  jet2_eta = jet_eta;
425  jet2_phi = jet_phi;
426  jet_pt =
427  jet_current_pt; // current highest jet gets et from the new highest
428  jet_eta = jet->eta();
429  jet_phi = jet->phi();
430  } else if (jet_current_pt > jet2_pt) {
431  jet2_pt = jet_current_pt;
432  jet2_eta = jet->eta();
433  jet2_phi = jet->phi();
434  }
435  }
437 
439  // Fill histograms if a jet found
440  // NOTE: if a photon was required to be found, but wasn't
441  // we wouldn't have made it to this point in the code
442  if (jet_pt > 0.0) {
443 
444  // Photon Plots
445  h_photon_et->Fill(photon_et);
446  h_photon_eta->Fill(photon_eta);
447  h_photon_count_bar->Fill(photon_count_bar);
448  h_photon_count_end->Fill(photon_count_end);
449 
450  // Photon Et hists for different orientations to the jet
451  if (fabs(photon_eta) < 1.45 &&
452  photon_passPhotonID) { // Lead photon is in barrel
453  if (fabs(jet_eta) < 1.45) { // jet is in barrel
454  if (photon_eta * jet_eta > 0) {
455  h_photon_et_jetcs->Fill(photon_et);
456  } else {
457  h_photon_et_jetco->Fill(photon_et);
458  }
459  } else if (jet_eta > 1.55 && jet_eta < 2.5) { // jet is in endcap
460  if (photon_eta * jet_eta > 0) {
461  h_photon_et_jetfs->Fill(photon_et);
462  } else {
463  h_photon_et_jetfo->Fill(photon_et);
464  }
465  }
466  } // END of Lead Photon is in Barrel
467 
468  // Jet Plots
469  h_jet_pt->Fill(jet_pt);
470  h_jet_eta->Fill(jet_eta);
471  h_jet_count->Fill(jet_count);
472  h_deltaPhi_photon_jet->Fill(abs(deltaPhi(photon_phi, jet_phi)));
473  if (abs(deltaPhi(photon_phi, jet_phi)) > 2.8)
474  h_deltaEt_photon_jet->Fill((photon_et - jet_pt) / photon_et);
475 
476  // 2nd Highest Jet Plots
477  if (jet2_pt > 0.0) {
478  h_jet2_pt->Fill(jet2_pt);
479  h_jet2_eta->Fill(jet2_eta);
480  h_jet2_ptOverPhotonEt->Fill(jet2_pt / photon_et);
481  h_deltaPhi_photon_jet2->Fill(abs(deltaPhi(photon_phi, jet2_phi)));
482  h_deltaPhi_jet_jet2->Fill(abs(deltaPhi(jet_phi, jet2_phi)));
483  h_deltaR_jet_jet2->Fill(deltaR(jet_eta, jet_phi, jet2_eta, jet2_phi));
484  h_deltaR_photon_jet2->Fill(
485  deltaR(photon_eta, photon_phi, jet2_eta, jet2_phi));
486  }
487  }
488  // End of Filling histograms
490 }
491 
492 // Local Variables:
493 // show-trailing-whitespace: t
494 // truncate-lines: t
495 // End:
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:215
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
Base class for all types of Jets.
Definition: Jet.h:20
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
virtual double phi() const final
momentum azimuthal angle
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
Strings::size_type size() const
Definition: TriggerNames.cc:39
tuple vertexCollection
QcdPhotonsDQM(const edm::ParameterSet &)
Constructor.
int iEvent
Definition: GenABIO.cc:230
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isValid() const
Definition: HandleBase.h:75
#define LogTrace(id)
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:74
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:276
T const * product() const
Definition: Handle.h:81
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:27
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
virtual ~QcdPhotonsDQM()
Destructor.
virtual double eta() const final
momentum pseudorapidity
void analyze(const edm::Event &, const edm::EventSetup &)
Get the analysis.
Definition: Run.h:43
virtual double pt() const final
transverse momentum