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

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 ()
 

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_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 theCaloJetCollectionLabel
 
DQMStoretheDbe
 
std::string theHltMenu
 
double theMinCaloJetPt
 
double theMinPhotonEt
 
edm::InputTag thePhotonCollectionLabel
 
double thePlotJetMaxEta
 
double thePlotPhotonMaxEt
 
double thePlotPhotonMaxEta
 
std::vector< std::string > thePlotTheseTriggersToo
 
bool theRequirePhotonFound
 
std::string theTriggerPathToPass
 
std::string theTriggerResultsCollection
 
edm::InputTag theVertexCollectionLabel
 

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)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

DQM offline for QCD-Photons

Date:
2010/06/16 15:53:52
Revision:
1.15
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().

65  {
66  // Get parameters from configuration file
67  theTriggerPathToPass = parameters.getParameter<string>("triggerPathToPass");
68  thePlotTheseTriggersToo = parameters.getParameter<vector<string> >("plotTheseTriggersToo");
69  theHltMenu = parameters.getParameter<string>("hltMenu");
70  theTriggerResultsCollection = parameters.getParameter<string>("triggerResultsCollection");
71  thePhotonCollectionLabel = parameters.getParameter<InputTag>("photonCollection");
72  theCaloJetCollectionLabel = parameters.getParameter<InputTag>("caloJetCollection");
73  theVertexCollectionLabel = parameters.getParameter<InputTag>("vertexCollection");
74  theMinCaloJetPt = parameters.getParameter<double>("minCaloJetPt");
75  theMinPhotonEt = parameters.getParameter<double>("minPhotonEt");
76  theRequirePhotonFound = parameters.getParameter<bool>("requirePhotonFound");
77  thePlotPhotonMaxEt = parameters.getParameter<double>("plotPhotonMaxEt");
78  thePlotPhotonMaxEta = parameters.getParameter<double>("plotPhotonMaxEta");
79  thePlotJetMaxEta = parameters.getParameter<double>("plotJetMaxEta");
80  // just to initialize
81  isValidHltConfig_ = false;
82 
83 }
double theMinCaloJetPt
Definition: QcdPhotonsDQM.h:71
T getParameter(std::string const &) const
double thePlotPhotonMaxEt
Definition: QcdPhotonsDQM.h:74
double thePlotPhotonMaxEta
Definition: QcdPhotonsDQM.h:75
std::string theTriggerPathToPass
Definition: QcdPhotonsDQM.h:64
double theMinPhotonEt
Definition: QcdPhotonsDQM.h:72
bool isValidHltConfig_
Definition: QcdPhotonsDQM.h:58
double thePlotJetMaxEta
Definition: QcdPhotonsDQM.h:76
bool theRequirePhotonFound
Definition: QcdPhotonsDQM.h:73
edm::InputTag thePhotonCollectionLabel
Definition: QcdPhotonsDQM.h:68
std::vector< std::string > thePlotTheseTriggersToo
Definition: QcdPhotonsDQM.h:65
std::string theTriggerResultsCollection
Definition: QcdPhotonsDQM.h:67
edm::InputTag theVertexCollectionLabel
Definition: QcdPhotonsDQM.h:70
edm::InputTag theCaloJetCollectionLabel
Definition: QcdPhotonsDQM.h:69
std::string theHltMenu
Definition: QcdPhotonsDQM.h:66
QcdPhotonsDQM::~QcdPhotonsDQM ( )
virtual

Destructor.

Definition at line 85 of file QcdPhotonsDQM.cc.

85  {
86 }

Member Function Documentation

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

Get the analysis.

Implements edm::EDAnalyzer.

Definition at line 171 of file QcdPhotonsDQM.cc.

References abs, Geom::deltaPhi(), deltaR(), edm::SortedCollection< T, SORT >::find(), flags, edm::EventSetup::get(), edm::Event::getByLabel(), EcalClusterLazyTools::getMaximum(), i, edm::HandleBase::isValid(), LogTrace, edm::Handle< T >::product(), findQualityFiles::size, cond::rpcobgas::time, and GoodVertex_cfg::vertexCollection.

