1 //********************************************************************
2 // Example of the reconstruction that generates the ESD
4 // a) file containing the ITS clusters
5 // (the AliITSFindClustersV2.C macro can be used to generate it)
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)
13 // AliESDs.root containing the ESD events
15 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
16 //********************************************************************
18 #if !defined(__CINT__) || defined(__MAKECINT__)
19 #include <Riostream.h>
22 #include "TStopwatch.h"
27 #include "AliRunLoader.h"
28 #include "AliLoader.h"
29 #include "AliHeader.h"
30 #include "AliGenEventHeader.h"
33 #include "AliESDpid.h"
36 #include "AliITSgeom.h"
37 #include "AliITStrackerV2.h"
38 #include "AliV0vertexer.h"
39 #include "AliCascadeVertexer.h"
40 #include "AliITSpidESD.h"
41 #include "AliITSLoader.h"
43 #include "AliTPCParam.h"
44 #include "AliTPCtracker.h"
45 #include "AliTPCpidESD.h"
46 #include "AliTPCLoader.h"
48 #include "AliTRDtracker.h"
49 #include "AliTRDPartID.h"
51 #include "AliTOFpidESD.h"
53 #include "AliTOFGeometry.h"
56 extern AliRun *gAlice;
59 Int_t AliESDtest(Int_t nev=1,Int_t run=0) {
61 /**** Initialization of the NewIO *******/
64 delete gAlice->GetRunLoader();
69 AliRunLoader *rl = AliRunLoader::Open("galice.root");
71 cerr<<"Can not open session"<<endl;
74 Int_t retval = rl->LoadgAlice();
76 cerr<<"AliESDtest.C : LoadgAlice returned error"<<endl;
80 retval = rl->LoadHeader();
82 cerr<<"AliESDtest.C : LoadHeader returned error"<<endl;
86 gAlice=rl->GetAliRun();
89 AliKalmanTrack::SetConvConst(
90 1000/0.299792458/gAlice->Field()->SolenoidField()
94 /**** The ITS corner ********************/
96 AliITSLoader* itsl = (AliITSLoader*)rl->GetLoader("ITSLoader");
98 cerr<<"AliESDtest.C : Can not get the ITS loader"<<endl;
101 itsl->LoadRecPoints("read");
103 AliITS *dITS = (AliITS*)gAlice->GetDetector("ITS");
105 cerr<<"AliESDtest.C : Can not find the ITS detector !"<<endl;
108 AliITSgeom *geom = dITS->GetITSgeom();
110 //An instance of the ITS tracker
111 AliITStrackerV2 itsTracker(geom);
113 //An instance of the ITS PID maker
114 Double_t parITS[]={35.5,0.11,10.};
115 AliITSpidESD itsPID(parITS);
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
126 AliV0vertexer vtxer(cuts);
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
137 AliCascadeVertexer cvtxer=AliCascadeVertexer(cts);
139 /**** The TPC corner ********************/
141 AliTPCLoader* tpcl = (AliTPCLoader*)rl->GetLoader("TPCLoader");
143 cerr<<"AliESDtest.C : can not get the TPC loader"<<endl;
146 tpcl->LoadRecPoints("read");
149 AliTPCParam *par=(AliTPCParam *)gDirectory->Get("75x40_100x60_150x60");
151 cerr<<"TPC parameters have not been found !\n";
155 //An instance of the TPC tracker
156 AliTPCtracker tpcTracker(par);
158 //An instance of the TPC PID maker
159 Double_t parTPC[]={45.0,0.08,10.};
160 AliTPCpidESD tpcPID(parTPC);
163 /**** The TRD corner ********************/
165 AliLoader* trdl = rl->GetLoader("TRDLoader");
167 cerr<<"AliESDtest.C : can not get the TRD loader"<<endl;
170 trdl->LoadRecPoints("read");
172 //An instance of the TRD tracker
174 AliTRDtracker trdTracker(gFile); //galice.root file
177 //An instance of the TRD PID maker
178 TFile* pidFile = TFile::Open("pid.root");
179 if (!pidFile->IsOpen()) {
180 cerr << "Can't get pid.root !\n";
183 AliTRDPartID* trdPID = (AliTRDPartID*) pidFile->Get("AliTRDPartID");
185 cerr << "Can't get PID object !\n";
191 /**** The TOF corner ********************/
192 AliTOF *dTOF = (AliTOF*)gAlice->GetDetector("TOF");
194 cerr<<"AliESDtest.C : Can not find the TOF detector !"<<endl;
197 AliTOFGeometry *tofGeo = dTOF->GetGeometry();
199 cerr<<"AliESDtest.C : Can not find the TOF geometry !"<<endl;
203 AliLoader* tofl = rl->GetLoader("TOFLoader");
205 cerr<<"AliESDtest.C : can not get the TOF loader"<<endl;
208 tofl->LoadDigits("read");
210 //Instance of the TOF PID maker
211 Double_t parTOF[]={130.,5.};
212 AliTOFpidESD tofPID(parTOF);
215 //rl->UnloadgAlice();
218 TFile *ef=TFile::Open("AliESDs.root","RECREATE");
219 if (!ef->IsOpen()) {cerr<<"Can't AliESDs.root !\n"; return 1;}
223 if (nev>rl->GetNumberOfEvents()) nev=rl->GetNumberOfEvents();
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();
228 event->SetRunNumber(run);
229 event->SetEventNumber(i);
233 //***** Primary vertex reconstruction (MC vertex position, for the moment)
235 rl->GetHeader()->GenEventHeader()->PrimaryVertex(v);
236 Double_t vtx[3]={v[0],v[1],v[2]};
242 event->SetVertex(vtx,cvtx);
243 cvtx[1]=cvtx[0]; cvtx[2]=cvtx[5]; //trackers use only the diag.elements
245 //***** Initial path towards the primary vertex
246 tpcTracker.SetVertex(vtx,cvtx);
247 TTree *tpcTree=tpcl->TreeR();
249 cerr<<"Can't get the TPC cluster tree !\n";
252 tpcTracker.LoadClusters(tpcTree);
253 rc+=tpcTracker.Clusters2Tracks(event);
254 tpcPID.MakePID(event); // preliminary PID
255 AliESDpid::MakePID(event); // for the ITS tracker
257 itsTracker.SetVertex(vtx,cvtx);
258 TTree *itsTree=itsl->TreeR();
260 cerr<<"Can't get the ITS cluster tree !\n";
263 itsTracker.LoadClusters(itsTree);
264 rc+=itsTracker.Clusters2Tracks(event);
267 //***** Back propagation towards the outer barrel detectors
268 rc+=itsTracker.PropagateBack(event);
270 rc+=tpcTracker.PropagateBack(event);
272 TTree *trdTree=trdl->TreeR();
274 cerr<<"Can't get the TRD cluster tree !\n";
277 trdTracker.LoadClusters(trdTree);
278 rc+=trdTracker.PropagateBack(event);
281 for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++) {
282 AliESDtrack* track = event->GetTrack(iTrack);
283 trdPID->MakePID(track);
287 TTree *tofTree=tofl->TreeD();
289 cerr<<"Can't get the TOF cluster tree !\n";
292 tofPID.LoadClusters(tofTree,tofGeo);
293 tofPID.MakePID(event);
294 tofPID.UnloadClusters();
297 //***** Now the final refit at the primary vertex...
298 rc+=trdTracker.RefitInward(event);
299 trdTracker.UnloadClusters();
301 rc+=tpcTracker.RefitInward(event);
302 tpcTracker.UnloadClusters();
303 tpcPID.MakePID(event);
305 rc+=itsTracker.RefitInward(event);
306 itsTracker.UnloadClusters();
307 itsPID.MakePID(event);
310 //***** Here is the combined PID
311 AliESDpid::MakePID(event);
314 //***** Hyperon reconstruction
315 vtxer.SetVertex(vtx);
316 rc+=vtxer.Tracks2V0vertices(event); // V0 finding
317 rc+=cvtxer.V0sTracks2CascadeVertices(event); // cascade finding
321 //***** Some final manipulations with this event
324 sprintf(ename,"%d",i);
326 if (!event->Write(ename)) rc++;
329 cerr<<"Something bad happened...\n";
333 timer.Stop(); timer.Print();