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