]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliESDtest.C
Modifications needed by the HBT analysis (P.Skowronski)
[u/mrichter/AliRoot.git] / STEER / AliESDtest.C
CommitLineData
ae982df3 1//********************************************************************
2// Example of the reconstruction that generates the ESD
3// Input files:
c630aafd 4// a) file containing the ITS clusters
ae982df3 5// (the AliITSFindClustersV2.C macro can be used to generate it)
c630aafd 6// b) file containing the TPC clusters
7// (the AliTPCFindClusters.C macro can be used to generate it)
8// c) file containing the TRD clusters
9// (the AliTRDdigits2cluster.C macro can be used to generate it)
10// d) file containing the TOF digits
11// (the AliTOFSDigits2Digits.C macro can be used to generate it)
ae982df3 12// Ouput file:
13// AliESDs.root containing the ESD events
14//
15// Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
16//********************************************************************
17
c630aafd 18#if !defined(__CINT__) || defined(__MAKECINT__)
ae982df3 19 #include <Riostream.h>
20 #include "TFile.h"
c630aafd 21 #include "TSystem.h"
ae982df3 22 #include "TStopwatch.h"
873f1f73 23 #include "TArrayF.h"
c630aafd 24
25 #include "AliMagF.h"
26 #include "AliRun.h"
27 #include "AliRunLoader.h"
28 #include "AliLoader.h"
873f1f73 29 #include "AliHeader.h"
30 #include "AliGenEventHeader.h"
ae982df3 31
32 #include "AliESD.h"
8c6a71ab 33 #include "AliESDpid.h"
c630aafd 34
35 #include "AliITS.h"
ae982df3 36 #include "AliITSgeom.h"
37 #include "AliITStrackerV2.h"
e23730c7 38 #include "AliV0vertexer.h"
39 #include "AliCascadeVertexer.h"
c630aafd 40 #include "AliITSpidESD.h"
41 #include "AliITSLoader.h"
42
43 #include "AliTPCParam.h"
44 #include "AliTPCtracker.h"
45 #include "AliTPCpidESD.h"
46 #include "AliTPCLoader.h"
47
79e94bf8 48 #include "AliTRDtracker.h"
49 #include "AliTRDPartID.h"
c630aafd 50
51 #include "AliTOFpidESD.h"
8c29135e 52 #include "AliTOF.h"
53 #include "AliTOFGeometry.h"
ae982df3 54#endif
55
c630aafd 56extern AliRun *gAlice;
57extern TFile *gFile;
58
a6efa2c6 59Int_t AliESDtest(Int_t nev=1,Int_t run=0) {
c630aafd 60
61/**** Initialization of the NewIO *******/
79e94bf8 62
c630aafd 63 if (gAlice) {
64 delete gAlice->GetRunLoader();
65 delete gAlice;
66 gAlice=0;
67 }
68
c630aafd 69 AliRunLoader *rl = AliRunLoader::Open("galice.root");
70 if (rl == 0x0) {
71 cerr<<"Can not open session"<<endl;
72 return 1;
73 }
74 Int_t retval = rl->LoadgAlice();
75 if (retval) {
76 cerr<<"AliESDtest.C : LoadgAlice returned error"<<endl;
77 delete rl;
78 return 1;
79 }
80 retval = rl->LoadHeader();
81 if (retval) {
82 cerr<<"AliESDtest.C : LoadHeader returned error"<<endl;
83 delete rl;
ae982df3 84 return 2;
85 }
c630aafd 86 gAlice=rl->GetAliRun();
87
ae982df3 88
c630aafd 89 AliKalmanTrack::SetConvConst(
90 1000/0.299792458/gAlice->Field()->SolenoidField()
91 );
ae982df3 92
93
c630aafd 94/**** The ITS corner ********************/
95
96 AliITSLoader* itsl = (AliITSLoader*)rl->GetLoader("ITSLoader");
97 if (itsl == 0x0) {
98 cerr<<"AliESDtest.C : Can not get the ITS loader"<<endl;
99 return 3;
100 }
101 itsl->LoadRecPoints("read");
102
103 AliITS *dITS = (AliITS*)gAlice->GetDetector("ITS");
104 if (!dITS) {
105 cerr<<"AliESDtest.C : Can not find the ITS detector !"<<endl;
ae982df3 106 return 4;
107 }
c630aafd 108 AliITSgeom *geom = dITS->GetITSgeom();
ae982df3 109
110 //An instance of the ITS tracker
111 AliITStrackerV2 itsTracker(geom);
8c6a71ab 112
113 //An instance of the ITS PID maker
431be10d 114 Double_t parITS[]={35.5,0.11,10.};
8c6a71ab 115 AliITSpidESD itsPID(parITS);
79e94bf8 116
e23730c7 117 //An instance of the V0 finder
118 Double_t cuts[]={33, // max. allowed chi2
119 0.16,// min. allowed negative daughter's impact parameter
120 0.05,// min. allowed positive daughter's impact parameter
121 0.080,// max. allowed DCA between the daughter tracks
122 0.998,// max. allowed cosine of V0's pointing angle
123 0.9, // min. radius of the fiducial volume
124 2.9 // max. radius of the fiducial volume
125 };
126 AliV0vertexer vtxer(cuts);
127
128 Double_t cts[]={33., // max. allowed chi2
129 0.05, // min. allowed V0 impact parameter
130 0.008, // window around the Lambda mass
131 0.035, // min. allowed bachelor's impact parameter
132 0.10, // max. allowed DCA between a V0 and a track
133 0.9985, //max. allowed cosine of the cascade pointing angle
134 0.9, // min. radius of the fiducial volume
135 2.9 // max. radius of the fiducial volume
136 };
137 AliCascadeVertexer cvtxer=AliCascadeVertexer(cts);
c630aafd 138
139/**** The TPC corner ********************/
140
141 AliTPCLoader* tpcl = (AliTPCLoader*)rl->GetLoader("TPCLoader");
142 if (tpcl == 0x0) {
143 cerr<<"AliESDtest.C : can not get the TPC loader"<<endl;
144 return 5;
145 }
146 tpcl->LoadRecPoints("read");
147
148 rl->CdGAFile();
149 AliTPCParam *par=(AliTPCParam *)gDirectory->Get("75x40_100x60_150x60");
150 if (!par) {
151 cerr<<"TPC parameters have not been found !\n";
79e94bf8 152 return 6;
153 }
c630aafd 154
155 //An instance of the TPC tracker
156 AliTPCtracker tpcTracker(par);
157
158 //An instance of the TPC PID maker
431be10d 159 Double_t parTPC[]={45.0,0.08,10.};
c630aafd 160 AliTPCpidESD tpcPID(parTPC);
161
162
163/**** The TRD corner ********************/
164
165 AliLoader* trdl = rl->GetLoader("TRDLoader");
166 if (trdl == 0x0) {
167 cerr<<"AliESDtest.C : can not get the TRD loader"<<endl;
168 return 5;
169 }
170 trdl->LoadRecPoints("read");
79e94bf8 171
172 //An instance of the TRD tracker
c630aafd 173 rl->CdGAFile();
174 AliTRDtracker trdTracker(gFile); //galice.root file
79e94bf8 175
c630aafd 176/*
79e94bf8 177 //An instance of the TRD PID maker
c630aafd 178 TFile* pidFile = TFile::Open("pid.root");
179 if (!pidFile->IsOpen()) {
180 cerr << "Can't get pid.root !\n";
181 return 7;
182 }
183 AliTRDPartID* trdPID = (AliTRDPartID*) pidFile->Get("AliTRDPartID");
184 if (!trdPID) {
185 cerr << "Can't get PID object !\n";
186 return 8;
187 }
188*/
189
190
191/**** The TOF corner ********************/
8c29135e 192 AliTOF *dTOF = (AliTOF*)gAlice->GetDetector("TOF");
193 if (!dTOF) {
194 cerr<<"AliESDtest.C : Can not find the TOF detector !"<<endl;
195 return 4;
196 }
197 AliTOFGeometry *tofGeo = dTOF->GetGeometry();
198 if (!tofGeo) {
199 cerr<<"AliESDtest.C : Can not find the TOF geometry !"<<endl;
200 return 4;
201 }
c630aafd 202
203 AliLoader* tofl = rl->GetLoader("TOFLoader");
204 if (tofl == 0x0) {
205 cerr<<"AliESDtest.C : can not get the TOF loader"<<endl;
206 return 5;
207 }
208 tofl->LoadDigits("read");
209
c630aafd 210 //Instance of the TOF PID maker
211 Double_t parTOF[]={130.,5.};
212 AliTOFpidESD tofPID(parTOF);
213
214
215 //rl->UnloadgAlice();
216
79e94bf8 217
218 TFile *ef=TFile::Open("AliESDs.root","RECREATE");
ae982df3 219 if (!ef->IsOpen()) {cerr<<"Can't AliESDs.root !\n"; return 1;}
220
221 TStopwatch timer;
222 Int_t rc=0;
c630aafd 223 if (nev>rl->GetNumberOfEvents()) nev=rl->GetNumberOfEvents();
ae982df3 224 //The loop over events
225 for (Int_t i=0; i<nev; i++) {
226 cerr<<"\n\nProcessing event number : "<<i<<endl;
227 AliESD *event=new AliESD();
a6efa2c6 228 event->SetRunNumber(run);
229 event->SetEventNumber(i);
a866ac60 230 event->SetMagneticField(gAlice->Field()->SolenoidField());
231
c630aafd 232 rl->GetEvent(i);
233
873f1f73 234//***** Primary vertex reconstruction (MC vertex position, for the moment)
235 TArrayF v(3);
236 rl->GetHeader()->GenEventHeader()->PrimaryVertex(v);
237 Double_t vtx[3]={v[0],v[1],v[2]};
238 Double_t cvtx[6]={
239 0.005,
240 0.000, 0.005,
241 0.000, 0.000, 0.010
242 };
243 event->SetVertex(vtx,cvtx);
244 cvtx[1]=cvtx[0]; cvtx[2]=cvtx[5]; //trackers use only the diag.elements
bb2ceb1f 245
246//***** Initial path towards the primary vertex
873f1f73 247 tpcTracker.SetVertex(vtx,cvtx);
c630aafd 248 TTree *tpcTree=tpcl->TreeR();
249 if (!tpcTree) {
250 cerr<<"Can't get the TPC cluster tree !\n";
251 return 4;
252 }
253 tpcTracker.LoadClusters(tpcTree);
ae982df3 254 rc+=tpcTracker.Clusters2Tracks(event);
f0ab2c5f 255 tpcPID.MakePID(event); // preliminary PID
256 AliESDpid::MakePID(event); // for the ITS tracker
ae982df3 257
873f1f73 258 itsTracker.SetVertex(vtx,cvtx);
c630aafd 259 TTree *itsTree=itsl->TreeR();
260 if (!itsTree) {
bb2ceb1f 261 cerr<<"Can't get the ITS cluster tree !\n";
c630aafd 262 return 4;
263 }
264 itsTracker.LoadClusters(itsTree);
ae982df3 265 rc+=itsTracker.Clusters2Tracks(event);
266
e23730c7 267
bb2ceb1f 268//***** Back propagation towards the outer barrel detectors
ae982df3 269 rc+=itsTracker.PropagateBack(event);
ae982df3 270
271 rc+=tpcTracker.PropagateBack(event);
8c6a71ab 272
c630aafd 273 TTree *trdTree=trdl->TreeR();
274 if (!trdTree) {
bb2ceb1f 275 cerr<<"Can't get the TRD cluster tree !\n";
c630aafd 276 return 4;
277 }
278 trdTracker.LoadClusters(trdTree);
79e94bf8 279 rc+=trdTracker.PropagateBack(event);
431be10d 280
c630aafd 281/*
79e94bf8 282 for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++) {
283 AliESDtrack* track = event->GetTrack(iTrack);
284 trdPID->MakePID(track);
285 }
c630aafd 286*/
431be10d 287
c630aafd 288 TTree *tofTree=tofl->TreeD();
289 if (!tofTree) {
290 cerr<<"Can't get the TOF cluster tree !\n";
291 return 4;
292 }
8c29135e 293 tofPID.LoadClusters(tofTree,tofGeo);
c630aafd 294 tofPID.MakePID(event);
295 tofPID.UnloadClusters();
296
79e94bf8 297
bb2ceb1f 298//***** Now the final refit at the primary vertex...
299 rc+=trdTracker.RefitInward(event);
300 trdTracker.UnloadClusters();
301
302 rc+=tpcTracker.RefitInward(event);
303 tpcTracker.UnloadClusters();
431be10d 304 tpcPID.MakePID(event);
bb2ceb1f 305
306 rc+=itsTracker.RefitInward(event);
307 itsTracker.UnloadClusters();
431be10d 308 itsPID.MakePID(event);
bb2ceb1f 309
310
431be10d 311//***** Here is the combined PID
312 AliESDpid::MakePID(event);
313
bb2ceb1f 314
315//***** Hyperon reconstruction
873f1f73 316 vtxer.SetVertex(vtx);
bb2ceb1f 317 rc+=vtxer.Tracks2V0vertices(event); // V0 finding
318 rc+=cvtxer.V0sTracks2CascadeVertices(event); // cascade finding
319
320
321
322//***** Some final manipulations with this event
ae982df3 323 if (rc==0) {
324 Char_t ename[100];
325 sprintf(ename,"%d",i);
79e94bf8 326 ef->cd();
ae982df3 327 if (!event->Write(ename)) rc++;
328 }
329 if (rc) {
330 cerr<<"Something bad happened...\n";
331 }
332 delete event;
333 }
334 timer.Stop(); timer.Print();
335
c630aafd 336 //pidFile->Close();
337
ae982df3 338 delete par;
c630aafd 339
ae982df3 340 ef->Close();
341
c630aafd 342 delete rl;
343
ae982df3 344 return rc;
345}