]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDtest.C
Removing some temporary solutions for the track parameters at the primary vertex...
[u/mrichter/AliRoot.git] / STEER / AliESDtest.C
1 //********************************************************************
2 //     Example of the reconstruction that generates the ESD
3 // Input files: 
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)
12 // Ouput file:
13 //      AliESDs.root containing the ESD events 
14 //
15 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
16 //********************************************************************
17
18 #if !defined(__CINT__) || defined(__MAKECINT__)
19   #include <Riostream.h>
20   #include "TFile.h"
21   #include "TSystem.h"
22   #include "TStopwatch.h"
23   #include "TGeant3.h"
24
25   #include "AliMagF.h"
26   #include "AliRun.h"
27   #include "AliRunLoader.h"
28   #include "AliLoader.h"
29
30   #include "AliESD.h"
31   #include "AliESDpid.h"
32
33   #include "AliITS.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"
40
41   #include "AliTPCParam.h"
42   #include "AliTPCtracker.h"
43   #include "AliTPCpidESD.h"
44   #include "AliTPCLoader.h"
45
46   #include "AliTRDtracker.h"
47   #include "AliTRDPartID.h"
48
49   #include "AliTOFpidESD.h"
50 #endif
51
52 extern TSystem *gSystem;
53 extern AliRun *gAlice;
54 extern TFile *gFile;
55
56 Int_t AliESDtest(Int_t nev=1) {
57
58 /**** Initialization of the NewIO *******/
59
60    if (gAlice) {
61       delete gAlice->GetRunLoader();
62       delete gAlice; 
63       gAlice=0;
64    }
65
66    gSystem->Load("libgeant321");     // needed for the PID in TOF 
67    new TGeant3("");                  // must be re-done !
68
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;
84       return 2;
85    }
86    gAlice=rl->GetAliRun();
87        
88
89    AliKalmanTrack::SetConvConst(
90       1000/0.299792458/gAlice->Field()->SolenoidField()
91    );
92
93
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;
106       return 4;
107    }
108    AliITSgeom *geom = dITS->GetITSgeom();
109
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);}
114    
115    //An instance of the ITS PID maker
116    Double_t parITS[]={34.,0.15,10.};
117    AliITSpidESD itsPID(parITS);
118
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
127                    };
128    AliV0vertexer vtxer(cuts);
129
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
138                    };
139    AliCascadeVertexer cvtxer=AliCascadeVertexer(cts);
140
141 /**** The TPC corner ********************/
142
143    AliTPCLoader* tpcl = (AliTPCLoader*)rl->GetLoader("TPCLoader");
144    if (tpcl == 0x0) {
145       cerr<<"AliESDtest.C : can not get the TPC loader"<<endl;
146       return 5;
147    }
148    tpcl->LoadRecPoints("read");
149
150    rl->CdGAFile();
151    AliTPCParam *par=(AliTPCParam *)gDirectory->Get("75x40_100x60_150x60");
152    if (!par) { 
153       cerr<<"TPC parameters have not been found !\n";
154       return 6;
155    }
156    
157    //An instance of the TPC tracker
158    AliTPCtracker tpcTracker(par);
159
160    //An instance of the TPC PID maker
161    Double_t parTPC[]={47.,0.10,10.};
162    AliTPCpidESD tpcPID(parTPC);
163
164
165 /**** The TRD corner ********************/
166
167    AliLoader* trdl = rl->GetLoader("TRDLoader");
168    if (trdl == 0x0) {
169       cerr<<"AliESDtest.C : can not get the TRD loader"<<endl;
170       return 5;
171    }
172    trdl->LoadRecPoints("read");
173
174    //An instance of the TRD tracker
175    rl->CdGAFile();
176    AliTRDtracker trdTracker(gFile);  //galice.root file
177
178 /*
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";
183      return 7;
184    }
185    AliTRDPartID* trdPID = (AliTRDPartID*) pidFile->Get("AliTRDPartID");
186    if (!trdPID) {
187      cerr << "Can't get PID object !\n";
188      return 8;
189    }
190 */
191
192
193 /**** The TOF corner ********************/
194
195    AliLoader* tofl = rl->GetLoader("TOFLoader");
196    if (tofl == 0x0) {
197       cerr<<"AliESDtest.C : can not get the TOF loader"<<endl;
198       return 5;
199    }
200    tofl->LoadDigits("read");
201
202
203    //Instance of the TOF PID maker
204    Double_t parTOF[]={130.,5.};
205    AliTOFpidESD tofPID(parTOF);
206
207
208    //rl->UnloadgAlice();
209
210
211    TFile *ef=TFile::Open("AliESDs.root","RECREATE");
212    if (!ef->IsOpen()) {cerr<<"Can't AliESDs.root !\n"; return 1;}
213
214    TStopwatch timer;
215    Int_t rc=0;
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(); 
221
222      rl->GetEvent(i);
223  
224
225 //***** Initial path towards the primary vertex
226      TTree *tpcTree=tpcl->TreeR();
227      if (!tpcTree) {
228         cerr<<"Can't get the TPC cluster tree !\n";
229         return 4;
230      }     
231      tpcTracker.LoadClusters(tpcTree);
232      rc+=tpcTracker.Clusters2Tracks(event);
233
234      TTree *itsTree=itsl->TreeR();
235      if (!itsTree) {
236         cerr<<"Can't get the ITS cluster tree !\n";
237         return 4;
238      }     
239      itsTracker.LoadClusters(itsTree);
240      rc+=itsTracker.Clusters2Tracks(event);
241
242
243 //***** Back propagation towards the outer barrel detectors
244      rc+=itsTracker.PropagateBack(event); 
245      //itsPID.MakePID(event);
246      
247      rc+=tpcTracker.PropagateBack(event);
248      tpcPID.MakePID(event);
249
250      TTree *trdTree=trdl->TreeR();
251      if (!trdTree) {
252         cerr<<"Can't get the TRD cluster tree !\n";
253         return 4;
254      } 
255      trdTracker.LoadClusters(trdTree);
256      rc+=trdTracker.PropagateBack(event);
257 /*
258      for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++) {
259        AliESDtrack* track = event->GetTrack(iTrack);
260        trdPID->MakePID(track);
261      }
262 */
263      TTree *tofTree=tofl->TreeD();
264      if (!tofTree) {
265         cerr<<"Can't get the TOF cluster tree !\n";
266         return 4;
267      } 
268      tofPID.LoadClusters(tofTree);
269      tofPID.MakePID(event);
270      tofPID.UnloadClusters();
271
272
273
274 //***** Here is the combined PID
275      AliESDpid::MakePID(event);
276
277
278
279 //***** Now the final refit at the primary vertex...
280      rc+=trdTracker.RefitInward(event);
281      trdTracker.UnloadClusters();
282
283      rc+=tpcTracker.RefitInward(event);
284      tpcTracker.UnloadClusters();
285
286      rc+=itsTracker.RefitInward(event); 
287      itsTracker.UnloadClusters();
288
289
290
291 //***** Hyperon reconstruction 
292      rc+=vtxer.Tracks2V0vertices(event);            // V0 finding
293      rc+=cvtxer.V0sTracks2CascadeVertices(event);   // cascade finding
294
295
296
297 //***** Some final manipulations with this event 
298      if (rc==0) {
299         Char_t ename[100]; 
300         sprintf(ename,"%d",i);
301         ef->cd();
302         if (!event->Write(ename)) rc++;
303      } 
304      if (rc) {
305         cerr<<"Something bad happened...\n";
306      }
307      delete event;
308    }
309    timer.Stop(); timer.Print();
310
311    //pidFile->Close();
312
313    delete par;
314
315    ef->Close();
316
317    delete rl;
318
319    return rc;
320 }