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),
147 // create reconstruction object with default parameters
151 //_____________________________________________________________________________
152 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
155 fRunLocalReconstruction(rec.fRunLocalReconstruction),
156 fRunVertexFinder(rec.fRunVertexFinder),
157 fRunTracking(rec.fRunTracking),
158 fFillESD(rec.fFillESD),
159 fGAliceFileName(rec.fGAliceFileName),
161 fStopOnError(rec.fStopOnError),
183 for (Int_t i = 0; i < fOptions.GetEntriesFast(); i++) {
184 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
188 //_____________________________________________________________________________
189 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
191 // assignment operator
193 this->~AliReconstruction();
194 new(this) AliReconstruction(rec);
198 //_____________________________________________________________________________
199 AliReconstruction::~AliReconstruction()
208 //_____________________________________________________________________________
209 void AliReconstruction::SetGAliceFile(const char* fileName)
211 // set the name of the galice file
213 fGAliceFileName = fileName;
216 //_____________________________________________________________________________
217 void AliReconstruction::SetOption(const char* detector, const char* option)
219 // set options for the reconstruction of a detector
221 TObject* obj = fOptions.FindObject(detector);
222 if (obj) fOptions.Remove(obj);
223 fOptions.Add(new TNamed(detector, option));
227 //_____________________________________________________________________________
228 Bool_t AliReconstruction::Run(const char* input)
230 // run the reconstruction
233 if (!input) input = fInput.Data();
234 TString fileName(input);
235 if (fileName.EndsWith("/")) {
236 fRawReader = new AliRawReaderFile(fileName);
237 } else if (fileName.EndsWith(".root")) {
238 fRawReader = new AliRawReaderRoot(fileName);
239 } else if (!fileName.IsNull()) {
240 fRawReader = new AliRawReaderDate(fileName);
241 fRawReader->SelectEvents(7);
244 // open the run loader
245 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
247 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
251 fRunLoader->LoadgAlice();
252 AliRun* aliRun = fRunLoader->GetAliRun();
254 AliError(Form("no gAlice object found in file %s",
255 fGAliceFileName.Data()));
260 AliTracker::SetFieldMap(gAlice->Field());
262 // load the reconstructor objects
263 TPluginManager* pluginManager = gROOT->GetPluginManager();
264 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
265 TString detName = fgkDetectorName[iDet];
266 TString recName = "Ali" + detName + "Reconstructor";
267 if (!gAlice->GetDetector(detName) && detName != "HLT") continue;
269 if(detName == "HLT") {
270 if (!gROOT->GetClass("AliLevel3")) {
271 gSystem->Load("libAliL3Src.so");
272 gSystem->Load("libAliL3Misc.so");
273 gSystem->Load("libAliL3Hough.so");
274 gSystem->Load("libAliL3Comp.so");
278 AliReconstructor* reconstructor = NULL;
279 // first check if a plugin is defined for the reconstructor
280 TPluginHandler* pluginHandler =
281 pluginManager->FindHandler("AliReconstructor", detName);
282 // if not, but the reconstructor class is implemented, add a plugin for it
283 if (!pluginHandler && gROOT->GetClass(recName.Data())) {
284 AliDebug(1, Form("defining plugin for %s", recName.Data()));
285 pluginManager->AddHandler("AliReconstructor", detName,
286 recName, detName, recName + "()");
287 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
289 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
290 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
292 // if there is no reconstructor class for the detector use the dummy one
293 if (!reconstructor && gAlice->GetDetector(detName)) {
294 AliDebug(1, Form("using dummy reconstructor for %s", detName.Data()));
295 reconstructor = new AliDummyReconstructor(gAlice->GetDetector(detName));
298 TObject* obj = fOptions.FindObject(detName.Data());
299 if (obj) reconstructor->SetOption(obj->GetTitle());
300 fReconstructors.Add(reconstructor);
304 // local reconstruction
305 if (!fRunLocalReconstruction.IsNull()) {
306 if (!RunLocalReconstruction(fRunLocalReconstruction)) {
307 if (fStopOnError) {CleanUp(); return kFALSE;}
310 if (!fRunVertexFinder && !fRunTracking && fFillESD.IsNull()) return kTRUE;
313 if (fRunVertexFinder && !CreateVertexer()) {
320 // get loaders and trackers
321 if (fRunTracking && !CreateTrackers()) {
328 // create the ESD output file and tree
329 TFile* file = TFile::Open("AliESDs.root", "RECREATE");
330 if (!file->IsOpen()) {
331 AliError("opening AliESDs.root failed");
332 if (fStopOnError) {CleanUp(file); return kFALSE;}
334 AliESD* esd = new AliESD;
335 TTree* tree = new TTree("esdTree", "Tree with ESD objects");
336 tree->Branch("ESD", "AliESD", &esd);
341 if (fRawReader) fRawReader->RewindEvents();
342 for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) {
343 AliInfo(Form("processing event %d", iEvent));
344 fRunLoader->GetEvent(iEvent);
345 if (fRawReader) fRawReader->NextEvent();
348 sprintf(fileName, "ESD_%d.%d_final.root",
349 aliRun->GetRunNumber(), aliRun->GetEvNumber());
350 if (!gSystem->AccessPathName(fileName)) continue;
353 esd->SetRunNumber(aliRun->GetRunNumber());
354 esd->SetEventNumber(aliRun->GetEvNumber());
355 esd->SetMagneticField(aliRun->Field()->SolenoidField());
358 if (fRunVertexFinder) {
359 if (!ReadESD(esd, "vertex")) {
360 if (!RunVertexFinder(esd)) {
361 if (fStopOnError) {CleanUp(file); return kFALSE;}
363 if (fCheckPointLevel > 0) WriteESD(esd, "vertex");
369 if (!ReadESD(esd, "tracking")) {
370 if (!RunTracking(esd)) {
371 if (fStopOnError) {CleanUp(file); return kFALSE;}
373 if (fCheckPointLevel > 0) WriteESD(esd, "tracking");
378 if (!fFillESD.IsNull()) {
379 if (!FillESD(esd, fFillESD)) {
380 if (fStopOnError) {CleanUp(file); return kFALSE;}
385 AliESDpid::MakePID(esd);
386 if (fCheckPointLevel > 1) WriteESD(esd, "PID");
391 if (fCheckPointLevel > 0) WriteESD(esd, "final");
403 //_____________________________________________________________________________
404 Bool_t AliReconstruction::RunLocalReconstruction(const TString& detectors)
406 // run the local reconstruction
408 TStopwatch stopwatch;
411 TString detStr = detectors;
412 for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
413 AliReconstructor* reconstructor =
414 (AliReconstructor*) fReconstructors[iDet];
415 TString detName = reconstructor->GetDetectorName();
416 if (IsSelected(detName, detStr)) {
417 AliInfo(Form("running reconstruction for %s", detName.Data()));
418 TStopwatch stopwatchDet;
419 stopwatchDet.Start();
421 fRawReader->RewindEvents();
422 reconstructor->Reconstruct(fRunLoader, fRawReader);
424 reconstructor->Reconstruct(fRunLoader);
426 AliInfo(Form("execution time for %s:", detName.Data()));
427 ToAliInfo(stopwatchDet.Print());
431 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
432 AliError(Form("the following detectors were not found: %s",
434 if (fStopOnError) return kFALSE;
437 AliInfo("execution time:");
438 ToAliInfo(stopwatch.Print());
443 //_____________________________________________________________________________
444 Bool_t AliReconstruction::RunVertexFinder(AliESD*& esd)
446 // run the barrel tracking
448 TStopwatch stopwatch;
451 AliESDVertex* vertex = NULL;
452 Double_t vtxPos[3] = {0, 0, 0};
453 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
455 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
456 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
457 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
461 AliInfo("running the ITS vertex finder");
462 vertex = fITSVertexer->FindVertexForCurrentEvent(fRunLoader->GetEventNumber());
464 AliWarning("Vertex not found");
465 vertex = new AliESDVertex();
468 vertex->SetTruePos(vtxPos); // store also the vertex from MC
472 AliInfo("getting the primary vertex from MC");
473 vertex = new AliESDVertex(vtxPos, vtxErr);
477 vertex->GetXYZ(vtxPos);
478 vertex->GetSigmaXYZ(vtxErr);
480 AliWarning("no vertex reconstructed");
481 vertex = new AliESDVertex(vtxPos, vtxErr);
483 esd->SetVertex(vertex);
484 if (fITSTracker) fITSTracker->SetVertex(vtxPos, vtxErr);
485 if (fTPCTracker) fTPCTracker->SetVertex(vtxPos, vtxErr);
486 if (fTRDTracker) fTRDTracker->SetVertex(vtxPos, vtxErr);
489 AliInfo("execution time:");
490 ToAliInfo(stopwatch.Print());
495 //_____________________________________________________________________________
496 Bool_t AliReconstruction::RunTracking(AliESD*& esd)
498 // run the barrel tracking
500 TStopwatch stopwatch;
504 AliError("no TPC tracker");
507 AliInfo("running tracking");
510 AliDebug(1, "TPC tracking");
511 fTPCLoader->LoadRecPoints("read");
512 TTree* tpcTree = fTPCLoader->TreeR();
514 AliError("Can't get the TPC cluster tree");
517 fTPCTracker->LoadClusters(tpcTree);
518 if (fTPCTracker->Clusters2Tracks(esd) != 0) {
519 AliError("TPC Clusters2Tracks failed");
522 if (fCheckPointLevel > 1) WriteESD(esd, "TPC.tracking");
525 AliWarning("no ITS tracker");
528 GetReconstructor("TPC")->FillESD(fRunLoader, esd); // preliminary
529 AliESDpid::MakePID(esd); // PID for the ITS tracker
532 AliDebug(1, "ITS tracking");
533 fITSLoader->LoadRecPoints("read");
534 TTree* itsTree = fITSLoader->TreeR();
536 Error("RunTracking", "Can't get the ITS cluster tree");
539 fITSTracker->LoadClusters(itsTree);
540 if (fITSTracker->Clusters2Tracks(esd) != 0) {
541 AliError("ITS Clusters2Tracks failed");
544 if (fCheckPointLevel > 1) WriteESD(esd, "ITS.tracking");
546 if (fTRDTracker || fTOFTracker) {
547 // ITS back propagation
548 AliDebug(1, "ITS back propagation");
549 if (fITSTracker->PropagateBack(esd) != 0) {
550 AliError("ITS backward propagation failed");
553 if (fCheckPointLevel > 1) WriteESD(esd, "ITS.back");
555 // TPC back propagation
556 AliDebug(1, "TPC back propagation");
557 if (fTPCTracker->PropagateBack(esd) != 0) {
558 AliError("TPC backward propagation failed");
561 if (fCheckPointLevel > 1) WriteESD(esd, "TPC.back");
564 AliWarning("no TRD tracker");
566 // TRD back propagation
567 AliDebug(1, "TRD back propagation");
568 fTRDLoader->LoadRecPoints("read");
569 TTree* trdTree = fTRDLoader->TreeR();
571 AliError("Can't get the TRD cluster tree");
574 fTRDTracker->LoadClusters(trdTree);
575 if (fTRDTracker->PropagateBack(esd) != 0) {
576 AliError("TRD backward propagation failed");
579 if (fCheckPointLevel > 1) WriteESD(esd, "TRD.back");
583 AliWarning("no TOF tracker");
585 // TOF back propagation
586 AliDebug(1, "TOF back propagation");
587 fTOFLoader->LoadDigits("read");
588 TTree* tofTree = fTOFLoader->TreeD();
590 AliError("Can't get the TOF digits tree");
593 fTOFTracker->LoadClusters(tofTree);
594 if (fTOFTracker->PropagateBack(esd) != 0) {
595 AliError("TOF backward propagation failed");
598 if (fCheckPointLevel > 1) WriteESD(esd, "TOF.back");
599 fTOFTracker->UnloadClusters();
600 fTOFLoader->UnloadDigits();
603 AliWarning("no RICH tracker");
605 // RICH back propagation
606 AliDebug(1, "RICH back propagation");
607 fRICHLoader->LoadRecPoints("read");
608 TTree* richTree = fRICHLoader->TreeR();
610 AliError("Can't get the RICH cluster tree");
613 fRICHTracker->LoadClusters(richTree);
614 if (fRICHTracker->PropagateBack(esd) != 0) {
615 AliError("RICH backward propagation failed");
618 if (fCheckPointLevel > 1) WriteESD(esd, "RICH.back");
619 fRICHTracker->UnloadClusters();
620 fRICHLoader->UnloadRecPoints();
626 AliDebug(1, "TRD inward refit");
627 if (fTRDTracker->RefitInward(esd) != 0) {
628 AliError("TRD inward refit failed");
631 if (fCheckPointLevel > 1) WriteESD(esd, "TRD.refit");
632 fTRDTracker->UnloadClusters();
633 fTRDLoader->UnloadRecPoints();
637 AliInfo("TPC inward refit");
638 if (fTPCTracker->RefitInward(esd) != 0) {
639 AliError("TPC inward refit failed");
642 if (fCheckPointLevel > 1) WriteESD(esd, "TPC.refit");
645 AliInfo("ITS inward refit");
646 if (fITSTracker->RefitInward(esd) != 0) {
647 AliError("ITS inward refit failed");
650 if (fCheckPointLevel > 1) WriteESD(esd, "ITS.refit");
653 fITSTracker->UnloadClusters();
654 fITSLoader->UnloadRecPoints();
657 fTPCTracker->UnloadClusters();
658 fTPCLoader->UnloadRecPoints();
660 AliInfo("execution time:");
661 ToAliInfo(stopwatch.Print());
666 //_____________________________________________________________________________
667 Bool_t AliReconstruction::FillESD(AliESD*& esd, const TString& detectors)
669 // fill the event summary data
671 TStopwatch stopwatch;
673 AliInfo("filling ESD");
675 TString detStr = detectors;
676 for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
677 AliReconstructor* reconstructor =
678 (AliReconstructor*) fReconstructors[iDet];
679 TString detName = reconstructor->GetDetectorName();
680 if (IsSelected(detName, detStr)) {
681 if (!ReadESD(esd, detName.Data())) {
682 AliDebug(1, Form("filling ESD for %s", detName.Data()));
684 reconstructor->FillESD(fRunLoader, fRawReader, esd);
686 reconstructor->FillESD(fRunLoader, esd);
688 if (fCheckPointLevel > 2) WriteESD(esd, detName.Data());
693 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
694 AliError(Form("the following detectors were not found: %s",
696 if (fStopOnError) return kFALSE;
699 AliInfo("execution time:");
700 ToAliInfo(stopwatch.Print());
706 //_____________________________________________________________________________
707 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
709 // check whether detName is contained in detectors
710 // if yes, it is removed from detectors
712 // check if all detectors are selected
713 if ((detectors.CompareTo("ALL") == 0) ||
714 detectors.BeginsWith("ALL ") ||
715 detectors.EndsWith(" ALL") ||
716 detectors.Contains(" ALL ")) {
721 // search for the given detector
722 Bool_t result = kFALSE;
723 if ((detectors.CompareTo(detName) == 0) ||
724 detectors.BeginsWith(detName+" ") ||
725 detectors.EndsWith(" "+detName) ||
726 detectors.Contains(" "+detName+" ")) {
727 detectors.ReplaceAll(detName, "");
731 // clean up the detectors string
732 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
733 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
734 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
739 //_____________________________________________________________________________
740 AliReconstructor* AliReconstruction::GetReconstructor(const char* detName) const
742 // get the reconstructor object for a detector
744 for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
745 AliReconstructor* reconstructor =
746 (AliReconstructor*) fReconstructors[iDet];
747 if (strcmp(reconstructor->GetDetectorName(), detName) == 0) {
748 return reconstructor;
754 //_____________________________________________________________________________
755 Bool_t AliReconstruction::CreateVertexer()
757 // create the vertexer
760 AliReconstructor* itsReconstructor = GetReconstructor("ITS");
761 if (itsReconstructor) {
762 fITSVertexer = itsReconstructor->CreateVertexer(fRunLoader);
765 AliWarning("couldn't create a vertexer for ITS");
766 if (fStopOnError) return kFALSE;
772 //_____________________________________________________________________________
773 Bool_t AliReconstruction::CreateTrackers()
775 // get the loaders and create the trackers
778 fITSLoader = fRunLoader->GetLoader("ITSLoader");
780 AliWarning("no ITS loader found");
781 if (fStopOnError) return kFALSE;
783 AliReconstructor* itsReconstructor = GetReconstructor("ITS");
784 if (itsReconstructor) {
785 fITSTracker = itsReconstructor->CreateTracker(fRunLoader);
788 AliWarning("couldn't create a tracker for ITS");
789 if (fStopOnError) return kFALSE;
794 fTPCLoader = fRunLoader->GetLoader("TPCLoader");
796 AliError("no TPC loader found");
797 if (fStopOnError) return kFALSE;
799 AliReconstructor* tpcReconstructor = GetReconstructor("TPC");
800 if (tpcReconstructor) {
801 fTPCTracker = tpcReconstructor->CreateTracker(fRunLoader);
804 AliError("couldn't create a tracker for TPC");
805 if (fStopOnError) return kFALSE;
810 fTRDLoader = fRunLoader->GetLoader("TRDLoader");
812 AliWarning("no TRD loader found");
813 if (fStopOnError) return kFALSE;
815 AliReconstructor* trdReconstructor = GetReconstructor("TRD");
816 if (trdReconstructor) {
817 fTRDTracker = trdReconstructor->CreateTracker(fRunLoader);
820 AliWarning("couldn't create a tracker for TRD");
821 if (fStopOnError) return kFALSE;
826 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
828 AliWarning("no TOF loader found");
829 if (fStopOnError) return kFALSE;
831 AliReconstructor* tofReconstructor = GetReconstructor("TOF");
832 if (tofReconstructor) {
833 fTOFTracker = tofReconstructor->CreateTracker(fRunLoader);
836 AliWarning("couldn't create a tracker for TOF");
837 if (fStopOnError) return kFALSE;
842 fRICHLoader = fRunLoader->GetLoader("RICHLoader");
844 AliWarning("no RICH loader found");
845 if (fStopOnError) return kFALSE;
847 AliReconstructor* tofReconstructor = GetReconstructor("RICH");
848 if (tofReconstructor) {
849 fRICHTracker = tofReconstructor->CreateTracker(fRunLoader);
852 AliWarning("couldn't create a tracker for RICH");
853 if (fStopOnError) return kFALSE;
860 //_____________________________________________________________________________
861 void AliReconstruction::CleanUp(TFile* file)
863 // delete trackers and the run loader and close and delete the file
865 fReconstructors.Delete();
892 //_____________________________________________________________________________
893 Bool_t AliReconstruction::ReadESD(AliESD*& esd, const char* recStep) const
895 // read the ESD event from a file
897 if (!esd) return kFALSE;
899 sprintf(fileName, "ESD_%d.%d_%s.root",
900 esd->GetRunNumber(), esd->GetEventNumber(), recStep);
901 if (gSystem->AccessPathName(fileName)) return kFALSE;
903 AliDebug(1, Form("reading ESD from file %s", fileName));
904 TFile* file = TFile::Open(fileName);
905 if (!file || !file->IsOpen()) {
906 AliError(Form("opening %s failed", fileName));
913 esd = (AliESD*) file->Get("ESD");
919 //_____________________________________________________________________________
920 void AliReconstruction::WriteESD(AliESD* esd, const char* recStep) const
922 // write the ESD event to a file
926 sprintf(fileName, "ESD_%d.%d_%s.root",
927 esd->GetRunNumber(), esd->GetEventNumber(), recStep);
929 AliDebug(1, Form("writing ESD to file %s", fileName));
930 TFile* file = TFile::Open(fileName, "recreate");
931 if (!file || !file->IsOpen()) {
932 AliError(Form("opening %s failed", fileName));