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
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) 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
 
DQMStoretheDbe
 
edm::InputTag theEndcapRecHitTag
 
edm::InputTag theJetCollectionLabel_
 
double theMinJetPt_
 
double theMinPhotonEt_
 
edm::InputTag thePhotonCollectionLabel_
 
double thePlotJetMaxEta_
 
double thePlotPhotonMaxEt_
 
double thePlotPhotonMaxEta_
 
std::vector< std::string > thePlotTheseTriggersToo_
 
bool theRequirePhotonFound_
 
std::string theTriggerPathToPass_
 
edm::InputTag theVertexCollectionLabel_
 
edm::InputTag trigTag_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- 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::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
CurrentProcessingContext const * currentContext () const
 
- 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

Date:
2012/10/10 04:00:00
Revision:
1.18
Author
Michael B. Anderson, University of Wisconsin Madison

Definition at line 27 of file QcdPhotonsDQM.h.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 65 of file QcdPhotonsDQM.cc.

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

65  {
66  // Get parameters from configuration file
67  theTriggerPathToPass_ = parameters.getParameter<string>("triggerPathToPass");
68  thePlotTheseTriggersToo_ = parameters.getParameter<vector<string> >("plotTheseTriggersToo");
69  trigTag_ = parameters.getUntrackedParameter<edm::InputTag>("trigTag");
70  thePhotonCollectionLabel_ = parameters.getParameter<InputTag>("photonCollection");
71  theJetCollectionLabel_ = parameters.getParameter<InputTag>("jetCollection");
72  theVertexCollectionLabel_ = 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  // just to initialize
82  isValidHltConfig_ = false;
83 
84  // coverity says...
91  h_jet2_eta = 0;
92  h_jet2_pt = 0;
94  h_jet_count = 0;
95  h_jet_eta = 0;
96  h_jet_pt = 0;
99  h_photon_et = 0;
101  h_photon_et_jetco = 0;
102  h_photon_et_jetcs = 0;
103  h_photon_et_jetfo = 0;
104  h_photon_et_jetfs = 0;
111  h_photon_eta = 0;
112  h_triggers_passed = 0;
113 
115 
116 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * h_photon_et
Definition: QcdPhotonsDQM.h:85
edm::InputTag theJetCollectionLabel_
Definition: QcdPhotonsDQM.h:68
double thePlotPhotonMaxEta_
Definition: QcdPhotonsDQM.h:74
edm::InputTag theEndcapRecHitTag
Definition: QcdPhotonsDQM.h:78
edm::InputTag trigTag_
Definition: QcdPhotonsDQM.h:66
MonitorElement * h_photon_eta
Definition: QcdPhotonsDQM.h:86
MonitorElement * h_photon_count_bar
Definition: QcdPhotonsDQM.h:87
double thePlotJetMaxEta_
Definition: QcdPhotonsDQM.h:75
MonitorElement * h_photon_et_jetco
MonitorElement * h_photon_et_ratio_co_fo
edm::InputTag theBarrelRecHitTag
Definition: QcdPhotonsDQM.h:77
std::vector< std::string > thePlotTheseTriggersToo_
Definition: QcdPhotonsDQM.h:65
double theMinPhotonEt_
Definition: QcdPhotonsDQM.h:71
double theMinJetPt_
Definition: QcdPhotonsDQM.h:70
MonitorElement * h_deltaR_jet_jet2
Definition: QcdPhotonsDQM.h:99
MonitorElement * h_deltaR_photon_jet2
MonitorElement * h_photon_et_jetfs
bool isValidHltConfig_
Definition: QcdPhotonsDQM.h:58
MonitorElement * h_deltaEt_photon_jet
Definition: QcdPhotonsDQM.h:94
MonitorElement * h_jet_count
Definition: QcdPhotonsDQM.h:91
MonitorElement * h_photon_et_ratio_cs_fs
MonitorElement * h_photon_et_ratio_cs_fo
MonitorElement * h_jet_eta
Definition: QcdPhotonsDQM.h:90
double thePlotPhotonMaxEt_
Definition: QcdPhotonsDQM.h:73
MonitorElement * h_photon_et_jetfo
edm::InputTag thePhotonCollectionLabel_
Definition: QcdPhotonsDQM.h:67
MonitorElement * h_jet2_ptOverPhotonEt
Definition: QcdPhotonsDQM.h:95
MonitorElement * h_photon_et_jetcs
MonitorElement * h_jet2_pt
Definition: QcdPhotonsDQM.h:96
MonitorElement * h_photon_et_beforeCuts
Definition: QcdPhotonsDQM.h:84
MonitorElement * h_jet_pt
Definition: QcdPhotonsDQM.h:89
MonitorElement * h_photon_et_ratio_fo_fs
DQMStore * theDbe
Definition: QcdPhotonsDQM.h:55
bool theRequirePhotonFound_
Definition: QcdPhotonsDQM.h:72
MonitorElement * h_photon_count_end
Definition: QcdPhotonsDQM.h:88
std::string theTriggerPathToPass_
Definition: QcdPhotonsDQM.h:64
MonitorElement * h_jet2_eta
Definition: QcdPhotonsDQM.h:97
MonitorElement * h_deltaPhi_photon_jet2
Definition: QcdPhotonsDQM.h:98
MonitorElement * h_photon_et_ratio_co_fs
edm::InputTag theVertexCollectionLabel_
Definition: QcdPhotonsDQM.h:69
MonitorElement * h_triggers_passed
Definition: QcdPhotonsDQM.h:83
MonitorElement * h_deltaPhi_jet_jet2
Definition: QcdPhotonsDQM.h:93
MonitorElement * h_deltaPhi_photon_jet
Definition: QcdPhotonsDQM.h:92
MonitorElement * h_photon_et_ratio_co_cs
QcdPhotonsDQM::~QcdPhotonsDQM ( )
virtual

Destructor.

Definition at line 118 of file QcdPhotonsDQM.cc.

118  {
119 }

Member Function Documentation

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

Get the analysis.

Implements edm::EDAnalyzer.

Definition at line 204 of file QcdPhotonsDQM.cc.

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

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

Inizialize parameters for histo binning.

Reimplemented from edm::EDAnalyzer.

Definition at line 122 of file QcdPhotonsDQM.cc.

References i, LogTrace, and AlCaHLTBitMon_QueryRunRegistry::string.

122  {
123 
124  logTraceName = "QcdPhotonAnalyzer";
125 
126  LogTrace(logTraceName)<<"Parameters initialization";
127 
128  theDbe->setCurrentFolder("Physics/QcdPhotons"); // Use folder with name of PAG
129 
130  std::stringstream aStringStream;
131  std::string aString;
132  aStringStream << theMinJetPt_;
133  aString = aStringStream.str();
134 
135  // Monitor of triggers passed
136  int numOfTriggersToMonitor = thePlotTheseTriggersToo_.size();
137  h_triggers_passed = theDbe->book1D("triggers_passed", "Events passing these trigger paths", numOfTriggersToMonitor, 0, numOfTriggersToMonitor);
138  for (int i=0; i<numOfTriggersToMonitor; i++) {
140  }
141 
142  // Keep the number of plots and number of bins to a minimum!
143  h_photon_et_beforeCuts = theDbe->book1D("photon_et_beforeCuts", "#gamma with highest E_{T};E_{T}(#gamma) (GeV)", 20, 0., thePlotPhotonMaxEt_);
144  h_photon_et = theDbe->book1D("photon_et", "#gamma with highest E_{T};E_{T}(#gamma) (GeV)", 20, 0., thePlotPhotonMaxEt_);
145  h_photon_eta = theDbe->book1D("photon_eta", "#gamma with highest E_{T};#eta(#gamma)", 40, -thePlotPhotonMaxEta_, thePlotPhotonMaxEta_);
146  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);
147  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);
148 
149  h_jet_pt = theDbe->book1D("jet_pt", "Jet with highest p_{T} (from "+theJetCollectionLabel_.label()+");p_{T}(1^{st} jet) (GeV)", 20, 0., thePlotPhotonMaxEt_);
150  h_jet_eta = theDbe->book1D("jet_eta", "Jet with highest p_{T} (from "+theJetCollectionLabel_.label()+");#eta(1^{st} jet)", 20, -thePlotJetMaxEta_, thePlotJetMaxEta_);
151  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);
152  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);
153  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);
154  h_jet_count = theDbe->book1D("jet_count", "Number of "+theJetCollectionLabel_.label()+" (p_{T} > "+aString+" GeV);Number of Jets", 8, -0.5, 7.5);
155  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_);
156  h_jet2_eta = theDbe->book1D("jet2_eta", "Jet with 2^{nd} highest p_{T} (from "+theJetCollectionLabel_.label()+");#eta(2^{nd} jet)", 20, -thePlotJetMaxEta_, thePlotJetMaxEta_);
157  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);
158  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);
159  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);
160  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);
161 
162  // Photon Et for different jet configurations
163  Float_t bins_et[] = {15,20,30,50,80};
164  int num_bins_et = 4;
165  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);
166  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);
167  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);
168  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);
169  h_photon_et_jetcs->getTH1F()->Sumw2();
170  h_photon_et_jetco->getTH1F()->Sumw2();
171  h_photon_et_jetfs->getTH1F()->Sumw2();
172  h_photon_et_jetfo->getTH1F()->Sumw2();
173  // Ratio of the above Photon Et distributions
174  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);
175  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);
176  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);
177  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);
178  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);
179  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);
180  h_photon_et_ratio_co_cs->getTH1F()->Sumw2();
181  h_photon_et_ratio_fo_fs->getTH1F()->Sumw2();
182  h_photon_et_ratio_cs_fs->getTH1F()->Sumw2();
183  h_photon_et_ratio_co_fs->getTH1F()->Sumw2();
184  h_photon_et_ratio_cs_fo->getTH1F()->Sumw2();
185  h_photon_et_ratio_co_fo->getTH1F()->Sumw2();
186 }
int i
Definition: DBlmapReader.cc:9
MonitorElement * h_photon_et
Definition: QcdPhotonsDQM.h:85
edm::InputTag theJetCollectionLabel_
Definition: QcdPhotonsDQM.h:68
double thePlotPhotonMaxEta_
Definition: QcdPhotonsDQM.h:74
MonitorElement * h_photon_eta
Definition: QcdPhotonsDQM.h:86
MonitorElement * h_photon_count_bar
Definition: QcdPhotonsDQM.h:87
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:722
double thePlotJetMaxEta_
Definition: QcdPhotonsDQM.h:75
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:70
MonitorElement * h_deltaR_jet_jet2
Definition: QcdPhotonsDQM.h:99
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:94
MonitorElement * h_jet_count
Definition: QcdPhotonsDQM.h:91
MonitorElement * h_photon_et_ratio_cs_fs
MonitorElement * h_photon_et_ratio_cs_fo
MonitorElement * h_jet_eta
Definition: QcdPhotonsDQM.h:90
double thePlotPhotonMaxEt_
Definition: QcdPhotonsDQM.h:73
MonitorElement * h_photon_et_jetfo
MonitorElement * h_jet2_ptOverPhotonEt
Definition: QcdPhotonsDQM.h:95
MonitorElement * h_photon_et_jetcs
MonitorElement * h_jet2_pt
Definition: QcdPhotonsDQM.h:96
MonitorElement * h_photon_et_beforeCuts
Definition: QcdPhotonsDQM.h:84
MonitorElement * h_jet_pt
Definition: QcdPhotonsDQM.h:89
MonitorElement * h_photon_et_ratio_fo_fs
DQMStore * theDbe
Definition: QcdPhotonsDQM.h:55
#define LogTrace(id)
MonitorElement * h_photon_count_end
Definition: QcdPhotonsDQM.h:88
MonitorElement * h_jet2_eta
Definition: QcdPhotonsDQM.h:97
MonitorElement * h_deltaPhi_photon_jet2
Definition: QcdPhotonsDQM.h:98
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:83
MonitorElement * h_deltaPhi_jet_jet2
Definition: QcdPhotonsDQM.h:93
MonitorElement * h_deltaPhi_photon_jet
Definition: QcdPhotonsDQM.h:92
MonitorElement * h_photon_et_ratio_co_cs
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:434
void QcdPhotonsDQM::beginRun ( const edm::Run iRun,
const edm::EventSetup iSet 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 192 of file QcdPhotonsDQM.cc.

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

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

Reimplemented from edm::EDAnalyzer.

Definition at line 478 of file QcdPhotonsDQM.cc.

478  {
479  if (num_events_in_run>0) {
481  }
488 }
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:83
MonitorElement * h_photon_et_ratio_co_cs

Member Data Documentation

MonitorElement* QcdPhotonsDQM::h_deltaEt_photon_jet
private

Definition at line 94 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaPhi_jet_jet2
private

Definition at line 93 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaPhi_photon_jet
private

Definition at line 92 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaPhi_photon_jet2
private

Definition at line 98 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaR_jet_jet2
private

Definition at line 99 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaR_photon_jet2
private

Definition at line 100 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet2_eta
private

Definition at line 97 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet2_pt
private

Definition at line 96 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet2_ptOverPhotonEt
private

Definition at line 95 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet_count
private

Definition at line 91 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet_eta
private

Definition at line 90 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet_pt
private

Definition at line 89 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_count_bar
private

Definition at line 87 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_count_end
private

Definition at line 88 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et
private

Definition at line 85 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_beforeCuts
private

Definition at line 84 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetco
private

Definition at line 103 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetcs
private

Definition at line 102 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetfo
private

Definition at line 105 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetfs
private

Definition at line 104 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_co_cs
private

Definition at line 107 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_co_fo
private

Definition at line 112 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_co_fs
private

Definition at line 110 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_cs_fo
private

Definition at line 111 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_cs_fs
private

Definition at line 109 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_fo_fs
private

Definition at line 108 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_eta
private

Definition at line 86 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_triggers_passed
private

Definition at line 83 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 80 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::theBarrelRecHitTag
private

Definition at line 77 of file QcdPhotonsDQM.h.

DQMStore* QcdPhotonsDQM::theDbe
private

Definition at line 55 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::theEndcapRecHitTag
private

Definition at line 78 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::theJetCollectionLabel_
private

Definition at line 68 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::theMinJetPt_
private

Definition at line 70 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::theMinPhotonEt_
private

Definition at line 71 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::thePhotonCollectionLabel_
private

Definition at line 67 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::thePlotJetMaxEta_
private

Definition at line 75 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::thePlotPhotonMaxEt_
private

Definition at line 73 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::thePlotPhotonMaxEta_
private

Definition at line 74 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 72 of file QcdPhotonsDQM.h.

std::string QcdPhotonsDQM::theTriggerPathToPass_
private

Definition at line 64 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::theVertexCollectionLabel_
private

Definition at line 69 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::trigTag_
private

Definition at line 66 of file QcdPhotonsDQM.h.