]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliReconstruction.cxx
create galice.root, headers and loaders if they don't exist (in case of raw data...
[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   }
474
475   // pass 2: ALL backwards
476   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
477     if (!fTracker[iDet]) continue;
478     AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
479
480     // load clusters
481     if (iDet > 1) {     // all except ITS, TPC
482       TTree* tree = NULL;
483       if (iDet == 3) {   // TOF
484         fLoader[iDet]->LoadDigits("read");
485         tree = fLoader[iDet]->TreeD();
486       } else {
487         fLoader[iDet]->LoadRecPoints("read");
488         tree = fLoader[iDet]->TreeR();
489       }
490       if (!tree) {
491         AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
492         return kFALSE;
493       }
494       fTracker[iDet]->LoadClusters(tree);
495     }
496
497     // run tracking
498     if (fTracker[iDet]->PropagateBack(esd) != 0) {
499       AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
500       return kFALSE;
501     }
502     if (fCheckPointLevel > 1) {
503       WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
504     }
505
506     // unload clusters
507     if (iDet > 2) {     // all except ITS, TPC, TRD
508       fTracker[iDet]->UnloadClusters();
509       if (iDet == 3) {   // TOF
510         fLoader[iDet]->UnloadDigits();
511       } else {
512         fLoader[iDet]->UnloadRecPoints();
513       }
514     }
515   }
516
517   // pass 3: TRD + TPC + ITS refit inwards
518   for (Int_t iDet = 2; iDet >= 0; iDet--) {
519     if (!fTracker[iDet]) continue;
520     AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
521
522     // run tracking
523     if (fTracker[iDet]->RefitInward(esd) != 0) {
524       AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
525       return kFALSE;
526     }
527     if (fCheckPointLevel > 1) {
528       WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
529     }
530
531     // unload clusters
532     fTracker[iDet]->UnloadClusters();
533     fLoader[iDet]->UnloadRecPoints();
534   }
535
536   AliInfo("execution time:");
537   ToAliInfo(stopwatch.Print());
538
539   return kTRUE;
540 }
541
542 //_____________________________________________________________________________
543 Bool_t AliReconstruction::FillESD(AliESD*& esd, const TString& detectors)
544 {
545 // fill the event summary data
546
547   TStopwatch stopwatch;
548   stopwatch.Start();
549   AliInfo("filling ESD");
550
551   TString detStr = detectors;
552   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
553     if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
554     AliReconstructor* reconstructor = GetReconstructor(iDet);
555     if (!reconstructor) continue;
556
557     if (!ReadESD(esd, fgkDetectorName[iDet])) {
558       AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
559       if (fRawReader) {
560         reconstructor->FillESD(fRunLoader, fRawReader, esd);
561       } else {
562         reconstructor->FillESD(fRunLoader, esd);
563       }
564       if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
565     }
566   }
567
568   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
569     AliError(Form("the following detectors were not found: %s", 
570                   detStr.Data()));
571     if (fStopOnError) return kFALSE;
572   }
573
574   AliInfo("execution time:");
575   ToAliInfo(stopwatch.Print());
576
577   return kTRUE;
578 }
579
580
581 //_____________________________________________________________________________
582 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
583 {
584 // check whether detName is contained in detectors
585 // if yes, it is removed from detectors
586
587   // check if all detectors are selected
588   if ((detectors.CompareTo("ALL") == 0) ||
589       detectors.BeginsWith("ALL ") ||
590       detectors.EndsWith(" ALL") ||
591       detectors.Contains(" ALL ")) {
592     detectors = "ALL";
593     return kTRUE;
594   }
595
596   // search for the given detector
597   Bool_t result = kFALSE;
598   if ((detectors.CompareTo(detName) == 0) ||
599       detectors.BeginsWith(detName+" ") ||
600       detectors.EndsWith(" "+detName) ||
601       detectors.Contains(" "+detName+" ")) {
602     detectors.ReplaceAll(detName, "");
603     result = kTRUE;
604   }
605
606   // clean up the detectors string
607   while (detectors.Contains("  ")) detectors.ReplaceAll("  ", " ");
608   while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
609   while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
610
611   return result;
612 }
613
614 //_____________________________________________________________________________
615 Bool_t AliReconstruction::InitRunLoader()
616 {
617 // get or create the run loader
618
619   if (gAlice) delete gAlice;
620   gAlice = NULL;
621
622   if (!gSystem->AccessPathName(fGAliceFileName.Data())) {
623     fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
624     if (!fRunLoader) {
625       AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
626       CleanUp();
627       return kFALSE;
628     }
629     if (fRunLoader->LoadgAlice() == 0) {
630       gAlice = fRunLoader->GetAliRun();
631       AliTracker::SetFieldMap(gAlice->Field());
632     }
633     if (!gAlice && !fRawReader) {
634       AliError(Form("no gAlice object found in file %s",
635                     fGAliceFileName.Data()));
636       CleanUp();
637       return kFALSE;
638     }
639
640   } else {               // galice.root does not exist
641     if (!fRawReader) {
642       AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
643       CleanUp();
644       return kFALSE;
645     }
646     fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
647                                     AliConfig::GetDefaultEventFolderName(),
648                                     "recreate");
649     if (!fRunLoader) {
650       AliError(Form("could not create run loader in file %s", 
651                     fGAliceFileName.Data()));
652       CleanUp();
653       return kFALSE;
654     }
655     fRunLoader->MakeTree("E");
656     Int_t iEvent = 0;
657     while (fRawReader->NextEvent()) {
658       fRunLoader->SetEventNumber(iEvent);
659       fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(), 
660                                      iEvent, iEvent);
661       fRunLoader->MakeTree("H");
662       fRunLoader->TreeE()->Fill();
663       iEvent++;
664     }
665     fRawReader->RewindEvents();
666     fRunLoader->WriteHeader("OVERWRITE");
667     fRunLoader->CdGAFile();
668     fRunLoader->Write(0, TObject::kOverwrite);
669 //    AliTracker::SetFieldMap(???);
670   }
671
672   return kTRUE;
673 }
674
675 //_____________________________________________________________________________
676 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
677 {
678 // get the reconstructor object and the loader for a detector
679
680   if (fReconstructor[iDet]) return fReconstructor[iDet];
681
682   // load the reconstructor object
683   TPluginManager* pluginManager = gROOT->GetPluginManager();
684   TString detName = fgkDetectorName[iDet];
685   TString recName = "Ali" + detName + "Reconstructor";
686   if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT")) return NULL;
687
688   if (detName == "HLT") {
689     if (!gROOT->GetClass("AliLevel3")) {
690       gSystem->Load("libAliL3Src.so");
691       gSystem->Load("libAliL3Misc.so");
692       gSystem->Load("libAliL3Hough.so");
693       gSystem->Load("libAliL3Comp.so");
694     }
695   }
696
697   AliReconstructor* reconstructor = NULL;
698   // first check if a plugin is defined for the reconstructor
699   TPluginHandler* pluginHandler = 
700     pluginManager->FindHandler("AliReconstructor", detName);
701   // if not, add a plugin for it
702   if (!pluginHandler) {
703     AliDebug(1, Form("defining plugin for %s", recName.Data()));
704     if (gSystem->Load("lib" + detName + "base.so") == 0) {
705       pluginManager->AddHandler("AliReconstructor", detName, 
706                                 recName, detName + "rec", recName + "()");
707     } else {
708       pluginManager->AddHandler("AliReconstructor", detName, 
709                                 recName, detName, recName + "()");
710     }
711     pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
712   }
713   if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
714     reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
715   }
716   if (reconstructor) {
717     TObject* obj = fOptions.FindObject(detName.Data());
718     if (obj) reconstructor->SetOption(obj->GetTitle());
719     fReconstructor[iDet] = reconstructor;
720   }
721
722   // get or create the loader
723   if (detName != "HLT") {
724     fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
725     if (!fLoader[iDet]) {
726       AliConfig::Instance()
727         ->CreateDetectorFolders(fRunLoader->GetEventFolder(), 
728                                 detName, detName);
729       // first check if a plugin is defined for the loader
730       TPluginHandler* pluginHandler = 
731         pluginManager->FindHandler("AliLoader", detName);
732       // if not, add a plugin for it
733       if (!pluginHandler) {
734         TString loaderName = "Ali" + detName + "Loader";
735         AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
736         pluginManager->AddHandler("AliLoader", detName, 
737                                   loaderName, detName + "base", 
738                                   loaderName + "(const char*, TFolder*)");
739         pluginHandler = pluginManager->FindHandler("AliLoader", detName);
740       }
741       if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
742         fLoader[iDet] = 
743           (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(), 
744                                                  fRunLoader->GetEventFolder());
745       }
746       if (!fLoader[iDet]) {   // use default loader
747         fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
748       }
749       if (!fLoader[iDet]) {
750         AliWarning(Form("couldn't get loader for %s", detName.Data()));
751         if (fStopOnError) return kFALSE;
752       } else {
753         fRunLoader->AddLoader(fLoader[iDet]);
754         fRunLoader->CdGAFile();
755         if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
756         fRunLoader->Write(0, TObject::kOverwrite);
757       }
758     }
759   }
760       
761   return reconstructor;
762 }
763
764 //_____________________________________________________________________________
765 Bool_t AliReconstruction::CreateVertexer()
766 {
767 // create the vertexer
768
769   fVertexer = NULL;
770   AliReconstructor* itsReconstructor = GetReconstructor(0);
771   if (itsReconstructor) {
772     fVertexer = itsReconstructor->CreateVertexer(fRunLoader);
773   }
774   if (!fVertexer) {
775     AliWarning("couldn't create a vertexer for ITS");
776     if (fStopOnError) return kFALSE;
777   }
778
779   return kTRUE;
780 }
781
782 //_____________________________________________________________________________
783 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
784 {
785 // create the trackers
786
787   TString detStr = detectors;
788   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
789     if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
790     AliReconstructor* reconstructor = GetReconstructor(iDet);
791     if (!reconstructor) continue;
792     TString detName = fgkDetectorName[iDet];
793     if (detName == "HLT") continue;
794
795     fTracker[iDet] = reconstructor->CreateTracker(fRunLoader);
796     if (!fTracker[iDet] && (iDet < 7)) {
797       AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
798       if (fStopOnError) return kFALSE;
799     }
800   }
801
802   return kTRUE;
803 }
804
805 //_____________________________________________________________________________
806 void AliReconstruction::CleanUp(TFile* file)
807 {
808 // delete trackers and the run loader and close and delete the file
809
810   for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
811     delete fReconstructor[iDet];
812     fReconstructor[iDet] = NULL;
813     fLoader[iDet] = NULL;
814     delete fTracker[iDet];
815     fTracker[iDet] = NULL;
816   }
817   delete fVertexer;
818   fVertexer = NULL;
819
820   delete fRunLoader;
821   fRunLoader = NULL;
822   delete fRawReader;
823   fRawReader = NULL;
824
825   if (file) {
826     file->Close();
827     delete file;
828   }
829 }
830
831
832 //_____________________________________________________________________________
833 Bool_t AliReconstruction::ReadESD(AliESD*& esd, const char* recStep) const
834 {
835 // read the ESD event from a file
836
837   if (!esd) return kFALSE;
838   char fileName[256];
839   sprintf(fileName, "ESD_%d.%d_%s.root", 
840           esd->GetRunNumber(), esd->GetEventNumber(), recStep);
841   if (gSystem->AccessPathName(fileName)) return kFALSE;
842
843   AliDebug(1, Form("reading ESD from file %s", fileName));
844   TFile* file = TFile::Open(fileName);
845   if (!file || !file->IsOpen()) {
846     AliError(Form("opening %s failed", fileName));
847     delete file;
848     return kFALSE;
849   }
850
851   gROOT->cd();
852   delete esd;
853   esd = (AliESD*) file->Get("ESD");
854   file->Close();
855   delete file;
856   return kTRUE;
857 }
858
859 //_____________________________________________________________________________
860 void AliReconstruction::WriteESD(AliESD* esd, const char* recStep) const
861 {
862 // write the ESD event to a file
863
864   if (!esd) return;
865   char fileName[256];
866   sprintf(fileName, "ESD_%d.%d_%s.root", 
867           esd->GetRunNumber(), esd->GetEventNumber(), recStep);
868
869   AliDebug(1, Form("writing ESD to file %s", fileName));
870   TFile* file = TFile::Open(fileName, "recreate");
871   if (!file || !file->IsOpen()) {
872     AliError(Form("opening %s failed", fileName));
873   } else {
874     esd->Write("ESD");
875     file->Close();
876   }
877   delete file;
878 }