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 reconstruction can be switched on or off for individual detectors by //
37 // rec.SetRunReconstruction("..."); //
39 // The argument is a (case sensitive) string with the names of the //
40 // detectors separated by a space. The special string "ALL" selects all //
41 // available detectors. This is the default. //
43 // The tracking in ITS, TPC and TRD and the creation of ESD tracks can be //
46 // rec.SetRunTracking(kFALSE); //
48 // The filling of additional ESD information can be steered by //
50 // rec.SetFillESD("..."); //
52 // Again, the string specifies the list of detectors. The default is "ALL". //
54 // The reconstruction requires digits as input. For the creation of digits //
55 // have a look at the class AliSimulation. //
57 // For debug purposes the method SetCheckPointLevel can be used. If the //
58 // argument is greater than 0, files with ESD events will be written after //
59 // selected steps of the reconstruction for each event: //
60 // level 1: after tracking and after filling of ESD (final) //
61 // level 2: in addition after each tracking step //
62 // level 3: in addition after the filling of ESD for each detector //
63 // If a final check point file exists for an event, this event will be //
64 // skipped in the reconstruction. The tracking and the filling of ESD for //
65 // a detector will be skipped as well, if the corresponding check point //
66 // file exists. The ESD event will then be loaded from the file instead. //
68 ///////////////////////////////////////////////////////////////////////////////
71 #include "AliReconstruction.h"
72 #include "AliRunLoader.h"
74 #include "AliModule.h"
75 #include "AliDetector.h"
76 #include "AliTracker.h"
78 #include "AliHeader.h"
79 #include "AliGenEventHeader.h"
80 #include "AliESDpid.h"
87 ClassImp(AliReconstruction)
90 //_____________________________________________________________________________
91 AliReconstruction::AliReconstruction(const char* gAliceFilename,
92 const char* name, const char* title) :
95 fRunReconstruction("ALL"),
98 fGAliceFileName(gAliceFilename),
112 // create reconstruction object with default parameters
116 //_____________________________________________________________________________
117 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
120 fRunReconstruction(rec.fRunReconstruction),
121 fRunTracking(rec.fRunTracking),
122 fFillESD(rec.fFillESD),
123 fGAliceFileName(rec.fGAliceFileName),
124 fStopOnError(rec.fStopOnError),
141 //_____________________________________________________________________________
142 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
144 // assignment operator
146 this->~AliReconstruction();
147 new(this) AliReconstruction(rec);
151 //_____________________________________________________________________________
152 AliReconstruction::~AliReconstruction()
160 //_____________________________________________________________________________
161 void AliReconstruction::SetGAliceFile(const char* fileName)
163 // set the name of the galice file
165 fGAliceFileName = fileName;
169 //_____________________________________________________________________________
170 Bool_t AliReconstruction::Run()
172 // run the reconstruction
174 // open the run loader
175 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
177 Error("Run", "no run loader found in file %s",
178 fGAliceFileName.Data());
182 fRunLoader->LoadgAlice();
183 AliRun* aliRun = fRunLoader->GetAliRun();
185 Error("Run", "no gAlice object found in file %s",
186 fGAliceFileName.Data());
192 // local reconstruction
193 if (!fRunReconstruction.IsNull()) {
194 if (!RunReconstruction(fRunReconstruction)) {
195 if (fStopOnError) {CleanUp(); return kFALSE;}
198 if (!fRunTracking && fFillESD.IsNull()) return kTRUE;
200 // get loaders and trackers
201 if (fRunTracking && !CreateTrackers()) {
208 // create the ESD output file
209 TFile* file = TFile::Open("AliESDs.root", "RECREATE");
210 if (!file->IsOpen()) {
211 Error("Run", "opening AliESDs.root failed");
212 if (fStopOnError) {CleanUp(file); return kFALSE;}
216 for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) {
217 Info("Run", "processing event %d", iEvent);
218 fRunLoader->GetEvent(iEvent);
221 sprintf(fileName, "ESD_%d.%d_final.root",
222 aliRun->GetRunNumber(), aliRun->GetEvNumber());
223 if (!gSystem->AccessPathName(fileName)) continue;
225 AliESD* esd = new AliESD;
226 esd->SetRunNumber(aliRun->GetRunNumber());
227 esd->SetEventNumber(aliRun->GetEvNumber());
228 esd->SetMagneticField(aliRun->Field()->SolenoidField());
232 if (!ReadESD(esd, "tracking")) {
233 if (!RunTracking(esd)) {
234 if (fStopOnError) {CleanUp(file); return kFALSE;}
236 if (fCheckPointLevel > 0) WriteESD(esd, "tracking");
241 if (!fFillESD.IsNull()) {
242 if (!FillESD(esd, fFillESD)) {
243 if (fStopOnError) {CleanUp(file); return kFALSE;}
248 AliESDpid::MakePID(esd);
249 if (fCheckPointLevel > 1) WriteESD(esd, "PID");
253 sprintf(name, "ESD%d", iEvent);
255 if (!esd->Write(name)) {
256 Error("Run", "writing ESD failed");
257 if (fStopOnError) {CleanUp(file); return kFALSE;}
261 if (fCheckPointLevel > 0) WriteESD(esd, "final");
271 //_____________________________________________________________________________
272 Bool_t AliReconstruction::RunReconstruction(const TString& detectors)
274 // run the reconstruction
276 TStopwatch stopwatch;
279 TString detStr = detectors;
280 TObjArray* detArray = fRunLoader->GetAliRun()->Detectors();
281 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
282 AliModule* det = (AliModule*) detArray->At(iDet);
283 if (!det || !det->IsActive()) continue;
284 if (IsSelected(det->GetName(), detStr)) {
285 Info("RunReconstruction", "running reconstruction for %s",
287 TStopwatch stopwatchDet;
288 stopwatchDet.Start();
290 Info("RunReconstruction", "execution time for %s:", det->GetName());
291 stopwatchDet.Print();
295 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
296 Error("RunReconstruction", "the following detectors were not found: %s",
298 if (fStopOnError) return kFALSE;
301 Info("RunReconstruction", "execution time:");
307 //_____________________________________________________________________________
308 Bool_t AliReconstruction::RunTracking(AliESD*& esd)
310 // run the barrel tracking
312 TStopwatch stopwatch;
315 // get the primary vertex (from MC for the moment)
317 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(vertex);
318 Double_t vtxPos[3] = {vertex[0], vertex[1], vertex[2]};
319 Double_t vtxCov[6] = {
324 Double_t vtxErr[3] = {vtxCov[0], vtxCov[2], vtxCov[5]}; // diag. elements
325 esd->SetVertex(vtxPos, vtxCov);
326 if (fITSTracker) fITSTracker->SetVertex(vtxPos, vtxErr);
327 if (fTPCTracker) fTPCTracker->SetVertex(vtxPos, vtxErr);
328 if (fTRDTracker) fTRDTracker->SetVertex(vtxPos, vtxErr);
329 if (fCheckPointLevel > 1) WriteESD(esd, "vertex");
332 Error("RunTracking", "no TPC tracker");
337 Info("RunTracking", "TPC tracking");
338 fTPCLoader->LoadRecPoints("read");
339 TTree* tpcTree = fTPCLoader->TreeR();
341 Error("RunTracking", "Can't get the TPC cluster tree");
344 fTPCTracker->LoadClusters(tpcTree);
345 if (fTPCTracker->Clusters2Tracks(esd) != 0) {
346 Error("RunTracking", "TPC Clusters2Tracks failed");
349 if (fCheckPointLevel > 1) WriteESD(esd, "TPC.tracking");
352 Warning("RunTracking", "no ITS tracker");
355 fRunLoader->GetAliRun()->GetDetector("TPC")->FillESD(esd); // preliminary
356 AliESDpid::MakePID(esd); // PID for the ITS tracker
359 Info("RunTracking", "ITS tracking");
360 fITSLoader->LoadRecPoints("read");
361 TTree* itsTree = fITSLoader->TreeR();
363 Error("RunTracking", "Can't get the ITS cluster tree");
366 fITSTracker->LoadClusters(itsTree);
367 if (fITSTracker->Clusters2Tracks(esd) != 0) {
368 Error("RunTracking", "ITS Clusters2Tracks failed");
371 if (fCheckPointLevel > 1) WriteESD(esd, "ITS.tracking");
374 Warning("RunTracking", "no TRD tracker");
376 // ITS back propagation
377 Info("RunTracking", "ITS back propagation");
378 if (fITSTracker->PropagateBack(esd) != 0) {
379 Error("RunTracking", "ITS backward propagation failed");
382 if (fCheckPointLevel > 1) WriteESD(esd, "ITS.back");
384 // TPC back propagation
385 Info("RunTracking", "TPC back propagation");
386 if (fTPCTracker->PropagateBack(esd) != 0) {
387 Error("RunTracking", "TPC backward propagation failed");
390 if (fCheckPointLevel > 1) WriteESD(esd, "TPC.back");
392 // TRD back propagation
393 Info("RunTracking", "TRD back propagation");
394 fTRDLoader->LoadRecPoints("read");
395 TTree* trdTree = fTRDLoader->TreeR();
397 Error("RunTracking", "Can't get the TRD cluster tree");
400 fTRDTracker->LoadClusters(trdTree);
401 if (fTRDTracker->PropagateBack(esd) != 0) {
402 Error("RunTracking", "TRD backward propagation failed");
405 if (fCheckPointLevel > 1) WriteESD(esd, "TRD.back");
408 Warning("RunTracking", "no TOF tracker");
410 // TOF back propagation
411 Info("RunTracking", "TOF back propagation");
412 fTOFLoader->LoadDigits("read");
413 TTree* tofTree = fTOFLoader->TreeD();
415 Error("RunTracking", "Can't get the TOF digits tree");
418 fTOFTracker->LoadClusters(tofTree);
419 if (fTOFTracker->PropagateBack(esd) != 0) {
420 Error("RunTracking", "TOF backward propagation failed");
423 if (fCheckPointLevel > 1) WriteESD(esd, "TOF.back");
424 fTOFTracker->UnloadClusters();
425 fTOFLoader->UnloadDigits();
429 Info("RunTracking", "TRD inward refit");
430 if (fTRDTracker->RefitInward(esd) != 0) {
431 Error("RunTracking", "TRD inward refit failed");
434 if (fCheckPointLevel > 1) WriteESD(esd, "TRD.refit");
435 fTRDTracker->UnloadClusters();
436 fTRDLoader->UnloadRecPoints();
439 Info("RunTracking", "TPC inward refit");
440 if (fTPCTracker->RefitInward(esd) != 0) {
441 Error("RunTracking", "TPC inward refit failed");
444 if (fCheckPointLevel > 1) WriteESD(esd, "TPC.refit");
447 Info("RunTracking", "ITS inward refit");
448 if (fITSTracker->RefitInward(esd) != 0) {
449 Error("RunTracking", "ITS inward refit failed");
452 if (fCheckPointLevel > 1) WriteESD(esd, "ITS.refit");
455 fITSTracker->UnloadClusters();
456 fITSLoader->UnloadRecPoints();
459 fTPCTracker->UnloadClusters();
460 fTPCLoader->UnloadRecPoints();
462 Info("RunTracking", "execution time:");
468 //_____________________________________________________________________________
469 Bool_t AliReconstruction::FillESD(AliESD*& esd, const TString& detectors)
471 // fill the event summary data
473 TStopwatch stopwatch;
476 TString detStr = detectors;
477 TObjArray* detArray = fRunLoader->GetAliRun()->Detectors();
478 for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
479 AliModule* det = (AliModule*) detArray->At(iDet);
480 if (!det || !det->IsActive()) continue;
481 if (IsSelected(det->GetName(), detStr)) {
482 if (!ReadESD(esd, det->GetName())) {
483 Info("FillESD", "filling ESD for %s",
486 if (fCheckPointLevel > 2) WriteESD(esd, det->GetName());
491 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
492 Error("FillESD", "the following detectors were not found: %s",
494 if (fStopOnError) return kFALSE;
497 Info("FillESD", "execution time:");
504 //_____________________________________________________________________________
505 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
507 // check whether detName is contained in detectors
508 // if yes, it is removed from detectors
510 // check if all detectors are selected
511 if ((detectors.CompareTo("ALL") == 0) ||
512 detectors.BeginsWith("ALL ") ||
513 detectors.EndsWith(" ALL") ||
514 detectors.Contains(" ALL ")) {
519 // search for the given detector
520 Bool_t result = kFALSE;
521 if ((detectors.CompareTo(detName) == 0) ||
522 detectors.BeginsWith(detName+" ") ||
523 detectors.EndsWith(" "+detName) ||
524 detectors.Contains(" "+detName+" ")) {
525 detectors.ReplaceAll(detName, "");
529 // clean up the detectors string
530 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
531 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
532 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
537 //_____________________________________________________________________________
538 Bool_t AliReconstruction::CreateTrackers()
540 // get the loaders and create the trackers
542 AliRun* aliRun = fRunLoader->GetAliRun();
545 fITSLoader = fRunLoader->GetLoader("ITSLoader");
547 Warning("CreateTrackers", "no ITS loader found");
548 if (fStopOnError) return kFALSE;
550 if (aliRun->GetDetector("ITS")) {
551 fITSTracker = aliRun->GetDetector("ITS")->CreateTracker();
554 Warning("CreateTrackers", "couldn't create a tracker for ITS");
555 if (fStopOnError) return kFALSE;
560 fTPCLoader = fRunLoader->GetLoader("TPCLoader");
562 Error("CreateTrackers", "no TPC loader found");
563 if (fStopOnError) return kFALSE;
565 if (aliRun->GetDetector("TPC")) {
566 fTPCTracker = aliRun->GetDetector("TPC")->CreateTracker();
569 Error("CreateTrackers", "couldn't create a tracker for TPC");
570 if (fStopOnError) return kFALSE;
575 fTRDLoader = fRunLoader->GetLoader("TRDLoader");
577 Warning("CreateTrackers", "no TRD loader found");
578 if (fStopOnError) return kFALSE;
580 if (aliRun->GetDetector("TRD")) {
581 fTRDTracker = aliRun->GetDetector("TRD")->CreateTracker();
584 Warning("CreateTrackers", "couldn't create a tracker for TRD");
585 if (fStopOnError) return kFALSE;
590 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
592 Warning("CreateTrackers", "no TOF loader found");
593 if (fStopOnError) return kFALSE;
595 if (aliRun->GetDetector("TOF")) {
596 fTOFTracker = aliRun->GetDetector("TOF")->CreateTracker();
599 Warning("CreateTrackers", "couldn't create a tracker for TOF");
600 if (fStopOnError) return kFALSE;
607 //_____________________________________________________________________________
608 void AliReconstruction::CleanUp(TFile* file)
610 // delete trackers and the run loader and close and delete the file
631 //_____________________________________________________________________________
632 Bool_t AliReconstruction::ReadESD(AliESD*& esd, const char* recStep) const
634 // read the ESD event from a file
636 if (!esd) return kFALSE;
638 sprintf(fileName, "ESD_%d.%d_%s.root",
639 esd->GetRunNumber(), esd->GetEventNumber(), recStep);
640 if (gSystem->AccessPathName(fileName)) return kFALSE;
642 Info("ReadESD", "reading ESD from file %s", fileName);
643 TFile* file = TFile::Open(fileName);
644 if (!file || !file->IsOpen()) {
645 Error("ReadESD", "opening %s failed", fileName);
652 esd = (AliESD*) file->Get("ESD");
658 //_____________________________________________________________________________
659 void AliReconstruction::WriteESD(AliESD* esd, const char* recStep) const
661 // write the ESD event to a file
665 sprintf(fileName, "ESD_%d.%d_%s.root",
666 esd->GetRunNumber(), esd->GetEventNumber(), recStep);
668 Info("WriteESD", "writing ESD to file %s", fileName);
669 TFile* file = TFile::Open(fileName, "recreate");
670 if (!file || !file->IsOpen()) {
671 Error("WriteESD", "opening %s failed", fileName);