CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
QcdPhotonsDQM Class Reference

#include <QcdPhotonsDQM.h>

Inheritance diagram for QcdPhotonsDQM:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &)
 Get the analysis. More...
 
void beginJob ()
 Inizialize parameters for histo binning. More...
 
void beginRun (const edm::Run &, const edm::EventSetup &)
 
void endJob (void)
 Save the histos. More...
 
void endRun (const edm::Run &, const edm::EventSetup &)
 
 QcdPhotonsDQM (const edm::ParameterSet &)
 Constructor. More...
 
virtual ~QcdPhotonsDQM ()
 Destructor. More...
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Attributes

MonitorElementh_deltaEt_photon_jet
 
MonitorElementh_deltaPhi_jet_jet2
 
MonitorElementh_deltaPhi_photon_jet
 
MonitorElementh_deltaPhi_photon_jet2
 
MonitorElementh_deltaR_jet_jet2
 
MonitorElementh_deltaR_photon_jet2
 
MonitorElementh_jet2_eta
 
MonitorElementh_jet2_pt
 
MonitorElementh_jet2_ptOverPhotonEt
 
MonitorElementh_jet_count
 
MonitorElementh_jet_eta
 
MonitorElementh_jet_pt
 
MonitorElementh_photon_count_bar
 
MonitorElementh_photon_count_end
 
MonitorElementh_photon_et
 
MonitorElementh_photon_et_beforeCuts
 
MonitorElementh_photon_et_jetco
 
MonitorElementh_photon_et_jetcs
 
MonitorElementh_photon_et_jetfo
 
MonitorElementh_photon_et_jetfs
 
MonitorElementh_photon_et_ratio_co_cs
 
MonitorElementh_photon_et_ratio_co_fo
 
MonitorElementh_photon_et_ratio_co_fs
 
MonitorElementh_photon_et_ratio_cs_fo
 
MonitorElementh_photon_et_ratio_cs_fs
 
MonitorElementh_photon_et_ratio_fo_fs
 
MonitorElementh_photon_eta
 
MonitorElementh_triggers_passed
 
HLTConfigProvider hltConfigProvider_
 
bool isValidHltConfig_
 
std::string logTraceName
 
int num_events_in_run
 
edm::InputTag theBarrelRecHitTag_
 
edm::EDGetTokenT
< EcalRecHitCollection
theBarrelRecHitToken_
 
DQMStoretheDbe
 
edm::InputTag theEndcapRecHitTag_
 
edm::EDGetTokenT
< EcalRecHitCollection
theEndcapRecHitToken_
 
edm::InputTag theJetCollectionLabel_
 
edm::EDGetTokenT< edm::View
< reco::Jet > > 
theJetCollectionToken_
 
double theMinJetPt_
 
double theMinPhotonEt_
 
edm::EDGetTokenT
< reco::PhotonCollection
thePhotonCollectionToken_
 
double thePlotJetMaxEta_
 
double thePlotPhotonMaxEt_
 
double thePlotPhotonMaxEta_
 
std::vector< std::string > thePlotTheseTriggersToo_
 
bool theRequirePhotonFound_
 
std::string theTriggerPathToPass_
 
edm::EDGetTokenT
< reco::VertexCollection
theVertexCollectionToken_
 
edm::EDGetTokenT
< edm::TriggerResults
trigTagToken_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

DQM offline for QCD-Photons

Author
Michael B. Anderson, University of Wisconsin Madison

Definition at line 29 of file QcdPhotonsDQM.h.

Constructor & Destructor Documentation

QcdPhotonsDQM::QcdPhotonsDQM ( const edm::ParameterSet parameters)

Constructor.

Definition at line 59 of file QcdPhotonsDQM.cc.

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), and cppFunctionSkipper::operator.

