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"
36 #include "AliKFParticle.h"
37 #include "AliKFVertex.h"
41 const AliESDV0Params AliESDv0::fgkParams;
43 AliESDv0::AliESDv0() :
47 fEffMass(TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass()),
65 //--------------------------------------------------------------------
66 // Default constructor (K0s)
67 //--------------------------------------------------------------------
69 for (Int_t i=0; i<3; i++) {
75 for (Int_t i=0; i<6; i++) {
79 for (Int_t i=0;i<6;i++){fClusters[0][i]=0; fClusters[1][i]=0;}
80 fNormDCAPrim[0]=fNormDCAPrim[1]=0;
81 for (Int_t i=0;i<3;i++){fAngle[i]=0;}
82 for (Int_t i=0;i<4;i++){fCausality[i]=0;}
85 AliESDv0::AliESDv0(const AliESDv0& v0) :
89 fEffMass(v0.fEffMass),
90 fDcaV0Daughters(v0.fDcaV0Daughters),
93 fDistSigma(v0.fDistSigma),
94 fChi2Before(v0.fChi2Before),
95 fChi2After(v0.fChi2After),
96 fPointAngleFi(v0.fPointAngleFi),
97 fPointAngleTh(v0.fPointAngleTh),
98 fPointAngle(v0.fPointAngle),
99 fPdgCode(v0.fPdgCode),
103 fNBefore(v0.fNBefore),
105 fOnFlyStatus(v0.fOnFlyStatus)
107 //--------------------------------------------------------------------
108 // The copy constructor
109 //--------------------------------------------------------------------
111 for (int i=0; i<3; i++) {
112 fPos[i] = v0.fPos[i];
113 fNmom[i] = v0.fNmom[i];
114 fPmom[i] = v0.fPmom[i];
116 for (int i=0; i<6; i++) {
117 fPosCov[i] = v0.fPosCov[i];
120 for (Int_t i=0; i<2; i++) {
121 fNormDCAPrim[i]=v0.fNormDCAPrim[i];
123 for (Int_t i=0;i<6;i++){
124 fClusters[0][i]=v0.fClusters[0][i];
125 fClusters[1][i]=v0.fClusters[1][i];
127 for (Int_t i=0;i<3;i++){
128 fAngle[i]=v0.fAngle[i];
130 for (Int_t i=0;i<4;i++){fCausality[i]=v0.fCausality[i];}
134 AliESDv0::AliESDv0(const AliExternalTrackParam &t1, Int_t i1,
135 const AliExternalTrackParam &t2, Int_t i2) :
139 fEffMass(TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass()),
157 //--------------------------------------------------------------------
158 // Main constructor (K0s)
159 //--------------------------------------------------------------------
161 //Make sure the daughters are ordered (needed for the on-the-fly V0s)
162 Short_t cN=t1.Charge(), cP=t2.Charge();
163 if ((cN>0) && (cN != cP)) {
164 fParamN.~AliExternalTrackParam();
165 new (&fParamN) AliExternalTrackParam(t2);
166 fParamP.~AliExternalTrackParam();
167 new (&fParamP) AliExternalTrackParam(t1);
174 for (Int_t i=0; i<6; i++) {
178 //Trivial estimation of the vertex parameters
179 Double_t alpha=t1.GetAlpha(), cs=TMath::Cos(alpha), sn=TMath::Sin(alpha);
182 Double_t px1=tmp[0], py1=tmp[1], pz1=tmp[2];
184 Double_t x1=tmp[0], y1=tmp[1], z1=tmp[2];
185 const Double_t ss=0.0005*0.0005;//a kind of a residual misalignment precision
186 Double_t sx1=sn*sn*t1.GetSigmaY2()+ss, sy1=cs*cs*t1.GetSigmaY2()+ss;
189 alpha=t2.GetAlpha(); cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
191 Double_t px2=tmp[0], py2=tmp[1], pz2=tmp[2];
193 Double_t x2=tmp[0], y2=tmp[1], z2=tmp[2];
194 Double_t sx2=sn*sn*t2.GetSigmaY2()+ss, sy2=cs*cs*t2.GetSigmaY2()+ss;
196 Double_t sz1=t1.GetSigmaZ2(), sz2=t2.GetSigmaZ2();
197 Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
198 Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
199 Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
200 fPos[0]=wx1*x1 + wx2*x2; fPos[1]=wy1*y1 + wy2*y2; fPos[2]=wz1*z1 + wz2*z2;
202 //fPos[0]=0.5*(x1+x2); fPos[1]=0.5*(y1+y2); fPos[2]=0.5*(z1+z2);
203 fNmom[0]=px1; fNmom[1]=py1; fNmom[2]=pz1;
204 fPmom[0]=px2; fPmom[1]=py2; fPmom[2]=pz2;
206 for (Int_t i=0;i<6;i++){fClusters[0][i]=0; fClusters[1][i]=0;}
207 fNormDCAPrim[0]=fNormDCAPrim[1]=0;
208 for (Int_t i=0;i<3;i++){fAngle[i]=0;}
209 for (Int_t i=0;i<4;i++){fCausality[i]=0;}
212 AliESDv0& AliESDv0::operator=(const AliESDv0 &v0)
214 //--------------------------------------------------------------------
215 // The assignment operator
216 //--------------------------------------------------------------------
218 if(this==&v0)return *this;
219 AliVParticle::operator=(v0);
220 fParamN = v0.fParamN;
221 fParamP = v0.fParamP;
222 fEffMass = v0.fEffMass;
223 fDcaV0Daughters = v0.fDcaV0Daughters;
224 fChi2V0 = v0.fChi2V0;
226 fDistSigma = v0.fDistSigma;
227 fChi2Before = v0.fChi2Before;
228 fChi2After = v0.fChi2After;
229 fPointAngleFi = v0.fPointAngleFi;
230 fPointAngleTh = v0.fPointAngleTh;
231 fPointAngle = v0.fPointAngle;
232 fPdgCode = v0.fPdgCode;
235 fStatus = v0.fStatus;
236 fNBefore = v0.fNBefore;
237 fNAfter = v0.fNAfter;
238 fOnFlyStatus = v0.fOnFlyStatus;
240 for (int i=0; i<3; i++) {
241 fPos[i] = v0.fPos[i];
242 fNmom[i] = v0.fNmom[i];
243 fPmom[i] = v0.fPmom[i];
245 for (int i=0; i<6; i++) {
246 fPosCov[i] = v0.fPosCov[i];
248 for (Int_t i=0; i<2; i++) {
249 fNormDCAPrim[i]=v0.fNormDCAPrim[i];
251 for (Int_t i=0;i<6;i++){
252 fClusters[0][i]=v0.fClusters[0][i];
253 fClusters[1][i]=v0.fClusters[1][i];
255 for (Int_t i=0;i<3;i++){
256 fAngle[i]=v0.fAngle[i];
258 for (Int_t i=0;i<4;i++){fCausality[i]=v0.fCausality[i];}
263 void AliESDv0::Copy(TObject& obj) const {
265 // this overwrites the virtual TOBject::Copy()
266 // to allow run time copying without casting
269 if(this==&obj)return;
270 AliESDv0 *robj = dynamic_cast<AliESDv0*>(&obj);
271 if(!robj)return; // not an aliesesv0
275 AliESDv0::~AliESDv0(){
276 //--------------------------------------------------------------------
278 //--------------------------------------------------------------------
281 // Start with AliVParticle functions
282 Double_t AliESDv0::E() const {
283 //--------------------------------------------------------------------
284 // This gives the energy assuming the ChangeMassHypothesis was called
285 //--------------------------------------------------------------------
289 Double_t AliESDv0::Y() const {
290 //--------------------------------------------------------------------
291 // This gives the energy assuming the ChangeMassHypothesis was called
292 //--------------------------------------------------------------------
296 // Then extend AliVParticle functions
297 Double_t AliESDv0::E(Int_t pdg) const {
298 //--------------------------------------------------------------------
299 // This gives the energy with the particle hypothesis as argument
300 //--------------------------------------------------------------------
301 Double_t mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
302 return TMath::Sqrt(mass*mass+P()*P());
305 Double_t AliESDv0::Y(Int_t pdg) const {
306 //--------------------------------------------------------------------
307 // This gives the rapidity with the particle hypothesis as argument
308 //--------------------------------------------------------------------
309 return 0.5*TMath::Log((E(pdg)+Pz())/(E(pdg)-Pz()+1.e-13));
312 // Now the functions for analysis consistency
313 Double_t AliESDv0::RapK0Short() const {
314 //--------------------------------------------------------------------
315 // This gives the pseudorapidity assuming a K0s particle
316 //--------------------------------------------------------------------
320 Double_t AliESDv0::RapLambda() const {
321 //--------------------------------------------------------------------
322 // This gives the pseudorapidity assuming a (Anti) Lambda particle
323 //--------------------------------------------------------------------
327 Double_t AliESDv0::AlphaV0() const {
328 //--------------------------------------------------------------------
329 // This gives the Armenteros-Podolanski alpha
330 //--------------------------------------------------------------------
331 TVector3 momNeg(fNmom[0],fNmom[1],fNmom[2]);
332 TVector3 momPos(fPmom[0],fPmom[1],fPmom[2]);
333 TVector3 momTot(Px(),Py(),Pz());
335 Double_t lQlNeg = momNeg.Dot(momTot)/momTot.Mag();
336 Double_t lQlPos = momPos.Dot(momTot)/momTot.Mag();
338 //return 1.-2./(1.+lQlNeg/lQlPos);
339 return (lQlPos - lQlNeg)/(lQlPos + lQlNeg);
342 Double_t AliESDv0::PtArmV0() const {
343 //--------------------------------------------------------------------
344 // This gives the Armenteros-Podolanski ptarm
345 //--------------------------------------------------------------------
346 TVector3 momNeg(fNmom[0],fNmom[1],fNmom[2]);
347 TVector3 momTot(Px(),Py(),Pz());
349 return momNeg.Perp(momTot);
352 // Eventually the older functions
353 Double_t AliESDv0::ChangeMassHypothesis(Int_t code) {
354 //--------------------------------------------------------------------
355 // This function changes the mass hypothesis for this V0
356 // and returns the "kinematical quality" of this hypothesis
357 //--------------------------------------------------------------------
359 Double_t piMass=TDatabasePDG::Instance()->GetParticle(kPiPlus)->Mass();
361 Double_t prMass=TDatabasePDG::Instance()->GetParticle(kProton)->Mass();
363 Double_t k0Mass=TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass();
365 Double_t l0Mass=TDatabasePDG::Instance()->GetParticle(kLambda0)->Mass();
367 Double_t nmass=piMass, pmass=piMass, mass=k0Mass, ps=0.206;
373 nmass=piMass; pmass=prMass; mass=l0Mass; ps=0.101; break;
375 pmass=piMass; nmass=prMass; mass=l0Mass; ps=0.101; break;
379 AliError("invalide PDG code ! Assuming K0s...");
384 Double_t pxn=fNmom[0], pyn=fNmom[1], pzn=fNmom[2];
385 Double_t pxp=fPmom[0], pyp=fPmom[1], pzp=fPmom[2];
387 Double_t en=TMath::Sqrt(nmass*nmass + pxn*pxn + pyn*pyn + pzn*pzn);
388 Double_t ep=TMath::Sqrt(pmass*pmass + pxp*pxp + pyp*pyp + pzp*pzp);
389 Double_t pxl=pxn+pxp, pyl=pyn+pyp, pzl=pzn+pzp;
390 Double_t pl=TMath::Sqrt(pxl*pxl + pyl*pyl + pzl*pzl);
392 fEffMass=TMath::Sqrt((en+ep)*(en+ep)-pl*pl);
394 Double_t beta=pl/(en+ep);
395 Double_t pln=(pxn*pxl + pyn*pyl + pzn*pzl)/pl;
396 Double_t plp=(pxp*pxl + pyp*pyl + pzp*pzl)/pl;
398 Double_t pt2=pxp*pxp + pyp*pyp + pzp*pzp - plp*plp;
400 Double_t a=(plp-pln)/(plp+pln);
401 a -= (pmass*pmass-nmass*nmass)/(mass*mass);
402 a = 0.25*beta*beta*mass*mass*a*a + pt2;
408 void AliESDv0::GetPxPyPz(Double_t &px, Double_t &py, Double_t &pz) const {
409 //--------------------------------------------------------------------
410 // This function returns V0's momentum (global)
411 //--------------------------------------------------------------------
412 px=fNmom[0]+fPmom[0];
413 py=fNmom[1]+fPmom[1];
414 pz=fNmom[2]+fPmom[2];
417 void AliESDv0::GetXYZ(Double_t &x, Double_t &y, Double_t &z) const {
418 //--------------------------------------------------------------------
419 // This function returns V0's position (global)
420 //--------------------------------------------------------------------
426 Float_t AliESDv0::GetD(Double_t x0, Double_t y0) const {
427 //--------------------------------------------------------------------
428 // This function returns V0's impact parameter calculated in 2D in XY plane
429 //--------------------------------------------------------------------
430 Double_t x=fPos[0],y=fPos[1];
431 Double_t px=fNmom[0]+fPmom[0];
432 Double_t py=fNmom[1]+fPmom[1];
434 Double_t dz=(x0-x)*py - (y0-y)*px;
435 Double_t d=TMath::Sqrt(dz*dz/(px*px+py*py));
439 Float_t AliESDv0::GetD(Double_t x0, Double_t y0, Double_t z0) const {
440 //--------------------------------------------------------------------
441 // This function returns V0's impact parameter calculated in 3D
442 //--------------------------------------------------------------------
443 Double_t x=fPos[0],y=fPos[1],z=fPos[2];
444 Double_t px=fNmom[0]+fPmom[0];
445 Double_t py=fNmom[1]+fPmom[1];
446 Double_t pz=fNmom[2]+fPmom[2];
448 Double_t dx=(y0-y)*pz - (z0-z)*py;
449 Double_t dy=(x0-x)*pz - (z0-z)*px;
450 Double_t dz=(x0-x)*py - (y0-y)*px;
451 Double_t d=TMath::Sqrt((dx*dx+dy*dy+dz*dz)/(px*px+py*py+pz*pz));
455 Float_t AliESDv0::GetV0CosineOfPointingAngle(Double_t refPointX, Double_t refPointY, Double_t refPointZ) const {
456 // calculates the pointing angle of the V0 wrt a reference point
458 Double_t momV0[3]; //momentum of the V0
459 GetPxPyPz(momV0[0],momV0[1],momV0[2]);
461 Double_t deltaPos[3]; //vector between the reference point and the V0 vertex
462 deltaPos[0] = fPos[0] - refPointX;
463 deltaPos[1] = fPos[1] - refPointY;
464 deltaPos[2] = fPos[2] - refPointZ;
466 Double_t momV02 = momV0[0]*momV0[0] + momV0[1]*momV0[1] + momV0[2]*momV0[2];
467 Double_t deltaPos2 = deltaPos[0]*deltaPos[0] + deltaPos[1]*deltaPos[1] + deltaPos[2]*deltaPos[2];
469 Double_t cosinePointingAngle = (deltaPos[0]*momV0[0] +
470 deltaPos[1]*momV0[1] +
471 deltaPos[2]*momV0[2] ) /
472 TMath::Sqrt(momV02 * deltaPos2);
474 return cosinePointingAngle;
478 // **** The following functions need to be revised
480 void AliESDv0::GetPosCov(Double_t cov[6]) const {
482 for (Int_t i=0; i<6; ++i) cov[i] = fPosCov[i];
486 Double_t AliESDv0::GetSigmaY(){
488 // return sigmay in y at vertex position using covariance matrix
490 const Double_t * cp = fParamP.GetCovariance();
491 const Double_t * cm = fParamN.GetCovariance();
492 Double_t sigmay = cp[0]+cm[0]+ cp[5]*(fParamP.GetX()-fRr)*(fParamP.GetX()-fRr)+ cm[5]*(fParamN.GetX()-fRr)*(fParamN.GetX()-fRr);
493 return (sigmay>0) ? TMath::Sqrt(sigmay):100;
496 Double_t AliESDv0::GetSigmaZ(){
498 // return sigmay in y at vertex position using covariance matrix
500 const Double_t * cp = fParamP.GetCovariance();
501 const Double_t * cm = fParamN.GetCovariance();
502 Double_t sigmaz = cp[2]+cm[2]+ cp[9]*(fParamP.GetX()-fRr)*(fParamP.GetX()-fRr)+ cm[9]*(fParamN.GetX()-fRr)*(fParamN.GetX()-fRr);
503 return (sigmaz>0) ? TMath::Sqrt(sigmaz):100;
506 Double_t AliESDv0::GetSigmaD0(){
508 // Sigma parameterization using covariance matrix
510 // sigma of distance between two tracks in vertex position
511 // sigma of DCA is proportianal to sigmaD0
512 // factor 2 difference is explained by the fact that the DCA is calculated at the position
513 // where the tracks as closest together ( not exact position of the vertex)
515 const Double_t * cp = fParamP.GetCovariance();
516 const Double_t * cm = fParamN.GetCovariance();
517 Double_t sigmaD0 = cp[0]+cm[0]+cp[2]+cm[2]+fgkParams.fPSigmaOffsetD0*fgkParams.fPSigmaOffsetD0;
518 sigmaD0 += ((fParamP.GetX()-fRr)*(fParamP.GetX()-fRr))*(cp[5]+cp[9]);
519 sigmaD0 += ((fParamN.GetX()-fRr)*(fParamN.GetX()-fRr))*(cm[5]+cm[9]);
520 return (sigmaD0>0)? TMath::Sqrt(sigmaD0):100;
524 Double_t AliESDv0::GetSigmaAP0(){
526 //Sigma parameterization using covariance matrices
528 Double_t prec = TMath::Sqrt((fNmom[0]+fPmom[0])*(fNmom[0]+fPmom[0])
529 +(fNmom[1]+fPmom[1])*(fNmom[1]+fPmom[1])
530 +(fNmom[2]+fPmom[2])*(fNmom[2]+fPmom[2]));
531 Double_t normp = TMath::Sqrt(fPmom[0]*fPmom[0]+fPmom[1]*fPmom[1]+fPmom[2]*fPmom[2])/prec; // fraction of the momenta
532 Double_t normm = TMath::Sqrt(fNmom[0]*fNmom[0]+fNmom[1]*fNmom[1]+fNmom[2]*fNmom[2])/prec;
533 const Double_t * cp = fParamP.GetCovariance();
534 const Double_t * cm = fParamN.GetCovariance();
535 Double_t sigmaAP0 = fgkParams.fPSigmaOffsetAP0*fgkParams.fPSigmaOffsetAP0; // minimal part
536 sigmaAP0 += (cp[5]+cp[9])*(normp*normp)+(cm[5]+cm[9])*(normm*normm); // angular resolution part
537 Double_t sigmaAP1 = GetSigmaD0()/(TMath::Abs(fRr)+0.01); // vertex position part
538 sigmaAP0 += 0.5*sigmaAP1*sigmaAP1;
539 return (sigmaAP0>0)? TMath::Sqrt(sigmaAP0):100;
542 Double_t AliESDv0::GetEffectiveSigmaD0(){
544 // minimax - effective Sigma parameterization
545 // p12 effective curvature and v0 radius postion used as parameters
547 Double_t p12 = TMath::Sqrt(fParamP.GetParameter()[4]*fParamP.GetParameter()[4]+
548 fParamN.GetParameter()[4]*fParamN.GetParameter()[4]);
549 Double_t sigmaED0= TMath::Max(TMath::Sqrt(fRr)-fgkParams.fPSigmaRminDE,0.0)*fgkParams.fPSigmaCoefDE*p12*p12;
552 sigmaED0 = TMath::Sqrt(sigmaED0+fgkParams.fPSigmaOffsetDE*fgkParams.fPSigmaOffsetDE);
553 return (sigmaED0<fgkParams.fPSigmaMaxDE) ? sigmaED0: fgkParams.fPSigmaMaxDE;
557 Double_t AliESDv0::GetEffectiveSigmaAP0(){
559 // effective Sigma parameterization of point angle resolution
561 Double_t p12 = TMath::Sqrt(fParamP.GetParameter()[4]*fParamP.GetParameter()[4]+
562 fParamN.GetParameter()[4]*fParamN.GetParameter()[4]);
563 Double_t sigmaAPE= fgkParams.fPSigmaBase0APE;
564 sigmaAPE+= fgkParams.fPSigmaR0APE/(fgkParams.fPSigmaR1APE+fRr);
565 sigmaAPE*= (fgkParams.fPSigmaP0APE+fgkParams.fPSigmaP1APE*p12);
566 sigmaAPE = TMath::Min(sigmaAPE,fgkParams.fPSigmaMaxAPE);
571 Double_t AliESDv0::GetMinimaxSigmaAP0(){
573 // calculate mini-max effective sigma of point angle resolution
575 //compv0->fTree->SetAlias("SigmaAP2","max(min((SigmaAP0+SigmaAPE0)*0.5,1.5*SigmaAPE0),0.5*SigmaAPE0+0.003)");
576 Double_t effectiveSigma = GetEffectiveSigmaAP0();
577 Double_t sigmaMMAP = 0.5*(GetSigmaAP0()+effectiveSigma);
578 sigmaMMAP = TMath::Min(sigmaMMAP, fgkParams.fPMaxFractionAP0*effectiveSigma);
579 sigmaMMAP = TMath::Max(sigmaMMAP, fgkParams.fPMinFractionAP0*effectiveSigma+fgkParams.fPMinAP0);
582 Double_t AliESDv0::GetMinimaxSigmaD0(){
584 // calculate mini-max sigma of dca resolution
586 //compv0->fTree->SetAlias("SigmaD2","max(min((SigmaD0+SigmaDE0)*0.5,1.5*SigmaDE0),0.5*SigmaDE0)");
587 Double_t effectiveSigma = GetEffectiveSigmaD0();
588 Double_t sigmaMMD0 = 0.5*(GetSigmaD0()+effectiveSigma);
589 sigmaMMD0 = TMath::Min(sigmaMMD0, fgkParams.fPMaxFractionD0*effectiveSigma);
590 sigmaMMD0 = TMath::Max(sigmaMMD0, fgkParams.fPMinFractionD0*effectiveSigma+fgkParams.fPMinD0);
595 Double_t AliESDv0::GetLikelihoodAP(Int_t mode0, Int_t mode1){
597 // get likelihood for point angle
599 Double_t sigmaAP = 0.007; //default sigma
602 sigmaAP = GetSigmaAP0(); // mode 0 - covariance matrix estimates used
605 sigmaAP = GetEffectiveSigmaAP0(); // mode 1 - effective sigma used
608 sigmaAP = GetMinimaxSigmaAP0(); // mode 2 - minimax sigma
611 Double_t apNorm = TMath::Min(TMath::ACos(fPointAngle)/sigmaAP,50.);
612 //normalized point angle, restricted - because of overflow problems in Exp
613 Double_t likelihood = 0;
616 likelihood = TMath::Exp(-0.5*apNorm*apNorm);
620 likelihood = (TMath::Exp(-0.5*apNorm*apNorm)+0.5* TMath::Exp(-0.25*apNorm*apNorm))/1.5;
624 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;
631 Double_t AliESDv0::GetLikelihoodD(Int_t mode0, Int_t mode1){
633 // get likelihood for DCA
635 Double_t sigmaD = 0.03; //default sigma
638 sigmaD = GetSigmaD0(); // mode 0 - covariance matrix estimates used
641 sigmaD = GetEffectiveSigmaD0(); // mode 1 - effective sigma used
644 sigmaD = GetMinimaxSigmaD0(); // mode 2 - minimax sigma
648 //Bo: Double_t dNorm = TMath::Min(fDist2/sigmaD,50.);
649 Double_t dNorm = TMath::Min(fDcaV0Daughters/sigmaD,50.);//Bo:
650 //normalized point angle, restricted - because of overflow problems in Exp
651 Double_t likelihood = 0;
654 likelihood = TMath::Exp(-2.*dNorm);
658 likelihood = (TMath::Exp(-2.*dNorm)+0.5* TMath::Exp(-dNorm))/1.5;
662 likelihood = (TMath::Exp(-2.*dNorm)+0.5* TMath::Exp(-dNorm)+0.25*TMath::Exp(-0.5*dNorm))/1.75;
670 Double_t AliESDv0::GetLikelihoodC(Int_t mode0, Int_t /*mode1*/) const {
672 // get likelihood for Causality
673 // !!! Causality variables defined in AliITStrackerMI !!!
674 // when more information was available
676 Double_t likelihood = 0.5;
677 Double_t minCausal = TMath::Min(fCausality[0],fCausality[1]);
678 Double_t maxCausal = TMath::Max(fCausality[0],fCausality[1]);
679 // minCausal = TMath::Max(minCausal,0.5*maxCausal);
680 //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");
685 likelihood = TMath::Power((1.05-2*(0.8-TMath::Exp(-maxCausal))),4.);
688 likelihood = TMath::Power(1.05-(2*(0.8-TMath::Exp(-maxCausal))+(2*(0.8-TMath::Exp(-minCausal))))*0.5,4.);
695 void AliESDv0::SetCausality(Float_t pb0, Float_t pb1, Float_t pa0, Float_t pa1)
700 fCausality[0] = pb0; // probability - track 0 exist before vertex
701 fCausality[1] = pb1; // probability - track 1 exist before vertex
702 fCausality[2] = pa0; // probability - track 0 exist close after vertex
703 fCausality[3] = pa1; // probability - track 1 exist close after vertex
705 void AliESDv0::SetClusters(const Int_t *clp, const Int_t *clm)
708 // Set its clusters indexes
710 for (Int_t i=0;i<6;i++) fClusters[0][i] = clp[i];
711 for (Int_t i=0;i<6;i++) fClusters[1][i] = clm[i];
714 Double_t AliESDv0::GetEffMass(UInt_t p1, UInt_t p2) const{
716 // calculate effective mass
718 const Double_t kpmass[5] = {TDatabasePDG::Instance()->GetParticle(kElectron)->Mass(),
719 TDatabasePDG::Instance()->GetParticle(kMuonMinus)->Mass(),
720 TDatabasePDG::Instance()->GetParticle(kPiPlus)->Mass(),
721 TDatabasePDG::Instance()->GetParticle(kKPlus)->Mass(),
722 TDatabasePDG::Instance()->GetParticle(kProton)->Mass()};
726 Double_t mass1 = kpmass[p1];
727 Double_t mass2 = kpmass[p2];
728 const Double_t *m1 = fPmom;
729 const Double_t *m2 = fNmom;
731 //if (fRP[p1]+fRM[p2]<fRP[p2]+fRM[p1]){
736 Double_t e1 = TMath::Sqrt(mass1*mass1+
740 Double_t e2 = TMath::Sqrt(mass2*mass2+
745 (m2[0]+m1[0])*(m2[0]+m1[0])+
746 (m2[1]+m1[1])*(m2[1]+m1[1])+
747 (m2[2]+m1[2])*(m2[2]+m1[2]);
749 mass = (e1+e2)*(e1+e2)-mass;
750 if (mass < 0.) mass = 0.;
751 return (TMath::Sqrt(mass));
753 if(p1>4 || p2>4) return -1;
754 const AliExternalTrackParam *paramP = GetParamP();
755 const AliExternalTrackParam *paramN = GetParamN();
756 if (paramP->GetParameter()[4]<0){
760 Double_t pmom[3]={0}, nmom[3]={0};
761 paramP->GetPxPyPz(pmom);
762 paramN->GetPxPyPz(nmom);
763 Double_t e12 = kpmass[p1]*kpmass[p1]+pmom[0]*pmom[0]+pmom[1]*pmom[1]+pmom[2]*pmom[2];
764 Double_t e22 = kpmass[p2]*kpmass[p2]+nmom[0]*nmom[0]+nmom[1]*nmom[1]+nmom[2]*nmom[2];
765 Double_t cmass = TMath::Sqrt(TMath::Max(kpmass[p1]*kpmass[p1]+kpmass[p2]*kpmass[p2]
766 +2.*(TMath::Sqrt(e12*e22)-pmom[0]*nmom[0]-pmom[1]*nmom[1]-pmom[2]*nmom[2]),0.));
773 Double_t AliESDv0::GetKFInfo(UInt_t p1, UInt_t p2, Int_t type) const{
777 // 1 - return err mass
780 const Int_t spdg[5]={kPositron,kMuonPlus,kPiPlus, kKPlus, kProton};
781 const AliExternalTrackParam *paramP = GetParamP();
782 const AliExternalTrackParam *paramN = GetParamN();
783 if (paramP->GetSign()<0){
787 AliKFParticle kfp1( *(paramP), spdg[p1] *TMath::Sign(1,p1) );
788 AliKFParticle kfp2( *(paramN), spdg[p2] *TMath::Sign(1,p2) );
789 AliKFParticle *v0KF = new AliKFParticle;
792 if (type==0) return v0KF->GetMass();
793 if (type==1) return v0KF->GetErrMass();
794 if (type==2) return v0KF->GetChi2();
799 Double_t AliESDv0::GetKFInfoScale(UInt_t p1, UInt_t p2, Int_t type, Double_t d1pt, Double_t s1pt) const{
803 // 1 - return err mass
806 // s1pt - scaling of 1/pt
807 // Important function to benchmark the pt resolution, and to find out systematic distortion
809 const Int_t spdg[5]={kPositron,kMuonPlus,kPiPlus, kKPlus, kProton};
810 const AliExternalTrackParam *paramP = GetParamP();
811 const AliExternalTrackParam *paramN = GetParamN();
812 if (paramP->GetSign()<0){
816 Double_t *pparam1 = (Double_t*)paramP->GetParameter();
817 Double_t *pparam2 = (Double_t*)paramN->GetParameter();
820 pparam1[4]*=(1+s1pt);
821 pparam2[4]*=(1+s1pt);
823 AliKFParticle kfp1( *paramP, spdg[p1] *TMath::Sign(1,p1) );
824 AliKFParticle kfp2( *paramN, spdg[p2] *TMath::Sign(1,p2) );
825 AliKFParticle *v0KF = new AliKFParticle;
828 if (type==0) return v0KF->GetMass();
829 if (type==1) return v0KF->GetErrMass();
830 if (type==2) return v0KF->GetChi2();