CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MBUEandQCDValidation.cc
Go to the documentation of this file.
1 /*Class MBUEandQCDValidation
2  *
3  * Class to fill dqm monitor elements from existing EDM file
4  *
5  * $Date: 2011/12/29 10:53:11 $
6  * $Revision: 1.4 $
7  */
8 
11 
12 #include "CLHEP/Units/defs.h"
13 #include "CLHEP/Units/PhysicalConstants.h"
14 #include "CLHEP/Units/SystemOfUnits.h"
15 
16 #include "fastjet/JetDefinition.hh"
17 #include "fastjet/ClusterSequence.hh"
18 
19 using namespace edm;
20 
22  _wmanager(iPSet),
23  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
24  genchjetCollection_(iPSet.getParameter<edm::InputTag>("genChjetsCollection")),
25  genjetCollection_(iPSet.getParameter<edm::InputTag>("genjetsCollection")),
26  verbosity_(iPSet.getUntrackedParameter<unsigned int>("verbosity",0))
27 {
28  dbe = 0;
30 
31  hepmcGPCollection.reserve(initSize);
32  hepmcCharge.reserve(initSize);
33 
35 
37 
38 }
39 
41 
42  delete theCalo;
43 
44 }
45 
47 {
48  if(dbe){
50  dbe->setCurrentFolder("Generator/MBUEandQCD");
51 
53 
54  // Number of analyzed events
55  nEvt = dbe->book1D("nEvt", "n analyzed Events", 1, 0., 1.);
56 
57  // Number of events with no forward trigger
58  nNoFwdTrig = dbe->book1D("nNoFwdTrig", "n Events no forward trigger", 1, 0., 1.);
59 
60  // Number of events with a single arm forward trigger
61  nSaFwdTrig = dbe->book1D("nSaFwdTrig", "n Events single arm forward trigger", 1, 0., 1.);
62 
63  // Number of events with b quark
64  nbquark = dbe->book1D("nbquark", "n Events with b quark", 1, 0., 1.);
65 
66  // Number of events with c and b quark
67  ncandbquark = dbe->book1D("ncandbquark", "n Events with c and b quark", 1, 0., 1.);
68 
69  // Number of events with c and no b quark
70  ncnobquark = dbe->book1D("ncnobquark", "n Events with c and no b quark", 1, 0., 1.);
71 
72 
73  // Number of selected events for QCD-09-010
74  nEvt1 = dbe->book1D("nEvt1", "n Events QCD-09-010", 1, 0., 1.);
75  // dNchdpt QCD-09-010
76  dNchdpt1 = dbe->book1D("dNchdpt1", "dNchdpt QCD-09-010", 30, 0., 6.);
77  // dNchdeta QCD-09-010
78  dNchdeta1 = dbe->book1D("dNchdeta1", "dNchdeta QCD-09-010", 10, -2.5, 2.5);
79  // Number of selected events for QCD-10-001
80 
81  nEvt2 = dbe->book1D("nEvt2", "n Events QCD-10-001", 1, 0., 1.);
82  // Leading track pt QCD-10-001
83  leadTrackpt = dbe->book1D("leadTrackpt", "leading track pt QCD-10-001", 200, 0., 100.);
84  // Leading track eta QCD-10-001
85  leadTracketa = dbe->book1D("leadTracketa", "leading track eta QCD-10-001", 50., -2.5,2.5);
86  // transverse charged particle density vs leading track pt
87  nChaDenLpt = dbe->bookProfile("nChaDenLpt", "charged density vs leading pt", 200, 0., 100., 0., 100., " ");
88  // transverse charged particle density vs leading track pt
89  sptDenLpt = dbe->bookProfile("sptDenLpt", "sum pt density vs leading pt", 200, 0., 100., 0., 300., " ");
90  // dNchdpt QCD-10-001 transverse
91  dNchdpt2 = dbe->book1D("dNchdpt2", "dNchdpt QCD-10-001", 200, 0., 100.);
92  // dNchdeta QCD-10-001 transverse
93  dNchdeta2 = dbe->book1D("dNchdeta2", "dNchdeta QCD-10-001", 50, -2.5, 2.5);
94  // nCha QCD-10-001 transverse
95  nCha = dbe->book1D("nCha", "n charged QCD-10-001", 100, 0., 100.);
96  // dNchdSpt transverse
97  dNchdSpt = dbe->book1D("dNchdSpt", "dNchdSpt QCD-10-001", 300, 0., 300.);
98  // dNchdphi
99  dNchdphi = dbe->bookProfile("dNchdphi", "dNchdphi QCD-10-001", nphiBin, -180., 180., 0., 30., " ");
100  // dSptdphi
101  dSptdphi = dbe->bookProfile("dSptdphi", "dSptdphi QCD-10-001", nphiBin, -180., 180., 0., 30., " ");
102 
103  // number of charged jets QCD-10-001
104  nChj = dbe->book1D("nChj", "n charged jets QCD-10-001", 30, 0, 30.);
105  // dNchjdeta QCD-10-001
106  dNchjdeta = dbe->book1D("dNchjdeta", "dNchjdeta QCD-10-001", 50, -2.5, 2.5);
107  // dNchjdpt QCD-10-001
108  dNchjdpt = dbe->book1D("dNchjdpt", "dNchjdpt QCD-10-001", 100, 0., 100.);
109  // leading charged jet pt QCD-10-001
110  leadChjpt = dbe->book1D("leadChjpt", "leadChjpt QCD-10-001", 100, 0., 100.);
111  // leading charged jet eta QCD-10-001
112  leadChjeta = dbe->book1D("leadChjeta", "leadChjeta QCD-10-001", 50, -2.5, 2.5);
113  // (pt1+pt2)/ptot
114  pt1pt2optotch = dbe->bookProfile("pt1pt2optotch", "sum 2 leading jets over ptot", 50, 0., 100., 0., 1., " ");
115 
116  // particle rates in tracker acceptance
117  nPPbar = dbe->book1D("nPPbar", "nPPbar QCD-10-001", 30, 0., 30.);
118  nKpm = dbe->book1D("nKpm", "nKpm QCD-10-001", 30, 0., 30.);
119  nK0s = dbe->book1D("nK0s", "nK0s QCD-10-001", 30, 0., 30.);
120  nL0 = dbe->book1D("nL0", "nL0 QCD-10-001", 30, 0., 30.);
121  nXim = dbe->book1D("nXim", "nXim QCD-10-001", 30, 0., 30.);
122  nOmega = dbe->book1D("nOmega", "nOmega QCD-10-001", 30, 0., 30.);
123 
124  pPPbar = dbe->book1D("pPPbar", "Log10(pt) PPbar QCD-10-001", 25, -2., 3.);
125  pKpm = dbe->book1D("pKpm", "Log10(pt) Kpm QCD-10-001", 25, -2., 3.);
126  pK0s = dbe->book1D("pK0s", "Log10(pt) K0s QCD-10-001", 25, -2., 3.);
127  pL0 = dbe->book1D("pL0", "Log10(pt) L0 QCD-10-001", 25, -2., 3.);
128  pXim = dbe->book1D("pXim", "Log10(pt) Xim QCD-10-001", 25, -2., 3.);
129  pOmega = dbe->book1D("pOmega", "Log10(pt) Omega QCD-10-001", 25, -2., 3.);
130 
131  // neutral rate in the barrel + HF acceptance
132  nNNbar = dbe->book1D("nNNbar", "nNNbar QCD-10-001", 30, 0., 30.);
133  nGamma = dbe->book1D("nGamma", "nGamma QCD-10-001", 50, 0., 200.);
134 
135  pNNbar = dbe->book1D("pNNbar", "Log10(pt) NNbar QCD-10-001", 25, -2., 3.);
136  pGamma = dbe->book1D("pGamma", "Log10(pt) Gamma QCD-10-001", 25, -2., 3.);
137 
138  // highest pt electron spectrum
139  elePt = dbe->book1D("elePt", "highest pt electron Log10(pt)", 30, -2., 4.);
140 
141  // highest pt muon spectrum
142  muoPt = dbe->book1D("muoPt", "highest pt muon Log10(pt)", 30, -2., 4.);
143 
144 
145  // number of selected di-jet events
146  nDijet = dbe->book1D("nDijet", "n Dijet Events", 1, 0., 1.);
147  // number of jets
148  nj = dbe->book1D("nj", "n jets ", 30, 0, 30.);
149  // dNjdeta
150  dNjdeta = dbe->book1D("dNjdeta", "dNjdeta ", 50, -5., 5.);
151  // dNjdpt
152  dNjdpt = dbe->book1D("dNjdpt", "dNjdpt ", 60, 0., 300.);
153  // (pt1+pt2)/ptot
154  pt1pt2optot = dbe->bookProfile("pt1pt2optot", "sum 2 leading jets over Et tot ", 60, 0., 300., 0., 1., " ");
155  // pt1-pt2
156  pt1pt2balance = dbe->book1D("pt1pt2balance", "2 leading jets pt difference ", 10, 0., 1.);
157  // pt1 pt2 Delta phi
158  pt1pt2Dphi = dbe->book1D("pt1pt2Dphi", "pt1 pt2 delta phi ", nphiBin, 0., 180.);
159  // pt1 pt2 invariant mass
160  pt1pt2InvM = dbe->book1D("pt1pt2InvM", "pt1 pt2 invariant mass ", 60, 0., 600.);
161  // pt3 fraction
162  pt3Frac = dbe->book1D("pt3Frac", "2 pt3 over pt1+pt2 ", 30, 0., 1.);
163  // sum of jets Et
164  sumJEt = dbe->book1D("sumJEt", "sum Jet Et ", 60, 0., 300.);
165  // fraction of missing Et over sum of jets Et
166  missEtosumJEt = dbe->book1D("missEtosumJEt", "missing Et over sumJet Et ", 30, 0., 1.);
167  // sum of final state particle Pt
168  sumPt = dbe->book1D("sumPt", "sum particle Pt ", 60, 0., 600.);
169  // sum of final state charged particle Pt
170  sumChPt = dbe->book1D("sumChPt", "sum charged particle Pt ", 60, 0., 300.);
171 
172  //Number of selected events for the HF energy flux analysis
173  nHFflow = dbe->book1D("nHFflow", "n HF flow events", 1, 0., 1.);
174  //Forward energy flow for MinBias BSC selection
175  dEdetaHFmb = dbe->bookProfile("dEdetaHFmb", "dEdeta HF MinBias", (int)CaloCellManager::nForwardEta, 0, (double)CaloCellManager::nForwardEta, 0., 300., " ");
176  //Forward energy flow for QCD dijet selection
177  dEdetaHFdj = dbe->bookProfile("dEdetaHFdj", "dEdeta HF QCD dijet", (int)CaloCellManager::nForwardEta, 0, (double)CaloCellManager::nForwardEta, 0., 300., " ");
178 
179  // FWD-10-001 like diffraction analysis
180  nHFSD = dbe->book1D("nHFSD","n single diffraction in HF", 1, 0., 1.);
181  // E-pz HF-
182  EmpzHFm = dbe->book1D("EmpzHFm", "E-pz HF- SD", 40, 0., 200.);
183  // Number of cells above threshold
184  ntHFm = dbe->book1D("ntHFm", "number of HF- tower SD", 20, 0., 20.);
185  // Energy in HF-
186  eneHFmSel = dbe->book1D("eneHFmSel", "energy in HF-", 40, 0., 200.);
187 
188  // number of jets accepted in the 'Jet-Multiplicity' analysis
189  _JM25njets = dbe->book1D("JM25njets", "n jets", 15, 0, 15.);
190  _JM25ht = dbe->book1D("JM25ht", "HT", 80, 0, 800.);
191  _JM25pt1 = dbe->book1D("JM25pt1", "pt", 40, 0, 200.);
192  _JM25pt2 = dbe->book1D("JM25pt2", "pt", 40, 0, 200.);
193  _JM25pt3 = dbe->book1D("JM25pt3", "pt", 40, 0, 200.);
194  _JM25pt4 = dbe->book1D("JM25pt4", "pt", 40, 0, 200.);
195 
196  _JM80njets = dbe->book1D("JM80njets", "n jets", 15, 0, 15.);
197  _JM80ht = dbe->book1D("JM80ht", "HT", 80, 300, 1100.);
198  _JM80pt1 = dbe->book1D("JM80pt1", "pt", 40, 60, 260.);
199  _JM80pt2 = dbe->book1D("JM80pt2", "pt", 40, 60, 260.);
200  _JM80pt3 = dbe->book1D("JM80pt3", "pt", 40, 60, 260.);
201  _JM80pt4 = dbe->book1D("JM80pt4", "pt", 40, 60, 260.);
202 
203 
204  // differential jet rates
205  djr10 = dbe->book1D("djr10", "Differential Jet Rate 1#rightarrow0", 60, -1., 5.);
206  djr21 = dbe->book1D("djr21", "Differential Jet Rate 2#rightarrow1", 60, -1., 5.);
207  djr32 = dbe->book1D("djr32", "Differential Jet Rate 3#rightarrow2", 60, -1., 5.);
208  djr43 = dbe->book1D("djr43", "Differential Jet Rate 4#rightarrow3", 60, -1., 5.);
209 
210  // sumET analysis
211  _sumEt = dbe->book1D("sumET", "Sum of stable particles Et", 150, 0, 600.);
212  _sumEt1 = dbe->book1D("sumET1", "Sum of stable particles Et (eta<0.5)", 150, 0, 200.);
213  _sumEt2 = dbe->book1D("sumET2", "Sum of stable particles Et (0.5<eta<1.0)", 150, 0, 200.);
214  _sumEt3 = dbe->book1D("sumET3", "Sum of stable particles Et (1.0<eta<1.5)", 150, 0, 200.);
215  _sumEt4 = dbe->book1D("sumET4", "Sum of stable particles Et (1.5<eta<2.0)", 150, 0, 200.);
216  _sumEt5 = dbe->book1D("sumET5", "Sum of stable particles Et (2.0<eta<5.0)", 150, 0, 200.);
217 
218 
219  }
220  return;
221 }
222 
224 void MBUEandQCDValidation::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup){
226  iSetup.getData( fPDGTable );
227  return;
228 }
229 void MBUEandQCDValidation::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){return;}
231 {
232 
235  iEvent.getByLabel(hepmcCollection_, evt);
236 
237  //Get HepMC EVENT
238  HepMC::GenEvent *myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
239 
240  double weight = _wmanager.weight(iEvent);
241 
242 
243  if ( verbosity_ > 0 ) { myGenEvent->print(); }
244 
245  double binW = 1.;
246 
247  hepmcGPCollection.clear();
248  hepmcCharge.clear();
249  for (unsigned int i = 0; i < eneInCell.size(); i++) { eneInCell[i] = 0.; }
250 
251  nEvt->Fill(0.5,weight);
252 
253  //Looping through HepMC::GenParticle collection to search for status 1 particles
254  double charge = 0.;
255  unsigned int nb = 0;
256  unsigned int nc = 0;
257  for (HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin(); iter != myGenEvent->particles_end(); ++iter){
258  if ( std::fabs((*iter)->pdg_id()) == 4 ) { nc++; }
259  if ( std::fabs((*iter)->pdg_id()) == 5 ) { nb++; }
260  if ( (*iter)->status() == 1) {
261  hepmcGPCollection.push_back(*iter);
262  const HepPDT::ParticleData* PData = fPDGTable->particle(HepPDT::ParticleID((*iter)->pdg_id()));
263  if(PData==0) { charge = -999.; }
264  else
265  charge = PData->charge();
266 
267  hepmcCharge.push_back(charge);
268 
269  if ( verbosity_ > 0 ) {
270  std::cout << "HepMC " << std::setw(14) << std::fixed << (*iter)->barcode()
271  << std::setw(14) << std::fixed << (*iter)->pdg_id()
272  << std::setw(14) << std::fixed << (*iter)->momentum().perp()
273  << std::setw(14) << std::fixed << (*iter)->momentum().eta()
274  << std::setw(14) << std::fixed << (*iter)->momentum().phi() << std::endl;
275  }
276 
277  }
278  }
279 
280  int nBSCp = 0; int nBSCm = 0; double eneHFp = 0.; double eneHFm = 0.; int nChapt05 = 0; int nChaVtx = 0;
281  for (unsigned int i = 0; i < hepmcGPCollection.size(); i++ ){
282  if ( !isNeutrino(i) ) {
283 
284  // BSC trigger
285 
286  if ( hepmcGPCollection[i]->momentum().eta() > 3.23 && hepmcGPCollection[i]->momentum().eta() < 4.65 ) { nBSCp++; }
287  if ( hepmcGPCollection[i]->momentum().eta() < -3.23 && hepmcGPCollection[i]->momentum().eta() > -4.65 ) { nBSCm++; }
288 
289  // number of charged particles in different selections
290 
291  if ( std::fabs(hepmcGPCollection[i]->momentum().eta()) < 2.5 && hepmcGPCollection[i]->momentum().perp() > 0.5 && isCharged(i) ) { nChapt05++; }
292  if ( std::fabs(hepmcGPCollection[i]->momentum().eta()) < 2.5 && hepmcGPCollection[i]->momentum().perp() > 0.1 && isCharged(i) ) { nChaVtx++; }
293  unsigned int theIndex = theCalo->getCellIndexFromAngle(hepmcGPCollection[i]->momentum().eta(),hepmcGPCollection[i]->momentum().phi());
294  if ( theIndex < CaloCellManager::nCaloCell ) eneInCell[theIndex] += hepmcGPCollection[i]->momentum().rho();
295  }
296  }
297 
298  // Forward calorimeters energy
299 
300  for (unsigned int icell = CaloCellManager::nBarrelCell+CaloCellManager::nEndcapCell; icell < CaloCellManager::nCaloCell; icell++ ) {
301  if ( theCalo->getCellFromIndex(icell)->getEtaMin() < 0. ) { eneHFm += eneInCell[icell]; }
302  else { eneHFp += eneInCell[icell]; }
303  }
304 
305  // QCD-09-010 selection
306  bool sel1 = false;
307  if ( (nBSCp > 0 || nBSCm > 0) && eneHFp >= 3. && eneHFm >= 3. ) { sel1 = true; }
308 
309  // QCD-10-001 selection
310  bool sel2 = false;
311  if ( (nBSCp >0 || nBSCm > 0) && nChaVtx >= 3 && nChapt05 > 1 ) { sel2 = true; }
312 
313  // no forward trigger selection
314  bool sel3 = false;
315  if ( nBSCp == 0 && nBSCm == 0 ) { sel3 = true; }
316 
317  // single arm forward trigger selection
318  bool sel4 = false;
319  if ( ( nBSCp>0 && nBSCm == 0 ) || ( nBSCm>0 && nBSCp == 0 ) ) { sel4 = true; }
320 
321  // BSC selection
322  bool sel5 = false;
323  if ( nBSCp > 0 && nBSCm > 0 ) { sel5 = true; }
324 
325  // basic JME-10-001, FWD-10-002 and Jet-Multiplicity selection
326  bool sel6 = false;
327  if ( sel5 && nChaVtx > 3 ) { sel6 = true; }
328 
329  // FWD-10-001 selection
330  bool sel7 = false;
331  if ( nChaVtx >= 3 && nBSCm > 0 && eneHFp < 8. ) { sel7 = true; }
332 
333  // Fill selection histograms
334  if ( sel1 ) nEvt1->Fill(0.5,weight);
335  if ( sel2 ) nEvt2->Fill(0.5,weight);
336  if ( sel3 ) nNoFwdTrig->Fill(0.5,weight);
337  if ( sel4 ) nSaFwdTrig->Fill(0.5,weight);
338  if ( sel6 ) nHFflow->Fill(0.5,weight);
339  if ( sel7 ) nHFSD->Fill(0.5,weight);
340 
341  if ( nb > 0 ) nbquark->Fill(0.5,weight);
342  if ( nb > 0 && nc > 0 ) ncandbquark->Fill(0.5,weight);
343  if ( nb == 0 && nc > 0 ) ncnobquark->Fill(0.5,weight);
344 
345  // track analyses
346  double ptMax = 0.;
347  unsigned int iMax = 0;
348  double ptot = 0.;
349  unsigned int ppbar = 0; unsigned int nnbar = 0; unsigned int kpm = 0; unsigned int k0s = 0; unsigned int l0 = 0; unsigned int gamma = 0;
350  unsigned int xim = 0; unsigned int omega = 0;
351  unsigned int ele = 0; unsigned int muo = 0;
352  unsigned int eleMax = 0;
353  unsigned int muoMax = 0;
354 
355  std::vector<double> hfMB (CaloCellManager::nForwardEta,0);
356  std::vector<double> hfDJ (CaloCellManager::nForwardEta,0);
357 
358  for (unsigned int i = 0; i < hepmcGPCollection.size(); i++ ){
359  double eta = hepmcGPCollection[i]->momentum().eta();
360  double pt = hepmcGPCollection[i]->momentum().perp();
361  int pdgId = hepmcGPCollection[i]->pdg_id();
362  if ( isCharged(i) && std::fabs(eta) < 2.5 ) {
363  if ( sel1 ) {
364  // QCD-09-010
365  binW = dNchdpt1->getTH1()->GetBinWidth(1);
366  dNchdpt1->Fill(pt,1./binW); // weight to account for the pt bin width
367  binW = dNchdeta1->getTH1()->GetBinWidth(1);
368  dNchdeta1->Fill(eta,1./binW); // weight to account for the eta bin width
369  }
370  // search for the leading track QCD-10-001
371  if ( sel2 ) {
372  if ( pt > ptMax ) { ptMax = pt; iMax = i; }
373  ptot += pt;
374 
375  // identified charged particle
376  if (std::abs(pdgId) == 2212) {
377  ppbar++;
378  pPPbar->Fill(std::log10(pt),weight);
379  }
380  else if (std::abs(pdgId) == 321) {
381  kpm++;
382  pKpm->Fill(std::log10(pt),weight);
383  }
384  else if (std::abs(pdgId) == 3312) {
385  xim++;
386  pXim->Fill(std::log10(pt),weight);
387  }
388  else if (std::abs(pdgId) == 3334) {
389  omega++;
390  pOmega->Fill(std::log10(pt),weight);
391  }
392  else if (std::abs(pdgId) == 11) {
393  ele++;
394  eleMax = i;
395  }
396  else if (std::abs(pdgId) == 13) {
397  muo++;
398  muoMax = i;
399  }
400  }
401  }
402  else if ( sel2 && isNeutral(i) && std::fabs(eta) < 2.5 ) {
403  if (std::abs(pdgId) == 310) {
404  k0s++;
405  pK0s->Fill(std::log10(pt),weight);
406  }
407  else if (std::abs(pdgId) == 3122) {
408  l0++;
409  pL0->Fill(std::log10(pt),weight);
410  }
411  }
412  else if ( sel2 && isNeutral(i) && std::fabs(eta) < 5.19 ) {
413  if (std::abs(pdgId) == 2112) {
414  nnbar++;
415  pNNbar->Fill(std::log10(pt),weight);
416  }
417  else if (std::abs(pdgId) == 22) {
418  gamma++;
419  pGamma->Fill(std::log10(pt),weight);
420  }
421  }
422  unsigned int iBin = getHFbin(eta);
423  if ( sel6 && !isNeutrino(i) && iBin < CaloCellManager::nForwardEta ) {
424  hfMB[iBin] += hepmcGPCollection[i]->momentum().rho();
425  }
426  }
427  nPPbar->Fill(ppbar,weight);
428  nNNbar->Fill(nnbar,weight);
429  nKpm->Fill(kpm,weight);
430  nK0s->Fill(k0s,weight);
431  nL0->Fill(l0,weight);
432  nXim->Fill(xim,weight);
433  nOmega->Fill(omega,weight);
434  nGamma->Fill(gamma,weight);
435 
436  if ( ele > 0 ) elePt->Fill(std::log10(hepmcGPCollection[eleMax]->momentum().perp()),weight);
437  if ( muo > 0 ) muoPt->Fill(std::log10(hepmcGPCollection[muoMax]->momentum().perp()),weight);
438 
439  leadTrackpt->Fill(hepmcGPCollection[iMax]->momentum().perp(),weight);
440  leadTracketa->Fill(hepmcGPCollection[iMax]->momentum().eta(),weight);
441 
442  std::vector<double> theEtaRanges(theCalo->getEtaRanges());
443 
444  for (unsigned int i = 0; i < CaloCellManager::nForwardEta; i++ ) {
446  dEdetaHFmb->Fill(i+0.5,hfMB[i]/binW);
447  }
448 
449  // FWD-10-001
450 
451  if ( sel7 ) {
452 
453  double empz = 0.;
454  unsigned int nCellOvTh = 0;
455  double threshold = 0.;
456 
457  for (unsigned int icell = 0; icell < eneInCell.size(); icell++ ) {
458 
459  if ( theCalo->getCellFromIndex(icell)->getSubSys() != CaloCellId::Forward ) { threshold = 3.; }
460  else { threshold = 4.; }
461 
462  if ( eneInCell[icell] > threshold ) {
463  if ( theCalo->getCellFromIndex(icell)->getSubSys() == CaloCellId::Forward ) { nCellOvTh++; }
464  empz += eneInCell[icell]*(1.-std::cos(theCalo->getCellFromIndex(icell)->getThetaCell()));
465  }
466 
467  }
468 
469  EmpzHFm->Fill(empz,weight);
470  ntHFm->Fill(nCellOvTh,weight);
471  eneHFmSel->Fill(eneHFm,weight);
472 
473  }
474 
475  // QCD-10-001
476  double phiMax = hepmcGPCollection[iMax]->momentum().phi();
477  std::vector<unsigned int> nchvsphi (nphiBin,0);
478  std::vector<double> sptvsphi (nphiBin,0.);
479  unsigned int nChaTra = 0;
480  double sptTra = 0.;
481 
482  double binPhiW = 360./nphiBin;
483  if ( sel2 ) {
484  for (unsigned int i = 0; i < hepmcGPCollection.size(); i++ ){
485  if ( isCharged(i) && std::fabs(hepmcGPCollection[i]->momentum().eta()) < 2. ) {
486  double thePhi = (hepmcGPCollection[i]->momentum().phi()-phiMax)/CLHEP::degree;
487  if ( thePhi < -180. ) { thePhi += 360.; }
488  else if ( thePhi > 180. ) { thePhi -= 360.; }
489  unsigned int thePhiBin = (int)((thePhi+180.)/binPhiW);
490  if ( thePhiBin == nphiBin ) { thePhiBin -= 1; }
491  nchvsphi[thePhiBin]++;
492  sptvsphi[thePhiBin] += hepmcGPCollection[i]->momentum().perp();
493  // analysis in the transverse region
494  if ( std::fabs(thePhi) > 60. && std::fabs(thePhi) < 120. ) {
495  nChaTra++;
496  sptTra += hepmcGPCollection[i]->momentum().perp();
497  binW = dNchdpt2->getTH1()->GetBinWidth(1);
498  dNchdpt2->Fill(hepmcGPCollection[i]->momentum().perp(),1./binW); // weight to account for the pt bin width
499  binW = dNchdeta2->getTH1()->GetBinWidth(1);
500  dNchdeta2->Fill(hepmcGPCollection[i]->momentum().eta(),1./binW); // weight to account for the eta bin width
501  }
502  }
503  }
504  nCha->Fill(nChaTra,weight);
505  binW = dNchdSpt->getTH1()->GetBinWidth(1);
506  dNchdSpt->Fill(sptTra,1.);
507  //how do one apply weights to a profile? MonitorElement doesn't allow to
508  nChaDenLpt->Fill(hepmcGPCollection[iMax]->momentum().perp(),nChaTra/4./CLHEP::twopi);
509  sptDenLpt->Fill(hepmcGPCollection[iMax]->momentum().perp(),sptTra/4./CLHEP::twopi);
510  for ( unsigned int i = 0; i < nphiBin; i++ ) {
511  double thisPhi = -180.+(i+0.5)*binPhiW;
512  dNchdphi->Fill(thisPhi,nchvsphi[i]/binPhiW/4.); // density in phi and eta
513  dSptdphi->Fill(thisPhi,sptvsphi[i]/binPhiW/4.); // density in phi and eta
514  }
515  }
516 
517  // Gather information in the charged GenJet collection
519  iEvent.getByLabel(genchjetCollection_, genChJets );
520 
521  unsigned int nJets = 0;
522  double pt1 = 0.;
523  double pt2 = 0.;
524  reco::GenJetCollection::const_iterator ij1 = genChJets->begin();
525  reco::GenJetCollection::const_iterator ij2 = genChJets->begin();
526  if ( sel2 ) {
527  for (reco::GenJetCollection::const_iterator iter=genChJets->begin();iter!=genChJets->end();++iter){
528  double eta = (*iter).eta();
529  double pt = (*iter).pt();
530  if ( verbosity_ > 0 ) {
531  std::cout << "GenJet " << std::setw(14) << std::fixed << (*iter).pt()
532  << std::setw(14) << std::fixed << (*iter).eta()
533  << std::setw(14) << std::fixed << (*iter).phi() << std::endl;
534  }
535  if ( std::fabs(eta) < 2. ) {
536  nJets++;
537  binW = dNchjdeta->getTH1()->GetBinWidth(1);
538  dNchjdeta->Fill(eta,1./binW);
539  binW = dNchjdpt->getTH1()->GetBinWidth(1);
540  dNchjdpt->Fill(pt,1./binW);
541  if ( pt >= pt1 ) { pt1 = pt; ij1 = iter; }
542  if ( pt < pt1 && pt >= pt2 ) { pt2 = pt; ij2 = iter; }
543  }
544  }
545 
546  nChj->Fill(nJets,weight);
547  if ( nJets > 0 && ij1 != genChJets->end() ) {
548  leadChjpt->Fill(pt1,weight);
549  leadChjeta->Fill((*ij1).eta(),weight);
550  if ( nJets > 1 && ij2 != genChJets->end() ) {
551  pt1pt2optotch->Fill(pt1+pt2,(pt1+pt2)/ptot);
552  }
553  }
554  }
555 
556 
557  // Gather information in the GenJet collection
559  iEvent.getByLabel(genjetCollection_, genJets );
560 
561  nJets = 0;
562  pt1 = 0.;
563  pt2 = 0.;
564  double pt3 = 0.;
565 
566 
567  // needed for Jet-Multiplicity Analysis
568  int jm25njets = 0;
569  double jm25HT = 0.;
570  double jm25pt1 = 0.;
571  double jm25pt2 = 0.;
572  double jm25pt3 = 0.;
573  double jm25pt4 = 0.;
574 
575  int jm80njets = 0;
576  double jm80HT = 0.;
577  double jm80pt1 = 0.;
578  double jm80pt2 = 0.;
579  double jm80pt3 = 0.;
580  double jm80pt4 = 0.;
581 
582 
583 
584  reco::GenJetCollection::const_iterator ij3 = genJets->begin();
585  if ( sel6 ) {
586  for (reco::GenJetCollection::const_iterator iter=genJets->begin();iter!=genJets->end();++iter){
587  double eta = (*iter).eta();
588  double pt = (*iter).pt();
589  if ( verbosity_ > 0 ) {
590  std::cout << "GenJet " << std::setw(14) << std::fixed << (*iter).pt()
591  << std::setw(14) << std::fixed << (*iter).eta()
592  << std::setw(14) << std::fixed << (*iter).phi() << std::endl;
593  }
594  if ( std::fabs(eta) < 5. ) {
595  nJets++;
596  if ( pt >= pt1 ) { pt1 = pt; ij1 = iter; }
597  if ( pt < pt1 && pt >= pt2 ) { pt2 = pt; ij2 = iter; }
598  if ( pt < pt2 && pt >= pt3 ) { pt3 = pt; ij3 = iter; }
599  }
600 
601  // find variables for Jet-Multiplicity Analysis
602  if(fabs(iter->eta()) < 3. && iter->pt()>25.) {
603  jm25njets++;
604  jm25HT += iter->pt();
605  if(iter->pt()>jm25pt1) {
606  jm25pt4 = jm25pt3;
607  jm25pt3 = jm25pt2;
608  jm25pt2 = jm25pt1;
609  jm25pt1 = iter->pt();
610  } else if(iter->pt()>jm25pt2) {
611  jm25pt4 = jm25pt3;
612  jm25pt3 = jm25pt2;
613  jm25pt2 = iter->pt();
614  } else if(iter->pt()>jm25pt3) {
615  jm25pt4 = jm25pt3;
616  jm25pt3 = iter->pt();
617  } else if(iter->pt()>jm25pt4) {
618  jm25pt4 = iter->pt();
619  }
620  // even harder jets...
621  if(iter->pt()>80.) {
622  jm80njets++;
623  jm80HT += iter->pt();
624  if(iter->pt()>jm80pt1) {
625  jm80pt4 = jm80pt3;
626  jm80pt3 = jm80pt2;
627  jm80pt2 = jm80pt1;
628  jm80pt1 = iter->pt();
629  } else if(iter->pt()>jm80pt2) {
630  jm80pt4 = jm80pt3;
631  jm80pt3 = jm80pt2;
632  jm80pt2 = iter->pt();
633  } else if(iter->pt()>jm80pt3) {
634  jm80pt4 = jm80pt3;
635  jm80pt3 = iter->pt();
636  } else if(iter->pt()>jm80pt4) {
637  jm80pt4 = iter->pt();
638  }
639  }
640 
641  }
642 
643  if(jm25njets>3) {
644  _JM25njets ->Fill(jm25njets,weight);
645  _JM25ht ->Fill(jm25HT,weight);
646  _JM25pt1 ->Fill(jm25pt1,weight);
647  _JM25pt2 ->Fill(jm25pt2,weight);
648  _JM25pt3 ->Fill(jm25pt3,weight);
649  _JM25pt4 ->Fill(jm25pt4,weight);
650  }
651  if(jm80njets>3) {
652  _JM80njets ->Fill(jm80njets,weight);
653  _JM80ht ->Fill(jm80HT,weight);
654  _JM80pt1 ->Fill(jm80pt1,weight);
655  _JM80pt2 ->Fill(jm80pt2,weight);
656  _JM80pt3 ->Fill(jm80pt3,weight);
657  _JM80pt4 ->Fill(jm80pt4,weight);
658  }
659  }
660 
661  // select a di-jet event JME-10-001 variant
662  double sumJetEt = 0; double sumPartPt = 0.; double sumChPartPt = 0.;
663  double jpx = 0; double jpy = 0;
664  if ( nJets >= 2 && ij1 != genJets->end() && ij2 != genJets->end() ) {
665  if ( (*ij1).pt() > 25. && (*ij1).pt() > 25. ) {
666  double deltaPhi = std::fabs((*ij1).phi()-(*ij2).phi())/CLHEP::degree;
667  if ( deltaPhi > 180. ) deltaPhi = 360.-deltaPhi;
668  pt1pt2Dphi->Fill(deltaPhi,weight);
669  if ( std::fabs(deltaPhi) > 2.5*CLHEP::degree ) {
670 
671  nDijet->Fill(0.5,weight);
672 
673  for (unsigned int i = 0; i < hepmcGPCollection.size(); i++ ){
674  double eta = hepmcGPCollection[i]->momentum().eta();
675  unsigned int iBin = getHFbin(eta);
676  if ( !isNeutrino(i) && iBin < CaloCellManager::nForwardEta ) {
677  hfDJ[iBin] += hepmcGPCollection[i]->momentum().rho();
678  }
679  if ( !isNeutrino(i) && std::fabs(eta) < 5. ) {
680  sumPartPt += hepmcGPCollection[i]->momentum().perp();
681  if ( isCharged(i) ) {
682  sumChPartPt += hepmcGPCollection[i]->momentum().perp();
683  }
684  }
685  }
686  for (unsigned int i = 0; i < CaloCellManager::nForwardEta; i++ ) {
688  dEdetaHFdj->Fill(i+0.5,hfDJ[i]/binW);
689  }
690 
691  double invMass = (*ij1).energy()*(*ij2).energy()-(*ij1).px()*(*ij2).px()-(*ij1).py()*(*ij2).py()-(*ij1).pz()*(*ij2).pz();
692  invMass = std::sqrt(invMass);
693  pt1pt2InvM->Fill(invMass,weight);
694 
695  sumPt->Fill(sumPartPt,weight);
696  sumChPt->Fill(sumChPartPt,weight);
697 
698  unsigned int nSelJets = 0;
699  for (reco::GenJetCollection::const_iterator iter=genJets->begin();iter!=genJets->end();++iter){
700  double pt = (*iter).pt();
701  double eta = (*iter).eta();
702  if ( std::fabs(eta) < 5. ) {
703  nSelJets++;
704  binW = dNjdeta->getTH1()->GetBinWidth(1);
705  dNjdeta->Fill(eta,1./binW*weight);
706  binW = dNjdpt->getTH1()->GetBinWidth(1);
707  dNjdpt->Fill(pt,1./binW*weight);
708  sumJetEt += (*iter).pt();
709  jpx += (*iter).px();
710  jpy += (*iter).py();
711  }
712  }
713 
714  nj->Fill(nSelJets,weight);
715  double mEt = std::sqrt(jpx*jpx+jpy*jpy);
716  sumJEt->Fill(sumJetEt,weight);
717  missEtosumJEt->Fill(mEt/sumJetEt,weight);
718 
719  if ( nSelJets >= 3 ) { pt3Frac->Fill((*ij3).pt()/(pt1+pt2),weight); }
720 
721  pt1pt2optot->Fill(pt1+pt2,(pt1+pt2)/sumJetEt);
722  pt1pt2balance->Fill((pt1-pt2)/(pt1+pt2),weight);
723  }
724  }
725  }
726  }
727 
728 
729  //compute differential jet rates
730  std::vector<const HepMC::GenParticle*> qcdActivity;
731  HepMCValidationHelper::removeIsolatedLeptons(myGenEvent, 0.2, 3., qcdActivity);
732  //HepMCValidationHelper::allStatus1(myGenEvent, qcdActivity);
733  //fill PseudoJets to use fastjet
734  std::vector<fastjet::PseudoJet> vecs;
735  int counterUser = 1;
736  std::vector<const HepMC::GenParticle*>::const_iterator iqcdact;
737  for (iqcdact = qcdActivity.begin(); iqcdact != qcdActivity.end(); ++iqcdact){
738  const HepMC::FourVector& fmom = (*iqcdact)->momentum();
739  fastjet::PseudoJet pseudoJet(fmom.px(), fmom.py(), fmom.pz(), fmom.e());
740  pseudoJet.set_user_index(counterUser);
741  vecs.push_back(pseudoJet);
742  ++counterUser;
743  }
744  //compute jets
745  fastjet::ClusterSequence cseq(vecs, fastjet::JetDefinition(fastjet::kt_algorithm, 1., fastjet::E_scheme));
746  //access the cluster sequence and get the relevant info
747  djr10->Fill(std::log10(sqrt(cseq.exclusive_dmerge(0))),weight);
748  djr21->Fill(std::log10(sqrt(cseq.exclusive_dmerge(1))),weight);
749  djr32->Fill(std::log10(sqrt(cseq.exclusive_dmerge(2))),weight);
750  djr43->Fill(std::log10(sqrt(cseq.exclusive_dmerge(3))),weight);
751 
752 
753  // compute sumEt for all stable particles
754  std::vector<const HepMC::GenParticle*> allStable;
755  HepMCValidationHelper::allStatus1(myGenEvent, allStable);
756 
757  double sumEt = 0.;
758  double sumEt1 = 0.;
759  double sumEt2 = 0.;
760  double sumEt3 = 0.;
761  double sumEt4 = 0.;
762  double sumEt5 = 0.;
763 
764  for(std::vector<const HepMC::GenParticle*>::const_iterator iter=allStable.begin();
765  iter != allStable.end(); ++iter) {
766 
767  double thisEta=fabs((*iter)->momentum().eta());
768 
769  if(thisEta < 5.) {
770  const HepMC::FourVector mom=(*iter)->momentum();
771  double px=mom.px();
772  double py=mom.py();
773  double pz=mom.pz();
774  double E=mom.e();
775  double thisSumEt = (
776  sqrt(px*px + py*py)*E /
777  sqrt(px*px + py*py + pz*pz)
778  );
779  sumEt += thisSumEt;
780  if(thisEta<1.0) sumEt1 += thisSumEt;
781  else if(thisEta<2.0) sumEt2 += thisSumEt;
782  else if(thisEta<3.0) sumEt3 += thisSumEt;
783  else if(thisEta<4.0) sumEt4 += thisSumEt;
784  else sumEt5 += thisSumEt;
785 
786  }
787  }
788 
789  if(sumEt>0.)
790  _sumEt->Fill(sumEt,weight);
791  if(sumEt1>0.)
792  _sumEt1->Fill(sumEt1,weight);
793  if(sumEt2>0.)
794  _sumEt2->Fill(sumEt2,weight);
795  if(sumEt3>0.)
796  _sumEt3->Fill(sumEt3,weight);
797  if(sumEt4>0.)
798  _sumEt4->Fill(sumEt4,weight);
799  if(sumEt5>0.)
800  _sumEt5->Fill(sumEt5,weight);
801 
802  delete myGenEvent;
803 }//analyze
804 
806 
807  bool status = false;
808  if ( hepmcGPCollection.size() < i+1 ) { return status; }
809  else { status = (hepmcCharge[i] != 0. && hepmcCharge[i] != -999.); }
810  return status;
811 
812 }
813 
815 
816  bool status = false;
817  int pdgId = std::abs(hepmcGPCollection[i]->pdg_id());
818  if ( hepmcGPCollection.size() < i+1 ) { return status; }
819  else { status = (hepmcCharge[i] == 0. && pdgId != 12 && pdgId != 14 && pdgId != 16) ; }
820  return status;
821 
822 }
823 
825 
826  bool status = false;
827  int pdgId = std::abs(hepmcGPCollection[i]->pdg_id());
828  if ( hepmcGPCollection.size() < i+1 ) { return status; }
829  else { status = (pdgId == 12 || pdgId == 14 || pdgId == 16) ; }
830  return status;
831 
832 }
833 
834 unsigned int MBUEandQCDValidation::getHFbin(double eta) {
835 
836  unsigned int iBin = 999;
837 
838  std::vector<double> theEtaRanges(theCalo->getEtaRanges());
839 
842  if ( std::fabs(eta) >= theEtaRanges[i] && std::fabs(eta) < theEtaRanges[i+1] )
844  }
845 
846  return iBin;
847 
848 }
bool isCharged(unsigned int i)
MonitorElement * ncnobquark
int i
Definition: DBlmapReader.cc:9
double getThetaCell()
Definition: CaloCellId.cc:49
MonitorElement * nbquark
static const unsigned int nphiBin
std::vector< const HepMC::GenParticle * > hepmcGPCollection
status 1 GenParticle collection
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:722
MonitorElement * sptDenLpt
MonitorElement * leadTracketa
static const unsigned int nForwardEta
MonitorElement * leadTrackpt
MBUEandQCDValidation(const edm::ParameterSet &)
static const unsigned int nBarrelCell
MonitorElement * dNchdpt1
#define abs(x)
Definition: mlp_lapack.h:159
MonitorElement * _JM80njets
DQMStore * dbe
ME&#39;s &quot;container&quot;.
void allStatus1(const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &status1)
MonitorElement * dNchjdeta
T eta() const
CaloCellManager * theCalo
manager of calorimetric cell structure
MonitorElement * dNchdeta1
static const unsigned int nEndcapEta
edm::InputTag genjetCollection_
double charge(const std::vector< uint8_t > &Ampls)
double getEtaMin() const
Definition: CaloCellId.h:26
MonitorElement * pt1pt2Dphi
MonitorElement * dNchdeta2
void getData(T &iHolder) const
Definition: EventSetup.h:67
void Fill(long long x)
MonitorElement * _JM25njets
std::vector< double > hepmcCharge
bool isNeutrino(unsigned int i)
MonitorElement * pt1pt2InvM
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.
std::vector< double > getEtaRanges()
CaloCellId * getCellFromIndex(unsigned int id)
int iEvent
Definition: GenABIO.cc:243
void removeIsolatedLeptons(const HepMC::GenEvent *all, double deltaR, double sumPt, std::vector< const HepMC::GenParticle * > &pruned)
MonitorElement * eneHFmSel
MonitorElement * nSaFwdTrig
MonitorElement * pt1pt2optot
MonitorElement * pt1pt2optotch
static const unsigned int nCaloCell
MonitorElement * leadChjpt
T sqrt(T t)
Definition: SSEVec.h:48
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
MonitorElement * nEvt1
QCD-09-010 analysis.
System getSubSys() const
Definition: CaloCellId.h:30
TH1 * getTH1(void) const
MonitorElement * bookProfile(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, const char *option="s")
Definition: DQMStore.cc:1036
MonitorElement * dEdetaHFdj
std::vector< double > eneInCell
static const unsigned int nEndcapCell
MonitorElement * pt1pt2balance
edm::InputTag hepmcCollection_
HepPDT::ParticleData ParticleData
static const unsigned int nBarrelEta
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
edm::InputTag genchjetCollection_
MonitorElement * leadChjeta
virtual void endRun(const edm::Run &, const edm::EventSetup &)
MonitorElement * nNoFwdTrig
virtual void analyze(const edm::Event &, const edm::EventSetup &)
MonitorElement * ncandbquark
unsigned int getHFbin(double eta)
bool isNeutral(unsigned int i)
T perp() const
Magnitude of transverse component.
MonitorElement * dEdetaHFmb
tuple cout
Definition: gather_cfg.py:121
int weight
Definition: histoStyle.py:50
tuple status
Definition: ntuplemaker.py:245
MonitorElement * nChaDenLpt
unsigned int getCellIndexFromAngle(double eta, double phi)
double weight(const edm::Event &)
virtual void beginRun(const edm::Run &, const edm::EventSetup &)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:434
static const unsigned int initSize
Definition: Run.h:36
Definition: DDAxes.h:10
MonitorElement * missEtosumJEt