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