]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliReconstruction.cxx
544c87cf1fc2f57edeaf9fb5393d2ac5dae90c0a
[u/mrichter/AliRoot.git] / STEER / AliReconstruction.cxx
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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 // class for running the reconstruction                                      //
21 //                                                                           //
22 // Clusters and tracks are created for all detectors and all events by       //
23 // typing:                                                                   //
24 //                                                                           //
25 //   AliReconstruction rec;                                                  //
26 //   rec.Run();                                                              //
27 //                                                                           //
28 // The Run method returns kTRUE in case of successful execution.             //
29 // The name of the galice file can be changed from the default               //
30 // "galice.root" by passing it as argument to the AliReconstruction          //
31 // constructor or by                                                         //
32 //                                                                           //
33 //   rec.SetGAliceFile("...");                                               //
34 //                                                                           //
35 // The local reconstruction can be switched on or off for individual         //
36 // detectors by                                                              //
37 //                                                                           //
38 //   rec.SetRunLocalReconstruction("...");                                   //
39 //                                                                           //
40 // The argument is a (case sensitive) string with the names of the           //
41 // detectors separated by a space. The special string "ALL" selects all      //
42 // available detectors. This is the default.                                 //
43 //                                                                           //
44 // The tracking in ITS, TPC and TRD and the creation of ESD tracks can be    //
45 // switched off by                                                           //
46 //                                                                           //
47 //   rec.SetRunTracking(kFALSE);                                             //
48 //                                                                           //
49 // The filling of additional ESD information can be steered by               //
50 //                                                                           //
51 //   rec.SetFillESD("...");                                                  //
52 //                                                                           //
53 // Again, the string specifies the list of detectors. The default is "ALL".  //
54 //                                                                           //
55 // The reconstruction requires digits as input. For the creation of digits   //
56 // have a look at the class AliSimulation.                                   //
57 //                                                                           //
58 // For debug purposes the method SetCheckPointLevel can be used. If the      //
59 // argument is greater than 0, files with ESD events will be written after   //
60 // selected steps of the reconstruction for each event:                      //
61 //   level 1: after tracking and after filling of ESD (final)                //
62 //   level 2: in addition after each tracking step                           //
63 //   level 3: in addition after the filling of ESD for each detector         //
64 // If a final check point file exists for an event, this event will be       //
65 // skipped in the reconstruction. The tracking and the filling of ESD for    //
66 // a detector will be skipped as well, if the corresponding check point      //
67 // file exists. The ESD event will then be loaded from the file instead.     //
68 //                                                                           //
69 ///////////////////////////////////////////////////////////////////////////////
70
71 #include <TArrayF.h>
72 #include <TFile.h>
73 #include <TSystem.h>
74 #include <TROOT.h>
75 #include <TPluginManager.h>
76 #include <TStopwatch.h>
77
78 #include "AliReconstruction.h"
79 #include "AliRunLoader.h"
80 #include "AliRun.h"
81 #include "AliTracker.h"
82 #include "AliESD.h"
83 #include "AliESDVertex.h"
84 #include "AliVertexer.h"
85 #include "AliHeader.h"
86 #include "AliGenEventHeader.h"
87 #include "AliESDpid.h"
88 #include "AliMagF.h"
89
90 ClassImp(AliReconstruction)
91
92
93 //_____________________________________________________________________________
94 const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "RICH", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "START", "VZERO", "CRT", "HLT"};
95
96 //_____________________________________________________________________________
97 AliReconstruction::AliReconstruction(const char* gAliceFilename,
98                                      const char* name, const char* title) :
99   TNamed(name, title),
100
101   fRunLocalReconstruction("ALL"),
102   fRunVertexFinder(kTRUE),
103   fRunTracking(kTRUE),
104   fFillESD("ALL"),
105   fGAliceFileName(gAliceFilename),
106   fStopOnError(kFALSE),
107   fCheckPointLevel(0),
108
109   fRunLoader(NULL),
110   fITSLoader(NULL),
111   fITSVertexer(NULL),
112   fITSTracker(NULL),
113   fTPCLoader(NULL),
114   fTPCTracker(NULL),
115   fTRDLoader(NULL),
116   fTRDTracker(NULL),
117   fTOFLoader(NULL),
118   fTOFTracker(NULL)
119 {
120 // create reconstruction object with default parameters
121
122 }
123
124 //_____________________________________________________________________________
125 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
126   TNamed(rec),
127
128   fRunLocalReconstruction(rec.fRunLocalReconstruction),
129   fRunVertexFinder(rec.fRunVertexFinder),
130   fRunTracking(rec.fRunTracking),
131   fFillESD(rec.fFillESD),
132   fGAliceFileName(rec.fGAliceFileName),
133   fStopOnError(rec.fStopOnError),
134   fCheckPointLevel(0),
135
136   fRunLoader(NULL),
137   fITSLoader(NULL),
138   fITSVertexer(NULL),
139   fITSTracker(NULL),
140   fTPCLoader(NULL),
141   fTPCTracker(NULL),
142   fTRDLoader(NULL),
143   fTRDTracker(NULL),
144   fTOFLoader(NULL),
145   fTOFTracker(NULL)
146 {
147 // copy constructor
148
149 }
150
151 //_____________________________________________________________________________
152 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
153 {
154 // assignment operator
155
156   this->~AliReconstruction();
157   new(this) AliReconstruction(rec);
158   return *this;
159 }
160
161 //_____________________________________________________________________________
162 AliReconstruction::~AliReconstruction()
163 {
164 // clean up
165
166   CleanUp();
167 }
168
169
170 //_____________________________________________________________________________
171 void AliReconstruction::SetGAliceFile(const char* fileName)
172 {
173 // set the name of the galice file
174
175   fGAliceFileName = fileName;
176 }
177
178
179 //_____________________________________________________________________________
180 Bool_t AliReconstruction::Run()
181 {
182 // run the reconstruction
183
184   // open the run loader
185   fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
186   if (!fRunLoader) {
187     Error("Run", "no run loader found in file %s", 
188           fGAliceFileName.Data());
189     CleanUp();
190     return kFALSE;
191   }
192   fRunLoader->LoadgAlice();
193   AliRun* aliRun = fRunLoader->GetAliRun();
194   if (!aliRun) {
195     Error("Run", "no gAlice object found in file %s", 
196           fGAliceFileName.Data());
197     CleanUp();
198     return kFALSE;
199   }
200   gAlice = aliRun;
201
202   // load the reconstructor objects
203   TPluginManager* pluginManager = gROOT->GetPluginManager();
204   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
205     TString detName = fgkDetectorName[iDet];
206     TString recName = "Ali" + detName + "Reconstructor";
207     if (!gAlice->GetDetector(detName) && detName != "HLT") continue;
208
209     if(detName == "HLT") {
210       if (!gROOT->GetClass("AliLevel3")) {
211         gSystem->Load("libAliL3Src.so");
212         gSystem->Load("libAliL3Misc.so");
213         gSystem->Load("libAliL3Hough.so");
214         gSystem->Load("libAliL3Comp.so");
215       }
216     }
217
218     AliReconstructor* reconstructor = NULL;
219     // first check if a plugin is defined for the reconstructor
220     TPluginHandler* pluginHandler = 
221       pluginManager->FindHandler("AliReconstructor", detName);
222     // if not, but the reconstructor class is implemented, add a plugin for it
223     if (!pluginHandler && gROOT->GetClass(recName.Data())) {
224       Info("Run", "defining plugin for %s", recName.Data());
225       pluginManager->AddHandler("AliReconstructor", detName, 
226                                 recName, detName, recName + "()");
227       pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
228     }
229     if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
230       reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
231     }
232     // if there is no reconstructor class for the detector use the dummy one
233     if (!reconstructor && gAlice->GetDetector(detName)) {
234       Info("Run", "using dummy reconstructor for %s", detName.Data());
235       reconstructor = new AliDummyReconstructor(gAlice->GetDetector(detName));
236     }
237     if (reconstructor) fReconstructors.Add(reconstructor);
238   }
239
240   // local reconstruction
241   if (!fRunLocalReconstruction.IsNull()) {
242     if (!RunLocalReconstruction(fRunLocalReconstruction)) {
243       if (fStopOnError) {CleanUp(); return kFALSE;}
244     }
245   }
246   if (!fRunVertexFinder && !fRunTracking && fFillESD.IsNull()) return kTRUE;
247
248   // get vertexer
249   if (fRunVertexFinder && !CreateVertexer()) {
250     if (fStopOnError) {
251       CleanUp(); 
252       return kFALSE;
253     }
254   }
255
256   // get loaders and trackers
257   if (fRunTracking && !CreateTrackers()) {
258     if (fStopOnError) {
259       CleanUp(); 
260       return kFALSE;
261     }      
262   }
263
264   // create the ESD output file and tree
265   TFile* file = TFile::Open("AliESDs.root", "RECREATE");
266   if (!file->IsOpen()) {
267     Error("Run", "opening AliESDs.root failed");
268     if (fStopOnError) {CleanUp(file); return kFALSE;}    
269   }
270   AliESD* esd = new AliESD;
271   TTree* tree = new TTree("esdTree", "Tree with ESD objects");
272   tree->Branch("ESD", "AliESD", &esd);
273   delete esd;
274   gROOT->cd();
275
276   // loop over events
277   for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) {
278     Info("Run", "processing event %d", iEvent);
279     fRunLoader->GetEvent(iEvent);
280
281     char fileName[256];
282     sprintf(fileName, "ESD_%d.%d_final.root", 
283             aliRun->GetRunNumber(), aliRun->GetEvNumber());
284     if (!gSystem->AccessPathName(fileName)) continue;
285
286     esd = new AliESD;
287     esd->SetRunNumber(aliRun->GetRunNumber());
288     esd->SetEventNumber(aliRun->GetEvNumber());
289     esd->SetMagneticField(aliRun->Field()->SolenoidField());
290
291     // vertex finder
292     if (fRunVertexFinder) {
293       if (!ReadESD(esd, "vertex")) {
294         if (!RunVertexFinder(esd)) {
295           if (fStopOnError) {CleanUp(file); return kFALSE;}
296         }
297         if (fCheckPointLevel > 0) WriteESD(esd, "vertex");
298       }
299     }
300
301     // barrel tracking
302     if (fRunTracking) {
303       if (!ReadESD(esd, "tracking")) {
304         if (!RunTracking(esd)) {
305           if (fStopOnError) {CleanUp(file); return kFALSE;}
306         }
307         if (fCheckPointLevel > 0) WriteESD(esd, "tracking");
308       }
309     }
310
311     // fill ESD
312     if (!fFillESD.IsNull()) {
313       if (!FillESD(esd, fFillESD)) {
314         if (fStopOnError) {CleanUp(file); return kFALSE;}
315       }
316     }
317
318     // combined PID
319     AliESDpid::MakePID(esd);
320     if (fCheckPointLevel > 1) WriteESD(esd, "PID");
321
322     // write ESD
323     tree->Fill();
324
325     if (fCheckPointLevel > 0) WriteESD(esd, "final");
326     delete esd;
327   }
328
329   file->cd();
330   tree->Write();
331   CleanUp(file);
332
333   return kTRUE;
334 }
335
336
337 //_____________________________________________________________________________
338 Bool_t AliReconstruction::RunLocalReconstruction(const TString& detectors)
339 {
340 // run the local reconstruction
341
342   TStopwatch stopwatch;
343   stopwatch.Start();
344
345   TString detStr = detectors;
346   for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
347     AliReconstructor* reconstructor = 
348       (AliReconstructor*) fReconstructors[iDet];
349     TString detName = reconstructor->GetDetectorName();
350     if (IsSelected(detName, detStr)) {
351       Info("RunReconstruction", "running reconstruction for %s", 
352            detName.Data());
353       TStopwatch stopwatchDet;
354       stopwatchDet.Start();
355       reconstructor->Reconstruct(fRunLoader);
356       Info("RunReconstruction", "execution time for %s:", detName.Data());
357       stopwatchDet.Print();
358     }
359   }
360
361   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
362     Error("RunReconstruction", "the following detectors were not found: %s", 
363           detStr.Data());
364     if (fStopOnError) return kFALSE;
365   }
366
367   Info("RunReconstruction", "execution time:");
368   stopwatch.Print();
369
370   return kTRUE;
371 }
372
373 //_____________________________________________________________________________
374 Bool_t AliReconstruction::RunVertexFinder(AliESD*& esd)
375 {
376 // run the barrel tracking
377
378   TStopwatch stopwatch;
379   stopwatch.Start();
380
381   AliESDVertex* vertex = NULL;
382   Double_t vtxPos[3] = {0, 0, 0};
383   Double_t vtxErr[3] = {0.07, 0.07, 0.1};
384   TArrayF mcVertex(3); 
385   if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
386     fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
387     for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
388   }
389
390   if (fITSVertexer) {
391     Info("RunVertexFinder", "running the ITS vertex finder");
392     fITSVertexer->SetDebug(1);
393     vertex = fITSVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
394     if(!vertex){
395       Warning("RunVertexFinder","Vertex not found \n");
396       vertex = new AliESDVertex();
397     }
398     else {
399       vertex->SetTruePos(vtxPos);  // store also the vertex from MC
400     }
401
402   } else {
403     Info("RunVertexFinder", "getting the primary vertex from MC");
404     vertex = new AliESDVertex(vtxPos, vtxErr);
405   }
406
407   if (vertex) {
408     vertex->GetXYZ(vtxPos);
409     vertex->GetSigmaXYZ(vtxErr);
410   } else {
411     Warning("RunVertexFinder", "no vertex reconstructed");
412     vertex = new AliESDVertex(vtxPos, vtxErr);
413   }
414   esd->SetVertex(vertex);
415   if (fITSTracker) fITSTracker->SetVertex(vtxPos, vtxErr);
416   if (fTPCTracker) fTPCTracker->SetVertex(vtxPos, vtxErr);
417   if (fTRDTracker) fTRDTracker->SetVertex(vtxPos, vtxErr);
418   delete vertex;
419
420   Info("RunVertexFinder", "execution time:");
421   stopwatch.Print();
422
423   return kTRUE;
424 }
425
426 //_____________________________________________________________________________
427 Bool_t AliReconstruction::RunTracking(AliESD*& esd)
428 {
429 // run the barrel tracking
430
431   TStopwatch stopwatch;
432   stopwatch.Start();
433
434   if (!fTPCTracker) {
435     Error("RunTracking", "no TPC tracker");
436     return kFALSE;
437   }
438
439   // TPC tracking
440   Info("RunTracking", "TPC tracking");
441   fTPCLoader->LoadRecPoints("read");
442   TTree* tpcTree = fTPCLoader->TreeR();
443   if (!tpcTree) {
444     Error("RunTracking", "Can't get the TPC cluster tree");
445     return kFALSE;
446   }
447   fTPCTracker->LoadClusters(tpcTree);
448   if (fTPCTracker->Clusters2Tracks(esd) != 0) {
449     Error("RunTracking", "TPC Clusters2Tracks failed");
450     return kFALSE;
451   }
452   if (fCheckPointLevel > 1) WriteESD(esd, "TPC.tracking");
453
454   if (!fITSTracker) {
455     Warning("RunTracking", "no ITS tracker");
456   } else {
457
458     fRunLoader->GetAliRun()->GetDetector("TPC")->FillESD(esd); // preliminary
459     AliESDpid::MakePID(esd);                  // PID for the ITS tracker
460
461     // ITS tracking
462     Info("RunTracking", "ITS tracking");
463     fITSLoader->LoadRecPoints("read");
464     TTree* itsTree = fITSLoader->TreeR();
465     if (!itsTree) {
466       Error("RunTracking", "Can't get the ITS cluster tree");
467       return kFALSE;
468     }
469     fITSTracker->LoadClusters(itsTree);
470     if (fITSTracker->Clusters2Tracks(esd) != 0) {
471       Error("RunTracking", "ITS Clusters2Tracks failed");
472       return kFALSE;
473     }
474     if (fCheckPointLevel > 1) WriteESD(esd, "ITS.tracking");
475
476     if (!fTRDTracker) {
477       Warning("RunTracking", "no TRD tracker");
478     } else {
479       // ITS back propagation
480       Info("RunTracking", "ITS back propagation");
481       if (fITSTracker->PropagateBack(esd) != 0) {
482         Error("RunTracking", "ITS backward propagation failed");
483         return kFALSE;
484       }
485       if (fCheckPointLevel > 1) WriteESD(esd, "ITS.back");
486
487       // TPC back propagation
488       Info("RunTracking", "TPC back propagation");
489       if (fTPCTracker->PropagateBack(esd) != 0) {
490         Error("RunTracking", "TPC backward propagation failed");
491         return kFALSE;
492       }
493       if (fCheckPointLevel > 1) WriteESD(esd, "TPC.back");
494
495       // TRD back propagation
496       Info("RunTracking", "TRD back propagation");
497       fTRDLoader->LoadRecPoints("read");
498       TTree* trdTree = fTRDLoader->TreeR();
499       if (!trdTree) {
500         Error("RunTracking", "Can't get the TRD cluster tree");
501         return kFALSE;
502       }
503       fTRDTracker->LoadClusters(trdTree);
504       if (fTRDTracker->PropagateBack(esd) != 0) {
505         Error("RunTracking", "TRD backward propagation failed");
506         return kFALSE;
507       }
508       if (fCheckPointLevel > 1) WriteESD(esd, "TRD.back");
509
510       if (!fTOFTracker) {
511         Warning("RunTracking", "no TOF tracker");
512       } else {
513         // TOF back propagation
514         Info("RunTracking", "TOF back propagation");
515         fTOFLoader->LoadDigits("read");
516         TTree* tofTree = fTOFLoader->TreeD();
517         if (!tofTree) {
518           Error("RunTracking", "Can't get the TOF digits tree");
519           return kFALSE;
520         }
521         fTOFTracker->LoadClusters(tofTree);
522         if (fTOFTracker->PropagateBack(esd) != 0) {
523           Error("RunTracking", "TOF backward propagation failed");
524           return kFALSE;
525         }
526         if (fCheckPointLevel > 1) WriteESD(esd, "TOF.back");
527         fTOFTracker->UnloadClusters();
528         fTOFLoader->UnloadDigits();
529       }
530
531       // TRD inward refit
532       Info("RunTracking", "TRD inward refit");
533       if (fTRDTracker->RefitInward(esd) != 0) {
534         Error("RunTracking", "TRD inward refit failed");
535         return kFALSE;
536       }
537       if (fCheckPointLevel > 1) WriteESD(esd, "TRD.refit");
538       fTRDTracker->UnloadClusters();
539       fTRDLoader->UnloadRecPoints();
540     
541       // TPC inward refit
542       Info("RunTracking", "TPC inward refit");
543       if (fTPCTracker->RefitInward(esd) != 0) {
544         Error("RunTracking", "TPC inward refit failed");
545         return kFALSE;
546       }
547       if (fCheckPointLevel > 1) WriteESD(esd, "TPC.refit");
548     
549       // ITS inward refit
550       Info("RunTracking", "ITS inward refit");
551       if (fITSTracker->RefitInward(esd) != 0) {
552         Error("RunTracking", "ITS inward refit failed");
553         return kFALSE;
554       }
555       if (fCheckPointLevel > 1) WriteESD(esd, "ITS.refit");
556
557     }  // if TRD tracker
558     fITSTracker->UnloadClusters();
559     fITSLoader->UnloadRecPoints();
560
561   }  // if ITS tracker
562   fTPCTracker->UnloadClusters();
563   fTPCLoader->UnloadRecPoints();
564
565   Info("RunTracking", "execution time:");
566   stopwatch.Print();
567
568   return kTRUE;
569 }
570
571 //_____________________________________________________________________________
572 Bool_t AliReconstruction::FillESD(AliESD*& esd, const TString& detectors)
573 {
574 // fill the event summary data
575
576   TStopwatch stopwatch;
577   stopwatch.Start();
578
579   TString detStr = detectors;
580   for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
581     AliReconstructor* reconstructor = 
582       (AliReconstructor*) fReconstructors[iDet];
583     TString detName = reconstructor->GetDetectorName();
584     if (IsSelected(detName, detStr)) {
585       if (!ReadESD(esd, detName.Data())) {
586         Info("FillESD", "filling ESD for %s", detName.Data());
587         reconstructor->FillESD(fRunLoader, esd);
588         if (fCheckPointLevel > 2) WriteESD(esd, detName.Data());
589       }
590     }
591   }
592
593   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
594     Error("FillESD", "the following detectors were not found: %s", 
595           detStr.Data());
596     if (fStopOnError) return kFALSE;
597   }
598
599   Info("FillESD", "execution time:");
600   stopwatch.Print();
601
602   return kTRUE;
603 }
604
605
606 //_____________________________________________________________________________
607 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
608 {
609 // check whether detName is contained in detectors
610 // if yes, it is removed from detectors
611
612   // check if all detectors are selected
613   if ((detectors.CompareTo("ALL") == 0) ||
614       detectors.BeginsWith("ALL ") ||
615       detectors.EndsWith(" ALL") ||
616       detectors.Contains(" ALL ")) {
617     detectors = "ALL";
618     return kTRUE;
619   }
620
621   // search for the given detector
622   Bool_t result = kFALSE;
623   if ((detectors.CompareTo(detName) == 0) ||
624       detectors.BeginsWith(detName+" ") ||
625       detectors.EndsWith(" "+detName) ||
626       detectors.Contains(" "+detName+" ")) {
627     detectors.ReplaceAll(detName, "");
628     result = kTRUE;
629   }
630
631   // clean up the detectors string
632   while (detectors.Contains("  ")) detectors.ReplaceAll("  ", " ");
633   while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
634   while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
635
636   return result;
637 }
638
639 //_____________________________________________________________________________
640 AliReconstructor* AliReconstruction::GetReconstructor(const char* detName) const
641 {
642 // get the reconstructor object for a detector
643
644   for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
645     AliReconstructor* reconstructor = 
646       (AliReconstructor*) fReconstructors[iDet];
647     if (strcmp(reconstructor->GetDetectorName(), detName) == 0) {
648       return reconstructor;
649     }
650   }
651   return NULL;
652 }
653
654 //_____________________________________________________________________________
655 Bool_t AliReconstruction::CreateVertexer()
656 {
657 // create the vertexer
658
659   fITSVertexer = NULL;
660   AliReconstructor* itsReconstructor = GetReconstructor("ITS");
661   if (itsReconstructor) {
662     fITSVertexer = itsReconstructor->CreateVertexer(fRunLoader);
663   }
664   if (!fITSVertexer) {
665     Warning("CreateVertexer", "couldn't create a vertexer for ITS");
666     if (fStopOnError) return kFALSE;
667   }
668
669   return kTRUE;
670 }
671
672 //_____________________________________________________________________________
673 Bool_t AliReconstruction::CreateTrackers()
674 {
675 // get the loaders and create the trackers
676
677   fITSTracker = NULL;
678   fITSLoader = fRunLoader->GetLoader("ITSLoader");
679   if (!fITSLoader) {
680     Warning("CreateTrackers", "no ITS loader found");
681     if (fStopOnError) return kFALSE;
682   } else {
683     AliReconstructor* itsReconstructor = GetReconstructor("ITS");
684     if (itsReconstructor) {
685       fITSTracker = itsReconstructor->CreateTracker(fRunLoader);
686     }
687     if (!fITSTracker) {
688       Warning("CreateTrackers", "couldn't create a tracker for ITS");
689       if (fStopOnError) return kFALSE;
690     }
691   }
692     
693   fTPCTracker = NULL;
694   fTPCLoader = fRunLoader->GetLoader("TPCLoader");
695   if (!fTPCLoader) {
696     Error("CreateTrackers", "no TPC loader found");
697     if (fStopOnError) return kFALSE;
698   } else {
699     AliReconstructor* tpcReconstructor = GetReconstructor("TPC");
700     if (tpcReconstructor) {
701       fTPCTracker = tpcReconstructor->CreateTracker(fRunLoader);
702     }
703     if (!fTPCTracker) {
704       Error("CreateTrackers", "couldn't create a tracker for TPC");
705       if (fStopOnError) return kFALSE;
706     }
707   }
708     
709   fTRDTracker = NULL;
710   fTRDLoader = fRunLoader->GetLoader("TRDLoader");
711   if (!fTRDLoader) {
712     Warning("CreateTrackers", "no TRD loader found");
713     if (fStopOnError) return kFALSE;
714   } else {
715     AliReconstructor* trdReconstructor = GetReconstructor("TRD");
716     if (trdReconstructor) {
717       fTRDTracker = trdReconstructor->CreateTracker(fRunLoader);
718     }
719     if (!fTRDTracker) {
720       Warning("CreateTrackers", "couldn't create a tracker for TRD");
721       if (fStopOnError) return kFALSE;
722     }
723   }
724     
725   fTOFTracker = NULL;
726   fTOFLoader = fRunLoader->GetLoader("TOFLoader");
727   if (!fTOFLoader) {
728     Warning("CreateTrackers", "no TOF loader found");
729     if (fStopOnError) return kFALSE;
730   } else {
731     AliReconstructor* tofReconstructor = GetReconstructor("TOF");
732     if (tofReconstructor) {
733       fTOFTracker = tofReconstructor->CreateTracker(fRunLoader);
734     }
735     if (!fTOFTracker) {
736       Warning("CreateTrackers", "couldn't create a tracker for TOF");
737       if (fStopOnError) return kFALSE;
738     }
739   }
740
741   return kTRUE;
742 }
743
744 //_____________________________________________________________________________
745 void AliReconstruction::CleanUp(TFile* file)
746 {
747 // delete trackers and the run loader and close and delete the file
748
749   fReconstructors.Delete();
750
751   delete fITSVertexer;
752   fITSVertexer = NULL;
753   delete fITSTracker;
754   fITSTracker = NULL;
755   delete fTPCTracker;
756   fTPCTracker = NULL;
757   delete fTRDTracker;
758   fTRDTracker = NULL;
759   delete fTOFTracker;
760   fTOFTracker = NULL;
761
762   delete fRunLoader;
763   fRunLoader = NULL;
764
765   if (file) {
766     file->Close();
767     delete file;
768   }
769 }
770
771
772 //_____________________________________________________________________________
773 Bool_t AliReconstruction::ReadESD(AliESD*& esd, const char* recStep) const
774 {
775 // read the ESD event from a file
776
777   if (!esd) return kFALSE;
778   char fileName[256];
779   sprintf(fileName, "ESD_%d.%d_%s.root", 
780           esd->GetRunNumber(), esd->GetEventNumber(), recStep);
781   if (gSystem->AccessPathName(fileName)) return kFALSE;
782
783   Info("ReadESD", "reading ESD from file %s", fileName);
784   TFile* file = TFile::Open(fileName);
785   if (!file || !file->IsOpen()) {
786     Error("ReadESD", "opening %s failed", fileName);
787     delete file;
788     return kFALSE;
789   }
790
791   gROOT->cd();
792   delete esd;
793   esd = (AliESD*) file->Get("ESD");
794   file->Close();
795   delete file;
796   return kTRUE;
797 }
798
799 //_____________________________________________________________________________
800 void AliReconstruction::WriteESD(AliESD* esd, const char* recStep) const
801 {
802 // write the ESD event to a file
803
804   if (!esd) return;
805   char fileName[256];
806   sprintf(fileName, "ESD_%d.%d_%s.root", 
807           esd->GetRunNumber(), esd->GetEventNumber(), recStep);
808
809   Info("WriteESD", "writing ESD to file %s", fileName);
810   TFile* file = TFile::Open(fileName, "recreate");
811   if (!file || !file->IsOpen()) {
812     Error("WriteESD", "opening %s failed", fileName);
813   } else {
814     esd->Write("ESD");
815     file->Close();
816   }
817   delete file;
818 }