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