]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDv0.cxx
Shuttle getters and setters for main OCDB/Reference, local OCDB/Reference, temp
[u/mrichter/AliRoot.git] / STEER / AliESDv0.cxx
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 <Riostream.h>
29 #include <TMath.h>
30 #include <TDatabasePDG.h>
31 #include <TPDGCode.h>
32 #include <TParticlePDG.h>
33
34 #include "AliLog.h"
35 #include "AliESDv0.h"
36 #include "AliExternalTrackParam.h"
37
38 ClassImp(AliESDv0)
39
40 AliESDV0Params  AliESDv0::fgkParams;
41
42 AliESDv0::AliESDv0() :
43   TObject(),
44   fOnFlyStatus(kFALSE),
45   fPdgCode(kK0Short),
46   fEffMass(TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass()),
47   fDcaV0Daughters(0),
48   fPointAngle(0),
49   fChi2V0(31.),
50   fNidx(-1),
51   fParamN(),
52   fPidx(-1),
53   fParamP(),
54   fRr(-1),
55   fStatus(0),
56   fDistSigma(0),
57   fChi2Before(0),
58   fNBefore(0),
59   fChi2After(0),
60   fNAfter(0),
61   fPointAngleFi(0),
62   fPointAngleTh(0)
63 {
64   //--------------------------------------------------------------------
65   // Default constructor  (K0s)
66   //--------------------------------------------------------------------
67
68   for (Int_t i=0; i<3; i++) {
69     fPos[i] = 0.;
70     fNmom[i] = 0.;
71     fPmom[i] = 0.;
72   }
73
74   for (Int_t i=0; i<6; i++) {
75     fPosCov[i]= 0.;
76   }
77
78   for (Int_t i=0;i<6;i++){fClusters[0][i]=0; fClusters[1][i]=0;}
79   fNormDCAPrim[0]=fNormDCAPrim[1]=0;
80   for (Int_t i=0;i<3;i++){fAngle[i]=0;}
81   for (Int_t i=0;i<4;i++){fCausality[i]=0;}
82 }
83
84 AliESDv0::AliESDv0(const AliESDv0& v0) :
85   TObject(v0),
86   fOnFlyStatus(v0.fOnFlyStatus),
87   fPdgCode(v0.fPdgCode),
88   fEffMass(v0.fEffMass),
89   fDcaV0Daughters(v0.fDcaV0Daughters),
90   fPointAngle(v0.fPointAngle),
91   fChi2V0(v0.fChi2V0),
92   fNidx(v0.fNidx),
93   fParamN(v0.fParamN),
94   fPidx(v0.fPidx),
95   fParamP(v0.fParamP),
96   fRr(v0.fRr),
97   fStatus(v0.fStatus),
98   fDistSigma(v0.fDistSigma),
99   fChi2Before(v0.fChi2Before),
100   fNBefore(v0.fNBefore),
101   fChi2After(v0.fChi2After),
102   fNAfter(v0.fNAfter),
103   fPointAngleFi(v0.fPointAngleFi),
104   fPointAngleTh(v0.fPointAngleTh)
105 {
106   //--------------------------------------------------------------------
107   // The copy constructor
108   //--------------------------------------------------------------------
109
110   for (int i=0; i<3; i++) {
111     fPos[i]  = v0.fPos[i];
112     fNmom[i] = v0.fNmom[i];
113     fPmom[i] = v0.fPmom[i];
114   }
115   for (int i=0; i<6; i++) {
116     fPosCov[i]  = v0.fPosCov[i];
117   }
118
119   for (Int_t i=0; i<2; i++) {
120     fNormDCAPrim[i]=v0.fNormDCAPrim[i];
121   }
122   for (Int_t i=0;i<6;i++){
123       fClusters[0][i]=v0.fClusters[0][i]; 
124       fClusters[1][i]=v0.fClusters[1][i];
125   }
126   for (Int_t i=0;i<3;i++){
127       fAngle[i]=v0.fAngle[i];
128   }
129   for (Int_t i=0;i<4;i++){fCausality[i]=v0.fCausality[i];}
130 }
131
132 AliESDv0::AliESDv0(const AliExternalTrackParam &t1, Int_t i1,
133                    const AliExternalTrackParam &t2, Int_t i2) :
134   TObject(),
135   fOnFlyStatus(kFALSE),
136   fPdgCode(kK0Short),
137   fEffMass(TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass()),
138   fDcaV0Daughters(0),
139   fPointAngle(0),
140   fChi2V0(31.),
141   fNidx(i1),
142   fParamN(t1),
143   fPidx(i2),
144   fParamP(t2),
145   fRr(-1),
146   fStatus(0),
147   fDistSigma(0),
148   fChi2Before(0),
149   fNBefore(0),
150   fChi2After(0),
151   fNAfter(0),
152   fPointAngleFi(0),
153   fPointAngleTh(0)
154 {
155   //--------------------------------------------------------------------
156   // Main constructor  (K0s)
157   //--------------------------------------------------------------------
158
159   for (Int_t i=0; i<6; i++) {
160     fPosCov[i]= 0.;
161   }
162
163   //Trivial estimation of the vertex parameters
164   Double_t alpha=t1.GetAlpha(), cs=TMath::Cos(alpha), sn=TMath::Sin(alpha);
165   Double_t tmp[3];
166   t1.GetPxPyPz(tmp);
167   Double_t px1=tmp[0], py1=tmp[1], pz1=tmp[2];
168   t1.GetXYZ(tmp);
169   Double_t  x1=tmp[0],  y1=tmp[1],  z1=tmp[2];
170   const Double_t ss=0.0005*0.0005;//a kind of a residual misalignment precision
171   Double_t sx1=sn*sn*t1.GetSigmaY2()+ss, sy1=cs*cs*t1.GetSigmaY2()+ss; 
172
173
174   alpha=t2.GetAlpha(); cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
175   t2.GetPxPyPz(tmp);
176   Double_t px2=tmp[0], py2=tmp[1], pz2=tmp[2];
177   t2.GetXYZ(tmp);
178   Double_t  x2=tmp[0],  y2=tmp[1],  z2=tmp[2];
179   Double_t sx2=sn*sn*t2.GetSigmaY2()+ss, sy2=cs*cs*t2.GetSigmaY2()+ss; 
180     
181   Double_t sz1=t1.GetSigmaZ2(), sz2=t2.GetSigmaZ2();
182   Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
183   Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
184   Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
185   fPos[0]=wx1*x1 + wx2*x2; fPos[1]=wy1*y1 + wy2*y2; fPos[2]=wz1*z1 + wz2*z2;
186
187   //fPos[0]=0.5*(x1+x2); fPos[1]=0.5*(y1+y2); fPos[2]=0.5*(z1+z2);
188   fNmom[0]=px1; fNmom[1]=py1; fNmom[2]=pz1; 
189   fPmom[0]=px2; fPmom[1]=py2; fPmom[2]=pz2;
190
191   for (Int_t i=0;i<6;i++){fClusters[0][i]=0; fClusters[1][i]=0;}
192   fNormDCAPrim[0]=fNormDCAPrim[1]=0;
193   for (Int_t i=0;i<3;i++){fAngle[i]=0;}
194   for (Int_t i=0;i<4;i++){fCausality[i]=0;}
195 }
196
197 AliESDv0::~AliESDv0(){
198   //--------------------------------------------------------------------
199   // Empty destructor
200   //--------------------------------------------------------------------
201 }
202
203
204
205 Double_t AliESDv0::ChangeMassHypothesis(Int_t code) {
206   //--------------------------------------------------------------------
207   // This function changes the mass hypothesis for this V0
208   // and returns the "kinematical quality" of this hypothesis 
209   //--------------------------------------------------------------------
210   static
211   Double_t piMass=TDatabasePDG::Instance()->GetParticle(kPiPlus)->Mass();
212   static
213   Double_t prMass=TDatabasePDG::Instance()->GetParticle(kProton)->Mass();
214   static
215   Double_t k0Mass=TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass();
216   static
217   Double_t l0Mass=TDatabasePDG::Instance()->GetParticle(kLambda0)->Mass();
218
219   Double_t nmass=piMass, pmass=piMass, mass=k0Mass, ps=0.206;
220
221   fPdgCode=code;
222
223   switch (code) {
224   case kLambda0:
225     nmass=piMass; pmass=prMass; mass=l0Mass; ps=0.101; break;
226   case kLambda0Bar:
227     pmass=piMass; nmass=prMass; mass=l0Mass; ps=0.101; break;
228   case kK0Short: 
229     break;
230   default:
231     AliError("invalide PDG code ! Assuming K0s...");
232     fPdgCode=kK0Short;
233     break;
234   }
235
236   Double_t pxn=fNmom[0], pyn=fNmom[1], pzn=fNmom[2]; 
237   Double_t pxp=fPmom[0], pyp=fPmom[1], pzp=fPmom[2];
238
239   Double_t en=TMath::Sqrt(nmass*nmass + pxn*pxn + pyn*pyn + pzn*pzn);
240   Double_t ep=TMath::Sqrt(pmass*pmass + pxp*pxp + pyp*pyp + pzp*pzp);
241   Double_t pxl=pxn+pxp, pyl=pyn+pyp, pzl=pzn+pzp;
242   Double_t pl=TMath::Sqrt(pxl*pxl + pyl*pyl + pzl*pzl);
243
244   fEffMass=TMath::Sqrt((en+ep)*(en+ep)-pl*pl);
245
246   Double_t beta=pl/(en+ep);
247   Double_t pln=(pxn*pxl + pyn*pyl + pzn*pzl)/pl;
248   Double_t plp=(pxp*pxl + pyp*pyl + pzp*pzl)/pl;
249
250   Double_t pt2=pxp*pxp + pyp*pyp + pzp*pzp - plp*plp;
251
252   Double_t a=(plp-pln)/(plp+pln);
253   a -= (pmass*pmass-nmass*nmass)/(mass*mass);
254   a = 0.25*beta*beta*mass*mass*a*a + pt2;
255
256   return (a - ps*ps);
257   
258 }
259
260 void AliESDv0::GetPxPyPz(Double_t &px, Double_t &py, Double_t &pz) const {
261   //--------------------------------------------------------------------
262   // This function returns V0's momentum (global)
263   //--------------------------------------------------------------------
264   px=fNmom[0]+fPmom[0]; 
265   py=fNmom[1]+fPmom[1]; 
266   pz=fNmom[2]+fPmom[2]; 
267 }
268
269 void AliESDv0::GetXYZ(Double_t &x, Double_t &y, Double_t &z) const {
270   //--------------------------------------------------------------------
271   // This function returns V0's position (global)
272   //--------------------------------------------------------------------
273   x=fPos[0]; 
274   y=fPos[1]; 
275   z=fPos[2]; 
276 }
277
278 Float_t AliESDv0::GetD(Double_t x0, Double_t y0, Double_t z0) const {
279   //--------------------------------------------------------------------
280   // This function returns V0's impact parameter
281   //--------------------------------------------------------------------
282   Double_t x=fPos[0],y=fPos[1],z=fPos[2];
283   Double_t px=fNmom[0]+fPmom[0];
284   Double_t py=fNmom[1]+fPmom[1];
285   Double_t pz=fNmom[2]+fPmom[2];
286
287   Double_t dx=(y0-y)*pz - (z0-z)*py; 
288   Double_t dy=(x0-x)*pz - (z0-z)*px;
289   Double_t dz=(x0-x)*py - (y0-y)*px;
290   Double_t d=TMath::Sqrt((dx*dx+dy*dy+dz*dz)/(px*px+py*py+pz*pz));
291   return d;
292 }
293
294
295 Float_t AliESDv0::GetV0CosineOfPointingAngle(Double_t& refPointX, Double_t& refPointY, Double_t& refPointZ) const {
296   // calculates the pointing angle of the V0 wrt a reference point
297
298   Double_t momV0[3]; //momentum of the V0
299   GetPxPyPz(momV0[0],momV0[1],momV0[2]);
300
301   Double_t deltaPos[3]; //vector between the reference point and the V0 vertex
302   deltaPos[0] = fPos[0] - refPointX;
303   deltaPos[1] = fPos[1] - refPointY;
304   deltaPos[2] = fPos[2] - refPointZ;
305
306   Double_t momV02    = momV0[0]*momV0[0] + momV0[1]*momV0[1] + momV0[2]*momV0[2];
307   Double_t deltaPos2 = deltaPos[0]*deltaPos[0] + deltaPos[1]*deltaPos[1] + deltaPos[2]*deltaPos[2];
308
309   Double_t cosinePointingAngle = (deltaPos[0]*momV0[0] +
310                                   deltaPos[1]*momV0[1] +
311                                   deltaPos[2]*momV0[2] ) /
312     TMath::Sqrt(momV02 * deltaPos2);
313   
314   return cosinePointingAngle;
315 }
316
317
318 // **** The following functions need to be revised
319
320 Double_t AliESDv0::GetSigmaY(){
321   //
322   // return sigmay in y  at vertex position  using covariance matrix 
323   //
324   const Double_t * cp  = fParamP.GetCovariance();
325   const Double_t * cm  = fParamN.GetCovariance();
326   Double_t sigmay = cp[0]+cm[0]+ cp[5]*(fParamP.GetX()-fRr)*(fParamP.GetX()-fRr)+ cm[5]*(fParamN.GetX()-fRr)*(fParamN.GetX()-fRr);
327   return (sigmay>0) ? TMath::Sqrt(sigmay):100;
328 }
329
330 Double_t AliESDv0::GetSigmaZ(){
331   //
332   // return sigmay in y  at vertex position  using covariance matrix 
333   //
334   const Double_t * cp  = fParamP.GetCovariance();
335   const Double_t * cm  = fParamN.GetCovariance();
336   Double_t sigmaz = cp[2]+cm[2]+ cp[9]*(fParamP.GetX()-fRr)*(fParamP.GetX()-fRr)+ cm[9]*(fParamN.GetX()-fRr)*(fParamN.GetX()-fRr);
337   return (sigmaz>0) ? TMath::Sqrt(sigmaz):100;
338 }
339
340 Double_t AliESDv0::GetSigmaD0(){
341   //
342   // Sigma parameterization using covariance matrix
343   //
344   // sigma of distance between two tracks in vertex position 
345   // sigma of DCA is proportianal to sigmaD0
346   // factor 2 difference is explained by the fact that the DCA is calculated at the position 
347   // where the tracks as closest together ( not exact position of the vertex)
348   //
349   const Double_t * cp      = fParamP.GetCovariance();
350   const Double_t * cm      = fParamN.GetCovariance();
351   Double_t sigmaD0   = cp[0]+cm[0]+cp[2]+cm[2]+fgkParams.fPSigmaOffsetD0*fgkParams.fPSigmaOffsetD0;
352   sigmaD0           += ((fParamP.GetX()-fRr)*(fParamP.GetX()-fRr))*(cp[5]+cp[9]);
353   sigmaD0           += ((fParamN.GetX()-fRr)*(fParamN.GetX()-fRr))*(cm[5]+cm[9]);
354   return (sigmaD0>0)? TMath::Sqrt(sigmaD0):100;
355 }
356
357
358 Double_t AliESDv0::GetSigmaAP0(){
359   //
360   //Sigma parameterization using covariance matrices
361   //
362   Double_t prec  = TMath::Sqrt((fNmom[0]+fPmom[0])*(fNmom[0]+fPmom[0])
363                               +(fNmom[1]+fPmom[1])*(fNmom[1]+fPmom[1])
364                               +(fNmom[2]+fPmom[2])*(fNmom[2]+fPmom[2]));
365   Double_t normp = TMath::Sqrt(fPmom[0]*fPmom[0]+fPmom[1]*fPmom[1]+fPmom[2]*fPmom[2])/prec;  // fraction of the momenta
366   Double_t normm = TMath::Sqrt(fNmom[0]*fNmom[0]+fNmom[1]*fNmom[1]+fNmom[2]*fNmom[2])/prec;  
367   const Double_t * cp      = fParamP.GetCovariance();
368   const Double_t * cm      = fParamN.GetCovariance();
369   Double_t sigmaAP0 = fgkParams.fPSigmaOffsetAP0*fgkParams.fPSigmaOffsetAP0;                           // minimal part
370   sigmaAP0 +=  (cp[5]+cp[9])*(normp*normp)+(cm[5]+cm[9])*(normm*normm);          // angular resolution part
371   Double_t sigmaAP1 = GetSigmaD0()/(TMath::Abs(fRr)+0.01);                       // vertex position part
372   sigmaAP0 +=  0.5*sigmaAP1*sigmaAP1;                              
373   return (sigmaAP0>0)? TMath::Sqrt(sigmaAP0):100;
374 }
375
376 Double_t AliESDv0::GetEffectiveSigmaD0(){
377   //
378   // minimax - effective Sigma parameterization 
379   // p12 effective curvature and v0 radius postion used as parameters  
380   //  
381   Double_t p12 = TMath::Sqrt(fParamP.GetParameter()[4]*fParamP.GetParameter()[4]+
382                              fParamN.GetParameter()[4]*fParamN.GetParameter()[4]);
383   Double_t sigmaED0= TMath::Max(TMath::Sqrt(fRr)-fgkParams.fPSigmaRminDE,0.0)*fgkParams.fPSigmaCoefDE*p12*p12;
384   sigmaED0*= sigmaED0;
385   sigmaED0*= sigmaED0;
386   sigmaED0 = TMath::Sqrt(sigmaED0+fgkParams.fPSigmaOffsetDE*fgkParams.fPSigmaOffsetDE);
387   return (sigmaED0<fgkParams.fPSigmaMaxDE) ? sigmaED0: fgkParams.fPSigmaMaxDE;
388 }
389
390
391 Double_t AliESDv0::GetEffectiveSigmaAP0(){
392   //
393   // effective Sigma parameterization of point angle resolution 
394   //
395   Double_t p12 = TMath::Sqrt(fParamP.GetParameter()[4]*fParamP.GetParameter()[4]+
396                              fParamN.GetParameter()[4]*fParamN.GetParameter()[4]);
397   Double_t sigmaAPE= fgkParams.fPSigmaBase0APE;
398   sigmaAPE+= fgkParams.fPSigmaR0APE/(fgkParams.fPSigmaR1APE+fRr);
399   sigmaAPE*= (fgkParams.fPSigmaP0APE+fgkParams.fPSigmaP1APE*p12);
400   sigmaAPE = TMath::Min(sigmaAPE,fgkParams.fPSigmaMaxAPE);
401   return sigmaAPE;
402 }
403
404
405 Double_t  AliESDv0::GetMinimaxSigmaAP0(){
406   //
407   // calculate mini-max effective sigma of point angle resolution
408   //
409   //compv0->fTree->SetAlias("SigmaAP2","max(min((SigmaAP0+SigmaAPE0)*0.5,1.5*SigmaAPE0),0.5*SigmaAPE0+0.003)");
410   Double_t    effectiveSigma = GetEffectiveSigmaAP0();
411   Double_t    sigmaMMAP = 0.5*(GetSigmaAP0()+effectiveSigma);
412   sigmaMMAP  = TMath::Min(sigmaMMAP, fgkParams.fPMaxFractionAP0*effectiveSigma);
413   sigmaMMAP  = TMath::Max(sigmaMMAP, fgkParams.fPMinFractionAP0*effectiveSigma+fgkParams.fPMinAP0);
414   return sigmaMMAP;
415 }
416 Double_t  AliESDv0::GetMinimaxSigmaD0(){
417   //
418   // calculate mini-max sigma of dca resolution
419   // 
420   //compv0->fTree->SetAlias("SigmaD2","max(min((SigmaD0+SigmaDE0)*0.5,1.5*SigmaDE0),0.5*SigmaDE0)");
421   Double_t    effectiveSigma = GetEffectiveSigmaD0();
422   Double_t    sigmaMMD0 = 0.5*(GetSigmaD0()+effectiveSigma);
423   sigmaMMD0  = TMath::Min(sigmaMMD0, fgkParams.fPMaxFractionD0*effectiveSigma);
424   sigmaMMD0  = TMath::Max(sigmaMMD0, fgkParams.fPMinFractionD0*effectiveSigma+fgkParams.fPMinD0);
425   return sigmaMMD0;
426 }
427
428
429 Double_t AliESDv0::GetLikelihoodAP(Int_t mode0, Int_t mode1){
430   //
431   // get likelihood for point angle
432   //
433   Double_t sigmaAP = 0.007;            //default sigma
434   switch (mode0){
435   case 0:
436     sigmaAP = GetSigmaAP0();           // mode 0  - covariance matrix estimates used 
437     break;
438   case 1:
439     sigmaAP = GetEffectiveSigmaAP0();  // mode 1 - effective sigma used
440     break;
441   case 2:
442     sigmaAP = GetMinimaxSigmaAP0();    // mode 2 - minimax sigma
443     break;
444   }
445   Double_t apNorm = TMath::Min(TMath::ACos(fPointAngle)/sigmaAP,50.);  
446   //normalized point angle, restricted - because of overflow problems in Exp
447   Double_t likelihood = 0;
448   switch(mode1){
449   case 0:
450     likelihood = TMath::Exp(-0.5*apNorm*apNorm);   
451     // one component
452     break;
453   case 1:
454     likelihood = (TMath::Exp(-0.5*apNorm*apNorm)+0.5* TMath::Exp(-0.25*apNorm*apNorm))/1.5;
455     // two components
456     break;
457   case 2:
458     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;
459     // three components
460     break;
461   }
462   return likelihood;
463 }
464
465 Double_t AliESDv0::GetLikelihoodD(Int_t mode0, Int_t mode1){
466   //
467   // get likelihood for DCA
468   //
469   Double_t sigmaD = 0.03;            //default sigma
470   switch (mode0){
471   case 0:
472     sigmaD = GetSigmaD0();           // mode 0  - covariance matrix estimates used 
473     break;
474   case 1:
475     sigmaD = GetEffectiveSigmaD0();  // mode 1 - effective sigma used
476     break;
477   case 2:
478     sigmaD = GetMinimaxSigmaD0();    // mode 2 - minimax sigma
479     break;
480   }
481
482   //Bo:  Double_t dNorm = TMath::Min(fDist2/sigmaD,50.);
483   Double_t dNorm = TMath::Min(fDcaV0Daughters/sigmaD,50.);//Bo:
484   //normalized point angle, restricted - because of overflow problems in Exp
485   Double_t likelihood = 0;
486   switch(mode1){
487   case 0:
488     likelihood = TMath::Exp(-2.*dNorm);   
489     // one component
490     break;
491   case 1:
492     likelihood = (TMath::Exp(-2.*dNorm)+0.5* TMath::Exp(-dNorm))/1.5;
493     // two components
494     break;
495   case 2:
496     likelihood = (TMath::Exp(-2.*dNorm)+0.5* TMath::Exp(-dNorm)+0.25*TMath::Exp(-0.5*dNorm))/1.75;
497     // three components
498     break;
499   }
500   return likelihood;
501
502 }
503
504 Double_t AliESDv0::GetLikelihoodC(Int_t mode0, Int_t /*mode1*/){
505   //
506   // get likelihood for Causality
507   // !!!  Causality variables defined in AliITStrackerMI !!! 
508   //      when more information was available
509   //  
510   Double_t likelihood = 0.5;
511   Double_t minCausal  = TMath::Min(fCausality[0],fCausality[1]);
512   Double_t maxCausal  = TMath::Max(fCausality[0],fCausality[1]);
513   //  minCausal           = TMath::Max(minCausal,0.5*maxCausal);
514   //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");
515   
516   switch(mode0){
517   case 0:
518     //normalization 
519     likelihood = TMath::Power((1.05-2*(0.8-TMath::Exp(-maxCausal))),4.);
520     break;
521   case 1:
522     likelihood = TMath::Power(1.05-(2*(0.8-TMath::Exp(-maxCausal))+(2*(0.8-TMath::Exp(-minCausal))))*0.5,4.);
523     break;
524   }
525   return likelihood;
526   
527 }
528
529 void AliESDv0::SetCausality(Float_t pb0, Float_t pb1, Float_t pa0, Float_t pa1)
530 {
531   //
532   // set probabilities
533   //
534   fCausality[0] = pb0;     // probability - track 0 exist before vertex
535   fCausality[1] = pb1;     // probability - track 1 exist before vertex
536   fCausality[2] = pa0;     // probability - track 0 exist close after vertex
537   fCausality[3] = pa1;     // probability - track 1 exist close after vertex
538 }
539 void  AliESDv0::SetClusters(Int_t *clp, Int_t *clm)
540 {
541   //
542   // Set its clusters indexes
543   //
544   for (Int_t i=0;i<6;i++) fClusters[0][i] = clp[i]; 
545   for (Int_t i=0;i<6;i++) fClusters[1][i] = clm[i]; 
546 }