CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
QcdPhotonsDQM Class Reference

#include <QcdPhotonsDQM.h>

Inheritance diagram for QcdPhotonsDQM:
DQMEDAnalyzer edm::one::EDProducer< edm::EndRunProducer, edm::one::WatchRuns, edm::EndLuminosityBlockProducer, edm::one::WatchLuminosityBlocks, edm::Accumulator > edm::one::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 Get the analysis. More...
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
 QcdPhotonsDQM (const edm::ParameterSet &)
 Constructor. More...
 
 ~QcdPhotonsDQM () override
 Destructor. More...
 
- Public Member Functions inherited from DQMEDAnalyzer
void accumulate (edm::Event const &event, edm::EventSetup const &setup) final
 
virtual void analyze (edm::Event const &, edm::EventSetup const &)
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
virtual void dqmBeginLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &)
 
virtual void dqmBeginRun (edm::Run const &, edm::EventSetup const &)
 
 DQMEDAnalyzer ()
 
virtual void dqmEndLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &)
 
virtual void dqmEndRun (edm::Run const &, edm::EventSetup const &)
 
void endLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &) final
 
void endLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &, edm::EventSetup const &) final
 
void endRunProduce (edm::Run &run, edm::EventSetup const &setup) final
 
virtual bool getCanSaveByLumi ()
 
- Public Member Functions inherited from edm::one::EDProducer< edm::EndRunProducer, edm::one::WatchRuns, edm::EndLuminosityBlockProducer, edm::one::WatchLuminosityBlocks, edm::Accumulator >
 EDProducer ()=default
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
- Public Member Functions inherited from edm::one::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () noexcept(false) override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

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_eta
 
MonitorElementh_triggers_passed
 
std::string logTraceName
 
edm::InputTag theBarrelRecHitTag_
 
edm::EDGetTokenT< EcalRecHitCollectiontheBarrelRecHitToken_
 
edm::InputTag theEndcapRecHitTag_
 
edm::EDGetTokenT< EcalRecHitCollectiontheEndcapRecHitToken_
 
edm::InputTag theJetCollectionLabel_
 
edm::EDGetTokenT< edm::View< reco::Jet > > theJetCollectionToken_
 
double theMinJetPt_
 
double theMinPhotonEt_
 
edm::EDGetTokenT< reco::PhotonCollectionthePhotonCollectionToken_
 
double thePlotJetMaxEta_
 
double thePlotPhotonMaxEt_
 
double thePlotPhotonMaxEta_
 
std::vector< std::string > thePlotTheseTriggersToo_
 
bool theRequirePhotonFound_
 
std::string theTriggerPathToPass_
 
edm::EDGetTokenT< reco::VertexCollectiontheVertexCollectionToken_
 
edm::EDGetTokenT< edm::TriggerResultstrigTagToken_
 

Additional Inherited Members

- Public Types inherited from DQMEDAnalyzer
typedef dqm::reco::DQMStore DQMStore
 
typedef dqm::reco::MonitorElement MonitorElement
 
- Public Types inherited from edm::one::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::ProducerBase
ProducesCollector producesCollector ()
 
- 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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
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)
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 

Detailed Description

DQM offline for QCD-Photons

Author
Michael B. Anderson, University of Wisconsin Madison

Definition at line 26 of file QcdPhotonsDQM.h.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 57 of file QcdPhotonsDQM.cc.

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

