]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - STEER/AliReconstruction.cxx
framework for reconstruction of raw data
[u/mrichter/AliRoot.git] / STEER / AliReconstruction.cxx
index 23683df8231a4abcce551af76abef4cece07e620..47f1f28510b83690fd116e594dd64e79bfc58d00 100644 (file)
 //                                                                           //
 ///////////////////////////////////////////////////////////////////////////////
 
+#include <TArrayF.h>
+#include <TFile.h>
+#include <TSystem.h>
+#include <TROOT.h>
+#include <TPluginManager.h>
+#include <TStopwatch.h>
 
 #include "AliReconstruction.h"
 #include "AliRunLoader.h"
 #include "AliRun.h"
-#include "AliModule.h"
-#include "AliDetector.h"
-#include "AliReconstructor.h"
+#include "AliRawReaderFile.h"
+#include "AliRawReaderDate.h"
+#include "AliRawReaderRoot.h"
 #include "AliTracker.h"
 #include "AliESD.h"
 #include "AliESDVertex.h"
 #include "AliGenEventHeader.h"
 #include "AliESDpid.h"
 #include "AliMagF.h"
-#include <TArrayF.h>
-#include <TSystem.h>
-#include <TROOT.h>
-#include <TPluginManager.h>
-
 
 ClassImp(AliReconstruction)
 
 
+//_____________________________________________________________________________
+const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "RICH", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "START", "VZERO", "CRT", "HLT"};
+
 //_____________________________________________________________________________
 AliReconstruction::AliReconstruction(const char* gAliceFilename,
                                     const char* name, const char* title) :
