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