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