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 // 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 //
45 // rec.SetGAliceFile("..."); //
47 // The local reconstruction can be switched on or off for individual //
50 // rec.SetRunLocalReconstruction("..."); //
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. //
56 // The reconstruction of the primary vertex position can be switched off by //
58 // rec.SetRunVertexFinder(kFALSE); //
60 // The tracking in ITS, TPC and TRD and the creation of ESD tracks can be //
63 // rec.SetRunTracking(kFALSE); //
65 // The filling of additional ESD information can be steered by //
67 // rec.SetFillESD("..."); //
69 // Again, the string specifies the list of detectors. The default is "ALL". //
71 // The reconstruction requires digits or raw data as input. For the creation //
72 // of digits and raw data have a look at the class AliSimulation. //
74 // For debug purposes the method SetCheckPointLevel can be used. If the //
75 // argument is greater than 0, files with ESD events will be written after //
76 // selected steps of the reconstruction for each event: //
77 // level 1: after tracking and after filling of ESD (final) //
78 // level 2: in addition after each tracking step //
79 // level 3: in addition after the filling of ESD for each detector //
80 // If a final check point file exists for an event, this event will be //
81 // skipped in the reconstruction. The tracking and the filling of ESD for //
82 // a detector will be skipped as well, if the corresponding check point //
83 // file exists. The ESD event will then be loaded from the file instead. //
85 ///////////////////////////////////////////////////////////////////////////////
91 #include <TPluginManager.h>
92 #include <TStopwatch.h>
94 #include "AliReconstruction.h"
96 #include "AliRunLoader.h"
98 #include "AliRawReaderFile.h"
99 #include "AliRawReaderDate.h"
100 #include "AliRawReaderRoot.h"
101 #include "AliTracker.h"
103 #include "AliESDVertex.h"
104 #include "AliVertexer.h"
105 #include "AliHeader.h"
106 #include "AliGenEventHeader.h"
107 #include "AliESDpid.h"
110 ClassImp(AliReconstruction)
113 //_____________________________________________________________________________
114 const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "RICH", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "START", "VZERO", "CRT", "HLT"};
116 //_____________________________________________________________________________
117 AliReconstruction::AliReconstruction(const char* gAliceFilename,
118 const char* name, const char* title) :
121 fRunLocalReconstruction("ALL"),
122 fRunVertexFinder(kTRUE),
125 fGAliceFileName(gAliceFilename),
127 fStopOnError(kFALSE),
145 // create reconstruction object with default parameters
149 //_____________________________________________________________________________
150 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
153 fRunLocalReconstruction(rec.fRunLocalReconstruction),
154 fRunVertexFinder(rec.fRunVertexFinder),
155 fRunTracking(rec.fRunTracking),
156 fFillESD(rec.fFillESD),
157 fGAliceFileName(rec.fGAliceFileName),
159 fStopOnError(rec.fStopOnError),
179 for (Int_t i = 0; i < fOptions.GetEntriesFast(); i++) {
180 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
184 //_____________________________________________________________________________
185 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
187 // assignment operator
189 this->~AliReconstruction();
190 new(this) AliReconstruction(rec);
194 //_____________________________________________________________________________
195 AliReconstruction::~AliReconstruction()
204 //_____________________________________________________________________________
205 void AliReconstruction::SetGAliceFile(const char* fileName)
207 // set the name of the galice file
209 fGAliceFileName = fileName;
212 //_____________________________________________________________________________
213 void AliReconstruction::SetOption(const char* detector, const char* option)
215 // set options for the reconstruction of a detector
217 TObject* obj = fOptions.FindObject(detector);
218 if (obj) fOptions.Remove(obj);
219 fOptions.Add(new TNamed(detector, option));
223 //_____________________________________________________________________________
224 Bool_t AliReconstruction::Run(const char* input)
226 // run the reconstruction
229 if (!input) input = fInput.Data();
230 TString fileName(input);
231 if (fileName.EndsWith("/")) {
232 fRawReader = new AliRawReaderFile(fileName);
233 } else if (fileName.EndsWith(".root")) {
234 fRawReader = new AliRawReaderRoot(fileName);
235 } else if (!fileName.IsNull()) {
236 fRawReader = new AliRawReaderDate(fileName);
237 fRawReader->SelectEvents(7);
240 // open the run loader
241 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
243 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
247 fRunLoader->LoadgAlice();
248 AliRun* aliRun = fRunLoader->GetAliRun();
250 AliError(Form("no gAlice object found in file %s",
251 fGAliceFileName.Data()));
256 AliTracker::SetFieldMap(gAlice->Field());
258 // load the reconstructor objects
259 TPluginManager* pluginManager = gROOT->GetPluginManager();
260 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
261 TString detName = fgkDetectorName[iDet];
262 TString recName = "Ali" + detName + "Reconstructor";
263 if (!gAlice->GetDetector(detName) && detName != "HLT") continue;
265 if(detName == "HLT") {
266 if (!gROOT->GetClass("AliLevel3")) {
267 gSystem->Load("libAliL3Src.so");
268 gSystem->Load("libAliL3Misc.so");
269 gSystem->Load("libAliL3Hough.so");
270 gSystem->Load("libAliL3Comp.so");
274 AliReconstructor* reconstructor = NULL;
275 // first check if a plugin is defined for the reconstructor
276 TPluginHandler* pluginHandler =
277 pluginManager->FindHandler("AliReconstructor", detName);
278 // if not, but the reconstructor class is implemented, add a plugin for it
279 if (!pluginHandler && gROOT->GetClass(recName.Data())) {
280 AliDebug(1, Form("defining plugin for %s", recName.Data()));
281 pluginManager->AddHandler("AliReconstructor", detName,
282 recName, detName, recName + "()");
283 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
285 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
286 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
288 // if there is no reconstructor class for the detector use the dummy one
289 if (!reconstructor && gAlice->GetDetector(detName)) {
290 AliDebug(1, Form("using dummy reconstructor for %s", detName.Data()));
291 reconstructor = new AliDummyReconstructor(gAlice->GetDetector(detName));
294 TObject* obj = fOptions.FindObject(detName.Data());
295 if (obj) reconstructor->SetOption(obj->GetTitle());
296 fReconstructors.Add(reconstructor);
300 // local reconstruction
301 if (!fRunLocalReconstruction.IsNull()) {
302 if (!RunLocalReconstruction(fRunLocalReconstruction)) {
303 if (fStopOnError) {CleanUp(); return kFALSE;}
306 if (!fRunVertexFinder && !fRunTracking && fFillESD.IsNull()) return kTRUE;
309 if (fRunVertexFinder && !CreateVertexer()) {
316 // get loaders and trackers
317 if (fRunTracking && !CreateTrackers()) {
324 // create the ESD output file and tree
325 TFile* file = TFile::Open("AliESDs.root", "RECREATE");
326 if (!file->IsOpen()) {
327 AliError("opening AliESDs.root failed");
328 if (fStopOnError) {CleanUp(file); return kFALSE;}
330 AliESD* esd = new AliESD;
331 TTree* tree = new TTree("esdTree", "Tree with ESD objects");
332 tree->Branch("ESD", "AliESD", &esd);
337 if (fRawReader) fRawReader->RewindEvents();
338 for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) {
339 AliInfo(Form("processing event %d", iEvent));
340 fRunLoader->GetEvent(iEvent);
341 if (fRawReader) fRawReader->NextEvent();
344 sprintf(fileName, "ESD_%d.%d_final.root",
345 aliRun->GetRunNumber(), aliRun->GetEvNumber());
346 if (!gSystem->AccessPathName(fileName)) continue;
349 esd->SetRunNumber(aliRun->GetRunNumber());
350 esd->SetEventNumber(aliRun->GetEvNumber());
351 esd->SetMagneticField(aliRun->Field()->SolenoidField());
354 if (fRunVertexFinder) {
355 if (!ReadESD(esd, "vertex")) {
356 if (!RunVertexFinder(esd)) {
357 if (fStopOnError) {CleanUp(file); return kFALSE;}
359 if (fCheckPointLevel > 0) WriteESD(esd, "vertex");
365 if (!ReadESD(esd, "tracking")) {
366 if (!RunTracking(esd)) {
367 if (fStopOnError) {CleanUp(file); return kFALSE;}
369 if (fCheckPointLevel > 0) WriteESD(esd, "tracking");
374 if (!fFillESD.IsNull()) {
375 if (!FillESD(esd, fFillESD)) {
376 if (fStopOnError) {CleanUp(file); return kFALSE;}
381 AliESDpid::MakePID(esd);
382 if (fCheckPointLevel > 1) WriteESD(esd, "PID");
387 if (fCheckPointLevel > 0) WriteESD(esd, "final");
399 //_____________________________________________________________________________
400 Bool_t AliReconstruction::RunLocalReconstruction(const TString& detectors)
402 // run the local reconstruction
404 TStopwatch stopwatch;
407 TString detStr = detectors;
408 for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
409 AliReconstructor* reconstructor =
410 (AliReconstructor*) fReconstructors[iDet];
411 TString detName = reconstructor->GetDetectorName();
412 if (IsSelected(detName, detStr)) {
413 AliInfo(Form("running reconstruction for %s", detName.Data()));
414 TStopwatch stopwatchDet;
415 stopwatchDet.Start();
417 fRawReader->RewindEvents();
418 reconstructor->Reconstruct(fRunLoader, fRawReader);
420 reconstructor->Reconstruct(fRunLoader);
422 AliInfo(Form("execution time for %s:", detName.Data()));
423 ToAliInfo(stopwatchDet.Print());
427 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
428 AliError(Form("the following detectors were not found: %s",
430 if (fStopOnError) return kFALSE;
433 AliInfo("execution time:");
434 ToAliInfo(stopwatch.Print());
439 //_____________________________________________________________________________
440 Bool_t AliReconstruction::RunVertexFinder(AliESD*& esd)
442 // run the barrel tracking
444 TStopwatch stopwatch;
447 AliESDVertex* vertex = NULL;
448 Double_t vtxPos[3] = {0, 0, 0};
449 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
451 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
452 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
453 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
457 AliInfo("running the ITS vertex finder");
458 fITSVertexer->SetDebug(1);
459 vertex = fITSVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
461 AliWarning("Vertex not found");
462 vertex = new AliESDVertex();
465 vertex->SetTruePos(vtxPos); // store also the vertex from MC
469 AliInfo("getting the primary vertex from MC");
470 vertex = new AliESDVertex(vtxPos, vtxErr);
474 vertex->GetXYZ(vtxPos);
475 vertex->GetSigmaXYZ(vtxErr);
477 AliWarning("no vertex reconstructed");
478 vertex = new AliESDVertex(vtxPos, vtxErr);
480 esd->SetVertex(vertex);
481 if (fITSTracker) fITSTracker->SetVertex(vtxPos, vtxErr);
482 if (fTPCTracker) fTPCTracker->SetVertex(vtxPos, vtxErr);
483 if (fTRDTracker) fTRDTracker->SetVertex(vtxPos, vtxErr);
486 AliInfo("execution time:");
487 ToAliInfo(stopwatch.Print());
492 //_____________________________________________________________________________
493 Bool_t AliReconstruction::RunTracking(AliESD*& esd)
495 // run the barrel tracking
497 TStopwatch stopwatch;
501 AliError("no TPC tracker");
504 AliInfo("running tracking");
507 AliDebug(1, "TPC tracking");
508 fTPCLoader->LoadRecPoints("read");
509 TTree* tpcTree = fTPCLoader->TreeR();
511 AliError("Can't get the TPC cluster tree");
514 fTPCTracker->LoadClusters(tpcTree);
515 if (fTPCTracker->Clusters2Tracks(esd) != 0) {
516 AliError("TPC Clusters2Tracks failed");
519 if (fCheckPointLevel > 1) WriteESD(esd, "TPC.tracking");
522 AliWarning("no ITS tracker");
525 GetReconstructor("TPC")->FillESD(fRunLoader, esd); // preliminary
526 AliESDpid::MakePID(esd); // PID for the ITS tracker
529 AliDebug(1, "ITS tracking");
530 fITSLoader->LoadRecPoints("read");
531 TTree* itsTree = fITSLoader->TreeR();
533 Error("RunTracking", "Can't get the ITS cluster tree");
536 fITSTracker->LoadClusters(itsTree);
537 if (fITSTracker->Clusters2Tracks(esd) != 0) {
538 AliError("ITS Clusters2Tracks failed");
541 if (fCheckPointLevel > 1) WriteESD(esd, "ITS.tracking");
544 AliWarning("no TRD tracker");
546 // ITS back propagation
547 AliDebug(1, "ITS back propagation");
548 if (fITSTracker->PropagateBack(esd) != 0) {
549 AliError("ITS backward propagation failed");
552 if (fCheckPointLevel > 1) WriteESD(esd, "ITS.back");
554 // TPC back propagation
555 AliDebug(1, "TPC back propagation");
556 if (fTPCTracker->PropagateBack(esd) != 0) {
557 AliError("TPC backward propagation failed");
560 if (fCheckPointLevel > 1) WriteESD(esd, "TPC.back");
562 // TRD back propagation
563 AliDebug(1, "TRD back propagation");
564 fTRDLoader->LoadRecPoints("read");
565 TTree* trdTree = fTRDLoader->TreeR();
567 AliError("Can't get the TRD cluster tree");
570 fTRDTracker->LoadClusters(trdTree);
571 if (fTRDTracker->PropagateBack(esd) != 0) {
572 AliError("TRD backward propagation failed");
575 if (fCheckPointLevel > 1) WriteESD(esd, "TRD.back");
578 AliWarning("no TOF tracker");
580 // TOF back propagation
581 AliDebug(1, "TOF back propagation");
582 fTOFLoader->LoadDigits("read");
583 TTree* tofTree = fTOFLoader->TreeD();
585 AliError("Can't get the TOF digits tree");
588 fTOFTracker->LoadClusters(tofTree);
589 if (fTOFTracker->PropagateBack(esd) != 0) {
590 AliError("TOF backward propagation failed");
593 if (fCheckPointLevel > 1) WriteESD(esd, "TOF.back");
594 fTOFTracker->UnloadClusters();
595 fTOFLoader->UnloadDigits();
599 AliDebug(1, "TRD inward refit");
600 if (fTRDTracker->RefitInward(esd) != 0) {
601 AliError("TRD inward refit failed");
604 if (fCheckPointLevel > 1) WriteESD(esd, "TRD.refit");
605 fTRDTracker->UnloadClusters();
606 fTRDLoader->UnloadRecPoints();
609 AliInfo("TPC inward refit");
610 if (fTPCTracker->RefitInward(esd) != 0) {
611 AliError("TPC inward refit failed");
614 if (fCheckPointLevel > 1) WriteESD(esd, "TPC.refit");
617 AliInfo("ITS inward refit");
618 if (fITSTracker->RefitInward(esd) != 0) {
619 AliError("ITS inward refit failed");
622 if (fCheckPointLevel > 1) WriteESD(esd, "ITS.refit");
625 fITSTracker->UnloadClusters();
626 fITSLoader->UnloadRecPoints();
629 fTPCTracker->UnloadClusters();
630 fTPCLoader->UnloadRecPoints();
632 AliInfo("execution time:");
633 ToAliInfo(stopwatch.Print());
638 //_____________________________________________________________________________
639 Bool_t AliReconstruction::FillESD(AliESD*& esd, const TString& detectors)
641 // fill the event summary data
643 TStopwatch stopwatch;
645 AliInfo("filling ESD");
647 TString detStr = detectors;
648 for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
649 AliReconstructor* reconstructor =
650 (AliReconstructor*) fReconstructors[iDet];
651 TString detName = reconstructor->GetDetectorName();
652 if (IsSelected(detName, detStr)) {
653 if (!ReadESD(esd, detName.Data())) {
654 AliDebug(1, Form("filling ESD for %s", detName.Data()));
656 reconstructor->FillESD(fRunLoader, fRawReader, esd);
658 reconstructor->FillESD(fRunLoader, esd);
660 if (fCheckPointLevel > 2) WriteESD(esd, detName.Data());
665 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
666 AliError(Form("the following detectors were not found: %s",
668 if (fStopOnError) return kFALSE;
671 AliInfo("execution time:");
672 ToAliInfo(stopwatch.Print());
678 //_____________________________________________________________________________
679 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
681 // check whether detName is contained in detectors
682 // if yes, it is removed from detectors
684 // check if all detectors are selected
685 if ((detectors.CompareTo("ALL") == 0) ||
686 detectors.BeginsWith("ALL ") ||
687 detectors.EndsWith(" ALL") ||
688 detectors.Contains(" ALL ")) {
693 // search for the given detector
694 Bool_t result = kFALSE;
695 if ((detectors.CompareTo(detName) == 0) ||
696 detectors.BeginsWith(detName+" ") ||
697 detectors.EndsWith(" "+detName) ||
698 detectors.Contains(" "+detName+" ")) {
699 detectors.ReplaceAll(detName, "");
703 // clean up the detectors string
704 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
705 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
706 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
711 //_____________________________________________________________________________
712 AliReconstructor* AliReconstruction::GetReconstructor(const char* detName) const
714 // get the reconstructor object for a detector
716 for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
717 AliReconstructor* reconstructor =
718 (AliReconstructor*) fReconstructors[iDet];
719 if (strcmp(reconstructor->GetDetectorName(), detName) == 0) {
720 return reconstructor;
726 //_____________________________________________________________________________
727 Bool_t AliReconstruction::CreateVertexer()
729 // create the vertexer
732 AliReconstructor* itsReconstructor = GetReconstructor("ITS");
733 if (itsReconstructor) {
734 fITSVertexer = itsReconstructor->CreateVertexer(fRunLoader);
737 AliWarning("couldn't create a vertexer for ITS");
738 if (fStopOnError) return kFALSE;
744 //_____________________________________________________________________________
745 Bool_t AliReconstruction::CreateTrackers()
747 // get the loaders and create the trackers
750 fITSLoader = fRunLoader->GetLoader("ITSLoader");
752 AliWarning("no ITS loader found");
753 if (fStopOnError) return kFALSE;
755 AliReconstructor* itsReconstructor = GetReconstructor("ITS");
756 if (itsReconstructor) {
757 fITSTracker = itsReconstructor->CreateTracker(fRunLoader);
760 AliWarning("couldn't create a tracker for ITS");
761 if (fStopOnError) return kFALSE;
766 fTPCLoader = fRunLoader->GetLoader("TPCLoader");
768 AliError("no TPC loader found");
769 if (fStopOnError) return kFALSE;
771 AliReconstructor* tpcReconstructor = GetReconstructor("TPC");
772 if (tpcReconstructor) {
773 fTPCTracker = tpcReconstructor->CreateTracker(fRunLoader);
776 AliError("couldn't create a tracker for TPC");
777 if (fStopOnError) return kFALSE;
782 fTRDLoader = fRunLoader->GetLoader("TRDLoader");
784 AliWarning("no TRD loader found");
785 if (fStopOnError) return kFALSE;
787 AliReconstructor* trdReconstructor = GetReconstructor("TRD");
788 if (trdReconstructor) {
789 fTRDTracker = trdReconstructor->CreateTracker(fRunLoader);
792 AliWarning("couldn't create a tracker for TRD");
793 if (fStopOnError) return kFALSE;
798 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
800 AliWarning("no TOF loader found");
801 if (fStopOnError) return kFALSE;
803 AliReconstructor* tofReconstructor = GetReconstructor("TOF");
804 if (tofReconstructor) {
805 fTOFTracker = tofReconstructor->CreateTracker(fRunLoader);
808 AliWarning("couldn't create a tracker for TOF");
809 if (fStopOnError) return kFALSE;
816 //_____________________________________________________________________________
817 void AliReconstruction::CleanUp(TFile* file)
819 // delete trackers and the run loader and close and delete the file
821 fReconstructors.Delete();
846 //_____________________________________________________________________________
847 Bool_t AliReconstruction::ReadESD(AliESD*& esd, const char* recStep) const
849 // read the ESD event from a file
851 if (!esd) return kFALSE;
853 sprintf(fileName, "ESD_%d.%d_%s.root",
854 esd->GetRunNumber(), esd->GetEventNumber(), recStep);
855 if (gSystem->AccessPathName(fileName)) return kFALSE;
857 AliDebug(1, Form("reading ESD from file %s", fileName));
858 TFile* file = TFile::Open(fileName);
859 if (!file || !file->IsOpen()) {
860 AliError(Form("opening %s failed", fileName));
867 esd = (AliESD*) file->Get("ESD");
873 //_____________________________________________________________________________
874 void AliReconstruction::WriteESD(AliESD* esd, const char* recStep) const
876 // write the ESD event to a file
880 sprintf(fileName, "ESD_%d.%d_%s.root",
881 esd->GetRunNumber(), esd->GetEventNumber(), recStep);
883 AliDebug(1, Form("writing ESD to file %s", fileName));
884 TFile* file = TFile::Open(fileName, "recreate");
885 if (!file || !file->IsOpen()) {
886 AliError(Form("opening %s failed", fileName));