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