59  {
60  // Get parameters from configuration file
61  theTriggerPathToPass_ = parameters.getParameter<string>("triggerPathToPass");
63  parameters.getParameter<vector<string> >("plotTheseTriggersToo");
64  theJetCollectionLabel_ = parameters.getParameter<InputTag>("jetCollection");
65  trigTagToken_ = consumes<edm::TriggerResults>(
66  parameters.getUntrackedParameter<edm::InputTag>("trigTag"));
67  thePhotonCollectionToken_ = consumes<reco::PhotonCollection>(
68  parameters.getParameter<InputTag>("photonCollection"));
69  theJetCollectionToken_ = consumes<edm::View<reco::Jet> >(
70  parameters.getParameter<InputTag>("jetCollection"));
71  theVertexCollectionToken_ = consumes<reco::VertexCollection>(
72  parameters.getParameter<InputTag>("vertexCollection"));
73  theMinJetPt_ = parameters.getParameter<double>("minJetPt");
74  theMinPhotonEt_ = parameters.getParameter<double>("minPhotonEt");
75  theRequirePhotonFound_ = parameters.getParameter<bool>("requirePhotonFound");
76  thePlotPhotonMaxEt_ = parameters.getParameter<double>("plotPhotonMaxEt");
77  thePlotPhotonMaxEta_ = parameters.getParameter<double>("plotPhotonMaxEta");
78  thePlotJetMaxEta_ = parameters.getParameter<double>("plotJetMaxEta");
79  theBarrelRecHitTag_ = parameters.getParameter<InputTag>("barrelRecHitTag");
80  theEndcapRecHitTag_ = parameters.getParameter<InputTag>("endcapRecHitTag");
81  theBarrelRecHitToken_ = consumes<EcalRecHitCollection>(
82  parameters.getParameter<InputTag>("barrelRecHitTag"));
83  theEndcapRecHitToken_ = consumes<EcalRecHitCollection>(
84  parameters.getParameter<InputTag>("endcapRecHitTag"));
85  // just to initialize
86  isValidHltConfig_ = false;
87 
88  // coverity says...
95  h_jet2_eta = 0;
96  h_jet2_pt = 0;
98  h_jet_count = 0;
99  h_jet_eta = 0;
100  h_jet_pt = 0;
101  h_photon_count_bar = 0;
102  h_photon_count_end = 0;
103  h_photon_et = 0;
105  h_photon_et_jetco = 0;
106  h_photon_et_jetcs = 0;
107  h_photon_et_jetfo = 0;
108  h_photon_et_jetfs = 0;
115  h_photon_eta = 0;
116  h_triggers_passed = 0;
117 
119 }
edm::InputTag theBarrelRecHitTag_
Definition: QcdPhotonsDQM.h:78
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< edm::View< reco::Jet > > theJetCollectionToken_
Definition: QcdPhotonsDQM.h:69
MonitorElement * h_photon_et
Definition: QcdPhotonsDQM.h:88
edm::InputTag theJetCollectionLabel_
Definition: QcdPhotonsDQM.h:66
double thePlotPhotonMaxEta_
Definition: QcdPhotonsDQM.h:75
MonitorElement * h_photon_eta
Definition: QcdPhotonsDQM.h:89
MonitorElement * h_photon_count_bar
Definition: QcdPhotonsDQM.h:90
edm::InputTag theEndcapRecHitTag_
Definition: QcdPhotonsDQM.h:79
edm::EDGetTokenT< edm::TriggerResults > trigTagToken_
Definition: QcdPhotonsDQM.h:67
double thePlotJetMaxEta_
Definition: QcdPhotonsDQM.h:76
MonitorElement * h_photon_et_jetco
MonitorElement * h_photon_et_ratio_co_fo
std::vector< std::string > thePlotTheseTriggersToo_
Definition: QcdPhotonsDQM.h:65
edm::EDGetTokenT< reco::VertexCollection > theVertexCollectionToken_
Definition: QcdPhotonsDQM.h:70
double theMinPhotonEt_
Definition: QcdPhotonsDQM.h:72
double theMinJetPt_
Definition: QcdPhotonsDQM.h:71
MonitorElement * h_deltaR_jet_jet2
MonitorElement * h_deltaR_photon_jet2
MonitorElement * h_photon_et_jetfs
bool isValidHltConfig_
Definition: QcdPhotonsDQM.h:58
MonitorElement * h_deltaEt_photon_jet
Definition: QcdPhotonsDQM.h:97
MonitorElement * h_jet_count
Definition: QcdPhotonsDQM.h:94
MonitorElement * h_photon_et_ratio_cs_fs
MonitorElement * h_photon_et_ratio_cs_fo
MonitorElement * h_jet_eta
Definition: QcdPhotonsDQM.h:93
double thePlotPhotonMaxEt_
Definition: QcdPhotonsDQM.h:74
MonitorElement * h_photon_et_jetfo
MonitorElement * h_jet2_ptOverPhotonEt
Definition: QcdPhotonsDQM.h:98
MonitorElement * h_photon_et_jetcs
MonitorElement * h_jet2_pt
Definition: QcdPhotonsDQM.h:99
MonitorElement * h_photon_et_beforeCuts
Definition: QcdPhotonsDQM.h:87
MonitorElement * h_jet_pt
Definition: QcdPhotonsDQM.h:92
MonitorElement * h_photon_et_ratio_fo_fs
DQMStore * theDbe
Definition: QcdPhotonsDQM.h:55
edm::EDGetTokenT< EcalRecHitCollection > theBarrelRecHitToken_
Definition: QcdPhotonsDQM.h:80
bool theRequirePhotonFound_
Definition: QcdPhotonsDQM.h:73
MonitorElement * h_photon_count_end
Definition: QcdPhotonsDQM.h:91
std::string theTriggerPathToPass_
Definition: QcdPhotonsDQM.h:64
MonitorElement * h_jet2_eta
MonitorElement * h_deltaPhi_photon_jet2
MonitorElement * h_photon_et_ratio_co_fs
MonitorElement * h_triggers_passed
Definition: QcdPhotonsDQM.h:86
edm::EDGetTokenT< reco::PhotonCollection > thePhotonCollectionToken_
Definition: QcdPhotonsDQM.h:68
MonitorElement * h_deltaPhi_jet_jet2
Definition: QcdPhotonsDQM.h:96
MonitorElement * h_deltaPhi_photon_jet
Definition: QcdPhotonsDQM.h:95
edm::EDGetTokenT< EcalRecHitCollection > theEndcapRecHitToken_
Definition: QcdPhotonsDQM.h:81
MonitorElement * h_photon_et_ratio_co_cs
QcdPhotonsDQM::~QcdPhotonsDQM ( )
virtual

Destructor.

Definition at line 121 of file QcdPhotonsDQM.cc.

121 {}

Member Function Documentation

void QcdPhotonsDQM::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
virtual

Get the analysis.

Implements edm::EDAnalyzer.

Definition at line 304 of file QcdPhotonsDQM.cc.

References a, funct::abs(), b, SiPixelRawToDigiRegional_cfi::deltaPhi, deltaR(), reco::LeafCandidate::eta(), edm::Event::getByToken(), i, edm::HandleBase::isValid(), metsig::jet, LogTrace, reco::LeafCandidate::phi(), edm::Handle< T >::product(), reco::LeafCandidate::pt(), edm::TriggerNames::size(), edm::TriggerNames::triggerName(), edm::Event::triggerNames(), trigNames, and GoodVertex_cfg::vertexCollection.

