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