CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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> inline bool GreaterByE (const T& a1, const T& a2) {return a1.E() > a2.E();}
10 
11  //sort by descending pt
12  inline bool sortByPt(const HepMC::GenParticle* a , const HepMC::GenParticle* b) {return a->momentum().perp() > b->momentum().perp(); }
13 
14  //sort by energy
15  inline bool sortByE(const HepMC::GenParticle* a , const HepMC::GenParticle* b) {return a->momentum().e() > b->momentum().e(); }
16 
17  //sort by rapidity
18  inline bool sortByRapidity(const HepMC::GenParticle* a , const HepMC::GenParticle* b) {
19  const HepMC::FourVector& amom = a->momentum();
20  const HepMC::FourVector& bmom = b->momentum();
21  double rapa = 0.5 * std::log( (amom.e() + amom.z()) / (amom.e() - amom.z()) );
22  double rapb = 0.5 * std::log( (bmom.e() + bmom.z()) / (bmom.e() - bmom.z()) );
23  return rapa > rapb;
24  }
25 
26  //sort by pseudorapidity
27  inline bool sortByPseudoRapidity(const HepMC::GenParticle* a , const HepMC::GenParticle* b) {return a->momentum().eta() > b->momentum().eta(); }
28 
29  void findDescendents(const HepMC::GenParticle* a, std::vector<const HepMC::GenParticle*>& descendents);
30 
31  //get all status 1 particles
32  void allStatus1(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& status1);
33 
34  //get all status 2 particles
35  void allStatus2(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& status2);
36 
37  //get all status 3 particles
38  void allStatus3(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& status3);
39 
40  //given a collection of leptons and the collection of all particles in the event,
41  //get the collection of photons closer than deltaR from any of the leptons
42  void findFSRPhotons(const std::vector<const HepMC::GenParticle*>& leptons,
43  const std::vector<const HepMC::GenParticle*>& all,
44  double deltaR,
45  std::vector<const HepMC::GenParticle*>& photons);
46 
47 
48  //given a collection of leptons and the collection of all particles in the event,
49  //get the collection of photons closer than deltaR from any of the leptons
50  void findFSRPhotons(const std::vector<const HepMC::GenParticle*>& leptons,
51  const HepMC::GenEvent* all,
52  double deltaR,
53  std::vector<const HepMC::GenParticle*>& photons);
54 
55  //returns true if a status 3 particle is a tau or if a status 1 particle is either an electron or a neutrino
57 
58  //returns true if a status 1 particle is a neutrino
59  bool isNeutrino(const HepMC::GenParticle* part );
60 
61  //returns true is status 3 particle is tau
62  bool isTau(const HepMC::GenParticle* part);
63 
64  //removes isolated leptons (their decay products in the case of taus) and possible fsr photons from the list of particles.
65  //this should result in a list of particles to be used for hadronic activity studies, such as construction of jets
66  //note: deltaR is both the cone in wich we compute "isolation" and the radius in which we look for fsr photons
67  void removeIsolatedLeptons(const HepMC::GenEvent* all, double deltaR, double sumPt, std::vector<const HepMC::GenParticle*>& pruned);
68 
69  //get all status 1 particles
70  void allStatus1(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& status1);
71 
72  //get all visible status1 particles
73  void allVisibleParticles(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& visible);
74 
75  //compute generated met
76  TLorentzVector genMet(const HepMC::GenEvent* all, double etamin = -9999., double etamax = 9999.);
77 
78 }
79 
80 #endif
bool isTau(const HepMC::GenParticle *part)
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)
void findDescendents(const HepMC::GenParticle *a, std::vector< const HepMC::GenParticle * > &descendents)
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
part
Definition: HCALResponse.h:20
TLorentzVector genMet(const HepMC::GenEvent *all, double etamin=-9999., double etamax=9999.)
double b
Definition: hdecay.h:120
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:121
long double T
tuple log
Definition: cmsBatch.py:341