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