]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliReconstruction.cxx
d36e33b8a7c241046eb1a1040ddcae9ab04d9e16
[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   // get the run loader
242   if (!InitRunLoader()) return kFALSE;
243
244   // local reconstruction
245   if (!fRunLocalReconstruction.IsNull()) {
246     if (!RunLocalReconstruction(fRunLocalReconstruction)) {
247       if (fStopOnError) {CleanUp(); return kFALSE;}
248     }
249   }
250   if (!fRunVertexFinder && fRunTracking.IsNull() && 
251       fFillESD.IsNull()) return kTRUE;
252
253   // get vertexer
254   if (fRunVertexFinder && !CreateVertexer()) {
255     if (fStopOnError) {
256       CleanUp(); 
257       return kFALSE;
258     }
259   }
260
261   // get trackers
262   if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
263     if (fStopOnError) {
264       CleanUp(); 
265       return kFALSE;
266     }      
267   }
268
269   // create the ESD output file and tree
270   TFile* file = TFile::Open("AliESDs.root", "RECREATE");
271   if (!file->IsOpen()) {
272     AliError("opening AliESDs.root failed");
273     if (fStopOnError) {CleanUp(file); return kFALSE;}    
274   }
275   AliESD* esd = new AliESD;
276   TTree* tree = new TTree("esdTree", "Tree with ESD objects");
277   tree->Branch("ESD", "AliESD", &esd);
278   delete esd;
279   gROOT->cd();
280
281   // loop over events
282   if (fRawReader) fRawReader->RewindEvents();
283   
284   for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) {
285     AliInfo(Form("processing event %d", iEvent));
286     fRunLoader->GetEvent(iEvent);
287     if (fRawReader) fRawReader->NextEvent();
288
289     char fileName[256];
290     sprintf(fileName, "ESD_%d.%d_final.root", 
291             fRunLoader->GetHeader()->GetRun(), 
292             fRunLoader->GetHeader()->GetEventNrInRun());
293     if (!gSystem->AccessPathName(fileName)) continue;
294
295     esd = new AliESD;
296     esd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
297     esd->SetEventNumber(fRunLoader->GetHeader()->GetEventNrInRun());
298     if (gAlice) {
299       esd->SetMagneticField(gAlice->Field()->SolenoidField());
300     } else {
301       // ???
302     }
303
304     // vertex finder
305     if (fRunVertexFinder) {
306       if (!ReadESD(esd, "vertex")) {
307         if (!RunVertexFinder(esd)) {
308           if (fStopOnError) {CleanUp(file); return kFALSE;}
309         }
310         if (fCheckPointLevel > 0) WriteESD(esd, "vertex");
311       }
312     }
313
314     // barrel tracking
315     if (!fRunTracking.IsNull()) {
316       if (!ReadESD(esd, "tracking")) {
317         if (!RunTracking(esd)) {
318           if (fStopOnError) {CleanUp(file); return kFALSE;}
319         }
320         if (fCheckPointLevel > 0) WriteESD(esd, "tracking");
321       }
322     }
323
324     // fill ESD
325     if (!fFillESD.IsNull()) {
326       if (!FillESD(esd, fFillESD)) {
327         if (fStopOnError) {CleanUp(file); return kFALSE;}
328       }
329     }
330
331     // combined PID
332     AliESDpid::MakePID(esd);
333     if (fCheckPointLevel > 1) WriteESD(esd, "PID");
334
335     // write ESD
336     tree->Fill();
337
338     if (fCheckPointLevel > 0) WriteESD(esd, "final");
339     delete esd;
340   }
341
342   file->cd();
343   tree->Write();
344   CleanUp(file);
345
346   return kTRUE;
347 }
348
349
350 //_____________________________________________________________________________
351 Bool_t AliReconstruction::RunLocalReconstruction(const TString& detectors)
352 {
353 // run the local reconstruction
354
355   TStopwatch stopwatch;
356   stopwatch.Start();
357
358   TString detStr = detectors;
359   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
360     if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
361     AliReconstructor* reconstructor = GetReconstructor(iDet);
362     if (!reconstructor) continue;
363
364     AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
365     TStopwatch stopwatchDet;
366     stopwatchDet.Start();
367     if (fRawReader) {
368       fRawReader->RewindEvents();
369       reconstructor->Reconstruct(fRunLoader, fRawReader);
370     } else {
371       reconstructor->Reconstruct(fRunLoader);
372     }
373     AliInfo(Form("execution time for %s:", fgkDetectorName[iDet]));
374     ToAliInfo(stopwatchDet.Print());
375   }
376
377   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
378     AliError(Form("the following detectors were not found: %s",
379                   detStr.Data()));
380     if (fStopOnError) return kFALSE;
381   }
382
383   AliInfo("execution time:");
384   ToAliInfo(stopwatch.Print());
385
386   return kTRUE;
387 }
388
389 //_____________________________________________________________________________
390 Bool_t AliReconstruction::RunVertexFinder(AliESD*& esd)
391 {
392 // run the barrel tracking
393
394   TStopwatch stopwatch;
395   stopwatch.Start();
396
397   AliESDVertex* vertex = NULL;
398   Double_t vtxPos[3] = {0, 0, 0};
399   Double_t vtxErr[3] = {0.07, 0.07, 0.1};
400   TArrayF mcVertex(3); 
401   if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
402     fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
403     for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
404   }
405
406   if (fVertexer) {
407     AliInfo("running the ITS vertex finder");
408     vertex = fVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
409     if(!vertex){
410       AliWarning("Vertex not found");
411       vertex = new AliESDVertex();
412     }
413     else {
414       vertex->SetTruePos(vtxPos);  // store also the vertex from MC
415     }
416
417   } else {
418     AliInfo("getting the primary vertex from MC");
419     vertex = new AliESDVertex(vtxPos, vtxErr);
420   }
421
422   if (vertex) {
423     vertex->GetXYZ(vtxPos);
424     vertex->GetSigmaXYZ(vtxErr);
425   } else {
426     AliWarning("no vertex reconstructed");
427     vertex = new AliESDVertex(vtxPos, vtxErr);
428   }
429   esd->SetVertex(vertex);
430   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
431     if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
432   }  
433   delete vertex;
434
435   AliInfo("execution time:");
436   ToAliInfo(stopwatch.Print());
437
438   return kTRUE;
439 }
440
441 //_____________________________________________________________________________
442 Bool_t AliReconstruction::RunTracking(AliESD*& esd)
443 {
444 // run the barrel tracking
445
446   TStopwatch stopwatch;
447   stopwatch.Start();
448
449   AliInfo("running tracking");
450
451   // pass 1: TPC + ITS inwards
452   for (Int_t iDet = 1; iDet >= 0; iDet--) {
453     if (!fTracker[iDet]) continue;
454     AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
455
456     // load clusters
457     fLoader[iDet]->LoadRecPoints("read");
458     TTree* tree = fLoader[iDet]->TreeR();
459     if (!tree) {
460       AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
461       return kFALSE;
462     }
463     fTracker[iDet]->LoadClusters(tree);
464
465     // run tracking
466     if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
467       AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
468       return kFALSE;
469     }
470     if (fCheckPointLevel > 1) {
471       WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
472     }
473     // preliminary PID in TPC needed by the ITS tracker
474     if (iDet == 1) {
475       GetReconstructor(1)->FillESD(fRunLoader, esd);
476       AliESDpid::MakePID(esd);
477     }
478   }
479
480   // pass 2: ALL backwards
481   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
482     if (!fTracker[iDet]) continue;
483     AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
484
485     // load clusters
486     if (iDet > 1) {     // all except ITS, TPC
487       TTree* tree = NULL;
488       if (iDet == 3) {   // TOF
489         fLoader[iDet]->LoadDigits("read");
490         tree = fLoader[iDet]->TreeD();
491       } else {
492         fLoader[iDet]->LoadRecPoints("read");
493         tree = fLoader[iDet]->TreeR();
494       }
495       if (!tree) {
496         AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
497         return kFALSE;
498       }
499       fTracker[iDet]->LoadClusters(tree);
500     }
501
502     // run tracking
503     if (fTracker[iDet]->PropagateBack(esd) != 0) {
504       AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
505       return kFALSE;
506     }
507     if (fCheckPointLevel > 1) {
508       WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
509     }
510
511     // unload clusters
512     if (iDet > 2) {     // all except ITS, TPC, TRD
513       fTracker[iDet]->UnloadClusters();
514       if (iDet == 3) {   // TOF
515         fLoader[iDet]->UnloadDigits();
516       } else {
517         fLoader[iDet]->UnloadRecPoints();
518       }
519     }
520   }
521
522   // pass 3: TRD + TPC + ITS refit inwards
523   for (Int_t iDet = 2; iDet >= 0; iDet--) {
524     if (!fTracker[iDet]) continue;
525     AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
526
527     // run tracking
528     if (fTracker[iDet]->RefitInward(esd) != 0) {
529       AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
530       return kFALSE;
531     }
532     if (fCheckPointLevel > 1) {
533       WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
534     }
535
536     // unload clusters
537     fTracker[iDet]->UnloadClusters();
538     fLoader[iDet]->UnloadRecPoints();
539   }
540
541   AliInfo("execution time:");
542   ToAliInfo(stopwatch.Print());
543
544   return kTRUE;
545 }
546
547 //_____________________________________________________________________________
548 Bool_t AliReconstruction::FillESD(AliESD*& esd, const TString& detectors)
549 {
550 // fill the event summary data
551
552   TStopwatch stopwatch;
553   stopwatch.Start();
554   AliInfo("filling ESD");
555
556   TString detStr = detectors;
557   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
558     if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
559     AliReconstructor* reconstructor = GetReconstructor(iDet);
560     if (!reconstructor) continue;
561
562     if (!ReadESD(esd, fgkDetectorName[iDet])) {
563       AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
564       if (fRawReader) {
565         reconstructor->FillESD(fRunLoader, fRawReader, esd);
566       } else {
567         reconstructor->FillESD(fRunLoader, esd);
568       }
569       if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
570     }
571   }
572
573   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
574     AliError(Form("the following detectors were not found: %s", 
575                   detStr.Data()));
576     if (fStopOnError) return kFALSE;
577   }
578
579   AliInfo("execution time:");
580   ToAliInfo(stopwatch.Print());
581
582   return kTRUE;
583 }
584
585
586 //_____________________________________________________________________________
587 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
588 {
589 // check whether detName is contained in detectors
590 // if yes, it is removed from detectors
591
592   // check if all detectors are selected
593   if ((detectors.CompareTo("ALL") == 0) ||
594       detectors.BeginsWith("ALL ") ||
595       detectors.EndsWith(" ALL") ||
596       detectors.Contains(" ALL ")) {
597     detectors = "ALL";
598     return kTRUE;
599   }
600
601   // search for the given detector
602   Bool_t result = kFALSE;
603   if ((detectors.CompareTo(detName) == 0) ||
604       detectors.BeginsWith(detName+" ") ||
605       detectors.EndsWith(" "+detName) ||
606       detectors.Contains(" "+detName+" ")) {
607     detectors.ReplaceAll(detName, "");
608     result = kTRUE;
609   }
610
611   // clean up the detectors string
612   while (detectors.Contains("  ")) detectors.ReplaceAll("  ", " ");
613   while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
614   while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
615
616   return result;
617 }
618
619 //_____________________________________________________________________________
620 Bool_t AliReconstruction::InitRunLoader()
621 {
622 // get or create the run loader
623
624   if (gAlice) delete gAlice;
625   gAlice = NULL;
626
627   if (!gSystem->AccessPathName(fGAliceFileName.Data())) {
628     fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
629     if (!fRunLoader) {
630       AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
631       CleanUp();
632       return kFALSE;
633     }
634     if (fRunLoader->LoadgAlice() == 0) {
635       gAlice = fRunLoader->GetAliRun();
636       AliTracker::SetFieldMap(gAlice->Field());
637     }
638     if (!gAlice && !fRawReader) {
639       AliError(Form("no gAlice object found in file %s",
640                     fGAliceFileName.Data()));
641       CleanUp();
642       return kFALSE;
643     }
644
645   } else {               // galice.root does not exist
646     if (!fRawReader) {
647       AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
648       CleanUp();
649       return kFALSE;
650     }
651     fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
652                                     AliConfig::GetDefaultEventFolderName(),
653                                     "recreate");
654     if (!fRunLoader) {
655       AliError(Form("could not create run loader in file %s", 
656                     fGAliceFileName.Data()));
657       CleanUp();
658       return kFALSE;
659     }
660     fRunLoader->MakeTree("E");
661     Int_t iEvent = 0;
662     while (fRawReader->NextEvent()) {
663       fRunLoader->SetEventNumber(iEvent);
664       fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(), 
665                                      iEvent, iEvent);
666       fRunLoader->MakeTree("H");
667       fRunLoader->TreeE()->Fill();
668       iEvent++;
669     }
670     fRawReader->RewindEvents();
671     fRunLoader->WriteHeader("OVERWRITE");
672     fRunLoader->CdGAFile();
673     fRunLoader->Write(0, TObject::kOverwrite);
674 //    AliTracker::SetFieldMap(???);
675   }
676
677   return kTRUE;
678 }
679
680 //_____________________________________________________________________________
681 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
682 {
683 // get the reconstructor object and the loader for a detector
684
685   if (fReconstructor[iDet]) return fReconstructor[iDet];
686
687   // load the reconstructor object
688   TPluginManager* pluginManager = gROOT->GetPluginManager();
689   TString detName = fgkDetectorName[iDet];
690   TString recName = "Ali" + detName + "Reconstructor";
691   if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT")) return NULL;
692
693   if (detName == "HLT") {
694     if (!gROOT->GetClass("AliLevel3")) {
695       gSystem->Load("libAliL3Src.so");
696       gSystem->Load("libAliL3Misc.so");
697       gSystem->Load("libAliL3Hough.so");
698       gSystem->Load("libAliL3Comp.so");
699     }
700   }
701
702   AliReconstructor* reconstructor = NULL;
703   // first check if a plugin is defined for the reconstructor
704   TPluginHandler* pluginHandler = 
705     pluginManager->FindHandler("AliReconstructor", detName);
706   // if not, add a plugin for it
707   if (!pluginHandler) {
708     AliDebug(1, Form("defining plugin for %s", recName.Data()));
709     if (gSystem->Load("lib" + detName + "base.so") == 0) {
710       pluginManager->AddHandler("AliReconstructor", detName, 
711                                 recName, detName + "rec", recName + "()");
712     } else {
713       pluginManager->AddHandler("AliReconstructor", detName, 
714                                 recName, detName, recName + "()");
715     }
716     pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
717   }
718   if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
719     reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
720   }
721   if (reconstructor) {
722     TObject* obj = fOptions.FindObject(detName.Data());
723     if (obj) reconstructor->SetOption(obj->GetTitle());
724     fReconstructor[iDet] = reconstructor;
725   }
726
727   // get or create the loader
728   if (detName != "HLT") {
729     fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
730     if (!fLoader[iDet]) {
731       AliConfig::Instance()
732         ->CreateDetectorFolders(fRunLoader->GetEventFolder(), 
733                                 detName, detName);
734       // first check if a plugin is defined for the loader
735       TPluginHandler* pluginHandler = 
736         pluginManager->FindHandler("AliLoader", detName);
737       // if not, add a plugin for it
738       if (!pluginHandler) {
739         TString loaderName = "Ali" + detName + "Loader";
740         AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
741         pluginManager->AddHandler("AliLoader", detName, 
742                                   loaderName, detName + "base", 
743                                   loaderName + "(const char*, TFolder*)");
744         pluginHandler = pluginManager->FindHandler("AliLoader", detName);
745       }
746       if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
747         fLoader[iDet] = 
748           (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(), 
749                                                  fRunLoader->GetEventFolder());
750       }
751       if (!fLoader[iDet]) {   // use default loader
752         fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
753       }
754       if (!fLoader[iDet]) {
755         AliWarning(Form("couldn't get loader for %s", detName.Data()));
756         if (fStopOnError) return NULL;
757       } else {
758         fRunLoader->AddLoader(fLoader[iDet]);
759         fRunLoader->CdGAFile();
760         if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
761         fRunLoader->Write(0, TObject::kOverwrite);
762       }
763     }
764   }
765       
766   return reconstructor;
767 }
768
769 //_____________________________________________________________________________
770 Bool_t AliReconstruction::CreateVertexer()
771 {
772 // create the vertexer
773
774   fVertexer = NULL;
775   AliReconstructor* itsReconstructor = GetReconstructor(0);
776   if (itsReconstructor) {
777     fVertexer = itsReconstructor->CreateVertexer(fRunLoader);
778   }
779   if (!fVertexer) {
780     AliWarning("couldn't create a vertexer for ITS");
781     if (fStopOnError) return kFALSE;
782   }
783
784   return kTRUE;
785 }
786
787 //_____________________________________________________________________________
788 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
789 {
790 // create the trackers
791
792   TString detStr = detectors;
793   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
794     if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
795     AliReconstructor* reconstructor = GetReconstructor(iDet);
796     if (!reconstructor) continue;
797     TString detName = fgkDetectorName[iDet];
798     if (detName == "HLT") continue;
799
800     fTracker[iDet] = reconstructor->CreateTracker(fRunLoader);
801     if (!fTracker[iDet] && (iDet < 7)) {
802       AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
803       if (fStopOnError) return kFALSE;
804     }
805   }
806
807   return kTRUE;
808 }
809
810 //_____________________________________________________________________________
811 void AliReconstruction::CleanUp(TFile* file)
812 {
813 // delete trackers and the run loader and close and delete the file
814
815   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
816     delete fReconstructor[iDet];
817     fReconstructor[iDet] = NULL;
818     fLoader[iDet] = NULL;
819     delete fTracker[iDet];
820     fTracker[iDet] = NULL;
821   }
822   delete fVertexer;
823   fVertexer = NULL;
824
825   delete fRunLoader;
826   fRunLoader = NULL;
827   delete fRawReader;
828   fRawReader = NULL;
829
830   if (file) {
831     file->Close();
832     delete file;
833   }
834 }
835
836
837 //_____________________________________________________________________________
838 Bool_t AliReconstruction::ReadESD(AliESD*& esd, const char* recStep) const
839 {
840 // read the ESD event from a file
841
842   if (!esd) return kFALSE;
843   char fileName[256];
844   sprintf(fileName, "ESD_%d.%d_%s.root", 
845           esd->GetRunNumber(), esd->GetEventNumber(), recStep);
846   if (gSystem->AccessPathName(fileName)) return kFALSE;
847
848   AliDebug(1, Form("reading ESD from file %s", fileName));
849   TFile* file = TFile::Open(fileName);
850   if (!file || !file->IsOpen()) {
851     AliError(Form("opening %s failed", fileName));
852     delete file;
853     return kFALSE;
854   }
855
856   gROOT->cd();
857   delete esd;
858   esd = (AliESD*) file->Get("ESD");
859   file->Close();
860   delete file;
861   return kTRUE;
862 }
863
864 //_____________________________________________________________________________
865 void AliReconstruction::WriteESD(AliESD* esd, const char* recStep) const
866 {
867 // write the ESD event to a file
868
869   if (!esd) return;
870   char fileName[256];
871   sprintf(fileName, "ESD_%d.%d_%s.root", 
872           esd->GetRunNumber(), esd->GetEventNumber(), recStep);
873
874   AliDebug(1, Form("writing ESD to file %s", fileName));
875   TFile* file = TFile::Open(fileName, "recreate");
876   if (!file || !file->IsOpen()) {
877     AliError(Form("opening %s failed", fileName));
878   } else {
879     esd->Write("ESD");
880     file->Close();
881   }
882   delete file;
883 }