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