171  {
173 
174  // short-circuit if hlt problems
175  if( ! isValidHltConfig_ ) return;
176 
177  LogTrace(logTraceName)<<"Analysis of event # ";
178 
180  // Did event pass HLT paths?
181  Handle<TriggerResults> HLTresults;
182  iEvent.getByLabel(InputTag(theTriggerResultsCollection, "", theHltMenu), HLTresults);
183  if (!HLTresults.isValid()) return;
184 
185  unsigned int triggerIndex; // index of trigger path
186  bool passed_HLT;
187 
188  // See if event passed trigger paths
189  // increment that bin in the trigger plot
190  for (unsigned int i=0; i<thePlotTheseTriggersToo.size(); i++) {
191  passed_HLT = false;
193  if (triggerIndex < HLTresults->size()) passed_HLT = HLTresults->accept(triggerIndex);
194  if (passed_HLT) h_triggers_passed->Fill(i);
195  }
196 
197  // Quit if the event did not pass the HLT path we care about
198  passed_HLT = false;
199  triggerIndex = hltConfigProvider_.triggerIndex(theTriggerPathToPass); // index of trigger path
200  if (triggerIndex < HLTresults->size()) passed_HLT = HLTresults->accept(triggerIndex);
201  if (!passed_HLT) return;
203 
204 
206  // Does event have valid vertex?
207  // Get the primary event vertex
208  Handle<VertexCollection> vertexHandle;
209  iEvent.getByLabel(theVertexCollectionLabel, vertexHandle);
210  VertexCollection vertexCollection = *(vertexHandle.product());
211  double vtx_ndof = -1.0;
212  double vtx_z = 0.0;
213  bool vtx_isFake = true;
214  if (vertexCollection.size()>0) {
215  vtx_ndof = vertexCollection.begin()->ndof();
216  vtx_z = vertexCollection.begin()->z();
217  vtx_isFake = false;
218  }
219  if (vtx_isFake || fabs(vtx_z)>15 || vtx_ndof<4) return;
221 
222 
224  // Did the event pass certain L1 Technical Trigger bits?
225  // It's probably beam halo
226  // TODO: ADD code
228 
229  // grab photons
230  Handle<PhotonCollection> photonCollection;
231  iEvent.getByLabel(thePhotonCollectionLabel, photonCollection);
232 
233  // If photon collection is empty, exit
234  if (!photonCollection.isValid()) return;
235 
236  // For finding spikes
237  Handle<EcalRecHitCollection> EBReducedRecHits;
238  iEvent.getByLabel("reducedEcalRecHitsEB", EBReducedRecHits);
239  Handle<EcalRecHitCollection> EEReducedRecHits;
240  iEvent.getByLabel("reducedEcalRecHitsEE", EEReducedRecHits);
241  EcalClusterLazyTools lazyTool(iEvent, iSetup, InputTag("reducedEcalRecHitsEB"), InputTag("reducedEcalRecHitsEE") );
242 
243 
244  // Find the highest et "decent" photon
245  float photon_et = -9.0;
246  float photon_eta = -9.0;
247  float photon_phi = -9.0;
248  bool photon_passPhotonID = false;
249  bool found_lead_pho = false;
250  int photon_count_bar = 0;
251  int photon_count_end = 0;
252  // Assumption: reco photons are ordered by Et
253  for (PhotonCollection::const_iterator recoPhoton = photonCollection->begin(); recoPhoton!=photonCollection->end(); recoPhoton++){
254 
255  // stop looping over photons once we get to too low Et
256  if ( recoPhoton->et() < theMinPhotonEt ) break;
257 
258  // Ignore ECAL Spikes
259  const reco::CaloClusterPtr seed = recoPhoton->superCluster()->seed();
260  DetId id = lazyTool.getMaximum(*seed).first; // Cluster shape variables
261  float time = -999., outOfTimeChi2 = -999., chi2 = -999.;
262  int flags=-1, severity = -1;
263  const EcalRecHitCollection & rechits = ( recoPhoton->isEB() ? *EBReducedRecHits : *EEReducedRecHits);
264  EcalRecHitCollection::const_iterator it = rechits.find( id );
265  if( it != rechits.end() ) {
266  time = it->time();
267  outOfTimeChi2 = it->outOfTimeChi2();
268  chi2 = it->chi2();
269  flags = it->recoFlag();
270 
272  iSetup.get<EcalSeverityLevelAlgoRcd>().get(sevlv);
273  severity = sevlv->severityLevel( id, rechits);
274  }
275  bool isNotSpike = ((recoPhoton->isEB() && (severity!=3 && severity!=4 ) && (flags != 2) ) || recoPhoton->isEE());
276  if (!isNotSpike) continue; // move on to next photon
277  // END of determining ECAL Spikes
278 
279  bool pho_current_passPhotonID = false;
280  bool pho_current_isEB = recoPhoton->isEB();
281  bool pho_current_isEE = recoPhoton->isEE();
282 
283  if ( pho_current_isEB && (recoPhoton->sigmaIetaIeta() < 0.01 || recoPhoton->hadronicOverEm() < 0.05) ) {
284  // Photon object in barrel passes photon ID
285  pho_current_passPhotonID = true;
286  photon_count_bar++;
287  } else if ( pho_current_isEE && (recoPhoton->hadronicOverEm() < 0.05) ) {
288  // Photon object in endcap passes photon ID
289  pho_current_passPhotonID = true;
290  photon_count_end++;
291  }
292 
293  if (!found_lead_pho) {
294  found_lead_pho = true;
295  photon_passPhotonID = pho_current_passPhotonID;
296  photon_et = recoPhoton->et();
297  photon_eta = recoPhoton->eta();
298  photon_phi = recoPhoton->phi();
299  }
300  }
301 
302  // If user requires a photon to be found, but none is, return.
303  // theRequirePhotonFound should pretty much always be set to 'True'
304  // except when running on qcd monte carlo just to see the jets.
305  if ( theRequirePhotonFound && (!photon_passPhotonID || photon_et<theMinPhotonEt) ) return;
306 
307 
309  // Find the highest et jet
310  Handle<CaloJetCollection> caloJetCollection;
311  iEvent.getByLabel (theCaloJetCollectionLabel,caloJetCollection);
312  if (!caloJetCollection.isValid()) return;
313 
314  float jet_pt = -8.0;
315  float jet_eta = -8.0;
316  float jet_phi = -8.0;
317  int jet_count = 0;
318  float jet2_pt = -9.0;
319  float jet2_eta = -9.0;
320  float jet2_phi = -9.0;
321  // Assumption: jets are ordered by Et
322  for (CaloJetCollection::const_iterator i_calojet = caloJetCollection->begin(); i_calojet != caloJetCollection->end(); i_calojet++) {
323 
324  float jet_current_pt = i_calojet->pt();
325 
326  // don't care about jets that overlap with the lead photon
327  if ( deltaR(i_calojet->eta(), i_calojet->phi(), photon_eta, photon_phi) < 0.5 ) continue;
328  // stop looping over jets once we get to too low Et
329  if (jet_current_pt < theMinCaloJetPt) break;
330 
331  jet_count++;
332  if (jet_current_pt > jet_pt) {
333  jet2_pt = jet_pt; // 2nd highest jet get's et from current highest
334  jet2_eta = jet_eta;
335  jet2_phi = jet_phi;
336  jet_pt = jet_current_pt; // current highest jet gets et from the new highest
337  jet_eta = i_calojet->eta();
338  jet_phi = i_calojet->phi();
339  } else if (jet_current_pt > jet2_pt) {
340  jet2_pt = jet_current_pt;
341  jet2_eta = i_calojet->eta();
342  jet2_phi = i_calojet->phi();
343  }
344  }
346 
347 
349  // Fill histograms if a jet found
350  // NOTE: if a photon was required to be found, but wasn't
351  // we wouldn't have made it to this point in the code
352  if ( jet_pt > 0.0 ) {
353 
354  // Photon Plots
355  h_photon_et ->Fill( photon_et );
356  h_photon_eta ->Fill( photon_eta );
357  h_photon_count_bar->Fill( photon_count_bar );
358  h_photon_count_end->Fill( photon_count_end );
359 
360  // Photon Et hists for different orientations to the jet
361  if ( fabs(photon_eta)<1.45 && photon_passPhotonID ) { // Lead photon is in barrel
362  if (fabs(jet_eta)<1.45){ // jet is in barrel
363  if (photon_eta*jet_eta>0) {
364  h_photon_et_jetcs->Fill(photon_et);
365  } else {
366  h_photon_et_jetco->Fill(photon_et);
367  }
368  } else if (jet_eta>1.55 && jet_eta<2.5) { // jet is in endcap
369  if (photon_eta*jet_eta>0) {
370  h_photon_et_jetfs->Fill(photon_et);
371  } else {
372  h_photon_et_jetfo->Fill(photon_et);
373  }
374  }
375  } // END of Lead Photon is in Barrel
376 
377  // Jet Plots
378  h_jet_pt ->Fill( jet_pt );
379  h_jet_eta ->Fill( jet_eta );
380  h_jet_count ->Fill( jet_count );
381  h_deltaPhi_photon_jet ->Fill( abs(deltaPhi(photon_phi, jet_phi)) );
382  if ( abs(deltaPhi(photon_phi,jet_phi))>2.8 ) h_deltaEt_photon_jet->Fill( (photon_et-jet_pt)/photon_et );
383 
384  // 2nd Highest Jet Plots
385  if ( jet2_pt > 0.0 ) {
386  h_jet2_pt ->Fill( jet2_pt );
387  h_jet2_eta ->Fill( jet2_eta );
388  h_jet2_ptOverPhotonEt ->Fill( jet2_pt/photon_et );
389  h_deltaPhi_photon_jet2->Fill( abs(deltaPhi(photon_phi, jet2_phi)) );
390  h_deltaPhi_jet_jet2 ->Fill( abs(deltaPhi( jet_phi, jet2_phi)) );
391  h_deltaR_jet_jet2 ->Fill( deltaR( jet_eta, jet_phi, jet2_eta, jet2_phi) );
392  h_deltaR_photon_jet2 ->Fill( deltaR(photon_eta, photon_phi, jet2_eta, jet2_phi) );
393  }
394  }
395  // End of Filling histograms
397 }
double theMinCaloJetPt
Definition: QcdPhotonsDQM.h:71
HLTConfigProvider hltConfigProvider_
Definition: QcdPhotonsDQM.h:57
int i
Definition: DBlmapReader.cc:9
MonitorElement * h_photon_et
Definition: QcdPhotonsDQM.h:82
MonitorElement * h_photon_eta
Definition: QcdPhotonsDQM.h:83
MonitorElement * h_photon_count_bar
Definition: QcdPhotonsDQM.h:84
MonitorElement * h_photon_et_jetco
double deltaPhi(float phi1, float phi2)
Definition: VectorUtil.h:30
MonitorElement * h_deltaR_jet_jet2
Definition: QcdPhotonsDQM.h:96
std::string theTriggerPathToPass
Definition: QcdPhotonsDQM.h:64
std::string logTraceName
Definition: QcdPhotonsDQM.h:61
double theMinPhotonEt
Definition: QcdPhotonsDQM.h:72
std::vector< T >::const_iterator const_iterator
MonitorElement * h_deltaR_photon_jet2
Definition: QcdPhotonsDQM.h:97
MonitorElement * h_photon_et_jetfs
bool isValidHltConfig_
Definition: QcdPhotonsDQM.h:58
#define abs(x)
Definition: mlp_lapack.h:159
std::vector< Variable::Flags > flags
Definition: MVATrainer.cc:135
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
MonitorElement * h_deltaEt_photon_jet
Definition: QcdPhotonsDQM.h:91
MonitorElement * h_jet_count
Definition: QcdPhotonsDQM.h:88
bool theRequirePhotonFound
Definition: QcdPhotonsDQM.h:73
MonitorElement * h_jet_eta
Definition: QcdPhotonsDQM.h:87
tuple vertexCollection
void Fill(long long x)
unsigned int triggerIndex(const std::string &triggerName) const
slot position of trigger path in trigger table (0 to size-1)
edm::InputTag thePhotonCollectionLabel
Definition: QcdPhotonsDQM.h:68
MonitorElement * h_photon_et_jetfo
MonitorElement * h_jet2_ptOverPhotonEt
Definition: QcdPhotonsDQM.h:92
MonitorElement * h_photon_et_jetcs
Definition: QcdPhotonsDQM.h:99
MonitorElement * h_jet2_pt
Definition: QcdPhotonsDQM.h:93
MonitorElement * h_jet_pt
Definition: QcdPhotonsDQM.h:86
std::vector< std::string > thePlotTheseTriggersToo
Definition: QcdPhotonsDQM.h:65
std::string theTriggerResultsCollection
Definition: QcdPhotonsDQM.h:67
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
#define LogTrace(id)
edm::InputTag theVertexCollectionLabel
Definition: QcdPhotonsDQM.h:70
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
MonitorElement * h_photon_count_end
Definition: QcdPhotonsDQM.h:85
MonitorElement * h_jet2_eta
Definition: QcdPhotonsDQM.h:94
Definition: DetId.h:20
MonitorElement * h_deltaPhi_photon_jet2
Definition: QcdPhotonsDQM.h:95
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: Handle.h:74
iterator find(key_type k)
edm::InputTag theCaloJetCollectionLabel
Definition: QcdPhotonsDQM.h:69
MonitorElement * h_triggers_passed
Definition: QcdPhotonsDQM.h:81
MonitorElement * h_deltaPhi_jet_jet2
Definition: QcdPhotonsDQM.h:90
MonitorElement * h_deltaPhi_photon_jet
Definition: QcdPhotonsDQM.h:89
std::string theHltMenu
Definition: QcdPhotonsDQM.h:66
tuple size
Write out results.
void QcdPhotonsDQM::beginJob ( void  )
virtual

