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