CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HepMCValidationHelper.h
Go to the documentation of this file.
1 #ifndef Validation_EventGenerator_HepMCValidationHelper
2 #define Validation_EventGenerator_HepMCValidationHelper
3 
4 #include <HepMC/GenEvent.h>
5 #include <vector>
6 #include "TLorentzVector.h"
7 
8 namespace HepMCValidationHelper {
9  template <class T>
10  inline bool GreaterByE(const T& a1, const T& a2) {
11  return a1.E() > a2.E();
12  }
13 
14  //sort by descending pt
15  inline bool sortByPt(const HepMC::GenParticle* a, const HepMC::GenParticle* b) {
16  return a->momentum().perp() > b->momentum().perp();
17  }
18 
19  //sort by energy
20  inline bool sortByE(const HepMC::GenParticle* a, const HepMC::GenParticle* b) {
21  return a->momentum().e() > b->momentum().e();
22  }
23 
24  //sort by rapidity
25  inline bool sortByRapidity(const HepMC::GenParticle* a, const HepMC::GenParticle* b) {
26  const HepMC::FourVector& amom = a->momentum();
27  const HepMC::FourVector& bmom = b->momentum();
28  double rapa = 0.5 * std::log((amom.e() + amom.z()) / (amom.e() - amom.z()));
29  double rapb = 0.5 * std::log((bmom.e() + bmom.z()) / (bmom.e() - bmom.z()));
30  return rapa > rapb;
31  }
32 
33  //sort by pseudorapidity
35  return a->momentum().eta() > b->momentum().eta();
36  }
37 
38  void findDescendents(const HepMC::GenParticle* a, std::vector<const HepMC::GenParticle*>& descendents);
39 
40  //get all status 1 particles
41  void allStatus1(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& status1);
42 
43  //get all status 2 particles
44  void allStatus2(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& status2);
45 
46  //get all status 3 particles
47  void allStatus3(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& status3);
48 
49  //given a collection of leptons and the collection of all particles in the event,
50  //get the collection of photons closer than deltaR from any of the leptons
51  void findFSRPhotons(const std::vector<const HepMC::GenParticle*>& leptons,
52  const std::vector<const HepMC::GenParticle*>& all,
53  double deltaR,
54  std::vector<const HepMC::GenParticle*>& photons);
55 
56  //given a collection of leptons and the collection of all particles in the event,
57  //get the collection of photons closer than deltaR from any of the leptons
58  void findFSRPhotons(const std::vector<const HepMC::GenParticle*>& leptons,
59  const HepMC::GenEvent* all,
60  double deltaR,
61  std::vector<const HepMC::GenParticle*>& photons);
62 
63  //returns true if a status 3 particle is a tau or if a status 1 particle is either an electron or a neutrino
65 
66  //returns true if a status 1 particle is a neutrino
67  bool isNeutrino(const HepMC::GenParticle* part);
68 
69  //returns true is status 3 particle is tau
70  bool isTau(const HepMC::GenParticle* part);
71 
72  //removes isolated leptons (their decay products in the case of taus) and possible fsr photons from the list of particles.
73  //this should result in a list of particles to be used for hadronic activity studies, such as construction of jets
74  //note: deltaR is both the cone in wich we compute "isolation" and the radius in which we look for fsr photons
76  double deltaR,
77  double sumPt,
78  std::vector<const HepMC::GenParticle*>& pruned);
79 
80  //get all status 1 particles
81  void allStatus1(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& status1);
82 
83  //get all visible status1 particles
84  void allVisibleParticles(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& visible);
85 
86  //compute generated met
87  TLorentzVector genMet(const HepMC::GenEvent* all, double etamin = -9999., double etamax = 9999.);
88 
89 } // namespace HepMCValidationHelper
90 
91 #endif
bool isTau(const HepMC::GenParticle *part)
static std::vector< std::string > checklist log
void allStatus3(const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &status3)
void findFSRPhotons(const std::vector< const HepMC::GenParticle * > &leptons, const std::vector< const HepMC::GenParticle * > &all, double deltaR, std::vector< const HepMC::GenParticle * > &photons)
bool sortByPt(const HepMC::GenParticle *a, const HepMC::GenParticle *b)
bool isChargedLepton(const HepMC::GenParticle *part)
void allStatus1(const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &status1)
bool sortByPseudoRapidity(const HepMC::GenParticle *a, const HepMC::GenParticle *b)
void removeIsolatedLeptons(const HepMC::GenEvent *all, double deltaR, double sumPt, std::vector< const HepMC::GenParticle * > &pruned)
void allStatus2(const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &status2)
bool sortByRapidity(const HepMC::GenParticle *a, const HepMC::GenParticle *b)
void allVisibleParticles(const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &visible)
def all
workaround iterator generators for ROOT classes
Definition: cmstools.py:25
void findDescendents(const HepMC::GenParticle *a, std::vector< const HepMC::GenParticle * > &descendents)
part
Definition: HCALResponse.h:20
TLorentzVector genMet(const HepMC::GenEvent *all, double etamin=-9999., double etamax=9999.)
double b
Definition: hdecay.h:118
bool isNeutrino(const HepMC::GenParticle *part)
bool sortByE(const HepMC::GenParticle *a, const HepMC::GenParticle *b)
bool GreaterByE(const T &a1, const T &a2)
double a
Definition: hdecay.h:119
long double T