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"
31 #include "AliESDpid.h"
34 #include "AliITSgeom.h"
35 #include "AliITStrackerV2.h"
36 #include "AliV0vertexer.h"
37 #include "AliCascadeVertexer.h"
38 #include "AliITSpidESD.h"
39 #include "AliITSLoader.h"
41 #include "AliTPCParam.h"
42 #include "AliTPCtracker.h"
43 #include "AliTPCpidESD.h"
44 #include "AliTPCLoader.h"
46 #include "AliTRDtracker.h"
47 #include "AliTRDPartID.h"
49 #include "AliTOFpidESD.h"
52 extern TSystem *gSystem;
53 extern AliRun *gAlice;
56 Int_t AliESDtest(Int_t nev=1) {
58 /**** Initialization of the NewIO *******/
61 delete gAlice->GetRunLoader();
66 gSystem->Load("libgeant321"); // needed for the PID in TOF
67 new TGeant3(""); // must be re-done !
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);
112 {Double_t xyz[]={0.,0.,0.}, ers[]={0.005, 0.005, 0.010};
113 itsTracker.SetVertex(xyz,ers);}
115 //An instance of the ITS PID maker
116 Double_t parITS[]={34.,0.15,10.};
117 AliITSpidESD itsPID(parITS);
119 //An instance of the V0 finder
120 Double_t cuts[]={33, // max. allowed chi2
121 0.16,// min. allowed negative daughter's impact parameter
122 0.05,// min. allowed positive daughter's impact parameter
123 0.080,// max. allowed DCA between the daughter tracks
124 0.998,// max. allowed cosine of V0's pointing angle
125 0.9, // min. radius of the fiducial volume
126 2.9 // max. radius of the fiducial volume
128 AliV0vertexer vtxer(cuts);
130 Double_t cts[]={33., // max. allowed chi2
131 0.05, // min. allowed V0 impact parameter
132 0.008, // window around the Lambda mass
133 0.035, // min. allowed bachelor's impact parameter
134 0.10, // max. allowed DCA between a V0 and a track
135 0.9985, //max. allowed cosine of the cascade pointing angle
136 0.9, // min. radius of the fiducial volume
137 2.9 // max. radius of the fiducial volume
139 AliCascadeVertexer cvtxer=AliCascadeVertexer(cts);
141 /**** The TPC corner ********************/
143 AliTPCLoader* tpcl = (AliTPCLoader*)rl->GetLoader("TPCLoader");
145 cerr<<"AliESDtest.C : can not get the TPC loader"<<endl;
148 tpcl->LoadRecPoints("read");
151 AliTPCParam *par=(AliTPCParam *)gDirectory->Get("75x40_100x60_150x60");
153 cerr<<"TPC parameters have not been found !\n";
157 //An instance of the TPC tracker
158 AliTPCtracker tpcTracker(par);
160 //An instance of the TPC PID maker
161 Double_t parTPC[]={47.,0.10,10.};
162 AliTPCpidESD tpcPID(parTPC);
165 /**** The TRD corner ********************/
167 AliLoader* trdl = rl->GetLoader("TRDLoader");
169 cerr<<"AliESDtest.C : can not get the TRD loader"<<endl;
172 trdl->LoadRecPoints("read");
174 //An instance of the TRD tracker
176 AliTRDtracker trdTracker(gFile); //galice.root file
179 //An instance of the TRD PID maker
180 TFile* pidFile = TFile::Open("pid.root");
181 if (!pidFile->IsOpen()) {
182 cerr << "Can't get pid.root !\n";
185 AliTRDPartID* trdPID = (AliTRDPartID*) pidFile->Get("AliTRDPartID");
187 cerr << "Can't get PID object !\n";
193 /**** The TOF corner ********************/
195 AliLoader* tofl = rl->GetLoader("TOFLoader");
197 cerr<<"AliESDtest.C : can not get the TOF loader"<<endl;
200 tofl->LoadDigits("read");
203 //Instance of the TOF PID maker
204 Double_t parTOF[]={130.,5.};
205 AliTOFpidESD tofPID(parTOF);
208 //rl->UnloadgAlice();
211 TFile *ef=TFile::Open("AliESDs.root","RECREATE");
212 if (!ef->IsOpen()) {cerr<<"Can't AliESDs.root !\n"; return 1;}
216 if (nev>rl->GetNumberOfEvents()) nev=rl->GetNumberOfEvents();
217 //The loop over events
218 for (Int_t i=0; i<nev; i++) {
219 cerr<<"\n\nProcessing event number : "<<i<<endl;
220 AliESD *event=new AliESD();
225 //***** Initial path towards the primary vertex
226 TTree *tpcTree=tpcl->TreeR();
228 cerr<<"Can't get the TPC cluster tree !\n";
231 tpcTracker.LoadClusters(tpcTree);
232 rc+=tpcTracker.Clusters2Tracks(event);
234 TTree *itsTree=itsl->TreeR();
236 cerr<<"Can't get the ITS cluster tree !\n";
239 itsTracker.LoadClusters(itsTree);
240 rc+=itsTracker.Clusters2Tracks(event);
243 //***** Back propagation towards the outer barrel detectors
244 rc+=itsTracker.PropagateBack(event);
245 //itsPID.MakePID(event);
247 rc+=tpcTracker.PropagateBack(event);
248 tpcPID.MakePID(event);
250 TTree *trdTree=trdl->TreeR();
252 cerr<<"Can't get the TRD cluster tree !\n";
255 trdTracker.LoadClusters(trdTree);
256 rc+=trdTracker.PropagateBack(event);
258 for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++) {
259 AliESDtrack* track = event->GetTrack(iTrack);
260 trdPID->MakePID(track);
263 TTree *tofTree=tofl->TreeD();
265 cerr<<"Can't get the TOF cluster tree !\n";
268 tofPID.LoadClusters(tofTree);
269 tofPID.MakePID(event);
270 tofPID.UnloadClusters();
274 //***** Here is the combined PID
275 AliESDpid::MakePID(event);
279 //***** Now the final refit at the primary vertex...
280 rc+=trdTracker.RefitInward(event);
281 trdTracker.UnloadClusters();
283 rc+=tpcTracker.RefitInward(event);
284 tpcTracker.UnloadClusters();
286 rc+=itsTracker.RefitInward(event);
287 itsTracker.UnloadClusters();
291 //***** Hyperon reconstruction
292 rc+=vtxer.Tracks2V0vertices(event); // V0 finding
293 rc+=cvtxer.V0sTracks2CascadeVertices(event); // cascade finding
297 //***** Some final manipulations with this event
300 sprintf(ename,"%d",i);
302 if (!event->Write(ename)) rc++;
305 cerr<<"Something bad happened...\n";
309 timer.Stop(); timer.Print();