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