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. //
29 // The name of the galice file can be changed from the default //
30 // "galice.root" by passing it as argument to the AliReconstruction //
31 // constructor or by //
33 // rec.SetGAliceFile("..."); //
35 // The local reconstruction can be switched on or off for individual //
38 // rec.SetRunLocalReconstruction("..."); //
40 // The argument is a (case sensitive) string with the names of the //
41 // detectors separated by a space. The special string "ALL" selects all //
42 // available detectors. This is the default. //
44 // The tracking in ITS, TPC and TRD and the creation of ESD tracks can be //
47 // rec.SetRunTracking(kFALSE); //
49 // The filling of additional ESD information can be steered by //
51 // rec.SetFillESD("..."); //
53 // Again, the string specifies the list of detectors. The default is "ALL". //
55 // The reconstruction requires digits as input. For the creation of digits //
56 // have a look at the class AliSimulation. //
58 // For debug purposes the method SetCheckPointLevel can be used. If the //
59 // argument is greater than 0, files with ESD events will be written after //
60 // selected steps of the reconstruction for each event: //
61 // level 1: after tracking and after filling of ESD (final) //
62 // level 2: in addition after each tracking step //
63 // level 3: in addition after the filling of ESD for each detector //
64 // If a final check point file exists for an event, this event will be //
65 // skipped in the reconstruction. The tracking and the filling of ESD for //
66 // a detector will be skipped as well, if the corresponding check point //
67 // file exists. The ESD event will then be loaded from the file instead. //
69 ///////////////////////////////////////////////////////////////////////////////
75 #include <TPluginManager.h>
77 #include "AliReconstruction.h"
78 #include "AliRunLoader.h"
80 #include "AliTracker.h"
82 #include "AliESDVertex.h"
83 #include "AliVertexer.h"
84 #include "AliHeader.h"
85 #include "AliGenEventHeader.h"
86 #include "AliESDpid.h"
89 ClassImp(AliReconstruction)
92 //_____________________________________________________________________________
93 const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "RICH", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "START", "VZERO", "CRT"};
95 //_____________________________________________________________________________
96 AliReconstruction::AliReconstruction(const char* gAliceFilename,
97 const char* name, const char* title) :
100 fRunLocalReconstruction("ALL"),
101 fRunVertexFinder(kTRUE),
104 fGAliceFileName(gAliceFilename),
105 fStopOnError(kFALSE),
119 // create reconstruction object with default parameters
123 //_____________________________________________________________________________
124 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
127 fRunLocalReconstruction(rec.fRunLocalReconstruction),
128 fRunVertexFinder(rec.fRunVertexFinder),
129 fRunTracking(rec.fRunTracking),
130 fFillESD(rec.fFillESD),
131 fGAliceFileName(rec.fGAliceFileName),
132 fStopOnError(rec.fStopOnError),
150 //_____________________________________________________________________________
151 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
153 // assignment operator
155 this->~AliReconstruction();
156 new(this) AliReconstruction(rec);
160 //_____________________________________________________________________________
161 AliReconstruction::~AliReconstruction()
169 //_____________________________________________________________________________
170 void AliReconstruction::SetGAliceFile(const char* fileName)
172 // set the name of the galice file
174 fGAliceFileName = fileName;
178 //_____________________________________________________________________________
179 Bool_t AliReconstruction::Run()
181 // run the reconstruction
183 // open the run loader
184 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
186 Error("Run", "no run loader found in file %s",
187 fGAliceFileName.Data());
191 fRunLoader->LoadgAlice();
192 AliRun* aliRun = fRunLoader->GetAliRun();
194 Error("Run", "no gAlice object found in file %s",
195 fGAliceFileName.Data());
201 // load the reconstructor objects
202 TPluginManager* pluginManager = gROOT->GetPluginManager();
203 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
204 TString detName = fgkDetectorName[iDet];
205 TString recName = "Ali" + detName + "Reconstructor";
206 if (!gAlice->GetDetector(detName)) continue;
208 AliReconstructor* reconstructor = NULL;
209 // first check if a plugin is defined for the reconstructor
210 TPluginHandler* pluginHandler =
211 pluginManager->FindHandler("AliReconstructor", detName);
212 // if not, but the reconstructor class is implemented, add a plugin for it
213 if (!pluginHandler && gROOT->GetClass(recName.Data())) {
214 Info("Run", "defining plugin for %s", recName.Data());
215 pluginManager->AddHandler("AliReconstructor", detName,
216 recName, detName, recName + "()");
217 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
219 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
220 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
222 // if there is no reconstructor class for the detector use the dummy one
223 if (!reconstructor) {
224 Info("Run", "using dummy reconstructor for %s", detName.Data());
225 reconstructor = new AliDummyReconstructor(gAlice->GetDetector(detName));
227 if (reconstructor) fReconstructors.Add(reconstructor);
230 // local reconstruction
231 if (!fRunLocalReconstruction.IsNull()) {
232 if (!RunLocalReconstruction(fRunLocalReconstruction)) {
233 if (fStopOnError) {CleanUp(); return kFALSE;}
236 if (!fRunVertexFinder && !fRunTracking && fFillESD.IsNull()) return kTRUE;
239 if (fRunVertexFinder && !CreateVertexer()) {
246 // get loaders and trackers
247 if (fRunTracking && !CreateTrackers()) {
254 // create the ESD output file
255 TFile* file = TFile::Open("AliESDs.root", "RECREATE");
256 if (!file->IsOpen()) {
257 Error("Run", "opening AliESDs.root failed");
258 if (fStopOnError) {CleanUp(file); return kFALSE;}
262 for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) {
263 Info("Run", "processing event %d", iEvent);
264 fRunLoader->GetEvent(iEvent);
267 sprintf(fileName, "ESD_%d.%d_final.root",
268 aliRun->GetRunNumber(), aliRun->GetEvNumber());
269 if (!gSystem->AccessPathName(fileName)) continue;
271 AliESD* esd = new AliESD;
272 esd->SetRunNumber(aliRun->GetRunNumber());
273 esd->SetEventNumber(aliRun->GetEvNumber());
274 esd->SetMagneticField(aliRun->Field()->SolenoidField());
277 if (fRunVertexFinder) {
278 if (!ReadESD(esd, "vertex")) {
279 if (!RunVertexFinder(esd)) {
280 if (fStopOnError) {CleanUp(file); return kFALSE;}
282 if (fCheckPointLevel > 0) WriteESD(esd, "vertex");
288 if (!ReadESD(esd, "tracking")) {
289 if (!RunTracking(esd)) {
290 if (fStopOnError) {CleanUp(file); return kFALSE;}
292 if (fCheckPointLevel > 0) WriteESD(esd, "tracking");
297 if (!fFillESD.IsNull()) {
298 if (!FillESD(esd, fFillESD)) {
299 if (fStopOnError) {CleanUp(file); return kFALSE;}
304 AliESDpid::MakePID(esd);
305 if (fCheckPointLevel > 1) WriteESD(esd, "PID");
309 sprintf(name, "ESD%d", iEvent);
311 if (!esd->Write(name)) {
312 Error("Run", "writing ESD failed");
313 if (fStopOnError) {CleanUp(file); return kFALSE;}
317 if (fCheckPointLevel > 0) WriteESD(esd, "final");
327 //_____________________________________________________________________________
328 Bool_t AliReconstruction::RunLocalReconstruction(const TString& detectors)
330 // run the local reconstruction
332 TStopwatch stopwatch;
335 TString detStr = detectors;
336 for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
337 AliReconstructor* reconstructor =
338 (AliReconstructor*) fReconstructors[iDet];
339 TString detName = reconstructor->GetDetectorName();
340 if (IsSelected(detName, detStr)) {
341 Info("RunReconstruction", "running reconstruction for %s",
343 TStopwatch stopwatchDet;
344 stopwatchDet.Start();
345 reconstructor->Reconstruct(fRunLoader);
346 Info("RunReconstruction", "execution time for %s:", detName.Data());
347 stopwatchDet.Print();
351 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
352 Error("RunReconstruction", "the following detectors were not found: %s",
354 if (fStopOnError) return kFALSE;
357 Info("RunReconstruction", "execution time:");
363 //_____________________________________________________________________________
364 Bool_t AliReconstruction::RunVertexFinder(AliESD*& esd)
366 // run the barrel tracking
368 TStopwatch stopwatch;
371 AliESDVertex* vertex = NULL;
372 Double_t vtxPos[3] = {0, 0, 0};
373 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
375 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
376 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
379 Info("RunVertexFinder", "running the ITS vertex finder");
380 fITSVertexer->SetDebug(1);
381 vertex = fITSVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
383 Warning("RunVertexFinder","Vertex not found \n");
384 vertex = new AliESDVertex();
387 vertex->SetTruePos(vtxPos); // store also the vertex from MC
391 Info("RunVertexFinder", "getting the primary vertex from MC");
392 vertex = new AliESDVertex(vtxPos, vtxErr);
396 vertex->GetXYZ(vtxPos);
397 vertex->GetSigmaXYZ(vtxErr);
399 Warning("RunVertexFinder", "no vertex reconstructed");
400 vertex = new AliESDVertex(vtxPos, vtxErr);
402 esd->SetVertex(vertex);
403 if (fITSTracker) fITSTracker->SetVertex(vtxPos, vtxErr);
404 if (fTPCTracker) fTPCTracker->SetVertex(vtxPos, vtxErr);
405 if (fTRDTracker) fTRDTracker->SetVertex(vtxPos, vtxErr);
408 Info("RunVertexFinder", "execution time:");
414 //_____________________________________________________________________________
415 Bool_t AliReconstruction::RunTracking(AliESD*& esd)
417 // run the barrel tracking
419 TStopwatch stopwatch;
423 Error("RunTracking", "no TPC tracker");
428 Info("RunTracking", "TPC tracking");
429 fTPCLoader->LoadRecPoints("read");
430 TTree* tpcTree = fTPCLoader->TreeR();
432 Error("RunTracking", "Can't get the TPC cluster tree");
435 fTPCTracker->LoadClusters(tpcTree);
436 if (fTPCTracker->Clusters2Tracks(esd) != 0) {
437 Error("RunTracking", "TPC Clusters2Tracks failed");
440 if (fCheckPointLevel > 1) WriteESD(esd, "TPC.tracking");
443 Warning("RunTracking", "no ITS tracker");
446 fRunLoader->GetAliRun()->GetDetector("TPC")->FillESD(esd); // preliminary
447 AliESDpid::MakePID(esd); // PID for the ITS tracker
450 Info("RunTracking", "ITS tracking");
451 fITSLoader->LoadRecPoints("read");
452 TTree* itsTree = fITSLoader->TreeR();
454 Error("RunTracking", "Can't get the ITS cluster tree");
457 fITSTracker->LoadClusters(itsTree);
458 if (fITSTracker->Clusters2Tracks(esd) != 0) {
459 Error("RunTracking", "ITS Clusters2Tracks failed");
462 if (fCheckPointLevel > 1) WriteESD(esd, "ITS.tracking");
465 Warning("RunTracking", "no TRD tracker");
467 // ITS back propagation
468 Info("RunTracking", "ITS back propagation");
469 if (fITSTracker->PropagateBack(esd) != 0) {
470 Error("RunTracking", "ITS backward propagation failed");
473 if (fCheckPointLevel > 1) WriteESD(esd, "ITS.back");
475 // TPC back propagation
476 Info("RunTracking", "TPC back propagation");
477 if (fTPCTracker->PropagateBack(esd) != 0) {
478 Error("RunTracking", "TPC backward propagation failed");
481 if (fCheckPointLevel > 1) WriteESD(esd, "TPC.back");
483 // TRD back propagation
484 Info("RunTracking", "TRD back propagation");
485 fTRDLoader->LoadRecPoints("read");
486 TTree* trdTree = fTRDLoader->TreeR();
488 Error("RunTracking", "Can't get the TRD cluster tree");
491 fTRDTracker->LoadClusters(trdTree);
492 if (fTRDTracker->PropagateBack(esd) != 0) {
493 Error("RunTracking", "TRD backward propagation failed");
496 if (fCheckPointLevel > 1) WriteESD(esd, "TRD.back");
499 Warning("RunTracking", "no TOF tracker");
501 // TOF back propagation
502 Info("RunTracking", "TOF back propagation");
503 fTOFLoader->LoadDigits("read");
504 TTree* tofTree = fTOFLoader->TreeD();
506 Error("RunTracking", "Can't get the TOF digits tree");
509 fTOFTracker->LoadClusters(tofTree);
510 if (fTOFTracker->PropagateBack(esd) != 0) {
511 Error("RunTracking", "TOF backward propagation failed");
514 if (fCheckPointLevel > 1) WriteESD(esd, "TOF.back");
515 fTOFTracker->UnloadClusters();
516 fTOFLoader->UnloadDigits();
520 Info("RunTracking", "TRD inward refit");
521 if (fTRDTracker->RefitInward(esd) != 0) {
522 Error("RunTracking", "TRD inward refit failed");
525 if (fCheckPointLevel > 1) WriteESD(esd, "TRD.refit");
526 fTRDTracker->UnloadClusters();
527 fTRDLoader->UnloadRecPoints();
530 Info("RunTracking", "TPC inward refit");
531 if (fTPCTracker->RefitInward(esd) != 0) {
532 Error("RunTracking", "TPC inward refit failed");
535 if (fCheckPointLevel > 1) WriteESD(esd, "TPC.refit");
538 Info("RunTracking", "ITS inward refit");
539 if (fITSTracker->RefitInward(esd) != 0) {
540 Error("RunTracking", "ITS inward refit failed");
543 if (fCheckPointLevel > 1) WriteESD(esd, "ITS.refit");
546 fITSTracker->UnloadClusters();
547 fITSLoader->UnloadRecPoints();
550 fTPCTracker->UnloadClusters();
551 fTPCLoader->UnloadRecPoints();
553 Info("RunTracking", "execution time:");
559 //_____________________________________________________________________________
560 Bool_t AliReconstruction::FillESD(AliESD*& esd, const TString& detectors)
562 // fill the event summary data
564 TStopwatch stopwatch;
567 TString detStr = detectors;
568 for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
569 AliReconstructor* reconstructor =
570 (AliReconstructor*) fReconstructors[iDet];
571 TString detName = reconstructor->GetDetectorName();
572 if (IsSelected(detName, detStr)) {
573 if (!ReadESD(esd, detName.Data())) {
574 Info("FillESD", "filling ESD for %s", detName.Data());
575 reconstructor->FillESD(fRunLoader, esd);
576 if (fCheckPointLevel > 2) WriteESD(esd, detName.Data());
581 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
582 Error("FillESD", "the following detectors were not found: %s",
584 if (fStopOnError) return kFALSE;
587 Info("FillESD", "execution time:");
594 //_____________________________________________________________________________
595 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
597 // check whether detName is contained in detectors
598 // if yes, it is removed from detectors
600 // check if all detectors are selected
601 if ((detectors.CompareTo("ALL") == 0) ||
602 detectors.BeginsWith("ALL ") ||
603 detectors.EndsWith(" ALL") ||
604 detectors.Contains(" ALL ")) {
609 // search for the given detector
610 Bool_t result = kFALSE;
611 if ((detectors.CompareTo(detName) == 0) ||
612 detectors.BeginsWith(detName+" ") ||
613 detectors.EndsWith(" "+detName) ||
614 detectors.Contains(" "+detName+" ")) {
615 detectors.ReplaceAll(detName, "");
619 // clean up the detectors string
620 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
621 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
622 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
627 //_____________________________________________________________________________
628 AliReconstructor* AliReconstruction::GetReconstructor(const char* detName) const
630 // get the reconstructor object for a detector
632 for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
633 AliReconstructor* reconstructor =
634 (AliReconstructor*) fReconstructors[iDet];
635 if (strcmp(reconstructor->GetDetectorName(), detName) == 0) {
636 return reconstructor;
642 //_____________________________________________________________________________
643 Bool_t AliReconstruction::CreateVertexer()
645 // create the vertexer
648 AliReconstructor* itsReconstructor = GetReconstructor("ITS");
649 if (itsReconstructor) {
650 fITSVertexer = itsReconstructor->CreateVertexer(fRunLoader);
653 Warning("CreateVertexer", "couldn't create a vertexer for ITS");
654 if (fStopOnError) return kFALSE;
660 //_____________________________________________________________________________
661 Bool_t AliReconstruction::CreateTrackers()
663 // get the loaders and create the trackers
666 fITSLoader = fRunLoader->GetLoader("ITSLoader");
668 Warning("CreateTrackers", "no ITS loader found");
669 if (fStopOnError) return kFALSE;
671 AliReconstructor* itsReconstructor = GetReconstructor("ITS");
672 if (itsReconstructor) {
673 fITSTracker = itsReconstructor->CreateTracker(fRunLoader);
676 Warning("CreateTrackers", "couldn't create a tracker for ITS");
677 if (fStopOnError) return kFALSE;
682 fTPCLoader = fRunLoader->GetLoader("TPCLoader");
684 Error("CreateTrackers", "no TPC loader found");
685 if (fStopOnError) return kFALSE;
687 AliReconstructor* tpcReconstructor = GetReconstructor("TPC");
688 if (tpcReconstructor) {
689 fTPCTracker = tpcReconstructor->CreateTracker(fRunLoader);
692 Error("CreateTrackers", "couldn't create a tracker for TPC");
693 if (fStopOnError) return kFALSE;
698 fTRDLoader = fRunLoader->GetLoader("TRDLoader");
700 Warning("CreateTrackers", "no TRD loader found");
701 if (fStopOnError) return kFALSE;
703 AliReconstructor* trdReconstructor = GetReconstructor("TRD");
704 if (trdReconstructor) {
705 fTRDTracker = trdReconstructor->CreateTracker(fRunLoader);
708 Warning("CreateTrackers", "couldn't create a tracker for TRD");
709 if (fStopOnError) return kFALSE;
714 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
716 Warning("CreateTrackers", "no TOF loader found");
717 if (fStopOnError) return kFALSE;
719 AliReconstructor* tofReconstructor = GetReconstructor("TOF");
720 if (tofReconstructor) {
721 fTOFTracker = tofReconstructor->CreateTracker(fRunLoader);
724 Warning("CreateTrackers", "couldn't create a tracker for TOF");
725 if (fStopOnError) return kFALSE;
732 //_____________________________________________________________________________
733 void AliReconstruction::CleanUp(TFile* file)
735 // delete trackers and the run loader and close and delete the file
737 fReconstructors.Delete();
760 //_____________________________________________________________________________
761 Bool_t AliReconstruction::ReadESD(AliESD*& esd, const char* recStep) const
763 // read the ESD event from a file
765 if (!esd) return kFALSE;
767 sprintf(fileName, "ESD_%d.%d_%s.root",
768 esd->GetRunNumber(), esd->GetEventNumber(), recStep);
769 if (gSystem->AccessPathName(fileName)) return kFALSE;
771 Info("ReadESD", "reading ESD from file %s", fileName);
772 TFile* file = TFile::Open(fileName);
773 if (!file || !file->IsOpen()) {
774 Error("ReadESD", "opening %s failed", fileName);
781 esd = (AliESD*) file->Get("ESD");
787 //_____________________________________________________________________________
788 void AliReconstruction::WriteESD(AliESD* esd, const char* recStep) const
790 // write the ESD event to a file
794 sprintf(fileName, "ESD_%d.%d_%s.root",
795 esd->GetRunNumber(), esd->GetEventNumber(), recStep);
797 Info("WriteESD", "writing ESD to file %s", fileName);
798 TFile* file = TFile::Open(fileName, "recreate");
799 if (!file || !file->IsOpen()) {
800 Error("WriteESD", "opening %s failed", fileName);