]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/ESD/AliKFParticleBase.h
Update timestamp for new data points simulation
[u/mrichter/AliRoot.git] / STEER / ESD / AliKFParticleBase.h
1 //---------------------------------------------------------------------------------
2 // The AliKFParticleBase 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 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
25 class 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
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 );
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   //* Set consruction method
83
84   void SetConstructMethod(Int_t m) {fConstructMethod = m;}
85
86   //* Set and get mass hypothesis of the particle
87   void SetMassHypo(Double_t m) { fMassHypo = m;}
88   const Double_t& GetMassHypo() const { return fMassHypo; }
89
90   //* Returns the sum of masses of the daughters
91   const Double_t& GetSumDaughterMass() const {return SumDaughterMass;}
92
93   //*
94   //*  ACCESSORS
95   //*
96
97   //* Simple accessors 
98
99   Double_t GetX    () const { return fP[0]; }
100   Double_t GetY    () const { return fP[1]; }
101   Double_t GetZ    () const { return fP[2]; }
102   Double_t GetPx   () const { return fP[3]; }
103   Double_t GetPy   () const { return fP[4]; }
104   Double_t GetPz   () const { return fP[5]; }
105   Double_t GetE    () const { return fP[6]; }
106   Double_t GetS    () const { return fP[7]; }
107   Int_t    GetQ    () const { return fQ;    }
108   Double_t GetChi2 () const { return fChi2; }
109   Int_t    GetNDF  () const { return fNDF;  }
110
111   const Double_t& X    () const { return fP[0]; }
112   const Double_t& Y    () const { return fP[1]; }
113   const Double_t& Z    () const { return fP[2]; }
114   const Double_t& Px   () const { return fP[3]; }
115   const Double_t& Py   () const { return fP[4]; }
116   const Double_t& Pz   () const { return fP[5]; }
117   const Double_t& E    () const { return fP[6]; }
118   const Double_t& S    () const { return fP[7]; }
119   const Int_t   & Q    () const { return fQ;    }
120   const Double_t& Chi2 () const { return fChi2; }
121   const Int_t   & NDF  () const { return fNDF;  }
122
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   
158
159   Double_t & Parameter ( Int_t i )        { return fP[i];       }
160   Double_t & Covariance( Int_t i )        { return fC[i];       }
161   Double_t & Covariance( Int_t i, Int_t j ) { return fC[IJ(i,j)]; }
162
163
164   //* 
165   //* CONSTRUCTION OF THE PARTICLE BY ITS DAUGHTERS AND MOTHER
166   //* USING THE KALMAN FILTER METHOD
167   //*
168
169
170   //* Simple way to add daughter ex. D0+= Pion; 
171
172   void operator +=( const AliKFParticleBase &Daughter );  
173
174   //* Add daughter track to the particle 
175
176   void AddDaughter( const AliKFParticleBase &Daughter );
177
178   void AddDaughterWithEnergyFit( const AliKFParticleBase &Daughter );
179   void AddDaughterWithEnergyCalc( const AliKFParticleBase &Daughter );
180   void AddDaughterWithEnergyFitMC( const AliKFParticleBase &Daughter ); //with mass constrained
181
182   //* Set production vertex 
183
184   void SetProductionVertex( const AliKFParticleBase &Vtx );
185
186   //* Set mass constraint 
187
188   void SetNonlinearMassConstraint( Double_t Mass );
189   void SetMassConstraint( Double_t Mass, Double_t SigmaMass = 0 );
190   
191   //* Set no decay length for resonances
192
193   void SetNoDecayLength();
194
195
196   //* Everything in one go  
197
198   void Construct( const AliKFParticleBase *vDaughters[], Int_t NDaughters, 
199                   const AliKFParticleBase *ProdVtx=0,   Double_t Mass=-1, Bool_t IsConstrained=0  );
200
201
202   //*
203   //*                   TRANSPORT
204   //* 
205   //*  ( main transportation parameter is S = SignedPath/Momentum )
206   //*  ( parameters of decay & production vertices are stored locally )
207   //*
208
209
210   //* Transport the particle to its decay vertex 
211
212   void TransportToDecayVertex();
213
214   //* Transport the particle to its production vertex 
215
216   void TransportToProductionVertex();
217
218   //* Transport the particle on dS parameter (SignedPath/Momentum) 
219
220   void TransportToDS( Double_t dS );
221
222   //* Particular extrapolators one can use 
223
224   Double_t GetDStoPointBz( Double_t Bz, const Double_t xyz[] ) const;
225   
226   void GetDStoParticleBz( Double_t Bz, const AliKFParticleBase &p, 
227                           Double_t &dS, Double_t &dS1       ) const ;
228  
229   // Double_t GetDStoPointCBM( const Double_t xyz[] ) 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 AliKFParticleBase &Vtx ) const;
243   Double_t GetDistanceFromParticle( const AliKFParticleBase &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 AliKFParticleBase &Vtx ) const;
251   Double_t GetDeviationFromParticle( const AliKFParticleBase &p ) const;  
252
253   //* Subtract the particle from the vertex  
254
255   void SubtractFromVertex( AliKFParticleBase &Vtx ) const;
256
257   //* Special method for creating gammas
258
259   void ConstructGammaBz( const AliKFParticleBase &daughter1,
260                          const AliKFParticleBase &daughter2, double Bz  );
261
262   //* return parameters for the Armenteros-Podolanski plot
263   static void GetArmenterosPodolanski(AliKFParticleBase& positive, AliKFParticleBase& 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  protected:
269
270   static Int_t IJ( Int_t i, Int_t j ){ 
271     return ( j<=i ) ? i*(i+1)/2+j :j*(j+1)/2+i;
272   }
273
274   Double_t & Cij( Int_t i, Int_t j ){ return fC[IJ(i,j)]; }
275
276   void Convert( bool ToProduction );
277   void TransportLine( Double_t S, Double_t P[], Double_t C[] ) const ;
278   Double_t GetDStoPointLine( const Double_t xyz[] ) const;
279
280   static Bool_t InvertSym3( const Double_t A[], Double_t Ainv[] );
281
282   static void MultQSQt( const Double_t Q[], const Double_t S[], 
283                         Double_t SOut[] );
284
285   static Double_t GetSCorrection( const Double_t Part[], const Double_t XYZ[] );
286
287   void GetMeasurement( const Double_t XYZ[], Double_t m[], Double_t V[] ) const ;
288
289   //* Mass constraint function. is needed for the nonlinear mass constraint and a fit with mass constraint
290   void SetMassConstraint( Double_t *mP, Double_t *mC, Double_t mJ[7][7], Double_t mass );
291
292   Double_t fP[8];  //* Main particle parameters {X,Y,Z,Px,Py,Pz,E,S[=DecayLength/P]}
293   Double_t fC[36]; //* Low-triangle covariance matrix of fP
294   Int_t    fQ;     //* Particle charge 
295   Int_t    fNDF;   //* Number of degrees of freedom 
296   Double_t fChi2;  //* Chi^2
297
298   Double_t fSFromDecay; //* Distance from decay vertex to current position
299
300   Bool_t fAtProductionVertex; //* Flag shows that the particle error along
301                               //* its trajectory is taken from production vertex    
302
303   Double_t fVtxGuess[3];  //* Guess for the position of the decay vertex 
304                           //* ( used for linearisation of equations )
305
306   Bool_t fIsLinearized;   //* Flag shows that the guess is present
307
308   Int_t fConstructMethod; //* Determines the method for the particle construction. 
309   //* 0 - Energy considered as an independent veriable, fitted independently from momentum, without any constraints on mass
310   //* 1 - Energy considered as a dependent variable, calculated from the momentum and mass hypothesis
311   //* 2 - Energy considered as an independent variable, fitted independently from momentum, with constraints on mass of daughter particle
312
313   Double_t SumDaughterMass;  //* sum of the daughter particles masses
314   Double_t fMassHypo;  //* sum of the daughter particles masses
315
316   ClassDef( AliKFParticleBase, 2);
317 };
318
319 #endif