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 cascade vertex class
20 // This is part of the Event Summary Data
21 // which contains the result of the reconstruction
22 // and is the main set of classes for analaysis
23 // Origin: Christian Kuhn, IReS, Strasbourg, christian.kuhn@ires.in2p3.fr
24 //-------------------------------------------------------------------------
26 #include <TDatabasePDG.h>
30 #include "AliExternalTrackParam.h"
32 #include "AliESDcascade.h"
34 ClassImp(AliESDcascade)
36 AliESDcascade::AliESDcascade() :
38 fEffMassXi(TDatabasePDG::Instance()->GetParticle(kXiMinus)->Mass()),
40 fDcaXiDaughters(1024),
44 //--------------------------------------------------------------------
45 // Default constructor (Xi-)
46 //--------------------------------------------------------------------
47 for (Int_t j=0; j<3; j++) {
53 fPosCovXi[1]=fPosCovXi[2]=0.;
59 fBachMomCov[1]=fBachMomCov[2]=0.;
65 AliESDcascade::~AliESDcascade() {
68 AliESDcascade::AliESDcascade(const AliESDv0 &v,
69 const AliExternalTrackParam &t, Int_t i) :
71 fEffMassXi(TDatabasePDG::Instance()->GetParticle(kXiMinus)->Mass()),
73 fDcaXiDaughters(1024),
77 //---------------------------------------------------------------------------------------------
78 // Main constructor (Xi-)
79 //---------------------------------------------------------------------------------------------
81 Double_t r[3]; t.GetXYZ(r);
82 Double_t x1=r[0], y1=r[1], z1=r[2]; // position of the bachelor
83 Double_t p[3]; t.GetPxPyPz(p);
84 Double_t px1=p[0], py1=p[1], pz1=p[2];// momentum of the bachelor track
86 Double_t x2,y2,z2; // position of the V0
88 Double_t px2,py2,pz2; // momentum of V0
89 v.GetPxPyPz(px2,py2,pz2);
91 Double_t a2=((x1-x2)*px2+(y1-y2)*py2+(z1-z2)*pz2)/(px2*px2+py2*py2+pz2*pz2);
93 Double_t xm=x2+a2*px2;
94 Double_t ym=y2+a2*py2;
95 Double_t zm=z2+a2*pz2;
97 //dca between V0 and bachelor
99 fDcaXiDaughters = TMath::Sqrt((x1-xm)*(x1-xm) + (y1-ym)*(y1-ym) + (z1-zm)*(z1-zm));
101 // position of the cascade decay
103 fPosXi[0]=0.5*(x1+xm); fPosXi[1]=0.5*(y1+ym); fPosXi[2]=0.5*(z1+zm);
106 // invariant mass of the cascade (default is Ximinus)
108 Double_t e1=TMath::Sqrt(0.13957*0.13957 + px1*px1 + py1*py1 + pz1*pz1);
109 Double_t e2=TMath::Sqrt(1.11568*1.11568 + px2*px2 + py2*py2 + pz2*pz2);
111 fEffMassXi=TMath::Sqrt((e1+e2)*(e1+e2)-
112 (px1+px2)*(px1+px2)-(py1+py2)*(py1+py2)-(pz1+pz2)*(pz1+pz2));
115 // momenta of the bachelor and the V0
117 fBachMom[0]=px1; fBachMom[1]=py1; fBachMom[2]=pz1;
119 //PH Covariance matrices: to be calculated correctly in the future
121 fPosCovXi[1]=fPosCovXi[2]=0.;
127 fBachMomCov[1]=fBachMomCov[2]=0.;
136 AliESDcascade::AliESDcascade(const AliESDcascade& cas) :
138 fEffMassXi(cas.fEffMassXi),
139 fChi2Xi(cas.fChi2Xi),
140 fDcaXiDaughters(cas.fDcaXiDaughters),
141 fPdgCodeXi(cas.fPdgCodeXi),
142 fBachIdx(cas.fBachIdx)
145 for (int i=0; i<3; i++) {
146 fPosXi[i] = cas.fPosXi[i];
147 fBachMom[i] = cas.fBachMom[i];
149 for (int i=0; i<6; i++) {
150 fPosCovXi[i] = cas.fPosCovXi[i];
151 fBachMomCov[i] = cas.fBachMomCov[i];
155 AliESDcascade& AliESDcascade::operator=(const AliESDcascade& cas){
158 // Assigmnet oeprator
161 if(this==&cas) return *this;
162 AliESDv0::operator=(cas);
164 fEffMassXi = cas.fEffMassXi;
165 fChi2Xi = cas.fChi2Xi;
166 fDcaXiDaughters = cas.fDcaXiDaughters;
167 fPdgCodeXi = cas.fPdgCodeXi;
168 fBachIdx = cas.fBachIdx;
169 for (int i=0; i<3; i++) {
170 fPosXi[i] = cas.fPosXi[i];
171 fBachMom[i] = cas.fBachMom[i];
173 for (int i=0; i<6; i++) {
174 fPosCovXi[i] = cas.fPosCovXi[i];
175 fBachMomCov[i] = cas.fBachMomCov[i];
180 void AliESDcascade::Copy(TObject &obj) const {
182 // this overwrites the virtual TOBject::Copy()
183 // to allow run time copying without casting
186 if(this==&obj)return;
187 AliESDcascade *robj = dynamic_cast<AliESDcascade*>(&obj);
188 if(!robj)return; // not an AliESDcascade
193 Double_t AliESDcascade::ChangeMassHypothesis(Double_t &v0q, Int_t code) {
194 //--------------------------------------------------------------------
195 // This function changes the mass hypothesis for this cascade
196 // and returns the "kinematical quality" of this hypothesis
197 // together with the "quality" of associated V0 (argument v0q)
198 //--------------------------------------------------------------------
199 Double_t nmass=0.13957, pmass=0.93827, ps0=0.101;
200 Double_t bmass=0.13957, mass =1.3213, ps =0.139;
211 nmass=0.93827; pmass=0.13957;
214 bmass=0.49368; mass=1.67245; ps=0.211;
217 nmass=0.93827; pmass=0.13957;
218 bmass=0.49368; mass=1.67245; ps=0.211;
221 AliError("Invalide PDG code ! Assuming XiMinus's...");
226 Double_t pxn=fNmom[0], pyn=fNmom[1], pzn=fNmom[2];
227 Double_t pxp=fPmom[0], pyp=fPmom[1], pzp=fPmom[2];
229 Double_t px0=pxn+pxp, py0=pyn+pyp, pz0=pzn+pzp;
230 Double_t p0=TMath::Sqrt(px0*px0 + py0*py0 + pz0*pz0);
232 Double_t e0=TMath::Sqrt(1.11568*1.11568 + p0*p0);
233 Double_t beta0=p0/e0;
234 Double_t pln=(pxn*px0 + pyn*py0 + pzn*pz0)/p0;
235 Double_t plp=(pxp*px0 + pyp*py0 + pzp*pz0)/p0;
236 Double_t pt2=pxp*pxp + pyp*pyp + pzp*pzp - plp*plp;
238 Double_t a=(plp-pln)/(plp+pln);
239 a -= (pmass*pmass-nmass*nmass)/(1.11568*1.11568);
240 a = 0.25*beta0*beta0*1.11568*1.11568*a*a + pt2;
246 Double_t pxb=fBachMom[0], pyb=fBachMom[1], pzb=fBachMom[2];
248 Double_t eb=TMath::Sqrt(bmass*bmass + pxb*pxb + pyb*pyb + pzb*pzb);
249 Double_t pxl=px0+pxb, pyl=py0+pyb, pzl=pz0+pzb;
250 Double_t pl=TMath::Sqrt(pxl*pxl + pyl*pyl + pzl*pzl);
252 fEffMassXi=TMath::Sqrt((e0+eb)*(e0+eb) - pl*pl);
254 Double_t beta=pl/(e0+eb);
255 Double_t pl0=(px0*pxl + py0*pyl + pz0*pzl)/pl;
256 Double_t plb=(pxb*pxl + pyb*pyl + pzb*pzl)/pl;
259 a=(pl0-plb)/(pl0+plb);
260 a -= (1.11568*1.11568-bmass*bmass)/(mass*mass);
261 a = 0.25*beta*beta*mass*mass*a*a + pt2;
267 AliESDcascade::GetPxPyPz(Double_t &px, Double_t &py, Double_t &pz) const {
268 //--------------------------------------------------------------------
269 // This function returns the cascade momentum (global)
270 //--------------------------------------------------------------------
271 px=fNmom[0]+fPmom[0]+fBachMom[0];
272 py=fNmom[1]+fPmom[1]+fBachMom[1];
273 pz=fNmom[2]+fPmom[2]+fBachMom[2];
276 void AliESDcascade::GetXYZcascade(Double_t &x, Double_t &y, Double_t &z) const {
277 //--------------------------------------------------------------------
278 // This function returns cascade position (global)
279 //--------------------------------------------------------------------
285 Double_t AliESDcascade::GetDcascade(Double_t x0, Double_t y0, Double_t z0) const {
286 //--------------------------------------------------------------------
287 // This function returns the cascade impact parameter
288 //--------------------------------------------------------------------
290 Double_t x=fPosXi[0],y=fPosXi[1],z=fPosXi[2];
291 Double_t px=fNmom[0]+fPmom[0]+fBachMom[0];
292 Double_t py=fNmom[1]+fPmom[1]+fBachMom[1];
293 Double_t pz=fNmom[2]+fPmom[2]+fBachMom[2];
295 Double_t dx=(y0-y)*pz - (z0-z)*py;
296 Double_t dy=(x0-x)*pz - (z0-z)*px;
297 Double_t dz=(x0-x)*py - (y0-y)*px;
298 Double_t d=TMath::Sqrt((dx*dx+dy*dy+dz*dz)/(px*px+py*py+pz*pz));
303 Double_t AliESDcascade::GetCascadeCosineOfPointingAngle(Double_t& refPointX, Double_t& refPointY, Double_t& refPointZ) const {
304 // calculates the pointing angle of the cascade wrt a reference point
306 Double_t momCas[3]; //momentum of the cascade
307 GetPxPyPz(momCas[0],momCas[1],momCas[2]);
309 Double_t deltaPos[3]; //vector between the reference point and the cascade vertex
310 deltaPos[0] = fPosXi[0] - refPointX;
311 deltaPos[1] = fPosXi[1] - refPointY;
312 deltaPos[2] = fPosXi[2] - refPointZ;
314 Double_t momCas2 = momCas[0]*momCas[0] + momCas[1]*momCas[1] + momCas[2]*momCas[2];
315 Double_t deltaPos2 = deltaPos[0]*deltaPos[0] + deltaPos[1]*deltaPos[1] + deltaPos[2]*deltaPos[2];
317 Double_t cosinePointingAngle = (deltaPos[0]*momCas[0] +
318 deltaPos[1]*momCas[1] +
319 deltaPos[2]*momCas[2] ) /
320 TMath::Sqrt(momCas2 * deltaPos2);
322 return cosinePointingAngle;
325 void AliESDcascade::GetPosCovXi(Double_t cov[6]) const {
327 for (Int_t i=0; i<6; ++i) cov[i] = fPosCovXi[i];