57  {
58  // Get parameters from configuration file
59  theTriggerPathToPass_ = parameters.getParameter<string>("triggerPathToPass");
60  thePlotTheseTriggersToo_ = parameters.getParameter<vector<string> >("plotTheseTriggersToo");
61  theJetCollectionLabel_ = parameters.getParameter<InputTag>("jetCollection");
62  trigTagToken_ = consumes<edm::TriggerResults>(parameters.getUntrackedParameter<edm::InputTag>("trigTag"));
63  thePhotonCollectionToken_ = consumes<reco::PhotonCollection>(parameters.getParameter<InputTag>("photonCollection"));
64  theJetCollectionToken_ = consumes<edm::View<reco::Jet> >(parameters.getParameter<InputTag>("jetCollection"));
65  theVertexCollectionToken_ = consumes<reco::VertexCollection>(parameters.getParameter<InputTag>("vertexCollection"));
66  theMinJetPt_ = parameters.getParameter<double>("minJetPt");
67  theMinPhotonEt_ = parameters.getParameter<double>("minPhotonEt");
68  theRequirePhotonFound_ = parameters.getParameter<bool>("requirePhotonFound");
69  thePlotPhotonMaxEt_ = parameters.getParameter<double>("plotPhotonMaxEt");
70  thePlotPhotonMaxEta_ = parameters.getParameter<double>("plotPhotonMaxEta");
71  thePlotJetMaxEta_ = parameters.getParameter<double>("plotJetMaxEta");
72  theBarrelRecHitTag_ = parameters.getParameter<InputTag>("barrelRecHitTag");
73  theEndcapRecHitTag_ = parameters.getParameter<InputTag>("endcapRecHitTag");
74  theBarrelRecHitToken_ = consumes<EcalRecHitCollection>(parameters.getParameter<InputTag>("barrelRecHitTag"));
75  theEndcapRecHitToken_ = consumes<EcalRecHitCollection>(parameters.getParameter<InputTag>("endcapRecHitTag"));
76 
77  // coverity says...
78  h_deltaEt_photon_jet = nullptr;
79  h_deltaPhi_jet_jet2 = nullptr;
80  h_deltaPhi_photon_jet = nullptr;
81  h_deltaPhi_photon_jet2 = nullptr;
82  h_deltaR_jet_jet2 = nullptr;
83  h_deltaR_photon_jet2 = nullptr;
84  h_jet2_eta = nullptr;
85  h_jet2_pt = nullptr;
86  h_jet2_ptOverPhotonEt = nullptr;
87  h_jet_count = nullptr;
88  h_jet_eta = nullptr;
89  h_jet_pt = nullptr;
90  h_photon_count_bar = nullptr;
91  h_photon_count_end = nullptr;
92  h_photon_et = nullptr;
93  h_photon_et_beforeCuts = nullptr;
94  h_photon_et_jetco = nullptr;
95  h_photon_et_jetcs = nullptr;
96  h_photon_et_jetfo = nullptr;
97  h_photon_et_jetfs = nullptr;
98  h_photon_eta = nullptr;
99  h_triggers_passed = nullptr;
100 }
edm::InputTag theBarrelRecHitTag_
Definition: QcdPhotonsDQM.h:61
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< edm::View< reco::Jet > > theJetCollectionToken_
Definition: QcdPhotonsDQM.h:52
MonitorElement * h_photon_et
Definition: QcdPhotonsDQM.h:69
edm::InputTag theJetCollectionLabel_
Definition: QcdPhotonsDQM.h:49
double thePlotPhotonMaxEta_
Definition: QcdPhotonsDQM.h:58
MonitorElement * h_photon_eta
Definition: QcdPhotonsDQM.h:70
MonitorElement * h_photon_count_bar
Definition: QcdPhotonsDQM.h:71
edm::InputTag theEndcapRecHitTag_
Definition: QcdPhotonsDQM.h:62
edm::EDGetTokenT< edm::TriggerResults > trigTagToken_
Definition: QcdPhotonsDQM.h:50
double thePlotJetMaxEta_
Definition: QcdPhotonsDQM.h:59
MonitorElement * h_photon_et_jetco
Definition: QcdPhotonsDQM.h:87
std::vector< std::string > thePlotTheseTriggersToo_
Definition: QcdPhotonsDQM.h:48
edm::EDGetTokenT< reco::VertexCollection > theVertexCollectionToken_
Definition: QcdPhotonsDQM.h:53
double theMinPhotonEt_
Definition: QcdPhotonsDQM.h:55
double theMinJetPt_
Definition: QcdPhotonsDQM.h:54
MonitorElement * h_deltaR_jet_jet2
Definition: QcdPhotonsDQM.h:83
MonitorElement * h_deltaR_photon_jet2
Definition: QcdPhotonsDQM.h:84
MonitorElement * h_photon_et_jetfs
Definition: QcdPhotonsDQM.h:88
MonitorElement * h_deltaEt_photon_jet
Definition: QcdPhotonsDQM.h:78
MonitorElement * h_jet_count
Definition: QcdPhotonsDQM.h:75
MonitorElement * h_jet_eta
Definition: QcdPhotonsDQM.h:74
double thePlotPhotonMaxEt_
Definition: QcdPhotonsDQM.h:57
MonitorElement * h_photon_et_jetfo
Definition: QcdPhotonsDQM.h:89
MonitorElement * h_jet2_ptOverPhotonEt
Definition: QcdPhotonsDQM.h:79
MonitorElement * h_photon_et_jetcs
Definition: QcdPhotonsDQM.h:86
MonitorElement * h_jet2_pt
Definition: QcdPhotonsDQM.h:80
MonitorElement * h_photon_et_beforeCuts
Definition: QcdPhotonsDQM.h:68
MonitorElement * h_jet_pt
Definition: QcdPhotonsDQM.h:73
edm::EDGetTokenT< EcalRecHitCollection > theBarrelRecHitToken_
Definition: QcdPhotonsDQM.h:63
bool theRequirePhotonFound_
Definition: QcdPhotonsDQM.h:56
MonitorElement * h_photon_count_end
Definition: QcdPhotonsDQM.h:72
std::string theTriggerPathToPass_
Definition: QcdPhotonsDQM.h:47
MonitorElement * h_jet2_eta
Definition: QcdPhotonsDQM.h:81
MonitorElement * h_deltaPhi_photon_jet2
Definition: QcdPhotonsDQM.h:82
MonitorElement * h_triggers_passed
Definition: QcdPhotonsDQM.h:67
edm::EDGetTokenT< reco::PhotonCollection > thePhotonCollectionToken_
Definition: QcdPhotonsDQM.h:51
MonitorElement * h_deltaPhi_jet_jet2
Definition: QcdPhotonsDQM.h:77
MonitorElement * h_deltaPhi_photon_jet
Definition: QcdPhotonsDQM.h:76
edm::EDGetTokenT< EcalRecHitCollection > theEndcapRecHitToken_
Definition: QcdPhotonsDQM.h:64
QcdPhotonsDQM::~QcdPhotonsDQM ( )
override

