FairRoot/PandaRoot
PndEventShape.h
Go to the documentation of this file.
1 #ifndef PNDEVENTSHAPE_H
2 #define PNDEVENTSHAPE_H 1
3 
4 #include "TLorentzVector.h"
5 #include "TVector3.h"
6 
7 #define FWMAX 6 // maximum Fox Wolfram moment
8 
9 class RhoCandList;
10 
12 public:
13  PndEventShape(RhoCandList &l, TLorentzVector cms, double neutMinE=0.0, double chrgMinP=0.0);
15 
16  // ******* multiplicities
17  int NParticles() const {return fN;} // number of particle candidates
18  int NCharged() const {return fnChrg;} // number of charged candidates
19  int NNeutral() const {return fnNeut;} // number of neutral candidates
20 
21  // ******* maxima of momenta
22  double PmaxLab() const {return fpmaxlab;} // max momentum in lab system
23  double PmaxCms() const {return fpmaxcms;} // max momentum im cms system
24  double PminLab() const {return fpminlab;} // min momentum in lab system
25  double PminCms() const {return fpmincms;} // min momentum im cms system
26  double Ptmax() const {return fptmax;} // max pt (same for lab and cms)
27  double Ptmin() const {return fptmin;} // min pt (same for lab and cms)
28  double PRapmax() const {return fprapmax;} // max pseudorapidity (lab)
29  double EmaxNeutLab() const {return femaxneutlab;} // max neutral energy (lab)
30  double EmaxNeutCms() const {return femaxneutcms;} // max neutral energy (cms)
31  double PmaxChrgLab() const {return fpmaxchlab;} // max charged momentum (lab)
32  double PmaxChrgCms() const {return fpmaxchcms;} // max charged momentum (cms)
33 
34  // ******* sum of energies/momenta (lab system)
35  double PtSumLab() const {return fptsumlab;} // sum of pt in (lab)
36  double NeutEtSumLab() const {return fneutetsumlab;} // sum of transvers energys of neutrals (lab)
37  double NeutESumLab() const {return fneutesumlab;} // sum of energys of neutrals (lab)
38  double ChrgPtSumLab() const {return fchrgptsumlab;} // sum of pt of charged (lab)
39  double ChrgPSumLab() const {return fchrgpsumlab;} // sum of momenta of charged (lab)
40 
41  // ******* sum of energies/momenta (cms system)
42  double PtSumCms() const {return fptsumcms;} // sum of pt in (cms)
43  double NeutEtSumCms() const {return fneutetsumcms;} // sum of transvers energys of neutrals (cms)
44  double NeutESumCms() const {return fneutesumcms;} // sum of energys of neutrals (cms)
45  double ChrgPtSumCms() const {return fchrgptsumcms;} // sum of pt of charged (cms)
46  double ChrgPSumCms() const {return fchrgpsumcms;} // sum of momenta of charged (cms)
47 
48  // ******* detector sepcific quantities
49  double DetEmcSum() const {return fdetemcsum;} // sum of EMC cluster energies
50  double DetEmcMax() const {return fdetemcmax;} // maximum of EMC cluster energies
51 
52  // ******* multiplicities with threshold
53  int MultPminLab(double pmin); // number of particles with p>pmin (lab frame)
54  int MultPmaxLab(double pmax); // number of particles with p<pmax (lab frame)
55  int MultPminCms(double pmin); // number of particles with p>pmin (cms frame)
56  int MultPmaxCms(double pmax); // number of particles with p<pmax (cms frame)
57 
58  // ******* multiplicities with threshold
59  int MultPtminLab(double ptmin); // number of particles with pt>pmin (lab frame)
60  int MultPtmaxLab(double ptmax); // number of particles with pt<pmax (lab frame)
61  int MultPtminCms(double ptmin); // number of particles with pt>pmin (cms frame)
62  int MultPtmaxCms(double ptmax); // number of particles with pt<pmax (cms frame)
63 
64  // ******* PID multiplicities with PID and momentum threshold
65  int MultElectronPminLab(double prob, double pmin=0); // number of electrons with p>pmin and PID prob>prob (lab frame)
66  int MultMuonPminLab(double prob, double pmin=0); // number of muons with p>pmin and PID prob>prob (lab frame)
67  int MultPionPminLab(double prob, double pmin=0); // number of pions with p>pmin and PID prob>prob (lab frame)
68  int MultKaonPminLab(double prob, double pmin=0); // number of kaons with p>pmin and PID prob>prob (lab frame)
69  int MultProtonPminLab(double prob, double pmin=0); // number of protons with p>pmin and PID prob>prob (lab frame)
70 
71  int MultElectronPminCms(double prob, double pmin=0); // number of electrons with p>pmin and PID prob>prob (cms frame)
72  int MultMuonPminCms(double prob, double pmin=0); // number of muons with p>pmin and PID prob>prob (cms frame)
73  int MultPionPminCms(double prob, double pmin=0); // number of pions with p>pmin and PID prob>prob (cms frame)
74  int MultKaonPminCms(double prob, double pmin=0); // number of kaons with p>pmin and PID prob>prob (cms frame)
75  int MultProtonPminCms(double prob, double pmin=0); // number of protons with p>pmin and PID prob>prob (cms frame)
76 
77  // ******* neutrals multiplicities with threshold
78  int MultNeutEminLab(double emin); // number of neutrals with E>emin (lab frame)
79  int MultNeutEmaxLab(double emax); // number of neutrals with E<emax (lab frame)
80  int MultNeutEminCms(double emin); // number of neutrals with E>emin (cms frame)
81  int MultNeutEmaxCms(double emax); // number of neutrals with E<emax (cms frame)
82 
83  // ******* charged multiplicities with threshold
84  int MultChrgPminLab(double pmin); // number of charged with p>pmin (lab frame)
85  int MultChrgPmaxLab(double pmax); // number of charged with p<pmax (lab frame)
86  int MultChrgPminCms(double pmin); // number of charged with p>pmin (cms frame)
87  int MultChrgPmaxCms(double pmax); // number of charged with p<pmax (cms frame)
88 
89  // ******* sums with threshold
90  double SumPminLab(double pmin); // sum of momenta with p>pmin (lab frame)
91  double SumPmaxLab(double pmax); // sum of momenta with p<pmax (lab frame)
92  double SumPminCms(double pmin); // sum of momenta with p>pmin (cms frame)
93  double SumPmaxCms(double pmax); // sum of momenta with p<pmax (cms frame)
94 
95  // ******* sums with threshold
96  double SumPtminLab(double ptmin); // sum of pt with pt>pmin (lab frame)
97  double SumPtmaxLab(double ptmax); // sum of pt with pt<pmax (lab frame)
98  double SumPtminCms(double ptmin); // sum of pt with pt>pmin (cms frame)
99  double SumPtmaxCms(double ptmax); // sum of pt with pt<pmax (cms frame)
100 
101  // ******* neutrals sums with threshold
102  double SumNeutEminLab(double emin); // sum of energies of neutrals with E>emin (lab frame)
103  double SumNeutEmaxLab(double emax); // sum of energies of neutrals with E<emax (lab frame)
104  double SumNeutEminCms(double emin); // sum of energies of neutrals with E>emin (cms frame)
105  double SumNeutEmaxCms(double emax); // sum of energies of neutrals with E<emax (cms frame)
106 
107  // ******* charged sums with threshold
108  double SumChrgPminLab(double pmin); // sum of momenta of charged with p>pmin (lab frame)
109  double SumChrgPmaxLab(double pmax); // sum of momenta of charged with p<pmax (lab frame)
110  double SumChrgPminCms(double pmin); // sum of momenta of charged with p>pmin (cms frame)
111  double SumChrgPmaxCms(double pmax); // sum of momenta of charged with p<pmax (cms frame)
112 
113  // ******* shape variables
114  double Sphericity(); // Sphericity
115  double Aplanarity(); // Aplanarity
116  double Planarity(); // Planarity
117  double Circularity(); // Cirularity
118 
119  double FoxWolfMomH(int order); // Fox Wolfram moment absolute H_i
120  double FoxWolfMomR(int order); // Fox Wolfram moment relative R_i = H_i/H_0
121 
122  double Thrust(int Nmax=4); // Thrust, with 0.5 < thr < 1.0
123  TVector3 ThrustVector(); // Direction of thrust vector
124 
125 private:
126 
127  void ComputeSphericity(); // compute sph, apl, pla
128  double Eps(const TVector3 v1, const TVector3 v2) {return (v1*v2)>0. ? 1. : -1.;} // aux for Thrust
129  double Legendre( int l, double x ); // Legendre function; auxilliary for Fox Wolfram moments
130  static bool CmpTVect3Mag(TVector3 v1, TVector3 v2) { return (v1.Mag()<v2.Mag()); }
131 
132  std::vector<TLorentzVector> fLabList; // List of 4-vectors in lab frame
133  std::vector<TLorentzVector> fCmsList; // List of 4-vectors in cms frame
134  std::vector<int> fCharge; // List of charges of particles
135  std::vector<double> fElProb; // List of electron probabilities
136  std::vector<double> fMuProb; // List of muon probabilities
137  std::vector<double> fPiProb; // List of pion probabilities
138  std::vector<double> fKaProb; // List of kaon probabilities
139  std::vector<double> fPrProb; // List of proton probabilities
140 
141  int fnChrg; // number of charged particles
142  int fnNeut; // number of neutral particles
143  int fN; // number of particles
144 
145  double fpmaxlab; // maximum momentum lab frame
146  double fpmaxcms; // maximum momentum cms frame
147  double fpminlab; // minimum momentum lab frame
148  double fpmincms; // minimum momentum cms frame
149  double fptmax; // maximum transvers momentum
150  double fptmin; // minimum transvers momentum
151  double fprapmax; // maximum pseudorapidity
152  double femaxneutlab; // max neutral energy (lab)
153  double femaxneutcms; // max neutral energy (cms)
154  double fpmaxchlab; // max charged momentum (lab)
155  double fpmaxchcms; // max charged momentum (cms)
156 
157  double fdetemcsum; // sum of EMC cluster energies
158  double fdetemcmax; // maximum of EMC cluster energies
159 
160  double fptsumlab; // sum of pt in (lab)
161  double fneutetsumlab; // sum of transvers energys of neutrals (lab)
162  double fneutesumlab; // sum of energys of neutrals (lab)
163  double fchrgptsumlab; // sum of pt of charged (lab)
164  double fchrgpsumlab; // sum of momenta of charged (lab)
165 
166  double fptsumcms; // sum of pt in (lab)
167  double fneutetsumcms; // sum of transvers energys of neutrals (lab)
168  double fneutesumcms; // sum of energys of neutrals (lab)
169  double fchrgptsumcms; // sum of pt of charged (lab)
170  double fchrgpsumcms; // sum of momenta of charged (lab)
171 
172  double fsph; // Sphericity
173  double fapl; // Aplanarity
174  double fpla; // Planarity
175  double fcir; // Circularity
176 
177  double fFWmom[FWMAX]; // Fox Wolfram moments up to FWMAX
178  bool fFWready; // did we compute FW moments?
179 
180  double fthr; // thrust
181  TVector3 fThrVect; // direction of thrust
182  TVector3 fBoost; // boost vector to go to requested frame
183 
184 
185 };
186 
187 #endif
double SumNeutEmaxCms(double emax)
int MultChrgPminCms(double pmin)
double SumNeutEminLab(double emin)
double femaxneutcms
double FoxWolfMomH(int order)
double ChrgPtSumCms() const
Definition: PndEventShape.h:45
double ChrgPSumCms() const
Definition: PndEventShape.h:46
int MultPmaxLab(double pmax)
int MultKaonPminCms(double prob, double pmin=0)
int MultProtonPminCms(double prob, double pmin=0)
double Sphericity()
double fneutetsumlab
std::vector< double > fKaProb
int MultPtmaxLab(double ptmax)
double ChrgPSumLab() const
Definition: PndEventShape.h:39
double NeutEtSumCms() const
Definition: PndEventShape.h:43
double DetEmcSum() const
Definition: PndEventShape.h:49
TVector3 ThrustVector()
std::vector< int > fCharge
double SumChrgPminCms(double pmin)
double Circularity()
double SumPminCms(double pmin)
double PtSumLab() const
Definition: PndEventShape.h:35
int MultPtminCms(double ptmin)
std::vector< TLorentzVector > fLabList
double fchrgptsumcms
double SumNeutEmaxLab(double emax)
int NParticles() const
Definition: PndEventShape.h:17
PndEventShape(RhoCandList &l, TLorentzVector cms, double neutMinE=0.0, double chrgMinP=0.0)
double fchrgpsumlab
double SumNeutEminCms(double emin)
static bool CmpTVect3Mag(TVector3 v1, TVector3 v2)
int MultKaonPminLab(double prob, double pmin=0)
double SumPtmaxCms(double ptmax)
double EmaxNeutLab() const
Definition: PndEventShape.h:29
int MultPtminLab(double ptmin)
double SumChrgPminLab(double pmin)
std::vector< double > fElProb
double EmaxNeutCms() const
Definition: PndEventShape.h:30
double femaxneutlab
int MultNeutEminCms(double emin)
int MultChrgPmaxCms(double pmax)
double DetEmcMax() const
Definition: PndEventShape.h:50
std::vector< double > fMuProb
int MultPtmaxCms(double ptmax)
double Planarity()
int MultNeutEmaxCms(double emax)
double PRapmax() const
Definition: PndEventShape.h:28
double Legendre(int l, double x)
double PmaxChrgCms() const
Definition: PndEventShape.h:32
int MultMuonPminLab(double prob, double pmin=0)
int MultPminLab(double pmin)
double fneutesumcms
int MultChrgPmaxLab(double pmax)
double NeutESumLab() const
Definition: PndEventShape.h:37
int MultNeutEmaxLab(double emax)
std::vector< double > fPiProb
int MultPionPminCms(double prob, double pmin=0)
double SumPtmaxLab(double ptmax)
double fneutetsumcms
double ChrgPtSumLab() const
Definition: PndEventShape.h:38
TVector3 fThrVect
double Eps(const TVector3 v1, const TVector3 v2)
double fFWmom[FWMAX]
void ComputeSphericity()
double SumPtminCms(double ptmin)
double PmaxLab() const
Definition: PndEventShape.h:22
double PminCms() const
Definition: PndEventShape.h:25
int MultPminCms(double pmin)
double PmaxChrgLab() const
Definition: PndEventShape.h:31
double fneutesumlab
double SumChrgPmaxLab(double pmax)
double Thrust(int Nmax=4)
double PtSumCms() const
Definition: PndEventShape.h:42
int MultPionPminLab(double prob, double pmin=0)
Double_t x
int MultMuonPminCms(double prob, double pmin=0)
int MultPmaxCms(double pmax)
double Ptmin() const
Definition: PndEventShape.h:27
double SumChrgPmaxCms(double pmax)
TVector3 v1
Definition: bump_analys.C:40
double SumPtminLab(double ptmin)
TVector3 fBoost
int MultNeutEminLab(double emin)
double SumPmaxCms(double pmax)
double SumPminLab(double pmin)
std::vector< double > fPrProb
int MultChrgPminLab(double pmin)
std::vector< TLorentzVector > fCmsList
double FoxWolfMomR(int order)
TVector3 v2
Definition: bump_analys.C:40
double NeutESumCms() const
Definition: PndEventShape.h:44
int MultElectronPminCms(double prob, double pmin=0)
int MultElectronPminLab(double prob, double pmin=0)
#define FWMAX
Definition: PndEventShape.h:7
int NNeutral() const
Definition: PndEventShape.h:19
int NCharged() const
Definition: PndEventShape.h:18
double PminLab() const
Definition: PndEventShape.h:24
int MultProtonPminLab(double prob, double pmin=0)
double fchrgptsumlab
double PmaxCms() const
Definition: PndEventShape.h:23
double Ptmax() const
Definition: PndEventShape.h:26
double SumPmaxLab(double pmax)
double fchrgpsumcms
double NeutEtSumLab() const
Definition: PndEventShape.h:36
double Aplanarity()