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