Destructor.

Definition at line 102 of file QcdPhotonsDQM.cc.

102 {}

Member Function Documentation

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

Get the analysis.

Definition at line 235 of file QcdPhotonsDQM.cc.

References a, funct::abs(), edm::HLTGlobalStatus::accept(), b, SiPixelRawToDigiRegional_cfi::deltaPhi, PbPb_ZMuSkimMuonDPG_cff::deltaR, reco::LeafCandidate::eta(), edm::Event::getByToken(), mps_fire::i, edm::HandleBase::isValid(), metsig::jet, jetfilter_cfi::jetCollection, LogTrace, reco::LeafCandidate::phi(), ExoticaDQM_cfi::photonCollection, edm::Handle< T >::product(), reco::LeafCandidate::pt(), edm::TriggerNames::size(), edm::TriggerNames::triggerName(), edm::Event::triggerNames(), trigNames, and spclusmultinvestigator_cfi::vertexCollection.

235  {
236  LogTrace(logTraceName) << "Analysis of event # ";
237 
239  // Did event pass HLT paths?
240  Handle<TriggerResults> HLTresults;
241  iEvent.getByToken(trigTagToken_, HLTresults);
242  if (!HLTresults.isValid()) {
243  // LogWarning("") << ">>> TRIGGER collection does not exist !!!";
244  return;
245  }
246  const edm::TriggerNames& trigNames = iEvent.triggerNames(*HLTresults);
247 
248  bool passed_HLT = false;
249 
250  // See if event passed trigger paths
251  // increment that bin in the trigger plot
252  for (unsigned int i = 0; i < thePlotTheseTriggersToo_.size(); i++) {
253  passed_HLT = false;
254  for (unsigned int ti = 0; (ti < trigNames.size()) && !passed_HLT; ++ti) {
255  size_t pos = trigNames.triggerName(ti).find(thePlotTheseTriggersToo_[i]);
256  if (pos == 0)
257  passed_HLT = HLTresults->accept(ti);
258  }
259  if (passed_HLT)
261  }
262 
263  // grab photons
265  iEvent.getByToken(thePhotonCollectionToken_, photonCollection);
266 
267  // If photon collection is empty, exit
268  if (!photonCollection.isValid())
269  return;
270 
271  // Quit if the event did not pass the HLT path we care about
272  passed_HLT = false;
273  {
274  // bool found=false;
275  for (unsigned int ti = 0; ti < trigNames.size(); ++ti) {
276  size_t pos = trigNames.triggerName(ti).find(theTriggerPathToPass_);
277  if (pos == 0) {
278  passed_HLT = HLTresults->accept(ti);
279  // found=true;
280  break;
281  }
282  }
283 
284  // Assumption: reco photons are ordered by Et
285  for (PhotonCollection::const_iterator recoPhoton = photonCollection->begin(); recoPhoton != photonCollection->end();
286  recoPhoton++) {
287  // stop looping over photons once we get to too low Et
288  if (recoPhoton->et() < theMinPhotonEt_)
289  break;
290 
291  h_photon_et_beforeCuts->Fill(recoPhoton->et());
292  break; // leading photon only
293  }
294 
295  if (!passed_HLT) {
296  return;
297  }
298  }
299 
301 
302  // std::cout << "\tpassed main trigger (" << theTriggerPathToPass_ << ")" <<
303  // std::endl;
304 
306  // Does event have valid vertex?
307  // Get the primary event vertex
308  Handle<VertexCollection> vertexHandle;
309  iEvent.getByToken(theVertexCollectionToken_, vertexHandle);
310  VertexCollection vertexCollection = *(vertexHandle.product());
311  // double vtx_ndof = -1.0;
312  // double vtx_z = 0.0;
313  // bool vtx_isFake = true;
314  // if (vertexCollection.size()>0) {
315  // vtx_ndof = vertexCollection.begin()->ndof();
316  // vtx_z = vertexCollection.begin()->z();
317  // vtx_isFake = false;
318  //}
319  // if (vtx_isFake || fabs(vtx_z)>15 || vtx_ndof<4) return;
320 
321  int nvvertex = 0;
322  for (unsigned int i = 0; i < vertexCollection.size(); ++i) {
323  if (vertexCollection[i].isValid())
324  nvvertex++;
325  }
326  if (nvvertex == 0)
327  return;
328 
330 
331  // std::cout << "\tpassed vertex selection" << std::endl;
332 
334  // Did the event pass certain L1 Technical Trigger bits?
335  // It's probably beam halo
336  // TODO: ADD code
338 
339  // For finding spikes
340  Handle<EcalRecHitCollection> EBReducedRecHits;
341  iEvent.getByToken(theBarrelRecHitToken_, EBReducedRecHits);
342  Handle<EcalRecHitCollection> EEReducedRecHits;
343  iEvent.getByToken(theEndcapRecHitToken_, EEReducedRecHits);
345 
346  // Find the highest et "decent" photon
347  float photon_et = -9.0;
348  float photon_eta = -9.0;
349  float photon_phi = -9.0;
350  bool photon_passPhotonID = false;
351  bool found_lead_pho = false;
352  int photon_count_bar = 0;
353  int photon_count_end = 0;
354  // False Assumption: reco photons are ordered by Et
355  // find the photon with highest et
356  auto pho_maxet = std::max_element(
357  photonCollection->begin(),
358  photonCollection->end(),
359  [](const PhotonCollection::value_type& a, const PhotonCollection::value_type& b) { return a.et() < b.et(); });
360  if (pho_maxet != photonCollection->end() && pho_maxet->et() >= theMinPhotonEt_) {
361  /*
362  // Ignore ECAL Spikes
363  const reco::CaloClusterPtr seed = pho_maxet->superCluster()->seed();
364  DetId id = lazyTool.getMaximum(*seed).first; // Cluster shape variables
365  // float time = -999., outOfTimeChi2 = -999., chi2 = -999.; // UNUSED
366  int flags=-1, severity = -1;
367  const EcalRecHitCollection & rechits = ( pho_maxet->isEB() ?
368  *EBReducedRecHits : *EEReducedRecHits);
369  EcalRecHitCollection::const_iterator it = rechits.find( id );
370  if( it != rechits.end() ) {
371  // time = it->time(); // UNUSED
372  // outOfTimeChi2 = it->outOfTimeChi2(); // UNUSED
373  // chi2 = it->chi2(); // UNUSED
374  flags = it->recoFlag();
375 
376  edm::ESHandle<EcalSeverityLevelAlgo> sevlv;
377  iSetup.get<EcalSeverityLevelAlgoRcd>().get(sevlv);
378  severity = sevlv->severityLevel( id, rechits);
379  }
380  bool isNotSpike = ((pho_maxet->isEB() && (severity!=3 && severity!=4 ) &&
381  (flags != 2) ) || pho_maxet->isEE());
382  if (!isNotSpike) continue; // move on to next photon
383  // END of determining ECAL Spikes
384  */
385 
386  bool pho_current_passPhotonID = false;
387  bool pho_current_isEB = pho_maxet->isEB();
388  bool pho_current_isEE = pho_maxet->isEE();
389 
390  if (pho_current_isEB && (pho_maxet->sigmaIetaIeta() < 0.01 || pho_maxet->hadronicOverEm() < 0.05)) {
391  // Photon object in barrel passes photon ID
392  pho_current_passPhotonID = true;
393  photon_count_bar++;
394  } else if (pho_current_isEE && (pho_maxet->hadronicOverEm() < 0.05)) {
395  // Photon object in endcap passes photon ID
396  pho_current_passPhotonID = true;
397  photon_count_end++;
398  }
399 
400  if (!found_lead_pho) {
401  found_lead_pho = true;
402  photon_passPhotonID = pho_current_passPhotonID;
403  photon_et = pho_maxet->et();
404  photon_eta = pho_maxet->eta();
405  photon_phi = pho_maxet->phi();
406  }
407  }
408 
409  // If user requires a photon to be found, but none is, return.
410  // theRequirePhotonFound should pretty much always be set to 'True'
411  // except when running on qcd monte carlo just to see the jets.
412  if (theRequirePhotonFound_ && (!photon_passPhotonID || photon_et < theMinPhotonEt_))
413  return;
414 
416  // Find the highest et jet
418  iEvent.getByToken(theJetCollectionToken_, jetCollection);
419  if (!jetCollection.isValid())
420  return;
421 
422  float jet_pt = -8.0;
423  float jet_eta = -8.0;
424  float jet_phi = -8.0;
425  int jet_count = 0;
426  float jet2_pt = -9.0;
427  float jet2_eta = -9.0;
428  float jet2_phi = -9.0;
429  // Assumption: jets are ordered by Et
430  for (unsigned int i_jet = 0; i_jet < jetCollection->size(); i_jet++) {
431  const Jet* jet = &jetCollection->at(i_jet);
432 
433  float jet_current_pt = jet->pt();
434 
435  // don't care about jets that overlap with the lead photon
436  if (deltaR(jet->eta(), jet->phi(), photon_eta, photon_phi) < 0.5)
437  continue;
438  // stop looping over jets once we get to too low Et
439  if (jet_current_pt < theMinJetPt_)
440  break;
441 
442  jet_count++;
443  if (jet_current_pt > jet_pt) {
444  jet2_pt = jet_pt; // 2nd highest jet get's et from current highest
445  jet2_eta = jet_eta;
446  jet2_phi = jet_phi;
447  jet_pt = jet_current_pt; // current highest jet gets et from the new highest
448  jet_eta = jet->eta();
449  jet_phi = jet->phi();
450  } else if (jet_current_pt > jet2_pt) {
451  jet2_pt = jet_current_pt;
452  jet2_eta = jet->eta();
453  jet2_phi = jet->phi();
454  }
455  }
457 
459  // Fill histograms if a jet found
460  // NOTE: if a photon was required to be found, but wasn't
461  // we wouldn't have made it to this point in the code
462  if (jet_pt > 0.0) {
463  // Photon Plots
464  h_photon_et->Fill(photon_et);
465  h_photon_eta->Fill(photon_eta);
466  h_photon_count_bar->Fill(photon_count_bar);
467  h_photon_count_end->Fill(photon_count_end);
468 
469  // Photon Et hists for different orientations to the jet
470  if (fabs(photon_eta) < 1.45 && photon_passPhotonID) { // Lead photon is in barrel
471  if (fabs(jet_eta) < 1.45) { // jet is in barrel
472  if (photon_eta * jet_eta > 0) {
473  h_photon_et_jetcs->Fill(photon_et);
474  } else {
475  h_photon_et_jetco->Fill(photon_et);
476  }
477  } else if (jet_eta > 1.55 && jet_eta < 2.5) { // jet is in endcap
478  if (photon_eta * jet_eta > 0) {
479  h_photon_et_jetfs->Fill(photon_et);
480  } else {
481  h_photon_et_jetfo->Fill(photon_et);
482  }
483  }
484  } // END of Lead Photon is in Barrel
485 
486  // Jet Plots
487  h_jet_pt->Fill(jet_pt);
488  h_jet_eta->Fill(jet_eta);
489  h_jet_count->Fill(jet_count);
490  h_deltaPhi_photon_jet->Fill(abs(deltaPhi(photon_phi, jet_phi)));
491  if (abs(deltaPhi(photon_phi, jet_phi)) > 2.8)
492  h_deltaEt_photon_jet->Fill((photon_et - jet_pt) / photon_et);
493 
494  // 2nd Highest Jet Plots
495  if (jet2_pt > 0.0) {
496  h_jet2_pt->Fill(jet2_pt);
497  h_jet2_eta->Fill(jet2_eta);
498  h_jet2_ptOverPhotonEt->Fill(jet2_pt / photon_et);
499  h_deltaPhi_photon_jet2->Fill(abs(deltaPhi(photon_phi, jet2_phi)));
500  h_deltaPhi_jet_jet2->Fill(abs(deltaPhi(jet_phi, jet2_phi)));
501  h_deltaR_jet_jet2->Fill(deltaR(jet_eta, jet_phi, jet2_eta, jet2_phi));
502  h_deltaR_photon_jet2->Fill(deltaR(photon_eta, photon_phi, jet2_eta, jet2_phi));
503  }
504  }
505  // End of Filling histograms
507 }
edm::EDGetTokenT< edm::View< reco::Jet > > theJetCollectionToken_
Definition: QcdPhotonsDQM.h:52
MonitorElement * h_photon_et
Definition: QcdPhotonsDQM.h:69
double eta() const final
momentum pseudorapidity
MonitorElement * h_photon_eta
Definition: QcdPhotonsDQM.h:70
MonitorElement * h_photon_count_bar
Definition: QcdPhotonsDQM.h:71
edm::EDGetTokenT< edm::TriggerResults > trigTagToken_
Definition: QcdPhotonsDQM.h:50
MonitorElement * h_photon_et_jetco
Definition: QcdPhotonsDQM.h:87
std::vector< std::string > thePlotTheseTriggersToo_
Definition: QcdPhotonsDQM.h:48
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
edm::EDGetTokenT< reco::VertexCollection > theVertexCollectionToken_
Definition: QcdPhotonsDQM.h:53
double theMinPhotonEt_
Definition: QcdPhotonsDQM.h:55
double theMinJetPt_
Definition: QcdPhotonsDQM.h:54
MonitorElement * h_deltaR_jet_jet2
Definition: QcdPhotonsDQM.h:83
bool accept() const
Has at least one path accepted the event?
std::string logTraceName
Definition: QcdPhotonsDQM.h:44
MonitorElement * h_deltaR_photon_jet2
Definition: QcdPhotonsDQM.h:84
MonitorElement * h_photon_et_jetfs
Definition: QcdPhotonsDQM.h:88
double pt() const final
transverse momentum
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
MonitorElement * h_deltaEt_photon_jet
Definition: QcdPhotonsDQM.h:78
Strings::size_type size() const
Definition: TriggerNames.cc:31
MonitorElement * h_jet_count
Definition: QcdPhotonsDQM.h:75
MonitorElement * h_jet_eta
Definition: QcdPhotonsDQM.h:74
void Fill(long long x)
MonitorElement * h_photon_et_jetfo
Definition: QcdPhotonsDQM.h:89
MonitorElement * h_jet2_ptOverPhotonEt
Definition: QcdPhotonsDQM.h:79
MonitorElement * h_photon_et_jetcs
Definition: QcdPhotonsDQM.h:86
Definition: Jet.py:1
MonitorElement * h_jet2_pt
Definition: QcdPhotonsDQM.h:80
MonitorElement * h_photon_et_beforeCuts
Definition: QcdPhotonsDQM.h:68
MonitorElement * h_jet_pt
Definition: QcdPhotonsDQM.h:73
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< EcalRecHitCollection > theBarrelRecHitToken_
Definition: QcdPhotonsDQM.h:63
bool isValid() const
Definition: HandleBase.h:70
#define LogTrace(id)
bool theRequirePhotonFound_
Definition: QcdPhotonsDQM.h:56
MonitorElement * h_photon_count_end
Definition: QcdPhotonsDQM.h:72
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:57
std::string theTriggerPathToPass_
Definition: QcdPhotonsDQM.h:47
MonitorElement * h_jet2_eta
Definition: QcdPhotonsDQM.h:81
MonitorElement * h_deltaPhi_photon_jet2
Definition: QcdPhotonsDQM.h:82
T const * product() const
Definition: Handle.h:69
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:22
double b
Definition: hdecay.h:118
double a
Definition: hdecay.h:119
MonitorElement * h_triggers_passed
Definition: QcdPhotonsDQM.h:67
edm::EDGetTokenT< reco::PhotonCollection > thePhotonCollectionToken_
Definition: QcdPhotonsDQM.h:51
MonitorElement * h_deltaPhi_jet_jet2
Definition: QcdPhotonsDQM.h:77
MonitorElement * h_deltaPhi_photon_jet
Definition: QcdPhotonsDQM.h:76
double phi() const final
momentum azimuthal angle
edm::EDGetTokenT< EcalRecHitCollection > theEndcapRecHitToken_
Definition: QcdPhotonsDQM.h:64
edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const override
Definition: Event.cc:265
void QcdPhotonsDQM::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  ,
edm::EventSetup const &   
)
overridevirtual

