]>
Commit | Line | Data |
---|---|---|
f826d409 | 1 | //--------------------------------------------------------------------------------- |
2 | // The AliKFParticle 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 is ALICE interface to general mathematics in AliKFParticleBase | |
14 | // | |
15 | // -= Copyright © ALICE HLT Group =- | |
16 | //_________________________________________________________________________________ | |
17 | ||
18 | #ifndef ALIKFPARTICLE_H | |
19 | #define ALIKFPARTICLE_H | |
20 | ||
21 | #include "AliKFParticleBase.h" | |
cfb5a2e1 | 22 | #include "TMath.h" |
f826d409 | 23 | |
706952f5 | 24 | class AliVTrack; |
cfb5a2e1 | 25 | class AliVVertex; |
706fa2e5 | 26 | class AliExternalTrackParam; |
f826d409 | 27 | |
f826d409 | 28 | class AliKFParticle :public AliKFParticleBase |
29 | { | |
30 | ||
31 | public: | |
32 | ||
33 | //* | |
34 | //* INITIALIZATION | |
35 | //* | |
36 | ||
4bbc290d | 37 | //* Set magnetic field for all particles |
38 | ||
39 | static void SetField( Double_t Bz ); | |
40 | ||
f826d409 | 41 | //* Constructor (empty) |
42 | ||
e7b09c95 | 43 | AliKFParticle():AliKFParticleBase(){ ; } |
f826d409 | 44 | |
45 | //* Destructor (empty) | |
46 | ||
e7b09c95 | 47 | ~AliKFParticle(){ ; } |
f826d409 | 48 | |
616ffc76 | 49 | //* Construction of mother particle by its 2-3-4 daughters |
50 | ||
537b06a4 | 51 | AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2, Bool_t gamma = kFALSE ); |
616ffc76 | 52 | |
53 | AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2, | |
54 | const AliKFParticle &d3 ); | |
55 | ||
56 | AliKFParticle( const AliKFParticle &d1, const AliKFParticle &d2, | |
57 | const AliKFParticle &d3, const AliKFParticle &d4 ); | |
58 | ||
e7b09c95 | 59 | //* Initialisation from "cartesian" coordinates ( X Y Z Px Py Pz ) |
60 | //* Parameters, covariance matrix, charge and PID hypothesis should be provided | |
61 | ||
7972a8db | 62 | void Create( const Double_t Param[], const Double_t Cov[], Int_t Charge, Int_t PID ); |
73599358 | 63 | void Create( const Double_t Param[], const Double_t Cov[], Int_t Charge, Double_t Mass ); |
e7b09c95 | 64 | |
65 | //* Initialisation from ALICE track, PID hypothesis shoould be provided | |
f826d409 | 66 | |
706952f5 | 67 | AliKFParticle( const AliVTrack &track, Int_t PID ); |
73599358 | 68 | AliKFParticle( const AliExternalTrackParam &track, Double_t Mass, Int_t Charge ); |
f826d409 | 69 | |
706952f5 | 70 | //* Initialisation from VVertex |
f826d409 | 71 | |
706952f5 | 72 | AliKFParticle( const AliVVertex &vertex ); |
f826d409 | 73 | |
f826d409 | 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 ; //* x of current position | |
89 | Double_t GetY () const ; //* y of current position | |
90 | Double_t GetZ () const ; //* z of current position | |
91 | Double_t GetPx () const ; //* x-compoment of 3-momentum | |
92 | Double_t GetPy () const ; //* y-compoment of 3-momentum | |
93 | Double_t GetPz () const ; //* z-compoment of 3-momentum | |
94 | Double_t GetE () const ; //* energy | |
95 | Double_t GetS () const ; //* decay length / momentum | |
96 | Int_t GetQ () const ; //* charge | |
97 | Double_t GetChi2 () const ; //* chi^2 | |
98 | Int_t GetNDF () const ; //* Number of Degrees of Freedom | |
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; } | |
f826d409 | 111 | |
112 | Double_t GetParameter ( int i ) const ; | |
113 | Double_t GetCovariance( int i ) const ; | |
114 | Double_t GetCovariance( int i, int j ) const ; | |
115 | ||
e7b09c95 | 116 | //* Accessors with calculations, value returned w/o error flag |
117 | ||
706952f5 | 118 | Double_t GetP () const; //* momentum |
119 | Double_t GetPt () const; //* transverse momentum | |
120 | Double_t GetEta () const; //* pseudorapidity | |
121 | Double_t GetPhi () const; //* phi | |
122 | Double_t GetMomentum () const; //* momentum (same as GetP() ) | |
123 | Double_t GetMass () const; //* mass | |
124 | Double_t GetDecayLength () const; //* decay length | |
446ce366 | 125 | Double_t GetDecayLengthXY () const; //* decay length in XY |
706952f5 | 126 | Double_t GetLifeTime () const; //* life time |
127 | Double_t GetR () const; //* distance to the origin | |
e7b09c95 | 128 | |
129 | //* Accessors to estimated errors | |
130 | ||
706952f5 | 131 | Double_t GetErrX () const ; //* x of current position |
e7b09c95 | 132 | Double_t GetErrY () const ; //* y of current position |
133 | Double_t GetErrZ () const ; //* z of current position | |
134 | Double_t GetErrPx () const ; //* x-compoment of 3-momentum | |
135 | Double_t GetErrPy () const ; //* y-compoment of 3-momentum | |
136 | Double_t GetErrPz () const ; //* z-compoment of 3-momentum | |
137 | Double_t GetErrE () const ; //* energy | |
138 | Double_t GetErrS () const ; //* decay length / momentum | |
706952f5 | 139 | Double_t GetErrP () const ; //* momentum |
140 | Double_t GetErrPt () const ; //* transverse momentum | |
141 | Double_t GetErrEta () const ; //* pseudorapidity | |
142 | Double_t GetErrPhi () const ; //* phi | |
143 | Double_t GetErrMomentum () const ; //* momentum | |
144 | Double_t GetErrMass () const ; //* mass | |
145 | Double_t GetErrDecayLength () const ; //* decay length | |
446ce366 | 146 | Double_t GetErrDecayLengthXY () const ; //* decay length in XY |
706952f5 | 147 | Double_t GetErrLifeTime () const ; //* life time |
148 | Double_t GetErrR () const ; //* distance to the origin | |
e7b09c95 | 149 | |
f826d409 | 150 | //* Accessors with calculations( &value, &estimated sigma ) |
151 | //* error flag returned (0 means no error during calculations) | |
152 | ||
706952f5 | 153 | int GetP ( Double_t &P, Double_t &SigmaP ) const ; //* momentum |
154 | int GetPt ( Double_t &Pt, Double_t &SigmaPt ) const ; //* transverse momentum | |
155 | int GetEta ( Double_t &Eta, Double_t &SigmaEta ) const ; //* pseudorapidity | |
156 | int GetPhi ( Double_t &Phi, Double_t &SigmaPhi ) const ; //* phi | |
157 | int GetMomentum ( Double_t &P, Double_t &SigmaP ) const ; //* momentum | |
158 | int GetMass ( Double_t &M, Double_t &SigmaM ) const ; //* mass | |
159 | int GetDecayLength ( Double_t &L, Double_t &SigmaL ) const ; //* decay length | |
446ce366 | 160 | int GetDecayLengthXY ( Double_t &L, Double_t &SigmaL ) const ; //* decay length in XY |
161 | int GetLifeTime ( Double_t &T, Double_t &SigmaT ) const ; //* life time | |
706952f5 | 162 | int GetR ( Double_t &R, Double_t &SigmaR ) const ; //* R |
f826d409 | 163 | |
e7b09c95 | 164 | |
f826d409 | 165 | //* |
166 | //* MODIFIERS | |
167 | //* | |
168 | ||
169 | Double_t & X () ; | |
170 | Double_t & Y () ; | |
171 | Double_t & Z () ; | |
172 | Double_t & Px () ; | |
173 | Double_t & Py () ; | |
174 | Double_t & Pz () ; | |
175 | Double_t & E () ; | |
176 | Double_t & S () ; | |
177 | Int_t & Q () ; | |
178 | Double_t & Chi2 () ; | |
179 | Int_t & NDF () ; | |
180 | ||
181 | Double_t & Parameter ( int i ) ; | |
182 | Double_t & Covariance( int i ) ; | |
183 | Double_t & Covariance( int i, int j ) ; | |
bb4bc052 | 184 | Double_t * Parameters () ; |
185 | Double_t * CovarianceMatrix() ; | |
f826d409 | 186 | |
187 | //* | |
188 | //* CONSTRUCTION OF THE PARTICLE BY ITS DAUGHTERS AND MOTHER | |
189 | //* USING THE KALMAN FILTER METHOD | |
190 | //* | |
191 | ||
192 | ||
616ffc76 | 193 | //* Add daughter to the particle |
f826d409 | 194 | |
616ffc76 | 195 | void AddDaughter( const AliKFParticle &Daughter ); |
f826d409 | 196 | |
616ffc76 | 197 | //* Add daughter via += operator: ex.{ D0; D0+=Pion; D0+= Kaon; } |
f826d409 | 198 | |
616ffc76 | 199 | void operator +=( const AliKFParticle &Daughter ); |
f826d409 | 200 | |
201 | //* Set production vertex | |
202 | ||
203 | void SetProductionVertex( const AliKFParticle &Vtx ); | |
204 | ||
e7b09c95 | 205 | //* Set mass constraint |
f826d409 | 206 | |
e7b09c95 | 207 | void SetMassConstraint( Double_t Mass, Double_t SigmaMass = 0 ); |
f826d409 | 208 | |
e7b09c95 | 209 | //* Set no decay length for resonances |
210 | ||
211 | void SetNoDecayLength(); | |
212 | ||
f826d409 | 213 | //* Everything in one go |
214 | ||
215 | void Construct( const AliKFParticle *vDaughters[], int NDaughters, | |
706952f5 | 216 | const AliKFParticle *ProdVtx=0, Double_t Mass=-1, Bool_t IsConstrained=0 ); |
f826d409 | 217 | |
218 | //* | |
219 | //* TRANSPORT | |
220 | //* | |
221 | //* ( main transportation parameter is S = SignedPath/Momentum ) | |
222 | //* ( parameters of decay & production vertices are stored locally ) | |
223 | //* | |
224 | ||
225 | //* Transport the particle to its decay vertex | |
226 | ||
227 | void TransportToDecayVertex(); | |
228 | ||
229 | //* Transport the particle to its production vertex | |
230 | ||
231 | void TransportToProductionVertex(); | |
232 | ||
233 | //* Transport the particle close to xyz[] point | |
234 | ||
235 | void TransportToPoint( const Double_t xyz[] ); | |
236 | ||
706952f5 | 237 | //* Transport the particle close to VVertex |
f826d409 | 238 | |
706952f5 | 239 | void TransportToVertex( const AliVVertex &v ); |
f826d409 | 240 | |
241 | //* Transport the particle close to another particle p | |
242 | ||
243 | void TransportToParticle( const AliKFParticle &p ); | |
244 | ||
245 | //* Transport the particle on dS parameter (SignedPath/Momentum) | |
246 | ||
247 | void TransportToDS( Double_t dS ); | |
248 | ||
249 | //* Get dS to a certain space point | |
250 | ||
251 | Double_t GetDStoPoint( const Double_t xyz[] ) const ; | |
252 | ||
253 | //* Get dS to other particle p (dSp for particle p also returned) | |
254 | ||
255 | void GetDStoParticle( const AliKFParticle &p, | |
256 | Double_t &DS, Double_t &DSp ) const ; | |
257 | ||
e7b09c95 | 258 | //* Get dS to other particle p in XY-plane |
259 | ||
260 | void GetDStoParticleXY( const AliKFParticleBase &p, | |
261 | Double_t &DS, Double_t &DSp ) const ; | |
262 | ||
f826d409 | 263 | //* |
264 | //* OTHER UTILITIES | |
265 | //* | |
266 | ||
267 | ||
268 | //* Calculate distance from another object [cm] | |
269 | ||
270 | Double_t GetDistanceFromVertex( const Double_t vtx[] ) const ; | |
271 | Double_t GetDistanceFromVertex( const AliKFParticle &Vtx ) const ; | |
706952f5 | 272 | Double_t GetDistanceFromVertex( const AliVVertex &Vtx ) const ; |
f826d409 | 273 | Double_t GetDistanceFromParticle( const AliKFParticle &p ) const ; |
274 | ||
275 | //* Calculate sqrt(Chi2/ndf) deviation from another object | |
276 | //* ( v = [xyz]-vertex, Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix ) | |
277 | ||
278 | Double_t GetDeviationFromVertex( const Double_t v[], const Double_t Cv[]=0 ) const ; | |
279 | Double_t GetDeviationFromVertex( const AliKFParticle &Vtx ) const ; | |
706952f5 | 280 | Double_t GetDeviationFromVertex( const AliVVertex &Vtx ) const ; |
f826d409 | 281 | Double_t GetDeviationFromParticle( const AliKFParticle &p ) const ; |
e7b09c95 | 282 | |
283 | //* Calculate distance from another object [cm] in XY-plane | |
284 | ||
446ce366 | 285 | Bool_t GetDistanceFromVertexXY( const Double_t vtx[], Double_t &val, Double_t &err ) const ; |
286 | Bool_t GetDistanceFromVertexXY( const Double_t vtx[], const Double_t Cv[], Double_t &val, Double_t &err ) const ; | |
287 | Bool_t GetDistanceFromVertexXY( const AliKFParticle &Vtx, Double_t &val, Double_t &err ) const ; | |
288 | Bool_t GetDistanceFromVertexXY( const AliVVertex &Vtx, Double_t &val, Double_t &err ) const ; | |
289 | ||
e7b09c95 | 290 | Double_t GetDistanceFromVertexXY( const Double_t vtx[] ) const ; |
291 | Double_t GetDistanceFromVertexXY( const AliKFParticle &Vtx ) const ; | |
706952f5 | 292 | Double_t GetDistanceFromVertexXY( const AliVVertex &Vtx ) const ; |
e7b09c95 | 293 | Double_t GetDistanceFromParticleXY( const AliKFParticle &p ) const ; |
294 | ||
295 | //* Calculate sqrt(Chi2/ndf) deviation from another object in XY plane | |
296 | //* ( v = [xyz]-vertex, Cv=[Cxx,Cxy,Cyy,Cxz,Cyz,Czz]-covariance matrix ) | |
297 | ||
298 | Double_t GetDeviationFromVertexXY( const Double_t v[], const Double_t Cv[]=0 ) const ; | |
299 | Double_t GetDeviationFromVertexXY( const AliKFParticle &Vtx ) const ; | |
706952f5 | 300 | Double_t GetDeviationFromVertexXY( const AliVVertex &Vtx ) const ; |
e7b09c95 | 301 | Double_t GetDeviationFromParticleXY( const AliKFParticle &p ) const ; |
302 | ||
303 | //* Calculate opennig angle between two particles | |
304 | ||
305 | Double_t GetAngle ( const AliKFParticle &p ) const ; | |
306 | Double_t GetAngleXY( const AliKFParticle &p ) const ; | |
307 | Double_t GetAngleRZ( const AliKFParticle &p ) const ; | |
f826d409 | 308 | |
309 | //* Subtract the particle from the vertex | |
310 | ||
311 | void SubtractFromVertex( AliKFParticle &v ) const ; | |
f826d409 | 312 | |
a65041d0 | 313 | //* Special method for creating gammas |
314 | ||
315 | void ConstructGamma( const AliKFParticle &daughter1, | |
316 | const AliKFParticle &daughter2 ); | |
91b25235 | 317 | |
318 | // * Pseudo Proper Time of decay = (r*pt) / |pt| * M/|pt| | |
319 | // @primVertex - primary vertex | |
320 | // @mass - mass of the mother particle (in the case of "Hb -> JPsi" it would be JPsi mass) | |
321 | // @*timeErr2 - squared error of the decay time. If timeErr2 = 0 it isn't calculated | |
322 | Double_t GetPseudoProperDecayTime( const AliKFParticle &primVertex, const Double_t& mass, Double_t* timeErr2 = 0 ) const; | |
323 | ||
f826d409 | 324 | protected: |
325 | ||
326 | //* | |
327 | //* INTERNAL STUFF | |
706952f5 | 328 | //* |
f826d409 | 329 | |
330 | //* Method to access ALICE field | |
331 | ||
4bbc290d | 332 | static Double_t GetFieldAlice(); |
f826d409 | 333 | |
334 | //* Other methods required by the abstract AliKFParticleBase class | |
335 | ||
336 | void GetFieldValue( const Double_t xyz[], Double_t B[] ) const ; | |
337 | void GetDStoParticle( const AliKFParticleBase &p, Double_t &DS, Double_t &DSp )const ; | |
338 | void Transport( Double_t dS, Double_t P[], Double_t C[] ) const ; | |
616ffc76 | 339 | static void GetExternalTrackParam( const AliKFParticleBase &p, Double_t &X, Double_t &Alpha, Double_t P[5] ) ; |
f826d409 | 340 | |
706952f5 | 341 | //void GetDStoParticleALICE( const AliKFParticleBase &p, Double_t &DS, Double_t &DS1 ) const; |
342 | ||
343 | ||
4bbc290d | 344 | private: |
345 | ||
346 | static Double_t fgBz; //* Bz compoment of the magnetic field | |
f826d409 | 347 | |
616ffc76 | 348 | ClassDef( AliKFParticle, 1 ); |
f826d409 | 349 | |
350 | }; | |
351 | ||
352 | ||
353 | ||
354 | //--------------------------------------------------------------------- | |
355 | // | |
356 | // Inline implementation of the AliKFParticle methods | |
357 | // | |
358 | //--------------------------------------------------------------------- | |
359 | ||
360 | ||
4bbc290d | 361 | inline void AliKFParticle::SetField( Double_t Bz ) |
362 | { | |
de0d0ceb | 363 | fgBz = Bz; |
4bbc290d | 364 | } |
f826d409 | 365 | |
616ffc76 | 366 | inline AliKFParticle::AliKFParticle( const AliKFParticle &d1, |
367 | const AliKFParticle &d2, | |
368 | const AliKFParticle &d3 ) | |
369 | { | |
370 | AliKFParticle mother; | |
371 | mother+= d1; | |
372 | mother+= d2; | |
373 | mother+= d3; | |
374 | *this = mother; | |
375 | } | |
376 | ||
377 | inline AliKFParticle::AliKFParticle( const AliKFParticle &d1, | |
378 | const AliKFParticle &d2, | |
379 | const AliKFParticle &d3, | |
380 | const AliKFParticle &d4 ) | |
381 | { | |
382 | AliKFParticle mother; | |
383 | mother+= d1; | |
384 | mother+= d2; | |
385 | mother+= d3; | |
386 | mother+= d4; | |
387 | *this = mother; | |
388 | } | |
389 | ||
390 | ||
f826d409 | 391 | inline void AliKFParticle::Initialize() |
392 | { | |
393 | AliKFParticleBase::Initialize(); | |
394 | } | |
395 | ||
396 | inline void AliKFParticle::SetVtxGuess( Double_t x, Double_t y, Double_t z ) | |
397 | { | |
398 | AliKFParticleBase::SetVtxGuess(x,y,z); | |
399 | } | |
400 | ||
401 | inline Double_t AliKFParticle::GetX () const | |
402 | { | |
403 | return AliKFParticleBase::GetX(); | |
404 | } | |
405 | ||
406 | inline Double_t AliKFParticle::GetY () const | |
407 | { | |
408 | return AliKFParticleBase::GetY(); | |
409 | } | |
410 | ||
411 | inline Double_t AliKFParticle::GetZ () const | |
412 | { | |
413 | return AliKFParticleBase::GetZ(); | |
414 | } | |
415 | ||
416 | inline Double_t AliKFParticle::GetPx () const | |
417 | { | |
418 | return AliKFParticleBase::GetPx(); | |
419 | } | |
420 | ||
421 | inline Double_t AliKFParticle::GetPy () const | |
422 | { | |
423 | return AliKFParticleBase::GetPy(); | |
424 | } | |
425 | ||
426 | inline Double_t AliKFParticle::GetPz () const | |
427 | { | |
428 | return AliKFParticleBase::GetPz(); | |
429 | } | |
430 | ||
431 | inline Double_t AliKFParticle::GetE () const | |
432 | { | |
433 | return AliKFParticleBase::GetE(); | |
434 | } | |
435 | ||
436 | inline Double_t AliKFParticle::GetS () const | |
437 | { | |
438 | return AliKFParticleBase::GetS(); | |
439 | } | |
440 | ||
441 | inline Int_t AliKFParticle::GetQ () const | |
442 | { | |
443 | return AliKFParticleBase::GetQ(); | |
444 | } | |
445 | ||
446 | inline Double_t AliKFParticle::GetChi2 () const | |
447 | { | |
448 | return AliKFParticleBase::GetChi2(); | |
449 | } | |
450 | ||
451 | inline Int_t AliKFParticle::GetNDF () const | |
452 | { | |
453 | return AliKFParticleBase::GetNDF(); | |
454 | } | |
455 | ||
456 | inline Double_t AliKFParticle::GetParameter ( int i ) const | |
457 | { | |
458 | return AliKFParticleBase::GetParameter(i); | |
459 | } | |
460 | ||
461 | inline Double_t AliKFParticle::GetCovariance( int i ) const | |
462 | { | |
463 | return AliKFParticleBase::GetCovariance(i); | |
464 | } | |
465 | ||
466 | inline Double_t AliKFParticle::GetCovariance( int i, int j ) const | |
467 | { | |
468 | return AliKFParticleBase::GetCovariance(i,j); | |
469 | } | |
470 | ||
e7b09c95 | 471 | |
706952f5 | 472 | inline Double_t AliKFParticle::GetP () const |
473 | { | |
474 | Double_t par, err; | |
475 | if( AliKFParticleBase::GetMomentum( par, err ) ) return 0; | |
476 | else return par; | |
477 | } | |
478 | ||
479 | inline Double_t AliKFParticle::GetPt () const | |
480 | { | |
481 | Double_t par, err; | |
482 | if( AliKFParticleBase::GetPt( par, err ) ) return 0; | |
483 | else return par; | |
484 | } | |
485 | ||
486 | inline Double_t AliKFParticle::GetEta () const | |
487 | { | |
488 | Double_t par, err; | |
489 | if( AliKFParticleBase::GetEta( par, err ) ) return 0; | |
490 | else return par; | |
491 | } | |
492 | ||
493 | inline Double_t AliKFParticle::GetPhi () const | |
494 | { | |
495 | Double_t par, err; | |
496 | if( AliKFParticleBase::GetPhi( par, err ) ) return 0; | |
497 | else return par; | |
498 | } | |
499 | ||
e7b09c95 | 500 | inline Double_t AliKFParticle::GetMomentum () const |
501 | { | |
502 | Double_t par, err; | |
503 | if( AliKFParticleBase::GetMomentum( par, err ) ) return 0; | |
504 | else return par; | |
505 | } | |
506 | ||
507 | inline Double_t AliKFParticle::GetMass () const | |
508 | { | |
509 | Double_t par, err; | |
510 | if( AliKFParticleBase::GetMass( par, err ) ) return 0; | |
511 | else return par; | |
512 | } | |
513 | ||
514 | inline Double_t AliKFParticle::GetDecayLength () const | |
515 | { | |
516 | Double_t par, err; | |
517 | if( AliKFParticleBase::GetDecayLength( par, err ) ) return 0; | |
518 | else return par; | |
519 | } | |
520 | ||
446ce366 | 521 | inline Double_t AliKFParticle::GetDecayLengthXY () const |
522 | { | |
523 | Double_t par, err; | |
524 | if( AliKFParticleBase::GetDecayLengthXY( par, err ) ) return 0; | |
525 | else return par; | |
526 | } | |
527 | ||
e7b09c95 | 528 | inline Double_t AliKFParticle::GetLifeTime () const |
529 | { | |
530 | Double_t par, err; | |
531 | if( AliKFParticleBase::GetLifeTime( par, err ) ) return 0; | |
532 | else return par; | |
533 | } | |
534 | ||
706952f5 | 535 | inline Double_t AliKFParticle::GetR () const |
536 | { | |
537 | Double_t par, err; | |
538 | if( AliKFParticleBase::GetR( par, err ) ) return 0; | |
539 | else return par; | |
540 | } | |
541 | ||
e7b09c95 | 542 | inline Double_t AliKFParticle::GetErrX () const |
543 | { | |
544 | return TMath::Sqrt(TMath::Abs( GetCovariance(0,0) )); | |
545 | } | |
546 | ||
547 | inline Double_t AliKFParticle::GetErrY () const | |
548 | { | |
549 | return TMath::Sqrt(TMath::Abs( GetCovariance(1,1) )); | |
550 | } | |
551 | ||
552 | inline Double_t AliKFParticle::GetErrZ () const | |
553 | { | |
554 | return TMath::Sqrt(TMath::Abs( GetCovariance(2,2) )); | |
555 | } | |
556 | ||
557 | inline Double_t AliKFParticle::GetErrPx () const | |
558 | { | |
559 | return TMath::Sqrt(TMath::Abs( GetCovariance(3,3) )); | |
560 | } | |
561 | ||
562 | inline Double_t AliKFParticle::GetErrPy () const | |
563 | { | |
564 | return TMath::Sqrt(TMath::Abs( GetCovariance(4,4) )); | |
565 | } | |
566 | ||
567 | inline Double_t AliKFParticle::GetErrPz () const | |
568 | { | |
569 | return TMath::Sqrt(TMath::Abs( GetCovariance(5,5) )); | |
570 | } | |
571 | ||
572 | inline Double_t AliKFParticle::GetErrE () const | |
573 | { | |
574 | return TMath::Sqrt(TMath::Abs( GetCovariance(6,6) )); | |
575 | } | |
576 | ||
577 | inline Double_t AliKFParticle::GetErrS () const | |
578 | { | |
579 | return TMath::Sqrt(TMath::Abs( GetCovariance(7,7) )); | |
580 | } | |
581 | ||
706952f5 | 582 | inline Double_t AliKFParticle::GetErrP () const |
583 | { | |
584 | Double_t par, err; | |
585 | if( AliKFParticleBase::GetMomentum( par, err ) ) return 1.e10; | |
586 | else return err; | |
587 | } | |
588 | ||
589 | inline Double_t AliKFParticle::GetErrPt () const | |
590 | { | |
591 | Double_t par, err; | |
592 | if( AliKFParticleBase::GetPt( par, err ) ) return 1.e10; | |
593 | else return err; | |
594 | } | |
595 | ||
596 | inline Double_t AliKFParticle::GetErrEta () const | |
597 | { | |
598 | Double_t par, err; | |
599 | if( AliKFParticleBase::GetEta( par, err ) ) return 1.e10; | |
600 | else return err; | |
601 | } | |
602 | ||
603 | inline Double_t AliKFParticle::GetErrPhi () const | |
604 | { | |
605 | Double_t par, err; | |
606 | if( AliKFParticleBase::GetPhi( par, err ) ) return 1.e10; | |
607 | else return err; | |
608 | } | |
609 | ||
e7b09c95 | 610 | inline Double_t AliKFParticle::GetErrMomentum () const |
611 | { | |
612 | Double_t par, err; | |
613 | if( AliKFParticleBase::GetMomentum( par, err ) ) return 1.e10; | |
614 | else return err; | |
615 | } | |
616 | ||
617 | inline Double_t AliKFParticle::GetErrMass () const | |
618 | { | |
619 | Double_t par, err; | |
620 | if( AliKFParticleBase::GetMass( par, err ) ) return 1.e10; | |
621 | else return err; | |
622 | } | |
623 | ||
624 | inline Double_t AliKFParticle::GetErrDecayLength () const | |
625 | { | |
626 | Double_t par, err; | |
627 | if( AliKFParticleBase::GetDecayLength( par, err ) ) return 1.e10; | |
628 | else return err; | |
629 | } | |
630 | ||
446ce366 | 631 | inline Double_t AliKFParticle::GetErrDecayLengthXY () const |
632 | { | |
633 | Double_t par, err; | |
634 | if( AliKFParticleBase::GetDecayLengthXY( par, err ) ) return 1.e10; | |
635 | else return err; | |
636 | } | |
637 | ||
e7b09c95 | 638 | inline Double_t AliKFParticle::GetErrLifeTime () const |
639 | { | |
640 | Double_t par, err; | |
641 | if( AliKFParticleBase::GetLifeTime( par, err ) ) return 1.e10; | |
642 | else return err; | |
643 | } | |
644 | ||
706952f5 | 645 | inline Double_t AliKFParticle::GetErrR () const |
646 | { | |
647 | Double_t par, err; | |
648 | if( AliKFParticleBase::GetR( par, err ) ) return 1.e10; | |
649 | else return err; | |
650 | } | |
651 | ||
652 | ||
653 | inline int AliKFParticle::GetP( Double_t &P, Double_t &SigmaP ) const | |
654 | { | |
655 | return AliKFParticleBase::GetMomentum( P, SigmaP ); | |
656 | } | |
657 | ||
658 | inline int AliKFParticle::GetPt( Double_t &Pt, Double_t &SigmaPt ) const | |
659 | { | |
660 | return AliKFParticleBase::GetPt( Pt, SigmaPt ); | |
661 | } | |
662 | ||
663 | inline int AliKFParticle::GetEta( Double_t &Eta, Double_t &SigmaEta ) const | |
664 | { | |
665 | return AliKFParticleBase::GetEta( Eta, SigmaEta ); | |
666 | } | |
667 | ||
668 | inline int AliKFParticle::GetPhi( Double_t &Phi, Double_t &SigmaPhi ) const | |
669 | { | |
670 | return AliKFParticleBase::GetPhi( Phi, SigmaPhi ); | |
671 | } | |
e7b09c95 | 672 | |
f826d409 | 673 | inline int AliKFParticle::GetMomentum( Double_t &P, Double_t &SigmaP ) const |
674 | { | |
675 | return AliKFParticleBase::GetMomentum( P, SigmaP ); | |
676 | } | |
677 | ||
678 | inline int AliKFParticle::GetMass( Double_t &M, Double_t &SigmaM ) const | |
679 | { | |
680 | return AliKFParticleBase::GetMass( M, SigmaM ); | |
681 | } | |
682 | ||
683 | inline int AliKFParticle::GetDecayLength( Double_t &L, Double_t &SigmaL ) const | |
684 | { | |
685 | return AliKFParticleBase::GetDecayLength( L, SigmaL ); | |
686 | } | |
687 | ||
446ce366 | 688 | inline int AliKFParticle::GetDecayLengthXY( Double_t &L, Double_t &SigmaL ) const |
689 | { | |
690 | return AliKFParticleBase::GetDecayLengthXY( L, SigmaL ); | |
691 | } | |
692 | ||
f826d409 | 693 | inline int AliKFParticle::GetLifeTime( Double_t &T, Double_t &SigmaT ) const |
694 | { | |
695 | return AliKFParticleBase::GetLifeTime( T, SigmaT ); | |
696 | } | |
697 | ||
706952f5 | 698 | inline int AliKFParticle::GetR( Double_t &R, Double_t &SigmaR ) const |
699 | { | |
700 | return AliKFParticleBase::GetR( R, SigmaR ); | |
701 | } | |
702 | ||
f826d409 | 703 | inline Double_t & AliKFParticle::X() |
704 | { | |
705 | return AliKFParticleBase::X(); | |
706 | } | |
707 | ||
708 | inline Double_t & AliKFParticle::Y() | |
709 | { | |
710 | return AliKFParticleBase::Y(); | |
711 | } | |
712 | ||
713 | inline Double_t & AliKFParticle::Z() | |
714 | { | |
715 | return AliKFParticleBase::Z(); | |
716 | } | |
717 | ||
718 | inline Double_t & AliKFParticle::Px() | |
719 | { | |
720 | return AliKFParticleBase::Px(); | |
721 | } | |
722 | ||
723 | inline Double_t & AliKFParticle::Py() | |
724 | { | |
725 | return AliKFParticleBase::Py(); | |
726 | } | |
727 | ||
728 | inline Double_t & AliKFParticle::Pz() | |
729 | { | |
730 | return AliKFParticleBase::Pz(); | |
731 | } | |
732 | ||
733 | inline Double_t & AliKFParticle::E() | |
734 | { | |
735 | return AliKFParticleBase::E(); | |
736 | } | |
737 | ||
738 | inline Double_t & AliKFParticle::S() | |
739 | { | |
740 | return AliKFParticleBase::S(); | |
741 | } | |
742 | ||
743 | inline Int_t & AliKFParticle::Q() | |
744 | { | |
745 | return AliKFParticleBase::Q(); | |
746 | } | |
747 | ||
748 | inline Double_t & AliKFParticle::Chi2() | |
749 | { | |
750 | return AliKFParticleBase::Chi2(); | |
751 | } | |
752 | ||
753 | inline Int_t & AliKFParticle::NDF() | |
754 | { | |
755 | return AliKFParticleBase::NDF(); | |
756 | } | |
757 | ||
758 | inline Double_t & AliKFParticle::Parameter ( int i ) | |
759 | { | |
760 | return AliKFParticleBase::Parameter(i); | |
761 | } | |
762 | ||
763 | inline Double_t & AliKFParticle::Covariance( int i ) | |
764 | { | |
765 | return AliKFParticleBase::Covariance(i); | |
766 | } | |
767 | ||
768 | inline Double_t & AliKFParticle::Covariance( int i, int j ) | |
769 | { | |
770 | return AliKFParticleBase::Covariance(i,j); | |
771 | } | |
772 | ||
bb4bc052 | 773 | inline Double_t * AliKFParticle::Parameters () |
774 | { | |
775 | return fP; | |
776 | } | |
777 | ||
778 | inline Double_t * AliKFParticle::CovarianceMatrix() | |
779 | { | |
780 | return fC; | |
781 | } | |
782 | ||
f826d409 | 783 | |
784 | inline void AliKFParticle::operator +=( const AliKFParticle &Daughter ) | |
785 | { | |
786 | AliKFParticleBase::operator +=( Daughter ); | |
787 | } | |
788 | ||
f826d409 | 789 | |
790 | inline void AliKFParticle::AddDaughter( const AliKFParticle &Daughter ) | |
791 | { | |
792 | AliKFParticleBase::AddDaughter( Daughter ); | |
793 | } | |
794 | ||
795 | inline void AliKFParticle::SetProductionVertex( const AliKFParticle &Vtx ) | |
796 | { | |
797 | AliKFParticleBase::SetProductionVertex( Vtx ); | |
798 | } | |
799 | ||
e7b09c95 | 800 | inline void AliKFParticle::SetMassConstraint( Double_t Mass, Double_t SigmaMass ) |
f826d409 | 801 | { |
e7b09c95 | 802 | AliKFParticleBase::SetMassConstraint( Mass, SigmaMass ); |
f826d409 | 803 | } |
e7b09c95 | 804 | |
805 | inline void AliKFParticle::SetNoDecayLength() | |
806 | { | |
807 | AliKFParticleBase::SetNoDecayLength(); | |
808 | } | |
809 | ||
f826d409 | 810 | inline void AliKFParticle::Construct( const AliKFParticle *vDaughters[], int NDaughters, |
706952f5 | 811 | const AliKFParticle *ProdVtx, Double_t Mass, Bool_t IsConstrained ) |
f826d409 | 812 | { |
813 | AliKFParticleBase::Construct( ( const AliKFParticleBase**)vDaughters, NDaughters, | |
706952f5 | 814 | ( const AliKFParticleBase*)ProdVtx, Mass, IsConstrained ); |
f826d409 | 815 | } |
816 | ||
817 | inline void AliKFParticle::TransportToDecayVertex() | |
818 | { | |
819 | AliKFParticleBase::TransportToDecayVertex(); | |
820 | } | |
821 | ||
822 | inline void AliKFParticle::TransportToProductionVertex() | |
823 | { | |
824 | AliKFParticleBase::TransportToProductionVertex(); | |
825 | } | |
826 | ||
827 | inline void AliKFParticle::TransportToPoint( const Double_t xyz[] ) | |
828 | { | |
829 | TransportToDS( GetDStoPoint(xyz) ); | |
830 | } | |
831 | ||
706952f5 | 832 | inline void AliKFParticle::TransportToVertex( const AliVVertex &v ) |
f826d409 | 833 | { |
4bbc290d | 834 | TransportToPoint( AliKFParticle(v).fP ); |
f826d409 | 835 | } |
836 | ||
837 | inline void AliKFParticle::TransportToParticle( const AliKFParticle &p ) | |
838 | { | |
839 | Double_t dS, dSp; | |
840 | GetDStoParticle( p, dS, dSp ); | |
841 | TransportToDS( dS ); | |
842 | } | |
843 | ||
844 | inline void AliKFParticle::TransportToDS( Double_t dS ) | |
845 | { | |
846 | AliKFParticleBase::TransportToDS( dS ); | |
847 | } | |
848 | ||
849 | inline Double_t AliKFParticle::GetDStoPoint( const Double_t xyz[] ) const | |
850 | { | |
851 | return AliKFParticleBase::GetDStoPointBz( GetFieldAlice(), xyz ); | |
852 | } | |
4bbc290d | 853 | |
f826d409 | 854 | |
855 | inline void AliKFParticle::GetDStoParticle( const AliKFParticle &p, | |
4bbc290d | 856 | Double_t &DS, Double_t &DSp ) const |
f826d409 | 857 | { |
e7b09c95 | 858 | GetDStoParticleXY( p, DS, DSp ); |
f826d409 | 859 | } |
860 | ||
861 | ||
862 | inline Double_t AliKFParticle::GetDistanceFromVertex( const Double_t vtx[] ) const | |
863 | { | |
864 | return AliKFParticleBase::GetDistanceFromVertex( vtx ); | |
865 | } | |
866 | ||
867 | inline Double_t AliKFParticle::GetDeviationFromVertex( const Double_t v[], | |
868 | const Double_t Cv[] ) const | |
869 | { | |
870 | return AliKFParticleBase::GetDeviationFromVertex( v, Cv); | |
871 | } | |
872 | ||
873 | inline Double_t AliKFParticle::GetDistanceFromVertex( const AliKFParticle &Vtx ) const | |
874 | { | |
875 | return AliKFParticleBase::GetDistanceFromVertex( Vtx ); | |
876 | } | |
877 | ||
878 | inline Double_t AliKFParticle::GetDeviationFromVertex( const AliKFParticle &Vtx ) const | |
879 | { | |
880 | return AliKFParticleBase::GetDeviationFromVertex( Vtx ); | |
881 | } | |
882 | ||
706952f5 | 883 | inline Double_t AliKFParticle::GetDistanceFromVertex( const AliVVertex &Vtx ) const |
f826d409 | 884 | { |
4bbc290d | 885 | return GetDistanceFromVertex( AliKFParticle(Vtx) ); |
f826d409 | 886 | } |
887 | ||
706952f5 | 888 | inline Double_t AliKFParticle::GetDeviationFromVertex( const AliVVertex &Vtx ) const |
f826d409 | 889 | { |
4bbc290d | 890 | return GetDeviationFromVertex( AliKFParticle(Vtx) ); |
f826d409 | 891 | } |
446ce366 | 892 | |
f826d409 | 893 | inline Double_t AliKFParticle::GetDistanceFromParticle( const AliKFParticle &p ) const |
894 | { | |
616ffc76 | 895 | return AliKFParticleBase::GetDistanceFromParticle( p ); |
f826d409 | 896 | } |
897 | ||
898 | inline Double_t AliKFParticle::GetDeviationFromParticle( const AliKFParticle &p ) const | |
899 | { | |
616ffc76 | 900 | return AliKFParticleBase::GetDeviationFromParticle( p ); |
f826d409 | 901 | } |
902 | ||
903 | inline void AliKFParticle::SubtractFromVertex( AliKFParticle &v ) const | |
904 | { | |
de0d0ceb | 905 | AliKFParticleBase::SubtractFromVertex( v ); |
f826d409 | 906 | } |
907 | ||
4bbc290d | 908 | inline Double_t AliKFParticle::GetFieldAlice() |
f826d409 | 909 | { |
4bbc290d | 910 | return fgBz; |
f826d409 | 911 | } |
912 | ||
913 | inline void AliKFParticle::GetFieldValue( const Double_t * /*xyz*/, Double_t B[] ) const | |
914 | { | |
915 | B[0] = B[1] = 0; | |
916 | B[2] = GetFieldAlice(); | |
917 | } | |
918 | ||
919 | inline void AliKFParticle::GetDStoParticle( const AliKFParticleBase &p, | |
4bbc290d | 920 | Double_t &DS, Double_t &DSp )const |
f826d409 | 921 | { |
e7b09c95 | 922 | GetDStoParticleXY( p, DS, DSp ); |
923 | } | |
924 | ||
925 | inline void AliKFParticle::GetDStoParticleXY( const AliKFParticleBase &p, | |
926 | Double_t &DS, Double_t &DSp ) const | |
927 | { | |
928 | AliKFParticleBase::GetDStoParticleBz( GetFieldAlice(), p, DS, DSp ) ; | |
706952f5 | 929 | //GetDStoParticleALICE( p, DS, DSp ) ; |
f826d409 | 930 | } |
931 | ||
932 | inline void AliKFParticle::Transport( Double_t dS, Double_t P[], Double_t C[] ) const | |
933 | { | |
934 | AliKFParticleBase::TransportBz( GetFieldAlice(), dS, P, C ); | |
935 | } | |
936 | ||
a65041d0 | 937 | inline void AliKFParticle::ConstructGamma( const AliKFParticle &daughter1, |
938 | const AliKFParticle &daughter2 ) | |
939 | { | |
940 | AliKFParticleBase::ConstructGammaBz( daughter1, daughter2, GetFieldAlice() ); | |
941 | } | |
942 | ||
f826d409 | 943 | #endif |