]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliCascadeVertex.cxx
Introducing Riostream.h
[u/mrichter/AliRoot.git] / ITS / AliCascadeVertex.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 //-------------------------------------------------------------------------
17 //               Implementation of the cascade vertex class
18 //
19 //    Origin: Christian Kuhn, IReS, Strasbourg, christian.kuhn@ires.in2p3.fr
20 //-------------------------------------------------------------------------
21 #include <Riostream.h>
22 #include <TMath.h>
23
24 #include "AliCascadeVertex.h"
25 #include "AliV0vertex.h"
26 #include "AliITStrackV2.h"
27
28 ClassImp(AliCascadeVertex)
29
30 AliCascadeVertex::AliCascadeVertex() : TObject() {
31   //--------------------------------------------------------------------
32   // Default constructor  (Xi-)
33   //--------------------------------------------------------------------
34   fPdgCode=kXiMinus;
35   fEffMass=1.32131;
36   fChi2=1.e+33;
37   fPos[0]=fPos[1]=fPos[2]=0.;
38   fPosCov[0]=fPosCov[1]=fPosCov[2]=fPosCov[3]=fPosCov[4]=fPosCov[5]=0.;
39 }
40
41
42
43 inline Double_t det(Double_t a00, Double_t a01, Double_t a10, Double_t a11){
44   // determinant 2x2
45   return a00*a11 - a01*a10;
46 }
47
48 inline Double_t det (Double_t a00,Double_t a01,Double_t a02,
49                      Double_t a10,Double_t a11,Double_t a12,
50                      Double_t a20,Double_t a21,Double_t a22) {
51   // determinant 3x3
52   return 
53   a00*det(a11,a12,a21,a22)-a01*det(a10,a12,a20,a22)+a02*det(a10,a11,a20,a21);
54 }
55
56
57
58 AliCascadeVertex::AliCascadeVertex(const AliV0vertex &v,const AliITStrackV2 &t) {
59   //--------------------------------------------------------------------
60   // Main constructor
61   //--------------------------------------------------------------------
62   fPdgCode=kXiMinus;
63
64   fV0lab[0]=v.GetNlabel(); fV0lab[1]=v.GetPlabel();
65   fBachLab=t.GetLabel(); 
66
67   //Trivial estimation of the vertex parameters
68   Double_t pt, phi, x, par[5];
69   Double_t alpha, cs, sn;
70
71   t.GetExternalParameters(x,par); alpha=t.GetAlpha();
72   pt=1./TMath::Abs(par[4]);
73   phi=TMath::ASin(par[2]) + alpha;  
74
75   // momentum of the bachelor track
76
77   Double_t px1=pt*TMath::Cos(phi), py1=pt*TMath::Sin(phi), pz1=pt*par[3];
78
79   cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
80
81   Double_t x1=x*cs - par[0]*sn; // position of the bachelor at dca (bachelor,V0)
82   Double_t y1=x*sn + par[0]*cs;
83   Double_t z1=par[1];
84
85   Double_t x2,y2,z2;          // position of the V0 
86   v.GetXYZ(x2,y2,z2);
87     
88   Double_t px2,py2,pz2;       // momentum of V0
89   v.GetPxPyPz(px2,py2,pz2);
90
91   Double_t a2=((x1-x2)*px2+(y1-y2)*py2+(z1-z2)*pz2)/(px2*px2+py2*py2+pz2*pz2);
92
93   Double_t xm=x2+a2*px2;
94   Double_t ym=y2+a2*py2;
95   Double_t zm=z2+a2*pz2;
96
97   // position of the cascade decay
98   
99   fPos[0]=0.5*(x1+xm); fPos[1]=0.5*(y1+ym); fPos[2]=0.5*(z1+zm);
100     
101
102   // invariant mass of the cascade (default is Ximinus)
103   
104   Double_t e1=TMath::Sqrt(0.13957*0.13957 + px1*px1 + py1*py1 + pz1*pz1);
105   Double_t e2=TMath::Sqrt(1.11568*1.11568 + px2*px2 + py2*py2 + pz2*pz2);
106   
107   fEffMass=TMath::Sqrt((e1+e2)*(e1+e2)-
108     (px1+px2)*(px1+px2)-(py1+py2)*(py1+py2)-(pz1+pz2)*(pz1+pz2));
109
110
111   // momenta of the bachelor and the V0
112   
113   fBachMom[0]=px1; fBachMom[1]=py1; fBachMom[2]=pz1; 
114   v.GetNPxPyPz(px2,py2,pz2);
115   fV0mom[0][0]=px2; fV0mom[0][1]=py2; fV0mom[0][2]=pz2;
116   v.GetPPxPyPz(px2,py2,pz2);
117   fV0mom[1][0]=px2; fV0mom[1][1]=py2; fV0mom[1][2]=pz2;
118
119
120   fChi2=7.;   
121
122 }
123
124 /*
125 Double_t AliCascadeVertex::ChangeMassHypothesis(Double_t &v0q, Int_t code) {
126   //--------------------------------------------------------------------
127   // This function changes the mass hypothesis for this cascade
128   // and returns the "kinematical quality" of this hypothesis
129   // together with the "quality" of associated V0 (argument v0q) 
130   //--------------------------------------------------------------------
131   Double_t nmass=0.13957, pmass=0.93827, des0=0.9437-0.1723; 
132   Double_t bmass=0.13957, mass =1.3213,  des =1.1243-0.1970;
133
134   fPdgCode=code;
135
136   switch (code) {
137   case 213: 
138        bmass=0.93827; 
139        break;
140   case kXiMinus:
141        break;
142   case kXiPlusBar:
143        nmass=0.93827; pmass=0.13957; des0=-des0; 
144        des=-des;
145        break;
146   case kOmegaMinus: 
147        bmass=0.49368; mass=1.67245; des=1.1355-0.5369;
148        break;
149   case kOmegaPlusBar: 
150        nmass=0.93827; pmass=0.13957; des0=-des0; 
151        bmass=0.49368; mass=1.67245; des=0.5369-1.1355;
152        break;
153   default:
154        cerr<<"AliCascadeVertex::ChangeMassHypothesis: ";
155        cerr<<"Invalide PDG code !  Assuming XiMinus's...\n";
156        fPdgCode=kXiMinus;
157     break;
158   }
159
160   Double_t pxn=fV0mom[0][0], pyn=fV0mom[0][1], pzn=fV0mom[0][2];
161   Double_t pxp=fV0mom[1][0], pyp=fV0mom[1][1], pzp=fV0mom[1][2];
162   Double_t en=TMath::Sqrt(nmass*nmass + pxn*pxn + pyn*pyn + pzn*pzn);
163   Double_t ep=TMath::Sqrt(pmass*pmass + pxp*pxp + pyp*pyp + pzp*pzp);
164   Double_t px0=pxn+pxp, py0=pyn+pyp, pz0=pzn+pzp;
165   Double_t p0=TMath::Sqrt(px0*px0 + py0*py0 + pz0*pz0);
166
167   Double_t gamma0=(en+ep)/1.11568, betagamma0=p0/1.11568;
168   Double_t pln=(pxn*px0 + pyn*py0 + pzn*pz0)/p0;
169   Double_t plp=(pxp*px0 + pyp*py0 + pzp*pz0)/p0;
170   Double_t plps=gamma0*plp - betagamma0*ep;
171
172   Double_t diff0=2*gamma0*plps + betagamma0*des0;
173
174
175   v0q=plp-pln-diff0;
176
177
178   Double_t pxb=fBachMom[0], pyb=fBachMom[1], pzb=fBachMom[2]; 
179
180   Double_t e0=TMath::Sqrt(1.11568*1.11568 + p0*p0);
181   Double_t eb=TMath::Sqrt(bmass*bmass + pxb*pxb + pyb*pyb + pzb*pzb);
182   Double_t pxl=px0+pxb, pyl=py0+pyb, pzl=pz0+pzb;
183   Double_t pl=TMath::Sqrt(pxl*pxl + pyl*pyl + pzl*pzl);
184   
185   fEffMass=TMath::Sqrt((e0+eb)*(e0+eb) - pl*pl);
186
187   Double_t gamma=(e0+eb)/mass, betagamma=pl/mass;
188   Double_t pl0=(px0*pxl + py0*pyl + pz0*pzl)/pl;
189   Double_t plb=(pxb*pxl + pyb*pyl + pzb*pzl)/pl;
190   Double_t pl0s=gamma*pl0 - betagamma*e0;
191
192   Double_t diff=2*gamma*pl0s + betagamma*des;
193
194   return (pl0-plb-diff);
195 }
196 */
197
198 Double_t AliCascadeVertex::ChangeMassHypothesis(Double_t &v0q, Int_t code) {
199   //--------------------------------------------------------------------
200   // This function changes the mass hypothesis for this cascade
201   // and returns the "kinematical quality" of this hypothesis
202   // together with the "quality" of associated V0 (argument v0q) 
203   //--------------------------------------------------------------------
204   Double_t nmass=0.13957, pmass=0.93827, ps0=0.101; 
205   Double_t bmass=0.13957, mass =1.3213,  ps =0.139;
206
207   fPdgCode=code;
208
209   switch (code) {
210   case 213: 
211        bmass=0.93827; 
212        break;
213   case kXiMinus:
214        break;
215   case kXiPlusBar:
216        nmass=0.93827; pmass=0.13957; 
217        break;
218   case kOmegaMinus: 
219        bmass=0.49368; mass=1.67245; ps=0.211;
220        break;
221   case kOmegaPlusBar: 
222        nmass=0.93827; pmass=0.13957; 
223        bmass=0.49368; mass=1.67245; ps=0.211;
224        break;
225   default:
226        cerr<<"AliCascadeVertex::ChangeMassHypothesis: ";
227        cerr<<"Invalide PDG code !  Assuming XiMinus's...\n";
228        fPdgCode=kXiMinus;
229     break;
230   }
231
232   Double_t pxn=fV0mom[0][0], pyn=fV0mom[0][1], pzn=fV0mom[0][2];
233   Double_t pxp=fV0mom[1][0], pyp=fV0mom[1][1], pzp=fV0mom[1][2];
234   Double_t px0=pxn+pxp, py0=pyn+pyp, pz0=pzn+pzp;
235   Double_t p0=TMath::Sqrt(px0*px0 + py0*py0 + pz0*pz0);
236
237   Double_t e0=TMath::Sqrt(1.11568*1.11568 + p0*p0);
238   Double_t beta0=p0/e0;
239   Double_t pln=(pxn*px0 + pyn*py0 + pzn*pz0)/p0;
240   Double_t plp=(pxp*px0 + pyp*py0 + pzp*pz0)/p0;
241   Double_t pt2=pxp*pxp + pyp*pyp + pzp*pzp - plp*plp;
242
243   Double_t a=(plp-pln)/(plp+pln);
244   a -= (pmass*pmass-nmass*nmass)/(1.11568*1.11568);
245   a = 0.25*beta0*beta0*1.11568*1.11568*a*a + pt2;
246
247
248   v0q=a - ps0*ps0;
249
250
251   Double_t pxb=fBachMom[0], pyb=fBachMom[1], pzb=fBachMom[2]; 
252
253   Double_t eb=TMath::Sqrt(bmass*bmass + pxb*pxb + pyb*pyb + pzb*pzb);
254   Double_t pxl=px0+pxb, pyl=py0+pyb, pzl=pz0+pzb;
255   Double_t pl=TMath::Sqrt(pxl*pxl + pyl*pyl + pzl*pzl);
256   
257   fEffMass=TMath::Sqrt((e0+eb)*(e0+eb) - pl*pl);
258
259   Double_t beta=pl/(e0+eb);
260   Double_t pl0=(px0*pxl + py0*pyl + pz0*pzl)/pl;
261   Double_t plb=(pxb*pxl + pyb*pyl + pzb*pzl)/pl;
262   pt2=p0*p0 - pl0*pl0;
263
264   a=(pl0-plb)/(pl0+plb);
265   a -= (1.11568*1.11568-bmass*bmass)/(mass*mass);
266   a = 0.25*beta*beta*mass*mass*a*a + pt2;
267
268   return (a - ps*ps);
269 }
270
271 void 
272 AliCascadeVertex::GetPxPyPz(Double_t &px, Double_t &py, Double_t &pz) const {
273   //--------------------------------------------------------------------
274   // This function returns the cascade momentum (global)
275   //--------------------------------------------------------------------
276   px=fV0mom[0][0]+fV0mom[1][0]+fBachMom[0]; 
277   py=fV0mom[0][1]+fV0mom[1][1]+fBachMom[1]; 
278   pz=fV0mom[0][2]+fV0mom[1][2]+fBachMom[2]; 
279 }
280
281 void AliCascadeVertex::GetXYZ(Double_t &x, Double_t &y, Double_t &z) const {
282   //--------------------------------------------------------------------
283   // This function returns cascade position (global)
284   //--------------------------------------------------------------------
285   x=fPos[0]; 
286   y=fPos[1]; 
287   z=fPos[2]; 
288 }
289
290 Double_t AliCascadeVertex::GetD(Double_t x0, Double_t y0, Double_t z0) const {
291   //--------------------------------------------------------------------
292   // This function returns the cascade impact parameter
293   //--------------------------------------------------------------------
294
295   Double_t x=fPos[0],y=fPos[1],z=fPos[2];
296   Double_t px=fV0mom[0][0]+fV0mom[1][0]+fBachMom[0];
297   Double_t py=fV0mom[0][1]+fV0mom[1][1]+fBachMom[1];
298   Double_t pz=fV0mom[0][2]+fV0mom[1][2]+fBachMom[2];
299
300   Double_t dx=(y0-y)*pz - (z0-z)*py; 
301   Double_t dy=(x0-x)*pz - (z0-z)*px;
302   Double_t dz=(x0-x)*py - (y0-y)*px;
303   Double_t d=TMath::Sqrt((dx*dx+dy*dy+dz*dz)/(px*px+py*py+pz*pz));
304
305   return d;
306 }
307