304  {
306 
308  // if( ! isValidHltConfig_ ) return;
309 
310  LogTrace(logTraceName) << "Analysis of event # ";
311 
313  // Did event pass HLT paths?
314  Handle<TriggerResults> HLTresults;
315  iEvent.getByToken(trigTagToken_, HLTresults);
316  if (!HLTresults.isValid()) {
317  // LogWarning("") << ">>> TRIGGER collection does not exist !!!";
318  return;
319  }
320  const edm::TriggerNames& trigNames = iEvent.triggerNames(*HLTresults);
321 
322  bool passed_HLT = false;
323 
324  // See if event passed trigger paths
325  // increment that bin in the trigger plot
326  for (unsigned int i = 0; i < thePlotTheseTriggersToo_.size(); i++) {
327  passed_HLT = false;
328  for (unsigned int ti = 0; (ti < trigNames.size()) && !passed_HLT; ++ti) {
329  size_t pos = trigNames.triggerName(ti).find(thePlotTheseTriggersToo_[i]);
330  if (pos == 0) passed_HLT = HLTresults->accept(ti);
331  }
332  if (passed_HLT) h_triggers_passed->Fill(i);
333  }
334 
335  // grab photons
336  Handle<PhotonCollection> photonCollection;
337  iEvent.getByToken(thePhotonCollectionToken_, photonCollection);
338 
339  // If photon collection is empty, exit
340  if (!photonCollection.isValid()) return;
341 
342  // Quit if the event did not pass the HLT path we care about
343  passed_HLT = false;
344  {
345  // bool found=false;
346  for (unsigned int ti = 0; ti < trigNames.size(); ++ti) {
347  size_t pos = trigNames.triggerName(ti).find(theTriggerPathToPass_);
348  if (pos == 0) {
349  passed_HLT = HLTresults->accept(ti);
350  // found=true;
351  break;
352  }
353  }
354 
355  // Assumption: reco photons are ordered by Et
356  for (PhotonCollection::const_iterator recoPhoton =
357  photonCollection->begin();
358  recoPhoton != photonCollection->end(); recoPhoton++) {
359  // stop looping over photons once we get to too low Et
360  if (recoPhoton->et() < theMinPhotonEt_) break;
361 
362  h_photon_et_beforeCuts->Fill(recoPhoton->et());
363  break; // leading photon only
364  }
365 
366  if (!passed_HLT) {
367  return;
368  }
369  }
370 
372 
373  // std::cout << "\tpassed main trigger (" << theTriggerPathToPass_ << ")" <<
374  // std::endl;
375 
377  // Does event have valid vertex?
378  // Get the primary event vertex
379  Handle<VertexCollection> vertexHandle;
380  iEvent.getByToken(theVertexCollectionToken_, vertexHandle);
381  VertexCollection vertexCollection = *(vertexHandle.product());
382  // double vtx_ndof = -1.0;
383  // double vtx_z = 0.0;
384  // bool vtx_isFake = true;
385  // if (vertexCollection.size()>0) {
386  // vtx_ndof = vertexCollection.begin()->ndof();
387  // vtx_z = vertexCollection.begin()->z();
388  // vtx_isFake = false;
389  //}
390  // if (vtx_isFake || fabs(vtx_z)>15 || vtx_ndof<4) return;
391 
392  int nvvertex = 0;
393  for (unsigned int i = 0; i < vertexCollection.size(); ++i) {
394  if (vertexCollection[i].isValid()) nvvertex++;
395  }
396  if (nvvertex == 0) return;
397 
399 
400  // std::cout << "\tpassed vertex selection" << std::endl;
401 
403  // Did the event pass certain L1 Technical Trigger bits?
404  // It's probably beam halo
405  // TODO: ADD code
407 
408  // For finding spikes
409  Handle<EcalRecHitCollection> EBReducedRecHits;
410  iEvent.getByToken(theBarrelRecHitToken_, EBReducedRecHits);
411  Handle<EcalRecHitCollection> EEReducedRecHits;
412  iEvent.getByToken(theEndcapRecHitToken_, EEReducedRecHits);
413  EcalClusterLazyTools lazyTool(iEvent, iSetup, theBarrelRecHitToken_,
415 
416  // Find the highest et "decent" photon
417  float photon_et = -9.0;
418  float photon_eta = -9.0;
419  float photon_phi = -9.0;
420  bool photon_passPhotonID = false;
421  bool found_lead_pho = false;
422  int photon_count_bar = 0;
423  int photon_count_end = 0;
424  // False Assumption: reco photons are ordered by Et
425  // find the photon with highest et
426  auto pho_maxet = std::max_element(
427  photonCollection->begin(), photonCollection->end(),
429  const PhotonCollection::value_type& b) { return a.et() < b.et(); });
430  if (pho_maxet != photonCollection->end() &&
431  pho_maxet->et() >= theMinPhotonEt_) {
432  /*
433  // Ignore ECAL Spikes
434  const reco::CaloClusterPtr seed = pho_maxet->superCluster()->seed();
435  DetId id = lazyTool.getMaximum(*seed).first; // Cluster shape variables
436  // float time = -999., outOfTimeChi2 = -999., chi2 = -999.; // UNUSED
437  int flags=-1, severity = -1;
438  const EcalRecHitCollection & rechits = ( pho_maxet->isEB() ?
439  *EBReducedRecHits : *EEReducedRecHits);
440  EcalRecHitCollection::const_iterator it = rechits.find( id );
441  if( it != rechits.end() ) {
442  // time = it->time(); // UNUSED
443  // outOfTimeChi2 = it->outOfTimeChi2(); // UNUSED
444  // chi2 = it->chi2(); // UNUSED
445  flags = it->recoFlag();
446 
447  edm::ESHandle<EcalSeverityLevelAlgo> sevlv;
448  iSetup.get<EcalSeverityLevelAlgoRcd>().get(sevlv);
449  severity = sevlv->severityLevel( id, rechits);
450  }
451  bool isNotSpike = ((pho_maxet->isEB() && (severity!=3 && severity!=4 ) &&
452  (flags != 2) ) || pho_maxet->isEE());
453  if (!isNotSpike) continue; // move on to next photon
454  // END of determining ECAL Spikes
455  */
456 
457  bool pho_current_passPhotonID = false;
458  bool pho_current_isEB = pho_maxet->isEB();
459  bool pho_current_isEE = pho_maxet->isEE();
460 
461  if (pho_current_isEB && (pho_maxet->sigmaIetaIeta() < 0.01 ||
462  pho_maxet->hadronicOverEm() < 0.05)) {
463  // Photon object in barrel passes photon ID
464  pho_current_passPhotonID = true;
465  photon_count_bar++;
466  } else if (pho_current_isEE && (pho_maxet->hadronicOverEm() < 0.05)) {
467  // Photon object in endcap passes photon ID
468  pho_current_passPhotonID = true;
469  photon_count_end++;
470  }
471 
472  if (!found_lead_pho) {
473  found_lead_pho = true;
474  photon_passPhotonID = pho_current_passPhotonID;
475  photon_et = pho_maxet->et();
476  photon_eta = pho_maxet->eta();
477  photon_phi = pho_maxet->phi();
478  }
479  }
480 
481  // If user requires a photon to be found, but none is, return.
482  // theRequirePhotonFound should pretty much always be set to 'True'
483  // except when running on qcd monte carlo just to see the jets.
485  (!photon_passPhotonID || photon_et < theMinPhotonEt_))
486  return;
487 
489  // Find the highest et jet
490  Handle<View<Jet> > jetCollection;
491  iEvent.getByToken(theJetCollectionToken_, jetCollection);
492  if (!jetCollection.isValid()) return;
493 
494  float jet_pt = -8.0;
495  float jet_eta = -8.0;
496  float jet_phi = -8.0;
497  int jet_count = 0;
498  float jet2_pt = -9.0;
499  float jet2_eta = -9.0;
500  float jet2_phi = -9.0;
501  // Assumption: jets are ordered by Et
502  for (unsigned int i_jet = 0; i_jet < jetCollection->size(); i_jet++) {
503  const Jet* jet = &jetCollection->at(i_jet);
504 
505  float jet_current_pt = jet->pt();
506 
507  // don't care about jets that overlap with the lead photon
508  if (deltaR(jet->eta(), jet->phi(), photon_eta, photon_phi) < 0.5) continue;
509  // stop looping over jets once we get to too low Et
510  if (jet_current_pt < theMinJetPt_) break;
511 
512  jet_count++;
513  if (jet_current_pt > jet_pt) {
514  jet2_pt = jet_pt; // 2nd highest jet get's et from current highest
515  jet2_eta = jet_eta;
516  jet2_phi = jet_phi;
517  jet_pt =
518  jet_current_pt; // current highest jet gets et from the new highest
519  jet_eta = jet->eta();
520  jet_phi = jet->phi();
521  } else if (jet_current_pt > jet2_pt) {
522  jet2_pt = jet_current_pt;
523  jet2_eta = jet->eta();
524  jet2_phi = jet->phi();
525  }
526  }
528 
530  // Fill histograms if a jet found
531  // NOTE: if a photon was required to be found, but wasn't
532  // we wouldn't have made it to this point in the code
533  if (jet_pt > 0.0) {
534 
535  // Photon Plots
536  h_photon_et->Fill(photon_et);
537  h_photon_eta->Fill(photon_eta);
538  h_photon_count_bar->Fill(photon_count_bar);
539  h_photon_count_end->Fill(photon_count_end);
540 
541  // Photon Et hists for different orientations to the jet
542  if (fabs(photon_eta) < 1.45 &&
543  photon_passPhotonID) { // Lead photon is in barrel
544  if (fabs(jet_eta) < 1.45) { // jet is in barrel
545  if (photon_eta * jet_eta > 0) {
546  h_photon_et_jetcs->Fill(photon_et);
547  } else {
548  h_photon_et_jetco->Fill(photon_et);
549  }
550  } else if (jet_eta > 1.55 && jet_eta < 2.5) { // jet is in endcap
551  if (photon_eta * jet_eta > 0) {
552  h_photon_et_jetfs->Fill(photon_et);
553  } else {
554  h_photon_et_jetfo->Fill(photon_et);
555  }
556  }
557  } // END of Lead Photon is in Barrel
558 
559  // Jet Plots
560  h_jet_pt->Fill(jet_pt);
561  h_jet_eta->Fill(jet_eta);
562  h_jet_count->Fill(jet_count);
563  h_deltaPhi_photon_jet->Fill(abs(deltaPhi(photon_phi, jet_phi)));
564  if (abs(deltaPhi(photon_phi, jet_phi)) > 2.8)
565  h_deltaEt_photon_jet->Fill((photon_et - jet_pt) / photon_et);
566 
567  // 2nd Highest Jet Plots
568  if (jet2_pt > 0.0) {
569  h_jet2_pt->Fill(jet2_pt);
570  h_jet2_eta->Fill(jet2_eta);
571  h_jet2_ptOverPhotonEt->Fill(jet2_pt / photon_et);
572  h_deltaPhi_photon_jet2->Fill(abs(deltaPhi(photon_phi, jet2_phi)));
573  h_deltaPhi_jet_jet2->Fill(abs(deltaPhi(jet_phi, jet2_phi)));
574  h_deltaR_jet_jet2->Fill(deltaR(jet_eta, jet_phi, jet2_eta, jet2_phi));
576  deltaR(photon_eta, photon_phi, jet2_eta, jet2_phi));
577  }
578  }
579  // End of Filling histograms
581 }
edm::EDGetTokenT< edm::View< reco::Jet > > theJetCollectionToken_
Definition: QcdPhotonsDQM.h:69
int i
Definition: DBlmapReader.cc:9
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:199
MonitorElement * h_photon_et
Definition: QcdPhotonsDQM.h:88
MonitorElement * h_photon_eta
Definition: QcdPhotonsDQM.h:89
MonitorElement * h_photon_count_bar
Definition: QcdPhotonsDQM.h:90
edm::EDGetTokenT< edm::TriggerResults > trigTagToken_
Definition: QcdPhotonsDQM.h:67
MonitorElement * h_photon_et_jetco
std::vector< std::string > thePlotTheseTriggersToo_
Definition: QcdPhotonsDQM.h:65
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
edm::EDGetTokenT< reco::VertexCollection > theVertexCollectionToken_
Definition: QcdPhotonsDQM.h:70
double theMinPhotonEt_
Definition: QcdPhotonsDQM.h:72
double theMinJetPt_
Definition: QcdPhotonsDQM.h:71
MonitorElement * h_deltaR_jet_jet2
Base class for all types of Jets.
Definition: Jet.h:20
std::string logTraceName
Definition: QcdPhotonsDQM.h:61
MonitorElement * h_deltaR_photon_jet2
MonitorElement * h_photon_et_jetfs
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
MonitorElement * h_deltaEt_photon_jet
Definition: QcdPhotonsDQM.h:97
Strings::size_type size() const
Definition: TriggerNames.cc:39
MonitorElement * h_jet_count
Definition: QcdPhotonsDQM.h:94
MonitorElement * h_jet_eta
Definition: QcdPhotonsDQM.h:93
tuple vertexCollection
void Fill(long long x)
MonitorElement * h_photon_et_jetfo
MonitorElement * h_jet2_ptOverPhotonEt
Definition: QcdPhotonsDQM.h:98
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
MonitorElement * h_photon_et_jetcs
MonitorElement * h_jet2_pt
Definition: QcdPhotonsDQM.h:99
MonitorElement * h_photon_et_beforeCuts
Definition: QcdPhotonsDQM.h:87
MonitorElement * h_jet_pt
Definition: QcdPhotonsDQM.h:92
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< EcalRecHitCollection > theBarrelRecHitToken_
Definition: QcdPhotonsDQM.h:80
bool isValid() const
Definition: HandleBase.h:76
Container::value_type value_type
#define LogTrace(id)
bool theRequirePhotonFound_
Definition: QcdPhotonsDQM.h:73
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
MonitorElement * h_photon_count_end
Definition: QcdPhotonsDQM.h:91
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:74
std::string theTriggerPathToPass_
Definition: QcdPhotonsDQM.h:64
MonitorElement * h_jet2_eta
MonitorElement * h_deltaPhi_photon_jet2
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
MonitorElement * h_triggers_passed
Definition: QcdPhotonsDQM.h:86
edm::EDGetTokenT< reco::PhotonCollection > thePhotonCollectionToken_
Definition: QcdPhotonsDQM.h:68
MonitorElement * h_deltaPhi_jet_jet2
Definition: QcdPhotonsDQM.h:96
MonitorElement * h_deltaPhi_photon_jet
Definition: QcdPhotonsDQM.h:95
virtual float pt() const GCC11_FINAL
transverse momentum
edm::EDGetTokenT< EcalRecHitCollection > theEndcapRecHitToken_
Definition: QcdPhotonsDQM.h:81
void QcdPhotonsDQM::beginJob ( void  )
virtual