@@ -102,10 +106,12 @@ AliReconstruction::AliReconstruction(const char* gAliceFilename,
   fRunTracking(kTRUE),
   fFillESD("ALL"),
   fGAliceFileName(gAliceFilename),
+  fInput(""),
   fStopOnError(kFALSE),
   fCheckPointLevel(0),
 
   fRunLoader(NULL),
+  fRawReader(NULL),
   fITSLoader(NULL),
   fITSVertexer(NULL),
   fITSTracker(NULL),
@@ -114,7 +120,10 @@ AliReconstruction::AliReconstruction(const char* gAliceFilename,
   fTRDLoader(NULL),
   fTRDTracker(NULL),
   fTOFLoader(NULL),
-  fTOFTracker(NULL)
+  fTOFTracker(NULL),
+
+  fReconstructors(),
+  fOptions()
 {
 // create reconstruction object with default parameters
 
@@ -129,10 +138,12 @@ AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
   fRunTracking(rec.fRunTracking),
   fFillESD(rec.fFillESD),
   fGAliceFileName(rec.fGAliceFileName),
+  fInput(rec.fInput),
   fStopOnError(rec.fStopOnError),
   fCheckPointLevel(0),
 
   fRunLoader(NULL),
+  fRawReader(NULL),
   fITSLoader(NULL),
   fITSVertexer(NULL),
   fITSTracker(NULL),
@@ -141,10 +152,16 @@ AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
   fTRDLoader(NULL),
   fTRDTracker(NULL),
   fTOFLoader(NULL),
-  fTOFTracker(NULL)
+  fTOFTracker(NULL),
+
+  fReconstructors(),
+  fOptions()
 {
 // copy constructor
 
+  for (Int_t i = 0; i < fOptions.GetEntriesFast(); i++) {
+    if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
+  }
 }
 
 //_____________________________________________________________________________
@@ -163,6 +180,7 @@ AliReconstruction::~AliReconstruction()
 // clean up
 
   CleanUp();
+  fOptions.Delete();
 }
 
 
@@ -174,12 +192,34 @@ void AliReconstruction::SetGAliceFile(const char* fileName)
   fGAliceFileName = fileName;
 }
 
+//_____________________________________________________________________________
+void AliReconstruction::SetOption(const char* detector, const char* option)
+{
+// set options for the reconstruction of a detector
+
+  TObject* obj = fOptions.FindObject(detector);
+  if (obj) fOptions.Remove(obj);
+  fOptions.Add(new TNamed(detector, option));
+}
+
 
 //_____________________________________________________________________________
-Bool_t AliReconstruction::Run()
+Bool_t AliReconstruction::Run(const char* input)
 {
 // run the reconstruction
 
+  // set the input
+  if (!input) input = fInput.Data();
+  TString fileName(input);
+  if (fileName.EndsWith("/")) {
+    fRawReader = new AliRawReaderFile(fileName);
+  } else if (fileName.EndsWith(".root")) {
+    fRawReader = new AliRawReaderRoot(fileName);
+  } else if (!fileName.IsNull()) {
+    fRawReader = new AliRawReaderDate(fileName);
+    fRawReader->SelectEvents(7);
+  }
+
   // open the run loader
   fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
   if (!fRunLoader) {
@@ -200,22 +240,44 @@ Bool_t AliReconstruction::Run()
 
   // load the reconstructor objects
   TPluginManager* pluginManager = gROOT->GetPluginManager();
-  if (!pluginManager->FindHandler("AliReconstructor", "TPC")) {
-    pluginManager->AddHandler("AliReconstructor", "TPC", 
-                             "AliTPCReconstructor", "TPC", 
-                             "AliTPCReconstructor()");
-  }
-  TObjArray* detArray = fRunLoader->GetAliRun()->Detectors();
-  for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
-    AliModule* det = (AliModule*) detArray->At(iDet);
-    if (!det || !det->IsActive()) continue;
+  for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
+    TString detName = fgkDetectorName[iDet];
+    TString recName = "Ali" + detName + "Reconstructor";
+    if (!gAlice->GetDetector(detName) && detName != "HLT") continue;
+
+    if(detName == "HLT") {
+      if (!gROOT->GetClass("AliLevel3")) {
+       gSystem->Load("libAliL3Src.so");
+       gSystem->Load("libAliL3Misc.so");
+       gSystem->Load("libAliL3Hough.so");
+       gSystem->Load("libAliL3Comp.so");
+      }
+    }
+
+    AliReconstructor* reconstructor = NULL;
+    // first check if a plugin is defined for the reconstructor
     TPluginHandler* pluginHandler = 
-      pluginManager->FindHandler("AliReconstructor", det->GetName());
-    if (!pluginHandler) continue;
-    if (pluginHandler->LoadPlugin() != 0) continue;
-    AliReconstructor* reconstructor = 
-      (AliReconstructor*) pluginHandler->ExecPlugin(0);
-    if (reconstructor) fReconstructors.Add(reconstructor);
+      pluginManager->FindHandler("AliReconstructor", detName);
+    // if not, but the reconstructor class is implemented, add a plugin for it
+    if (!pluginHandler && gROOT->GetClass(recName.Data())) {
+      Info("Run", "defining plugin for %s", recName.Data());
+      pluginManager->AddHandler("AliReconstructor", detName, 
+                               recName, detName, recName + "()");
+      pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
+    }
+    if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
+      reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
+    }
+    // if there is no reconstructor class for the detector use the dummy one
+    if (!reconstructor && gAlice->GetDetector(detName)) {
+      Info("Run", "using dummy reconstructor for %s", detName.Data());
+      reconstructor = new AliDummyReconstructor(gAlice->GetDetector(detName));
+    }
+    if (reconstructor) {
+      TObject* obj = fOptions.FindObject(detName.Data());
+      if (obj) reconstructor->SetOption(obj->GetTitle());
+      fReconstructors.Add(reconstructor);
+    }
   }
 
   // local reconstruction
@@ -242,24 +304,31 @@ Bool_t AliReconstruction::Run()
     }      
   }
 
-  // create the ESD output file
+  // create the ESD output file and tree
   TFile* file = TFile::Open("AliESDs.root", "RECREATE");
   if (!file->IsOpen()) {
     Error("Run", "opening AliESDs.root failed");
     if (fStopOnError) {CleanUp(file); return kFALSE;}    
   }
