HackAnalysis  2
heputil.h
1 
2 // Following code contains some code from Event.h in HEPUtils:
3 // -*- C++ -*-
4 //
5 // This file is part of HEPUtils -- https://bitbucket.org/andybuckley/heputils
6 // Copyright (C) 2013-2018 Andy Buckley <andy.buckley@cern.ch>
7 //
8 // Embedding of HEPUtils code in other projects is permitted provided this
9 // notice is retained and the HEPUtils namespace and include path are changed.
10 //
11 
12 
25 #pragma once
26 
27 #include "ExHEPUtils/BinnedFn.h"
28 #include "ExHEPUtils/Jet.h"
29 #include "ExHEPUtils/MathUtils.h"
30 #include "ExHEPUtils/Vectors.h"
31 #include "ExHEPUtils/Particle.h"
32 #include "pdg.h"
33 
34 //#include "useYoda.h"
35 
36 
37 
38 
39 
40 
41 namespace HEP {
42 
43  // Need a function to calculate the ET from a vector (defined as E sin theta). Honestly I don't know why that is not there
44  // can therefore define as E * pT/p = E*pT/sqrt(pT^2 + pz^2) = E/sqrt(1 + pz^2/pT^2)
45  // in fastjet this is inline double Et() const {return (_kt2==0) ? 0.0 : _E/sqrt(1.0+_pz*_pz/_kt2);}
46  // in HEPUtils we store the mass of the four-vector instead of the energy (go figure) so we use instead
47  // pT * sqrt((1 + m^2/p^2)) provided that p^2 is not zero
48 
49 
50 
51 
52  inline double ET2(P4 p4)
53  {
55  double p2 = p4.p2();
56  if(p2 == 0)
57  {
58  return 0.0;
59  };
60 
61  return p4.pT2()*(1.0+p4.m2()/p2);
62  };
63 
64  inline double ET(P4 p4)
65  {
66  return(sqrt(ET2(p4)));
67  };
68 
69  inline bool compare_pT(const HEP::Particle* p1, const HEP::Particle* p2)
70  {
71  return (p1->pT() > p2->pT()); // to sort with highest pt first
72  };
73 
74 /*
75  template<typename T1> bool compare_pT(T1* p1, T1* p2)
76  {
77  return (p1->pT() > p2->pT()); // to sort with highest pt first
78  };
79 */
80 
81  inline bool compare_pTJ(const HEP::Jet* p1, const HEP::Jet* p2)
82  {
83  return (p1->pT() > p2->pT()); // to sort with highest pt first
84  };
85 
86 
87  inline bool vertex_match(const HEP::P4 &v1, const HEP::P4 &v2)
88  {
89  HEP::P4 vdiff = v1-v2;
90  if ((v1-v2).p() > 1.0) return false; // vertices closer than a mm are combined
91  if ((v1-v2).T() > 1E-3) return false;
92  return true;
93 
94  }
95 
96 
97  // glorified container for final-state particles + jet constituents, with a method to translate production/decay vertices
98  class PileupEvent {
99  private:
100  std::vector<P4*> _jet_constituents;
101  std::vector<Particle*> _particles;
102  std::vector<P4*> _neutralhadrons;
103  std::vector<Particle*> _chargedhadrons;
104  public:
105 
106  void clear() {
107 
108  #define DELCLEAR(v) do { if (!v.empty()) for (size_t i = 0; i < v.size(); ++i) delete v[i]; v.clear(); } while (0)
109  DELCLEAR(_jet_constituents);
110  DELCLEAR(_particles);
111  DELCLEAR(_neutralhadrons);
112  DELCLEAR(_chargedhadrons);
113  #undef DELCLEAR
114  }
116  PileupEvent() { clear(); }
117 
120  clear();
121  }
122 
123 
124  void add_neutral_hadron(double px, double py, double pz, double e) {
125  P4* thad= new P4(px, py,pz, e);
126  _neutralhadrons.push_back(thad);
127  }
128 
129  void add_neutral_hadron(const P4 &p4)
130  {
131  P4* thad = new P4(p4);
132  _neutralhadrons.push_back(thad);
133 
134  }
135 
136 
137  void add_charged_hadron(double px, double py, double pz, double e, int pid) {
138  Particle* thad= new Particle(px, py, pz, e, pid);
139  _chargedhadrons.push_back(thad);
140  }
141 
142 
143  void add_charged_hadron(Particle* p)
144  {
145 
146  _chargedhadrons.push_back(p);
147 
148  }
149  /*
150  void add_charged_hadron(const Particle &p4)
151  {
152  Particle* thad = new Particle(p4);
153  _chargedhadrons.push_back(thad);
154 
155  }
156  */
157 
158  void add_jet_constituent(double px, double py, double pz, double e) {
159  P4* thad= new P4(px, py,pz, e);
160  _jet_constituents.push_back(thad);
161  }
162 
163  void add_jet_constituent(const P4 &p4)
164  {
165  P4* thad = new P4(p4);
166  _jet_constituents.push_back(thad);
167 
168  }
169 
170  void add_particle(Particle* p) {
171  _particles.push_back(p);
172 
173  }
174 
176  const std::vector<P4*>& get_jet_constituents() const {
177  return _jet_constituents;
178  }
179 
181  const std::vector<Particle*>& get_charged_hadrons() const {
182  return _chargedhadrons;
183  }
184 
186  const std::vector<P4*>& get_neutral_hadrons() const {
187  return _neutralhadrons;
188  }
189 
191  const std::vector<Particle*>& get_particles() const {
192  return _particles;
193  }
194 
195  std::vector<Particle*> translate_particles(double dz, double dt) {
196  std::vector<Particle*> out_particles;
197  for(auto part : _particles)
198  {
199  Particle* new_part = new Particle(part);
200  new_part->translate_zt(dz,dt);
201  out_particles.push_back(new_part);
202  }
203  return out_particles;
204  }
205 
206  std::vector<Particle*> translate_charged_hadrons(double dz, double dt) {
207  std::vector<Particle*> out_particles;
208  for(auto part : _chargedhadrons)
209  {
210  Particle* new_part = new Particle(part);
211  new_part->translate_zt(dz,dt);
212  out_particles.push_back(new_part);
213  }
214  return out_particles;
215  }
216 
217 
218 
219  };
220 
221 
223  // I've had to replace the HEPUtils event because it was insufficient
224 
225 
226 
227  class Event {
228 
229  private:
230 
232 
233 
235  double _centralweight;
236 
238  std::vector<double> _weights;
239 
241 
242  std::vector<Particle*> _photons, _electrons, _muons, _taus, _invisibles, _HSCPs, _chargedhadrons, _neutralRhadrons;
244 
246  mutable std::vector<Jet*> _jets;
247 
249  P4 _pmiss;
250 
251  std::vector<P4*> _neutralhadrons;
252 
253  // For pileup
254  std::vector<P4*> _HCALs, _ECALs; // These should be pointers to momenta owned by the event
255 
256  std::vector<Particle*> _charged_tracks; // these should be pointers to other particles owned by the event
258 
261  void operator = (const Event& e) {
262  clear(); //< Delete current particles
263  _weights = e._weights;
264  _centralweight = e._centralweight;
265  _photons = e._photons;
266  _electrons = e._electrons;
267  _muons = e._muons;
268  _taus = e._taus;
269  _invisibles = e._invisibles;
270  _HSCPs=e._HSCPs;
271  _jets = e._jets;
272  _pmiss = e._pmiss;
273  _neutralhadrons = e._neutralhadrons;
274  _chargedhadrons = e._chargedhadrons;
275  _HCALs=e._HCALs;
276  _ECALs=e._ECALs;
277  _charged_tracks=e._charged_tracks;
278  }
279 
280 
281  public:
282  double Average_rho;
284  void clear() {
285  _weights.clear();
286  _centralweight=0.0;
287  // TODO: indexed loop -> for (Particle* p : particles()) delete p;
288  #define DELCLEAR(v) do { if (!v.empty()) for (size_t i = 0; i < v.size(); ++i) delete v[i]; v.clear(); } while (0)
289  DELCLEAR(_photons);
290  DELCLEAR(_electrons);
291  DELCLEAR(_muons);
292  DELCLEAR(_taus);
293  DELCLEAR(_invisibles);
294  DELCLEAR(_HSCPs);
295  DELCLEAR(_jets);
296  DELCLEAR(_neutralhadrons);
297  DELCLEAR(_chargedhadrons);
298  //DELCLEAR(_ECALs);
299  //DELCLEAR(_HCALs);
300  //DELCLEAR(_charged_tracks);
301  #undef DELCLEAR
302 
303  _pmiss.clear();
304  }
306  Event() { clear(); }
307 
309  Event(const std::vector<Particle*>& ps, const std::vector<double>& weights=std::vector<double>()) {
310  clear();
311  _weights = weights;
312  add_particles(ps);
313  }
314 
316  ~Event() {
317  clear();
318  }
319 
320 
321  // ECAL, HCAL, charged track pointers are not owned by the event: these are pointers to vectors/particles that are
322  // also held elsewhere and will be de-allocated elsewhere.
323 
324  void add_ECAL(P4* p) {
325  _ECALs.push_back(p);
326 
327  }
328 
329  void add_HCAL(P4* p) {
330  _HCALs.push_back(p);
331 
332  }
333 
334  void add_charged_track(HEP::Particle* p) {
335  _charged_tracks.push_back(p);
336 
337  }
338 
339  const std::vector<P4*>& get_ECALs() const {
340  return _ECALs;
341  }
342 
343  const std::vector<P4*>& get_HCALs() const {
344  return _HCALs;
345  }
346 
347  const std::vector<Particle*>& get_charged_tracks() const {
348  return _charged_tracks;
349  }
350 
351 
352  void add_invisible(Particle* p) {
353  _invisibles.push_back(p);
354  }
355 
356  void add_HSCP(Particle* p) {
357  _HSCPs.push_back(p);
358  }
359 
360 
361  void add_neutral_hadron(double px, double py, double pz, double e) {
362  P4* thad= new P4(px, py,pz, e);
363  _neutralhadrons.push_back(thad);
364  }
365 
366  void add_neutral_hadron(const P4 &p4)
367  {
368  P4* thad = new P4(p4);
369  _neutralhadrons.push_back(thad);
370 
371  }
372  void add_neutral_hadrons(const std::vector<P4*> &hadrs)
373  {
374  for (size_t i = 0; i < hadrs.size(); ++i) add_neutral_hadron(*hadrs[i]);
375  }
376 
377 
378 
379  void add_charged_hadron(double px, double py, double pz, double e, int pid) {
380  Particle* thad= new Particle(px, py, pz, e, pid);
381  _chargedhadrons.push_back(thad);
382  }
383 
384  void add_charged_hadron(const P4 &p4, int pid)
385  {
386  Particle* thad = new Particle(p4, pid);
387  _chargedhadrons.push_back(thad);
388 
389  }
390 
391  void add_charged_hadron(Particle* p)
392  {
393 
394  _chargedhadrons.push_back(p);
395 
396  }
397 
398  void add_charged_hadrons(const std::vector<Particle*> &hadrs)
399  {
400  for (size_t i = 0; i < hadrs.size(); ++i) add_charged_hadron(hadrs[i]);
401  }
402  /*
403  void add_charged_hadron(const Particle &p4)
404  {
405  Particle* thad = new Particle(p4);
406  _chargedhadrons.push_back(thad);
407 
408  }
409 
410 
411  void add_charged_hadrons(const std::vector<Particle*> &hadrs)
412  {
413  for (size_t i = 0; i < hadrs.size(); ++i) add_charged_hadron(*hadrs[i]);
414  }
415  */
416 
417  const std::vector<Particle*>& get_charged_hadrons() const {
418  return _chargedhadrons;
419  }
420 
422  const std::vector<P4*>& get_neutral_hadrons() const {
423  return _neutralhadrons;
424  }
425 
427  Event* clone() const {
428  Event* rtn = new Event();
429  cloneTo(rtn);
430  return rtn;
431  }
432 
434  void cloneTo(Event* e) const {
435  assert(e != NULL);
436  cloneTo(*e);
437  }
438 
440  void cloneTo(Event& e) const {
441  e.set_weights(_weights);
442  e.set_weight(_centralweight); // new for central weight
443  const std::vector<Particle*> ps = particles();
444  for (size_t i = 0; i < ps.size(); ++i) {
445  e.add_particle(new Particle(*ps[i]));
446  }
447  const std::vector<Jet*> js = jets();
448  for (size_t i = 0; i < js.size(); ++i) {
449  e.add_jet(new Jet(*js[i]));
450  }
451  e._pmiss = _pmiss;
452  }
453 
455 
456 
458 
459 
461  void set_weights(const std::vector<double>& ws) {
462  _weights = ws;
463  }
464 
465 
466  // don't like this
467  /*
469  void set_weight(double w) {
470  _weights.clear();
471  _weights.push_back(w);
472  }
473  */
474 
475  void set_weight(double w) {
476  _centralweight =w;
477  }
478 
480  const std::vector<double>& weights() const {
481  return _weights;
482  }
483 
485  std::vector<double>& weights() {
486  return _weights;
487  }
488 
489  /*
491  double weight(size_t i=0) const {
492  if (_weights.empty()) {
493  if (i == 0) return 1;
494  throw std::runtime_error("Trying to access non-default weight from empty weight vector");
495  }
496  return _weights[i];
497  }
498 
499  */
501  double weight(size_t i) const {
502  if (_weights.empty()) {
503  throw std::runtime_error("Trying to access non-default weight from empty weight vector");
504  }
505  return _weights[i];
506  }
507  double weight() const {
508  return _centralweight;
509  }
511 
512 
524 
525 
526  if (!p->is_prompt()) // This is misleading since we do want dispaced vertices
527  {
528  if(p->abspid() == 1000024)
529  {
530  _HSCPs.push_back(p);
531  }
532  else if(p->pid() == 1000022)
533  {
534  _invisibles.push_back(p);
535  }
536  else
537  {
538  delete p;
539  }
540  /*
541  if( (p->abspid()) < 1000000 || (p->apspid() > 9000000))
542  {
543  delete p;
544  }
545  else if(p->abspid() == 1000024)
546  {
547  _HSCPs.push_back(p);
548  }
549  else if(p->pid() == 1000022)
550  {
551 
552  }
553  else
554  {
555  delete p;
556  }
557  */
558 
559  }
560  else if (p->pid() == 22)
561  _photons.push_back(p);
562  else if (p->abspid() == 11)
563  _electrons.push_back(p);
564  else if (p->abspid() == 13)
565  _muons.push_back(p);
566  else if (p->abspid() == 15)
567  _taus.push_back(p);
568  //else if (p->abspid() == 12 || p->abspid() == 14 || p->abspid() == 16 ||
569  // p->pid() == 1000022 || p->pid() == 1000039 ||
570  // in_range(p->pid(), 50, 60)) //< invert definition to specify all *visibles*?
571  else if (PDG::isInvisible(p->abspid()))
572  _invisibles.push_back(p);
573  //else if(p->abspid() == 1000024)
574  else if(PDG::isCharged(p->abspid()))
575  _HSCPs.push_back(p);
576  else
577  delete p; //
578  }
579 
580 
584  void add_particles(const std::vector<Particle*>& ps) {
585  for (size_t i = 0; i < ps.size(); ++i) add_particle(ps[i]);
586  }
587 
588 
592  std::vector<Particle*> particles() const {
593  // Add together all the vectors of the different particle types
594  std::vector<Particle*> rtn;
595  // rtn.reserve(visible_particles().size() + _invisibles.size());
596  rtn.reserve(_photons.size() + _electrons.size() + _muons.size() + _taus.size() + _invisibles.size());
597  #define APPEND_VEC(vec) rtn.insert(rtn.end(), vec.begin(), vec.end())
598  // APPEND_VEC(visible_particles());
599  APPEND_VEC(_photons);
600  APPEND_VEC(_electrons);
601  APPEND_VEC(_muons);
602  APPEND_VEC(_taus);
603  APPEND_VEC(_invisibles);
604  #undef APPEND_VEC
605  return rtn;
608  }
609 
610 
614  std::vector<Particle*> visible_particles() const {
615  // Add together all the vectors of the different particle types
616  std::vector<Particle*> rtn;
617  rtn.reserve(_photons.size() + _electrons.size() + _muons.size() + _taus.size());
618  #define APPEND_VEC(vec) rtn.insert(rtn.end(), vec.begin(), vec.end() )
619  APPEND_VEC(_photons);
620  APPEND_VEC(_electrons);
621  APPEND_VEC(_muons);
622  APPEND_VEC(_taus);
623  #undef APPEND_VEC
624  return rtn;
627  }
628 
629 
631  const std::vector<Particle*>& HSCPs() const {
632  return _HSCPs;
633  }
635  std::vector<Particle*>& HSCPs() {
636  return _HSCPs;
637  }
638 
639 
641  const std::vector<Particle*>& invisible_particles() const {
642  return _invisibles;
643  }
645  std::vector<Particle*>& invisible_particles() {
646  return _invisibles;
647  }
648 
649 
651  const std::vector<Particle*>& electrons() const {
652  return _electrons;
653  }
655  std::vector<Particle*>& electrons() {
656  return _electrons;
657  }
658 
659 
661  const std::vector<Particle*>& muons() const {
662  return _muons;
663  }
665  std::vector<Particle*>& muons() {
666  return _muons;
667  }
668 
669 
671  const std::vector<Particle*>& taus() const {
672  return _taus;
673  }
675  std::vector<Particle*>& taus() {
676  return _taus;
677  }
678 
679 
681  const std::vector<Particle*>& photons() const {
682  return _photons;
683  }
685  std::vector<Particle*>& photons() {
686  return _photons;
687  }
688 
689 
692 
693 
695  const std::vector<Jet*>& jets() const {
696  std::sort(_jets.begin(), _jets.end(), _cmpPtDesc);
697  return _jets;
698  }
699 
701  std::vector<Jet*>& jets() {
702  std::sort(_jets.begin(), _jets.end(), _cmpPtDesc);
703  return _jets;
704  }
705 
709  void set_jets(const std::vector<Jet*>& jets) {
710  _jets = jets;
711  }
712 
716  void add_jet(Jet* j) {
717  _jets.push_back(j);
718  }
719 
721 
722 
724 
725 
729  const P4& missingmom() const {
730  return _pmiss;
731  }
732 
736  void set_missingmom(const P4& pmiss) {
737  _pmiss = pmiss;
738  }
739 
741  double met() const {
742  return missingmom().pT();
743  }
744 
746 
747  };
748 
749 
750 }
751 
752 
753 
754 
Simple event class, separating particles into classes.
Definition: heputil.h:227
void add_particles(const std::vector< Particle * > &ps)
Definition: heputil.h:584
double met() const
Get the missing transverse momentum in GeV.
Definition: heputil.h:741
const std::vector< Particle * > & HSCPs() const
Get invisible final state particles.
Definition: heputil.h:631
std::vector< Particle * > visible_particles() const
Definition: heputil.h:614
void set_missingmom(const P4 &pmiss)
Set the missing momentum vector.
Definition: heputil.h:736
Event()
Default constructor.
Definition: heputil.h:306
const std::vector< Particle * > & invisible_particles() const
Get invisible final state particles.
Definition: heputil.h:641
const std::vector< Particle * > & muons() const
Get prompt muons.
Definition: heputil.h:661
std::vector< Particle * > & photons()
Get prompt photons (non-const)
Definition: heputil.h:685
const std::vector< Particle * > & photons() const
Get prompt photons.
Definition: heputil.h:681
std::vector< Particle * > & HSCPs()
Get invisible final state particles (non-const)
Definition: heputil.h:635
const std::vector< Particle * > & electrons() const
Get prompt electrons.
Definition: heputil.h:651
const std::vector< double > & weights() const
Get the event weights (const)
Definition: heputil.h:480
std::vector< Particle * > & electrons()
Get prompt electrons (non-const)
Definition: heputil.h:655
std::vector< double > & weights()
Get the event weights (non-const)
Definition: heputil.h:485
~Event()
Destructor (cleans up all passed Particles and calculated Jets)
Definition: heputil.h:316
const std::vector< Particle * > & taus() const
Get prompt (hadronic) taus.
Definition: heputil.h:671
void add_particle(Particle *p)
Definition: heputil.h:523
std::vector< Particle * > & invisible_particles()
Get invisible final state particles (non-const)
Definition: heputil.h:645
void cloneTo(Event &e) const
Clone a deep copy (new Particles and Jets allocated) into the provided event object.
Definition: heputil.h:440
Event(const std::vector< Particle * > &ps, const std::vector< double > &weights=std::vector< double >())
Constructor from a list of Particles, I don't use it.
Definition: heputil.h:309
const std::vector< Jet * > & jets() const
Get anti-kT 0.4 jets (not including charged leptons or photons)
Definition: heputil.h:695
double weight(size_t i) const
Get a single event weight – the nominal, by default.
Definition: heputil.h:501
void set_weights(const std::vector< double > &ws)
Set the event weights (also possible directly via non-const reference)
Definition: heputil.h:461
std::vector< Particle * > particles() const
Definition: heputil.h:592
void add_jet(Jet *j)
Add a jet to the jets collection.
Definition: heputil.h:716
void clear()
Empty the event's particle, jet and MET collections.
Definition: heputil.h:284
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
void cloneTo(Event *e) const
Clone a deep copy (new Particles and Jets allocated) into the provided event pointer.
Definition: heputil.h:434
Event * clone() const
Clone a copy on the heap.
Definition: heputil.h:427
std::vector< Particle * > & taus()
Get prompt (hadronic) taus (non-const)
Definition: heputil.h:675
std::vector< Particle * > & muons()
Get prompt muons (non-const)
Definition: heputil.h:665
const P4 & missingmom() const
Get the missing momentum vector.
Definition: heputil.h:729
std::vector< Jet * > & jets()
Get anti-kT 0.4 jets (not including charged leptons or photons) (non-const)
Definition: heputil.h:701
Definition: Jet.h:21
double pT() const
Get the squared transverse momentum.
Definition: Jet.h:96
A 4-momentum class for vectors.
Definition: Vectors.h:45
double p2() const
Get E^2.
Definition: Vectors.h:441
double m2() const
Get m^2.
Definition: Vectors.h:408
double pT() const
Get the transverse momentum (same as rho)
Definition: Vectors.h:453
double pT2() const
Get the transverse momentum squared (same as rho2)
Definition: Vectors.h:451
void clear()
Set the components to zero.
Definition: Vectors.h:82
Definition: Particle.h:24
int pid() const
Get PDG ID code.
Definition: Particle.h:187
int abspid() const
Get abs PDG ID code.
Definition: Particle.h:189
double pT() const
Get the squared transverse momentum.
Definition: Particle.h:165
bool is_prompt() const
Is this particle connected to the hard process or from a hadron/tau decay?
Definition: Particle.h:176
Definition: heputil.h:98
~PileupEvent()
Destructor (cleans up all passed Particles and calculated Jets)
Definition: heputil.h:119
const std::vector< Particle * > & get_charged_hadrons() const
Get charged hadrons.
Definition: heputil.h:181
PileupEvent()
Default constructor.
Definition: heputil.h:116
const std::vector< P4 * > & get_jet_constituents() const
Get jet constituents.
Definition: heputil.h:176
const std::vector< Particle * > & get_particles() const
Get particles.
Definition: heputil.h:191
const std::vector< P4 * > & get_neutral_hadrons() const
Get neutral hadrons.
Definition: heputil.h:186
Modified by Mark Goodsell goodsell@lpthe.jussieu.fr
Definition: ATLAS_SUSY_2018_16.cc:27
double ET2(P4 p4)
Definition: heputil.h:52
bool _cmpPtDesc(const Jet *a, const Jet *b)
Function/functor for container<const Jet*> sorting (cf. std::less)
Definition: Jet.h:123