Inizialize parameters for histo binning.

Reimplemented from edm::EDAnalyzer.

Definition at line 123 of file QcdPhotonsDQM.cc.

References i, LogTrace, and AlCaHLTBitMon_QueryRunRegistry::string.

123  {
124 
125  logTraceName = "QcdPhotonAnalyzer";
126 
127  LogTrace(logTraceName) << "Parameters initialization";
128 
129  theDbe->setCurrentFolder(
130  "Physics/QcdPhotons"); // Use folder with name of PAG
131 
132  std::stringstream aStringStream;
133  std::string aString;
134  aStringStream << theMinJetPt_;
135  aString = aStringStream.str();
136 
137  // Monitor of triggers passed
138  int numOfTriggersToMonitor = thePlotTheseTriggersToo_.size();
140  theDbe->book1D("triggers_passed", "Events passing these trigger paths",
141  numOfTriggersToMonitor, 0, numOfTriggersToMonitor);
142  for (int i = 0; i < numOfTriggersToMonitor; i++) {
144  }
145 
146  // Keep the number of plots and number of bins to a minimum!
147  h_photon_et_beforeCuts = theDbe->book1D(
148  "photon_et_beforeCuts", "#gamma with highest E_{T};E_{T}(#gamma) (GeV)",
149  20, 0., thePlotPhotonMaxEt_);
150  h_photon_et = theDbe->book1D("photon_et",
151  "#gamma with highest E_{T};E_{T}(#gamma) (GeV)",
152  20, 0., thePlotPhotonMaxEt_);
153  h_photon_eta =
154  theDbe->book1D("photon_eta", "#gamma with highest E_{T};#eta(#gamma)", 40,
156  h_photon_count_bar = theDbe->book1D(
157  "photon_count_bar",
158  "Number of #gamma's passing selection (Barrel);Number of #gamma's", 8,
159  -0.5, 7.5);
160  h_photon_count_end = theDbe->book1D(
161  "photon_count_end",
162  "Number of #gamma's passing selection (Endcap);Number of #gamma's", 8,
163  -0.5, 7.5);
164 
165  h_jet_pt = theDbe->book1D("jet_pt", "Jet with highest p_{T} (from " +
167  ");p_{T}(1^{st} jet) (GeV)",
168  20, 0., thePlotPhotonMaxEt_);
169  h_jet_eta = theDbe->book1D("jet_eta", "Jet with highest p_{T} (from " +
171  ");#eta(1^{st} jet)",
174  theDbe->book1D("deltaPhi_photon_jet",
175  "#Delta#phi between Highest E_{T} #gamma and "
176  "jet;#Delta#phi(#gamma,1^{st} jet)",
177  20, 0, 3.1415926);
179  theDbe->book1D("deltaPhi_jet_jet2",
180  "#Delta#phi between Highest E_{T} jet and 2^{nd} "
181  "jet;#Delta#phi(1^{st} jet,2^{nd} jet)",
182  20, 0, 3.1415926);
183  h_deltaEt_photon_jet = theDbe->book1D(
184  "deltaEt_photon_jet",
185  "(E_{T}(#gamma)-p_{T}(jet))/E_{T}(#gamma) when #Delta#phi(#gamma,1^{st} "
186  "jet) > 2.8;#DeltaE_{T}(#gamma,1^{st} jet)/E_{T}(#gamma)",
187  20, -1.0, 1.0);
188  h_jet_count = theDbe->book1D(
189  "jet_count", "Number of " + theJetCollectionLabel_.label() +
190  " (p_{T} > " + aString + " GeV);Number of Jets",
191  8, -0.5, 7.5);
192  h_jet2_pt = theDbe->book1D("jet2_pt", "Jet with 2^{nd} highest p_{T} (from " +
194  ");p_{T}(2^{nd} jet) (GeV)",
195  20, 0., thePlotPhotonMaxEt_);
196  h_jet2_eta = theDbe->book1D(
197  "jet2_eta", "Jet with 2^{nd} highest p_{T} (from " +
198  theJetCollectionLabel_.label() + ");#eta(2^{nd} jet)",
201  theDbe->book1D("jet2_ptOverPhotonEt",
202  "p_{T}(2^{nd} highest jet) / E_{T}(#gamma);p_{T}(2^{nd} "
203  "Jet)/E_{T}(#gamma)",
204  20, 0.0, 4.0);
206  theDbe->book1D("deltaPhi_photon_jet2",
207  "#Delta#phi between Highest E_{T} #gamma and 2^{nd} "
208  "highest jet;#Delta#phi(#gamma,2^{nd} jet)",
209  20, 0, 3.1415926);
210  h_deltaR_jet_jet2 = theDbe->book1D("deltaR_jet_jet2",
211  "#DeltaR between Highest Jet and 2^{nd} "
212  "Highest;#DeltaR(1^{st} jet,2^{nd} jet)",
213  30, 0, 6.0);
215  theDbe->book1D("deltaR_photon_jet2",
216  "#DeltaR between Highest E_{T} #gamma and 2^{nd} "
217  "jet;#DeltaR(#gamma, 2^{nd} jet)",
218  30, 0, 6.0);
219 
220  // Photon Et for different jet configurations
221  Float_t bins_et[] = {15, 20, 30, 50, 80};
222  int num_bins_et = 4;
224  theDbe->book1D("photon_et_jetcs",
225  "#gamma with highest E_{T} (#eta(jet)<1.45, "
226  "#eta(#gamma)#eta(jet)>0);E_{T}(#gamma) (GeV)",
227  num_bins_et, bins_et);
229  theDbe->book1D("photon_et_jetco",
230  "#gamma with highest E_{T} (#eta(jet)<1.45, "
231  "#eta(#gamma)#eta(jet)<0);E_{T}(#gamma) (GeV)",
232  num_bins_et, bins_et);
234  theDbe->book1D("photon_et_jetfs",
235  "#gamma with highest E_{T} (1.55<#eta(jet)<2.5, "
236  "#eta(#gamma)#eta(jet)>0);E_{T}(#gamma) (GeV)",
237  num_bins_et, bins_et);
239  theDbe->book1D("photon_et_jetfo",
240  "#gamma with highest E_{T} (1.55<#eta(jet)<2.5, "
241  "#eta(#gamma)#eta(jet)<0);E_{T}(#gamma) (GeV)",
242  num_bins_et, bins_et);
243  h_photon_et_jetcs->getTH1F()->Sumw2();
244  h_photon_et_jetco->getTH1F()->Sumw2();
245  h_photon_et_jetfs->getTH1F()->Sumw2();
246  h_photon_et_jetfo->getTH1F()->Sumw2();
247  // Ratio of the above Photon Et distributions
249  "photon_et_ratio_00_co_cs",
250  "D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)<0) / D(|#eta(jet)|<1.45, "
251  "#eta(jet)*#eta(#gamma)>0);E_{T}(#gamma) (GeV); ratio",
252  num_bins_et, bins_et);
254  theDbe->book1D("photon_et_ratio_01_fo_fs",
255  "D(1.55<|#eta(jet)|<2.6, #eta(jet)*#eta(#gamma)<0) / "
256  "D(1.55<|#eta(jet)|<2.6, "
257  "#eta(jet)*#eta(#gamma)>0);E_{T}(#gamma) (GeV); ratio",
258  num_bins_et, bins_et);
260  "photon_et_ratio_02_cs_fs",
261  "D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)>0) / D(1.55<|#eta(jet)|<2.6, "
262  "#eta(jet)*#eta(#gamma)>0);E_{T}(#gamma) (GeV); ratio",
263  num_bins_et, bins_et);
265  "photon_et_ratio_03_co_fs",
266  "D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)<0) / D(1.55<|#eta(jet)|<2.6, "
267  "#eta(jet)*#eta(#gamma)>0);E_{T}(#gamma) (GeV); ratio",
268  num_bins_et, bins_et);
270  "photon_et_ratio_04_cs_fo",
271  "D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)>0) / D(1.55<|#eta(jet)|<2.6, "
272  "#eta(jet)*#eta(#gamma)<0);E_{T}(#gamma) (GeV); ratio",
273  num_bins_et, bins_et);
275  "photon_et_ratio_05_co_fo",
276  "D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)<0) / D(1.55<|#eta(jet)|<2.6, "
277  "#eta(jet)*#eta(#gamma)<0);E_{T}(#gamma) (GeV); ratio",
278  num_bins_et, bins_et);
279  h_photon_et_ratio_co_cs->getTH1F()->Sumw2();
280  h_photon_et_ratio_fo_fs->getTH1F()->Sumw2();
281  h_photon_et_ratio_cs_fs->getTH1F()->Sumw2();
282  h_photon_et_ratio_co_fs->getTH1F()->Sumw2();
283  h_photon_et_ratio_cs_fo->getTH1F()->Sumw2();
284  h_photon_et_ratio_co_fo->getTH1F()->Sumw2();
285 }
int i
Definition: DBlmapReader.cc:9
MonitorElement * h_photon_et
Definition: QcdPhotonsDQM.h:88
edm::InputTag theJetCollectionLabel_
Definition: QcdPhotonsDQM.h:66
double thePlotPhotonMaxEta_
Definition: QcdPhotonsDQM.h:75
MonitorElement * h_photon_eta
Definition: QcdPhotonsDQM.h:89
MonitorElement * h_photon_count_bar
Definition: QcdPhotonsDQM.h:90
double thePlotJetMaxEta_
Definition: QcdPhotonsDQM.h:76
MonitorElement * h_photon_et_jetco
MonitorElement * h_photon_et_ratio_co_fo
std::vector< std::string > thePlotTheseTriggersToo_
Definition: QcdPhotonsDQM.h:65
double theMinJetPt_
Definition: QcdPhotonsDQM.h:71
MonitorElement * h_deltaR_jet_jet2
std::string logTraceName
Definition: QcdPhotonsDQM.h:61
MonitorElement * h_deltaR_photon_jet2
MonitorElement * h_photon_et_jetfs
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)
MonitorElement * h_deltaEt_photon_jet
Definition: QcdPhotonsDQM.h:97
MonitorElement * h_jet_count
Definition: QcdPhotonsDQM.h:94
MonitorElement * h_photon_et_ratio_cs_fs
MonitorElement * h_photon_et_ratio_cs_fo
MonitorElement * h_jet_eta
Definition: QcdPhotonsDQM.h:93
double thePlotPhotonMaxEt_
Definition: QcdPhotonsDQM.h:74
MonitorElement * h_photon_et_jetfo
MonitorElement * h_jet2_ptOverPhotonEt
Definition: QcdPhotonsDQM.h:98
MonitorElement * h_photon_et_jetcs
MonitorElement * h_jet2_pt
Definition: QcdPhotonsDQM.h:99
MonitorElement * h_photon_et_beforeCuts
Definition: QcdPhotonsDQM.h:87
MonitorElement * h_jet_pt
Definition: QcdPhotonsDQM.h:92
MonitorElement * h_photon_et_ratio_fo_fs
DQMStore * theDbe
Definition: QcdPhotonsDQM.h:55
#define LogTrace(id)
MonitorElement * h_photon_count_end
Definition: QcdPhotonsDQM.h:91
MonitorElement * h_jet2_eta
MonitorElement * h_deltaPhi_photon_jet2
TH1F * getTH1F(void) const
std::string const & label() const
Definition: InputTag.h:42
MonitorElement * h_photon_et_ratio_co_fs
MonitorElement * h_triggers_passed
Definition: QcdPhotonsDQM.h:86
MonitorElement * h_deltaPhi_jet_jet2
Definition: QcdPhotonsDQM.h:96
MonitorElement * h_deltaPhi_photon_jet
Definition: QcdPhotonsDQM.h:95
MonitorElement * h_photon_et_ratio_co_cs
void QcdPhotonsDQM::beginRun ( const edm::Run iRun,
const edm::EventSetup iSet 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 290 of file QcdPhotonsDQM.cc.

291  {
292 
293  // passed as parameter to HLTConfigProvider::init(), not yet used
294  bool isConfigChanged = false;
295 
296  // isValidHltConfig_ could be used to short-circuit analyze() in case of
297  // problems
299  hltConfigProvider_.init(iRun, iSet, "HLT", isConfigChanged);
300 
301  num_events_in_run = 0;
302 }
HLTConfigProvider hltConfigProvider_
Definition: QcdPhotonsDQM.h:57
bool isValidHltConfig_
Definition: QcdPhotonsDQM.h:58
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
void QcdPhotonsDQM::endJob ( void  )
virtual

Save the histos.

Reimplemented from edm::EDAnalyzer.

Definition at line 583 of file QcdPhotonsDQM.cc.

583 {}
void QcdPhotonsDQM::endRun ( const edm::Run run,
const edm::EventSetup es 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 585 of file QcdPhotonsDQM.cc.

585  {
586  if (num_events_in_run > 0) {
588  }
601 }
MonitorElement * h_photon_et_jetco
MonitorElement * h_photon_et_ratio_co_fo
MonitorElement * h_photon_et_jetfs
MonitorElement * h_photon_et_ratio_cs_fs
MonitorElement * h_photon_et_ratio_cs_fo
MonitorElement * h_photon_et_jetfo
MonitorElement * h_photon_et_jetcs
MonitorElement * h_photon_et_ratio_fo_fs
TH1F * getTH1F(void) const
MonitorElement * h_photon_et_ratio_co_fs
MonitorElement * h_triggers_passed
Definition: QcdPhotonsDQM.h:86
MonitorElement * h_photon_et_ratio_co_cs

Member Data Documentation

MonitorElement* QcdPhotonsDQM::h_deltaEt_photon_jet
private

Definition at line 97 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaPhi_jet_jet2
private

Definition at line 96 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaPhi_photon_jet
private

Definition at line 95 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaPhi_photon_jet2
private

Definition at line 101 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaR_jet_jet2
private

Definition at line 102 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaR_photon_jet2
private

Definition at line 103 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet2_eta
private

Definition at line 100 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet2_pt
private

Definition at line 99 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet2_ptOverPhotonEt
private

Definition at line 98 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet_count
private

Definition at line 94 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet_eta
private

Definition at line 93 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet_pt
private

Definition at line 92 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_count_bar
private

Definition at line 90 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_count_end
private

Definition at line 91 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et
private

Definition at line 88 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_beforeCuts
private

Definition at line 87 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetco
private

Definition at line 106 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetcs
private

Definition at line 105 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetfo
private

Definition at line 108 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetfs
private

Definition at line 107 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_co_cs
private

Definition at line 110 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_co_fo
private

Definition at line 115 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_co_fs
private

Definition at line 113 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_cs_fo
private

Definition at line 114 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_cs_fs
private

Definition at line 112 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_fo_fs
private

Definition at line 111 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_eta
private

Definition at line 89 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_triggers_passed
private

Definition at line 86 of file QcdPhotonsDQM.h.

HLTConfigProvider QcdPhotonsDQM::hltConfigProvider_
private

Definition at line 57 of file QcdPhotonsDQM.h.

bool QcdPhotonsDQM::isValidHltConfig_
private

Definition at line 58 of file QcdPhotonsDQM.h.

std::string QcdPhotonsDQM::logTraceName
private

Definition at line 61 of file QcdPhotonsDQM.h.

int QcdPhotonsDQM::num_events_in_run
private

Definition at line 83 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::theBarrelRecHitTag_
private

Definition at line 78 of file QcdPhotonsDQM.h.

edm::EDGetTokenT<EcalRecHitCollection> QcdPhotonsDQM::theBarrelRecHitToken_
private

Definition at line 80 of file QcdPhotonsDQM.h.

DQMStore* QcdPhotonsDQM::theDbe
private

Definition at line 55 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::theEndcapRecHitTag_
private

Definition at line 79 of file QcdPhotonsDQM.h.

edm::EDGetTokenT<EcalRecHitCollection> QcdPhotonsDQM::theEndcapRecHitToken_
private

Definition at line 81 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::theJetCollectionLabel_
private

Definition at line 66 of file QcdPhotonsDQM.h.

edm::EDGetTokenT<edm::View<reco::Jet> > QcdPhotonsDQM::theJetCollectionToken_
private

Definition at line 69 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::theMinJetPt_
private

Definition at line 71 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::theMinPhotonEt_
private

Definition at line 72 of file QcdPhotonsDQM.h.

edm::EDGetTokenT<reco::PhotonCollection> QcdPhotonsDQM::thePhotonCollectionToken_
private

Definition at line 68 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::thePlotJetMaxEta_
private

Definition at line 76 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::thePlotPhotonMaxEt_
private

Definition at line 74 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::thePlotPhotonMaxEta_
private

Definition at line 75 of file QcdPhotonsDQM.h.

std::vector<std::string> QcdPhotonsDQM::thePlotTheseTriggersToo_
private

Definition at line 65 of file QcdPhotonsDQM.h.

bool QcdPhotonsDQM::theRequirePhotonFound_
private

Definition at line 73 of file QcdPhotonsDQM.h.

std::string QcdPhotonsDQM::theTriggerPathToPass_
private

Definition at line 64 of file QcdPhotonsDQM.h.

edm::EDGetTokenT<reco::VertexCollection> QcdPhotonsDQM::theVertexCollectionToken_
private

Definition at line 70 of file QcdPhotonsDQM.h.

edm::EDGetTokenT<edm::TriggerResults> QcdPhotonsDQM::trigTagToken_
private

Definition at line 67 of file QcdPhotonsDQM.h.