HackAnalysis  2
Default_HEPMC.h
1 #pragma once
2 
3 
4 void SortEvent_DefaultHEPMC(HepMC::GenEvent *event, HEP::Event &OutEvt, int npileups, std::vector<HEP::PileupEvent*> &pileups,double &Rjet, std::mt19937 &engine){
5 
6 fastjet::JetDefinition jetDef(fastjet::antikt_algorithm, Rjet);
7 std::vector <HEP::Jet*> outjets;
8 //std::vector <HEPUtils::Jet> outjets;
9 std::vector <fastjet::PseudoJet> inclusiveJets, jets;
10 std::vector <fastjet::PseudoJet> pileupJets, pJInputs; // for calculating rho
11 
12 std::vector <fastjet::PseudoJet> fjInputs, bquark_momenta; //particles for applying jet clustering
13  // Loop over particles in the event: kinematic distribution
14 fjInputs.resize(0);
15 pJInputs.resize(0);
16 bquark_momenta.resize(0);
17 //cout << "Starting search for jets" << endl;
18 HepMC::GenEvent::particle_iterator p = event->particles_begin();
19 
20 //HEP::P4 p4miss(0.0,0.0,0.0,0.0);
21 
22 while( p != event->particles_end())
23  {
24 
25  int pid=(*p)->pdg_id();
26  int apid=abs(pid);
27 
28  //cout << pid << endl;
29  if( isfinal(*p))
30  {
31  HEP::Particle* newpart = new HEP::Particle((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e(),pid);
32  // All particles other than invisibles and muons are jet constituents. Bollecks, really we only want photons
33  //if( pid == 12 || pid == 14 || pid == 16 || pid == 18 || pid == 13 || pid == 11)
34  if ( (*p)->production_vertex() )
35  {
36 
37  HepMC::FourVector pos = (*p)->production_vertex()->position();
38  if (pos.rho() > 0.0001)
39  {
40  // check if it's actually produced inside the detector
41  if((pos.perp() > 12e3) || (abs(pos.z()) > 23e3 )) { p++; delete newpart; continue;}
42  newpart->set_prod(HEP::P4(pos.x(),pos.y(),pos.z(),pos.t()));
43  }
44  }
45 
46  if(PDG::isHadron(pid))
47  {
48  //p4miss-=newpart->mom();
49  int pcharge=PDG::get3Q(pid);
50  //if(PDG::isCharged(pid))
51  if(pcharge != 0)
52  {
53  //OutEvt.add_charged_hadron((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e());
54 
55 
56  newpart->set_3Q(pcharge);
57  OutEvt.add_charged_hadron(newpart);
58  }
59  else
60  {
61  OutEvt.add_neutral_hadron((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e());
62  delete newpart;
63  }
64  }
65  else
66  {
67  newpart->set_prompt(); // sets it to be prompt, otherwise the routine deletes it!!!
68 
69  if(apid == 22)
70  {
71  newpart->set_3Q(0);
72 
73  }
74 
75  if((apid == 11) || (apid == 13))
76  {
77  // Some debug code when I was looking at displaced vertex decays
78  // Turned out you can get what looks like a four-electron vertex when an emitted photon
79  // produced an electron-positron pair in the same place as the main decay
80  /*
81  if((newpart->prod_vertex().pT()) > 1530 && (newpart->prod_vertex().pT() < 1532))
82  {
83  cout << "Dodgy electron " << newpart->mom() << endl;
84  DebugMyElectron(event,*p);
85  }
86  */
87  if(pid > 0)
88  {
89  newpart->set_3Q(-3);
90 
91  }
92  else
93  {
94  newpart->set_3Q(3);
95  }
96  }
97  else
98  {
99  if(apid > 13) // care about heavy stable charged particles
100  {
101  newpart->set_3Q(PDG::get3Q(pid));
102  }
103  }
104 
105  OutEvt.add_particle(newpart);
106  }
107 
108 
109  if( pid == 12 || pid == 14 || pid == 16 || pid == 18 ) // include electrons and muons in jets
110  {
111 
112 
113  p++;
114  continue;
115  }
116  /*
117  if(PDG::get3Q(pid) != 0)
118  {
119  p4miss-=newpart->mom();
120  }
121  */
122  if((pid > 1000000) && (pid < 9000000))
123  {
124 
125  p++;
126  continue;
127  }
128  HepMC::FourVector v4((*p)->momentum());
129  //cout << "v4x: " << v4.px() << endl;
130 
131 
132  fastjet::PseudoJet particleTemp(v4.px(),v4.py(),v4.pz(),v4.e());
133  //particleTemp.set_user_index(njetparts);
134  //njetparts++;
135  fjInputs.push_back(particleTemp);
136 
137  p++;
138  //constituent_pids.push_back(apid);
139  //if(njetparts != constituent_pids.size())
140  //{
141  // cout << "Error in size of pids in jets!" << endl;
142  //}
143  }
144  else
145  {
146  //cout << "check here" << endl;
147  if(!((*p)->end_vertex())) { p++; continue;}
148 
149  if(isDT(*p,event))
150  {
151  // cout << "DT: " << pid << std::endl;
152  HEP::Particle* newpart = new HEP::Particle((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e(),pid);
153  HepMC::FourVector pos = (*p)->production_vertex()->position();
154  newpart->set_prod(HEP::P4(pos.x(),pos.y(),pos.z(),pos.t()));
155  HepMC::FourVector dec = (*p)->end_vertex()->position();
156  newpart->set_decay(HEP::P4(dec.x(),dec.y(),dec.z(),dec.t()));
157  newpart->set_3Q(PDG::get3Q(pid));
158  newpart->set_prompt();
159  OutEvt.add_particle(newpart);
160  p++; continue;
161  }
162 
163  if(isMetaStableHadron(*p,event))
164  {
165  HEP::Particle* newpart = new HEP::Particle((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e(),pid);
166  HepMC::FourVector pos = (*p)->production_vertex()->position();
167  newpart->set_prod(HEP::P4(pos.x(),pos.y(),pos.z(),pos.t()));
168  HepMC::FourVector dec = (*p)->end_vertex()->position();
169 
170  newpart->set_decay(HEP::P4(dec.x(),dec.y(),dec.z(),dec.t()));
171  int ch = PDG::get3Q(pid);
172  newpart->set_3Q(ch);
173  newpart->set_prompt();
174  if(ch !=0)
175  {
176 
177  OutEvt.add_charged_hadron(newpart);
178  //OutEvt.add_HSCP(newpart);
179 
180  }
181  else
182  {
183  delete newpart;
184  }
185 
186  }
187 
188  // don't care about b-tagging
189  /*
190  if(((*p)->end_vertex()) && (apid == 5))// use quarks as proxy, , (ought to use b hadrons, but whatever)
191  {
192  // ought to check to make sure they are part of the final process and not the incoming one, but whatever, I doubt it matters!
193  //
194  //Make sure it's the last b quark
195  bool islast = true;
196  for ( HepMC::GenVertex::particle_iterator des =(*p)->end_vertex()->particles_begin(HepMC::descendants); des != (*p)->end_vertex()->particles_end(HepMC::descendants);
197  ++des ) {
198  if (abs((*des)->pdg_id()) == 5)
199  {
200  islast=false;
201  break;
202  }
203  }
204  if(islast)
205  {
206  HepMC::FourVector v4((*p)->momentum());
207 
208  fastjet::PseudoJet tpjet(v4.px(),v4.py(),v4.pz(),v4.e());
209 
210  bquark_momenta.push_back(tpjet);
211  }
212 
213  //bquark_momenta.push_back(fastjet::PseudoJet(v4.px(),v4.py(),v4.pz(),v4.e()));
214  }
215  */
216 
217 
218  p++;
219 
220  //cout << "check done" << endl;
221  }
222 
223  };
224 
225 /*std::cout << "Found " << bquark_momenta.size() << " b quarks" << endl;
226 for( auto bm : bquark_momenta)
227 {
228  cout << bm.perp() << ", ";
229 
230 }
231 cout << endl;
232 */
233 
234 /*
235 for (int i = 0; i < event.size(); ++i) {
236  // Require visible particles inside detector.
237 if (! isfinal(event[i]) || event[i].isLepton() || abs(event[i].id()) == 1000024 || abs(event[i].id()) == 1000022
238  || event[i].id() == 23
239  || abs(event[i].id()) == 24){continue;}
240  fastjet::PseudoJet particleTemp(event[i].px(),event[i].py(),event[i].pz(),event[i].e());
241  particleTemp.set_user_index(i);
242  fjInputs.push_back(particleTemp);
243  }
244 
245 */
246 fastjet::ClusterSequence clustSeq(fjInputs, jetDef);
247 inclusiveJets = clustSeq.inclusive_jets(1.0); // argument is pTjet
248  jets = sorted_by_pt(inclusiveJets);
249 
250 std::uniform_int_distribution<> rd_el(0,pileups.size()-1);
251 
252 static std::normal_distribution<double> zdist(0,53.0);
253 static std::normal_distribution<double> tdist(0,160e-12);
254 
255 //auto starttime = std::chrono::high_resolution_clock::now();
256 //
257 //time_t starttime = time(NULL);
258 
259 for(int jp=0; jp< npileups; jp++)
260 {
261  int pileup_number=rd_el(engine);
262  //cout << "selecting pileup number " << pileup_number << std::endl;
263  //for(auto nhad : pileups[pileup_number]->get_neutral_hadrons())
264  // {
265  // OutEvt.add_neutral_hadron(nhad)
266  // }
267  double dz=zdist(engine);
268  double dt=tdist(engine);
269  OutEvt.add_particles(pileups[pileup_number]->translate_particles(dz,dt));
270  OutEvt.add_neutral_hadrons(pileups[pileup_number]->get_neutral_hadrons());
271  //OutEvt.add_charged_hadrons(pileups[pileup_number]->get_charged_hadrons());
272  OutEvt.add_charged_hadrons(pileups[pileup_number]->translate_charged_hadrons(dz,dt));
273  // OutEvt.add_hadron_vertices(pileups[pileup_number]->translate_hadron_vertices(dz,dt));
274 }
275 
276 // now to compute Jets including pileup
277 pJInputs.clear();
278 for(auto const chad : OutEvt.get_charged_hadrons())
279 {
280  fastjet::PseudoJet particleTemp(chad->mom().px(),chad->mom().py(),chad->mom().pz(),chad->E());
281  //particleTemp.set_user_index(njetparts);
282  //njetparts++;
283  pJInputs.push_back(particleTemp);
284 
285 }
286 
287 for(auto const neut : OutEvt.get_neutral_hadrons())
288 {
289  fastjet::PseudoJet particleTemp(neut->px(),neut->py(),neut->pz(),neut->E());
290  //particleTemp.set_user_index(njetparts);
291  //njetparts++;
292  pJInputs.push_back(particleTemp);
293 
294 }
295 
296 for(auto const vis : OutEvt.visible_particles())
297 {
298  fastjet::PseudoJet particleTemp(vis->mom().px(),vis->mom().py(),vis->mom().pz(),vis->mom().E());
299  //particleTemp.set_user_index(njetparts);
300  //njetparts++;
301  pJInputs.push_back(particleTemp);
302 
303 }
304 
305 //fastjet::ClusterSequence PJSeq(pJInputs, jetDef);
306 //pileupJets = PJSeq.inclusive_jets(1.0); // argument is pTjet
307 //GridMedianBackgroundEstimator bge(double max_rapidity,double requested_grid_spacing);
311 //rho = bge.rho();
313 //sigma = bge.sigma();
314 
315 // 0.6 according to https://arxiv.org/pdf/1607.03663.pdf
316 // I'll use 0.5 since it's the cone size
317 fastjet::GridMedianBackgroundEstimator bge(4.7,0.5);
318 //bge.set_particles(pileupJets);
319 bge.set_particles(pJInputs);
320 OutEvt.Average_rho=bge.rho();
321 
322 pileupJets.clear();
323 
324 //auto endtime = std::chrono::high_resolution_clock::now();
325 //std::chrono::duration<double> timediff = std::chrono::duration_cast<std::chrono::duration<double>>(endtime-starttime);
326 //double timediff = difftime(endtime, starttime) * 1000.0;
327 
328 /*
329 cout.precision(10);
330 cout << "Filling pileup took " << timediff.count() << std::endl;
331 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;
332 cout << "Now filling event of n particles: " << event.size() << std::endl ;
333 */
334 
335 
336 
337 
338 // now do b tagging via method of looking if the delta_R of each jet to each initial b_quark is large enough
339 //int btagged=0;
340 //std::uniform_real_distribution<double> rd(0.0,1.0);
341 /*
342 for (auto n : electron_numbers)
343 {
344  cout << n << ", " ;
345 }
346 cout << endl;
347 */
348 
349 //cout << "Found " << bquark_momenta.size() << "bquarks" << endl;
350 //cout << "Found " << candidate_electrons.size() << " electron candidates" << endl;
351 //bool btag;
352 
353 
354 
355 
356 //cout << "smearing jets" << endl;
357 // now smear the jets
358 //smearJets(outjets,engine);
359 
360 for(auto jet : jets)
361 {
362  HEP::Jet* tjet = new HEP::Jet(jet.px(),jet.py(),jet.pz(),jet.E());
363  outjets.push_back(tjet);
364 
365 }
366 
367 OutEvt.set_jets(outjets);
368 // now return them
369 
370 //return outjets;
371 
372 
373 }
374 
375 
376 void filleventHEPMC(HEP::Event &OutEvt, HepMC::GenEvent* evt, int npileups, std::vector<HEP::PileupEvent*> &pileups, std::mt19937 &engine)
377 //void filleventHEPMC(HEP::Event &OutEvt, HepMC::GenEvent* evt, SettingsContainer &runsettings, std::vector<HEP::PileupEvent*> &pileups, std::mt19937 &engine)
378 {
379 // Function for filling a heputils event from the HepMC one, doing the smearing etc.
380 
381 // First let's set the weight
382 // this needs correcting
383 //OutEvt.set_weights(evt->weights());
384 if(evt->weights().size() > 0)
385  {
386 
387  OutEvt.set_weight(evt->weights()[0]);
388  }
389  else
390  {
391  OutEvt.set_weight(1.0);
392  }
393 
394 //HEPUtils::Event OutEvt;
395 
396 //int npileups = runsettings.npileup;
397 //double dR = runsettings.Rjet;
398 double dR=0.4;
399 
400 //std::vector<HEP::Jet*> tjets;
401 
402 //cout << "Getting jets " << endl;
403 
404 //cout << "sorting event" << endl;
405 SortEvent_DefaultHEPMC(evt,OutEvt,npileups,pileups,dR,engine);
406 // fill up the jets
407 //OutEvt.set_jets(tjets);
408 
409 // fill up final state particles, compute MET
410 
411 //cout << "filling leptons" << endl;
412 
413 //betterCalcMissingMomentum(evt, OutEvt);
414 anotherCalcMissingMomentum(evt, OutEvt);
415 //calcMissingMomentum(OutEvt, engine);
416 
417 //cout << "Before returning we have " << OutEvt.electrons().size() << " electrons and " << OutEvt.muons().size() << " muons" << endl;
418 //cout << "Before smearing we have " << OutEvt.electrons().size() << " electrons and " << OutEvt.muons().size() << " muons" << endl;
419 
420 
421 /*
422 filterEfficiency(OutEvt.jets(),JetEff,engine);
423 
424 
425 filterEfficiency(OutEvt.electrons(),ElectronEff,engine);
426 filterEfficiency(OutEvt.muons(),MuonEff,engine);
427 
428 smearElectronEnergy(OutEvt.electrons(),engine);
429 smearMuonMomentum(OutEvt.muons(),engine);
430 */
431 
432 
433 //cout << "After smearing we have " << OutEvt.electrons().size() << " electrons and " << OutEvt.muons().size() << " muons" << endl;
434 //return OutEvt;
435 
436 
437 
438 
439 }
440 
441 
442 
443 
444 
445 void filleventHEPMC_FatJet(HEP::Event &OutEvt, HepMC::GenEvent* evt, int npileups, std::vector<HEP::PileupEvent*> &pileups, std::mt19937 &engine)
446 //void filleventHEPMC(HEP::Event &OutEvt, HepMC::GenEvent* evt, SettingsContainer &runsettings, std::vector<HEP::PileupEvent*> &pileups, std::mt19937 &engine)
447 {
448 
449 if(evt->weights().size() > 0)
450  {
451 
452  OutEvt.set_weight(evt->weights()[0]);
453  }
454  else
455  {
456  OutEvt.set_weight(1.0);
457  }
458 
459 double dR=1.0;
460 SortEvent_DefaultHEPMC(evt,OutEvt,npileups,pileups,dR,engine);
461 anotherCalcMissingMomentum(evt, OutEvt);
462 
463 
464 }
Simple event class, separating particles into classes.
Definition: heputil.h:227
void add_particles(const std::vector< Particle * > &ps)
Definition: heputil.h:584
std::vector< Particle * > visible_particles() const
Definition: heputil.h:614
void add_particle(Particle *p)
Definition: heputil.h:523
const std::vector< P4 * > & get_neutral_hadrons() const
Get prompt electrons.
Definition: heputil.h:422
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
void set_prompt(bool isprompt=true)
Set promptness.
Definition: Particle.h:178