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