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
 
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 60 of file QcdPhotonsDQM.cc.

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

60  {
61  // Get parameters from configuration file
62  theTriggerPathToPass_ = parameters.getParameter<string>("triggerPathToPass");
63  thePlotTheseTriggersToo_ = 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 
120 }
edm::InputTag theBarrelRecHitTag_
Definition: QcdPhotonsDQM.h:80
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< edm::View< reco::Jet > > theJetCollectionToken_
Definition: QcdPhotonsDQM.h:71
MonitorElement * h_photon_et
Definition: QcdPhotonsDQM.h:90
edm::InputTag theJetCollectionLabel_
Definition: QcdPhotonsDQM.h:68
double thePlotPhotonMaxEta_
Definition: QcdPhotonsDQM.h:77
MonitorElement * h_photon_eta
Definition: QcdPhotonsDQM.h:91
MonitorElement * h_photon_count_bar
Definition: QcdPhotonsDQM.h:92
edm::InputTag theEndcapRecHitTag_
Definition: QcdPhotonsDQM.h:81
edm::EDGetTokenT< edm::TriggerResults > trigTagToken_
Definition: QcdPhotonsDQM.h:69
double thePlotJetMaxEta_
Definition: QcdPhotonsDQM.h:78
MonitorElement * h_photon_et_jetco
MonitorElement * h_photon_et_ratio_co_fo
std::vector< std::string > thePlotTheseTriggersToo_
Definition: QcdPhotonsDQM.h:67
edm::EDGetTokenT< reco::VertexCollection > theVertexCollectionToken_
Definition: QcdPhotonsDQM.h:72
double theMinPhotonEt_
Definition: QcdPhotonsDQM.h:74
double theMinJetPt_
Definition: QcdPhotonsDQM.h:73
MonitorElement * h_deltaR_jet_jet2
MonitorElement * h_deltaR_photon_jet2
MonitorElement * h_photon_et_jetfs
bool isValidHltConfig_
Definition: QcdPhotonsDQM.h:60
MonitorElement * h_deltaEt_photon_jet
Definition: QcdPhotonsDQM.h:99
MonitorElement * h_jet_count
Definition: QcdPhotonsDQM.h:96
MonitorElement * h_photon_et_ratio_cs_fs
MonitorElement * h_photon_et_ratio_cs_fo
MonitorElement * h_jet_eta
Definition: QcdPhotonsDQM.h:95
double thePlotPhotonMaxEt_
Definition: QcdPhotonsDQM.h:76
MonitorElement * h_photon_et_jetfo
MonitorElement * h_jet2_ptOverPhotonEt
MonitorElement * h_photon_et_jetcs
MonitorElement * h_jet2_pt
MonitorElement * h_photon_et_beforeCuts
Definition: QcdPhotonsDQM.h:89
MonitorElement * h_jet_pt
Definition: QcdPhotonsDQM.h:94
MonitorElement * h_photon_et_ratio_fo_fs
DQMStore * theDbe
Definition: QcdPhotonsDQM.h:57
edm::EDGetTokenT< EcalRecHitCollection > theBarrelRecHitToken_
Definition: QcdPhotonsDQM.h:82
bool theRequirePhotonFound_
Definition: QcdPhotonsDQM.h:75
MonitorElement * h_photon_count_end
Definition: QcdPhotonsDQM.h:93
std::string theTriggerPathToPass_
Definition: QcdPhotonsDQM.h:66
MonitorElement * h_jet2_eta
MonitorElement * h_deltaPhi_photon_jet2
MonitorElement * h_photon_et_ratio_co_fs
MonitorElement * h_triggers_passed
Definition: QcdPhotonsDQM.h:88
edm::EDGetTokenT< reco::PhotonCollection > thePhotonCollectionToken_
Definition: QcdPhotonsDQM.h:70
MonitorElement * h_deltaPhi_jet_jet2
Definition: QcdPhotonsDQM.h:98
MonitorElement * h_deltaPhi_photon_jet
Definition: QcdPhotonsDQM.h:97
edm::EDGetTokenT< EcalRecHitCollection > theEndcapRecHitToken_
Definition: QcdPhotonsDQM.h:83
MonitorElement * h_photon_et_ratio_co_cs
QcdPhotonsDQM::~QcdPhotonsDQM ( )
virtual

Destructor.

Definition at line 122 of file QcdPhotonsDQM.cc.

122  {
123 }

