]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliCascadeVertexer.cxx
Deletion of PID root files included
[u/mrichter/AliRoot.git] / ITS / AliCascadeVertexer.cxx
CommitLineData
a9a2d814 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 vertexer class
18//
19// Origin: Christian Kuhn, IReS, Strasbourg, christian.kuhn@ires.in2p3.fr
20//-------------------------------------------------------------------------
21#include <TFile.h>
22#include <TTree.h>
23#include <TObjArray.h>
24#include <iostream.h>
25
26#include "AliCascadeVertex.h"
27#include "AliCascadeVertexer.h"
28#include "AliV0vertex.h"
29#include "AliITStrackV2.h"
30
31ClassImp(AliCascadeVertexer)
32
33Int_t
34AliCascadeVertexer::V0sTracks2CascadeVertices(const TFile *inp, TFile *out) {
35
36 //--------------------------------------------------------------------
37 //This function reconstructs cascade vertices
38 //--------------------------------------------------------------------
39 TFile *in=(TFile*)inp;
40 TDirectory *savedir=gDirectory;
41
42 // Tree for vertices(V0's)
43
44 TTree *vtxTree=(TTree*)gDirectory->Get("TreeV");
45 TBranch *branch=vtxTree->GetBranch("vertices");
46 Int_t nentrV0=(Int_t)vtxTree->GetEntries();
47
48 TObjArray vtxV0(nentrV0);
49
50 // fill TObjArray vtxV0 with vertices
51
52 Int_t i, nV0=0;
53 for (i=0; i<nentrV0; i++) {
54
55 AliV0vertex *ioVertex=new AliV0vertex;
56 branch->SetAddress(&ioVertex);
57 vtxTree->GetEvent(i);
58 nV0++;
59 vtxV0.AddLast(ioVertex);
60
61 }
62
63
64 in->cd();
65
66 // Tree for tracks
67
68 TTree *trkTree=(TTree*)in->Get("TreeT_ITS_0");
69 branch=trkTree->GetBranch("tracks");
70 Int_t nentr=(Int_t)trkTree->GetEntries();
71
72 TObjArray trks(nentr);
73
74 // fill TObjArray trks with tracks
75
76 Int_t ntrack=0;
77
78 for (i=0; i<nentr; i++) {
79
80 AliITStrackV2 *iotrack=new AliITStrackV2;
81 branch->SetAddress(&iotrack);
82 trkTree->GetEvent(i);
4ab260d6 83
84 iotrack->PropagateTo(3.,0.0023,65.19); iotrack->PropagateTo(2.5,0.,0.);
85
a9a2d814 86 ntrack++; trks.AddLast(iotrack);
87
88 }
89
90 // create Tree for cascades
91
92 out->cd();
93
94 TTree cascTree("TreeCasc","Tree with cascades");
95 AliCascadeVertex *ioCascade=0;
96 cascTree.Branch("cascades","AliCascadeVertex",&ioCascade,32000,0);
97
98 // loop on all vertices
99
100 Double_t massLambda=1.11568;
101 Int_t ncasc=0;
102
103 for (i=0; i<nV0; i++) {
104
105 AliV0vertex *V0ver=(AliV0vertex *)vtxV0.UncheckedAt(i);
106
107 V0ver->ChangeMassHypothesis(kLambda0); //I.B.
108
109 if (V0ver->GetEffMass()<massLambda-fMassWin || // condition of the V0 mass window (cut fMassWin)
110 V0ver->GetEffMass()>massLambda+fMassWin) continue;
111
112 if (V0ver->GetD(0,0,0)<fDV0min) continue; // condition of minimum impact parameter of the V0 (cut fDV0min)
113 // here why not cuting on pointing angle ???
114
115 // for each vertex in the good mass range, loop on all tracks (= bachelor candidates)
116
117 for (Int_t k=0; k<ntrack; k++) {
118
119 AliITStrackV2 *bachtrk=(AliITStrackV2 *)trks.UncheckedAt(k);
120
121 if (TMath::Abs(bachtrk->GetD())<fDBachMin) continue; // eliminate to small impact parameters
122
123 if (V0ver->GetPdgCode()==kLambda0 && bachtrk->Get1Pt()<0.) continue; // condition on V0 label
124 if (V0ver->GetPdgCode()==kLambda0Bar && bachtrk->Get1Pt()>0.) continue; // + good sign for bachelor
125
126 AliV0vertex V0(*V0ver), *pV0=&V0;
127 AliITStrackV2 bt(*bachtrk), *pbt=&bt;
128
129 // calculation of the distance of closest approach between the V0 and the bachelor
130
131 Double_t dca=PropagateToDCA(pV0,pbt);
132 if (dca > fDCAmax) continue; // cut on dca
133
134 // construction of a cascade object
135
136 AliCascadeVertex cascade(*pV0,*pbt);
137 if (cascade.GetChi2() > fChi2max) continue;
138
139 // get cascade decay position (V0, bachelor)
140
141 Double_t x,y,z; cascade.GetXYZ(x,y,z);
142 Double_t r2=x*x + y*y;
143 if (r2 > fRmax*fRmax) continue; // condition on fiducial zone
144 if (r2 < fRmin*fRmin) continue;
145
4ab260d6 146 {
147 //I.B.
148 Double_t x1,y1,z1; V0ver->GetXYZ(x1,y1,z1);
149 if (r2 > (x1*x1+y1*y1)) continue;
150 if (z*z > z1*z1) continue;
151 }
152
a9a2d814 153 // get cascade momentum
154
155 Double_t px,py,pz; cascade.GetPxPyPz(px,py,pz);
156 Double_t p2=px*px+py*py+pz*pz;
157 Double_t cost=(x*px+y*py+z*pz)/TMath::Sqrt(p2*(r2+z*z));
158
159 if (cost<fCPAmax) continue; // condition on the cascade pointing angle
160
161 //cascade.ChangeMassHypothesis(); //default is Xi
162
163
164 ioCascade=&cascade; cascTree.Fill();
165
166 ncasc++;
167
168 }
169 }
170
171 cerr<<"Number of reconstructed cascades: "<<ncasc<<endl;
172
173 cascTree.Write();
174
175 trks.Delete();
176 vtxV0.Delete();
177
178 savedir->cd();
179
180 return 0;
181}
182
183inline Double_t det(Double_t a00, Double_t a01, Double_t a10, Double_t a11){
184 // determinant 2x2
185 return a00*a11 - a01*a10;
186}
187
188inline Double_t det (Double_t a00,Double_t a01,Double_t a02,
189 Double_t a10,Double_t a11,Double_t a12,
190 Double_t a20,Double_t a21,Double_t a22) {
191 // determinant 3x3
192 return
193 a00*det(a11,a12,a21,a22)-a01*det(a10,a12,a20,a22)+a02*det(a10,a11,a20,a21);
194}
195
196
197
198
199Double_t
200AliCascadeVertexer::PropagateToDCA(AliV0vertex *v, AliITStrackV2 *t) {
201 //--------------------------------------------------------------------
202 // This function returns the DCA between the V0 and the track
203 //--------------------------------------------------------------------
204
205 Double_t phi, x, par[5];
206 Double_t alpha, cs1, sn1;
207
208 t->GetExternalParameters(x,par); alpha=t->GetAlpha();
209 phi=TMath::ASin(par[2]) + alpha;
210 Double_t px1=TMath::Cos(phi), py1=TMath::Sin(phi), pz1=par[3];
211
212
213 cs1=TMath::Cos(alpha); sn1=TMath::Sin(alpha);
214 Double_t x1=x*cs1 - par[0]*sn1;
215 Double_t y1=x*sn1 + par[0]*cs1;
216 Double_t z1=par[1];
217
218
219 Double_t x2,y2,z2; // position and momentum of V0
220 Double_t px2,py2,pz2;
221
222 v->GetXYZ(x2,y2,z2);
223 v->GetPxPyPz(px2,py2,pz2);
224
225// calculation dca
226
227 Double_t dd= det(x2-x1,y2-y1,z2-z1,px1,py1,pz1,px2,py2,pz2);
228 Double_t ax= det(py1,pz1,py2,pz2);
229 Double_t ay=-det(px1,pz1,px2,pz2);
230 Double_t az= det(px1,py1,px2,py2);
231
232 Double_t dca=TMath::Abs(dd)/TMath::Sqrt(ax*ax + ay*ay + az*az);
233
234//points of the DCA
235 Double_t t1 = det(x2-x1,y2-y1,z2-z1,px2,py2,pz2,ax,ay,az)/
236 det(px1,py1,pz1,px2,py2,pz2,ax,ay,az);
237
238 x1 += px1*t1; y1 += py1*t1; //z1 += pz1*t1;
239
240
241 //propagate track to the points of DCA
242
243 x1=x1*cs1 + y1*sn1;
244 if (!t->PropagateTo(x1,0.,0.)) {
245 cerr<<"AliV0vertexer::PropagateToDCA: propagation failed !\n";
246 return 1.e+33;
247 }
248
249 return dca;
250}
251
252
253
254
255
256
257
258
259
260
261