HackAnalysis  2
TEST_Pythia.h
1 
2 
3 Pythia8::Vec4 getMissingMomentum_TEST(Pythia8::Event &event)
4 {
5  Pythia8::Vec4 missingETvec;
6  std::vector<int> metaStableMothers;
7  //First check for all final particles:
8  for (int i = 0; i < event.size(); ++i) {
9  Pythia8::Particle ptc = event[i];
10  if (!ptc.isFinal()) continue;
11  //If particle has been produced inside or after the calorimeter, skip it and include its parent
12  // These for ATLAS
13  //if (abs(ptc.vProd().pz()) > 6100. || ptc.vProd().pT() > 2300.){
14  // These for CMS
15  //if (abs(ptc.vProd().pz()) > 5500. || ptc.vProd().pT() > 2950.){
16  // CMS muon chamber. Assuming here something invisible produced by something charged
17  if (fabs(ptc.vProd().pz()) > 11000. || ptc.vProd().pT() > 7000.){
18  int mom = ptc.mother1();
19  //cout << "Adding mom " << mom << " from " << ptc.name() << " and index " << i << endl;
20  if(std::find(metaStableMothers.begin(), metaStableMothers.end(), mom) == metaStableMothers.end()) {
21  metaStableMothers.push_back(mom);
22  }
23  continue;
24  }
25 
26 
27  //Neutrinos are included in MET:
28  if (ptc.idAbs() == 12 || ptc.idAbs() == 14 || ptc.idAbs() == 16){
29  //cout << "Including: " << i << " " << ptc.name() << endl;
30  missingETvec += ptc.p();
31  continue;
32  }
33  if (ptc.isHadron()) continue; //Ignore usual hadrons
34  if (ptc.idAbs() == 22) continue; //Ignore photons
35  if (ptc.chargeType() !=0) continue; // Ignore anything charged!!
36  if (ptc.isVisible() && ptc.m() < 10.) continue; //Ignore light visible particles
37 
38  //Add to missing momentum
39  //cout << "Including: " << i << " " << ptc.name() << endl;
40  missingETvec += ptc.p();
41  }
42 
43  //Now check all meta-stable mothers:
44  for (int imom = 0; imom < metaStableMothers.size(); ++imom) {
45  Pythia8::Particle ptc = event[metaStableMothers[imom]];
46 
47 
48  //Neutrinos are included in MET:
49  if (ptc.idAbs() == 12 || ptc.idAbs() == 14 || ptc.idAbs() == 16){
50  //cout << "Including mom: " << metaStableMothers[imom] << " " << ptc.id() << endl;
51  missingETvec += ptc.p();
52  continue;
53  }
54  if (ptc.isHadron()) continue; //Ignore usual hadrons
55  if (ptc.idAbs() == 22) continue; //Ignore photons
56  if (ptc.chargeType() !=0) continue; // Ignore anything charged!!
57  if (ptc.isVisible() && ptc.m() < 10.) continue; //Ignore light visible particles
58 
59  //Add to missing momentum
60  //cout << "Including mom: " << metaStableMothers[imom] << " " << ptc.id() << endl;
61  missingETvec += ptc.p();
62  }
63 
64  return missingETvec;
65 
66 
67 };
68 
69 
70 /*
71 bool isECAL_ATLAS(HEP::P4* vprod, HEP::P4* vdec)
72 {
73 
74 
75 }
76 */
77 
78 
79 
80 
81 
82 void SortEvent_TEST(Pythia8::Event &event, HEP::Event &OutEvt, int npileups, std::vector<HEP::PileupEvent*> &pileups,double &Rjet, std::mt19937 &engine){
83 
84 fastjet::JetDefinition jetDef(fastjet::antikt_algorithm, Rjet);
85 std::vector <HEP::Jet*> outjets;
86 //std::vector <HEPUtils::Jet> outjets;
87 std::vector <fastjet::PseudoJet> inclusiveJets, jets;
88 
89 std::vector <fastjet::PseudoJet> pileupJets, pJInputs; // for calculating rho
90 pJInputs.resize(0);
91 
92 
93 
94 std::vector <fastjet::PseudoJet> fjInputs, bquark_momenta; //particles for applying jet clustering
95  // Loop over particles in the event: kinematic distribution
96 fjInputs.resize(0);
97 bquark_momenta.resize(0);
98 //cout << "Starting search for jets" << endl;
99 
100 
101 //for(int mpart = 0; mpart< event.size(); mpart++)
102 for(int mpart = event.size()-1; mpart>=0 ; mpart--)
103  {
104 
105  Pythia8::Particle* prt = &event[mpart];
106 
107  int pid=prt->id();
108  int apid=abs(pid);
109 
110 
111  if(prt->isFinal())
112  {
113  // This needs to either be added or deleted -- otherwise have a memory leak!
114  HEP::Particle* newpart = new HEP::Particle(prt->px(),prt->py(),prt->pz(),prt->e(),pid);
115 
116 
117  if(prt->isHadron() && prt->isVisible())
118  {
119  Pythia8::Vec4 pos = prt->vProd();
120 
121  if (pos.pAbs() > 0.0001)
122  {
123  newpart->set_prod(HEP::P4(pos.px(),pos.py(),pos.pz(),pos.e()));
124  }
125 
126  if(prt->isCharged())
127  {
128  //OutEvt.add_charged_hadron(newpart->mom());
129 
130  newpart->set_3Q(prt->chargeType());
131  OutEvt.add_charged_hadron(newpart);
132 
133  }
134  else
135  {
136  /*
137  if (pos.pAbs() > 0.0001)
138  {
139 
140  }
141  */
142  OutEvt.add_neutral_hadron(newpart->mom());
143  delete newpart;
144  }
145  }
146  else // not hadron -> either lepton, photon, hscp or invisible
147  {
148 
149  // look at the vertices
150  //if ( (*p)->production_vertex() ) { // i.e. does it have one: we want the position
151 
152  Pythia8::Vec4 pos = prt->vProd();
153 
154 
155  if (pos.pAbs() > 0.0001)
156  {
157  newpart->set_prod(HEP::P4(pos.px(),pos.py(),pos.pz(),pos.e()));
158  }
159 
160 
161  // }
162  //newpart->set_prompt(); // sets it to not be prompt, otherwise the routine deletes it!!!
163  if(apid == 11 || apid == 13 || apid == 22 ) // leptons or photons that *may* be "signal," or "background" and go into jets
164  {
165  //if(!fromHadron(npart,event))
166  //{
167  newpart->set_prompt(); // won't delete particle!
168 
169  if(apid == 11 || apid ==13 )
170  {
171  if(pid > 0) {newpart->set_3Q(-3);} else {newpart->set_3Q(3);}
172 
173  }
174 
175 
176  OutEvt.add_particle(newpart);
177 
178  continue; // don't include in the jet
179  /*}
180  else
181  {
182  if(apid == 22)
183  {
184  OutEvt.add_neutral_hadron(newpart->mom()); // photons contribute to calo energy like 'neutral hadrons'
185  }
186  else
187  {
188  OutEvt.add_charged_hadron(newpart->mom());
189  }
190 
191  }
192  */
193 
194 
195  }
196  else
197  {
198  // Stable particle, if it's a BSM particle let's find the charge
199  if(apid > 13)
200  {
201  newpart->set_3Q(prt->chargeType());
202  }
203  newpart->set_prompt();
204  }
205 
206 
207 
208  OutEvt.add_particle(newpart); // particle is not prompt by default, so is deleted if it comes from a hadron
209  }
210 
211 
212  if( pid == 12 || pid == 14 || pid == 16 || pid == 18 ) // include electrons and muons, not neutrinos
213  {
214  continue;
215  }
216  if((apid > 1000000) && (apid < 9000000))
217  {
218  continue;
219  }
220  //HepMC::FourVector v4((*p)->momentum());
221  Pythia8::Vec4 v4=prt->p();
222  fastjet::PseudoJet particleTemp(v4.px(),v4.py(),v4.pz(),v4.e());
223  fjInputs.push_back(particleTemp);
224  }
225  else
226  {
227  // cout << "is not final ";
228  if(ischargedTrackHadron(prt,event))
229  {
230 
231  Pythia8::Vec4 pos = prt->vProd();
232  Pythia8::Vec4 vdec = prt->vDec();
233  if(vdec.pT() < 33.0) continue; // check to make sure it doesn't decay too soon
234 
235  // Now must allocate or delete
236  HEP::Particle* newpart = new HEP::Particle(prt->px(),prt->py(),prt->pz(),prt->e(),pid);
237  if (pos.pAbs() > 0.0001)
238  {
239  newpart->set_prod(HEP::P4(pos.px(),pos.py(),pos.pz(),pos.e()));
240  }
241  newpart->set_decay(HEP::P4(vdec.px(),vdec.py(),vdec.pz(),vdec.e()));
242  newpart->set_3Q(prt->chargeType());
243 
244  OutEvt.add_charged_hadron(newpart);
245  continue;
246 
247  }
248 
249 
250 
251  if(apid == 15) // tau
252  {
253  bool islast = true;
254 
255  std::vector<int> daughters = prt->daughterList();
256  for(int di = 0; di < daughters.size(); di ++)
257  {
258  if (event[daughters[di]].idAbs() == 15)
259  {
260  islast=false;
261  break;
262  }
263 
264  }
265 
266  if(islast)
267  {
268  HEP::Particle* newpart = new HEP::Particle(prt->px(),prt->py(),prt->pz(),prt->e(),pid);
269  newpart->set_prompt();
270  Pythia8::Vec4 pos = prt->vProd();
271  if (pos.pAbs() > 0.001)
272  {
273  newpart->set_prod(HEP::P4(pos.px(),pos.py(),pos.pz(),pos.e()));
274  }
275  pos = prt->vDec();
276  if (pos.pAbs() > 0.001)
277  {
278  newpart->set_decay(HEP::P4(pos.px(),pos.py(),pos.pz(),pos.e()));
279  }
280 
281  OutEvt.add_particle(newpart);
282  continue;
283  }
284  }
285 
286 
287  if(isHSCP(prt,event))
288  {
289  //cout << "it is!" << prt->vDec().pT() << endl;
290  //cout << prt->id() << " it is!" << prt->vDec().pAbs() << " : " << prt->tau() << ", " << prt->tau0() << endl;
291 
292  HEP::Particle* newpart = new HEP::Particle(prt->px(),prt->py(),prt->pz(),prt->e(),pid);
293  Pythia8::Vec4 pos = prt->vProd();
294 
295  if (pos.pAbs() > 0.001)
296  {
297  newpart->set_prod(HEP::P4(pos.px(),pos.py(),pos.pz(),pos.e()));
298  }
299  else
300  {
301  newpart->set_prompt();
302  }
303  //pos.clear();
304  pos = prt->vDec(); // We have a particle that is
305 
306  if (pos.pAbs() > 0.0001)
307  {
308  newpart->set_decay(HEP::P4(pos.px(),pos.py(),pos.pz(),pos.e()));
309  }
310  newpart->set_3Q(prt->chargeType());
311  OutEvt.add_HSCP(newpart);
312  continue;
313  }
314 
315 
316 
317  if(apid == 5)// use quarks as proxy, , (ought to use b hadrons, but whatever)
318  {
319  // ought to check to make sure they are part of the final process and not the incoming one, but whatever, I doubt it matters!
320  //
321  //Make sure it's the last b quark
322  bool islast = true;
323 
324  std::vector<int> daughters = prt->daughterList();
325  for(int di = 0; di < daughters.size(); di ++)
326  {
327  if (event[daughters[di]].idAbs() == 5)
328  {
329  islast=false;
330  break;
331  }
332 
333  }
334 
335  if(islast)
336  {
337  Pythia8::Vec4 v4 = prt->p();
338 
339  fastjet::PseudoJet tpjet(v4.px(),v4.py(),v4.pz(),v4.e());
340 
341  bquark_momenta.push_back(tpjet);
342  }
343  }
344 
345  //bquark_momenta.push_back(fastjet::PseudoJet(v4.px(),v4.py(),v4.pz(),v4.e()));
346 
347  //cout << "check done" << endl;
348  }
349  };
350 
351 std::uniform_int_distribution<> rd_el(0,pileups.size()-1);
352 
353 static std::normal_distribution<double> zdist(0,53.0);
354 static std::normal_distribution<double> tdist(0,160e-12);
355 
356 
357 // Now we
358 
359 
360 
361 
362 
363 fastjet::ClusterSequence clustSeq(fjInputs, jetDef);
364 inclusiveJets = clustSeq.inclusive_jets(1.0); // argument is pTjet
365 jets = sorted_by_pt(inclusiveJets);
366 
367 
368 
369 
370 pJInputs.clear();
371 
372 /* Enable for pileup subtraction */
373 /*
374 for (auto const pJ : fjInputs)
375 {
376  pJInputs.push_back(pJ);
377 }
378 */
379 
380 for(int jp=0; jp< npileups; jp++)
381 {
382  int pileup_number=rd_el(engine);
383 
384  double dz=zdist(engine);
385  double dt=tdist(engine);
386 
387  // The translate routines create *NEW* vectors -> this is important because "add_particles/hadrons" takes ownership
388  OutEvt.add_particles(pileups[pileup_number]->translate_particles(dz,dt));
389  OutEvt.add_charged_hadrons(pileups[pileup_number]->translate_charged_hadrons(dz,dt));
390 
391  // The add_neutral_hadrons routine creates new'd copies
392  OutEvt.add_neutral_hadrons(pileups[pileup_number]->get_neutral_hadrons());
393 
394 
395 
396 
397  /* Enable for pileup subtraction */
398  /*
399  // Add pileup to the computation of the median rho
400  for(auto const chad : pileups[pileup_number]->get_charged_hadrons())
401  {
402  fastjet::PseudoJet particleTemp(chad->mom().px(),chad->mom().py(),chad->mom().pz(),chad->E());
403  //particleTemp.set_user_index(njetparts);
404  //njetparts++;
405  pJInputs.push_back(particleTemp);
406 
407  }
408 
409  for(auto const neut : pileups[pileup_number]->get_neutral_hadrons())
410  {
411  fastjet::PseudoJet particleTemp(neut->px(),neut->py(),neut->pz(),neut->E());
412  //particleTemp.set_user_index(njetparts);
413  //njetparts++;
414  pJInputs.push_back(particleTemp);
415 
416  }
417 
418  for(auto const vis : pileups[pileup_number]->get_particles())
419  {
420  fastjet::PseudoJet particleTemp(vis->mom().px(),vis->mom().py(),vis->mom().pz(),vis->mom().E());
421  //particleTemp.set_user_index(njetparts);
422  //njetparts++;
423  pJInputs.push_back(particleTemp);
424 
425  }
426  */
427 
428 
429 }
430 
431 
432 
433 
434 
435 //fastjet::ClusterSequence PJSeq(pJInputs, jetDef);
436 //pileupJets = PJSeq.inclusive_jets(1.0); // argument is pTjet
437 //GridMedianBackgroundEstimator bge(double max_rapidity,double requested_grid_spacing);
441 //rho = bge.rho();
443 //sigma = bge.sigma();
444 
445 // 0.6 according to https://arxiv.org/pdf/1607.03663.pdf
446 // I'll use 0.5 since it's the cone size
447 
448 /*
449 fastjet::GridMedianBackgroundEstimator bge(4.7,0.5);
450 //bge.set_particles(pileupJets);
451 bge.set_particles(pJInputs);
452 OutEvt.Average_rho=bge.rho();
453 */
454 pileupJets.clear();
455 
456 
457 
458 
459 
460 // now do b tagging via method of looking if the delta_R of each jet to each initial b_quark is large enough
461 bool btag;
462 for(auto const& jjet: jets)
463 {
464  btag=false;
465  for( auto const& bquark : bquark_momenta)
466  {
467  if(bquark.delta_R(jjet) < Rjet) { btag = true; break;}
468  }
469  // now put into a ExHEPUtils jet with the appropriate tag
470 
471  // but apply the 77% efficiency in b-tagging:
472  //smearbtag(&btag,jjet.perp(),engine);
473 
474  HEP::Jet* tjet = new HEP::Jet(jjet.px(),jjet.py(),jjet.pz(),jjet.E(),btag);
475 
476  outjets.push_back(tjet);
477 }
478 
479 
480 
481 // now smear the jets once we have efficiencies for CMS
482 //smearJets(outjets,engine);
483 
484 
485 
486 OutEvt.set_jets(outjets);
487 
488 
489 
490 
491 
492 }
493 
494 
495 
496 
497 
498 
499 
500 void filleventPYTHIA_TEST(HEP::Event &OutEvt, Pythia8::Event &evt, double weight, int npileups, std::vector<HEP::PileupEvent*> &pileups, std::mt19937 &engine)
501 {
502 // Function for filling a heputils event from the HepMC one, doing the smearing etc.
503 
504 // First let's set the weight
505 
506 OutEvt.set_weight(weight);
507 
508 double dR=0.4;
509 //std::vector<HEP::Jet*> tjets;
510 
511 SortEvent_TEST(evt,OutEvt,npileups,pileups,dR,engine);
512 // fill up the jets
513 //OutEvt.set_jets(tjets);
514 
515 // fill up final state particles, compute MET
516 
517 
518 
519 OutEvt.set_missingmom(PVtoH(getMissingMomentum_TEST(evt)));
520 
521 //OutEvt.set_missingmom(smearMET(OutEvt.missingmom(),engine));
522 
523 // add smearing once we have efficiencies for CMS
524 /*
525 filterEfficiency(OutEvt.jets(),JetEff,engine);
526 
527 
528 filterEfficiency(OutEvt.electrons(),ElectronEff,engine);
529 filterEfficiency(OutEvt.muons(),MuonEff,engine);
530 
531 smearElectronEnergy(OutEvt.electrons(),engine);
532 smearMuonMomentum(OutEvt.muons(),engine);
533 */
534 
535 
536 
537 }
Simple event class, separating particles into classes.
Definition: heputil.h:227
void add_particles(const std::vector< Particle * > &ps)
Definition: heputil.h:584
void set_missingmom(const P4 &pmiss)
Set the missing momentum vector.
Definition: heputil.h:736
void add_particle(Particle *p)
Definition: heputil.h:523
void set_jets(const std::vector< Jet * > &jets)
Set the jets collection.
Definition: heputil.h:709
Definition: Jet.h:21
A 4-momentum class for vectors.
Definition: Vectors.h:45
Definition: Particle.h:24
const P4 & mom() const
Get the 4 vector.
Definition: Particle.h:109
void set_prompt(bool isprompt=true)
Set promptness.
Definition: Particle.h:178