]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDcascade.cxx
Check the status of the command pipe, it might be 0x0 if the process cannot allocate...
[u/mrichter/AliRoot.git] / STEER / AliESDcascade.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 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 //     Modified by: Antonin Maire,IPHC, Antonin.Maire@ires.in2p3.fr
25 //            and  Boris Hippolyte,IPHC, hippolyt@in2p3.fr 
26 //-------------------------------------------------------------------------
27
28 #include <TDatabasePDG.h>
29 #include <TMath.h>
30 #include <TVector3.h>
31
32 #include "AliLog.h"
33 #include "AliExternalTrackParam.h"
34 #include "AliESDv0.h"
35 #include "AliESDcascade.h"
36
37 ClassImp(AliESDcascade)
38
39 AliESDcascade::AliESDcascade() : 
40   AliESDv0(),
41   fEffMassXi(TDatabasePDG::Instance()->GetParticle(kXiMinus)->Mass()),
42   fChi2Xi(1024),
43   fDcaXiDaughters(1024),
44   fPdgCodeXi(kXiMinus),
45   fBachIdx(-1)
46 {
47   //--------------------------------------------------------------------
48   // Default constructor  (Xi-)
49   //--------------------------------------------------------------------
50   for (Int_t j=0; j<3; j++) {
51     fPosXi[j]=0.;
52     fBachMom[j]=0.;
53   }
54
55   fPosCovXi[0]=1024;
56   fPosCovXi[1]=fPosCovXi[2]=0.;
57   fPosCovXi[3]=1024;
58   fPosCovXi[4]=0.;
59   fPosCovXi[5]=1024;
60
61   fBachMomCov[0]=1024;
62   fBachMomCov[1]=fBachMomCov[2]=0.;
63   fBachMomCov[3]=1024;
64   fBachMomCov[4]=0.;
65   fBachMomCov[5]=1024;
66 }
67
68 AliESDcascade::AliESDcascade(const AliESDcascade& cas) :
69   AliESDv0(cas),
70   fEffMassXi(cas.fEffMassXi),
71   fChi2Xi(cas.fChi2Xi),
72   fDcaXiDaughters(cas.fDcaXiDaughters),
73   fPdgCodeXi(cas.fPdgCodeXi),
74   fBachIdx(cas.fBachIdx)
75 {
76   //--------------------------------------------------------------------
77   // The copy constructor
78   //--------------------------------------------------------------------
79   for (int i=0; i<3; i++) {
80     fPosXi[i]     = cas.fPosXi[i];
81     fBachMom[i] = cas.fBachMom[i];
82   }
83   for (int i=0; i<6; i++) {
84     fPosCovXi[i]   = cas.fPosCovXi[i];
85     fBachMomCov[i] = cas.fBachMomCov[i];
86   }
87 }
88
89 AliESDcascade::AliESDcascade(const AliESDv0 &v,
90                              const AliExternalTrackParam &t, Int_t i) : 
91   AliESDv0(v),
92   fEffMassXi(TDatabasePDG::Instance()->GetParticle(kXiMinus)->Mass()),
93   fChi2Xi(1024),
94   fDcaXiDaughters(1024),
95   fPdgCodeXi(kXiMinus),
96   fBachIdx(i)
97 {
98   //--------------------------------------------------------------------
99   // Main constructor  (Xi-)
100   //--------------------------------------------------------------------
101
102   Double_t r[3]; t.GetXYZ(r);
103   Double_t x1=r[0], y1=r[1], z1=r[2]; // position of the bachelor
104   Double_t p[3]; t.GetPxPyPz(p);
105   Double_t px1=p[0], py1=p[1], pz1=p[2];// momentum of the bachelor track
106
107   Double_t x2,y2,z2;          // position of the V0 
108   v.GetXYZ(x2,y2,z2);    
109   Double_t px2,py2,pz2;       // momentum of V0
110   v.GetPxPyPz(px2,py2,pz2);
111
112   Double_t a2=((x1-x2)*px2+(y1-y2)*py2+(z1-z2)*pz2)/(px2*px2+py2*py2+pz2*pz2);
113
114   Double_t xm=x2+a2*px2;
115   Double_t ym=y2+a2*py2;
116   Double_t zm=z2+a2*pz2;
117
118   // position of the cascade decay
119   
120   fPosXi[0]=0.5*(x1+xm); fPosXi[1]=0.5*(y1+ym); fPosXi[2]=0.5*(z1+zm);
121     
122
123   // invariant mass of the cascade (default is Ximinus)
124   
125   Double_t e1=TMath::Sqrt(0.13957*0.13957 + px1*px1 + py1*py1 + pz1*pz1);
126   Double_t e2=TMath::Sqrt(1.11568*1.11568 + px2*px2 + py2*py2 + pz2*pz2);
127   
128   fEffMassXi=TMath::Sqrt((e1+e2)*(e1+e2)-
129     (px1+px2)*(px1+px2)-(py1+py2)*(py1+py2)-(pz1+pz2)*(pz1+pz2));
130
131
132   // momenta of the bachelor and the V0
133   
134   fBachMom[0]=px1; fBachMom[1]=py1; fBachMom[2]=pz1; 
135
136   // Setting pdg code and fixing charge
137   if (t.Charge()<0)
138     fPdgCodeXi = kXiMinus;
139   else
140     fPdgCodeXi = kXiPlusBar;
141
142   //PH Covariance matrices: to be calculated correctly in the future
143   fPosCovXi[0]=1024;
144   fPosCovXi[1]=fPosCovXi[2]=0.;
145   fPosCovXi[3]=1024;
146   fPosCovXi[4]=0.;
147   fPosCovXi[5]=1024;
148
149   fBachMomCov[0]=1024;
150   fBachMomCov[1]=fBachMomCov[2]=0.;
151   fBachMomCov[3]=1024;
152   fBachMomCov[4]=0.;
153   fBachMomCov[5]=1024;
154
155   fChi2Xi=1024.; 
156
157 }
158
159 AliESDcascade& AliESDcascade::operator=(const AliESDcascade& cas)
160 {
161   //--------------------------------------------------------------------
162   // The assignment operator
163   //--------------------------------------------------------------------
164
165   if(this==&cas) return *this;
166   AliESDv0::operator=(cas);
167
168   fEffMassXi = cas.fEffMassXi;
169   fChi2Xi    = cas.fChi2Xi;
170   fDcaXiDaughters = cas.fDcaXiDaughters;
171   fPdgCodeXi      = cas.fPdgCodeXi;
172   fBachIdx        = cas.fBachIdx;
173   for (int i=0; i<3; i++) {
174     fPosXi[i]     = cas.fPosXi[i];
175     fBachMom[i]   = cas.fBachMom[i];
176   }
177   for (int i=0; i<6; i++) {
178     fPosCovXi[i]   = cas.fPosCovXi[i];
179     fBachMomCov[i] = cas.fBachMomCov[i];
180   }
181   return *this;
182 }
183
184 void AliESDcascade::Copy(TObject &obj) const {
185   
186   // this overwrites the virtual TOBject::Copy()
187   // to allow run time copying without casting
188   // in AliESDEvent
189
190   if(this==&obj)return;
191   AliESDcascade *robj = dynamic_cast<AliESDcascade*>(&obj);
192   if(!robj)return; // not an AliESDcascade
193   *robj = *this;
194 }
195
196 AliESDcascade::~AliESDcascade() {
197   //--------------------------------------------------------------------
198   // Empty destructor
199   //--------------------------------------------------------------------
200 }
201
202 // Start with AliVParticle functions
203 Double_t AliESDcascade::E() const {
204   //--------------------------------------------------------------------
205   // This gives the energy assuming the ChangeMassHypothesis was called
206   //--------------------------------------------------------------------
207   return E(fPdgCodeXi);
208 }
209
210 Double_t AliESDcascade::Y() const {
211   //--------------------------------------------------------------------
212   // This gives the energy assuming the ChangeMassHypothesis was called
213   //--------------------------------------------------------------------
214   return Y(fPdgCodeXi);
215 }
216
217 // Then extend AliVParticle functions
218 Double_t AliESDcascade::E(Int_t pdg) const {
219   //--------------------------------------------------------------------
220   // This gives the energy with the particle hypothesis as argument 
221   //--------------------------------------------------------------------
222   Double_t mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
223   return TMath::Sqrt(mass*mass+P()*P());
224 }
225
226 Double_t AliESDcascade::Y(Int_t pdg) const {
227   //--------------------------------------------------------------------
228   // This gives the rapidity with the particle hypothesis as argument 
229   //--------------------------------------------------------------------
230   return 0.5*TMath::Log((E(pdg)+Pz())/(E(pdg)-Pz()+1.e-13));
231 }
232
233 // Now the functions for analysis consistency
234 Double_t AliESDcascade::RapXi() const {
235   //--------------------------------------------------------------------
236   // This gives the pseudorapidity assuming a (Anti) Xi particle
237   //--------------------------------------------------------------------
238   return Y(kXiMinus);
239 }
240
241 Double_t AliESDcascade::RapOmega() const {
242   //--------------------------------------------------------------------
243   // This gives the pseudorapidity assuming a (Anti) Omega particle
244   //--------------------------------------------------------------------
245   return Y(kOmegaMinus);
246 }
247
248 Double_t AliESDcascade::AlphaXi() const {
249   //--------------------------------------------------------------------
250   // This gives the Armenteros-Podolanski alpha
251   //--------------------------------------------------------------------
252   TVector3 momBach(fBachMom[0],fBachMom[1],fBachMom[2]);
253   TVector3 momV0(fNmom[0]+fPmom[0],fNmom[1]+fPmom[1],fNmom[2]+fPmom[2]);
254   TVector3 momTot(Px(),Py(),Pz());
255
256   Double_t QlBach = momBach.Dot(momTot)/momTot.Mag();
257   Double_t QlV0 = momV0.Dot(momTot)/momTot.Mag();
258
259   return 1.-2./(1.+QlBach/QlV0);
260 }
261
262 Double_t AliESDcascade::PtArmXi() const {
263   //--------------------------------------------------------------------
264   // This gives the Armenteros-Podolanski ptarm
265   //--------------------------------------------------------------------
266   TVector3 momBach(fBachMom[0],fBachMom[1],fBachMom[2]);
267   TVector3 momTot(Px(),Py(),Pz());
268
269   return momBach.Perp(momTot);
270 }
271
272 // Then the older functions
273 Double_t AliESDcascade::ChangeMassHypothesis(Double_t &v0q, Int_t code) {
274   //--------------------------------------------------------------------
275   // This function changes the mass hypothesis for this cascade
276   // and returns the "kinematical quality" of this hypothesis
277   // together with the "quality" of associated V0 (argument v0q) 
278   //--------------------------------------------------------------------
279   Double_t nmass=0.13957, pmass=0.93827, ps0=0.101; 
280   Double_t bmass=0.13957, mass =1.3213,  ps =0.139;
281
282   if (Charge()*code<0)
283     fPdgCodeXi = code;
284   else {
285     AliWarning("Chosen PDG code does not match the sign of the bachelor... Corrected !!");
286     fPdgCodeXi = -code;
287   }
288
289   switch (fPdgCodeXi) {
290   case kXiMinus:
291        break;
292   case kXiPlusBar:
293        nmass=0.93827; pmass=0.13957; 
294        break;
295   case kOmegaMinus: 
296        bmass=0.49368; mass=1.67245; ps=0.211;
297        break;
298   case kOmegaPlusBar: 
299        nmass=0.93827; pmass=0.13957; 
300        bmass=0.49368; mass=1.67245; ps=0.211;
301        break;
302   default:
303        AliError("Invalide PDG code !  Assuming a Xi particle...");
304        if (Charge()<0) {
305          fPdgCodeXi=kXiMinus;
306        }
307        else {
308          fPdgCodeXi=kXiPlusBar;
309          nmass=0.93827; pmass=0.13957; 
310        }
311     break;
312   }
313
314   Double_t pxn=fNmom[0], pyn=fNmom[1], pzn=fNmom[2];
315   Double_t pxp=fPmom[0], pyp=fPmom[1], pzp=fPmom[2];
316
317   Double_t px0=pxn+pxp, py0=pyn+pyp, pz0=pzn+pzp;
318   Double_t p0=TMath::Sqrt(px0*px0 + py0*py0 + pz0*pz0);
319
320   Double_t e0=TMath::Sqrt(1.11568*1.11568 + p0*p0);
321   Double_t beta0=p0/e0;
322   Double_t pln=(pxn*px0 + pyn*py0 + pzn*pz0)/p0;
323   Double_t plp=(pxp*px0 + pyp*py0 + pzp*pz0)/p0;
324   Double_t pt2=pxp*pxp + pyp*pyp + pzp*pzp - plp*plp;
325
326   Double_t a=(plp-pln)/(plp+pln);
327   a -= (pmass*pmass-nmass*nmass)/(1.11568*1.11568);
328   a = 0.25*beta0*beta0*1.11568*1.11568*a*a + pt2;
329
330
331   v0q=a - ps0*ps0;
332
333
334   Double_t pxb=fBachMom[0], pyb=fBachMom[1], pzb=fBachMom[2]; 
335
336   Double_t eb=TMath::Sqrt(bmass*bmass + pxb*pxb + pyb*pyb + pzb*pzb);
337   Double_t pxl=px0+pxb, pyl=py0+pyb, pzl=pz0+pzb;
338   Double_t pl=TMath::Sqrt(pxl*pxl + pyl*pyl + pzl*pzl);
339   
340   fEffMassXi=TMath::Sqrt((e0+eb)*(e0+eb) - pl*pl);
341
342   Double_t beta=pl/(e0+eb);
343   Double_t pl0=(px0*pxl + py0*pyl + pz0*pzl)/pl;
344   Double_t plb=(pxb*pxl + pyb*pyl + pzb*pzl)/pl;
345   pt2=p0*p0 - pl0*pl0;
346
347   a=(pl0-plb)/(pl0+plb);
348   a -= (1.11568*1.11568-bmass*bmass)/(mass*mass);
349   a = 0.25*beta*beta*mass*mass*a*a + pt2;
350
351   return (a - ps*ps);
352 }
353
354 void 
355 AliESDcascade::GetPxPyPz(Double_t &px, Double_t &py, Double_t &pz) const {
356   //--------------------------------------------------------------------
357   // This function returns the cascade momentum (global)
358   //--------------------------------------------------------------------
359   px=fNmom[0]+fPmom[0]+fBachMom[0];
360   py=fNmom[1]+fPmom[1]+fBachMom[1]; 
361   pz=fNmom[2]+fPmom[2]+fBachMom[2]; 
362 }
363
364 void AliESDcascade::GetXYZcascade(Double_t &x, Double_t &y, Double_t &z) const {
365   //--------------------------------------------------------------------
366   // This function returns cascade position (global)
367   //--------------------------------------------------------------------
368   x=fPosXi[0];
369   y=fPosXi[1];
370   z=fPosXi[2];
371 }
372
373 Double_t AliESDcascade::GetDcascade(Double_t x0, Double_t y0, Double_t z0) const {
374   //--------------------------------------------------------------------
375   // This function returns the cascade impact parameter
376   //--------------------------------------------------------------------
377
378   Double_t x=fPosXi[0],y=fPosXi[1],z=fPosXi[2];
379   Double_t px=fNmom[0]+fPmom[0]+fBachMom[0];
380   Double_t py=fNmom[1]+fPmom[1]+fBachMom[1];
381   Double_t pz=fNmom[2]+fPmom[2]+fBachMom[2];
382
383   Double_t dx=(y0-y)*pz - (z0-z)*py; 
384   Double_t dy=(x0-x)*pz - (z0-z)*px;
385   Double_t dz=(x0-x)*py - (y0-y)*px;
386   Double_t d=TMath::Sqrt((dx*dx+dy*dy+dz*dz)/(px*px+py*py+pz*pz));
387
388   return d;
389 }
390
391 Double_t AliESDcascade::GetCascadeCosineOfPointingAngle(Double_t& refPointX, Double_t& refPointY, Double_t& refPointZ) const {
392   // calculates the pointing angle of the cascade wrt a reference point
393
394   Double_t momCas[3]; //momentum of the cascade
395   GetPxPyPz(momCas[0],momCas[1],momCas[2]);
396
397   Double_t deltaPos[3]; //vector between the reference point and the cascade vertex
398   deltaPos[0] = fPosXi[0] - refPointX;
399   deltaPos[1] = fPosXi[1] - refPointY;
400   deltaPos[2] = fPosXi[2] - refPointZ;
401
402   Double_t momCas2    = momCas[0]*momCas[0] + momCas[1]*momCas[1] + momCas[2]*momCas[2];
403   Double_t deltaPos2 = deltaPos[0]*deltaPos[0] + deltaPos[1]*deltaPos[1] + deltaPos[2]*deltaPos[2];
404
405   Double_t cosinePointingAngle = (deltaPos[0]*momCas[0] +
406                                   deltaPos[1]*momCas[1] +
407                                   deltaPos[2]*momCas[2] ) /
408     TMath::Sqrt(momCas2 * deltaPos2);
409   
410   return cosinePointingAngle;
411 }
412
413 void AliESDcascade::GetPosCovXi(Double_t cov[6]) const {
414
415   for (Int_t i=0; i<6; ++i) cov[i] = fPosCovXi[i];
416 }