+  AliESD* esd = new AliESD;
+  TTree* tree = new TTree("esdTree", "Tree with ESD objects");
+  tree->Branch("ESD", "AliESD", &esd);
+  delete esd;
+  gROOT->cd();
 
   // loop over events
+  if (fRawReader) fRawReader->RewindEvents();
   for (Int_t iEvent = 0; iEvent < fRunLoader->GetNumberOfEvents(); iEvent++) {
     Info("Run", "processing event %d", iEvent);
     fRunLoader->GetEvent(iEvent);
+    if (fRawReader) fRawReader->NextEvent();
 
     char fileName[256];
     sprintf(fileName, "ESD_%d.%d_final.root", 
            aliRun->GetRunNumber(), aliRun->GetEvNumber());
     if (!gSystem->AccessPathName(fileName)) continue;
 
-    AliESD* esd = new AliESD;
+    esd = new AliESD;
     esd->SetRunNumber(aliRun->GetRunNumber());
     esd->SetEventNumber(aliRun->GetEvNumber());
     esd->SetMagneticField(aliRun->Field()->SolenoidField());
@@ -296,19 +365,14 @@ Bool_t AliReconstruction::Run()
     if (fCheckPointLevel > 1) WriteESD(esd, "PID");
 
     // write ESD
-    char name[100]; 
-    sprintf(name, "ESD%d", iEvent);
-    file->cd();
-    if (!esd->Write(name)) {
-      Error("Run", "writing ESD failed");
-      if (fStopOnError) {CleanUp(file); return kFALSE;}
-    }
-    file->Flush();
+    tree->Fill();
 
     if (fCheckPointLevel > 0) WriteESD(esd, "final");
     delete esd;
   }
 
+  file->cd();
+  tree->Write();
   CleanUp(file);
 
   return kTRUE;
@@ -324,35 +388,33 @@ Bool_t AliReconstruction::RunLocalReconstruction(const TString& detectors)
   stopwatch.Start();
 
   TString detStr = detectors;
-  TObjArray* detArray = fRunLoader->GetAliRun()->Detectors();
-  for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
-    AliModule* det = (AliModule*) detArray->At(iDet);
-    if (!det || !det->IsActive()) continue;
-    if (IsSelected(det->GetName(), detStr)) {
-      Info("RunReconstruction", "running reconstruction for %s", 
-          det->GetName());
+  for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
+    AliReconstructor* reconstructor = 
+      (AliReconstructor*) fReconstructors[iDet];
+    TString detName = reconstructor->GetDetectorName();
+    if (IsSelected(detName, detStr)) {
+      Info("RunLocalReconstruction", "running reconstruction for %s", 
+          detName.Data());
       TStopwatch stopwatchDet;
       stopwatchDet.Start();
-      AliReconstructor* reconstructor = (AliReconstructor*) 
-       fReconstructors.FindObject("Ali" + TString(det->GetName()) + 
-                                  "Reconstructor");
-      if (reconstructor) {
-       reconstructor->Reconstruct(fRunLoader);
+      if (fRawReader) {
+       fRawReader->RewindEvents();
+       reconstructor->Reconstruct(fRunLoader, fRawReader);
       } else {
-       det->Reconstruct();
+       reconstructor->Reconstruct(fRunLoader);
       }
-      Info("RunReconstruction", "execution time for %s:", det->GetName());
+      Info("RunLocalReconstruction", "execution time for %s:", detName.Data());
       stopwatchDet.Print();
     }
   }
 
   if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
-    Error("RunReconstruction", "the following detectors were not found: %s", 
-         detStr.Data());
+    Error("RunLocalReconstruction", 
+         "the following detectors were not found: %s", detStr.Data());
     if (fStopOnError) return kFALSE;
   }
 
-  Info("RunReconstruction", "execution time:");
+  Info("RunLocalReconstruction", "execution time:");
   stopwatch.Print();
 
   return kTRUE;
