1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // class for running the reconstruction //
22 // Clusters and tracks are created for all detectors and all events by //
25 // AliReconstruction rec; //
28 // The Run method returns kTRUE in case of successful execution. //
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 //
33 // rec.SetInput("..."); //
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 //
41 // By default all events are reconstructed. The reconstruction can be //
42 // limited to a range of events by giving the index of the first and the //
43 // last event as an argument to the Run method or by calling //
45 // rec.SetEventRange(..., ...); //
47 // The index -1 (default) can be used for the last event to indicate no //
48 // upper limit of the event range. //
50 // The name of the galice file can be changed from the default //
51 // "galice.root" by passing it as argument to the AliReconstruction //
52 // constructor or by //
54 // rec.SetGAliceFile("..."); //
56 // The local reconstruction can be switched on or off for individual //
59 // rec.SetRunLocalReconstruction("..."); //
61 // The argument is a (case sensitive) string with the names of the //
62 // detectors separated by a space. The special string "ALL" selects all //
63 // available detectors. This is the default. //
65 // The reconstruction of the primary vertex position can be switched off by //
67 // rec.SetRunVertexFinder(kFALSE); //
69 // The tracking and the creation of ESD tracks can be switched on for //
70 // selected detectors by //
72 // rec.SetRunTracking("..."); //
74 // The filling of additional ESD information can be steered by //
76 // rec.SetFillESD("..."); //
78 // Again, for both methods the string specifies the list of detectors. //
79 // The default is "ALL". //
81 // The call of the shortcut method //
83 // rec.SetRunReconstruction("..."); //
85 // is equivalent to calling SetRunLocalReconstruction, SetRunTracking and //
86 // SetFillESD with the same detector selecting string as argument. //
88 // The reconstruction requires digits or raw data as input. For the creation //
89 // of digits and raw data have a look at the class AliSimulation. //
91 // For debug purposes the method SetCheckPointLevel can be used. If the //
92 // argument is greater than 0, files with ESD events will be written after //
93 // selected steps of the reconstruction for each event: //
94 // level 1: after tracking and after filling of ESD (final) //
95 // level 2: in addition after each tracking step //
96 // level 3: in addition after the filling of ESD for each detector //
97 // If a final check point file exists for an event, this event will be //
98 // skipped in the reconstruction. The tracking and the filling of ESD for //
99 // a detector will be skipped as well, if the corresponding check point //
100 // file exists. The ESD event will then be loaded from the file instead. //
102 ///////////////////////////////////////////////////////////////////////////////
108 #include <TPluginManager.h>
109 #include <TStopwatch.h>
111 #include "AliReconstruction.h"
112 #include "AliReconstructor.h"
114 #include "AliRunLoader.h"
116 #include "AliRawReaderFile.h"
117 #include "AliRawReaderDate.h"
118 #include "AliRawReaderRoot.h"
119 #include "AliTracker.h"
121 #include "AliESDVertex.h"
122 #include "AliVertexer.h"
123 #include "AliHeader.h"
124 #include "AliGenEventHeader.h"
126 #include "AliESDpid.h"
129 ClassImp(AliReconstruction)
132 //_____________________________________________________________________________
133 const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "RICH", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "START", "VZERO", "CRT", "HLT"};
135 //_____________________________________________________________________________
136 AliReconstruction::AliReconstruction(const char* gAliceFilename,
137 const char* name, const char* title) :
140 fRunLocalReconstruction("ALL"),
141 fRunVertexFinder(kTRUE),
142 fRunHLTTracking(kFALSE),
145 fGAliceFileName(gAliceFilename),
149 fStopOnError(kFALSE),
158 // create reconstruction object with default parameters
160 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
161 fReconstructor[iDet] = NULL;
162 fLoader[iDet] = NULL;
163 fTracker[iDet] = NULL;
168 //_____________________________________________________________________________
169 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
172 fRunLocalReconstruction(rec.fRunLocalReconstruction),
173 fRunVertexFinder(rec.fRunVertexFinder),
174 fRunHLTTracking(rec.fRunHLTTracking),
175 fRunTracking(rec.fRunTracking),
176 fFillESD(rec.fFillESD),
177 fGAliceFileName(rec.fGAliceFileName),
179 fFirstEvent(rec.fFirstEvent),
180 fLastEvent(rec.fLastEvent),
181 fStopOnError(rec.fStopOnError),
192 for (Int_t i = 0; i < fOptions.GetEntriesFast(); i++) {
193 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
195 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
196 fReconstructor[iDet] = NULL;
197 fLoader[iDet] = NULL;
198 fTracker[iDet] = NULL;
202 //_____________________________________________________________________________
203 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
205 // assignment operator
207 this->~AliReconstruction();
208 new(this) AliReconstruction(rec);
212 //_____________________________________________________________________________
213 AliReconstruction::~AliReconstruction()
222 //_____________________________________________________________________________
223 void AliReconstruction::SetGAliceFile(const char* fileName)
225 // set the name of the galice file
227 fGAliceFileName = fileName;
230 //_____________________________________________________________________________
231 void AliReconstruction::SetOption(const char* detector, const char* option)
233 // set options for the reconstruction of a detector
235 TObject* obj = fOptions.FindObject(detector);
236 if (obj) fOptions.Remove(obj);
237 fOptions.Add(new TNamed(detector, option));
241 //_____________________________________________________________________________
242 Bool_t AliReconstruction::Run(const char* input,
243 Int_t firstEvent, Int_t lastEvent)
245 // run the reconstruction
248 if (!input) input = fInput.Data();
249 TString fileName(input);
250 if (fileName.EndsWith("/")) {
251 fRawReader = new AliRawReaderFile(fileName);
252 } else if (fileName.EndsWith(".root")) {
253 fRawReader = new AliRawReaderRoot(fileName);
254 } else if (!fileName.IsNull()) {
255 fRawReader = new AliRawReaderDate(fileName);
256 fRawReader->SelectEvents(7);
259 // get the run loader
260 if (!InitRunLoader()) return kFALSE;
262 // local reconstruction
263 if (!fRunLocalReconstruction.IsNull()) {
264 if (!RunLocalReconstruction(fRunLocalReconstruction)) {
265 if (fStopOnError) {CleanUp(); return kFALSE;}
268 // if (!fRunVertexFinder && fRunTracking.IsNull() &&
269 // fFillESD.IsNull()) return kTRUE;
272 if (fRunVertexFinder && !CreateVertexer()) {
280 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
287 // get the possibly already existing ESD file and tree
288 AliESD* esd = new AliESD; AliESD* hltesd = new AliESD;
289 TFile* fileOld = NULL;
290 TTree* treeOld = NULL; TTree *hlttreeOld = NULL;
291 if (!gSystem->AccessPathName("AliESDs.root")){
292 gSystem->CopyFile("AliESDs.root", "AliESDs.old.root", kTRUE);
293 fileOld = TFile::Open("AliESDs.old.root");
294 if (fileOld && fileOld->IsOpen()) {
295 treeOld = (TTree*) fileOld->Get("esdTree");
296 if (treeOld) treeOld->SetBranchAddress("ESD", &esd);
297 hlttreeOld = (TTree*) fileOld->Get("HLTesdTree");
298 if (hlttreeOld) hlttreeOld->SetBranchAddress("ESD", &hltesd);
302 // create the ESD output file and tree
303 TFile* file = TFile::Open("AliESDs.root", "RECREATE");
304 if (!file->IsOpen()) {
305 AliError("opening AliESDs.root failed");
306 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
308 TTree* tree = new TTree("esdTree", "Tree with ESD objects");
309 tree->Branch("ESD", "AliESD", &esd);
310 TTree* hlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
311 hlttree->Branch("ESD", "AliESD", &hltesd);
312 delete esd; delete hltesd;
313 esd = NULL; hltesd = NULL;
317 if (fRawReader) fRawReader->RewindEvents();
319 for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) {
320 if (fRawReader) fRawReader->NextEvent();
321 if ((iEvent < firstEvent) || ((lastEvent >= 0) && (iEvent > lastEvent))) {
322 // copy old ESD to the new one
324 treeOld->SetBranchAddress("ESD", &esd);
325 treeOld->GetEntry(iEvent);
329 hlttreeOld->SetBranchAddress("ESD", &hltesd);
330 hlttreeOld->GetEntry(iEvent);
336 AliInfo(Form("processing event %d", iEvent));
337 fRunLoader->GetEvent(iEvent);
340 sprintf(fileName, "ESD_%d.%d_final.root",
341 fRunLoader->GetHeader()->GetRun(),
342 fRunLoader->GetHeader()->GetEventNrInRun());
343 if (!gSystem->AccessPathName(fileName)) continue;
345 // local reconstruction
346 if (!fRunLocalReconstruction.IsNull()) {
347 if (!RunLocalEventReconstruction(fRunLocalReconstruction)) {
348 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
352 esd = new AliESD; hltesd = new AliESD;
353 esd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
354 hltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
355 esd->SetEventNumber(fRunLoader->GetHeader()->GetEventNrInRun());
356 hltesd->SetEventNumber(fRunLoader->GetHeader()->GetEventNrInRun());
358 esd->SetMagneticField(gAlice->Field()->SolenoidField());
359 hltesd->SetMagneticField(gAlice->Field()->SolenoidField());
365 if (fRunVertexFinder) {
366 if (!ReadESD(esd, "vertex")) {
367 if (!RunVertexFinder(esd)) {
368 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
370 if (fCheckPointLevel > 0) WriteESD(esd, "vertex");
375 if (!fRunTracking.IsNull()) {
376 if (fRunHLTTracking) {
377 hltesd->SetVertex(esd->GetVertex());
378 if (!RunHLTTracking(hltesd)) {
379 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
385 if (!fRunTracking.IsNull()) {
386 if (!ReadESD(esd, "tracking")) {
387 if (!RunTracking(esd)) {
388 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
390 if (fCheckPointLevel > 0) WriteESD(esd, "tracking");
395 if (!fFillESD.IsNull()) {
396 if (!FillESD(esd, fFillESD)) {
397 if (fStopOnError) {CleanUp(file, fileOld); return kFALSE;}
402 AliESDpid::MakePID(esd);
403 if (fCheckPointLevel > 1) WriteESD(esd, "PID");
410 if (fCheckPointLevel > 0) WriteESD(esd, "final");
411 delete esd; delete hltesd;
412 esd = NULL; hltesd = NULL;
418 CleanUp(file, fileOld);
424 //_____________________________________________________________________________
425 Bool_t AliReconstruction::RunLocalReconstruction(const TString& detectors)
427 // run the local reconstruction
429 TStopwatch stopwatch;
432 TString detStr = detectors;
433 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
434 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
435 AliReconstructor* reconstructor = GetReconstructor(iDet);
436 if (!reconstructor) continue;
437 if (reconstructor->HasLocalReconstruction()) continue;
439 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
440 TStopwatch stopwatchDet;
441 stopwatchDet.Start();
443 fRawReader->RewindEvents();
444 reconstructor->Reconstruct(fRunLoader, fRawReader);
446 reconstructor->Reconstruct(fRunLoader);
448 AliInfo(Form("execution time for %s:", fgkDetectorName[iDet]));
449 ToAliInfo(stopwatchDet.Print());
452 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
453 AliError(Form("the following detectors were not found: %s",
455 if (fStopOnError) return kFALSE;
458 AliInfo("execution time:");
459 ToAliInfo(stopwatch.Print());
464 //_____________________________________________________________________________
465 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
467 // run the local reconstruction
469 TStopwatch stopwatch;
472 TString detStr = detectors;
473 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
474 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
475 AliReconstructor* reconstructor = GetReconstructor(iDet);
476 if (!reconstructor) continue;
477 AliLoader* loader = fLoader[iDet];
479 // conversion of digits
480 if (fRawReader && reconstructor->HasDigitConversion()) {
481 AliInfo(Form("converting raw data digits into root objects for %s",
482 fgkDetectorName[iDet]));
483 TStopwatch stopwatchDet;
484 stopwatchDet.Start();
485 loader->LoadDigits("update");
486 loader->CleanDigits();
487 loader->MakeDigitsContainer();
488 TTree* digitsTree = loader->TreeD();
489 reconstructor->ConvertDigits(fRawReader, digitsTree);
490 loader->WriteDigits("OVERWRITE");
491 loader->UnloadDigits();
492 AliDebug(1, Form("execution time for %s:", fgkDetectorName[iDet]));
493 ToAliDebug(1, stopwatchDet.Print());
496 // local reconstruction
497 if (!reconstructor->HasLocalReconstruction()) continue;
498 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
499 TStopwatch stopwatchDet;
500 stopwatchDet.Start();
501 loader->LoadRecPoints("update");
502 loader->CleanRecPoints();
503 loader->MakeRecPointsContainer();
504 TTree* clustersTree = loader->TreeR();
505 if (fRawReader && !reconstructor->HasDigitConversion()) {
506 reconstructor->Reconstruct(fRawReader, clustersTree);
508 loader->LoadDigits("read");
509 TTree* digitsTree = loader->TreeD();
511 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
512 if (fStopOnError) return kFALSE;
514 reconstructor->Reconstruct(digitsTree, clustersTree);
516 loader->UnloadDigits();
518 loader->WriteRecPoints("OVERWRITE");
519 loader->UnloadRecPoints();
520 AliDebug(1, Form("execution time for %s:", fgkDetectorName[iDet]));
521 ToAliDebug(1, stopwatchDet.Print());
524 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
525 AliError(Form("the following detectors were not found: %s",
527 if (fStopOnError) return kFALSE;
530 AliInfo("execution time:");
531 ToAliInfo(stopwatch.Print());
536 //_____________________________________________________________________________
537 Bool_t AliReconstruction::RunVertexFinder(AliESD*& esd)
539 // run the barrel tracking
541 TStopwatch stopwatch;
544 AliESDVertex* vertex = NULL;
545 Double_t vtxPos[3] = {0, 0, 0};
546 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
548 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
549 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
550 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
554 AliInfo("running the ITS vertex finder");
555 if (fLoader[0]) fLoader[0]->LoadRecPoints();
556 vertex = fVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
557 if (fLoader[0]) fLoader[0]->UnloadRecPoints();
559 AliWarning("Vertex not found");
560 vertex = new AliESDVertex();
563 vertex->SetTruePos(vtxPos); // store also the vertex from MC
567 AliInfo("getting the primary vertex from MC");
568 vertex = new AliESDVertex(vtxPos, vtxErr);
572 vertex->GetXYZ(vtxPos);
573 vertex->GetSigmaXYZ(vtxErr);
575 AliWarning("no vertex reconstructed");
576 vertex = new AliESDVertex(vtxPos, vtxErr);
578 esd->SetVertex(vertex);
579 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
580 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
584 AliInfo("execution time:");
585 ToAliInfo(stopwatch.Print());
590 //_____________________________________________________________________________
591 Bool_t AliReconstruction::RunHLTTracking(AliESD*& esd)
593 // run the HLT barrel tracking
595 TStopwatch stopwatch;
599 AliError("Missing runLoader!");
603 AliInfo("running HLT tracking");
605 // Get a pointer to the HLT reconstructor
606 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
607 if (!reconstructor) return kFALSE;
610 for (Int_t iDet = 1; iDet >= 0; iDet--) {
611 TString detName = fgkDetectorName[iDet];
612 AliDebug(1, Form("%s HLT tracking", detName.Data()));
613 reconstructor->SetOption(detName.Data());
614 AliTracker *tracker = reconstructor->CreateTracker(fRunLoader);
616 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
617 if (fStopOnError) return kFALSE;
620 Double_t vtxErr[3]={0.005,0.005,0.010};
621 const AliESDVertex *vertex = esd->GetVertex();
622 vertex->GetXYZ(vtxPos);
623 tracker->SetVertex(vtxPos,vtxErr);
625 fLoader[iDet]->LoadRecPoints("read");
626 TTree* tree = fLoader[iDet]->TreeR();
628 AliError(Form("Can't get the %s cluster tree", detName.Data()));
631 tracker->LoadClusters(tree);
633 if (tracker->Clusters2Tracks(esd) != 0) {
634 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
638 tracker->UnloadClusters();
643 AliInfo("execution time:");
644 ToAliInfo(stopwatch.Print());
649 //_____________________________________________________________________________
650 Bool_t AliReconstruction::RunTracking(AliESD*& esd)
652 // run the barrel tracking
654 TStopwatch stopwatch;
657 AliInfo("running tracking");
659 // pass 1: TPC + ITS inwards
660 for (Int_t iDet = 1; iDet >= 0; iDet--) {
661 if (!fTracker[iDet]) continue;
662 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
665 fLoader[iDet]->LoadRecPoints("read");
666 TTree* tree = fLoader[iDet]->TreeR();
668 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
671 fTracker[iDet]->LoadClusters(tree);
674 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
675 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
678 if (fCheckPointLevel > 1) {
679 WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
681 // preliminary PID in TPC needed by the ITS tracker
683 GetReconstructor(1)->FillESD(fRunLoader, esd);
684 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
685 AliESDpid::MakePID(esd);
689 // pass 2: ALL backwards
690 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
691 if (!fTracker[iDet]) continue;
692 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
695 if (iDet > 1) { // all except ITS, TPC
697 if (iDet == 3) { // TOF
698 fLoader[iDet]->LoadDigits("read");
699 tree = fLoader[iDet]->TreeD();
701 fLoader[iDet]->LoadRecPoints("read");
702 tree = fLoader[iDet]->TreeR();
705 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
708 fTracker[iDet]->LoadClusters(tree);
712 if (fTracker[iDet]->PropagateBack(esd) != 0) {
713 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
716 if (fCheckPointLevel > 1) {
717 WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
721 if (iDet > 2) { // all except ITS, TPC, TRD
722 fTracker[iDet]->UnloadClusters();
723 if (iDet == 3) { // TOF
724 fLoader[iDet]->UnloadDigits();
726 fLoader[iDet]->UnloadRecPoints();
731 // pass 3: TRD + TPC + ITS refit inwards
732 for (Int_t iDet = 2; iDet >= 0; iDet--) {
733 if (!fTracker[iDet]) continue;
734 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
737 if (fTracker[iDet]->RefitInward(esd) != 0) {
738 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
741 if (fCheckPointLevel > 1) {
742 WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
746 fTracker[iDet]->UnloadClusters();
747 fLoader[iDet]->UnloadRecPoints();
750 AliInfo("execution time:");
751 ToAliInfo(stopwatch.Print());
756 //_____________________________________________________________________________
757 Bool_t AliReconstruction::FillESD(AliESD*& esd, const TString& detectors)
759 // fill the event summary data
761 TStopwatch stopwatch;
763 AliInfo("filling ESD");
765 TString detStr = detectors;
766 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
767 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
768 AliReconstructor* reconstructor = GetReconstructor(iDet);
769 if (!reconstructor) continue;
771 if (!ReadESD(esd, fgkDetectorName[iDet])) {
772 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
773 TTree* clustersTree = NULL;
774 if (reconstructor->HasLocalReconstruction() && fLoader[iDet]) {
775 fLoader[iDet]->LoadRecPoints("read");
776 clustersTree = fLoader[iDet]->TreeR();
778 AliError(Form("Can't get the %s clusters tree",
779 fgkDetectorName[iDet]));
780 if (fStopOnError) return kFALSE;
783 if (fRawReader && !reconstructor->HasDigitConversion()) {
784 reconstructor->FillESD(fRawReader, clustersTree, esd);
786 TTree* digitsTree = NULL;
788 fLoader[iDet]->LoadDigits("read");
789 digitsTree = fLoader[iDet]->TreeD();
791 AliError(Form("Can't get the %s digits tree",
792 fgkDetectorName[iDet]));
793 if (fStopOnError) return kFALSE;
796 reconstructor->FillESD(digitsTree, clustersTree, esd);
797 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
799 if (reconstructor->HasLocalReconstruction() && fLoader[iDet]) {
800 fLoader[iDet]->UnloadRecPoints();
804 reconstructor->FillESD(fRunLoader, fRawReader, esd);
806 reconstructor->FillESD(fRunLoader, esd);
808 if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
812 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
813 AliError(Form("the following detectors were not found: %s",
815 if (fStopOnError) return kFALSE;
818 AliInfo("execution time:");
819 ToAliInfo(stopwatch.Print());
825 //_____________________________________________________________________________
826 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
828 // check whether detName is contained in detectors
829 // if yes, it is removed from detectors
831 // check if all detectors are selected
832 if ((detectors.CompareTo("ALL") == 0) ||
833 detectors.BeginsWith("ALL ") ||
834 detectors.EndsWith(" ALL") ||
835 detectors.Contains(" ALL ")) {
840 // search for the given detector
841 Bool_t result = kFALSE;
842 if ((detectors.CompareTo(detName) == 0) ||
843 detectors.BeginsWith(detName+" ") ||
844 detectors.EndsWith(" "+detName) ||
845 detectors.Contains(" "+detName+" ")) {
846 detectors.ReplaceAll(detName, "");
850 // clean up the detectors string
851 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
852 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
853 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
858 //_____________________________________________________________________________
859 Bool_t AliReconstruction::InitRunLoader()
861 // get or create the run loader
863 if (gAlice) delete gAlice;
866 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
867 // load all base libraries to get the loader classes
868 TString libs = gSystem->GetLibraries();
869 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
870 TString detName = fgkDetectorName[iDet];
871 if (detName == "HLT") continue;
872 if (libs.Contains("lib" + detName + "base.so")) continue;
873 gSystem->Load("lib" + detName + "base.so");
875 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
877 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
881 fRunLoader->CdGAFile();
882 if (gFile->GetKey(AliRunLoader::GetGAliceName())) {
883 if (fRunLoader->LoadgAlice() == 0) {
884 gAlice = fRunLoader->GetAliRun();
885 AliTracker::SetFieldMap(gAlice->Field());
888 if (!gAlice && !fRawReader) {
889 AliError(Form("no gAlice object found in file %s",
890 fGAliceFileName.Data()));
895 } else { // galice.root does not exist
897 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
901 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
902 AliConfig::GetDefaultEventFolderName(),
905 AliError(Form("could not create run loader in file %s",
906 fGAliceFileName.Data()));
910 fRunLoader->MakeTree("E");
912 while (fRawReader->NextEvent()) {
913 fRunLoader->SetEventNumber(iEvent);
914 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
916 fRunLoader->MakeTree("H");
917 fRunLoader->TreeE()->Fill();
920 fRawReader->RewindEvents();
921 fRunLoader->WriteHeader("OVERWRITE");
922 fRunLoader->CdGAFile();
923 fRunLoader->Write(0, TObject::kOverwrite);
924 // AliTracker::SetFieldMap(???);
930 //_____________________________________________________________________________
931 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
933 // get the reconstructor object and the loader for a detector
935 if (fReconstructor[iDet]) return fReconstructor[iDet];
937 // load the reconstructor object
938 TPluginManager* pluginManager = gROOT->GetPluginManager();
939 TString detName = fgkDetectorName[iDet];
940 TString recName = "Ali" + detName + "Reconstructor";
941 if (gAlice && !gAlice->GetDetector(detName) && (detName != "HLT")) return NULL;
943 if (detName == "HLT") {
944 if (!gROOT->GetClass("AliLevel3")) {
945 gSystem->Load("libAliL3Src.so");
946 gSystem->Load("libAliL3Misc.so");
947 gSystem->Load("libAliL3Hough.so");
948 gSystem->Load("libAliL3Comp.so");
952 AliReconstructor* reconstructor = NULL;
953 // first check if a plugin is defined for the reconstructor
954 TPluginHandler* pluginHandler =
955 pluginManager->FindHandler("AliReconstructor", detName);
956 // if not, add a plugin for it
957 if (!pluginHandler) {
958 AliDebug(1, Form("defining plugin for %s", recName.Data()));
959 TString libs = gSystem->GetLibraries();
960 if (libs.Contains("lib" + detName + "base.so") ||
961 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
962 pluginManager->AddHandler("AliReconstructor", detName,
963 recName, detName + "rec", recName + "()");
965 pluginManager->AddHandler("AliReconstructor", detName,
966 recName, detName, recName + "()");
968 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
970 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
971 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
974 TObject* obj = fOptions.FindObject(detName.Data());
975 if (obj) reconstructor->SetOption(obj->GetTitle());
976 reconstructor->Init(fRunLoader);
977 fReconstructor[iDet] = reconstructor;
980 // get or create the loader
981 if (detName != "HLT") {
982 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
983 if (!fLoader[iDet]) {
984 AliConfig::Instance()
985 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
987 // first check if a plugin is defined for the loader
988 TPluginHandler* pluginHandler =
989 pluginManager->FindHandler("AliLoader", detName);
990 // if not, add a plugin for it
991 if (!pluginHandler) {
992 TString loaderName = "Ali" + detName + "Loader";
993 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
994 pluginManager->AddHandler("AliLoader", detName,
995 loaderName, detName + "base",
996 loaderName + "(const char*, TFolder*)");
997 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
999 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
1001 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
1002 fRunLoader->GetEventFolder());
1004 if (!fLoader[iDet]) { // use default loader
1005 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
1007 if (!fLoader[iDet]) {
1008 AliWarning(Form("couldn't get loader for %s", detName.Data()));
1009 if (fStopOnError) return NULL;
1011 fRunLoader->AddLoader(fLoader[iDet]);
1012 fRunLoader->CdGAFile();
1013 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
1014 fRunLoader->Write(0, TObject::kOverwrite);
1019 return reconstructor;
1022 //_____________________________________________________________________________
1023 Bool_t AliReconstruction::CreateVertexer()
1025 // create the vertexer
1028 AliReconstructor* itsReconstructor = GetReconstructor(0);
1029 if (itsReconstructor) {
1030 fVertexer = itsReconstructor->CreateVertexer(fRunLoader);
1033 AliWarning("couldn't create a vertexer for ITS");
1034 if (fStopOnError) return kFALSE;
1040 //_____________________________________________________________________________
1041 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
1043 // create the trackers
1045 TString detStr = detectors;
1046 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1047 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1048 AliReconstructor* reconstructor = GetReconstructor(iDet);
1049 if (!reconstructor) continue;
1050 TString detName = fgkDetectorName[iDet];
1051 if (detName == "HLT") {
1052 fRunHLTTracking = kTRUE;
1056 fTracker[iDet] = reconstructor->CreateTracker(fRunLoader);
1057 if (!fTracker[iDet] && (iDet < 7)) {
1058 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1059 if (fStopOnError) return kFALSE;
1066 //_____________________________________________________________________________
1067 void AliReconstruction::CleanUp(TFile* file, TFile* fileOld)
1069 // delete trackers and the run loader and close and delete the file
1071 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1072 delete fReconstructor[iDet];
1073 fReconstructor[iDet] = NULL;
1074 fLoader[iDet] = NULL;
1075 delete fTracker[iDet];
1076 fTracker[iDet] = NULL;
1094 gSystem->Unlink("AliESDs.old.root");
1099 //_____________________________________________________________________________
1100 Bool_t AliReconstruction::ReadESD(AliESD*& esd, const char* recStep) const
1102 // read the ESD event from a file
1104 if (!esd) return kFALSE;
1106 sprintf(fileName, "ESD_%d.%d_%s.root",
1107 esd->GetRunNumber(), esd->GetEventNumber(), recStep);
1108 if (gSystem->AccessPathName(fileName)) return kFALSE;
1110 AliDebug(1, Form("reading ESD from file %s", fileName));
1111 TFile* file = TFile::Open(fileName);
1112 if (!file || !file->IsOpen()) {
1113 AliError(Form("opening %s failed", fileName));
1120 esd = (AliESD*) file->Get("ESD");
1126 //_____________________________________________________________________________
1127 void AliReconstruction::WriteESD(AliESD* esd, const char* recStep) const
1129 // write the ESD event to a file
1133 sprintf(fileName, "ESD_%d.%d_%s.root",
1134 esd->GetRunNumber(), esd->GetEventNumber(), recStep);
1136 AliDebug(1, Form("writing ESD to file %s", fileName));
1137 TFile* file = TFile::Open(fileName, "recreate");
1138 if (!file || !file->IsOpen()) {
1139 AliError(Form("opening %s failed", fileName));