Implements DQMEDAnalyzer.

Definition at line 104 of file QcdPhotonsDQM.cc.

References dqm::dqmstoreimpl::DQMStore::IBooker::book1D(), mps_fire::i, LogTrace, hlt_dqm_clientPB-live_cfg::me, dqm::impl::MonitorElement::setBinLabel(), dqm::dqmstoreimpl::DQMStore::IBooker::setCurrentFolder(), and AlCaHLTBitMon_QueryRunRegistry::string.

104  {
105  logTraceName = "QcdPhotonAnalyzer";
106 
107  LogTrace(logTraceName) << "Parameters initialization";
108 
109  ibooker.setCurrentFolder("Physics/QcdPhotons"); // Use folder with name of PAG
110 
111  std::stringstream aStringStream;
112  std::string aString;
113  aStringStream << theMinJetPt_;
114  aString = aStringStream.str();
115 
116  // Monitor of triggers passed
117  int numOfTriggersToMonitor = thePlotTheseTriggersToo_.size();
118  h_triggers_passed = ibooker.book1D(
119  "triggers_passed", "Events passing these trigger paths", numOfTriggersToMonitor, 0, numOfTriggersToMonitor);
120  for (int i = 0; i < numOfTriggersToMonitor; i++) {
122  }
123 
124  // Keep the number of plots and number of bins to a minimum!
125  h_photon_et_beforeCuts = ibooker.book1D(
126  "photon_et_beforeCuts", "#gamma with highest E_{T};E_{T}(#gamma) (GeV)", 20, 0., thePlotPhotonMaxEt_);
127  h_photon_et =
128  ibooker.book1D("photon_et", "#gamma with highest E_{T};E_{T}(#gamma) (GeV)", 20, 0., thePlotPhotonMaxEt_);
129  h_photon_eta = ibooker.book1D(
130  "photon_eta", "#gamma with highest E_{T};#eta(#gamma)", 40, -thePlotPhotonMaxEta_, thePlotPhotonMaxEta_);
131  h_photon_count_bar = ibooker.book1D(
132  "photon_count_bar", "Number of #gamma's passing selection (Barrel);Number of #gamma's", 8, -0.5, 7.5);
133  h_photon_count_end = ibooker.book1D(
134  "photon_count_end", "Number of #gamma's passing selection (Endcap);Number of #gamma's", 8, -0.5, 7.5);
135  h_jet_pt =
136  ibooker.book1D("jet_pt",
137  "Jet with highest p_{T} (from " + theJetCollectionLabel_.label() + ");p_{T}(1^{st} jet) (GeV)",
138  20,
139  0.,
141  h_jet_eta = ibooker.book1D("jet_eta",
142  "Jet with highest p_{T} (from " + theJetCollectionLabel_.label() + ");#eta(1^{st} jet)",
143  20,
147  ibooker.book1D("deltaPhi_photon_jet",
148  "#Delta#phi between Highest E_{T} #gamma and jet;#Delta#phi(#gamma,1^{st} jet)",
149  20,
150  0,
151  3.1415926);
152  h_deltaPhi_jet_jet2 = ibooker.book1D("deltaPhi_jet_jet2",
153  "#Delta#phi between Highest E_{T} jet and 2^{nd} "
154  "jet;#Delta#phi(1^{st} jet,2^{nd} jet)",
155  20,
156  0,
157  3.1415926);
158  h_deltaEt_photon_jet = ibooker.book1D("deltaEt_photon_jet",
159  "(E_{T}(#gamma)-p_{T}(jet))/E_{T}(#gamma) when #Delta#phi(#gamma,1^{st} "
160  "jet) > 2.8;#DeltaE_{T}(#gamma,1^{st} jet)/E_{T}(#gamma)",
161  20,
162  -1.0,
163  1.0);
164  h_jet_count =
165  ibooker.book1D("jet_count",
166  "Number of " + theJetCollectionLabel_.label() + " (p_{T} > " + aString + " GeV);Number of Jets",
167  8,
168  -0.5,
169  7.5);
170  h_jet2_pt = ibooker.book1D(
171  "jet2_pt",
172  "Jet with 2^{nd} highest p_{T} (from " + theJetCollectionLabel_.label() + ");p_{T}(2^{nd} jet) (GeV)",
173  20,
174  0.,
176  h_jet2_eta =
177  ibooker.book1D("jet2_eta",
178  "Jet with 2^{nd} highest p_{T} (from " + theJetCollectionLabel_.label() + ");#eta(2^{nd} jet)",
179  20,
182  h_jet2_ptOverPhotonEt = ibooker.book1D(
183  "jet2_ptOverPhotonEt", "p_{T}(2^{nd} highest jet) / E_{T}(#gamma);p_{T}(2^{nd} Jet)/E_{T}(#gamma)", 20, 0.0, 4.0);
184  h_deltaPhi_photon_jet2 = ibooker.book1D("deltaPhi_photon_jet2",
185  "#Delta#phi between Highest E_{T} #gamma and 2^{nd} "
186  "highest jet;#Delta#phi(#gamma,2^{nd} jet)",
187  20,
188  0,
189  3.1415926);
190  h_deltaR_jet_jet2 = ibooker.book1D(
191  "deltaR_jet_jet2", "#DeltaR between Highest Jet and 2^{nd} Highest;#DeltaR(1^{st} jet,2^{nd} jet)", 30, 0, 6.0);
192  h_deltaR_photon_jet2 = ibooker.book1D("deltaR_photon_jet2",
193  "#DeltaR between Highest E_{T} #gamma and 2^{nd} "
194  "jet;#DeltaR(#gamma, 2^{nd} jet)",
195  30,
196  0,
197  6.0);
198 
199  // Photon Et for different jet configurations
200  Float_t bins_et[] = {15, 20, 30, 50, 80};
201  int num_bins_et = 4;
202  h_photon_et_jetcs = ibooker.book1D("photon_et_jetcs",
203  "#gamma with highest E_{T} (#eta(jet)<1.45, "
204  "#eta(#gamma)#eta(jet)>0);E_{T}(#gamma) (GeV)",
205  num_bins_et,
206  bins_et);
207  h_photon_et_jetco = ibooker.book1D("photon_et_jetco",
208  "#gamma with highest E_{T} (#eta(jet)<1.45, "
209  "#eta(#gamma)#eta(jet)<0);E_{T}(#gamma) (GeV)",
210  num_bins_et,
211  bins_et);
212  h_photon_et_jetfs = ibooker.book1D("photon_et_jetfs",
213  "#gamma with highest E_{T} (1.55<#eta(jet)<2.5, "
214  "#eta(#gamma)#eta(jet)>0);E_{T}(#gamma) (GeV)",
215  num_bins_et,
216  bins_et);
217  h_photon_et_jetfo = ibooker.book1D("photon_et_jetfo",
218  "#gamma with highest E_{T} (1.55<#eta(jet)<2.5, "
219  "#eta(#gamma)#eta(jet)<0);E_{T}(#gamma) (GeV)",
220  num_bins_et,
221  bins_et);
222 
223  auto setSumw2 = [](MonitorElement* me) {
224  if (me->getTH1F()->GetSumw2N() == 0) {
225  me->enableSumw2();
226  }
227  };
228 
229  setSumw2(h_photon_et_jetcs);
230  setSumw2(h_photon_et_jetco);
231  setSumw2(h_photon_et_jetfs);
232  setSumw2(h_photon_et_jetfo);
233 }
MonitorElement * h_photon_et
Definition: QcdPhotonsDQM.h:69
edm::InputTag theJetCollectionLabel_
Definition: QcdPhotonsDQM.h:49
double thePlotPhotonMaxEta_
Definition: QcdPhotonsDQM.h:58
MonitorElement * h_photon_eta
Definition: QcdPhotonsDQM.h:70
MonitorElement * h_photon_count_bar
Definition: QcdPhotonsDQM.h:71
double thePlotJetMaxEta_
Definition: QcdPhotonsDQM.h:59
MonitorElement * h_photon_et_jetco
Definition: QcdPhotonsDQM.h:87
std::vector< std::string > thePlotTheseTriggersToo_
Definition: QcdPhotonsDQM.h:48
double theMinJetPt_
Definition: QcdPhotonsDQM.h:54
MonitorElement * h_deltaR_jet_jet2
Definition: QcdPhotonsDQM.h:83
std::string logTraceName
Definition: QcdPhotonsDQM.h:44
MonitorElement * h_deltaR_photon_jet2
Definition: QcdPhotonsDQM.h:84
MonitorElement * h_photon_et_jetfs
Definition: QcdPhotonsDQM.h:88
MonitorElement * h_deltaEt_photon_jet
Definition: QcdPhotonsDQM.h:78
MonitorElement * h_jet_count
Definition: QcdPhotonsDQM.h:75
MonitorElement * h_jet_eta
Definition: QcdPhotonsDQM.h:74
double thePlotPhotonMaxEt_
Definition: QcdPhotonsDQM.h:57
MonitorElement * h_photon_et_jetfo
Definition: QcdPhotonsDQM.h:89
MonitorElement * h_jet2_ptOverPhotonEt
Definition: QcdPhotonsDQM.h:79
MonitorElement * h_photon_et_jetcs
Definition: QcdPhotonsDQM.h:86
MonitorElement * h_jet2_pt
Definition: QcdPhotonsDQM.h:80
MonitorElement * h_photon_et_beforeCuts
Definition: QcdPhotonsDQM.h:68
MonitorElement * h_jet_pt
Definition: QcdPhotonsDQM.h:73
virtual 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)
#define LogTrace(id)
MonitorElement * h_photon_count_end
Definition: QcdPhotonsDQM.h:72
MonitorElement * h_jet2_eta
Definition: QcdPhotonsDQM.h:81
MonitorElement * h_deltaPhi_photon_jet2
Definition: QcdPhotonsDQM.h:82
std::string const & label() const
Definition: InputTag.h:36
MonitorElement * h_triggers_passed
Definition: QcdPhotonsDQM.h:67
MonitorElement * h_deltaPhi_jet_jet2
Definition: QcdPhotonsDQM.h:77
MonitorElement * h_deltaPhi_photon_jet
Definition: QcdPhotonsDQM.h:76

Member Data Documentation

MonitorElement* QcdPhotonsDQM::h_deltaEt_photon_jet
private

Definition at line 78 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaPhi_jet_jet2
private

Definition at line 77 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaPhi_photon_jet
private

Definition at line 76 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaPhi_photon_jet2
private

Definition at line 82 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaR_jet_jet2
private

Definition at line 83 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaR_photon_jet2
private

Definition at line 84 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet2_eta
private

Definition at line 81 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet2_pt
private

Definition at line 80 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet2_ptOverPhotonEt
private

Definition at line 79 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet_count
private

Definition at line 75 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet_eta
private

Definition at line 74 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet_pt
private

Definition at line 73 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_count_bar
private

Definition at line 71 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_count_end
private

Definition at line 72 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et
private

Definition at line 69 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_beforeCuts
private

Definition at line 68 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetco
private

Definition at line 87 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetcs
private

Definition at line 86 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetfo
private

Definition at line 89 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetfs
private

Definition at line 88 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_eta
private

Definition at line 70 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_triggers_passed
private

Definition at line 67 of file QcdPhotonsDQM.h.

std::string QcdPhotonsDQM::logTraceName
private

Definition at line 44 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::theBarrelRecHitTag_
private

Definition at line 61 of file QcdPhotonsDQM.h.

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

Definition at line 63 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::theEndcapRecHitTag_
private

Definition at line 62 of file QcdPhotonsDQM.h.

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

Definition at line 64 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::theJetCollectionLabel_
private

Definition at line 49 of file QcdPhotonsDQM.h.

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

Definition at line 52 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::theMinJetPt_
private

Definition at line 54 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::theMinPhotonEt_
private

Definition at line 55 of file QcdPhotonsDQM.h.

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

Definition at line 51 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::thePlotJetMaxEta_
private

Definition at line 59 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::thePlotPhotonMaxEt_
private

Definition at line 57 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::thePlotPhotonMaxEta_
private

Definition at line 58 of file QcdPhotonsDQM.h.

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

Definition at line 48 of file QcdPhotonsDQM.h.

bool QcdPhotonsDQM::theRequirePhotonFound_
private

Definition at line 56 of file QcdPhotonsDQM.h.

std::string QcdPhotonsDQM::theTriggerPathToPass_
private

Definition at line 47 of file QcdPhotonsDQM.h.

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

Definition at line 53 of file QcdPhotonsDQM.h.

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

Definition at line 50 of file QcdPhotonsDQM.h.