@@ -370,8 +432,10 @@ Bool_t AliReconstruction::RunVertexFinder(AliESD*& esd)
   Double_t vtxPos[3] = {0, 0, 0};
   Double_t vtxErr[3] = {0.07, 0.07, 0.1};
   TArrayF mcVertex(3); 
-  fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
-  for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
+  if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
+    fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
+    for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
+  }
 
   if (fITSVertexer) {
     Info("RunVertexFinder", "running the ITS vertex finder");
@@ -441,7 +505,7 @@ Bool_t AliReconstruction::RunTracking(AliESD*& esd)
     Warning("RunTracking", "no ITS tracker");
   } else {
 
-    fRunLoader->GetAliRun()->GetDetector("TPC")->FillESD(esd); // preliminary
+    GetReconstructor("TPC")->FillESD(fRunLoader, esd); // preliminary
     AliESDpid::MakePID(esd);                  // PID for the ITS tracker
 
     // ITS tracking
@@ -563,23 +627,19 @@ Bool_t AliReconstruction::FillESD(AliESD*& esd, const TString& detectors)
   stopwatch.Start();
 
   TString detStr = detectors;
-  TObjArray* detArray = fRunLoader->GetAliRun()->Detectors();
-  for (Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++) {
-    AliModule* det = (AliModule*) detArray->At(iDet);
-    if (!det || !det->IsActive()) continue;
-    if (IsSelected(det->GetName(), detStr)) {
-      if (!ReadESD(esd, det->GetName())) {
-       Info("FillESD", "filling ESD for %s", 
-            det->GetName());
-       AliReconstructor* reconstructor = (AliReconstructor*) 
-         fReconstructors.FindObject("Ali" + TString(det->GetName()) +
-                                    "Reconstructor");
-       if (reconstructor) {
-         reconstructor->FillESD(fRunLoader, esd);
+  for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
+    AliReconstructor* reconstructor = 
+      (AliReconstructor*) fReconstructors[iDet];
+    TString detName = reconstructor->GetDetectorName();
+    if (IsSelected(detName, detStr)) {
+      if (!ReadESD(esd, detName.Data())) {
+       Info("FillESD", "filling ESD for %s", detName.Data());
+       if (fRawReader) {
+         reconstructor->FillESD(fRunLoader, fRawReader, esd);
        } else {
-         det->FillESD(esd);
+         reconstructor->FillESD(fRunLoader, esd);
        }
-       if (fCheckPointLevel > 2) WriteESD(esd, det->GetName());
+       if (fCheckPointLevel > 2) WriteESD(esd, detName.Data());
       }
     }
   }
@@ -630,21 +690,30 @@ Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
   return result;
 }
 
+//_____________________________________________________________________________
+AliReconstructor* AliReconstruction::GetReconstructor(const char* detName) const
+{
+// get the reconstructor object for a detector
+
+  for (Int_t iDet = 0; iDet < fReconstructors.GetEntriesFast(); iDet++) {
+    AliReconstructor* reconstructor = 
+      (AliReconstructor*) fReconstructors[iDet];
+    if (strcmp(reconstructor->GetDetectorName(), detName) == 0) {
+      return reconstructor;
+    }
+  }
+  return NULL;
+}
+
 //_____________________________________________________________________________
 Bool_t AliReconstruction::CreateVertexer()
 {
 // create the vertexer
 
   fITSVertexer = NULL;
-  AliReconstructor* itsReconstructor = (AliReconstructor*) 
-    fReconstructors.FindObject("AliITSReconstructor");
+  AliReconstructor* itsReconstructor = GetReconstructor("ITS");
   if (itsReconstructor) {
     fITSVertexer = itsReconstructor->CreateVertexer(fRunLoader);
-  } else {
-    AliRun* aliRun = fRunLoader->GetAliRun();
-    if (aliRun->GetDetector("ITS")) {
-      fITSVertexer = aliRun->GetDetector("ITS")->CreateVertexer();
-    }
   }
   if (!fITSVertexer) {
     Warning("CreateVertexer", "couldn't create a vertexer for ITS");
@@ -659,22 +728,15 @@ Bool_t AliReconstruction::CreateTrackers()
 {
 // get the loaders and create the trackers
 
-  AliRun* aliRun = fRunLoader->GetAliRun();
-
   fITSTracker = NULL;
   fITSLoader = fRunLoader->GetLoader("ITSLoader");
   if (!fITSLoader) {
     Warning("CreateTrackers", "no ITS loader found");
     if (fStopOnError) return kFALSE;
   } else {
-    AliReconstructor* itsReconstructor = (AliReconstructor*) 
-      fReconstructors.FindObject("AliITSReconstructor");
+    AliReconstructor* itsReconstructor = GetReconstructor("ITS");
     if (itsReconstructor) {
       fITSTracker = itsReconstructor->CreateTracker(fRunLoader);
-    } else {
-      if (aliRun->GetDetector("ITS")) {
-       fITSTracker = aliRun->GetDetector("ITS")->CreateTracker();
-      }
     }
     if (!fITSTracker) {
       Warning("CreateTrackers", "couldn't create a tracker for ITS");
@@ -688,14 +750,9 @@ Bool_t AliReconstruction::CreateTrackers()
     Error("CreateTrackers", "no TPC loader found");
     if (fStopOnError) return kFALSE;
   } else {
-    AliReconstructor* tpcReconstructor = (AliReconstructor*) 
-      fReconstructors.FindObject("AliTPCReconstructor");
+    AliReconstructor* tpcReconstructor = GetReconstructor("TPC");
     if (tpcReconstructor) {
       fTPCTracker = tpcReconstructor->CreateTracker(fRunLoader);
-    } else {
-      if (aliRun->GetDetector("TPC")) {
-       fTPCTracker = aliRun->GetDetector("TPC")->CreateTracker();
-      }
     }
     if (!fTPCTracker) {
       Error("CreateTrackers", "couldn't create a tracker for TPC");
@@ -709,14 +766,9 @@ Bool_t AliReconstruction::CreateTrackers()
     Warning("CreateTrackers", "no TRD loader found");
     if (fStopOnError) return kFALSE;
   } else {
-    AliReconstructor* trdReconstructor = (AliReconstructor*) 
-      fReconstructors.FindObject("AliTRDReconstructor");
+    AliReconstructor* trdReconstructor = GetReconstructor("TRD");
     if (trdReconstructor) {
       fTRDTracker = trdReconstructor->CreateTracker(fRunLoader);
-    } else {
-      if (aliRun->GetDetector("TRD")) {
-       fTRDTracker = aliRun->GetDetector("TRD")->CreateTracker();
-      }
     }
     if (!fTRDTracker) {
       Warning("CreateTrackers", "couldn't create a tracker for TRD");
@@ -730,14 +782,9 @@ Bool_t AliReconstruction::CreateTrackers()
     Warning("CreateTrackers", "no TOF loader found");
     if (fStopOnError) return kFALSE;
   } else {
-    AliReconstructor* tofReconstructor = (AliReconstructor*) 
-      fReconstructors.FindObject("AliTOFReconstructor");
+    AliReconstructor* tofReconstructor = GetReconstructor("TOF");
     if (tofReconstructor) {
       fTOFTracker = tofReconstructor->CreateTracker(fRunLoader);
-    } else {
-      if (aliRun->GetDetector("TOF")) {
-       fTOFTracker = aliRun->GetDetector("TOF")->CreateTracker();
-      }
     }
     if (!fTOFTracker) {
       Warning("CreateTrackers", "couldn't create a tracker for TOF");
@@ -768,6 +815,8 @@ void AliReconstruction::CleanUp(TFile* file)
 
   delete fRunLoader;
   fRunLoader = NULL;
+  delete fRawReader;
+  fRawReader = NULL;
 
   if (file) {
     file->Close();