CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DiJetVarAnalyzer.cc
Go to the documentation of this file.
2 
5 
8 
11 
14 
16 
20 
21 #include <cmath>
22 
23 //------------------------------------------------------------------------------
24 // A simple constructor which takes as inoput only the name of the PF jet collection
27  jetCollectionTag_ (conf.getUntrackedParameter<edm::InputTag>("jetCollectionTag")),
28  //dijetVarCollectionTag_ (conf.getUntrackedParameter<edm::InputTag>("dijetVarCollectionTag")),
29  widejetsCollectionTag_ (conf.getUntrackedParameter<edm::InputTag>("widejetsCollectionTag")),
30  metCollectionTag_ (conf.getUntrackedParameter<edm::InputTag>("metCollectionTag")),
31  metCleanCollectionTag_ (conf.getUntrackedParameter<edm::InputTag>("metCleanCollectionTag")),
32  numwidejets_ (conf.getParameter<unsigned int>("numwidejets")),
33  etawidejets_ (conf.getParameter<double>("etawidejets")),
34  ptwidejets_ (conf.getParameter<double>("ptwidejets")),
35  detawidejets_ (conf.getParameter<double>("detawidejets")),
36  dphiwidejets_ (conf.getParameter<double>("dphiwidejets")),
37  maxEMfraction_ (conf.getParameter<double>("maxEMfraction")),
38  maxHADfraction_ (conf.getParameter<double>("maxHADfraction")),
39  HLTpathMain_ (triggerExpression::parse( conf.getParameter<std::string>("HLTpathMain") )),
40  HLTpathMonitor_ (triggerExpression::parse( conf.getParameter<std::string>("HLTpathMonitor") )),
41  triggerConfiguration_ (conf.getParameterSet("triggerConfiguration"))
42 {
43 }
44 
45 //------------------------------------------------------------------------------
46 // Nothing to destroy: the DQM service thinks about everything
48  delete HLTpathMain_;
49  delete HLTpathMonitor_;
50 }
51 
52 //------------------------------------------------------------------------------
53 void DiJetVarAnalyzer::beginRun( const edm::Run & iRun, const edm::EventSetup & iSetup){
54 }
55 
56 //------------------------------------------------------------------------------
57 // Usual analyze method
59 
60  using namespace std;
61  using namespace edm;
62  using namespace reco;
63 
64  // ## Get jet collection
66  iEvent.getByLabel(jetCollectionTag_,calojets_handle);
67 
68  // Loop over all the selected jets ( defined at DQM/DataScouting/python/dijetScouting_cff.py )
69  double thisHT = 0.;
70  for(reco::CaloJetCollection::const_iterator it = calojets_handle->begin(); it != calojets_handle->end(); ++it)
71  {
72  //cout << "== jet: " << it->pt() << " " << it->eta() << " " << it->phi() << endl;
73  m_selJets_pt->Fill( it->pt() );
74  m_selJets_eta->Fill( it->eta() );
75  m_selJets_phi->Fill( it->phi() );
76  m_selJets_hadEnergyFraction->Fill( it->energyFractionHadronic() );
77  m_selJets_emEnergyFraction->Fill( it->emEnergyFraction() );
78  m_selJets_towersArea->Fill( it->towersArea() );
79 
80  thisHT += it->pt();
81  }
82 
83  // HT
84  m_HT_inclusive->Fill(thisHT);
85 
86  // ## Get widejets
88  iEvent.getByLabel (widejetsCollectionTag_,widejets_handle);
89 
90  TLorentzVector wj1;
91  TLorentzVector wj2;
92  TLorentzVector wdijet;
93 
94  double MJJWide = -1;
95  double DeltaEtaJJWide = -1;
96  double DeltaPhiJJWide = -1;
97 
98  if( widejets_handle->size() >= 2 )
99  {
100  wj1.SetPtEtaPhiM(widejets_handle->at(0).pt(),
101  widejets_handle->at(0).eta(),
102  widejets_handle->at(0).phi(),
103  widejets_handle->at(0).mass()
104  );
105  wj2.SetPtEtaPhiM(widejets_handle->at(1).pt(),
106  widejets_handle->at(1).eta(),
107  widejets_handle->at(1).phi(),
108  widejets_handle->at(1).mass()
109  );
110 
111  wdijet = wj1 + wj2;
112 
113  MJJWide = wdijet.M();
114  DeltaEtaJJWide = fabs(wj1.Eta()-wj2.Eta());
115  DeltaPhiJJWide = fabs(wj1.DeltaPhi(wj2));
116 
117  // cout << "== j1 wide: " << wj1.Pt() << " " << wj1.Eta() << " " << wj1.Phi() << " " << wj1.M() << endl;
118  // cout << "== j2 wide: " << wj2.Pt() << " " << wj2.Eta() << " " << wj2.Phi() << " " << wj2.M() << endl;
119  // cout << "== MJJWide: " << MJJWide << endl;
120  // cout << "== DeltaEtaJJWide: " << DeltaEtaJJWide << endl;
121  // cout << "== DeltaPhiJJWide: " << DeltaPhiJJWide << endl;
122  }
123 
124  // ## Get met collection
125 
126  // met
128  iEvent.getByLabel(metCollectionTag_,calomet_handle);
129 
130  // met cleaned
131  edm::Handle<reco::CaloMETCollection> calometClean_handle;
132  iEvent.getByLabel(metCleanCollectionTag_,calometClean_handle);
133 
134  if( calomet_handle.isValid() && calometClean_handle.isValid() )
135  {
136  // std::cout << "---" << std::endl;
137  // std::cout << "== calomet: " << (calomet_handle->front()).pt() << " " << (calomet_handle->front()).phi() << std::endl;
138  // std::cout << "== calometClean: " << (calometClean_handle->front()).pt() << " " << (calometClean_handle->front()).phi() << std::endl;
139  // std::cout << "== calomet - calometClean: " << (calomet_handle->front()).pt() - (calometClean_handle->front()).pt() << std::endl;
140  // std::cout << "---" << std::endl;
141  m_metCases->Fill(0);
142  m_metDiff->Fill( (calomet_handle->front()).pt() - (calometClean_handle->front()).pt() );
143  m_metVSmetclean->Fill( (calomet_handle->front()).pt() , (calometClean_handle->front()).pt() );
144  }
145  else if( calomet_handle.isValid() && !calometClean_handle.isValid() )
146  {
147  m_metCases->Fill(1);
148  m_metCaseNoMetClean->Fill((calomet_handle->front()).pt());
149  }
150  else if( !calomet_handle.isValid() && calometClean_handle.isValid() )
151  {
152  m_metCases->Fill(2);
153  }
154 
155  // ## Event Selection
156  bool pass_nocut=false;
157  bool pass_twowidejets=false;
158  bool pass_etaptcuts=false;
159  bool pass_deta=false;
160  bool pass_JetIDtwojets=true;
161  bool pass_dphi=false;
162  bool pass_metFilter=true;
163  //--
164  bool pass_deta_L4=false;
165  bool pass_deta_L3=false;
166  bool pass_deta_L2=false;
167 
168  bool pass_fullsel_NOdeta=false;
169  bool pass_fullsel_detaL4=false;
170  bool pass_fullsel_detaL3=false;
171  bool pass_fullsel_detaL2=false;
172  bool pass_fullsel=false;
173 
174  // No cut
175  pass_nocut=true;
176 
177  // Two wide jets
178  if( widejets_handle->size() >= numwidejets_ )
179  {
180  // Two wide jets
181  pass_twowidejets=true;
182 
183  // Eta/pt cuts
184  if( fabs(wj1.Eta())<etawidejets_ && wj1.Pt()>ptwidejets_
185  &&
186  fabs(wj2.Eta())<etawidejets_ && wj2.Pt()>ptwidejets_
187  )
188  {
189  pass_etaptcuts=true;
190  }
191 
192  // Deta cut
193  if( DeltaEtaJJWide < detawidejets_ )
194  pass_deta=true;
195 
196  // Dphi cut
197  if( DeltaPhiJJWide > dphiwidejets_ )
198  pass_dphi=true;
199 
200  // Other Deta cuts
201  if( DeltaEtaJJWide < 4.0 )
202  pass_deta_L4=true;
203 
204  if( DeltaEtaJJWide < 3.0 )
205  pass_deta_L3=true;
206 
207  if( DeltaEtaJJWide < 2.0 )
208  pass_deta_L2=true;
209 
210  }
211 
212  // Jet id two leading jets
213 
214  if( calojets_handle->size() >= numwidejets_ )
215  {
216  // first jet
217  reco::CaloJetCollection::const_iterator thisJet = calojets_handle->begin();
218  //cout << "== thisJet1: " << thisJet->pt() << " " << thisJet->eta() << " " << thisJet->phi() << endl;
219  if( thisJet->energyFractionHadronic()>maxHADfraction_ || thisJet->emEnergyFraction()>maxEMfraction_ )
220  pass_JetIDtwojets=false;
221 
222  // second jet
223  ++thisJet;
224  //cout << "== thisJet2: " << thisJet->pt() << " " << thisJet->eta() << " " << thisJet->phi() << endl;
225  if( thisJet->energyFractionHadronic()>maxHADfraction_ || thisJet->emEnergyFraction()>maxEMfraction_ )
226  pass_JetIDtwojets=false;
227  }
228 
229  // Met filter
230  if( calomet_handle.isValid() && calometClean_handle.isValid() )
231  {
232  if( fabs ( (calomet_handle->front()).pt() - (calometClean_handle->front()).pt() ) > 0.1 )
233  pass_metFilter=false;
234  }
235 
236  // Full selection (no deta cut)
237  if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_JetIDtwojets && pass_dphi && pass_metFilter )
238  pass_fullsel_NOdeta=true;
239 
240  // Full selection (various deta cuts)
241  if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_JetIDtwojets && pass_dphi && pass_metFilter && pass_deta_L4 )
242  pass_fullsel_detaL4=true;
243  if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_JetIDtwojets && pass_dphi && pass_metFilter && pass_deta_L3 )
244  pass_fullsel_detaL3=true;
245  if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_JetIDtwojets && pass_dphi && pass_metFilter && pass_deta_L2 )
246  pass_fullsel_detaL2=true;
247 
248  // Full selection (default deta cut)
249  if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_deta && pass_JetIDtwojets && pass_dphi && pass_metFilter )
250  pass_fullsel=true;
251 
252  // ## Fill Histograms
253 
254  // Cut-flow plot
255  if( pass_nocut )
256  m_cutFlow->Fill(0);
257  if( pass_nocut && pass_twowidejets )
258  m_cutFlow->Fill(1);
259  if( pass_nocut && pass_twowidejets && pass_etaptcuts )
260  m_cutFlow->Fill(2);
261  if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_deta )
262  m_cutFlow->Fill(3);
263  if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_deta && pass_JetIDtwojets )
264  m_cutFlow->Fill(4);
265  if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_deta && pass_JetIDtwojets && pass_dphi )
266  m_cutFlow->Fill(5);
267  if( pass_fullsel )
268  m_cutFlow->Fill(6);
269 
270  // After full selection
271  if( pass_fullsel )
272  {
273  // 1D histograms
274  m_MjjWide_finalSel->Fill(MJJWide);
276  m_DetajjWide_finalSel->Fill(DeltaEtaJJWide);
277  m_DphijjWide_finalSel->Fill(DeltaPhiJJWide);
278  m_HT_finalSel->Fill(thisHT);
279  }
280 
281  // After full selection (without "noise" filters)
282  if( pass_nocut && pass_twowidejets && pass_etaptcuts && pass_deta )
283  {
286  }
287 
288  // After full selection (except DeltaEta cut)
289  if( pass_fullsel_NOdeta )
290  {
291  // 1D histograms
292  m_DetajjWide->Fill(DeltaEtaJJWide);
293 
294  if( DeltaEtaJJWide >= 0.0 && DeltaEtaJJWide < 0.5 )
295  m_MjjWide_deta_0p0_0p5->Fill(MJJWide);
296  if( DeltaEtaJJWide >= 0.5 && DeltaEtaJJWide < 1.0 )
297  m_MjjWide_deta_0p5_1p0->Fill(MJJWide);
298  if( DeltaEtaJJWide >= 1.0 && DeltaEtaJJWide < 1.5 )
299  m_MjjWide_deta_1p0_1p5->Fill(MJJWide);
300  if( DeltaEtaJJWide >= 1.5 && DeltaEtaJJWide < 2.0 )
301  m_MjjWide_deta_1p5_2p0->Fill(MJJWide);
302  if( DeltaEtaJJWide >= 2.0 && DeltaEtaJJWide < 2.5 )
303  m_MjjWide_deta_2p0_2p5->Fill(MJJWide);
304  if( DeltaEtaJJWide >= 2.5 && DeltaEtaJJWide < 3.0 )
305  m_MjjWide_deta_2p5_3p0->Fill(MJJWide);
306  if( DeltaEtaJJWide >= 3.0 )
307  m_MjjWide_deta_3p0_inf->Fill(MJJWide);
308 
309  // 2D histograms
310  m_DetajjVsMjjWide->Fill(MJJWide,DeltaEtaJJWide);
311  m_DetajjVsMjjWide_rebin->Fill(MJJWide,DeltaEtaJJWide);
312  }
313 
314 
315  // ## Get Trigger Info
316 
317  // HLT paths for DataScouting
318  // DST_HT250_v1
319  // DST_L1HTT_Or_L1MultiJet_v1
320  // DST_Mu5_HT250_v1
321  // DST_Ele8_CaloIdL_CaloIsoVL_TrkIdVL_TrkIsoVL_HT250_v1
322 
323  int HLTpathMain_fired = -1;
324  int HLTpathMonitor_fired = -1;
325 
326 
328  // invalid HLT configuration, skip the processing
329 
330  // if the L1 or HLT configurations have changed, (re)initialize the filters (including during the first event)
334 
335  // log the expanded configuration
336  // std::cout << "HLT selector configurations updated" << std::endl;
337  // std::cout << "HLTpathMain: " << *HLTpathMain_ << std::endl;
338  // std::cout << "HLTpathMonitor: " << *HLTpathMonitor_ << std::endl;
339  }
340 
341  HLTpathMain_fired = (*HLTpathMain_)(triggerConfiguration_);
342  HLTpathMonitor_fired = (*HLTpathMonitor_)(triggerConfiguration_);
343 
344  // The OR of the two should always be "1"
345  // std::cout << *HLTpathMain_ << ": " << HLTpathMain_fired << " -- " << *HLTpathMonitor_ << ": " << HLTpathMonitor_fired << std::endl;
346  }
347 
348  // ## Trigger Efficiency Curves
349 
350  //denominator - full sel NO deta cut
351  if( pass_fullsel_NOdeta && HLTpathMonitor_fired == 1 )
352  {
353  m_MjjWide_den_NOdeta->Fill(MJJWide);
354 
355  //numerator
356  if( HLTpathMain_fired == 1)
357  {
358  m_MjjWide_num_NOdeta->Fill(MJJWide);
359  }
360  }
361 
362  //denominator - full sel deta < 4.0
363  if( pass_fullsel_detaL4 && HLTpathMonitor_fired == 1 )
364  {
365  m_MjjWide_den_detaL4->Fill(MJJWide);
366 
367  //numerator
368  if( HLTpathMain_fired == 1)
369  {
370  m_MjjWide_num_detaL4->Fill(MJJWide);
371  }
372  }
373 
374  //denominator - full sel deta < 3.0
375  if( pass_fullsel_detaL3 && HLTpathMonitor_fired == 1 )
376  {
377  m_MjjWide_den_detaL3->Fill(MJJWide);
378 
379  //numerator
380  if( HLTpathMain_fired == 1)
381  {
382  m_MjjWide_num_detaL3->Fill(MJJWide);
383  }
384  }
385 
386  //denominator - full sel deta < 2.0
387  if( pass_fullsel_detaL2 && HLTpathMonitor_fired == 1 )
388  {
389  m_MjjWide_den_detaL2->Fill(MJJWide);
390 
391  //numerator
392  if( HLTpathMain_fired == 1)
393  {
394  m_MjjWide_num_detaL2->Fill(MJJWide);
395  }
396  }
397 
398  //denominator - full sel default deta cut (typically 1.3)
399  if( pass_fullsel && HLTpathMonitor_fired == 1 )
400  {
401  m_MjjWide_den->Fill(MJJWide);
402 
403  //numerator
404  if( HLTpathMain_fired == 1)
405  {
406  m_MjjWide_num->Fill(MJJWide);
407  }
408  }
409 
410 
411 }
412 
414 }
415 
416 //------------------------------------------------------------------------------
417 // Function to book the Monitoring Elements.
419 
420  // ==> TO BE UPDATED FOR sqrt(s)=8 TeV
421  const int N_mass_bins=83;
422  float massBins[N_mass_bins+1] = {1, 3, 6, 10, 16, 23, 31, 40, 50, 61, 74, 88, 103, 119, 137, 156, 176, 197, 220, 244, 270, 296, 325, 354, 386, 419, 453, 489, 526, 565, 606, 649, 693, 740, 788, 838, 890, 944, 1000, 1058, 1118, 1181, 1246, 1313, 1383, 1455, 1530, 1607, 1687, 1770, 1856, 1945, 2037, 2132, 2231, 2332, 2438, 2546, 2659, 2775, 2895, 3019, 3147, 3279, 3416, 3558, 3704, 3854, 4010, 4171, 4337, 4509, 4686, 4869, 5058, 5253, 5455, 5663, 5877, 6099, 6328, 6564, 6808, 7000};
423 
424  // 1D histograms
425  m_cutFlow = bookH1withSumw2( "h1_cutFlow",
426  "Cut Flow",
427  7,0.,7.,
428  "Cut",
429  "Number of events"
430  );
431  m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(1,"No cut");
432  m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(2,"N(WideJets)>=2");
433  m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(3,"|#eta|<2.5 , pT>30");
434  m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(4,"|#Delta#eta|<1.3");
435  m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(5,"JetID");
436  m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(6,"|#Delta#phi|>#pi/3");
437  m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(7,"|met-metClean|>0.1");
438 
439  m_MjjWide_finalSel = bookH1withSumw2( "h1_MjjWide_finalSel",
440  "M_{jj} WideJets (final selection)",
441  8000,0.,8000.,
442  "M_{jj} WideJets [GeV]",
443  "Number of events"
444  );
445 
446  m_MjjWide_finalSel_varbin = bookH1withSumw2BinArray( "h1_MjjWide_finalSel_varbin",
447  "M_{jj} WideJets (final selection)",
448  N_mass_bins, massBins,
449  "M_{jj} WideJets [GeV]",
450  "Number of events"
451  );
452 
453  m_MjjWide_finalSel_WithoutNoiseFilter = bookH1withSumw2( "h1_MjjWide_finalSel_WithoutNoiseFilter",
454  "M_{jj} WideJets (final selection, without noise filters)",
455  8000,0.,8000.,
456  "M_{jj} WideJets [GeV]",
457  "Number of events"
458  );
459 
460  m_MjjWide_finalSel_WithoutNoiseFilter_varbin = bookH1withSumw2BinArray( "h1_MjjWide_finalSel_WithoutNoiseFilter_varbin",
461  "M_{jj} WideJets (final selection, without noise filters)",
462  N_mass_bins, massBins,
463  "M_{jj} WideJets [GeV]",
464  "Number of events"
465  );
466 
467  m_MjjWide_deta_0p0_0p5 = bookH1withSumw2( "h1_MjjWide_deta_0p0_0p5",
468  "M_{jj} WideJets (0.0<=#Delta#eta<0.5)",
469  8000,0.,8000.,
470  "M_{jj} WideJets [GeV]",
471  "Number of events"
472  );
473 
474  m_MjjWide_deta_0p5_1p0 = bookH1withSumw2( "h1_MjjWide_deta_0p5_1p0",
475  "M_{jj} WideJets (0.5<=#Delta#eta<1.0)",
476  8000,0.,8000.,
477  "M_{jj} WideJets [GeV]",
478  "Number of events"
479  );
480 
481  m_MjjWide_deta_1p0_1p5 = bookH1withSumw2( "h1_MjjWide_deta_1p0_1p5",
482  "M_{jj} WideJets (1.0<=#Delta#eta<1.5)",
483  8000,0.,8000.,
484  "M_{jj} WideJets [GeV]",
485  "Number of events"
486  );
487 
488  m_MjjWide_deta_1p5_2p0 = bookH1withSumw2( "h1_MjjWide_deta_1p5_2p0",
489  "M_{jj} WideJets (1.5<=#Delta#eta<2.0)",
490  8000,0.,8000.,
491  "M_{jj} WideJets [GeV]",
492  "Number of events"
493  );
494 
495  m_MjjWide_deta_2p0_2p5 = bookH1withSumw2( "h1_MjjWide_deta_2p0_2p5",
496  "M_{jj} WideJets (2.0<=#Delta#eta<2.5)",
497  8000,0.,8000.,
498  "M_{jj} WideJets [GeV]",
499  "Number of events"
500  );
501 
502  m_MjjWide_deta_2p5_3p0 = bookH1withSumw2( "h1_MjjWide_deta_2p5_3p0",
503  "M_{jj} WideJets (2.5<=#Delta#eta<3.0)",
504  8000,0.,8000.,
505  "M_{jj} WideJets [GeV]",
506  "Number of events"
507  );
508 
509  m_MjjWide_deta_3p0_inf = bookH1withSumw2( "h1_MjjWide_deta_3p0_inf",
510  "M_{jj} WideJets (#Delta#eta>=3.0)",
511  8000,0.,8000.,
512  "M_{jj} WideJets [GeV]",
513  "Number of events"
514  );
515 
516  m_MjjWide_den_NOdeta = bookH1withSumw2( "h1_MjjWide_den_NOdeta",
517  "HLT Efficiency Studies (no deta cut)",
518  400,0.,2000.,
519  "M_{jj} WideJets [GeV]",
520  "Number of events"
521  );
522 
523  m_MjjWide_num_NOdeta = bookH1withSumw2( "h1_MjjWide_num_NOdeta",
524  "HLT Efficiency Studies (no deta cut)",
525  400,0.,2000.,
526  "M_{jj} WideJets [GeV]",
527  "Number of events"
528  );
529 
530  m_MjjWide_den_detaL4 = bookH1withSumw2( "h1_MjjWide_den_detaL4",
531  "HLT Efficiency Studies (deta cut < 4.0)",
532  400,0.,2000.,
533  "M_{jj} WideJets [GeV]",
534  "Number of events"
535  );
536 
537  m_MjjWide_num_detaL4 = bookH1withSumw2( "h1_MjjWide_num_detaL4",
538  "HLT Efficiency Studies (deta cut < 4.0)",
539  400,0.,2000.,
540  "M_{jj} WideJets [GeV]",
541  "Number of events"
542  );
543 
544  m_MjjWide_den_detaL3 = bookH1withSumw2( "h1_MjjWide_den_detaL3",
545  "HLT Efficiency Studies (deta cut < 3.0)",
546  400,0.,2000.,
547  "M_{jj} WideJets [GeV]",
548  "Number of events"
549  );
550 
551  m_MjjWide_num_detaL3 = bookH1withSumw2( "h1_MjjWide_num_detaL3",
552  "HLT Efficiency Studies (deta cut < 3.0)",
553  400,0.,2000.,
554  "M_{jj} WideJets [GeV]",
555  "Number of events"
556  );
557 
558  m_MjjWide_den_detaL2 = bookH1withSumw2( "h1_MjjWide_den_detaL2",
559  "HLT Efficiency Studies (deta cut < 2.0)",
560  400,0.,2000.,
561  "M_{jj} WideJets [GeV]",
562  "Number of events"
563  );
564 
565  m_MjjWide_num_detaL2 = bookH1withSumw2( "h1_MjjWide_num_detaL2",
566  "HLT Efficiency Studies (deta cut < 2.0)",
567  400,0.,2000.,
568  "M_{jj} WideJets [GeV]",
569  "Number of events"
570  );
571 
572  m_MjjWide_den = bookH1withSumw2( "h1_MjjWide_den",
573  "HLT Efficiency Studies (default deta cut)",
574  400,0.,2000.,
575  "M_{jj} WideJets [GeV]",
576  "Number of events"
577  );
578 
579  m_MjjWide_num = bookH1withSumw2( "h1_MjjWide_num",
580  "HLT Efficiency Studies (default deta cut)",
581  400,0.,2000.,
582  "M_{jj} WideJets [GeV]",
583  "Number of events"
584  );
585 
586  m_DetajjWide_finalSel = bookH1withSumw2( "h1_DetajjWide_finalSel",
587  "#Delta#eta_{jj} WideJets (final selection)",
588  100,0.,5.,
589  "#Delta#eta_{jj} WideJets",
590  "Number of events"
591  );
592 
593  m_DetajjWide = bookH1withSumw2( "h1_DetajjWide",
594  "#Delta#eta_{jj} WideJets (final selection except #Delta#eta cut)",
595  100,0.,5.,
596  "#Delta#eta_{jj} WideJets",
597  "Number of events"
598  );
599 
600  m_DphijjWide_finalSel = bookH1withSumw2( "h1_DphijjWide_finalSel",
601  "#Delta#phi_{jj} WideJets (final selection)",
602  100,0.,TMath::Pi()+0.0001,
603  "#Delta#phi_{jj} WideJets [rad.]",
604  "Number of events"
605  );
606 
607 
608  m_selJets_pt = bookH1withSumw2( "h1_selJets_pt",
609  "Selected CaloJets",
610  500,0.,5000.,
611  "Jet Pt [GeV]",
612  "Number of events"
613  );
614 
615  m_selJets_eta = bookH1withSumw2( "h1_selJets_eta",
616  "Selected CaloJets",
617  100,-5.,5.,
618  "#eta",
619  "Number of events"
620  );
621 
622  m_selJets_phi = bookH1withSumw2( "h1_selJets_phi",
623  "Selected CaloJets",
624  100,-TMath::Pi(),TMath::Pi(),
625  "#phi (rad.)",
626  "Number of events"
627  );
628 
629  m_selJets_hadEnergyFraction = bookH1withSumw2( "h1_selJets_hadEnergyFraction",
630  "Selected CaloJets",
631  110,0.,1.1,
632  "HAD Energy Fraction",
633  "Number of events"
634  );
635 
636  m_selJets_emEnergyFraction = bookH1withSumw2( "h1_selJets_emEnergyFraction",
637  "Selected CaloJets",
638  110,0.,1.1,
639  "EM Energy Fraction",
640  "Number of events"
641  );
642 
643  m_selJets_towersArea = bookH1withSumw2( "h1_selJets_towersArea",
644  "Selected CaloJets",
645  200,0.,2.,
646  "towers area",
647  "Number of events"
648  );
649 
650  m_metDiff = bookH1withSumw2( "h1_metDiff",
651  "Met - MetCleaned",
652  500,-1000.,1000.,
653  "met - metcleaned [GeV]",
654  "Number of events"
655  );
656 
657  m_metCases = bookH1withSumw2( "h1_metCases",
658  "Met cases",
659  3,0.,3.,
660  "case",
661  "Number of events"
662  );
663  m_metCases->getTH1()->GetXaxis()->SetBinLabel(1,"met , metclean");
664  m_metCases->getTH1()->GetXaxis()->SetBinLabel(2,"met , !metclean");
665  m_metCases->getTH1()->GetXaxis()->SetBinLabel(3,"!met , metclean");
666 
667 
668  m_metCaseNoMetClean = bookH1withSumw2( "h1_metCaseNoMetClean",
669  "Met - MetCleaned",
670  1000,0.,2000.,
671  "MET [GeV]",
672  "Number of events"
673  );
674 
675  m_HT_inclusive = bookH1withSumw2( "h1_HT_inclusive",
676  "HT (inclusive)",
677  150,0.,15000.,
678  "HT [GeV]",
679  "Number of events"
680  );
681 
682  m_HT_finalSel = bookH1withSumw2( "h1_HT_finalSel",
683  "HT (final selection)",
684  150,0.,15000.,
685  "HT [GeV]",
686  "Number of events"
687  );
688 
689 
690  // 2D histograms
691  m_DetajjVsMjjWide = bookH2withSumw2("h2_DetajjVsMjjWide",
692  "#Delta#eta_{jj} vs M_{jj} WideJets",
693  8000,0.,8000.,
694  100,0.,5.,
695  "M_{jj} WideJets [GeV]",
696  "#Delta#eta_{jj} WideJets");
697 
698  m_DetajjVsMjjWide_rebin = bookH2withSumw2("h2_DetajjVsMjjWide_rebin",
699  "#Delta#eta_{jj} vs M_{jj} WideJets",
700  400,0.,8000.,
701  50,0.,5.,
702  "M_{jj} WideJets [GeV]",
703  "#Delta#eta_{jj} WideJets");
704 
705  m_metVSmetclean = bookH2withSumw2("h2_metVSmetclean",
706  "MET clean vs MET",
707  100,0.,2000.,
708  100,0.,2000.,
709  "MET [GeV]",
710  "MET clean [GeV]");
711 
712 
713 }
714 
MonitorElement * m_HT_finalSel
const double Pi
MonitorElement * m_MjjWide_num_NOdeta
MonitorElement * m_MjjWide_num_detaL4
MonitorElement * m_MjjWide_den_detaL3
edm::InputTag metCleanCollectionTag_
virtual void bookMEs()
virtual void endRun(edm::Run const &, edm::EventSetup const &)
MonitorElement * bookH2withSumw2(const std::string &name, const std::string &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, const std::string &titleX="", const std::string &titleY="", Option_t *option="COLZ")
MonitorElement * m_HT_inclusive
MonitorElement * m_MjjWide_den_NOdeta
MonitorElement * m_MjjWide_num_detaL2
unsigned int numwidejets_
Evaluator * parse(const T &text)
MonitorElement * m_selJets_towersArea
MonitorElement * bookH1withSumw2(const std::string &name, const std::string &title, int nchX, double lowX, double highX, const std::string &titleX="", const std::string &titleY="Events", Option_t *option="E1 P")
edm::InputTag metCollectionTag_
MonitorElement * m_metCaseNoMetClean
ParameterSet const & getParameterSet(ParameterSetID const &id)
edm::InputTag jetCollectionTag_
MonitorElement * m_MjjWide_deta_0p0_0p5
triggerExpression::Evaluator * HLTpathMonitor_
triggerExpression::Data triggerConfiguration_
MonitorElement * m_metDiff
MonitorElement * m_MjjWide_num_detaL3
MonitorElement * m_MjjWide_deta_3p0_inf
MonitorElement * m_MjjWide_deta_1p0_1p5
MonitorElement * m_selJets_emEnergyFraction
DiJetVarAnalyzer(const edm::ParameterSet &)
MonitorElement * m_selJets_phi
MonitorElement * bookH1withSumw2BinArray(const std::string &name, const std::string &title, int nchX, float *xbinsize, const std::string &titleX="", const std::string &titleY="Events", Option_t *option="E1 P")
MonitorElement * m_MjjWide_num
void Fill(long long x)
virtual void analyze(const edm::Event &, const edm::EventSetup &)
MonitorElement * m_selJets_pt
virtual void beginRun(const edm::Run &, const edm::EventSetup &)
int iEvent
Definition: GenABIO.cc:243
MonitorElement * m_MjjWide_finalSel_varbin
MonitorElement * m_cutFlow
virtual void init(const Data &data)
MonitorElement * m_MjjWide_finalSel_WithoutNoiseFilter
MonitorElement * m_metVSmetclean
MonitorElement * m_MjjWide_den
edm::InputTag widejetsCollectionTag_
TH1 * getTH1(void) const
triggerExpression::Evaluator * HLTpathMain_
MonitorElement * m_DphijjWide_finalSel
MonitorElement * m_MjjWide_deta_0p5_1p0
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
tuple conf
Definition: dbtoconf.py:185
MonitorElement * m_MjjWide_deta_2p0_2p5
MonitorElement * m_MjjWide_den_detaL2
MonitorElement * m_DetajjVsMjjWide_rebin
MonitorElement * m_selJets_eta
MonitorElement * m_MjjWide_deta_2p5_3p0
MonitorElement * m_DetajjVsMjjWide
bool setEvent(const edm::Event &event, const edm::EventSetup &setup)
MonitorElement * m_MjjWide_deta_1p5_2p0
MonitorElement * m_MjjWide_finalSel_WithoutNoiseFilter_varbin
MonitorElement * m_DetajjWide_finalSel
MonitorElement * m_MjjWide_finalSel
MonitorElement * m_MjjWide_den_detaL4
MonitorElement * m_DetajjWide
virtual ~DiJetVarAnalyzer()
MonitorElement * m_selJets_hadEnergyFraction
MonitorElement * m_metCases
Definition: Run.h:33