PandaRoot
KFParticleBase.h
Go to the documentation of this file.
1 //---------------------------------------------------------------------------------
2 // The KFParticleBase class
3 // .
4 // @author S.Gorbunov, I.Kisel, I.Kulakov, M.Zyzak
5 // @version 1.0
6 // @since 13.05.07
7 //
8 // Class to reconstruct and store the decayed particle parameters.
9 // The method is described in CBM-SOFT note 2007-003,
10 // ``Reconstruction of decayed particles based on the Kalman filter'',
11 // http://www.gsi.de/documents/DOC-2007-May-14-1.pdf
12 //
13 // This class describes general mathematics which is used by KFParticle class
14 //
15 // -= Copyright &copy ALICE HLT and CBM L1 Groups =-
16 //_________________________________________________________________________________
17 
18 #ifndef ALIKFPARTICLEBASE_H
19 #define ALIKFPARTICLEBASE_H
20 
21 //#include "TObject.h"
22 #include "RootTypesDef.h"
23 
24 #include <vector>
25 
26 class KFParticleBase : public TObject {
27 
28  public:
29  //*
30  //* ABSTRACT METHODS HAVE TO BE DEFINED IN USER CLASS
31  //*
32 
33  //* Virtual method to access the magnetic field
34 
35  virtual void GetFieldValue(const Double_t xyz[], Double_t B[]) const = 0;
36 
37  //* Virtual methods needed for particle transportation
38  //* One can use particular implementations for collider (only Bz component)
39  //* geometry and for fixed-target (CBM-like) geometry which are provided below
40  //* in TRANSPORT section
41 
42  //* Get dS to xyz[] space point
43 
44  virtual Double_t GetDStoPoint(const Double_t xyz[]) const = 0;
45 
46  //* Get dS to other particle p (dSp for particle p also returned)
47 
48  virtual void GetDStoParticle(const KFParticleBase &p, Double_t &DS, Double_t &DSp) const = 0;
49 
50  //* Transport on dS value along trajectory, output to P,C
51 
52  virtual void Transport(Double_t dS, Double_t P[], Double_t C[]) const = 0;
53 
54  //*
55  //* INITIALIZATION
56  //*
57 
58  //* Constructor
59 
61 
62  //* Destructor
63 
64  virtual ~KFParticleBase() { ; }
65 
66  //* Initialisation from "cartesian" coordinates ( X Y Z Px Py Pz )
67  //* Parameters, covariance matrix, charge, and mass hypothesis should be provided
68 
69  void Initialize(const Double_t Param[], const Double_t Cov[], Int_t Charge, Double_t Mass);
70 
71  //* Initialise covariance matrix and set current parameters to 0.0
72 
73  void Initialize();
74 
75  //* Set decay vertex parameters for linearisation
76 
77  void SetVtxGuess(Double_t x, Double_t y, Double_t z);
78 
79  //* Set consruction method
80 
82 
83  //* Set and get mass hypothesis of the particle
84  void SetMassHypo(Double_t m) { fMassHypo = m; }
85  const Double_t &GetMassHypo() const { return fMassHypo; }
86 
87  //* Returns the sum of masses of the daughters
88  const Double_t &GetSumDaughterMass() const { return SumDaughterMass; }
89 
90  //*
91  //* ACCESSORS
92  //*
93 
94  //* Simple accessors
95 
96  Double_t GetX() const { return fP[0]; }
97  Double_t GetY() const { return fP[1]; }
98  Double_t GetZ() const { return fP[2]; }
99  Double_t GetPx() const { return fP[3]; }
100  Double_t GetPy() const { return fP[4]; }
101  Double_t GetPz() const { return fP[5]; }
102  Double_t GetE() const { return fP[6]; }
103  Double_t GetS() const { return fP[7]; }
104  Int_t GetQ() const { return fQ; }
105  Double_t GetChi2() const { return fChi2; }
106  Int_t GetNDF() const { return fNDF; }
107 
108  const Double_t &X() const { return fP[0]; }
109  const Double_t &Y() const { return fP[1]; }
110  const Double_t &Z() const { return fP[2]; }
111  const Double_t &Px() const { return fP[3]; }
112  const Double_t &Py() const { return fP[4]; }
113  const Double_t &Pz() const { return fP[5]; }
114  const Double_t &E() const { return fP[6]; }
115  const Double_t &S() const { return fP[7]; }
116  const Int_t &Q() const { return fQ; }
117  const Double_t &Chi2() const { return fChi2; }
118  const Int_t &NDF() const { return fNDF; }
119 
120  Double_t GetParameter(Int_t i) const { return fP[i]; }
121  Double_t GetCovariance(Int_t i) const { return fC[i]; }
122  Double_t GetCovariance(Int_t i, Int_t j) const { return fC[IJ(i, j)]; }
123 
124  //* Accessors with calculations( &value, &estimated sigma )
125  //* error flag returned (0 means no error during calculations)
126 
127  Int_t GetMomentum(Double_t &P, Double_t &SigmaP) const;
128  Int_t GetPt(Double_t &Pt, Double_t &SigmaPt) const;
129  Int_t GetEta(Double_t &Eta, Double_t &SigmaEta) const;
130  Int_t GetPhi(Double_t &Phi, Double_t &SigmaPhi) const;
131  Int_t GetMass(Double_t &M, Double_t &SigmaM) const;
132  Int_t GetDecayLength(Double_t &L, Double_t &SigmaL) const;
133  Int_t GetDecayLengthXY(Double_t &L, Double_t &SigmaL) const;
134  Int_t GetLifeTime(Double_t &T, Double_t &SigmaT) const;
135  Int_t GetR(Double_t &R, Double_t &SigmaR) const;
136 
137  //*
138  //* MODIFIERS
139  //*
140 
141  Double_t &X() { return fP[0]; }
142  Double_t &Y() { return fP[1]; }
143  Double_t &Z() { return fP[2]; }
144  Double_t &Px() { return fP[3]; }
145  Double_t &Py() { return fP[4]; }
146  Double_t &Pz() { return fP[5]; }
147  Double_t &E() { return fP[6]; }
148  Double_t &S() { return fP[7]; }
149  Int_t &Q() { return fQ; }
150  Double_t &Chi2() { return fChi2; }
151  Int_t &NDF() { return fNDF; }
152 
153  Double_t &Parameter(Int_t i) { return fP[i]; }
154  Double_t &Covariance(Int_t i) { return fC[i]; }
155  Double_t &Covariance(Int_t i, Int_t j) { return fC[IJ(i, j)]; }
156 
157  //*
158  //* CONSTRUCTION OF THE PARTICLE BY ITS DAUGHTERS AND MOTHER
159  //* USING THE KALMAN FILTER METHOD
160  //*
161 
162  //* Simple way to add daughter ex. D0+= Pion;
163 
164  void operator+=(const KFParticleBase &Daughter);
165 
166  //* Add daughter track to the particle
167 
168  void AddDaughter(const KFParticleBase &Daughter);
169 
170  void AddDaughterWithEnergyFit(const KFParticleBase &Daughter);
171  void AddDaughterWithEnergyCalc(const KFParticleBase &Daughter);
172  void AddDaughterWithEnergyFitMC(const KFParticleBase &Daughter); // with mass constrained
173 
174  //* Set production vertex
175 
176  void SetProductionVertex(const KFParticleBase &Vtx);
177 
178  //* Set mass constraint
179 
180  void SetNonlinearMassConstraint(Double_t Mass);
181  void SetMassConstraint(Double_t Mass, Double_t SigmaMass = 0);
182 
183  //* Set no decay length for resonances
184 
185  void SetNoDecayLength();
186 
187  //* Everything in one go
188 
189  void Construct(const KFParticleBase *vDaughters[], Int_t nDaughters, const KFParticleBase *ProdVtx = nullptr, Double_t Mass = -1, Bool_t IsConstrained = 0);
190 
191  //*
192  //* TRANSPORT
193  //*
194  //* ( main transportation parameter is S = SignedPath/Momentum )
195  //* ( parameters of decay & production vertices are stored locally )
196  //*
197 
198  //* Transport the particle to its decay vertex
199 
200  void TransportToDecayVertex();
201 
202  //* Transport the particle to its production vertex
203 
205 
206  //* Transport the particle on dS parameter (SignedPath/Momentum)
207 
208  void TransportToDS(Double_t dS);
209 
210  //* Particular extrapolators one can use
211 
212  Double_t GetDStoPointBz(Double_t Bz, const Double_t xyz[]) const;
213  Double_t GetDStoPointBy(Double_t By, const Double_t xyz[]) const;
214 
215  void GetDStoParticleBz(Double_t Bz, const KFParticleBase &p, Double_t &dS, Double_t &dS1) const;
216  void GetDStoParticleBy(Double_t B, const KFParticleBase &p, Double_t &dS, Double_t &dS1) const;
217 
218  Double_t GetDStoPointCBM(const Double_t xyz[]) const;
219  void GetDStoParticleCBM(const KFParticleBase &p, Double_t &dS, Double_t &dS1) const;
220 
221  void TransportBz(Double_t Bz, Double_t dS, Double_t P[], Double_t C[]) const;
222  void TransportCBM(Double_t dS, Double_t P[], Double_t C[]) const;
223 
224  //*
225  //* OTHER UTILITIES
226  //*
227 
228  //* Calculate distance from another object [cm]
229 
230  Double_t GetDistanceFromVertex(const Double_t vtx[]) const;
231  Double_t GetDistanceFromVertex(const KFParticleBase &Vtx) const;
232  Double_t GetDistanceFromParticle(const KFParticleBase &p) const;
233 
234  //* Calculate sqrt(Chi2/ndf) deviation from vertex
235  //* v = [xyz], Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix
236 
237  Double_t GetDeviationFromVertex(const Double_t v[], const Double_t Cv[] = 0) const;
238  Double_t GetDeviationFromVertex(const KFParticleBase &Vtx) const;
239  Double_t GetDeviationFromParticle(const KFParticleBase &p) const;
240 
241  //* Subtract the particle from the vertex
242 
243  void SubtractFromVertex(KFParticleBase &Vtx) const;
244 
245  //* Special method for creating gammas
246 
247  void ConstructGammaBz(const KFParticleBase &daughter1, const KFParticleBase &daughter2, double Bz);
248 
249  //* return parameters for the Armenteros-Podolanski plot
250  static void GetArmenterosPodolanski(KFParticleBase &positive, KFParticleBase &negative, Double_t QtAlfa[2]);
251 
252  //* Rotates the KFParticle object around OZ axis, OZ axis is set by the vertex position
253  void RotateXY(Double_t angle, Double_t Vtx[3]);
254 
255  int Id() const { return fId; };
256  int NDaughters() const { return fDaughtersIds.size(); };
257  const std::vector<int> &DaughterIds() const { return fDaughtersIds; };
258 
259  void SetId(int id) { fId = id; }; // should be always used (manualy)
260  void AddDaughterId(int id) { fDaughtersIds.push_back(id); };
261 
262  void SetPDG(int pdg) { fPDG = pdg; }
263  int GetPDG() const { return fPDG; }
264 
265  protected:
266  static Int_t IJ(Int_t i, Int_t j) { return (j <= i) ? i * (i + 1) / 2 + j : j * (j + 1) / 2 + i; }
267 
268  Double_t &Cij(Int_t i, Int_t j) { return fC[IJ(i, j)]; }
269 
270  void Convert(bool ToProduction);
271  void TransportLine(Double_t S, Double_t P[], Double_t C[]) const;
272  Double_t GetDStoPointLine(const Double_t xyz[]) const;
273  void GetDStoParticleLine(const KFParticleBase &p, Double_t &dS, Double_t &dS1) const;
274 
275  void GetDSIter(const KFParticleBase &p, Double_t const &dS, Double_t x[3], Double_t dx[3], Double_t ddx[3]) const;
276 
277  static Bool_t InvertSym3(const Double_t A[], Double_t Ainv[]);
278 
279  static void MultQSQt(const Double_t Q[], const Double_t S[], Double_t SOut[]);
280 
281  static Double_t GetSCorrection(const Double_t Part[], const Double_t XYZ[]);
282 
283  void GetMeasurement(const Double_t XYZ[], Double_t m[], Double_t V[]) const;
284 
285  //* Mass constraint function. is needed for the nonlinear mass constraint and a fit with mass constraint
286  void SetMassConstraint(Double_t *mP, Double_t *mC, Double_t mJ[7][7], Double_t mass);
287 
288  Double_t fP[8]; //* Main particle parameters {X,Y,Z,Px,Py,Pz,E,S[=DecayLength/P]}
289  Double_t fC[36]; //* Low-triangle covariance matrix of fP
290  Int_t fQ; //* Particle charge
291  Int_t fNDF; //* Number of degrees of freedom
292  Double_t fChi2; //* Chi^2
293 
294  Double_t fSFromDecay; //* Distance from decay vertex to current position
295 
296  Bool_t fAtProductionVertex; //* Flag shows that the particle error along
297  //* its trajectory is taken from production vertex
298 
299  Double_t fVtxGuess[3]; //* Guess for the position of the decay vertex
300  //* ( used for linearisation of equations )
301 
302  Bool_t fIsLinearized; //* Flag shows that the guess is present
303 
304  Int_t fConstructMethod; //* Determines the method for the particle construction.
305  //* 0 - Energy considered as an independent veriable, fitted independently from momentum, without any constraints on mass
306  //* 1 - Energy considered as a dependent variable, calculated from the momentum and mass hypothesis
307  //* 2 - Energy considered as an independent variable, fitted independently from momentum, with constraints on mass of daughter particle
308 
309  Double_t SumDaughterMass; //* sum of the daughter particles masses
310  Double_t fMassHypo; //* sum of the daughter particles masse
311 
312  int fId; // id of particle
313  std::vector<int> fDaughtersIds; // id of particles it created from. if size == 1 then this is id of track. TODO use in functions. why unsigned short int doesn't work???
314 
315  int fPDG; // pdg hypothesis
316 
317  // ClassDef( KFParticleBase, 1 );
318 };
319 
320 #endif
Double_t GetPy() const
Double_t GetDistanceFromVertex(const Double_t vtx[]) const
Double_t GetZ() const
Double_t & Y()
void SetConstructMethod(Int_t m)
void TransportToProductionVertex()
void TransportToDecayVertex()
static void GetArmenterosPodolanski(KFParticleBase &positive, KFParticleBase &negative, Double_t QtAlfa[2])
void GetDStoParticleLine(const KFParticleBase &p, Double_t &dS, Double_t &dS1) const
void SetNonlinearMassConstraint(Double_t Mass)
__m128 m
Definition: P4_F32vec4.h:26
Int_t GetQ() const
Double_t GetPx() const
Double_t fVtxGuess[3]
Int_t GetMomentum(Double_t &P, Double_t &SigmaP) const
Double_t GetPz() const
Double_t fP[8]
void TransportCBM(Double_t dS, Double_t P[], Double_t C[]) const
Int_t GetPhi(Double_t &Phi, Double_t &SigmaPhi) const
Double_t & S()
Double_t & Chi2()
const Int_t & Q() const
Double_t GetDStoPointBz(Double_t Bz, const Double_t xyz[]) const
const Double_t & Pz() const
void GetDStoParticleCBM(const KFParticleBase &p, Double_t &dS, Double_t &dS1) const
Double_t GetDStoPointBy(Double_t By, const Double_t xyz[]) const
Int_t GetMass(Double_t &M, Double_t &SigmaM) const
Double_t GetY() const
const Double_t & Py() const
void SubtractFromVertex(KFParticleBase &Vtx) const
void SetVtxGuess(Double_t x, Double_t y, Double_t z)
Double_t GetParameter(Int_t i) const
void operator+=(const KFParticleBase &Daughter)
const Double_t & Chi2() const
void AddDaughterWithEnergyCalc(const KFParticleBase &Daughter)
static Int_t IJ(Int_t i, Int_t j)
const Double_t & Y() const
__m128 v
Definition: P4_F32vec4.h:3
static Double_t GetSCorrection(const Double_t Part[], const Double_t XYZ[])
unsigned int i
Definition: P4_F32vec4.h:21
Double_t & X()
Int_t GetNDF() const
int Id() const
void SetNoDecayLength()
virtual void GetDStoParticle(const KFParticleBase &p, Double_t &DS, Double_t &DSp) const =0
Double_t & Cij(Int_t i, Int_t j)
int GetPDG() const
Int_t GetEta(Double_t &Eta, Double_t &SigmaEta) const
Int_t GetPt(Double_t &Pt, Double_t &SigmaPt) const
void GetDStoParticleBy(Double_t B, const KFParticleBase &p, Double_t &dS, Double_t &dS1) const
void SetMassHypo(Double_t m)
void TransportLine(Double_t S, Double_t P[], Double_t C[]) const
void SetProductionVertex(const KFParticleBase &Vtx)
virtual void Transport(Double_t dS, Double_t P[], Double_t C[]) const =0
Double_t & E()
Double_t SumDaughterMass
virtual void GetFieldValue(const Double_t xyz[], Double_t B[]) const =0
const Double_t & GetMassHypo() const
void GetDStoParticleBz(Double_t Bz, const KFParticleBase &p, Double_t &dS, Double_t &dS1) const
Double_t GetCovariance(Int_t i, Int_t j) const
void SetMassConstraint(Double_t Mass, Double_t SigmaMass=0)
Double_t fSFromDecay
Int_t GetDecayLengthXY(Double_t &L, Double_t &SigmaL) const
const Int_t & NDF() const
void RotateXY(Double_t angle, Double_t Vtx[3])
Double_t GetDeviationFromVertex(const Double_t v[], const Double_t Cv[]=0) const
Double_t fC[36]
virtual ~KFParticleBase()
const std::vector< int > & DaughterIds() const
int NDaughters() const
static void MultQSQt(const Double_t Q[], const Double_t S[], Double_t SOut[])
Double_t GetCovariance(Int_t i) const
Int_t GetR(Double_t &R, Double_t &SigmaR) const
void SetPDG(int pdg)
Double_t GetChi2() const
void AddDaughterWithEnergyFit(const KFParticleBase &Daughter)
const Double_t & Px() const
Double_t GetDStoPointCBM(const Double_t xyz[]) const
const Double_t & GetSumDaughterMass() const
void AddDaughterId(int id)
void Convert(bool ToProduction)
Double_t GetX() const
void AddDaughterWithEnergyFitMC(const KFParticleBase &Daughter)
Double_t GetE() const
void Construct(const KFParticleBase *vDaughters[], Int_t nDaughters, const KFParticleBase *ProdVtx=nullptr, Double_t Mass=-1, Bool_t IsConstrained=0)
void GetMeasurement(const Double_t XYZ[], Double_t m[], Double_t V[]) const
void TransportToDS(Double_t dS)
Double_t GetS() const
void GetDSIter(const KFParticleBase &p, Double_t const &dS, Double_t x[3], Double_t dx[3], Double_t ddx[3]) const
const Double_t & X() const
Double_t & Px()
Double_t GetDeviationFromParticle(const KFParticleBase &p) const
Double_t & Z()
static Bool_t InvertSym3(const Double_t A[], Double_t Ainv[])
Bool_t fAtProductionVertex
Double_t & Covariance(Int_t i, Int_t j)
Double_t fMassHypo
Double_t GetDStoPointLine(const Double_t xyz[]) const
void AddDaughter(const KFParticleBase &Daughter)
Double_t & Parameter(Int_t i)
Double_t & Pz()
std::vector< int > fDaughtersIds
Double_t GetDistanceFromParticle(const KFParticleBase &p) const
void SetId(int id)
const Double_t & S() const
Int_t GetDecayLength(Double_t &L, Double_t &SigmaL) const
Double_t & Py()
const Double_t & E() const
const Double_t & Z() const
void ConstructGammaBz(const KFParticleBase &daughter1, const KFParticleBase &daughter2, double Bz)
virtual Double_t GetDStoPoint(const Double_t xyz[]) const =0
Double_t & Covariance(Int_t i)
void TransportBz(Double_t Bz, Double_t dS, Double_t P[], Double_t C[]) const
Int_t GetLifeTime(Double_t &T, Double_t &SigmaT) const