]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/ESD/AliKFParticleBase.h
restore threshold getters after parameter dynamics update (fw v. >= A012)
[u/mrichter/AliRoot.git] / STEER / ESD / AliKFParticleBase.h
CommitLineData
f826d409 1//---------------------------------------------------------------------------------
2// The AliKFParticleBase class
3// .
4// @author S.Gorbunov, I.Kisel
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 AliKFParticle class
14//
15// -= Copyright &copy ALICE HLT Group =-
16//_________________________________________________________________________________
17
18
19
20#ifndef ALIKFPARTICLEBASE_H
21#define ALIKFPARTICLEBASE_H
22
23#include "TObject.h"
24
25class AliKFParticleBase :public TObject {
26
27 public:
28
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 AliKFParticleBase &p,
49 Double_t &DS, Double_t &DSp ) const = 0;
50
51 //* Transport on dS value along trajectory, output to P,C
52
53 virtual void Transport( Double_t dS, Double_t P[], Double_t C[] ) const = 0;
54
55
56
57 //*
58 //* INITIALIZATION
59 //*
60
61 //* Constructor
62
63 AliKFParticleBase();
64
65 //* Destructor
66
e7b09c95 67 virtual ~AliKFParticleBase() { ; }
68
69 //* Initialisation from "cartesian" coordinates ( X Y Z Px Py Pz )
70 //* Parameters, covariance matrix, charge, and mass hypothesis should be provided
71
72 void Initialize( const Double_t Param[], const Double_t Cov[], Int_t Charge, Double_t Mass );
f826d409 73
74 //* Initialise covariance matrix and set current parameters to 0.0
75
76 void Initialize();
77
78 //* Set decay vertex parameters for linearisation
79
80 void SetVtxGuess( Double_t x, Double_t y, Double_t z );
81
82 //*
83 //* ACCESSORS
84 //*
85
86 //* Simple accessors
87
88 Double_t GetX () const { return fP[0]; }
89 Double_t GetY () const { return fP[1]; }
90 Double_t GetZ () const { return fP[2]; }
91 Double_t GetPx () const { return fP[3]; }
92 Double_t GetPy () const { return fP[4]; }
93 Double_t GetPz () const { return fP[5]; }
94 Double_t GetE () const { return fP[6]; }
95 Double_t GetS () const { return fP[7]; }
96 Int_t GetQ () const { return fQ; }
97 Double_t GetChi2 () const { return fChi2; }
98 Int_t GetNDF () const { return fNDF; }
91b25235 99
100 const Double_t& X () const { return fP[0]; }
101 const Double_t& Y () const { return fP[1]; }
102 const Double_t& Z () const { return fP[2]; }
103 const Double_t& Px () const { return fP[3]; }
104 const Double_t& Py () const { return fP[4]; }
105 const Double_t& Pz () const { return fP[5]; }
106 const Double_t& E () const { return fP[6]; }
107 const Double_t& S () const { return fP[7]; }
108 const Int_t & Q () const { return fQ; }
109 const Double_t& Chi2 () const { return fChi2; }
110 const Int_t & NDF () const { return fNDF; }
111
f826d409 112
113 Double_t GetParameter ( Int_t i ) const { return fP[i]; }
114 Double_t GetCovariance( Int_t i ) const { return fC[i]; }
115 Double_t GetCovariance( Int_t i, Int_t j ) const { return fC[IJ(i,j)]; }
116
117 //* Accessors with calculations( &value, &estimated sigma )
118 //* error flag returned (0 means no error during calculations)
119
706952f5 120 Int_t GetMomentum ( Double_t &P, Double_t &SigmaP ) const ;
121 Int_t GetPt ( Double_t &Pt, Double_t &SigmaPt ) const ;
122 Int_t GetEta ( Double_t &Eta, Double_t &SigmaEta ) const ;
123 Int_t GetPhi ( Double_t &Phi, Double_t &SigmaPhi ) const ;
124 Int_t GetMass ( Double_t &M, Double_t &SigmaM ) const ;
125 Int_t GetDecayLength ( Double_t &L, Double_t &SigmaL ) const ;
446ce366 126 Int_t GetDecayLengthXY ( Double_t &L, Double_t &SigmaL ) const ;
706952f5 127 Int_t GetLifeTime ( Double_t &T, Double_t &SigmaT ) const ;
128 Int_t GetR ( Double_t &R, Double_t &SigmaR ) const ;
f826d409 129
130 //*
131 //* MODIFIERS
132 //*
133
134 Double_t & X () { return fP[0]; }
135 Double_t & Y () { return fP[1]; }
136 Double_t & Z () { return fP[2]; }
137 Double_t & Px () { return fP[3]; }
138 Double_t & Py () { return fP[4]; }
139 Double_t & Pz () { return fP[5]; }
140 Double_t & E () { return fP[6]; }
141 Double_t & S () { return fP[7]; }
142 Int_t & Q () { return fQ; }
143 Double_t & Chi2 () { return fChi2; }
144 Int_t & NDF () { return fNDF; }
145
91b25235 146
147
f826d409 148 Double_t & Parameter ( Int_t i ) { return fP[i]; }
149 Double_t & Covariance( Int_t i ) { return fC[i]; }
150 Double_t & Covariance( Int_t i, Int_t j ) { return fC[IJ(i,j)]; }
151
152
153 //*
154 //* CONSTRUCTION OF THE PARTICLE BY ITS DAUGHTERS AND MOTHER
155 //* USING THE KALMAN FILTER METHOD
156 //*
157
158
159 //* Simple way to add daughter ex. D0+= Pion;
160
161 void operator +=( const AliKFParticleBase &Daughter );
162
163 //* Add daughter track to the particle
164
165 void AddDaughter( const AliKFParticleBase &Daughter );
166
167 //* Set production vertex
168
169 void SetProductionVertex( const AliKFParticleBase &Vtx );
170
e7b09c95 171 //* Set mass constraint
f826d409 172
e7b09c95 173 void SetMassConstraint( Double_t Mass, Double_t SigmaMass = 0 );
f826d409 174
e7b09c95 175 //* Set no decay length for resonances
176
177 void SetNoDecayLength();
178
179
f826d409 180 //* Everything in one go
181
182 void Construct( const AliKFParticleBase *vDaughters[], Int_t NDaughters,
706952f5 183 const AliKFParticleBase *ProdVtx=0, Double_t Mass=-1, Bool_t IsConstrained=0 );
f826d409 184
185
186 //*
187 //* TRANSPORT
188 //*
189 //* ( main transportation parameter is S = SignedPath/Momentum )
190 //* ( parameters of decay & production vertices are stored locally )
191 //*
192
193
194 //* Transport the particle to its decay vertex
195
196 void TransportToDecayVertex();
197
198 //* Transport the particle to its production vertex
199
200 void TransportToProductionVertex();
201
202 //* Transport the particle on dS parameter (SignedPath/Momentum)
203
204 void TransportToDS( Double_t dS );
205
206 //* Particular extrapolators one can use
207
208 Double_t GetDStoPointBz( Double_t Bz, const Double_t xyz[] ) const;
209
210 void GetDStoParticleBz( Double_t Bz, const AliKFParticleBase &p,
211 Double_t &dS, Double_t &dS1 ) const ;
212
213 // Double_t GetDStoPointCBM( const Double_t xyz[] ) const;
214
215 void TransportBz( Double_t Bz, Double_t dS, Double_t P[], Double_t C[] ) const;
216 void TransportCBM( Double_t dS, Double_t P[], Double_t C[] ) const;
217
218
219 //*
220 //* OTHER UTILITIES
221 //*
222
223 //* Calculate distance from another object [cm]
224
225 Double_t GetDistanceFromVertex( const Double_t vtx[] ) const;
226 Double_t GetDistanceFromVertex( const AliKFParticleBase &Vtx ) const;
616ffc76 227 Double_t GetDistanceFromParticle( const AliKFParticleBase &p ) const;
f826d409 228
229 //* Calculate sqrt(Chi2/ndf) deviation from vertex
230 //* v = [xyz], Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix
231
232 Double_t GetDeviationFromVertex( const Double_t v[],
233 const Double_t Cv[]=0 ) const;
234 Double_t GetDeviationFromVertex( const AliKFParticleBase &Vtx ) const;
616ffc76 235 Double_t GetDeviationFromParticle( const AliKFParticleBase &p ) const;
f826d409 236
237 //* Subtract the particle from the vertex
238
de0d0ceb 239 void SubtractFromVertex( AliKFParticleBase &Vtx ) const;
240
a65041d0 241 //* Special method for creating gammas
242
243 void ConstructGammaBz( const AliKFParticleBase &daughter1,
244 const AliKFParticleBase &daughter2, double Bz );
245
f826d409 246 protected:
247
248 static Int_t IJ( Int_t i, Int_t j ){
249 return ( j<=i ) ? i*(i+1)/2+j :j*(j+1)/2+i;
250 }
251
252 Double_t & Cij( Int_t i, Int_t j ){ return fC[IJ(i,j)]; }
253
254 void Convert( bool ToProduction );
255 void TransportLine( Double_t S, Double_t P[], Double_t C[] ) const ;
256 Double_t GetDStoPointLine( const Double_t xyz[] ) const;
257
a65041d0 258 static Bool_t InvertSym3( const Double_t A[], Double_t Ainv[] );
259
f826d409 260 static void MultQSQt( const Double_t Q[], const Double_t S[],
261 Double_t SOut[] );
262
706952f5 263 static Double_t GetSCorrection( const Double_t Part[], const Double_t XYZ[] );
264
616ffc76 265 void GetMeasurement( const Double_t XYZ[], Double_t m[], Double_t V[] ) const ;
f826d409 266
267 Double_t fP[8]; //* Main particle parameters {X,Y,Z,Px,Py,Pz,E,S[=DecayLength/P]}
268 Double_t fC[36]; //* Low-triangle covariance matrix of fP
269 Int_t fQ; //* Particle charge
270 Int_t fNDF; //* Number of degrees of freedom
271 Double_t fChi2; //* Chi^2
272
273 Double_t fSFromDecay; //* Distance from decay vertex to current position
274
275 Bool_t fAtProductionVertex; //* Flag shows that the particle error along
276 //* its trajectory is taken from production vertex
277
278 Double_t fVtxGuess[3]; //* Guess for the position of the decay vertex
279 //* ( used for linearisation of equations )
280
281 Bool_t fIsLinearized; //* Flag shows that the guess is present
282
283 ClassDef( AliKFParticleBase, 1 );
284};
285
286#endif