FairRoot/PandaRoot
RhoEventShapes.h
Go to the documentation of this file.
1 #ifndef RHOEVENTSHAPES_H
2 #define RHOEVENTSHAPES_H
3 // //
5 // RhoEventShapes //
6 // //
7 // Functions to calculate the event & event shape variables //
8 // //
9 // Author List: //
10 // Klaus Goetzen, GSI, May 2013 //
11 // Ralf Kliemt, HIM/GSI May.2013 //
12 // //
14 
15 
16 #include "TLorentzVector.h"
17 #include "TVector3.h"
18 #include "TMatrixD.h"
19 #include "TMatrixDEigen.h"
20 #include "RhoCandList.h"
21 
22 #define FWMAX 6 // maximum Fox Wolfram moment
23 
25 public:
26  RhoEventShapes(RhoCandList &l, TLorentzVector cms);
27  virtual ~RhoEventShapes(){};
28 
29  // ******* multiplicities
30  int NParticles() const {return fN;} // number of particle candidates
31  int NCharged() const {return fnChrg;} // number of charged candidates
32  int NNeutral() const {return fnNeut;} // number of neutral candidates
33 
34  // ******* maxima of momenta
35  double PmaxLab() const {return fpmaxlab;} // max momentum in lab system
36  double PmaxCms() const {return fpmaxcms;} // max momentum im cms system
37  double Ptmax() const {return fptmax;} // max pt (same for lab and cms)
38 
39  // ******* sum of energies/momenta (lab system)
40  double PtSumLab() const {return fptsumlab;} // sum of pt in (lab)
41  double NeutEtSumLab() const {return fneutetsumlab;} // sum of transvers energys of neutrals (lab)
42  double NeutESumLab() const {return fneutesumlab;} // sum of energys of neutrals (lab)
43  double ChrgPtSumLab() const {return fchrgptsumlab;} // sum of pt of charged (lab)
44  double ChrgPSumLab() const {return fchrgpsumlab;} // sum of momenta of charged (lab)
45 
46  // ******* sum of energies/momenta (cms system)
47  double PtSumCms() const {return fptsumcms;} // sum of pt in (lab)
48  double NeutEtSumCms() const {return fneutetsumcms;} // sum of transvers energys of neutrals (lab)
49  double NeutESumCms() const {return fneutesumcms;} // sum of energys of neutrals (lab)
50  double ChrgPtSumCms() const {return fchrgptsumcms;} // sum of pt of charged (lab)
51  double ChrgPSumCms() const {return fchrgpsumcms;} // sum of momenta of charged (lab)
52 
53  // ******* multiplicities with threshold
54  int MultPminLab(double pmin); // number of particles with p>pmin (lab frame)
55  int MultPmaxLab(double pmax); // number of particles with p<pmax (lab frame)
56  int MultPminCms(double pmin); // number of particles with p>pmin (cms frame)
57  int MultPmaxCms(double pmax); // number of particles with p<pmax (cms frame)
58 
59  // ******* shape variables
60  double Sphericity(); // Sphericity
61  double Aplanarity(); // Aplanarity
62  double Planarity(); // Planarity
63 
64  double FoxWolfMomH(int order); // Fox Wolfram moment absolute H_i
65  double FoxWolfMomR(int order); // Fox Wolfram moment relative R_i = H_i/H_0
66 
67  double Thrust();
68  TVector3 ThrustVector();
69 
70 private:
71 
72  void ComputeSphericity(); // compute sph, apl, pla
73  double Eps(const TVector3 v1, const TVector3 v2) {return (v1*v2)>0. ? 1. : -1.;} // aux for Thrust
74  double Legendre( int l, double x ); // Legendre function; auxilliary for Fox Wolfram moments
75 
76  std::vector<TLorentzVector> fLabList;
77  std::vector<TLorentzVector> fCmsList;
78  std::vector<int> fCharge;
79 
80  int fnChrg; // number of charged particles
81  int fnNeut; // number of neutral particles
82  int fN; // number of particles
83 
84  double fpmaxlab; // maximum momentum lab frame
85  double fpmaxcms; // maximum momentum cms frame
86  double fptmax; // maximum transvers momentum
87 
88  double fptsumlab; // sum of pt in (lab)
89  double fneutetsumlab; // sum of transvers energys of neutrals (lab)
90  double fneutesumlab; // sum of energys of neutrals (lab)
91  double fchrgptsumlab; // sum of pt of charged (lab)
92  double fchrgpsumlab; // sum of momenta of charged (lab)
93 
94  double fptsumcms; // sum of pt in (lab)
95  double fneutetsumcms; // sum of transvers energys of neutrals (lab)
96  double fneutesumcms; // sum of energys of neutrals (lab)
97  double fchrgptsumcms; // sum of pt of charged (lab)
98  double fchrgpsumcms; // sum of momenta of charged (lab)
99 
100  double fsph; // Sphericity
101  double fapl; // Aplanarity
102  double fpla; // Planarity
103 
104  double fFWmom[FWMAX]; // Fox Wolfram moments up to FWMAX
105  bool fFWready; // did we compute FW moments?
106 
107  double fthr; // thrust
108  TVector3 fThrVect; // direction of thrust
109  TVector3 fBoost; // boost vector to go to requested frame
110 
111  bool fNeutralOnly; // only compute for neutrals
112  bool fChargedOnly; // only compute for charged
113 
114  public:
116 
117 };
118 #endif
double fchrgptsumlab
#define FWMAX
double PtSumCms() const
double PtSumLab() const
int NNeutral() const
double fchrgptsumcms
double Ptmax() const
int MultPminLab(double pmin)
double NeutEtSumLab() const
double ChrgPtSumCms() const
double NeutESumCms() const
std::vector< TLorentzVector > fCmsList
! List of 4-vectors in cms frame
RhoEventShapes(RhoCandList &l, TLorentzVector cms)
ClassDef(RhoEventShapes, 1)
double FoxWolfMomH(int order)
double ChrgPtSumLab() const
double fneutetsumcms
std::vector< int > fCharge
! List of charges of particles
double NeutEtSumCms() const
double NeutESumLab() const
void ComputeSphericity()
double PmaxLab() const
int MultPmaxCms(double pmax)
double fneutetsumlab
double Legendre(int l, double x)
Double_t x
int MultPmaxLab(double pmax)
TVector3 ThrustVector()
double Eps(const TVector3 v1, const TVector3 v2)
TVector3 v1
Definition: bump_analys.C:40
double ChrgPSumCms() const
TVector3 v2
Definition: bump_analys.C:40
virtual ~RhoEventShapes()
double PmaxCms() const
int NParticles() const
double fFWmom[FWMAX]
double ChrgPSumLab() const
std::vector< TLorentzVector > fLabList
! List of 4-vectors in lab frame
int MultPminCms(double pmin)
int NCharged() const
double FoxWolfMomR(int order)