AliRoot Core  3dc7879 (3dc7879)
AliESDv0.cxx
Go to the documentation of this file.
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  * *
4  * Author: The ALICE Off-line Project. *
5  * Contributors are mentioned in the code where appropriate. *
6  * *
7  * Permission to use, copy, modify and distribute this software and its *
8  * documentation strictly for non-commercial purposes is hereby granted *
9  * without fee, provided that the above copyright notice appears in all *
10  * copies and that both the copyright notice and this permission notice *
11  * appear in the supporting documentation. The authors make no claims *
12  * about the suitability of this software for any purpose. It is *
13  * provided "as is" without express or implied warranty. *
14  **************************************************************************/
15 
16 /* $Id$ */
17 
18 //-------------------------------------------------------------------------
19 // Implementation of the ESD V0 vertex class
20 // This class is part of the Event Data Summary
21 // set of classes and contains information about
22 // V0 kind vertexes generated by a neutral particle
23 // Origin: Iouri Belikov, IReS, Strasbourg, Jouri.Belikov@cern.ch
24 // Modified by: Marian Ivanov, CERN, Marian.Ivanov@cern.ch
25 // and Boris Hippolyte,IPHC, hippolyt@in2p3.fr
26 //-------------------------------------------------------------------------
27 
28 #include <TMath.h>
29 #include <TDatabasePDG.h>
30 #include <TParticlePDG.h>
31 #include <TVector3.h>
32 #include <TMatrixD.h>
33 
34 #include "AliLog.h"
35 #include "AliESDv0.h"
36 #include "AliESDV0Params.h"
37 #include "AliKFParticle.h"
38 #include "AliKFVertex.h"
39 #include "AliESDVertex.h"
40 
41 ClassImp(AliESDv0)
42 
43 const AliESDV0Params AliESDv0::fgkParams;
44 
46  AliVParticle(),
47  fParamN(),
48  fParamP(),
49  fEffMass(TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass()),
50  fDcaV0Daughters(0),
51  fChi2V0(0.),
52  fRr(0),
53  fDistSigma(0),
54  fChi2Before(0),
55  fChi2After(0),
56  fPointAngleFi(0),
57  fPointAngleTh(0),
58  fPointAngle(0),
59  fPdgCode(kK0Short),
60  fNidx(0),
61  fPidx(0),
62  fStatus(0),
63  fNBefore(0),
64  fNAfter(0),
65  fOnFlyStatus(kFALSE)
66 {
67  //--------------------------------------------------------------------
68  // Default constructor (K0s)
69  //--------------------------------------------------------------------
70 
71  for (Int_t i=0; i<3; i++) {
72  fPos[i] = 0.;
73  fNmom[i] = 0.;
74  fPmom[i] = 0.;
75  }
76 
77  for (Int_t i=0; i<6; i++) {
78  fPosCov[i]= 0.;
79  }
80 
81  for (Int_t i=0;i<6;i++){fClusters[0][i]=0; fClusters[1][i]=0;}
82  fNormDCAPrim[0]=fNormDCAPrim[1]=0;
83  for (Int_t i=0;i<3;i++){fAngle[i]=0;}
84  for (Int_t i=0;i<4;i++){fCausality[i]=0;}
85 }
86 
88  AliVParticle(v0),
89  fParamN(v0.fParamN),
90  fParamP(v0.fParamP),
91  fEffMass(v0.fEffMass),
92  fDcaV0Daughters(v0.fDcaV0Daughters),
93  fChi2V0(v0.fChi2V0),
94  fRr(v0.fRr),
95  fDistSigma(v0.fDistSigma),
96  fChi2Before(v0.fChi2Before),
97  fChi2After(v0.fChi2After),
98  fPointAngleFi(v0.fPointAngleFi),
99  fPointAngleTh(v0.fPointAngleTh),
100  fPointAngle(v0.fPointAngle),
101  fPdgCode(v0.fPdgCode),
102  fNidx(v0.fNidx),
103  fPidx(v0.fPidx),
104  fStatus(v0.fStatus),
105  fNBefore(v0.fNBefore),
106  fNAfter(v0.fNAfter),
107  fOnFlyStatus(v0.fOnFlyStatus)
108 {
109  //--------------------------------------------------------------------
110  // The copy constructor
111  //--------------------------------------------------------------------
112 
113  Short_t cN=fParamN.Charge(), cP=fParamP.Charge();
114  Bool_t swp = cN>0 && (cN != cP);
115  if (swp) {
116  fParamN = v0.fParamP;
117  fParamP = v0.fParamN;
118  fNidx=v0.fPidx;
119  fPidx=v0.fNidx;
120  }
121 
122  for (int i=0; i<3; i++) {
123  fPos[i] = v0.fPos[i];
124  fNmom[i] = v0.fNmom[i];
125  fPmom[i] = v0.fPmom[i];
126  }
127  for (int i=0; i<6; i++) {
128  fPosCov[i] = v0.fPosCov[i];
129  }
130 
131  for (Int_t i=0; i<2; i++) {
132  fNormDCAPrim[i]=v0.fNormDCAPrim[i];
133  }
134  if (swp) {
135  for (Int_t i=0;i<6;i++){
136  fClusters[0][i]=v0.fClusters[1][i];
137  fClusters[1][i]=v0.fClusters[0][i];
138  }
139  }
140  else {
141  for (Int_t i=0;i<6;i++){
142  fClusters[0][i]=v0.fClusters[0][i];
143  fClusters[1][i]=v0.fClusters[1][i];
144  }
145  }
146  for (Int_t i=0;i<3;i++){
147  fAngle[i]=v0.fAngle[i];
148  }
149  for (Int_t i=0;i<4;i++){fCausality[i]=v0.fCausality[i];}
150 }
151 
152 
154  const AliExternalTrackParam &t2, Int_t i2) :
155  AliVParticle(),
156  fParamN(t1),
157  fParamP(t2),
158  fEffMass(TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass()),
159  fDcaV0Daughters(0),
160  fChi2V0(0.),
161  fRr(0),
162  fDistSigma(0),
163  fChi2Before(0),
164  fChi2After(0),
165  fPointAngleFi(0),
166  fPointAngleTh(0),
167  fPointAngle(0),
168  fPdgCode(kK0Short),
169  fNidx(i1),
170  fPidx(i2),
171  fStatus(0),
172  fNBefore(0),
173  fNAfter(0),
174  fOnFlyStatus(kFALSE)
175 {
176  //--------------------------------------------------------------------
177  // Main constructor (K0s)
178  //--------------------------------------------------------------------
179 
180  //Make sure the daughters are ordered (needed for the on-the-fly V0s)
181  Short_t cN=t1.Charge(), cP=t2.Charge();
182  Bool_t swp = cN>0 && (cN != cP);
183  if (swp) {
184  fParamN = t2;
185  fParamP = t1;
186  fNidx = i2;
187  fPidx = i1;
188  }
189 
190  for (Int_t i=0; i<6; i++) {
191  fPosCov[i]= 0.;
192  }
193 
194  //Trivial estimation of the vertex parameters
195  Double_t alpha=fParamN.GetAlpha(), cs=TMath::Cos(alpha), sn=TMath::Sin(alpha);
196  Double_t tmp[3];
197  fParamN.GetPxPyPz(tmp);
198  Double_t px1=tmp[0], py1=tmp[1], pz1=tmp[2];
199  fParamN.GetXYZ(tmp);
200  Double_t x1=tmp[0], y1=tmp[1], z1=tmp[2];
201  const Double_t ss=0.0005*0.0005;//a kind of a residual misalignment precision
202  Double_t sx1=sn*sn*fParamN.GetSigmaY2()+ss, sy1=cs*cs*fParamN.GetSigmaY2()+ss;
203 
204 
205  alpha=fParamP.GetAlpha(); cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
206  fParamP.GetPxPyPz(tmp);
207  Double_t px2=tmp[0], py2=tmp[1], pz2=tmp[2];
208  fParamP.GetXYZ(tmp);
209  Double_t x2=tmp[0], y2=tmp[1], z2=tmp[2];
210  Double_t sx2=sn*sn*fParamP.GetSigmaY2()+ss, sy2=cs*cs*fParamP.GetSigmaY2()+ss;
211 
212  Double_t sz1=fParamN.GetSigmaZ2(), sz2=fParamP.GetSigmaZ2();
213  Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
214  Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
215  Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
216  fPos[0]=wx1*x1 + wx2*x2; fPos[1]=wy1*y1 + wy2*y2; fPos[2]=wz1*z1 + wz2*z2;
217 
218  //fPos[0]=0.5*(x1+x2); fPos[1]=0.5*(y1+y2); fPos[2]=0.5*(z1+z2);
219  fNmom[0]=px1; fNmom[1]=py1; fNmom[2]=pz1;
220  fPmom[0]=px2; fPmom[1]=py2; fPmom[2]=pz2;
221 
222  for (Int_t i=0;i<6;i++){fClusters[0][i]=0; fClusters[1][i]=0;}
223  fNormDCAPrim[0]=fNormDCAPrim[1]=0;
224  for (Int_t i=0;i<3;i++){fAngle[i]=0;}
225  for (Int_t i=0;i<4;i++){fCausality[i]=0;}
226 }
227 
228 static Bool_t GetWeight(TMatrixD &w, const AliExternalTrackParam &t) {
229  //
230  // Returns the global weight matrix w = Transpose[G2P]*Inverse[Cpar]*G2P ,
231  // where the matrix Cpar is the transverse part of the t covariance
232  // in "parallel" system (i.e. the system with X axis parallel to momentum).
233  // The matrix G2P performs the transformation global -> "parallel".
234  //
235  Double_t phi=t.GetAlpha() + TMath::ASin(t.GetSnp());
236  Double_t sp=TMath::Sin(phi);
237  Double_t cp=TMath::Cos(phi);
238 
239  Double_t tgl=t.GetTgl();
240  Double_t cl=1/TMath::Sqrt(1.+ tgl*tgl);
241  Double_t sl=tgl*cl;
242 
243  TMatrixD g2p(3,3); //global --> parallel
244  g2p(0,0)= cp*cl; g2p(0,1)= sp*cl; g2p(0,2)=sl;
245  g2p(1,0)=-sp; g2p(1,1)= cp; g2p(1,2)=0.;
246  g2p(2,0)=-sl*cp; g2p(2,1)=-sl*sp; g2p(2,2)=cl;
247 
248  Double_t alpha=t.GetAlpha();
249  Double_t c=TMath::Cos(alpha), s=TMath::Sin(alpha);
250  TMatrixD l2g(3,3); //local --> global
251  l2g(0,0)= c; l2g(0,1)=-s; l2g(0,2)= 0;
252  l2g(1,0)= s; l2g(1,1)= c; l2g(1,2)= 0;
253  l2g(2,0)= 0; l2g(2,1)= 0; l2g(2,2)= 1;
254 
255  Double_t sy2=t.GetSigmaY2(), syz=t.GetSigmaZY(), sz2=t.GetSigmaZ2();
256  TMatrixD cvl(3,3); //local covariance
257  cvl(0,0)=0; cvl(0,1)=0; cvl(0,2)=0;
258  cvl(1,0)=0; cvl(1,1)=sy2; cvl(1,2)=syz;
259  cvl(2,0)=0; cvl(2,1)=syz; cvl(2,2)=sz2;
260 
261  TMatrixD l2p(g2p, TMatrixD::kMult, l2g);
262  TMatrixD cvp(3,3); //parallel covariance
263  cvp=l2p*cvl*TMatrixD(TMatrixD::kTransposed,l2p);
264 
265  Double_t det=cvp(1,1)*cvp(2,2) - cvp(1,2)*cvp(2,1);
266  if (TMath::Abs(det)<kAlmost0) return kFALSE;
267 
268  const Double_t m=100*100; //A large uncertainty in the momentum direction
269  const Double_t eps=1/m;
270  TMatrixD u(3,3); //Inverse of the transverse part of the parallel covariance
271  u(0,0)=eps; u(0,1)=0; u(0,2)=0;
272  u(1,0)=0; u(1,1)= cvp(2,2)/det; u(1,2)=-cvp(2,1)/det;
273  u(2,0)=0; u(2,1)=-cvp(1,2)/det; u(2,2)= cvp(1,1)/det;
274 
275  w=TMatrixD(TMatrixD::kTransposed,g2p)*u*g2p;
276 
277  return kTRUE;
278 }
279 
281  //--------------------------------------------------------------------
282  // Refit this vertex
283  //--------------------------------------------------------------------
284  fStatus=0;
285 
286  //Save the daughters' momenta
289 
290  //Trivial estimation of the V0 vertex parameters
291  Double_t r1[3]; fParamN.GetXYZ(r1);
292  Double_t r2[3]; fParamP.GetXYZ(r2);
293  for (Int_t i=0; i<3; i++) fPos[i]=0.5*(r1[i]+r2[i]);
294  fPosCov[1]=fPosCov[3]=fPosCov[4]=0.;
295  fPosCov[0]=(fPos[0]-r1[0])*(fPos[0]-r1[0])/12.;
296  fPosCov[2]=(fPos[1]-r1[1])*(fPos[1]-r1[1])/12.;
297  fPosCov[5]=(fPos[2]-r1[2])*(fPos[2]-r1[2])/12.;
298  fChi2V0=12.;
299 
300 
301  //Try to improve the V0 vertex parameters
302  TMatrixD w1(3,3);
303  if (!GetWeight(w1,fParamN)) return fStatus;
304  TMatrixD w2(3,3);
305  if (!GetWeight(w2,fParamP)) return fStatus;
306  TMatrixD cv(w1); cv+=w2;
307  cv.Invert();
308  if (!cv.IsValid()) return fStatus;
309 
310  //Covariance of the V0 vertex
311  fPosCov[0]=cv(0,0);
312  fPosCov[1]=cv(1,0); fPosCov[2]=cv(1,1);
313  fPosCov[3]=cv(2,0); fPosCov[4]=cv(2,1); fPosCov[5]=cv(2,2);
314 
315  //Position of the V0 vertex
316  TMatrixD cw1(cv,TMatrixD::kMult,w1);
317  for (Int_t i=0; i<3; i++) {
318  fPos[i]=r2[i];
319  for (Int_t j=0; j<3; j++) fPos[i] += cw1(i,j)*(r1[j] - r2[j]);
320  }
321 
322  //Chi2 of the V0 vertex
323  fChi2V0=0.;
324  Double_t res1[3]={r1[0]-fPos[0],r1[1]-fPos[1],r1[2]-fPos[2]};
325  Double_t res2[3]={r2[0]-fPos[0],r2[1]-fPos[1],r2[2]-fPos[2]};
326  for (Int_t i=0; i<3; i++)
327  for (Int_t j=0; j<3; j++)
328  fChi2V0 += res1[i]*res1[j]*w1(i,j) + res2[i]*res2[j]*w2(i,j);
329 
330  //Successful fit
331  fStatus=1;
332 
333  return fStatus;
334 }
335 
337 {
338  //--------------------------------------------------------------------
339  // The assignment operator
340  //--------------------------------------------------------------------
341 
342  if(this==&v0)return *this;
344  fParamN = v0.fParamN;
345  fParamP = v0.fParamP;
346  fEffMass = v0.fEffMass;
348  fChi2V0 = v0.fChi2V0;
349  fRr = v0.fRr;
350  fDistSigma = v0.fDistSigma;
352  fChi2After = v0.fChi2After;
356  fPdgCode = v0.fPdgCode;
357  fNidx = v0.fNidx;
358  fPidx = v0.fPidx;
359  fStatus = v0.fStatus;
360  fNBefore = v0.fNBefore;
361  fNAfter = v0.fNAfter;
363 
364  for (int i=0; i<3; i++) {
365  fPos[i] = v0.fPos[i];
366  fNmom[i] = v0.fNmom[i];
367  fPmom[i] = v0.fPmom[i];
368  }
369  for (int i=0; i<6; i++) {
370  fPosCov[i] = v0.fPosCov[i];
371  }
372  for (Int_t i=0; i<2; i++) {
373  fNormDCAPrim[i]=v0.fNormDCAPrim[i];
374  }
375  for (Int_t i=0;i<6;i++){
376  fClusters[0][i]=v0.fClusters[0][i];
377  fClusters[1][i]=v0.fClusters[1][i];
378  }
379  for (Int_t i=0;i<3;i++){
380  fAngle[i]=v0.fAngle[i];
381  }
382  for (Int_t i=0;i<4;i++){fCausality[i]=v0.fCausality[i];}
383 
384  return *this;
385 }
386 
387 void AliESDv0::Copy(TObject& obj) const {
388 
389  // this overwrites the virtual TOBject::Copy()
390  // to allow run time copying without casting
391  // in AliESDEvent
392 
393  if(this==&obj)return;
394  AliESDv0 *robj = dynamic_cast<AliESDv0*>(&obj);
395  if(!robj)return; // not an aliesesv0
396  *robj = *this;
397 }
398 
400  //--------------------------------------------------------------------
401  // Empty destructor
402  //--------------------------------------------------------------------
403 }
404 
405 // Start with AliVParticle functions
406 Double_t AliESDv0::E() const {
407  //--------------------------------------------------------------------
408  // This gives the energy assuming the ChangeMassHypothesis was called
409  //--------------------------------------------------------------------
410  return E(fPdgCode);
411 }
412 
413 Double_t AliESDv0::Y() const {
414  //--------------------------------------------------------------------
415  // This gives the energy assuming the ChangeMassHypothesis was called
416  //--------------------------------------------------------------------
417  return Y(fPdgCode);
418 }
419 
420 // Then extend AliVParticle functions
421 Double_t AliESDv0::E(Int_t pdg) const {
422  //--------------------------------------------------------------------
423  // This gives the energy with the particle hypothesis as argument
424  //--------------------------------------------------------------------
425  Double_t mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
426  return TMath::Sqrt(mass*mass+P()*P());
427 }
428 
429 Double_t AliESDv0::Y(Int_t pdg) const {
430  //--------------------------------------------------------------------
431  // This gives the rapidity with the particle hypothesis as argument
432  //--------------------------------------------------------------------
433  return 0.5*TMath::Log((E(pdg)+Pz())/(E(pdg)-Pz()+1.e-13));
434 }
435 
436 // Now the functions for analysis consistency
437 Double_t AliESDv0::RapK0Short() const {
438  //--------------------------------------------------------------------
439  // This gives the pseudorapidity assuming a K0s particle
440  //--------------------------------------------------------------------
441  return Y(kK0Short);
442 }
443 
444 Double_t AliESDv0::RapLambda() const {
445  //--------------------------------------------------------------------
446  // This gives the pseudorapidity assuming a (Anti) Lambda particle
447  //--------------------------------------------------------------------
448  return Y(kLambda0);
449 }
450 
451 Double_t AliESDv0::AlphaV0() const {
452  //--------------------------------------------------------------------
453  // This gives the Armenteros-Podolanski alpha
454  //--------------------------------------------------------------------
455  TVector3 momNeg(fNmom[0],fNmom[1],fNmom[2]);
456  TVector3 momPos(fPmom[0],fPmom[1],fPmom[2]);
457  TVector3 momTot(Px(),Py(),Pz());
458 
459  Double_t lQlNeg = momNeg.Dot(momTot)/momTot.Mag();
460  Double_t lQlPos = momPos.Dot(momTot)/momTot.Mag();
461 
462  //return 1.-2./(1.+lQlNeg/lQlPos);
463  return (lQlPos - lQlNeg)/(lQlPos + lQlNeg);
464 }
465 
466 Double_t AliESDv0::PtArmV0() const {
467  //--------------------------------------------------------------------
468  // This gives the Armenteros-Podolanski ptarm
469  //--------------------------------------------------------------------
470  TVector3 momNeg(fNmom[0],fNmom[1],fNmom[2]);
471  TVector3 momTot(Px(),Py(),Pz());
472 
473  return momNeg.Perp(momTot);
474 }
475 
476 // Eventually the older functions
477 Double_t AliESDv0::ChangeMassHypothesis(Int_t code) {
478  //--------------------------------------------------------------------
479  // This function changes the mass hypothesis for this V0
480  // and returns the "kinematical quality" of this hypothesis
481  //--------------------------------------------------------------------
482  static
483  Double_t piMass=TDatabasePDG::Instance()->GetParticle(kPiPlus)->Mass();
484  static
485  Double_t prMass=TDatabasePDG::Instance()->GetParticle(kProton)->Mass();
486  static
487  Double_t k0Mass=TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass();
488  static
489  Double_t l0Mass=TDatabasePDG::Instance()->GetParticle(kLambda0)->Mass();
490 
491  Double_t nmass=piMass, pmass=piMass, mass=k0Mass, ps=0.206;
492 
493  fPdgCode=code;
494 
495  switch (code) {
496  case kLambda0:
497  nmass=piMass; pmass=prMass; mass=l0Mass; ps=0.101; break;
498  case kLambda0Bar:
499  pmass=piMass; nmass=prMass; mass=l0Mass; ps=0.101; break;
500  case kK0Short:
501  break;
502  default:
503  AliError("invalide PDG code ! Assuming K0s...");
504  fPdgCode=kK0Short;
505  break;
506  }
507 
508  Double_t pxn=fNmom[0], pyn=fNmom[1], pzn=fNmom[2];
509  Double_t pxp=fPmom[0], pyp=fPmom[1], pzp=fPmom[2];
510 
511  Double_t en=TMath::Sqrt(nmass*nmass + pxn*pxn + pyn*pyn + pzn*pzn);
512  Double_t ep=TMath::Sqrt(pmass*pmass + pxp*pxp + pyp*pyp + pzp*pzp);
513  Double_t pxl=pxn+pxp, pyl=pyn+pyp, pzl=pzn+pzp;
514  Double_t pl=TMath::Sqrt(pxl*pxl + pyl*pyl + pzl*pzl);
515 
516  fEffMass=TMath::Sqrt((en+ep)*(en+ep)-pl*pl);
517 
518  Double_t beta=pl/(en+ep);
519  Double_t pln=(pxn*pxl + pyn*pyl + pzn*pzl)/pl;
520  Double_t plp=(pxp*pxl + pyp*pyl + pzp*pzl)/pl;
521 
522  Double_t pt2=pxp*pxp + pyp*pyp + pzp*pzp - plp*plp;
523 
524  Double_t a=(plp-pln)/(plp+pln);
525  a -= (pmass*pmass-nmass*nmass)/(mass*mass);
526  a = 0.25*beta*beta*mass*mass*a*a + pt2;
527 
528  return (a - ps*ps);
529 
530 }
531 
532 void AliESDv0::GetPxPyPz(Double_t &px, Double_t &py, Double_t &pz) const {
533  //--------------------------------------------------------------------
534  // This function returns V0's momentum (global)
535  //--------------------------------------------------------------------
536  px=fNmom[0]+fPmom[0];
537  py=fNmom[1]+fPmom[1];
538  pz=fNmom[2]+fPmom[2];
539 }
540 
541 void AliESDv0::GetXYZ(Double_t &x, Double_t &y, Double_t &z) const {
542  //--------------------------------------------------------------------
543  // This function returns V0's position (global)
544  //--------------------------------------------------------------------
545  x=fPos[0];
546  y=fPos[1];
547  z=fPos[2];
548 }
549 
550 Float_t AliESDv0::GetD(Double_t x0, Double_t y0) const {
551  //--------------------------------------------------------------------
552  // This function returns V0's impact parameter calculated in 2D in XY plane
553  //--------------------------------------------------------------------
554  Double_t x=fPos[0],y=fPos[1];
555  Double_t px=fNmom[0]+fPmom[0];
556  Double_t py=fNmom[1]+fPmom[1];
557 
558  Double_t dz=(x0-x)*py - (y0-y)*px;
559  Double_t d=TMath::Sqrt(dz*dz/(px*px+py*py));
560  return d;
561 }
562 
563 Float_t AliESDv0::GetD(Double_t x0, Double_t y0, Double_t z0) const {
564  //--------------------------------------------------------------------
565  // This function returns V0's impact parameter calculated in 3D
566  //--------------------------------------------------------------------
567  Double_t x=fPos[0],y=fPos[1],z=fPos[2];
568  Double_t px=fNmom[0]+fPmom[0];
569  Double_t py=fNmom[1]+fPmom[1];
570  Double_t pz=fNmom[2]+fPmom[2];
571 
572  Double_t dx=(y0-y)*pz - (z0-z)*py;
573  Double_t dy=(x0-x)*pz - (z0-z)*px;
574  Double_t dz=(x0-x)*py - (y0-y)*px;
575  Double_t d=TMath::Sqrt((dx*dx+dy*dy+dz*dz)/(px*px+py*py+pz*pz));
576  return d;
577 }
578 
579 Float_t AliESDv0::GetV0CosineOfPointingAngle(Double_t refPointX, Double_t refPointY, Double_t refPointZ) const {
580  // calculates the pointing angle of the V0 wrt a reference point
581 
582  Double_t momV0[3]; //momentum of the V0
583  GetPxPyPz(momV0[0],momV0[1],momV0[2]);
584 
585  Double_t deltaPos[3]; //vector between the reference point and the V0 vertex
586  deltaPos[0] = fPos[0] - refPointX;
587  deltaPos[1] = fPos[1] - refPointY;
588  deltaPos[2] = fPos[2] - refPointZ;
589 
590  Double_t momV02 = momV0[0]*momV0[0] + momV0[1]*momV0[1] + momV0[2]*momV0[2];
591  Double_t deltaPos2 = deltaPos[0]*deltaPos[0] + deltaPos[1]*deltaPos[1] + deltaPos[2]*deltaPos[2];
592 
593  Double_t cosinePointingAngle = (deltaPos[0]*momV0[0] +
594  deltaPos[1]*momV0[1] +
595  deltaPos[2]*momV0[2] ) /
596  TMath::Sqrt(momV02 * deltaPos2);
597 
598  return cosinePointingAngle;
599 }
600 
601 
602 // **** The following functions need to be revised
603 
604 void AliESDv0::GetPosCov(Double_t cov[6]) const {
605 
606  for (Int_t i=0; i<6; ++i) cov[i] = fPosCov[i];
607 
608 }
609 
611  //
612  // return sigmay in y at vertex position using covariance matrix
613  //
614  const Double_t * cp = fParamP.GetCovariance();
615  const Double_t * cm = fParamN.GetCovariance();
616  Double_t sigmay = cp[0]+cm[0]+ cp[5]*(fParamP.GetX()-fRr)*(fParamP.GetX()-fRr)+ cm[5]*(fParamN.GetX()-fRr)*(fParamN.GetX()-fRr);
617  return (sigmay>0) ? TMath::Sqrt(sigmay):100;
618 }
619 
621  //
622  // return sigmay in y at vertex position using covariance matrix
623  //
624  const Double_t * cp = fParamP.GetCovariance();
625  const Double_t * cm = fParamN.GetCovariance();
626  Double_t sigmaz = cp[2]+cm[2]+ cp[9]*(fParamP.GetX()-fRr)*(fParamP.GetX()-fRr)+ cm[9]*(fParamN.GetX()-fRr)*(fParamN.GetX()-fRr);
627  return (sigmaz>0) ? TMath::Sqrt(sigmaz):100;
628 }
629 
631  //
632  // Sigma parameterization using covariance matrix
633  //
634  // sigma of distance between two tracks in vertex position
635  // sigma of DCA is proportianal to sigmaD0
636  // factor 2 difference is explained by the fact that the DCA is calculated at the position
637  // where the tracks as closest together ( not exact position of the vertex)
638  //
639  const Double_t * cp = fParamP.GetCovariance();
640  const Double_t * cm = fParamN.GetCovariance();
641  Double_t sigmaD0 = cp[0]+cm[0]+cp[2]+cm[2]+fgkParams.fPSigmaOffsetD0*fgkParams.fPSigmaOffsetD0;
642  sigmaD0 += ((fParamP.GetX()-fRr)*(fParamP.GetX()-fRr))*(cp[5]+cp[9]);
643  sigmaD0 += ((fParamN.GetX()-fRr)*(fParamN.GetX()-fRr))*(cm[5]+cm[9]);
644  return (sigmaD0>0)? TMath::Sqrt(sigmaD0):100;
645 }
646 
647 
649  //
650  //Sigma parameterization using covariance matrices
651  //
652  Double_t prec = TMath::Sqrt((fNmom[0]+fPmom[0])*(fNmom[0]+fPmom[0])
653  +(fNmom[1]+fPmom[1])*(fNmom[1]+fPmom[1])
654  +(fNmom[2]+fPmom[2])*(fNmom[2]+fPmom[2]));
655  Double_t normp = TMath::Sqrt(fPmom[0]*fPmom[0]+fPmom[1]*fPmom[1]+fPmom[2]*fPmom[2])/prec; // fraction of the momenta
656  Double_t normm = TMath::Sqrt(fNmom[0]*fNmom[0]+fNmom[1]*fNmom[1]+fNmom[2]*fNmom[2])/prec;
657  const Double_t * cp = fParamP.GetCovariance();
658  const Double_t * cm = fParamN.GetCovariance();
659  Double_t sigmaAP0 = fgkParams.fPSigmaOffsetAP0*fgkParams.fPSigmaOffsetAP0; // minimal part
660  sigmaAP0 += (cp[5]+cp[9])*(normp*normp)+(cm[5]+cm[9])*(normm*normm); // angular resolution part
661  Double_t sigmaAP1 = GetSigmaD0()/(TMath::Abs(fRr)+0.01); // vertex position part
662  sigmaAP0 += 0.5*sigmaAP1*sigmaAP1;
663  return (sigmaAP0>0)? TMath::Sqrt(sigmaAP0):100;
664 }
665 
667  //
668  // minimax - effective Sigma parameterization
669  // p12 effective curvature and v0 radius postion used as parameters
670  //
671  Double_t p12 = TMath::Sqrt(fParamP.GetParameter()[4]*fParamP.GetParameter()[4]+
673  Double_t sigmaED0= TMath::Max(TMath::Sqrt(fRr)-fgkParams.fPSigmaRminDE,0.0)*fgkParams.fPSigmaCoefDE*p12*p12;
674  sigmaED0*= sigmaED0;
675  sigmaED0*= sigmaED0;
676  sigmaED0 = TMath::Sqrt(sigmaED0+fgkParams.fPSigmaOffsetDE*fgkParams.fPSigmaOffsetDE);
677  return (sigmaED0<fgkParams.fPSigmaMaxDE) ? sigmaED0: fgkParams.fPSigmaMaxDE;
678 }
679 
680 
682  //
683  // effective Sigma parameterization of point angle resolution
684  //
685  Double_t p12 = TMath::Sqrt(fParamP.GetParameter()[4]*fParamP.GetParameter()[4]+
687  Double_t sigmaAPE= fgkParams.fPSigmaBase0APE;
689  sigmaAPE*= (fgkParams.fPSigmaP0APE+fgkParams.fPSigmaP1APE*p12);
690  sigmaAPE = TMath::Min(sigmaAPE,fgkParams.fPSigmaMaxAPE);
691  return sigmaAPE;
692 }
693 
694 
696  //
697  // calculate mini-max effective sigma of point angle resolution
698  //
699  //compv0->fTree->SetAlias("SigmaAP2","max(min((SigmaAP0+SigmaAPE0)*0.5,1.5*SigmaAPE0),0.5*SigmaAPE0+0.003)");
700  Double_t effectiveSigma = GetEffectiveSigmaAP0();
701  Double_t sigmaMMAP = 0.5*(GetSigmaAP0()+effectiveSigma);
702  sigmaMMAP = TMath::Min(sigmaMMAP, fgkParams.fPMaxFractionAP0*effectiveSigma);
703  sigmaMMAP = TMath::Max(sigmaMMAP, fgkParams.fPMinFractionAP0*effectiveSigma+fgkParams.fPMinAP0);
704  return sigmaMMAP;
705 }
707  //
708  // calculate mini-max sigma of dca resolution
709  //
710  //compv0->fTree->SetAlias("SigmaD2","max(min((SigmaD0+SigmaDE0)*0.5,1.5*SigmaDE0),0.5*SigmaDE0)");
711  Double_t effectiveSigma = GetEffectiveSigmaD0();
712  Double_t sigmaMMD0 = 0.5*(GetSigmaD0()+effectiveSigma);
713  sigmaMMD0 = TMath::Min(sigmaMMD0, fgkParams.fPMaxFractionD0*effectiveSigma);
714  sigmaMMD0 = TMath::Max(sigmaMMD0, fgkParams.fPMinFractionD0*effectiveSigma+fgkParams.fPMinD0);
715  return sigmaMMD0;
716 }
717 
718 
719 Double_t AliESDv0::GetLikelihoodAP(Int_t mode0, Int_t mode1){
720  //
721  // get likelihood for point angle
722  //
723  Double_t sigmaAP = 0.007; //default sigma
724  switch (mode0){
725  case 0:
726  sigmaAP = GetSigmaAP0(); // mode 0 - covariance matrix estimates used
727  break;
728  case 1:
729  sigmaAP = GetEffectiveSigmaAP0(); // mode 1 - effective sigma used
730  break;
731  case 2:
732  sigmaAP = GetMinimaxSigmaAP0(); // mode 2 - minimax sigma
733  break;
734  }
735  Double_t apNorm = TMath::Min(TMath::ACos(fPointAngle)/sigmaAP,50.);
736  //normalized point angle, restricted - because of overflow problems in Exp
737  Double_t likelihood = 0;
738  switch(mode1){
739  case 0:
740  likelihood = TMath::Exp(-0.5*apNorm*apNorm);
741  // one component
742  break;
743  case 1:
744  likelihood = (TMath::Exp(-0.5*apNorm*apNorm)+0.5* TMath::Exp(-0.25*apNorm*apNorm))/1.5;
745  // two components
746  break;
747  case 2:
748  likelihood = (TMath::Exp(-0.5*apNorm*apNorm)+0.5* TMath::Exp(-0.25*apNorm*apNorm)+0.25*TMath::Exp(-0.125*apNorm*apNorm))/1.75;
749  // three components
750  break;
751  }
752  return likelihood;
753 }
754 
755 Double_t AliESDv0::GetLikelihoodD(Int_t mode0, Int_t mode1){
756  //
757  // get likelihood for DCA
758  //
759  Double_t sigmaD = 0.03; //default sigma
760  switch (mode0){
761  case 0:
762  sigmaD = GetSigmaD0(); // mode 0 - covariance matrix estimates used
763  break;
764  case 1:
765  sigmaD = GetEffectiveSigmaD0(); // mode 1 - effective sigma used
766  break;
767  case 2:
768  sigmaD = GetMinimaxSigmaD0(); // mode 2 - minimax sigma
769  break;
770  }
771 
772  //Bo: Double_t dNorm = TMath::Min(fDist2/sigmaD,50.);
773  Double_t dNorm = TMath::Min(fDcaV0Daughters/sigmaD,50.);//Bo:
774  //normalized point angle, restricted - because of overflow problems in Exp
775  Double_t likelihood = 0;
776  switch(mode1){
777  case 0:
778  likelihood = TMath::Exp(-2.*dNorm);
779  // one component
780  break;
781  case 1:
782  likelihood = (TMath::Exp(-2.*dNorm)+0.5* TMath::Exp(-dNorm))/1.5;
783  // two components
784  break;
785  case 2:
786  likelihood = (TMath::Exp(-2.*dNorm)+0.5* TMath::Exp(-dNorm)+0.25*TMath::Exp(-0.5*dNorm))/1.75;
787  // three components
788  break;
789  }
790  return likelihood;
791 
792 }
793 
794 Double_t AliESDv0::GetLikelihoodC(Int_t mode0, Int_t /*mode1*/) const {
795  //
796  // get likelihood for Causality
797  // !!! Causality variables defined in AliITStrackerMI !!!
798  // when more information was available
799  //
800  Double_t likelihood = 0.5;
801  Double_t minCausal = TMath::Min(fCausality[0],fCausality[1]);
802  Double_t maxCausal = TMath::Max(fCausality[0],fCausality[1]);
803  // minCausal = TMath::Max(minCausal,0.5*maxCausal);
804  //compv0->fTree->SetAlias("LCausal","(1.05-(2*(0.8-exp(-max(RC.fV0rec.fCausality[0],RC.fV0rec.fCausality[1])))+2*(0.8-exp(-min(RC.fV0rec.fCausality[0],RC.fV0rec.fCausality[1]))))/2)**4");
805 
806  switch(mode0){
807  case 0:
808  //normalization
809  likelihood = TMath::Power((1.05-2*(0.8-TMath::Exp(-maxCausal))),4.);
810  break;
811  case 1:
812  likelihood = TMath::Power(1.05-(2*(0.8-TMath::Exp(-maxCausal))+(2*(0.8-TMath::Exp(-minCausal))))*0.5,4.);
813  break;
814  }
815  return likelihood;
816 
817 }
818 
819 void AliESDv0::SetCausality(Float_t pb0, Float_t pb1, Float_t pa0, Float_t pa1)
820 {
821  //
822  // set probabilities
823  //
824  fCausality[0] = pb0; // probability - track 0 exist before vertex
825  fCausality[1] = pb1; // probability - track 1 exist before vertex
826  fCausality[2] = pa0; // probability - track 0 exist close after vertex
827  fCausality[3] = pa1; // probability - track 1 exist close after vertex
828 }
829 void AliESDv0::SetClusters(const Int_t *clp, const Int_t *clm)
830 {
831  //
832  // Set its clusters indexes
833  //
834  for (Int_t i=0;i<6;i++) fClusters[0][i] = clp[i];
835  for (Int_t i=0;i<6;i++) fClusters[1][i] = clm[i];
836 }
837 
838 Double_t AliESDv0::GetEffMass(UInt_t p1, UInt_t p2) const{
839  //
840  // calculate effective mass
841  //
842  const Double_t kpmass[5] = {TDatabasePDG::Instance()->GetParticle(kElectron)->Mass(),
843  TDatabasePDG::Instance()->GetParticle(kMuonMinus)->Mass(),
844  TDatabasePDG::Instance()->GetParticle(kPiPlus)->Mass(),
845  TDatabasePDG::Instance()->GetParticle(kKPlus)->Mass(),
846  TDatabasePDG::Instance()->GetParticle(kProton)->Mass()};
847  /*
848  if (p1>4) return -1;
849  if (p2>4) return -1;
850  Double_t mass1 = kpmass[p1];
851  Double_t mass2 = kpmass[p2];
852  const Double_t *m1 = fPmom;
853  const Double_t *m2 = fNmom;
854  //
855  //if (fRP[p1]+fRM[p2]<fRP[p2]+fRM[p1]){
856  // m1 = fPM;
857  // m2 = fPP;
858  //}
859  //
860  Double_t e1 = TMath::Sqrt(mass1*mass1+
861  m1[0]*m1[0]+
862  m1[1]*m1[1]+
863  m1[2]*m1[2]);
864  Double_t e2 = TMath::Sqrt(mass2*mass2+
865  m2[0]*m2[0]+
866  m2[1]*m2[1]+
867  m2[2]*m2[2]);
868  Double_t mass =
869  (m2[0]+m1[0])*(m2[0]+m1[0])+
870  (m2[1]+m1[1])*(m2[1]+m1[1])+
871  (m2[2]+m1[2])*(m2[2]+m1[2]);
872 
873  mass = (e1+e2)*(e1+e2)-mass;
874  if (mass < 0.) mass = 0.;
875  return (TMath::Sqrt(mass));
876  */
877  if(p1>4 || p2>4) return -1;
878  return GetEffMassExplicit(kpmass[p1],kpmass[p2]);
879 }
880 
881 Double_t AliESDv0::GetEffMassExplicit(Double_t m1, Double_t m2) const {
882  //
883  // calculate effective mass with given masses of decay products
884  //
885  const AliExternalTrackParam *paramP = GetParamP();
886  const AliExternalTrackParam *paramN = GetParamN();
887  if (paramP->GetParameter()[4]<0){
888  paramP=GetParamN();
889  paramN=GetParamP();
890  }
891  Double_t pmom[3]={0}, nmom[3]={0};
892  paramP->GetPxPyPz(pmom);
893  paramN->GetPxPyPz(nmom);
894  Double_t e12 = m1*m1+pmom[0]*pmom[0]+pmom[1]*pmom[1]+pmom[2]*pmom[2];
895  Double_t e22 = m2*m2+nmom[0]*nmom[0]+nmom[1]*nmom[1]+nmom[2]*nmom[2];
896  Double_t cmass = TMath::Sqrt(TMath::Max(m1*m1+m2*m2
897  +2.*(TMath::Sqrt(e12*e22)-pmom[0]*nmom[0]-pmom[1]*nmom[1]-pmom[2]*nmom[2]),0.));
898  return cmass;
899 
900 }
901 
902 
903 
904 Double_t AliESDv0::GetKFInfo(UInt_t p1, UInt_t p2, Int_t type) const{
905  //
906  // type:
907  // 0 - return mass
908  // 1 - return err mass
909  // 2 - return chi2
910  //
911  const Int_t spdg[5]={kPositron,kMuonPlus,kPiPlus, kKPlus, kProton};
912  const AliExternalTrackParam *paramP = GetParamP();
913  const AliExternalTrackParam *paramN = GetParamN();
914  if (paramP->GetParameter()[4]<0){
915  paramP=GetParamN();
916  paramN=GetParamP();
917  }
918  AliKFParticle kfp1( *(paramP), spdg[p1] *TMath::Sign(1,p1) );
919  AliKFParticle kfp2( *(paramN), spdg[p2] *TMath::Sign(1,p2) );
920  AliKFParticle v0KF;
921  v0KF+=kfp1;
922  v0KF+=kfp2;
923  if (type==0) return v0KF.GetMass();
924  if (type==1) return v0KF.GetErrMass();
925  if (type==2) return v0KF.GetChi2();
926  if (type==3) return v0KF.GetPt();
927 
928  return 0;
929 }
930 
940 Double_t AliESDv0::GetKFInfoScale(UInt_t p1, UInt_t p2, Int_t type, Double_t d1pt, Double_t s1pt, Double_t eLoss, Int_t flag) const{
941  // Used for fast numerical studies of impact of residual miscalibration
942  // Input:
943  // p1,p2 - particle type (0-el,1-mu, 2-pi, 3-K, 4-Proton)
944  // type
945  // 0 - return mass
946  // 1 - return err mass
947  // 2 - return chi2
948  // 3 - return Pt of the V0 ()
949  // distrotion
950  // d1pt - 1/pt shift
951  // s1pt - scaling of 1/pt
952  // Important function to benchmark the pt resolution, and to find out systematic distortion
953  // Used for fast numerical emulation of impact of the resifual misscalibration
954  //
955  /*
956  Example analysis using filtered trees:
957  // ratio of pt: 0.0015 1/GeV shift -> 0 mean shift for the K0 but ~ 1.1 % shift at 10 GeV for the Lambda
958  treeLambda->Draw("v0.GetKFInfoScale(4,2,3,0.0015,0.0)/v0.GetKFInfoScale(4,2,3,0.0,0)-1:v0.Pt()>>ptratioLambda(50,0,20)","v0.Pt()<20","prof",100000);
959  treeALambda->Draw("v0.GetKFInfoScale(2,4,3,0.0015,0.0)/v0.GetKFInfoScale(2,4,3,0.0,0)-1:v0.Pt()>>ptratioALambda(50,0,20)","v0.Pt()<20","profsame",100000);
960  treeK0->Draw("v0.GetKFInfoScale(2,2,3,0.0015,0.0)/v0.GetKFInfoScale(2,2,3,0.0,0)-1:v0.Pt()>>ptratioK0(50,0,20)","v0.Pt()<20","profsame",100000);
961  // pt spectra ratios:
962  // ratio of spectra - e.g q/pt shift ~ 0.0015 1/GeV -> Lambda - 10 % difference in pt yeald at 10 GeV - K0 - unaffected
963  treeLambda->Draw("v0.GetKFInfoScale(4,2,3,0.0015,0.0)>>hisLambdaShifted0015(50,0,20)","","",1000000);
964  treeLambda->Draw("v0.GetKFInfoScale(4,2,3,0.00,0.0)>>hisLambdaShifted0(50,0,20)","","",1000000);
965  TH1F ratioLambdaShifted0015 = *hisLambdaShifted0015;
966  ratioLambdaShifted0015.Sumw2();
967  ratioLambdaShifted0015.Divide(hisLambdaShifted0);
968  */
969  const Int_t spdg[5]={kPositron,kMuonPlus,kPiPlus, kKPlus, kProton};
970  if (p1>4 || p2>4)return 0;
971  AliExternalTrackParam paramP;
972  AliExternalTrackParam paramN;
973  if (GetParamP()->GetParameter()[4]<0){
974  paramP=*(GetParamN());
975  paramN=*(GetParamP());
976  }else{
977  paramP=*(GetParamP());
978  paramN=*(GetParamN());
979  }
980  // calculate delta of energy loss
981  Double_t bg1=paramP.P() /TDatabasePDG::Instance()->GetParticle(spdg[p1])->Mass();
982  Double_t bg2=paramN.P() /TDatabasePDG::Instance()->GetParticle(spdg[p2])->Mass();
985 
986  Double_t *pparam1 = (Double_t*)paramP.GetParameter();
987  Double_t *pparam2 = (Double_t*)paramN.GetParameter();
988  if (flag&0x1) pparam1[4]+=d1pt;
989  if (flag&0x2) pparam2[4]+=d1pt;
990  if (flag&0x1) pparam1[4]*=(1+s1pt)*(1.+eLoss*dP1/paramP.P());
991  if (flag&0x2) pparam2[4]*=(1+s1pt)*(1.+eLoss*dP2/paramN.P());
992 
993  //
994  AliKFParticle kfp1( paramP, spdg[p1] *TMath::Sign(1,p1) );
995  AliKFParticle kfp2( paramN, spdg[p2] *TMath::Sign(1,p2) );
996  AliKFParticle *v0KF = new AliKFParticle;
997  *(v0KF)+=kfp1;
998  *(v0KF)+=kfp2;
999  if (type==0) return v0KF->GetMass();
1000  if (type==1) return v0KF->GetErrMass();
1001  if (type==2) return v0KF->GetChi2();
1002  if (type==3) return v0KF->GetPt();
1003  return 0;
1004 }
1005 
1006 
1008  // Build an AliESDVertex from this AliESDv0
1010  return vtx;
1011 }
1012 
Double_t fPMaxFractionD0
Double_t GetSigmaZY() const
Double32_t fAngle[3]
Definition: AliESDv0.h:172
Int_t fNidx
its clusters CKBrev
Definition: AliESDv0.h:180
Double_t fPSigmaP0APE
Double32_t fNormDCAPrim[2]
Definition: AliESDv0.h:164
Double32_t fEffMass
Definition: AliESDv0.h:157
Double32_t fNmom[3]
Definition: AliESDv0.h:162
Double_t GetSigmaZ2() const
Double_t GetLikelihoodD(Int_t mode0, Int_t mode1)
Definition: AliESDv0.cxx:755
virtual ~AliESDv0()
Definition: AliESDv0.cxx:399
Double_t fPSigmaOffsetAP0
Short_t fStatus
Definition: AliESDv0.h:185
virtual Double_t E() const
Definition: AliESDv0.cxx:406
Double_t fPMaxFractionAP0
virtual Double_t Y() const
Definition: AliESDv0.cxx:413
Double_t GetChi2() const
Short_t fNBefore
Definition: AliESDv0.h:186
Double_t GetErrMass() const
Double32_t fDcaV0Daughters
Definition: AliESDv0.h:158
void SetCausality(Float_t pb0, Float_t pb1, Float_t pa0, Float_t pa1)
Definition: AliESDv0.cxx:819
Bool_t GetXYZ(Double_t *p) const
Double32_t fDistSigma
Definition: AliESDv0.h:166
AliVParticle & operator=(const AliVParticle &vPart)
Int_t fClusters[2][6]
Definition: AliESDv0.h:179
static Double_t BetheBlochGeant(Double_t bg, Double_t kp0=2.33, Double_t kp1=0.20, Double_t kp2=3.00, Double_t kp3=173e-9, Double_t kp4=0.49848)
Double_t GetSigmaD0()
Definition: AliESDv0.cxx:630
Double_t fPMinFractionAP0
Double_t GetKFInfoScale(UInt_t p1, UInt_t p2, Int_t type, Double_t d1pt, Double_t s1pt, Double_t eLoss=0, Int_t flag=0x3) const
Definition: AliESDv0.cxx:940
AliExternalTrackParam fParamN
Definition: AliESDv0.h:152
virtual void Copy(TObject &obj) const
Definition: AliESDv0.cxx:387
Double32_t fChi2Before
Definition: AliESDv0.h:167
void GetXYZ(Double_t &x, Double_t &y, Double_t &z) const
Definition: AliESDv0.cxx:541
Double_t GetAlpha() const
Double_t AlphaV0() const
Definition: AliESDv0.cxx:451
Int_t fPidx
Definition: AliESDv0.h:181
Double_t GetSigmaZ()
Definition: AliESDv0.cxx:620
virtual Double_t GetTgl() const
Double32_t fPointAngleTh
Definition: AliESDv0.h:174
Double32_t fRr
Definition: AliESDv0.h:165
const Double_t * GetParameter() const
const Float_t kAlmost0
Double_t GetLikelihoodAP(Int_t mode0, Int_t mode1)
Definition: AliESDv0.cxx:719
Double_t GetEffMass() const
Definition: AliESDv0.h:77
Double32_t fChi2After
Definition: AliESDv0.h:168
Double_t fPSigmaCoefDE
Double_t RapK0Short() const
Definition: AliESDv0.cxx:437
const AliExternalTrackParam * GetParamP() const
Definition: AliESDv0.h:94
virtual Short_t Charge() const
Double_t GetMinimaxSigmaD0()
Definition: AliESDv0.cxx:706
Double_t GetSigmaY2() const
Double_t fPSigmaMaxDE
Double_t GetV0CosineOfPointingAngle() const
Definition: AliESDv0.h:90
AliESDVertex GetVertex() const
Definition: AliESDv0.cxx:1007
Double_t ChangeMassHypothesis(Int_t code=kK0Short)
Definition: AliESDv0.cxx:477
Double32_t fPos[3]
Definition: AliESDv0.h:160
Double32_t fCausality[4]
Definition: AliESDv0.h:171
const Double_t * GetCovariance() const
Double32_t fPointAngle
Definition: AliESDv0.h:175
Int_t fPdgCode
Definition: AliESDv0.h:178
AliExternalTrackParam fParamP
Definition: AliESDv0.h:153
Bool_t GetPxPyPz(Double_t *p) const
Bool_t fOnFlyStatus
Definition: AliESDv0.h:189
void SetClusters(const Int_t *clp, const Int_t *clm)
Definition: AliESDv0.cxx:829
AliTPCcalibV0 v0
void GetPxPyPz(Double_t &px, Double_t &py, Double_t &pz) const
Definition: AliESDv0.cxx:532
Double_t fPSigmaP1APE
AliESDv0 & operator=(const AliESDv0 &v0)
Definition: AliESDv0.cxx:336
Double32_t fPmom[3]
Definition: AliESDv0.h:163
Double_t fPSigmaR1APE
Double_t GetPt() const
Double_t GetEffMassExplicit(Double_t m1, Double_t m2) const
Definition: AliESDv0.cxx:881
Int_t Refit()
Definition: AliESDv0.cxx:280
void GetPosCov(Double_t cov[6]) const
Definition: AliESDv0.cxx:604
Double_t PtArmV0() const
Definition: AliESDv0.cxx:466
Double_t fPMinFractionD0
Double_t GetKFInfo(UInt_t p1, UInt_t p2, Int_t type) const
Definition: AliESDv0.cxx:904
static const AliESDV0Params fgkParams
Definition: AliESDv0.h:194
Double_t GetMinimaxSigmaAP0()
Definition: AliESDv0.cxx:695
virtual Double_t Px() const
Definition: AliESDv0.h:38
Double_t fPSigmaOffsetDE
Double_t GetSigmaAP0()
Definition: AliESDv0.cxx:648
Double_t GetSigmaY()
Definition: AliESDv0.cxx:610
Double_t GetEffectiveSigmaAP0()
Definition: AliESDv0.cxx:681
#define AliError(message)
Definition: AliLog.h:591
virtual Double_t Py() const
Definition: AliESDv0.h:39
Double_t fPSigmaR0APE
static Bool_t GetWeight(TMatrixD &w, const AliExternalTrackParam &t)
Definition: AliESDv0.cxx:228
AliESDv0()
Definition: AliESDv0.cxx:45
Double_t GetMass() const
virtual Double_t P() const
Definition: AliESDv0.h:42
Double_t fPSigmaMaxAPE
Double_t fPSigmaBase0APE
Float_t GetD(Double_t x0, Double_t y0) const
Definition: AliESDv0.cxx:550
Double_t GetEffectiveSigmaD0()
Definition: AliESDv0.cxx:666
virtual Double_t Pz() const
Definition: AliESDv0.h:40
Double32_t fPointAngleFi
Definition: AliESDv0.h:173
Double_t fPSigmaRminDE
Double_t GetLikelihoodC(Int_t mode0, Int_t mode1) const
Definition: AliESDv0.cxx:794
Double_t fPSigmaOffsetD0
const AliExternalTrackParam * GetParamN() const
Definition: AliESDv0.h:95
Double32_t fPosCov[6]
Definition: AliESDv0.h:161
Double32_t fChi2V0
Definition: AliESDv0.h:159
Short_t fNAfter
Definition: AliESDv0.h:187
Double_t RapLambda() const
Definition: AliESDv0.cxx:444