Member Function Documentation

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

Get the analysis.

Implements edm::EDAnalyzer.

Definition at line 208 of file QcdPhotonsDQM.cc.

References funct::abs(), 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.

208  {
210 
212  //if( ! isValidHltConfig_ ) return;
213 
214  LogTrace(logTraceName)<<"Analysis of event # ";
215 
217  // Did event pass HLT paths?
218  Handle<TriggerResults> HLTresults;
219  iEvent.getByToken(trigTagToken_, HLTresults);
220  if (!HLTresults.isValid()) {
221  //LogWarning("") << ">>> TRIGGER collection does not exist !!!";
222  return;
223  }
224  const edm::TriggerNames & trigNames = iEvent.triggerNames(*HLTresults);
225 
226  bool passed_HLT=false;
227 
228  // See if event passed trigger paths
229  // increment that bin in the trigger plot
230  for (unsigned int i=0; i<thePlotTheseTriggersToo_.size(); i++) {
231  passed_HLT = false;
232  for (unsigned int ti=0; (ti<trigNames.size()) && !passed_HLT; ++ti) {
233  size_t pos = trigNames.triggerName(ti).find(thePlotTheseTriggersToo_[i]);
234  if (pos==0) passed_HLT = HLTresults->accept(ti);
235  }
236  if (passed_HLT) h_triggers_passed->Fill(i);
237  }
238 
239 
240  // grab photons
241  Handle<PhotonCollection> photonCollection;
242  iEvent.getByToken(thePhotonCollectionToken_, photonCollection);
243 
244  // If photon collection is empty, exit
245  if (!photonCollection.isValid()) return;
246 
247 
248  // Quit if the event did not pass the HLT path we care about
249  passed_HLT = false;
250  {
251  //bool found=false;
252  for (unsigned int ti=0; ti<trigNames.size(); ++ti) {
253  size_t pos = trigNames.triggerName(ti).find(theTriggerPathToPass_);
254  if (pos==0) {
255  passed_HLT = HLTresults->accept(ti);
256  //found=true;
257  break;
258  }
259  }
260 
261  // Assumption: reco photons are ordered by Et
262  for (PhotonCollection::const_iterator recoPhoton = photonCollection->begin(); recoPhoton!=photonCollection->end(); recoPhoton++){
263  // stop looping over photons once we get to too low Et
264  if ( recoPhoton->et() < theMinPhotonEt_ ) break;
265 
266  h_photon_et_beforeCuts->Fill(recoPhoton->et());
267  break; // leading photon only
268  }
269 
270  if (!passed_HLT) {
271  return;
272  }
273  }
274 
276 
277  //std::cout << "\tpassed main trigger (" << theTriggerPathToPass_ << ")" << std::endl;
278 
280  // Does event have valid vertex?
281  // Get the primary event vertex
282  Handle<VertexCollection> vertexHandle;
283  iEvent.getByToken(theVertexCollectionToken_, vertexHandle);
284  VertexCollection vertexCollection = *(vertexHandle.product());
285  //double vtx_ndof = -1.0;
286  //double vtx_z = 0.0;
287  //bool vtx_isFake = true;
288  //if (vertexCollection.size()>0) {
289  // vtx_ndof = vertexCollection.begin()->ndof();
290  // vtx_z = vertexCollection.begin()->z();
291  // vtx_isFake = false;
292  //}
293  //if (vtx_isFake || fabs(vtx_z)>15 || vtx_ndof<4) return;
294 
295  int nvvertex = 0;
296  for (unsigned int i=0; i<vertexCollection.size(); ++i) {
297  if (vertexCollection[i].isValid()) nvvertex++;
298  }
299  if (nvvertex==0) return;
300 
302 
303  //std::cout << "\tpassed vertex selection" << std::endl;
304 
305 
307  // Did the event pass certain L1 Technical Trigger bits?
308  // It's probably beam halo
309  // TODO: ADD code
311 
312 
313  // For finding spikes
314  Handle<EcalRecHitCollection> EBReducedRecHits;
315  iEvent.getByToken(theBarrelRecHitToken_, EBReducedRecHits);
316  Handle<EcalRecHitCollection> EEReducedRecHits;
317  iEvent.getByToken(theEndcapRecHitToken_, EEReducedRecHits);
319 
320 
321  // Find the highest et "decent" photon
322  float photon_et = -9.0;
323  float photon_eta = -9.0;
324  float photon_phi = -9.0;
325  bool photon_passPhotonID = false;
326  bool found_lead_pho = false;
327  int photon_count_bar = 0;
328  int photon_count_end = 0;
329  // Assumption: reco photons are ordered by Et
330  for (PhotonCollection::const_iterator recoPhoton = photonCollection->begin(); recoPhoton!=photonCollection->end(); recoPhoton++){
331 
332  // stop looping over photons once we get to too low Et
333  if ( recoPhoton->et() < theMinPhotonEt_ ) break;
334 
335  /*
336  // Ignore ECAL Spikes
337  const reco::CaloClusterPtr seed = recoPhoton->superCluster()->seed();
338  DetId id = lazyTool.getMaximum(*seed).first; // Cluster shape variables
339  // float time = -999., outOfTimeChi2 = -999., chi2 = -999.; // UNUSED
340  int flags=-1, severity = -1;
341  const EcalRecHitCollection & rechits = ( recoPhoton->isEB() ? *EBReducedRecHits : *EEReducedRecHits);
342  EcalRecHitCollection::const_iterator it = rechits.find( id );
343  if( it != rechits.end() ) {
344  // time = it->time(); // UNUSED
345  // outOfTimeChi2 = it->outOfTimeChi2(); // UNUSED
346  // chi2 = it->chi2(); // UNUSED
347  flags = it->recoFlag();
348 
349  edm::ESHandle<EcalSeverityLevelAlgo> sevlv;
350  iSetup.get<EcalSeverityLevelAlgoRcd>().get(sevlv);
351  severity = sevlv->severityLevel( id, rechits);
352  }
353  bool isNotSpike = ((recoPhoton->isEB() && (severity!=3 && severity!=4 ) && (flags != 2) ) || recoPhoton->isEE());
354  if (!isNotSpike) continue; // move on to next photon
355  // END of determining ECAL Spikes
356  */
357 
358  bool pho_current_passPhotonID = false;
359  bool pho_current_isEB = recoPhoton->isEB();
360  bool pho_current_isEE = recoPhoton->isEE();
361 
362  if ( pho_current_isEB && (recoPhoton->sigmaIetaIeta() < 0.01 || recoPhoton->hadronicOverEm() < 0.05) ) {
363  // Photon object in barrel passes photon ID
364  pho_current_passPhotonID = true;
365  photon_count_bar++;
366  } else if ( pho_current_isEE && (recoPhoton->hadronicOverEm() < 0.05) ) {
367  // Photon object in endcap passes photon ID
368  pho_current_passPhotonID = true;
369  photon_count_end++;
370  }
371 
372  if (!found_lead_pho) {
373  found_lead_pho = true;
374  photon_passPhotonID = pho_current_passPhotonID;
375  photon_et = recoPhoton->et();
376  photon_eta = recoPhoton->eta();
377  photon_phi = recoPhoton->phi();
378  }
379  }
380 
381  // If user requires a photon to be found, but none is, return.
382  // theRequirePhotonFound should pretty much always be set to 'True'
383  // except when running on qcd monte carlo just to see the jets.
384  if ( theRequirePhotonFound_ && (!photon_passPhotonID || photon_et<theMinPhotonEt_) ) return;
385 
386 
388  // Find the highest et jet
389  Handle<View<Jet> > jetCollection;
390  iEvent.getByToken (theJetCollectionToken_,jetCollection);
391  if (!jetCollection.isValid()) return;
392 
393  float jet_pt = -8.0;
394  float jet_eta = -8.0;
395  float jet_phi = -8.0;
396  int jet_count = 0;
397  float jet2_pt = -9.0;
398  float jet2_eta = -9.0;
399  float jet2_phi = -9.0;
400  // Assumption: jets are ordered by Et
401  for (unsigned int i_jet = 0; i_jet < jetCollection->size(); i_jet++) {
402  const Jet* jet = & jetCollection->at(i_jet);
403 
404  float jet_current_pt = jet->pt();
405 
406  // don't care about jets that overlap with the lead photon
407  if ( deltaR(jet->eta(), jet->phi(), photon_eta, photon_phi) < 0.5 ) continue;
408  // stop looping over jets once we get to too low Et
409  if (jet_current_pt < theMinJetPt_) break;
410 
411  jet_count++;
412  if (jet_current_pt > jet_pt) {
413  jet2_pt = jet_pt; // 2nd highest jet get's et from current highest
414  jet2_eta = jet_eta;
415  jet2_phi = jet_phi;
416  jet_pt = jet_current_pt; // current highest jet gets et from the new highest
417  jet_eta = jet->eta();
418  jet_phi = jet->phi();
419  } else if (jet_current_pt > jet2_pt) {
420  jet2_pt = jet_current_pt;
421  jet2_eta = jet->eta();
422  jet2_phi = jet->phi();
423  }
424  }
426 
427 
429  // Fill histograms if a jet found
430  // NOTE: if a photon was required to be found, but wasn't
431  // we wouldn't have made it to this point in the code
432  if ( jet_pt > 0.0 ) {
433 
434  // Photon Plots
435  h_photon_et ->Fill( photon_et );
436  h_photon_eta ->Fill( photon_eta );
437  h_photon_count_bar->Fill( photon_count_bar );
438  h_photon_count_end->Fill( photon_count_end );
439 
440  // Photon Et hists for different orientations to the jet
441  if ( fabs(photon_eta)<1.45 && photon_passPhotonID ) { // Lead photon is in barrel
442  if (fabs(jet_eta)<1.45){ // jet is in barrel
443  if (photon_eta*jet_eta>0) {
444  h_photon_et_jetcs->Fill(photon_et);
445  } else {
446  h_photon_et_jetco->Fill(photon_et);
447  }
448  } else if (jet_eta>1.55 && jet_eta<2.5) { // jet is in endcap
449  if (photon_eta*jet_eta>0) {
450  h_photon_et_jetfs->Fill(photon_et);
451  } else {
452  h_photon_et_jetfo->Fill(photon_et);
453  }
454  }
455  } // END of Lead Photon is in Barrel
456 
457  // Jet Plots
458  h_jet_pt ->Fill( jet_pt );
459  h_jet_eta ->Fill( jet_eta );
460  h_jet_count ->Fill( jet_count );
461  h_deltaPhi_photon_jet ->Fill( abs(deltaPhi(photon_phi, jet_phi)) );
462  if ( abs(deltaPhi(photon_phi,jet_phi))>2.8 ) h_deltaEt_photon_jet->Fill( (photon_et-jet_pt)/photon_et );
463 
464  // 2nd Highest Jet Plots
465  if ( jet2_pt > 0.0 ) {
466  h_jet2_pt ->Fill( jet2_pt );
467  h_jet2_eta ->Fill( jet2_eta );
468  h_jet2_ptOverPhotonEt ->Fill( jet2_pt/photon_et );
469  h_deltaPhi_photon_jet2->Fill( abs(deltaPhi(photon_phi, jet2_phi)) );
470  h_deltaPhi_jet_jet2 ->Fill( abs(deltaPhi( jet_phi, jet2_phi)) );
471  h_deltaR_jet_jet2 ->Fill( deltaR( jet_eta, jet_phi, jet2_eta, jet2_phi) );
472  h_deltaR_photon_jet2 ->Fill( deltaR(photon_eta, photon_phi, jet2_eta, jet2_phi) );
473  }
474  }
475  // End of Filling histograms
477 }
edm::InputTag theBarrelRecHitTag_
Definition: QcdPhotonsDQM.h:80
edm::EDGetTokenT< edm::View< reco::Jet > > theJetCollectionToken_
Definition: QcdPhotonsDQM.h:71
int i
Definition: DBlmapReader.cc:9
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:204
MonitorElement * h_photon_et
Definition: QcdPhotonsDQM.h:90
MonitorElement * h_photon_eta
Definition: QcdPhotonsDQM.h:91
MonitorElement * h_photon_count_bar
Definition: QcdPhotonsDQM.h:92
edm::InputTag theEndcapRecHitTag_
Definition: QcdPhotonsDQM.h:81
edm::EDGetTokenT< edm::TriggerResults > trigTagToken_
Definition: QcdPhotonsDQM.h:69
MonitorElement * h_photon_et_jetco
std::vector< std::string > thePlotTheseTriggersToo_
Definition: QcdPhotonsDQM.h:67
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
edm::EDGetTokenT< reco::VertexCollection > theVertexCollectionToken_
Definition: QcdPhotonsDQM.h:72
double theMinPhotonEt_
Definition: QcdPhotonsDQM.h:74
double theMinJetPt_
Definition: QcdPhotonsDQM.h:73
MonitorElement * h_deltaR_jet_jet2
Base class for all types of Jets.
Definition: Jet.h:20
std::string logTraceName
Definition: QcdPhotonsDQM.h:63
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:99
Strings::size_type size() const
Definition: TriggerNames.cc:39
MonitorElement * h_jet_count
Definition: QcdPhotonsDQM.h:96
MonitorElement * h_jet_eta
Definition: QcdPhotonsDQM.h:95
tuple vertexCollection
void Fill(long long x)
MonitorElement * h_photon_et_jetfo
MonitorElement * h_jet2_ptOverPhotonEt
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
MonitorElement * h_photon_et_jetcs
MonitorElement * h_jet2_pt
MonitorElement * h_photon_et_beforeCuts
Definition: QcdPhotonsDQM.h:89
MonitorElement * h_jet_pt
Definition: QcdPhotonsDQM.h:94
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< EcalRecHitCollection > theBarrelRecHitToken_
Definition: QcdPhotonsDQM.h:82
bool isValid() const
Definition: HandleBase.h:76
#define LogTrace(id)
bool theRequirePhotonFound_
Definition: QcdPhotonsDQM.h:75
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:93
static const char *const trigNames[]
Definition: EcalDumpRaw.cc:74
std::string theTriggerPathToPass_
Definition: QcdPhotonsDQM.h:66
MonitorElement * h_jet2_eta
MonitorElement * h_deltaPhi_photon_jet2
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:27
T const * product() const
Definition: Handle.h:81
MonitorElement * h_triggers_passed
Definition: QcdPhotonsDQM.h:88
edm::EDGetTokenT< reco::PhotonCollection > thePhotonCollectionToken_
Definition: QcdPhotonsDQM.h:70
MonitorElement * h_deltaPhi_jet_jet2
Definition: QcdPhotonsDQM.h:98
MonitorElement * h_deltaPhi_photon_jet
Definition: QcdPhotonsDQM.h:97
virtual float pt() const GCC11_FINAL
transverse momentum
edm::EDGetTokenT< EcalRecHitCollection > theEndcapRecHitToken_
Definition: QcdPhotonsDQM.h:83
void QcdPhotonsDQM::beginJob ( void  )
virtual

