HackAnalysis  2
General_HEPMC.h
1 #pragma once
2 bool isfinal(const HepMC::GenParticle* p )
3 {
4 //if ( (!(p->end_vertex())) && (p->status()==1) ) return true;
5 //if(p->status() ==3) return false;
6 
7 if ( (!(p->end_vertex())))
8 {
9  if(p->status() != 1)
10  {
11  cout << "Particle with no end vertex and status" << p->status() << endl;
12  return false;
13  }
14  return true;
15 }
16 
17 return false;
18 
19 };
20 
21 
22 /*
23 
24 HepMC: perp() gives pT, rho() gives 3-vector length
25 
26 
27 
28 */
29 
30 
31 bool isDT(const HepMC::GenParticle* p, const HepMC::GenEvent *event)
32 //Return True (False) if the particle is (is not) a HSCP
33 {
34 int pid = p->pdg_id();
35  int pch = PDG::get3Q(pid);
36 
37  //Skip neutral particles
38  if (pch == 0) return false;
39 
40 
41  // Care if it's a hadron! The metastable charged hadrons will be included elsewhere ...
42  if(PDG::isHadron(pid)) return false;
43 
44 
45  //Skip light particles
46  //if (fabs(particle->m()) < 20.) return false;
47  if(p->momentum().m() < 20.) return false;
48 
49 
50 
51  //Check if it is not an intermediate state (has itself as daughter)
52  for ( HepMC::GenVertex::particle_iterator des = p->end_vertex()->particles_begin(HepMC::descendants); des != p->end_vertex()->particles_end(HepMC::descendants);
53  ++des ) {
54  if ((*des)->pdg_id() == pid)
55  {
56  return false;
57 
58  }
59  }
60 
61 
62  // check produced and decays inside
63  if(!p->end_vertex()) return false;
64  if(!p->production_vertex()) return false;
65 
66  HepMC::FourVector pos = p->production_vertex()->position();
67  HepMC::FourVector dec = p->end_vertex()->position();
68 
69  double pperp = pos.perp();
70  double dperp = dec.perp();
71  double pz = pos.pz();
72  double dz = dec.pz();
73 
74  // check if track is long enough
75  if(((dperp - pperp) < 3.0 ) && (fabs(dz-pz) < 3.0)) return false;
76 
77  // check whether it shows up as a track
78  //if((pperp > 1100.0) || (dperp < 160.0) || (fabs(pz) > 2700.0) || (fabs(dz) > 2700.0)) return false;
79  if((pperp > 1100.0) || (dperp < 160.0) || (fabs(pz) > 2700.0)) return false;
80  //cout << "so far so good " << pid << " " << pz << ", " << dz << "," << dperp << std::endl;
81  //std::vector<int> daughters = particle->daughterListRecursive();
82  //for (int i = 0; i < daughters.size(); ++i){
83  // int idaughter = daughters[i];
84  // Pythia8::Particle daughterParticle = event[idaughter];
85  // if (daughterParticle.id() == particle->id()) return false;
86  //}
87 
88  return true;
89 }
90 
91 bool isMetaStableHadron(const HepMC::GenParticle* p, const HepMC::GenEvent *event)
92 //Return True (False) if the particle is (is not) a HSCP
93 {
94 int pid = p->pdg_id();
95 // int pch = PDG::get3Q(pid);
96  //Skip neutral particles
97 // if (pch == 0) return false;
98  // skip hadrons
99 if(!(PDG::isHadron(pid))) return false;
100 
101 
102 
103 
104  //Skip light particles
105  //if (fabs(particle->m()) < 20.) return false;
106  if(p->momentum().m() < 20.) return false;
107 
108  // We're interested in the tracker
109  //if( (p->production_vertex()->position().perp() > 1100.) || (fabs(p->production_vertex()->position().pz()) > 3100.)) {return false;}
110  // Let's say must be produced in the pixel detector, decays somewhere in the tracker or outside. That way, we include it in
111  // the isolation calculation only once
112 
113  double prodT=p->production_vertex()->position().perp();
114  double prodZ=p->production_vertex()->position().pz();
115 
116  if( (prodT > 200.) || (fabs(prodZ) > 270.)) {return false;}
117 
118 
119 
120 
121  //Check if it is not an intermediate state (has itself as daughter)
122  for ( HepMC::GenVertex::particle_iterator des = p->end_vertex()->particles_begin(HepMC::descendants); des != p->end_vertex()->particles_end(HepMC::descendants);
123  ++des ) {
124  if ((*des)->pdg_id() == pid)
125  {
126  return false;
127 
128  }
129  }
130  double decT=p->end_vertex()->position().perp();
131  double decZ=p->end_vertex()->position().pz();
132 
133  //if( (p->end_vertex()->position().perp() < 33.0)){return false;}
134  if( (decT < 300.0) ){return false;}
135 
136 
137 
138  //std::vector<int> daughters = particle->daughterListRecursive();
139  //for (int i = 0; i < daughters.size(); ++i){
140  // int idaughter = daughters[i];
141  // Pythia8::Particle daughterParticle = event[idaughter];
142  // if (daughterParticle.id() == particle->id()) return false;
143  //}
144 
145  return true;
146 }
147 
148 
149 
150 
151 
152 bool isHSCP(const HepMC::GenParticle* p, const HepMC::GenEvent *event)
153 //Return True (False) if the particle is (is not) a HSCP
154 {
155 int pid = p->pdg_id();
156  int pch = PDG::get3Q(pid);
157  //Skip neutral particles
158  if (pch == 0) return false;
159  // skip hadrons
160 if(PDG::isHadron(pid)) return false;
161 
162 
163 
164 
165  //Skip light particles
166  //if (fabs(particle->m()) < 20.) return false;
167  if(p->momentum().m() < 20.) return false;
168 
169  //if( (prodVertex.pT() > 1100.) || (abs(prodVertex.pz()) > 3100.)) {return false;}
170 
171  //Check if it is not an intermediate state (has itself as daughter)
172  for ( HepMC::GenVertex::particle_iterator des = p->end_vertex()->particles_begin(HepMC::descendants); des != p->end_vertex()->particles_end(HepMC::descendants);
173  ++des ) {
174  if ((*des)->pdg_id() == pid)
175  {
176  return false;
177 
178  }
179  }
180 
181  //std::vector<int> daughters = particle->daughterListRecursive();
182  //for (int i = 0; i < daughters.size(); ++i){
183  // int idaughter = daughters[i];
184  // Pythia8::Particle daughterParticle = event[idaughter];
185  // if (daughterParticle.id() == particle->id()) return false;
186  //}
187 
188  return true;
189 }
190 
191 /*
192 inline bool fromHadron(int n, const Pythia8::Event& evt) {
193  // Root particle is invalid
194  if (n == 0) return false;
195  const Pythia8::Particle& p = evt[n];
196  if (p.isHadron()) return true;
197  if (p.isParton()) return false; // stop the walking at the end of the hadron level
198  for (int m : p.motherList()) {
199  if (fromHadron(m, evt)) return true;
200  }
201  return false;
202  }
203 
204 bool isfinal(const HepMC::GenParticle* p )
205 {
206 if ( !p->end_vertex() && p->status()==1 ) return true;
207 return false;
208 
209 };
210 */
211 
212 /*
213 double d0Calc(const HepMC::FourVector xv, const HepMC::FourVector p)
214 //Compute the transverse impact parameter for a particle with momentum p created at vertex xv
215 //assuming no magnetic field
216 // Look at x-y plane and find clsost distance: this is |xv| sin varphi, where the angle is between the vertex
217 // and the momentum
218 // this is just given by xv cross phat, where each is represented as a 2d vector
219 {
220  double calpha = p.px()/p.perp();
221  double salpha = p.py()/p.perp();
222 
223  double a0 = xv.px()*salpha - xv.py()*calpha;
224 
225  return fabs(a0);
226 }
227 
228 double dzCalc(const HepMC::FourVector xv, const HepMC::FourVector p)
229 //Compute the longitudinal impact parameter for a particle with momentum p created at vertex xv
230 //assuming no magnetic field
231 // dz is the z value at the point of closest d0.
232 // v_min = xv - lambda \hat{p}
233 // lambda = |xv| cos varphi = xv dot \hat{p} (in 2d)
234 {
235  //double calpha = p.px()/p.pT();
236  //double salpha = p.py()/p.pT();
237 
238  //double a0 = xv.px()*salpha - xv.py()*calpha;
239  double lambda= (p.px()*xv.px()+p.py()*xv.py())/p.perp();
240  double zmin = xv.pz()-lambda*p.pz()/p.perp();
241  return fabs(zmin);
242 }
243 
244 double dz_sintheta(const HepMC::FourVector xv, const HepMC::FourVector p)
245 //Compute the longitudinal impact parameter for a particle with momentum p created at vertex xv
246 {
247  //double calpha = p.px()/p.pT();
248  //double salpha = p.py()/p.pT();
249 
250  //double a0 = xv.px()*salpha - xv.py()*calpha;
251  if(p.rho() == 0.0) return 0.0;
252  double dz = dzCalc(xv,p);
253 
254  return fabs(dz*p.perp()/p.rho());
255 }
256 */
257 
258 void DebugMyElectron(HepMC::GenEvent *event,const HepMC::GenParticle* p )
259 {
260 
261 const HepMC::GenParticle* q = new HepMC::GenParticle;
262 q=p;
263 cout << "Diagnosing" << endl;
264 int pid = p->pdg_id();
265 bool foundGenuineMother=false;
266 while((!foundGenuineMother) && (q->production_vertex()))
267 {
268  for ( HepMC::GenVertex::particle_iterator mother= q->production_vertex()->particles_begin(HepMC::parents);
269  mother != q->production_vertex()->particles_end(HepMC::parents); ++mother)
270  {
271  cout << "mother id: " << (*mother)->pdg_id() << ", ";
272  if((*mother)->pdg_id() == pid)
273  {
274  q=*mother;
275  break;
276  }
277  else
278  {
279  foundGenuineMother=true;
280  break;
281  }
282  }
283 }
284 cout << endl;
285 }
286 
287 
288 
289 
290 
291 bool SAFEisFromHadron(HepMC::GenEvent *event,const HepMC::GenParticle* p )
292 {
293  // follow the chain up until we get to a vertex where there are no mothers equal to p
294  //if ( !(*p)->production_vertex() ) return false; // weird case here, the particle has no production vertex???
295 
296  bool foundGenuineMother=false;
297  int pid = p->pdg_id();
298  const HepMC::GenParticle* q = new HepMC::GenParticle;
299 
300  int mpid;
301  bool ismotherhadron;
302  q = p;
303  while((!foundGenuineMother) && (q->production_vertex()))
304  {
305  foundGenuineMother= true;
306  for ( HepMC::GenVertex::particle_iterator mother= q->production_vertex()->particles_begin(HepMC::parents);
307  mother != q->production_vertex()->particles_end(HepMC::parents);
308  ++mother )
309  {
310  if((*mother)->pdg_id() == pid)
311  {
312  q = *mother; // the particle iterator points to the pointer to the particle
313  foundGenuineMother=false;
314  break;
315  }
316  else
317  {
318  mpid=(*mother)->pdg_id();
319  ismotherhadron=PDG::isHadron(mpid);
320  if(ismotherhadron)
321  {
322  foundGenuineMother=true;
323  q = *mother;
324  break;
325  }
326  }
327 
328  }
329  };
330 
331  // here we've found the genuine mother
332  if(ismotherhadron)
333  {
334  // cout << " Hadronic mother " << mpid << endl;
335  return true;
336  }
337 
338  /*
339  else
340  {
341  cout << " non-hadronic mother(s): " ;
342  for ( HepMC::GenVertex::particle_iterator mother= q->production_vertex()->particles_begin(HepMC::parents);
343  mother != q->production_vertex()->particles_end(HepMC::parents);
344  ++mother )
345  {
346  cout << (*mother)->pdg_id() << ", ";
347  }
348 
349  cout << endl;
350  }
351  */
352 
353  return false; //ismotherhadron;
354 
355 }
356 
357 bool isFromHadron(HepMC::GenEvent *event,const HepMC::GenParticle* p )
358 {
359  // follow the chain up until we get to a vertex where there are no mothers equal to p
360  //if ( !(*p)->production_vertex() ) return false; // weird case here, the particle has no production vertex???
361 
362  bool foundGenuineMother=false;
363  int pid = p->pdg_id();
364  int apid = abs(pid);
365 
366  if (apid < 10) return false; // quarks
367  if (apid == 21) return false; // gluon
368  if (PDG::isHadron(pid)) return true;
369 
370  // really we're only interested in electrons or photons that come from hadrons
371  // So ignore BSM particles
372  // This basically means anything with pdg code > 24 after we have already checked if it's a hadron
373 
374  if(apid > 24) return false;
375 
376  if(!(p->production_vertex())) return false;
377 
378  for ( HepMC::GenVertex::particle_iterator mother= p->production_vertex()->particles_begin(HepMC::parents);
379  mother != p->production_vertex()->particles_end(HepMC::parents);
380  ++mother )
381  {
382  if (isFromHadron(event, *mother)) return true;
383 
384  };
385 
386 
387  return false; //ismotherhadron;
388 
389 }
390 
391 
392 void calcMissingMomentum(HEP::Event &OutEvt, std::mt19937 &engine)
393 {
394  // look at the invisibles. Bit stupid I know but whatever.
395  HEP::P4 pmiss;
396  pmiss.clear();
397 
398  /*
399  if (abs(ptc.vProd().pz()) > 6100. || ptc.vProd().pT() > 2300.)
400 
401  */
402 
403  for(HEP::Particle* invisible : OutEvt.invisible_particles())
404  {
405  // std::cout << fabs(invisible->prod_vertex().pz()) << ", " << invisible->prod_vertex().pT() << std::endl;
406  //if((fabs(invisible->prod_vertex().pz()) > 6100.) || (invisible->prod_vertex().pT() > 2300.)) continue;
407  //if((fabs(invisible->prod_vertex().pz()) > 3150.) || (invisible->prod_vertex().pT() > 1200.)) continue;
408  pmiss+=invisible->mom();
409  }
410 
411  /*
412  for(HEP::Particle* invisible : OutEvt.HSCPs())
413  {
414  //if((fabs(invisible->prod_vertex().pz()) > 6100.) || (invisible->prod_vertex().pT() > 2300.)) continue;
415  // if((fabs(invisible->decay_vertex().pz()) < 6100.) || (invisible->decay_vertex().pT() < 2300.)) continue;
416  pmiss+=invisible->mom();
417  }
418  */
419  OutEvt.set_missingmom(pmiss);
420 
421  // smear here ...
422 }
423 
424 
425 
426 void anotherCalcMissingMomentum(HepMC::GenEvent *event,HEP::Event &OutEvt)
427 {
428 
429  HEP::P4 pmiss;
430  pmiss.clear();
431 
432  HepMC::GenEvent::particle_iterator p = event->particles_begin();
433  //HepMC::FourVector HMCpmiss;
434 //HEP::P4 p4miss(0.0,0.0,0.0,0.0);
435 
436 while( p != event->particles_end())
437  {
438 
439  int pid=(*p)->pdg_id();
440  int apid=abs(pid);
441  //HEP::P4 vprod(0.0,0.0,0.0,0.0);
442  double prodZ=0.0;
443  double prodT=0.0;
444  if ( (*p)->production_vertex() )
445  {
446  HepMC::FourVector pos = (*p)->production_vertex()->position();
447  if (pos.rho() > 0.0001)
448  {
449  prodT=pos.perp();
450  prodZ=fabs(pos.z());
451 
452  // If produced outside the detector then it is no use to us
453  if((prodT > 7e3) || (prodZ > 11e3 )) { p++; continue;}
454 
455 
456 
457  }
458 
459  }
460 
461  if( isfinal(*p))
462  {
463  // Need it to be produced before end of HCAL, unless it is produced
464  // with something like a muon track (ignore that case)
465  //if((prodT > 2.95e3) || (prodZ > 5.5e3 )) { p++; continue;}
466  if(PDG::isInvisible(pid))
467  {
468  pmiss+=HEP::P4((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e());
469 
470  }
471 
472 
473 
474 
475  }
476  else
477  {
478  // want the last *visible* particle that decays to an invisible outside the HCAL (but still inside the detector)
479  // NO! Actually we want the last *invisible* metastable particle that decays outside
480  if(!PDG::isInvisible(pid)) {p++; continue;}
481  // see if it decays to an invisible particle beyond the calo
482 
483  double decT = (*p)->end_vertex()->position().perp();
484  double decz = fabs((*p)->end_vertex()->position().pz());
485  // Decay must be outside the detector, we already know it's produced inside
486 
487  // If it decays outside the muon system then it doesn't contribute to pmiss!!
488  //if((decT > 7e3 ) || (decz > 11e3)) { p++; continue;}
489 
490  // If it decays in the HCAL then it doesn't contribute, but we've already checked it's in the detector
491  //if((decT < 2.95e3 ) && (decz < 5.5e3)) { p++; continue;}
492  // change to the muon system; if it decays inside then it is no good, it has to match above
493  // Here we are only interested in invisibles produced in the detector that escape
494  if((decT < 7.0e3 ) && (decz < 11.0e3)) { p++; continue;}
495 
496 
497 
498  //
499  bool islast = true;
500  //bool invisibleproduct = false;
501 
502  for ( HepMC::GenVertex::particle_iterator des =(*p)->end_vertex()->particles_begin(HepMC::descendants); des != (*p)->end_vertex()->particles_end(HepMC::descendants);
503  ++des ) {
504  if ((*des)->pdg_id() == pid)
505  {
506  islast=false;
507  break;
508  }
509  /*
510  else
511  {
512  if(PDG::isInvisible((*des)->pdg_id()))
513  {
514  invisibleproduct=true;
515 
516  }
517  }
518  */
519  }
520 
521  if(!islast) {p++; continue;}
522  //if(invisibleproduct)
523  //{
524  pmiss+=HEP::P4((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e());
525  //}
526 
527  }
528 
529  p++;
530 
531 
532 
533  }
534  OutEvt.set_missingmom(pmiss);
535 }
536 
537 
538 void betterCalcMissingMomentum(HepMC::GenEvent *event,HEP::Event &OutEvt)
539 {
540  HEP::P4 pmiss;
541  pmiss.clear();
542 
543  HepMC::GenEvent::particle_iterator p = event->particles_begin();
544  //HepMC::FourVector HMCpmiss;
545 //HEP::P4 p4miss(0.0,0.0,0.0,0.0);
546 
547 while( p != event->particles_end())
548  {
549 
550  int pid=(*p)->pdg_id();
551  int apid=abs(pid);
552  //HEP::P4 vprod(0.0,0.0,0.0,0.0);
553  double prodZ=0.0;
554  double prodT=0.0;
555  if ( (*p)->production_vertex() )
556  {
557  HepMC::FourVector pos = (*p)->production_vertex()->position();
558  if (pos.rho() > 0.0001)
559  {
560  prodT=pos.perp();
561  prodZ=fabs(pos.z());
562  if((prodT > 12e3) || (prodZ > 23e3 )) { p++; continue;}
563 
564 
565 
566  }
567 
568  }
569 
570 
571 
572  //cout << pid << endl;
573  if( isfinal(*p))
574  {
575 
576  if(PDG::isHadron(pid))
577  {
578  // make sure produced before end of HCAL
579  if((prodT > 3000.0)|| (prodZ > 5000.0)) {p++; continue;}
580  pmiss-=HEP::P4((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e());
581  }
582  else
583  {
584  if( (apid == 22) || (PDG::get3Q(pid) !=0 ) )
585  {
586  if((prodT > 1700.0)|| (prodZ > 4000.0)) {p++; continue;}
587  pmiss-=HEP::P4((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e());
588  }
589  }
590 
591  }
592  else
593  {
594  // check last of its kind
595  if((!(*p)->end_vertex()) || (!(*p)->production_vertex())) {p++; continue;}
596 
597  bool islast = true;
598  for ( HepMC::GenVertex::particle_iterator des =(*p)->end_vertex()->particles_begin(HepMC::descendants); des != (*p)->end_vertex()->particles_end(HepMC::descendants);
599  ++des ) {
600  if ((*des)->pdg_id() == pid)
601  {
602  islast=false;
603  break;
604  }
605  }
606 
607  if(!islast) { p++; continue;}
608 
609  double decZ= fabs((*p)->end_vertex()->position().pz());
610  double decT=(*p)->end_vertex()->position().perp();
611 
612  if(PDG::isHadron(pid))
613  {
614  if((prodT > 3000.0)|| (prodZ > 5000.0) || (decZ < 4e3 ) ||(decT < 1.77e3)) {p++; continue;}
615  pmiss-=HEP::P4((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e());
616  }
617  else
618  {
619  //if( (apid == 22) || (apid == 11) )
620  if(PDG::isCharged(pid))
621  {
622  // treat as a muon
623  //if((prodT > 1700.0)|| (prodZ > 4000.0) || (decZ < 3.15e3 ) ||(decT < 1.2e3)) {p++; continue;}
624  if( (prodT > 1.77e3) || (prodZ > 4000.0) || (decZ < 1.2e4) || (decT < 7.0e3)){ p++; continue;}
625  pmiss-=HEP::P4((*p)->momentum().px(),(*p)->momentum().py(),(*p)->momentum().pz(),(*p)->momentum().e());
626  }
627  }
628 
629  }
630 
631 
632  p++;
633  }
634 
635  //HEP::P4 pmiss(HMCpmiss.px(),HMCpmiss.py(),HMCpmiss.pz(),HMCpmiss.E());
636  OutEvt.set_missingmom(pmiss);
637 
638 
639 }
640 
641 
642 
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
const std::vector< Particle * > & invisible_particles() const
Get invisible final state particles.
Definition: heputil.h:641
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