Inizialize parameters for histo binning.

Reimplemented from edm::EDAnalyzer.

Definition at line 89 of file QcdPhotonsDQM.cc.

References i, LogTrace, and cmsCodeRules.cppFunctionSkipper::operator.

89  {
90 
91  logTraceName = "QcdPhotonAnalyzer";
92 
93  LogTrace(logTraceName)<<"Parameters initialization";
95 
96  theDbe->setCurrentFolder("Physics/QcdPhotons"); // Use folder with name of PAG
97 
98  std::stringstream aStringStream;
99  std::string aString;
100  aStringStream << theMinCaloJetPt;
101  aString = aStringStream.str();
102 
103  // Monitor of triggers passed
104  int numOfTriggersToMonitor = thePlotTheseTriggersToo.size();
105  h_triggers_passed = theDbe->book1D("h_triggers_passed", "Events passing these trigger paths", numOfTriggersToMonitor, 0, numOfTriggersToMonitor);
106  for (int i=0; i<numOfTriggersToMonitor; i++) {
108  }
109 
110  // Keep the number of plots and number of bins to a minimum!
111  h_photon_et = theDbe->book1D("h_photon_et", "#gamma with highest E_{T};E_{T}(#gamma) (GeV)", 20, 0., thePlotPhotonMaxEt);
112  h_photon_eta = theDbe->book1D("h_photon_eta", "#gamma with highest E_{T};#eta(#gamma)", 40, -thePlotPhotonMaxEta, thePlotPhotonMaxEta);
113  h_photon_count_bar = theDbe->book1D("h_photon_count_bar","Number of #gamma's passing selection (Barrel);Number of #gamma's", 8, -0.5, 7.5);
114  h_photon_count_end = theDbe->book1D("h_photon_count_end","Number of #gamma's passing selection (Endcap);Number of #gamma's", 8, -0.5, 7.5);
115 
116  h_jet_pt = theDbe->book1D("h_jet_pt", "Jet with highest p_{T} (from "+theCaloJetCollectionLabel.label()+");p_{T}(1^{st} jet) (GeV)", 20, 0., thePlotPhotonMaxEt);
117  h_jet_eta = theDbe->book1D("h_jet_eta", "Jet with highest p_{T} (from "+theCaloJetCollectionLabel.label()+");#eta(1^{st} jet)", 20, -thePlotJetMaxEta, thePlotJetMaxEta);
118  h_deltaPhi_photon_jet = theDbe->book1D("h_deltaPhi_photon_jet", "#Delta#phi between Highest E_{T} #gamma and jet;#Delta#phi(#gamma,1^{st} jet)", 20, 0, 3.1415926);
119  h_deltaPhi_jet_jet2 = theDbe->book1D("h_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);
120  h_deltaEt_photon_jet = theDbe->book1D("h_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);
121  h_jet_count = theDbe->book1D("h_jet_count", "Number of "+theCaloJetCollectionLabel.label()+" (p_{T} > "+aString+" GeV);Number of Jets", 8, -0.5, 7.5);
122  h_jet2_pt = theDbe->book1D("h_jet2_pt", "Jet with 2^{nd} highest p_{T} (from "+theCaloJetCollectionLabel.label()+");p_{T}(2^{nd} jet) (GeV)", 20, 0., thePlotPhotonMaxEt);
123  h_jet2_eta = theDbe->book1D("h_jet2_eta", "Jet with 2^{nd} highest p_{T} (from "+theCaloJetCollectionLabel.label()+");#eta(2^{nd} jet)", 20, -thePlotJetMaxEta, thePlotJetMaxEta);
124  h_jet2_ptOverPhotonEt = theDbe->book1D("h_jet2_ptOverPhotonEt", "p_{T}(2^{nd} highest jet) / E_{T}(#gamma);p_{T}(2^{nd} Jet)/E_{T}(#gamma)", 20, 0.0, 4.0);
125  h_deltaPhi_photon_jet2 = theDbe->book1D("h_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);
126  h_deltaR_jet_jet2 = theDbe->book1D("h_deltaR_jet_jet2", "#DeltaR between Highest Jet and 2^{nd} Highest;#DeltaR(1^{st} jet,2^{nd} jet)", 30, 0, 6.0);
127  h_deltaR_photon_jet2 = theDbe->book1D("h_deltaR_photon_jet2", "#DeltaR between Highest E_{T} #gamma and 2^{nd} jet;#DeltaR(#gamma, 2^{nd} jet)", 30, 0, 6.0);
128 
129  // Photon Et for different jet configurations
130  Float_t bins_et[] = {15,20,30,50,80};
131  int num_bins_et = 4;
132  h_photon_et_jetcs = theDbe->book1D("h_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);
133  h_photon_et_jetco = theDbe->book1D("h_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);
134  h_photon_et_jetfs = theDbe->book1D("h_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);
135  h_photon_et_jetfo = theDbe->book1D("h_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);
136  h_photon_et_jetcs->getTH1F()->Sumw2();
137  h_photon_et_jetco->getTH1F()->Sumw2();
138  h_photon_et_jetfs->getTH1F()->Sumw2();
139  h_photon_et_jetfo->getTH1F()->Sumw2();
140  // Ratio of the above Photon Et distributions
141  h_photon_et_ratio_co_cs = theDbe->book1D("h_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);
142  h_photon_et_ratio_fo_fs = theDbe->book1D("h_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);
143  h_photon_et_ratio_cs_fs = theDbe->book1D("h_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);
144  h_photon_et_ratio_co_fs = theDbe->book1D("h_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);
145  h_photon_et_ratio_cs_fo = theDbe->book1D("h_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);
146  h_photon_et_ratio_co_fo = theDbe->book1D("h_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);
147  h_photon_et_ratio_co_cs->getTH1F()->Sumw2();
148  h_photon_et_ratio_fo_fs->getTH1F()->Sumw2();
149  h_photon_et_ratio_cs_fs->getTH1F()->Sumw2();
150  h_photon_et_ratio_co_fs->getTH1F()->Sumw2();
151  h_photon_et_ratio_cs_fo->getTH1F()->Sumw2();
152  h_photon_et_ratio_co_fo->getTH1F()->Sumw2();
153 }
double theMinCaloJetPt
Definition: QcdPhotonsDQM.h:71
int i
Definition: DBlmapReader.cc:9
double thePlotPhotonMaxEt
Definition: QcdPhotonsDQM.h:74
MonitorElement * h_photon_et
Definition: QcdPhotonsDQM.h:82
MonitorElement * h_photon_eta
Definition: QcdPhotonsDQM.h:83
MonitorElement * h_photon_count_bar
Definition: QcdPhotonsDQM.h:84
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:514
MonitorElement * h_photon_et_jetco
MonitorElement * h_photon_et_ratio_co_fo
MonitorElement * h_deltaR_jet_jet2
Definition: QcdPhotonsDQM.h:96
double thePlotPhotonMaxEta
Definition: QcdPhotonsDQM.h:75
std::string logTraceName
Definition: QcdPhotonsDQM.h:61
MonitorElement * h_deltaR_photon_jet2
Definition: QcdPhotonsDQM.h:97
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:91
MonitorElement * h_jet_count
Definition: QcdPhotonsDQM.h:88
double thePlotJetMaxEta
Definition: QcdPhotonsDQM.h:76
MonitorElement * h_photon_et_ratio_cs_fs
MonitorElement * h_photon_et_ratio_cs_fo
MonitorElement * h_jet_eta
Definition: QcdPhotonsDQM.h:87
MonitorElement * h_photon_et_jetfo
MonitorElement * h_jet2_ptOverPhotonEt
Definition: QcdPhotonsDQM.h:92
MonitorElement * h_photon_et_jetcs
Definition: QcdPhotonsDQM.h:99
MonitorElement * h_jet2_pt
Definition: QcdPhotonsDQM.h:93
MonitorElement * h_jet_pt
Definition: QcdPhotonsDQM.h:86
MonitorElement * h_photon_et_ratio_fo_fs
DQMStore * theDbe
Definition: QcdPhotonsDQM.h:55
std::vector< std::string > thePlotTheseTriggersToo
Definition: QcdPhotonsDQM.h:65
#define LogTrace(id)
MonitorElement * h_photon_count_end
Definition: QcdPhotonsDQM.h:85
MonitorElement * h_jet2_eta
Definition: QcdPhotonsDQM.h:94
MonitorElement * h_deltaPhi_photon_jet2
Definition: QcdPhotonsDQM.h:95
TH1F * getTH1F(void) const
std::string const & label() const
Definition: InputTag.h:25
edm::InputTag theCaloJetCollectionLabel
Definition: QcdPhotonsDQM.h:69
MonitorElement * h_photon_et_ratio_co_fs
MonitorElement * h_triggers_passed
Definition: QcdPhotonsDQM.h:81
MonitorElement * h_deltaPhi_jet_jet2
Definition: QcdPhotonsDQM.h:90
MonitorElement * h_deltaPhi_photon_jet
Definition: QcdPhotonsDQM.h:89
MonitorElement * h_photon_et_ratio_co_cs
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:232
void QcdPhotonsDQM::beginRun ( const edm::Run r,
const edm::EventSetup iSetup 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 159 of file QcdPhotonsDQM.cc.

159  {
160 
161  // passed as parameter to HLTConfigProvider::init(), not yet used
162  bool isConfigChanged = false;
163 
164  // isValidHltConfig_ used to short-circuit analyze() in case of problems
165  isValidHltConfig_ = hltConfigProvider_.init( r, iSetup, theHltMenu, isConfigChanged );
166 
167  num_events_in_run = 0;
168 }
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)
std::string theHltMenu
Definition: QcdPhotonsDQM.h:66
void QcdPhotonsDQM::endJob ( void  )
virtual

Save the histos.

Reimplemented from edm::EDAnalyzer.

Definition at line 400 of file QcdPhotonsDQM.cc.

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

Reimplemented from edm::EDAnalyzer.

Definition at line 402 of file QcdPhotonsDQM.cc.

402  {
403  if (num_events_in_run>0) {
405  }
412 }
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
Definition: QcdPhotonsDQM.h:99
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:81
MonitorElement * h_photon_et_ratio_co_cs

Member Data Documentation

MonitorElement* QcdPhotonsDQM::h_deltaEt_photon_jet
private

Definition at line 91 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaPhi_jet_jet2
private

Definition at line 90 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaPhi_photon_jet
private

Definition at line 89 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaPhi_photon_jet2
private

Definition at line 95 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaR_jet_jet2
private

Definition at line 96 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_deltaR_photon_jet2
private

Definition at line 97 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet2_eta
private

Definition at line 94 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet2_pt
private

Definition at line 93 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet2_ptOverPhotonEt
private

Definition at line 92 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet_count
private

Definition at line 88 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet_eta
private

Definition at line 87 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_jet_pt
private

Definition at line 86 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_count_bar
private

Definition at line 84 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_count_end
private

Definition at line 85 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et
private

Definition at line 82 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetco
private

Definition at line 100 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetcs
private

Definition at line 99 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetfo
private

Definition at line 102 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_jetfs
private

Definition at line 101 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_co_cs
private

Definition at line 104 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_co_fo
private

Definition at line 109 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_co_fs
private

Definition at line 107 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_cs_fo
private

Definition at line 108 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_cs_fs
private

Definition at line 106 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_et_ratio_fo_fs
private

Definition at line 105 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_photon_eta
private

Definition at line 83 of file QcdPhotonsDQM.h.

MonitorElement* QcdPhotonsDQM::h_triggers_passed
private

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

edm::InputTag QcdPhotonsDQM::theCaloJetCollectionLabel
private

Definition at line 69 of file QcdPhotonsDQM.h.

DQMStore* QcdPhotonsDQM::theDbe
private

Definition at line 55 of file QcdPhotonsDQM.h.

std::string QcdPhotonsDQM::theHltMenu
private

Definition at line 66 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::theMinCaloJetPt
private

Definition at line 71 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::theMinPhotonEt
private

Definition at line 72 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::thePhotonCollectionLabel
private

Definition at line 68 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::thePlotJetMaxEta
private

Definition at line 76 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::thePlotPhotonMaxEt
private

Definition at line 74 of file QcdPhotonsDQM.h.

double QcdPhotonsDQM::thePlotPhotonMaxEta
private

Definition at line 75 of file QcdPhotonsDQM.h.

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

Definition at line 65 of file QcdPhotonsDQM.h.

bool QcdPhotonsDQM::theRequirePhotonFound
private

Definition at line 73 of file QcdPhotonsDQM.h.

std::string QcdPhotonsDQM::theTriggerPathToPass
private

Definition at line 64 of file QcdPhotonsDQM.h.

std::string QcdPhotonsDQM::theTriggerResultsCollection
private

Definition at line 67 of file QcdPhotonsDQM.h.

edm::InputTag QcdPhotonsDQM::theVertexCollectionLabel
private

Definition at line 70 of file QcdPhotonsDQM.h.