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