1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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 //-------------------------------------------------------------------------
29 #include <TDatabasePDG.h>
30 #include <TParticlePDG.h>
35 #include "AliESDV0Params.h"
39 const AliESDV0Params AliESDv0::fgkParams;
41 AliESDv0::AliESDv0() :
45 fEffMass(TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass()),
63 //--------------------------------------------------------------------
64 // Default constructor (K0s)
65 //--------------------------------------------------------------------
67 for (Int_t i=0; i<3; i++) {
73 for (Int_t i=0; i<6; i++) {
77 for (Int_t i=0;i<6;i++){fClusters[0][i]=0; fClusters[1][i]=0;}
78 fNormDCAPrim[0]=fNormDCAPrim[1]=0;
79 for (Int_t i=0;i<3;i++){fAngle[i]=0;}
80 for (Int_t i=0;i<4;i++){fCausality[i]=0;}
83 AliESDv0::AliESDv0(const AliESDv0& v0) :
87 fEffMass(v0.fEffMass),
88 fDcaV0Daughters(v0.fDcaV0Daughters),
91 fDistSigma(v0.fDistSigma),
92 fChi2Before(v0.fChi2Before),
93 fChi2After(v0.fChi2After),
94 fPointAngleFi(v0.fPointAngleFi),
95 fPointAngleTh(v0.fPointAngleTh),
96 fPointAngle(v0.fPointAngle),
97 fPdgCode(v0.fPdgCode),
101 fNBefore(v0.fNBefore),
103 fOnFlyStatus(v0.fOnFlyStatus)
105 //--------------------------------------------------------------------
106 // The copy constructor
107 //--------------------------------------------------------------------
109 for (int i=0; i<3; i++) {
110 fPos[i] = v0.fPos[i];
111 fNmom[i] = v0.fNmom[i];
112 fPmom[i] = v0.fPmom[i];
114 for (int i=0; i<6; i++) {
115 fPosCov[i] = v0.fPosCov[i];
118 for (Int_t i=0; i<2; i++) {
119 fNormDCAPrim[i]=v0.fNormDCAPrim[i];
121 for (Int_t i=0;i<6;i++){
122 fClusters[0][i]=v0.fClusters[0][i];
123 fClusters[1][i]=v0.fClusters[1][i];
125 for (Int_t i=0;i<3;i++){
126 fAngle[i]=v0.fAngle[i];
128 for (Int_t i=0;i<4;i++){fCausality[i]=v0.fCausality[i];}
132 AliESDv0::AliESDv0(const AliExternalTrackParam &t1, Int_t i1,
133 const AliExternalTrackParam &t2, Int_t i2) :
137 fEffMass(TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass()),
155 //--------------------------------------------------------------------
156 // Main constructor (K0s)
157 //--------------------------------------------------------------------
159 //Make sure the daughters are ordered (needed for the on-the-fly V0s)
160 Short_t cN=t1.Charge(), cP=t2.Charge();
161 if ((cN>0) && (cN != cP)) {
162 fParamN.~AliExternalTrackParam();
163 new (&fParamN) AliExternalTrackParam(t2);
164 fParamP.~AliExternalTrackParam();
165 new (&fParamP) AliExternalTrackParam(t1);
168 for (Int_t i=0; i<6; i++) {
172 //Trivial estimation of the vertex parameters
173 Double_t alpha=t1.GetAlpha(), cs=TMath::Cos(alpha), sn=TMath::Sin(alpha);
176 Double_t px1=tmp[0], py1=tmp[1], pz1=tmp[2];
178 Double_t x1=tmp[0], y1=tmp[1], z1=tmp[2];
179 const Double_t ss=0.0005*0.0005;//a kind of a residual misalignment precision
180 Double_t sx1=sn*sn*t1.GetSigmaY2()+ss, sy1=cs*cs*t1.GetSigmaY2()+ss;
183 alpha=t2.GetAlpha(); cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
185 Double_t px2=tmp[0], py2=tmp[1], pz2=tmp[2];
187 Double_t x2=tmp[0], y2=tmp[1], z2=tmp[2];
188 Double_t sx2=sn*sn*t2.GetSigmaY2()+ss, sy2=cs*cs*t2.GetSigmaY2()+ss;
190 Double_t sz1=t1.GetSigmaZ2(), sz2=t2.GetSigmaZ2();
191 Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
192 Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
193 Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
194 fPos[0]=wx1*x1 + wx2*x2; fPos[1]=wy1*y1 + wy2*y2; fPos[2]=wz1*z1 + wz2*z2;
196 //fPos[0]=0.5*(x1+x2); fPos[1]=0.5*(y1+y2); fPos[2]=0.5*(z1+z2);
197 fNmom[0]=px1; fNmom[1]=py1; fNmom[2]=pz1;
198 fPmom[0]=px2; fPmom[1]=py2; fPmom[2]=pz2;
200 for (Int_t i=0;i<6;i++){fClusters[0][i]=0; fClusters[1][i]=0;}
201 fNormDCAPrim[0]=fNormDCAPrim[1]=0;
202 for (Int_t i=0;i<3;i++){fAngle[i]=0;}
203 for (Int_t i=0;i<4;i++){fCausality[i]=0;}
206 AliESDv0& AliESDv0::operator=(const AliESDv0 &v0)
208 //--------------------------------------------------------------------
209 // The assignment operator
210 //--------------------------------------------------------------------
212 if(this==&v0)return *this;
213 AliVParticle::operator=(v0);
214 fParamN = v0.fParamN;
215 fParamP = v0.fParamP;
216 fEffMass = v0.fEffMass;
217 fDcaV0Daughters = v0.fDcaV0Daughters;
218 fChi2V0 = v0.fChi2V0;
220 fDistSigma = v0.fDistSigma;
221 fChi2Before = v0.fChi2Before;
222 fChi2After = v0.fChi2After;
223 fPointAngleFi = v0.fPointAngleFi;
224 fPointAngleTh = v0.fPointAngleTh;
225 fPointAngle = v0.fPointAngle;
226 fPdgCode = v0.fPdgCode;
229 fStatus = v0.fStatus;
230 fNBefore = v0.fNBefore;
231 fNAfter = v0.fNAfter;
232 fOnFlyStatus = v0.fOnFlyStatus;
234 for (int i=0; i<3; i++) {
235 fPos[i] = v0.fPos[i];
236 fNmom[i] = v0.fNmom[i];
237 fPmom[i] = v0.fPmom[i];
239 for (int i=0; i<6; i++) {
240 fPosCov[i] = v0.fPosCov[i];
242 for (Int_t i=0; i<2; i++) {
243 fNormDCAPrim[i]=v0.fNormDCAPrim[i];
245 for (Int_t i=0;i<6;i++){
246 fClusters[0][i]=v0.fClusters[0][i];
247 fClusters[1][i]=v0.fClusters[1][i];
249 for (Int_t i=0;i<3;i++){
250 fAngle[i]=v0.fAngle[i];
252 for (Int_t i=0;i<4;i++){fCausality[i]=v0.fCausality[i];}
257 void AliESDv0::Copy(TObject& obj) const {
259 // this overwrites the virtual TOBject::Copy()
260 // to allow run time copying without casting
263 if(this==&obj)return;
264 AliESDv0 *robj = dynamic_cast<AliESDv0*>(&obj);
265 if(!robj)return; // not an aliesesv0
269 AliESDv0::~AliESDv0(){
270 //--------------------------------------------------------------------
272 //--------------------------------------------------------------------
275 // Start with AliVParticle functions
276 Double_t AliESDv0::E() const {
277 //--------------------------------------------------------------------
278 // This gives the energy assuming the ChangeMassHypothesis was called
279 //--------------------------------------------------------------------
283 Double_t AliESDv0::Y() const {
284 //--------------------------------------------------------------------
285 // This gives the energy assuming the ChangeMassHypothesis was called
286 //--------------------------------------------------------------------
290 // Then extend AliVParticle functions
291 Double_t AliESDv0::E(Int_t pdg) const {
292 //--------------------------------------------------------------------
293 // This gives the energy with the particle hypothesis as argument
294 //--------------------------------------------------------------------
295 Double_t mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
296 return TMath::Sqrt(mass*mass+P()*P());
299 Double_t AliESDv0::Y(Int_t pdg) const {
300 //--------------------------------------------------------------------
301 // This gives the rapidity with the particle hypothesis as argument
302 //--------------------------------------------------------------------
303 return 0.5*TMath::Log((E(pdg)+Pz())/(E(pdg)-Pz()+1.e-13));
306 // Now the functions for analysis consistency
307 Double_t AliESDv0::RapK0Short() const {
308 //--------------------------------------------------------------------
309 // This gives the pseudorapidity assuming a K0s particle
310 //--------------------------------------------------------------------
314 Double_t AliESDv0::RapLambda() const {
315 //--------------------------------------------------------------------
316 // This gives the pseudorapidity assuming a (Anti) Lambda particle
317 //--------------------------------------------------------------------
321 Double_t AliESDv0::AlphaV0() const {
322 //--------------------------------------------------------------------
323 // This gives the Armenteros-Podolanski alpha
324 //--------------------------------------------------------------------
325 TVector3 momNeg(fNmom[0],fNmom[1],fNmom[2]);
326 TVector3 momPos(fPmom[0],fPmom[1],fPmom[2]);
327 TVector3 momTot(Px(),Py(),Pz());
329 Double_t lQlNeg = momNeg.Dot(momTot)/momTot.Mag();
330 Double_t lQlPos = momPos.Dot(momTot)/momTot.Mag();
332 //return 1.-2./(1.+lQlNeg/lQlPos);
333 return (lQlPos - lQlNeg)/(lQlPos + lQlNeg);
336 Double_t AliESDv0::PtArmV0() const {
337 //--------------------------------------------------------------------
338 // This gives the Armenteros-Podolanski ptarm
339 //--------------------------------------------------------------------
340 TVector3 momNeg(fNmom[0],fNmom[1],fNmom[2]);
341 TVector3 momTot(Px(),Py(),Pz());
343 return momNeg.Perp(momTot);
346 // Eventually the older functions
347 Double_t AliESDv0::ChangeMassHypothesis(Int_t code) {
348 //--------------------------------------------------------------------
349 // This function changes the mass hypothesis for this V0
350 // and returns the "kinematical quality" of this hypothesis
351 //--------------------------------------------------------------------
353 Double_t piMass=TDatabasePDG::Instance()->GetParticle(kPiPlus)->Mass();
355 Double_t prMass=TDatabasePDG::Instance()->GetParticle(kProton)->Mass();
357 Double_t k0Mass=TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass();
359 Double_t l0Mass=TDatabasePDG::Instance()->GetParticle(kLambda0)->Mass();
361 Double_t nmass=piMass, pmass=piMass, mass=k0Mass, ps=0.206;
367 nmass=piMass; pmass=prMass; mass=l0Mass; ps=0.101; break;
369 pmass=piMass; nmass=prMass; mass=l0Mass; ps=0.101; break;
373 AliError("invalide PDG code ! Assuming K0s...");
378 Double_t pxn=fNmom[0], pyn=fNmom[1], pzn=fNmom[2];
379 Double_t pxp=fPmom[0], pyp=fPmom[1], pzp=fPmom[2];
381 Double_t en=TMath::Sqrt(nmass*nmass + pxn*pxn + pyn*pyn + pzn*pzn);
382 Double_t ep=TMath::Sqrt(pmass*pmass + pxp*pxp + pyp*pyp + pzp*pzp);
383 Double_t pxl=pxn+pxp, pyl=pyn+pyp, pzl=pzn+pzp;
384 Double_t pl=TMath::Sqrt(pxl*pxl + pyl*pyl + pzl*pzl);
386 fEffMass=TMath::Sqrt((en+ep)*(en+ep)-pl*pl);
388 Double_t beta=pl/(en+ep);
389 Double_t pln=(pxn*pxl + pyn*pyl + pzn*pzl)/pl;
390 Double_t plp=(pxp*pxl + pyp*pyl + pzp*pzl)/pl;
392 Double_t pt2=pxp*pxp + pyp*pyp + pzp*pzp - plp*plp;
394 Double_t a=(plp-pln)/(plp+pln);
395 a -= (pmass*pmass-nmass*nmass)/(mass*mass);
396 a = 0.25*beta*beta*mass*mass*a*a + pt2;
402 void AliESDv0::GetPxPyPz(Double_t &px, Double_t &py, Double_t &pz) const {
403 //--------------------------------------------------------------------
404 // This function returns V0's momentum (global)
405 //--------------------------------------------------------------------
406 px=fNmom[0]+fPmom[0];
407 py=fNmom[1]+fPmom[1];
408 pz=fNmom[2]+fPmom[2];
411 void AliESDv0::GetXYZ(Double_t &x, Double_t &y, Double_t &z) const {
412 //--------------------------------------------------------------------
413 // This function returns V0's position (global)
414 //--------------------------------------------------------------------
420 Float_t AliESDv0::GetD(Double_t x0, Double_t y0, Double_t z0) const {
421 //--------------------------------------------------------------------
422 // This function returns V0's impact parameter
423 //--------------------------------------------------------------------
424 Double_t x=fPos[0],y=fPos[1],z=fPos[2];
425 Double_t px=fNmom[0]+fPmom[0];
426 Double_t py=fNmom[1]+fPmom[1];
427 Double_t pz=fNmom[2]+fPmom[2];
429 Double_t dx=(y0-y)*pz - (z0-z)*py;
430 Double_t dy=(x0-x)*pz - (z0-z)*px;
431 Double_t dz=(x0-x)*py - (y0-y)*px;
432 Double_t d=TMath::Sqrt((dx*dx+dy*dy+dz*dz)/(px*px+py*py+pz*pz));
436 Float_t AliESDv0::GetV0CosineOfPointingAngle(Double_t refPointX, Double_t refPointY, Double_t refPointZ) const {
437 // calculates the pointing angle of the V0 wrt a reference point
439 Double_t momV0[3]; //momentum of the V0
440 GetPxPyPz(momV0[0],momV0[1],momV0[2]);
442 Double_t deltaPos[3]; //vector between the reference point and the V0 vertex
443 deltaPos[0] = fPos[0] - refPointX;
444 deltaPos[1] = fPos[1] - refPointY;
445 deltaPos[2] = fPos[2] - refPointZ;
447 Double_t momV02 = momV0[0]*momV0[0] + momV0[1]*momV0[1] + momV0[2]*momV0[2];
448 Double_t deltaPos2 = deltaPos[0]*deltaPos[0] + deltaPos[1]*deltaPos[1] + deltaPos[2]*deltaPos[2];
450 Double_t cosinePointingAngle = (deltaPos[0]*momV0[0] +
451 deltaPos[1]*momV0[1] +
452 deltaPos[2]*momV0[2] ) /
453 TMath::Sqrt(momV02 * deltaPos2);
455 return cosinePointingAngle;
459 // **** The following functions need to be revised
461 void AliESDv0::GetPosCov(Double_t cov[6]) const {
463 for (Int_t i=0; i<6; ++i) cov[i] = fPosCov[i];
467 Double_t AliESDv0::GetSigmaY(){
469 // return sigmay in y at vertex position using covariance matrix
471 const Double_t * cp = fParamP.GetCovariance();
472 const Double_t * cm = fParamN.GetCovariance();
473 Double_t sigmay = cp[0]+cm[0]+ cp[5]*(fParamP.GetX()-fRr)*(fParamP.GetX()-fRr)+ cm[5]*(fParamN.GetX()-fRr)*(fParamN.GetX()-fRr);
474 return (sigmay>0) ? TMath::Sqrt(sigmay):100;
477 Double_t AliESDv0::GetSigmaZ(){
479 // return sigmay in y at vertex position using covariance matrix
481 const Double_t * cp = fParamP.GetCovariance();
482 const Double_t * cm = fParamN.GetCovariance();
483 Double_t sigmaz = cp[2]+cm[2]+ cp[9]*(fParamP.GetX()-fRr)*(fParamP.GetX()-fRr)+ cm[9]*(fParamN.GetX()-fRr)*(fParamN.GetX()-fRr);
484 return (sigmaz>0) ? TMath::Sqrt(sigmaz):100;
487 Double_t AliESDv0::GetSigmaD0(){
489 // Sigma parameterization using covariance matrix
491 // sigma of distance between two tracks in vertex position
492 // sigma of DCA is proportianal to sigmaD0
493 // factor 2 difference is explained by the fact that the DCA is calculated at the position
494 // where the tracks as closest together ( not exact position of the vertex)
496 const Double_t * cp = fParamP.GetCovariance();
497 const Double_t * cm = fParamN.GetCovariance();
498 Double_t sigmaD0 = cp[0]+cm[0]+cp[2]+cm[2]+fgkParams.fPSigmaOffsetD0*fgkParams.fPSigmaOffsetD0;
499 sigmaD0 += ((fParamP.GetX()-fRr)*(fParamP.GetX()-fRr))*(cp[5]+cp[9]);
500 sigmaD0 += ((fParamN.GetX()-fRr)*(fParamN.GetX()-fRr))*(cm[5]+cm[9]);
501 return (sigmaD0>0)? TMath::Sqrt(sigmaD0):100;
505 Double_t AliESDv0::GetSigmaAP0(){
507 //Sigma parameterization using covariance matrices
509 Double_t prec = TMath::Sqrt((fNmom[0]+fPmom[0])*(fNmom[0]+fPmom[0])
510 +(fNmom[1]+fPmom[1])*(fNmom[1]+fPmom[1])
511 +(fNmom[2]+fPmom[2])*(fNmom[2]+fPmom[2]));
512 Double_t normp = TMath::Sqrt(fPmom[0]*fPmom[0]+fPmom[1]*fPmom[1]+fPmom[2]*fPmom[2])/prec; // fraction of the momenta
513 Double_t normm = TMath::Sqrt(fNmom[0]*fNmom[0]+fNmom[1]*fNmom[1]+fNmom[2]*fNmom[2])/prec;
514 const Double_t * cp = fParamP.GetCovariance();
515 const Double_t * cm = fParamN.GetCovariance();
516 Double_t sigmaAP0 = fgkParams.fPSigmaOffsetAP0*fgkParams.fPSigmaOffsetAP0; // minimal part
517 sigmaAP0 += (cp[5]+cp[9])*(normp*normp)+(cm[5]+cm[9])*(normm*normm); // angular resolution part
518 Double_t sigmaAP1 = GetSigmaD0()/(TMath::Abs(fRr)+0.01); // vertex position part
519 sigmaAP0 += 0.5*sigmaAP1*sigmaAP1;
520 return (sigmaAP0>0)? TMath::Sqrt(sigmaAP0):100;
523 Double_t AliESDv0::GetEffectiveSigmaD0(){
525 // minimax - effective Sigma parameterization
526 // p12 effective curvature and v0 radius postion used as parameters
528 Double_t p12 = TMath::Sqrt(fParamP.GetParameter()[4]*fParamP.GetParameter()[4]+
529 fParamN.GetParameter()[4]*fParamN.GetParameter()[4]);
530 Double_t sigmaED0= TMath::Max(TMath::Sqrt(fRr)-fgkParams.fPSigmaRminDE,0.0)*fgkParams.fPSigmaCoefDE*p12*p12;
533 sigmaED0 = TMath::Sqrt(sigmaED0+fgkParams.fPSigmaOffsetDE*fgkParams.fPSigmaOffsetDE);
534 return (sigmaED0<fgkParams.fPSigmaMaxDE) ? sigmaED0: fgkParams.fPSigmaMaxDE;
538 Double_t AliESDv0::GetEffectiveSigmaAP0(){
540 // effective Sigma parameterization of point angle resolution
542 Double_t p12 = TMath::Sqrt(fParamP.GetParameter()[4]*fParamP.GetParameter()[4]+
543 fParamN.GetParameter()[4]*fParamN.GetParameter()[4]);
544 Double_t sigmaAPE= fgkParams.fPSigmaBase0APE;
545 sigmaAPE+= fgkParams.fPSigmaR0APE/(fgkParams.fPSigmaR1APE+fRr);
546 sigmaAPE*= (fgkParams.fPSigmaP0APE+fgkParams.fPSigmaP1APE*p12);
547 sigmaAPE = TMath::Min(sigmaAPE,fgkParams.fPSigmaMaxAPE);
552 Double_t AliESDv0::GetMinimaxSigmaAP0(){
554 // calculate mini-max effective sigma of point angle resolution
556 //compv0->fTree->SetAlias("SigmaAP2","max(min((SigmaAP0+SigmaAPE0)*0.5,1.5*SigmaAPE0),0.5*SigmaAPE0+0.003)");
557 Double_t effectiveSigma = GetEffectiveSigmaAP0();
558 Double_t sigmaMMAP = 0.5*(GetSigmaAP0()+effectiveSigma);
559 sigmaMMAP = TMath::Min(sigmaMMAP, fgkParams.fPMaxFractionAP0*effectiveSigma);
560 sigmaMMAP = TMath::Max(sigmaMMAP, fgkParams.fPMinFractionAP0*effectiveSigma+fgkParams.fPMinAP0);
563 Double_t AliESDv0::GetMinimaxSigmaD0(){
565 // calculate mini-max sigma of dca resolution
567 //compv0->fTree->SetAlias("SigmaD2","max(min((SigmaD0+SigmaDE0)*0.5,1.5*SigmaDE0),0.5*SigmaDE0)");
568 Double_t effectiveSigma = GetEffectiveSigmaD0();
569 Double_t sigmaMMD0 = 0.5*(GetSigmaD0()+effectiveSigma);
570 sigmaMMD0 = TMath::Min(sigmaMMD0, fgkParams.fPMaxFractionD0*effectiveSigma);
571 sigmaMMD0 = TMath::Max(sigmaMMD0, fgkParams.fPMinFractionD0*effectiveSigma+fgkParams.fPMinD0);
576 Double_t AliESDv0::GetLikelihoodAP(Int_t mode0, Int_t mode1){
578 // get likelihood for point angle
580 Double_t sigmaAP = 0.007; //default sigma
583 sigmaAP = GetSigmaAP0(); // mode 0 - covariance matrix estimates used
586 sigmaAP = GetEffectiveSigmaAP0(); // mode 1 - effective sigma used
589 sigmaAP = GetMinimaxSigmaAP0(); // mode 2 - minimax sigma
592 Double_t apNorm = TMath::Min(TMath::ACos(fPointAngle)/sigmaAP,50.);
593 //normalized point angle, restricted - because of overflow problems in Exp
594 Double_t likelihood = 0;
597 likelihood = TMath::Exp(-0.5*apNorm*apNorm);
601 likelihood = (TMath::Exp(-0.5*apNorm*apNorm)+0.5* TMath::Exp(-0.25*apNorm*apNorm))/1.5;
605 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;
612 Double_t AliESDv0::GetLikelihoodD(Int_t mode0, Int_t mode1){
614 // get likelihood for DCA
616 Double_t sigmaD = 0.03; //default sigma
619 sigmaD = GetSigmaD0(); // mode 0 - covariance matrix estimates used
622 sigmaD = GetEffectiveSigmaD0(); // mode 1 - effective sigma used
625 sigmaD = GetMinimaxSigmaD0(); // mode 2 - minimax sigma
629 //Bo: Double_t dNorm = TMath::Min(fDist2/sigmaD,50.);
630 Double_t dNorm = TMath::Min(fDcaV0Daughters/sigmaD,50.);//Bo:
631 //normalized point angle, restricted - because of overflow problems in Exp
632 Double_t likelihood = 0;
635 likelihood = TMath::Exp(-2.*dNorm);
639 likelihood = (TMath::Exp(-2.*dNorm)+0.5* TMath::Exp(-dNorm))/1.5;
643 likelihood = (TMath::Exp(-2.*dNorm)+0.5* TMath::Exp(-dNorm)+0.25*TMath::Exp(-0.5*dNorm))/1.75;
651 Double_t AliESDv0::GetLikelihoodC(Int_t mode0, Int_t /*mode1*/) const {
653 // get likelihood for Causality
654 // !!! Causality variables defined in AliITStrackerMI !!!
655 // when more information was available
657 Double_t likelihood = 0.5;
658 Double_t minCausal = TMath::Min(fCausality[0],fCausality[1]);
659 Double_t maxCausal = TMath::Max(fCausality[0],fCausality[1]);
660 // minCausal = TMath::Max(minCausal,0.5*maxCausal);
661 //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");
666 likelihood = TMath::Power((1.05-2*(0.8-TMath::Exp(-maxCausal))),4.);
669 likelihood = TMath::Power(1.05-(2*(0.8-TMath::Exp(-maxCausal))+(2*(0.8-TMath::Exp(-minCausal))))*0.5,4.);
676 void AliESDv0::SetCausality(Float_t pb0, Float_t pb1, Float_t pa0, Float_t pa1)
681 fCausality[0] = pb0; // probability - track 0 exist before vertex
682 fCausality[1] = pb1; // probability - track 1 exist before vertex
683 fCausality[2] = pa0; // probability - track 0 exist close after vertex
684 fCausality[3] = pa1; // probability - track 1 exist close after vertex
686 void AliESDv0::SetClusters(const Int_t *clp, const Int_t *clm)
689 // Set its clusters indexes
691 for (Int_t i=0;i<6;i++) fClusters[0][i] = clp[i];
692 for (Int_t i=0;i<6;i++) fClusters[1][i] = clm[i];
695 Double_t AliESDv0::GetEffMass(UInt_t p1, UInt_t p2) const{
697 // calculate effective mass
699 const Double_t kpmass[5] = {TDatabasePDG::Instance()->GetParticle(kElectron)->Mass(),
700 TDatabasePDG::Instance()->GetParticle(kMuonMinus)->Mass(),
701 TDatabasePDG::Instance()->GetParticle(kPiPlus)->Mass(),
702 TDatabasePDG::Instance()->GetParticle(kKPlus)->Mass(),
703 TDatabasePDG::Instance()->GetParticle(kProton)->Mass()};
707 Double_t mass1 = kpmass[p1];
708 Double_t mass2 = kpmass[p2];
709 const Double_t *m1 = fPmom;
710 const Double_t *m2 = fNmom;
712 //if (fRP[p1]+fRM[p2]<fRP[p2]+fRM[p1]){
717 Double_t e1 = TMath::Sqrt(mass1*mass1+
721 Double_t e2 = TMath::Sqrt(mass2*mass2+
726 (m2[0]+m1[0])*(m2[0]+m1[0])+
727 (m2[1]+m1[1])*(m2[1]+m1[1])+
728 (m2[2]+m1[2])*(m2[2]+m1[2]);
730 mass = (e1+e2)*(e1+e2)-mass;
731 if (mass < 0.) mass = 0.;
732 return (TMath::Sqrt(mass));
734 if(p1>4 || p2>4) return -1;
735 Double_t e12 = kpmass[p1]*kpmass[p1]+fPmom[0]*fPmom[0]+fPmom[1]*fPmom[1]+fPmom[2]*fPmom[2];
736 Double_t e22 = kpmass[p2]*kpmass[p2]+fNmom[0]*fNmom[0]+fNmom[1]*fNmom[1]+fNmom[2]*fNmom[2];
737 Double_t cmass = TMath::Sqrt(TMath::Max(kpmass[p1]*kpmass[p1]+kpmass[p2]*kpmass[p2]
738 +2.*(TMath::Sqrt(e12*e22)-fPmom[0]*fNmom[0]-fPmom[1]*fNmom[1]-fPmom[2]*fNmom[2]),0.));