HackAnalysis  2
ATLAS_HEPMC.h
1 #pragma once
2 
3 
4 void anotherCalcMissingMomentum_ATLAS(HepMC::GenEvent *event,HEP::Event &OutEvt)
5 {
6 
7  HEP::P4 pmiss;
8  pmiss.clear();
9 
10  HepMC::GenEvent::particle_iterator p = event->particles_begin();
11  //HepMC::FourVector HMCpmiss;
12 //HEP::P4 p4miss(0.0,0.0,0.0,0.0);
13 
14 while( p != event->particles_end())
15  {
16 
17  int pid=(*p)->pdg_id();
18  int apid=abs(pid);
19  //HEP::P4 vprod(0.0,0.0,0.0,0.0);
20  double prodZ=0.0;
21  double prodT=0.0;
22  if ( (*p)->production_vertex() )
23  {
24  HepMC::FourVector pos = (*p)->production_vertex()->position();
25  if (pos.rho() > 0.0001)
26  {
27  prodT=pos.perp();
28  prodZ=fabs(pos.z());
29 
30  // If produced outside the detector then it is no use to us
31  if((prodT > 7e3) || (prodZ > 11e3 )) { p++; continue;}
32 
33 
34 
35  }
36 
37  }
38 
39  if( isfinal(*p))
40  {
41  // Need it to be produced before end of HCAL, unless it is produced
42  // with something like a muon track (ignore that case)
43  //if((prodT > 2.95e3) || (prodZ > 5.5e3 )) { p++; continue;}
44  if(PDG::isInvisible(pid))
45  {
46  pmiss+=HEP::P4((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e());
47 
48  }
49 
50 
51 
52 
53  }
54  else
55  {
56  // want the last *visible* particle that decays to an invisible outside the HCAL (but still inside the detector)
57  // NO! Actually we want the last *invisible* metastable particle that decays outside
58  if(!PDG::isInvisible(pid)) {p++; continue;}
59  // see if it decays to an invisible particle beyond the calo
60 
61  double decT = (*p)->end_vertex()->position().perp();
62  double decz = fabs((*p)->end_vertex()->position().pz());
63  // Decay must be outside the detector, we already know it's produced inside
64 
65  // If it decays outside the muon system then it doesn't contribute to pmiss!!
66  //if((decT > 7e3 ) || (decz > 11e3)) { p++; continue;}
67 
68  // If it decays in the HCAL then it doesn't contribute, but we've already checked it's in the detector
69  //if((decT < 2.95e3 ) && (decz < 5.5e3)) { p++; continue;}
70  // change to the muon system; if it decays inside then it is no good, it has to match above
71  // Here we are only interested in invisibles produced in the detector that escape
72  if((decT < 10.0e3 ) && (decz < 13.0e3)) { p++; continue;}
73 
74 
75 
76  //
77  bool islast = true;
78  //bool invisibleproduct = false;
79 
80  for ( HepMC::GenVertex::particle_iterator des =(*p)->end_vertex()->particles_begin(HepMC::descendants); des != (*p)->end_vertex()->particles_end(HepMC::descendants);
81  ++des ) {
82  if ((*des)->pdg_id() == pid)
83  {
84  islast=false;
85  break;
86  }
87  /*
88  else
89  {
90  if(PDG::isInvisible((*des)->pdg_id()))
91  {
92  invisibleproduct=true;
93 
94  }
95  }
96  */
97  }
98 
99  if(!islast) {p++; continue;}
100  //if(invisibleproduct)
101  //{
102  pmiss+=HEP::P4((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e());
103  //}
104 
105  }
106 
107  p++;
108 
109 
110 
111  }
112  OutEvt.set_missingmom(pmiss);
113 }
114 
115 
116 
117 
118 void SortEvent_ATLASHEPMC(HepMC::GenEvent *event, HEP::Event &OutEvt, int npileups, std::vector<HEP::PileupEvent*> &pileups,double &Rjet, std::mt19937 &engine){
119 
120 fastjet::JetDefinition jetDef(fastjet::antikt_algorithm, Rjet);
121 std::vector <HEP::Jet*> outjets;
122 //std::vector <HEPUtils::Jet> outjets;
123 std::vector <fastjet::PseudoJet> inclusiveJets, jets;
124 std::vector <fastjet::PseudoJet> pileupJets, pJInputs; // for calculating rho
125 
126 std::vector <fastjet::PseudoJet> fjInputs, bquark_momenta; //particles for applying jet clustering
127  // Loop over particles in the event: kinematic distribution
128 fjInputs.resize(0);
129 pJInputs.resize(0);
130 bquark_momenta.resize(0);
131 //cout << "Starting search for jets" << endl;
132 HepMC::GenEvent::particle_iterator p = event->particles_begin();
133 
134 //HEP::P4 p4miss(0.0,0.0,0.0,0.0);
135 
136 while( p != event->particles_end())
137  {
138 
139  int pid=(*p)->pdg_id();
140  int apid=abs(pid);
141 
142  //cout << pid << endl;
143  if( isfinal(*p))
144  {
145  HEP::Particle* newpart = new HEP::Particle((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e(),pid);
146  // All particles other than invisibles and muons are jet constituents. Bollecks, really we only want photons
147  //if( pid == 12 || pid == 14 || pid == 16 || pid == 18 || pid == 13 || pid == 11)
148  if ( (*p)->production_vertex() )
149  {
150 
151  HepMC::FourVector pos = (*p)->production_vertex()->position();
152  if (pos.rho() > 0.0001)
153  {
154  // check if it's actually produced inside the detector
155  if((pos.perp() > 12e3) || (abs(pos.z()) > 23e3 )) { p++; delete newpart; continue;}
156  newpart->set_prod(HEP::P4(pos.x(),pos.y(),pos.z(),pos.t()));
157  }
158  }
159 
160  if(PDG::isHadron(pid))
161  {
162  //p4miss-=newpart->mom();
163  int pcharge=PDG::get3Q(pid);
164  //if(PDG::isCharged(pid))
165  if(pcharge != 0)
166  {
167  //OutEvt.add_charged_hadron((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e());
168  newpart->set_3Q(pcharge);
169  OutEvt.add_charged_hadron(newpart);
170  }
171  else
172  {
173  OutEvt.add_neutral_hadron((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e());
174  delete newpart;
175  }
176  }
177  else
178  {
179  newpart->set_prompt(); // sets it to be prompt, otherwise the routine deletes it!!!
180 
181 
182 
183  if(apid == 22)
184  {
185  newpart->set_3Q(0);
186 
187  }
188 
189  if((apid == 11) || (apid == 13) || (apid == 15))
190  {
191  // Some debug code when I was looking at displaced vertex decays
192  // Turned out you can get what looks like a four-electron vertex when an emitted photon
193  // produced an electron-positron pair in the same place as the main decay
194  /*
195  if((newpart->prod_vertex().pT()) > 1530 && (newpart->prod_vertex().pT() < 1532))
196  {
197  cout << "Dodgy electron " << newpart->mom() << endl;
198  DebugMyElectron(event,*p);
199  }
200  */
201  if(pid > 0)
202  {
203  newpart->set_3Q(-3);
204 
205  }
206  else
207  {
208  newpart->set_3Q(3);
209  }
210  }
211  else
212  {
213  if(apid > 25) // care about heavy stable charged particles
214  {
215  newpart->set_3Q(PDG::get3Q(pid));
216  }
217  }
218 
219  if(apid == 11 || apid == 13 || apid == 22 ) // leptons or photons that *may* be "signal," or "background" and go into jets
220  {
221  if(!isFromHadron(event, *p ))
222  {
223  OutEvt.add_particle(newpart);
224  p++;
225  continue; // don't include in jet ...
226  }
227  else
228  {
229  if(apid ==22) // photon like a neutral hadron for isolation purposes
230  {
231  OutEvt.add_neutral_hadron(newpart->mom());
232  delete newpart; // delete the photon but keep it in jet.
233  }
234  else
235  {
236  OutEvt.add_charged_hadron(newpart); // event now owns the charged 'hadron' -> don't delete it or try to add it again!
237  }
238 
239  }
240  }
241  else // care about stable BSM particles etc ...
242  {
243  OutEvt.add_particle(newpart);
244  }
245  /*
246  if( (apid != 11 && apid != 13 && apid != 22 ) || (! isFromHadron(event, *p )))
247  {
248  // keep particles except for leptons and photons that are clustered into jet from hadrons
249  OutEvt.add_particle(newpart);
250 
251  }
252  else
253  {
254  delete newpart;
255  }
256  */
257 
258  } // end if(isHadron)
259 
260 
261 
262  if( pid == 12 || pid == 14 || pid == 16 || pid == 18 ) // don't want neutrinos in jets
263  {
264 
265 
266  p++;
267  continue;
268  }
269  /*
270  if(PDG::get3Q(pid) != 0)
271  {
272  p4miss-=newpart->mom();
273  }
274  */
275  if((pid > 1000000) && (pid < 9000000)) // no bsm particles in jet ...
276  {
277 
278  p++;
279  continue;
280  }
281  HepMC::FourVector v4((*p)->momentum());
282  //cout << "v4x: " << v4.px() << endl;
283 
284 
285  fastjet::PseudoJet particleTemp(v4.px(),v4.py(),v4.pz(),v4.e());
286  //particleTemp.set_user_index(njetparts);
287  //njetparts++;
288  fjInputs.push_back(particleTemp);
289 
290  p++;
291  //constituent_pids.push_back(apid);
292  //if(njetparts != constituent_pids.size())
293  //{
294  // cout << "Error in size of pids in jets!" << endl;
295  //}
296  }
297  else
298  {
299  //cout << "check here" << endl;
300  if(!((*p)->end_vertex())) { p++; continue;}
301 
302  // Care about b-tagging and tau tagging!
303 
304  if(((*p)->end_vertex()) && (apid == 5))// use quarks as proxy, , (ought to use b hadrons, but whatever)
305  {
306  // ought to check to make sure they are part of the final process and not the incoming one, but whatever, I doubt it matters!
307  //
308  //Make sure it's the last b quark
309  bool islast = true;
310  for ( HepMC::GenVertex::particle_iterator des =(*p)->end_vertex()->particles_begin(HepMC::descendants); des != (*p)->end_vertex()->particles_end(HepMC::descendants);
311  ++des ) {
312  if (abs((*des)->pdg_id()) == 5)
313  {
314  islast=false;
315  break;
316  }
317  }
318  if(islast)
319  {
320  HepMC::FourVector v4((*p)->momentum());
321 
322  fastjet::PseudoJet tpjet(v4.px(),v4.py(),v4.pz(),v4.e());
323 
324  bquark_momenta.push_back(tpjet);
325  // cout << "adding b quark!" << endl;
326  }
327 
328  //bquark_momenta.push_back(fastjet::PseudoJet(v4.px(),v4.py(),v4.pz(),v4.e()));
329  }
330 
331  if((apid == 15) && ((*p)->end_vertex()))//taus!
332  {
333  // ought to check to make sure they are part of the final process and not the incoming one, but whatever, I doubt it matters!
334  //
335  //Make sure it's the last b quark
336  bool islast = true;
337  for ( HepMC::GenVertex::particle_iterator des =(*p)->end_vertex()->particles_begin(HepMC::descendants); des != (*p)->end_vertex()->particles_end(HepMC::descendants);
338  ++des ) {
339  if (abs((*des)->pdg_id()) == 15)
340  {
341  islast=false;
342  break;
343  }
344  }
345  if(islast)
346  {
347  HEP::Particle* newpart = new HEP::Particle((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e(),pid);
348  HepMC::FourVector pos = (*p)->production_vertex()->position();
349  newpart->set_prod(HEP::P4(pos.x(),pos.y(),pos.z(),pos.t()));
350  HepMC::FourVector dec = (*p)->end_vertex()->position();
351  newpart->set_decay(HEP::P4(dec.x(),dec.y(),dec.z(),dec.t()));
352  newpart->set_3Q(PDG::get3Q(pid));
353  newpart->set_prompt();
354  OutEvt.add_particle(newpart);
355  p++; continue;
356  }
357 
358  // end adding taus
359  }
360 
361 
362  if(isDT(*p,event))
363  {
364  // cout << "DT: " << pid << std::endl;
365  HEP::Particle* newpart = new HEP::Particle((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e(),pid);
366  HepMC::FourVector pos = (*p)->production_vertex()->position();
367  newpart->set_prod(HEP::P4(pos.x(),pos.y(),pos.z(),pos.t()));
368  HepMC::FourVector dec = (*p)->end_vertex()->position();
369  newpart->set_decay(HEP::P4(dec.x(),dec.y(),dec.z(),dec.t()));
370  newpart->set_3Q(PDG::get3Q(pid));
371  newpart->set_prompt();
372  OutEvt.add_particle(newpart);
373  p++; continue;
374  }
375 
376  if(isMetaStableHadron(*p,event))
377  {
378  HEP::Particle* newpart = new HEP::Particle((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e(),pid);
379  HepMC::FourVector pos = (*p)->production_vertex()->position();
380  newpart->set_prod(HEP::P4(pos.x(),pos.y(),pos.z(),pos.t()));
381  HepMC::FourVector dec = (*p)->end_vertex()->position();
382 
383  newpart->set_decay(HEP::P4(dec.x(),dec.y(),dec.z(),dec.t()));
384  int ch = PDG::get3Q(pid);
385  newpart->set_3Q(ch);
386  newpart->set_prompt();
387  if(ch !=0)
388  {
389 
390  OutEvt.add_charged_hadron(newpart);
391  //OutEvt.add_HSCP(newpart);
392 
393  }
394  else
395  {
396  delete newpart;
397  }
398 
399  }
400 
401 
402 
403  p++;
404 
405  //cout << "check done" << endl;
406  }
407 
408  };
409 
410 /*std::cout << "Found " << bquark_momenta.size() << " b quarks" << endl;
411 for( auto bm : bquark_momenta)
412 {
413  cout << bm.perp() << ", ";
414 
415 }
416 cout << endl;
417 */
418 
419 /*
420 for (int i = 0; i < event.size(); ++i) {
421  // Require visible particles inside detector.
422 if (! isfinal(event[i]) || event[i].isLepton() || abs(event[i].id()) == 1000024 || abs(event[i].id()) == 1000022
423  || event[i].id() == 23
424  || abs(event[i].id()) == 24){continue;}
425  fastjet::PseudoJet particleTemp(event[i].px(),event[i].py(),event[i].pz(),event[i].e());
426  particleTemp.set_user_index(i);
427  fjInputs.push_back(particleTemp);
428  }
429 
430 */
431 fastjet::ClusterSequence clustSeq(fjInputs, jetDef);
432 inclusiveJets = clustSeq.inclusive_jets(1.0); // argument is pTjet
433  jets = sorted_by_pt(inclusiveJets);
434 
435 std::uniform_int_distribution<> rd_el(0,pileups.size()-1);
436 
437 static std::normal_distribution<double> zdist(0,53.0);
438 static std::normal_distribution<double> tdist(0,160e-12);
439 
440 //auto starttime = std::chrono::high_resolution_clock::now();
441 //
442 //time_t starttime = time(NULL);
443 
444 
445 /* Enable for pileup subtraction */
446  /*
447 
448 for(int jp=0; jp< npileups; jp++)
449 {
450  int pileup_number=rd_el(engine);
451  //cout << "selecting pileup number " << pileup_number << std::endl;
452  //for(auto nhad : pileups[pileup_number]->get_neutral_hadrons())
453  // {
454  // OutEvt.add_neutral_hadron(nhad)
455  // }
456  double dz=zdist(engine);
457  double dt=tdist(engine);
458  OutEvt.add_particles(pileups[pileup_number]->translate_particles(dz,dt));
459  OutEvt.add_neutral_hadrons(pileups[pileup_number]->get_neutral_hadrons());
460  //OutEvt.add_charged_hadrons(pileups[pileup_number]->get_charged_hadrons());
461  OutEvt.add_charged_hadrons(pileups[pileup_number]->translate_charged_hadrons(dz,dt));
462  // OutEvt.add_hadron_vertices(pileups[pileup_number]->translate_hadron_vertices(dz,dt));
463 }
464 
465 
466 
467 // now to compute Jets including pileup
468 pJInputs.clear();
469 for(auto const chad : OutEvt.get_charged_hadrons())
470 {
471  fastjet::PseudoJet particleTemp(chad->mom().px(),chad->mom().py(),chad->mom().pz(),chad->E());
472  //particleTemp.set_user_index(njetparts);
473  //njetparts++;
474  pJInputs.push_back(particleTemp);
475 
476 }
477 
478 for(auto const neut : OutEvt.get_neutral_hadrons())
479 {
480  fastjet::PseudoJet particleTemp(neut->px(),neut->py(),neut->pz(),neut->E());
481  //particleTemp.set_user_index(njetparts);
482  //njetparts++;
483  pJInputs.push_back(particleTemp);
484 
485 }
486 
487 for(auto const vis : OutEvt.visible_particles())
488 {
489  fastjet::PseudoJet particleTemp(vis->mom().px(),vis->mom().py(),vis->mom().pz(),vis->mom().E());
490  //particleTemp.set_user_index(njetparts);
491  //njetparts++;
492  pJInputs.push_back(particleTemp);
493 
494 }
495 */
496 
497 
498 //fastjet::ClusterSequence PJSeq(pJInputs, jetDef);
499 //pileupJets = PJSeq.inclusive_jets(1.0); // argument is pTjet
500 //GridMedianBackgroundEstimator bge(double max_rapidity,double requested_grid_spacing);
504 //rho = bge.rho();
506 //sigma = bge.sigma();
507 
508 // 0.6 according to https://arxiv.org/pdf/1607.03663.pdf
509 // I'll use 0.5 since it's the cone size
510 
511 
512 /* // enable this if you want it!!
513 fastjet::GridMedianBackgroundEstimator bge(4.7,0.5);
514 //bge.set_particles(pileupJets);
515 bge.set_particles(pJInputs);
516 OutEvt.Average_rho=bge.rho();
517 */
518 
519 
520 
521 pileupJets.clear();
522 
523 //auto endtime = std::chrono::high_resolution_clock::now();
524 //std::chrono::duration<double> timediff = std::chrono::duration_cast<std::chrono::duration<double>>(endtime-starttime);
525 //double timediff = difftime(endtime, starttime) * 1000.0;
526 
527 /*
528 cout.precision(10);
529 cout << "Filling pileup took " << timediff.count() << std::endl;
530 cout << "Event now has " << OutEvt.get_charged_hadrons().size() << " charged hadrons, " << OutEvt.get_neutral_hadrons().size() << " neutral hadrons, " << OutEvt.visible_particles().size() << " particles " << std::endl;
531 cout << "Now filling event of n particles: " << event.size() << std::endl ;
532 */
533 
534 
535 
536 
537 // now do b tagging via method of looking if the delta_R of each jet to each initial b_quark is large enough
538 //int btagged=0;
539 //std::uniform_real_distribution<double> rd(0.0,1.0);
540 /*
541 for (auto n : electron_numbers)
542 {
543  cout << n << ", " ;
544 }
545 cout << endl;
546 */
547 
548 //cout << "Found " << bquark_momenta.size() << "bquarks" << endl;
549 //cout << "Found " << candidate_electrons.size() << " electron candidates" << endl;
550 //bool btag;
551 
552 
553 
554 
555 //cout << "smearing jets" << endl;
556 // now smear the jets
557 //smearJets(outjets,engine);
558  // now do b tagging via method of looking if the delta_R of each jet to each initial b_quark is large enough
559 bool btag;
560 
561 for(auto jet : jets)
562 {
563  btag=false;
564  for( auto const& bquark : bquark_momenta)
565  {
566  if(bquark.delta_R(jet) < Rjet) { btag = true; break;}
567  }
568 
569  HEP::Jet* tjet = new HEP::Jet(jet.px(),jet.py(),jet.pz(),jet.E(),btag);
570  outjets.push_back(tjet);
571 
572 }
573 
574 OutEvt.set_jets(outjets);
575 // now return them
576 
577 //return outjets;
578 
579 
580 }
581 
582 
583 void filleventHEPMC_ATLAS(HEP::Event &OutEvt, HepMC::GenEvent* evt, int npileups, std::vector<HEP::PileupEvent*> &pileups, std::mt19937 &engine)
584 //void filleventHEPMC_ATLAS(HEP::Event &OutEvt, HepMC::GenEvent* evt, SettingsContainer &runsettings, std::vector<HEP::PileupEvent*> &pileups, std::mt19937 &engine)
585 {
586 // Function for filling a heputils event from the HepMC one, doing the smearing etc.
587 
588 // First let's set the weight
589 // this needs correcting
590 //OutEvt.set_weights(evt->weights());
591 if(evt->weights().size() > 0)
592  {
593 
594  OutEvt.set_weight(evt->weights()[0]);
595  }
596  else
597  {
598  OutEvt.set_weight(1.0);
599  }
600 
601 
602 //int npileups = runsettings.npileup;
603 //double dR = runsettings.Rjet;
604 //HEPUtils::Event OutEvt;
605 double dR=0.4;
606 
607 SortEvent_ATLASHEPMC(evt,OutEvt,npileups,pileups,dR,engine);
608 
609 anotherCalcMissingMomentum_ATLAS(evt, OutEvt);
610 
611 
612 
613 
614 
615 }
Simple event class, separating particles into classes.
Definition: heputil.h:227
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
void clear()
Set the components to zero.
Definition: Vectors.h:82
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