Inizialize parameters for histo binning.

Reimplemented from edm::EDAnalyzer.

Definition at line 126 of file QcdPhotonsDQM.cc.

References i, LogTrace, and AlCaHLTBitMon_QueryRunRegistry::string.

126  {
127 
128  logTraceName = "QcdPhotonAnalyzer";
129 
130  LogTrace(logTraceName)<<"Parameters initialization";
131 
132  theDbe->setCurrentFolder("Physics/QcdPhotons"); // Use folder with name of PAG
133 
134  std::stringstream aStringStream;
135  std::string aString;
136  aStringStream << theMinJetPt_;
137  aString = aStringStream.str();
138 
139  // Monitor of triggers passed
140  int numOfTriggersToMonitor = thePlotTheseTriggersToo_.size();
141  h_triggers_passed = theDbe->book1D("triggers_passed", "Events passing these trigger paths", 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("photon_et_beforeCuts", "#gamma with highest E_{T};E_{T}(#gamma) (GeV)", 20, 0., thePlotPhotonMaxEt_);
148  h_photon_et = theDbe->book1D("photon_et", "#gamma with highest E_{T};E_{T}(#gamma) (GeV)", 20, 0., thePlotPhotonMaxEt_);
149  h_photon_eta = theDbe->book1D("photon_eta", "#gamma with highest E_{T};#eta(#gamma)", 40, -thePlotPhotonMaxEta_, thePlotPhotonMaxEta_);
150  h_photon_count_bar = theDbe->book1D("photon_count_bar","Number of #gamma's passing selection (Barrel);Number of #gamma's", 8, -0.5, 7.5);
151  h_photon_count_end = theDbe->book1D("photon_count_end","Number of #gamma's passing selection (Endcap);Number of #gamma's", 8, -0.5, 7.5);
152 
153  h_jet_pt = theDbe->book1D("jet_pt", "Jet with highest p_{T} (from "+theJetCollectionLabel_.label()+");p_{T}(1^{st} jet) (GeV)", 20, 0., thePlotPhotonMaxEt_);
154  h_jet_eta = theDbe->book1D("jet_eta", "Jet with highest p_{T} (from "+theJetCollectionLabel_.label()+");#eta(1^{st} jet)", 20, -thePlotJetMaxEta_, thePlotJetMaxEta_);
155  h_deltaPhi_photon_jet = theDbe->book1D("deltaPhi_photon_jet", "#Delta#phi between Highest E_{T} #gamma and jet;#Delta#phi(#gamma,1^{st} jet)", 20, 0, 3.1415926);
156  h_deltaPhi_jet_jet2 = theDbe->book1D("deltaPhi_jet_jet2", "#Delta#phi between Highest E_{T} jet and 2^{nd} jet;#Delta#phi(1^{st} jet,2^{nd} jet)", 20, 0, 3.1415926);
157  h_deltaEt_photon_jet = theDbe->book1D("deltaEt_photon_jet", "(E_{T}(#gamma)-p_{T}(jet))/E_{T}(#gamma) when #Delta#phi(#gamma,1^{st} jet) > 2.8;#DeltaE_{T}(#gamma,1^{st} jet)/E_{T}(#gamma)", 20, -1.0, 1.0);
158  h_jet_count = theDbe->book1D("jet_count", "Number of "+theJetCollectionLabel_.label()+" (p_{T} > "+aString+" GeV);Number of Jets", 8, -0.5, 7.5);
159  h_jet2_pt = theDbe->book1D("jet2_pt", "Jet with 2^{nd} highest p_{T} (from "+theJetCollectionLabel_.label()+");p_{T}(2^{nd} jet) (GeV)", 20, 0., thePlotPhotonMaxEt_);
160  h_jet2_eta = theDbe->book1D("jet2_eta", "Jet with 2^{nd} highest p_{T} (from "+theJetCollectionLabel_.label()+");#eta(2^{nd} jet)", 20, -thePlotJetMaxEta_, thePlotJetMaxEta_);
161  h_jet2_ptOverPhotonEt = theDbe->book1D("jet2_ptOverPhotonEt", "p_{T}(2^{nd} highest jet) / E_{T}(#gamma);p_{T}(2^{nd} Jet)/E_{T}(#gamma)", 20, 0.0, 4.0);
162  h_deltaPhi_photon_jet2 = theDbe->book1D("deltaPhi_photon_jet2","#Delta#phi between Highest E_{T} #gamma and 2^{nd} highest jet;#Delta#phi(#gamma,2^{nd} jet)", 20, 0, 3.1415926);
163  h_deltaR_jet_jet2 = theDbe->book1D("deltaR_jet_jet2", "#DeltaR between Highest Jet and 2^{nd} Highest;#DeltaR(1^{st} jet,2^{nd} jet)", 30, 0, 6.0);
164  h_deltaR_photon_jet2 = theDbe->book1D("deltaR_photon_jet2", "#DeltaR between Highest E_{T} #gamma and 2^{nd} jet;#DeltaR(#gamma, 2^{nd} jet)", 30, 0, 6.0);
165 
166  // Photon Et for different jet configurations
167  Float_t bins_et[] = {15,20,30,50,80};
168  int num_bins_et = 4;
169  h_photon_et_jetcs = theDbe->book1D("photon_et_jetcs", "#gamma with highest E_{T} (#eta(jet)<1.45, #eta(#gamma)#eta(jet)>0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
170  h_photon_et_jetco = theDbe->book1D("photon_et_jetco", "#gamma with highest E_{T} (#eta(jet)<1.45, #eta(#gamma)#eta(jet)<0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
171  h_photon_et_jetfs = theDbe->book1D("photon_et_jetfs", "#gamma with highest E_{T} (1.55<#eta(jet)<2.5, #eta(#gamma)#eta(jet)>0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
172  h_photon_et_jetfo = theDbe->book1D("photon_et_jetfo", "#gamma with highest E_{T} (1.55<#eta(jet)<2.5, #eta(#gamma)#eta(jet)<0);E_{T}(#gamma) (GeV)", num_bins_et, bins_et);
173  h_photon_et_jetcs->getTH1F()->Sumw2();
174  h_photon_et_jetco->getTH1F()->Sumw2();
175  h_photon_et_jetfs->getTH1F()->Sumw2();
176  h_photon_et_jetfo->getTH1F()->Sumw2();
177  // Ratio of the above Photon Et distributions
178  h_photon_et_ratio_co_cs = theDbe->book1D("photon_et_ratio_00_co_cs", "D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)<0) / D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)>0);E_{T}(#gamma) (GeV); ratio", num_bins_et, bins_et);
179  h_photon_et_ratio_fo_fs = theDbe->book1D("photon_et_ratio_01_fo_fs", "D(1.55<|#eta(jet)|<2.6, #eta(jet)*#eta(#gamma)<0) / D(1.55<|#eta(jet)|<2.6, #eta(jet)*#eta(#gamma)>0);E_{T}(#gamma) (GeV); ratio", num_bins_et, bins_et);
180  h_photon_et_ratio_cs_fs = theDbe->book1D("photon_et_ratio_02_cs_fs", "D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)>0) / D(1.55<|#eta(jet)|<2.6, #eta(jet)*#eta(#gamma)>0);E_{T}(#gamma) (GeV); ratio", num_bins_et, bins_et);
181  h_photon_et_ratio_co_fs = theDbe->book1D("photon_et_ratio_03_co_fs", "D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)<0) / D(1.55<|#eta(jet)|<2.6, #eta(jet)*#eta(#gamma)>0);E_{T}(#gamma) (GeV); ratio", num_bins_et, bins_et);
182  h_photon_et_ratio_cs_fo = theDbe->book1D("photon_et_ratio_04_cs_fo", "D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)>0) / D(1.55<|#eta(jet)|<2.6, #eta(jet)*#eta(#gamma)<0);E_{T}(#gamma) (GeV); ratio", num_bins_et, bins_et);
183  h_photon_et_ratio_co_fo = theDbe->book1D("photon_et_ratio_05_co_fo", "D(|#eta(jet)|<1.45, #eta(jet)*#eta(#gamma)<0) / D(1.55<|#eta(jet)|<2.6, #eta(jet)*#eta(#gamma)<0);E_{T}(#gamma) (GeV); ratio", num_bins_et, bins_et);
184  h_photon_et_ratio_co_cs->getTH1F()->Sumw2();
185  h_photon_et_ratio_fo_fs->getTH1F()->Sumw2();
186  h_photon_et_ratio_cs_fs->getTH1F()->Sumw2();
187  h_photon_et_ratio_co_fs->getTH1F()->Sumw2();
188  h_photon_et_ratio_cs_fo->getTH1F()->Sumw2();
189  h_photon_et_ratio_co_fo->getTH1F()->Sumw2();
190 }
int i
Definition: DBlmapReader.cc:9
MonitorElement * h_photon_et
Definition: QcdPhotonsDQM.h:90
edm::InputTag theJetCollectionLabel_
Definition: QcdPhotonsDQM.h:68
double thePlotPhotonMaxEta_
Definition: QcdPhotonsDQM.h:77
MonitorElement * h_photon_eta
Definition: QcdPhotonsDQM.h:91
MonitorElement * h_photon_count_bar
Definition: QcdPhotonsDQM.h:92
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:872
double thePlotJetMaxEta_
Definition: QcdPhotonsDQM.h:78
MonitorElement * h_photon_et_jetco
MonitorElement * h_photon_et_ratio_co_fo
std::vector< std::string > thePlotTheseTriggersToo_
Definition: QcdPhotonsDQM.h:67
double theMinJetPt_
Definition: QcdPhotonsDQM.h:73
MonitorElement * h_deltaR_jet_jet2
std::string logTraceName
Definition: QcdPhotonsDQM.h:63
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:99
MonitorElement * h_jet_count
Definition: QcdPhotonsDQM.h:96
MonitorElement * h_photon_et_ratio_cs_fs
MonitorElement * h_photon_et_ratio_cs_fo
MonitorElement * h_jet_eta
Definition: QcdPhotonsDQM.h:95
double thePlotPhotonMaxEt_
Definition: QcdPhotonsDQM.h:76
MonitorElement * h_photon_et_jetfo
MonitorElement * h_jet2_ptOverPhotonEt
MonitorElement * h_photon_et_jetcs
MonitorElement * h_jet2_pt
MonitorElement * h_photon_et_beforeCuts
Definition: QcdPhotonsDQM.h:89
MonitorElement * h_jet_pt
Definition: QcdPhotonsDQM.h:94
MonitorElement * h_photon_et_ratio_fo_fs
DQMStore * theDbe
Definition: QcdPhotonsDQM.h:57
#define LogTrace(id)
MonitorElement * h_photon_count_end
Definition: QcdPhotonsDQM.h:93
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:88
MonitorElement * h_deltaPhi_jet_jet2
Definition: QcdPhotonsDQM.h:98
MonitorElement * h_deltaPhi_photon_jet
Definition: QcdPhotonsDQM.h:97
MonitorElement * h_photon_et_ratio_co_cs
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:584
void QcdPhotonsDQM::beginRun ( const edm::Run iRun,
const edm::EventSetup iSet 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 196 of file QcdPhotonsDQM.cc.

196  {
197 
198  // passed as parameter to HLTConfigProvider::init(), not yet used
199  bool isConfigChanged = false;
200 
201  // isValidHltConfig_ could be used to short-circuit analyze() in case of problems
202  isValidHltConfig_ = hltConfigProvider_.init( iRun, iSet, "HLT", isConfigChanged );
203 
204  num_events_in_run = 0;
205 }
HLTConfigProvider hltConfigProvider_
Definition: QcdPhotonsDQM.h:59
bool isValidHltConfig_
Definition: QcdPhotonsDQM.h:60
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 480 of file QcdPhotonsDQM.cc.

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

Reimplemented from edm::EDAnalyzer.

Definition at line 482 of file QcdPhotonsDQM.cc.

482  {
483  if (num_events_in_run>0) {
485  }
492 }
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:88
MonitorElement * h_photon_et_ratio_co_cs

Member Data Documentation

MonitorElement* QcdPhotonsDQM::h_deltaEt_photon_jet
private

Definition at line 99 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaPhi_jet_jet2
private

Definition at line 98 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaPhi_photon_jet
private

Definition at line 97 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaPhi_photon_jet2
private

Definition at line 103 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaR_jet_jet2
private

Definition at line 104 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaR_photon_jet2
private

Definition at line 105 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet2_eta
private

Definition at line 102 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet2_pt
private

Definition at line 101 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet2_ptOverPhotonEt
private

Definition at line 100 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet_count
private

Definition at line 96 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet_eta
private

Definition at line 95 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet_pt
private

Definition at line 94 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_count_bar
private

Definition at line 92 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_count_end
private

Definition at line 93 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et
private

Definition at line 90 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_beforeCuts
private

Definition at line 89 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetco
private

Definition at line 108 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetcs
private

Definition at line 107 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetfo
private

Definition at line 110 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetfs
private

Definition at line 109 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_co_cs
private

Definition at line 112 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_co_fo
private

Definition at line 117 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_co_fs
private

Definition at line 115 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_cs_fo
private

Definition at line 116 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_cs_fs
private

Definition at line 114 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_fo_fs
private

Definition at line 113 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_eta
private

Definition at line 91 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_triggers_passed
private

Definition at line 88 of file QcdPhotonsDQM.h.

HLTConfigProvider QcdPhotonsDQM::hltConfigProvider_
private

Definition at line 59 of file QcdPhotonsDQM.h.

bool QcdPhotonsDQM::isValidHltConfig_
private

Definition at line 60 of file QcdPhotonsDQM.h.

std::string QcdPhotonsDQM::logTraceName
private

Definition at line 63 of file QcdPhotonsDQM.h.

int QcdPhotonsDQM::num_events_in_run
private

Definition at line 85 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::theBarrelRecHitTag_
private

Definition at line 80 of file QcdPhotonsDQM.h.

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

Definition at line 82 of file QcdPhotonsDQM.h.

DQMStore* QcdPhotonsDQM::theDbe
private

Definition at line 57 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::theEndcapRecHitTag_
private

Definition at line 81 of file QcdPhotonsDQM.h.

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

Definition at line 83 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::theJetCollectionLabel_
private

Definition at line 68 of file QcdPhotonsDQM.h.

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

Definition at line 71 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::theMinJetPt_
private

Definition at line 73 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::theMinPhotonEt_
private

Definition at line 74 of file QcdPhotonsDQM.h.

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

Definition at line 70 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::thePlotJetMaxEta_
private

Definition at line 78 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::thePlotPhotonMaxEt_
private

Definition at line 76 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::thePlotPhotonMaxEta_
private

Definition at line 77 of file QcdPhotonsDQM.h.

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

Definition at line 67 of file QcdPhotonsDQM.h.

bool QcdPhotonsDQM::theRequirePhotonFound_
private

Definition at line 75 of file QcdPhotonsDQM.h.

std::string QcdPhotonsDQM::theTriggerPathToPass_
private

Definition at line 66 of file QcdPhotonsDQM.h.

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

Definition at line 72 of file QcdPhotonsDQM.h.

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

Definition at line 69 of file QcdPhotonsDQM.h.