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