]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
SetRunReconstruction renamed to SetRunLoaclReconstruction. First implementation of...
authorhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 29 Mar 2004 16:02:02 +0000 (16:02 +0000)
committerhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 29 Mar 2004 16:02:02 +0000 (16:02 +0000)
12 files changed:
STEER/AliReconstruction.cxx
STEER/AliReconstruction.h
STEER/AliReconstructor.cxx [new file with mode: 0644]
STEER/AliReconstructor.h [new file with mode: 0644]
STEER/STEERLinkDef.h
STEER/libSTEER.pkg
TPC/AliTPC.cxx
TPC/AliTPC.h
TPC/AliTPCReconstructor.cxx [new file with mode: 0644]
TPC/AliTPCReconstructor.h [new file with mode: 0644]
TPC/TPCLinkDef.h
TPC/libTPC.pkg

index 8af343d8db4914af4f2317f397523354f0d7b67b..23683df8231a4abcce551af76abef4cece07e620 100644 (file)
 //                                                                           //
 //   rec.SetGAliceFile("...");                                               //
 //                                                                           //
-// The reconstruction can be switched on or off for individual detectors by  //
+// The local reconstruction can be switched on or off for individual         //
+// detectors by                                                              //
 //                                                                           //
-//   rec.SetRunReconstruction("...");                                        //
+//   rec.SetRunLocalReconstruction("...");                                   //
 //                                                                           //
 // The argument is a (case sensitive) string with the names of the           //
 // detectors separated by a space. The special string "ALL" selects all      //
@@ -73,6 +74,7 @@
 #include "AliRun.h"
 #include "AliModule.h"
 #include "AliDetector.h"
+#include "AliReconstructor.h"
 #include "AliTracker.h"
 #include "AliESD.h"
 #include "AliESDVertex.h"
@@ -84,6 +86,7 @@
 #include <TArrayF.h>
 #include <TSystem.h>
 #include <TROOT.h>
+#include <TPluginManager.h>
 
 
 ClassImp(AliReconstruction)
@@ -94,7 +97,7 @@ AliReconstruction::AliReconstruction(const char* gAliceFilename,
                                     const char* name, const char* title) :
   TNamed(name, title),
 
-  fRunReconstruction("ALL"),
+  fRunLocalReconstruction("ALL"),
   fRunVertexFinder(kTRUE),
   fRunTracking(kTRUE),
   fFillESD("ALL"),
@@ -121,7 +124,7 @@ AliReconstruction::AliReconstruction(const char* gAliceFilename,
 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
   TNamed(rec),
 
-  fRunReconstruction(rec.fRunReconstruction),
+  fRunLocalReconstruction(rec.fRunLocalReconstruction),
   fRunVertexFinder(rec.fRunVertexFinder),
   fRunTracking(rec.fRunTracking),
   fFillESD(rec.fFillESD),
@@ -195,9 +198,29 @@ Bool_t AliReconstruction::Run()
   }
   gAlice = aliRun;
 
+  // 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;
+    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);
+  }
+
   // local reconstruction
-  if (!fRunReconstruction.IsNull()) {
-    if (!RunReconstruction(fRunReconstruction)) {
+  if (!fRunLocalReconstruction.IsNull()) {
+    if (!RunLocalReconstruction(fRunLocalReconstruction)) {
       if (fStopOnError) {CleanUp(); return kFALSE;}
     }
   }
@@ -293,9 +316,9 @@ Bool_t AliReconstruction::Run()
 
 
 //_____________________________________________________________________________
-Bool_t AliReconstruction::RunReconstruction(const TString& detectors)
+Bool_t AliReconstruction::RunLocalReconstruction(const TString& detectors)
 {
-// run the reconstruction
+// run the local reconstruction
 
   TStopwatch stopwatch;
   stopwatch.Start();
@@ -310,7 +333,14 @@ Bool_t AliReconstruction::RunReconstruction(const TString& detectors)
           det->GetName());
       TStopwatch stopwatchDet;
       stopwatchDet.Start();
-      det->Reconstruct();
+      AliReconstructor* reconstructor = (AliReconstructor*) 
+       fReconstructors.FindObject("Ali" + TString(det->GetName()) + 
+                                  "Reconstructor");
+      if (reconstructor) {
+       reconstructor->Reconstruct(fRunLoader);
+      } else {
+       det->Reconstruct();
+      }
       Info("RunReconstruction", "execution time for %s:", det->GetName());
       stopwatchDet.Print();
     }
@@ -541,7 +571,14 @@ Bool_t AliReconstruction::FillESD(AliESD*& esd, const TString& detectors)
       if (!ReadESD(esd, det->GetName())) {
        Info("FillESD", "filling ESD for %s", 
             det->GetName());
-       det->FillESD(esd);
+       AliReconstructor* reconstructor = (AliReconstructor*) 
+         fReconstructors.FindObject("Ali" + TString(det->GetName()) +
+                                    "Reconstructor");
+       if (reconstructor) {
+         reconstructor->FillESD(fRunLoader, esd);
+       } else {
+         det->FillESD(esd);
+       }
        if (fCheckPointLevel > 2) WriteESD(esd, det->GetName());
       }
     }
@@ -599,9 +636,15 @@ Bool_t AliReconstruction::CreateVertexer()
 // create the vertexer
 
   fITSVertexer = NULL;
-  AliRun* aliRun = fRunLoader->GetAliRun();
-  if (aliRun->GetDetector("ITS")) {
-    fITSVertexer = aliRun->GetDetector("ITS")->CreateVertexer();
+  AliReconstructor* itsReconstructor = (AliReconstructor*) 
+    fReconstructors.FindObject("AliITSReconstructor");
+  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");
@@ -624,8 +667,14 @@ Bool_t AliReconstruction::CreateTrackers()
     Warning("CreateTrackers", "no ITS loader found");
     if (fStopOnError) return kFALSE;
   } else {
-    if (aliRun->GetDetector("ITS")) {
-      fITSTracker = aliRun->GetDetector("ITS")->CreateTracker();
+    AliReconstructor* itsReconstructor = (AliReconstructor*) 
+      fReconstructors.FindObject("AliITSReconstructor");
+    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");
@@ -639,8 +688,14 @@ Bool_t AliReconstruction::CreateTrackers()
     Error("CreateTrackers", "no TPC loader found");
     if (fStopOnError) return kFALSE;
   } else {
-    if (aliRun->GetDetector("TPC")) {
-      fTPCTracker = aliRun->GetDetector("TPC")->CreateTracker();
+    AliReconstructor* tpcReconstructor = (AliReconstructor*) 
+      fReconstructors.FindObject("AliTPCReconstructor");
+    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");
@@ -654,8 +709,14 @@ Bool_t AliReconstruction::CreateTrackers()
     Warning("CreateTrackers", "no TRD loader found");
     if (fStopOnError) return kFALSE;
   } else {
-    if (aliRun->GetDetector("TRD")) {
-      fTRDTracker = aliRun->GetDetector("TRD")->CreateTracker();
+    AliReconstructor* trdReconstructor = (AliReconstructor*) 
+      fReconstructors.FindObject("AliTRDReconstructor");
+    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");
@@ -669,8 +730,14 @@ Bool_t AliReconstruction::CreateTrackers()
     Warning("CreateTrackers", "no TOF loader found");
     if (fStopOnError) return kFALSE;
   } else {
-    if (aliRun->GetDetector("TOF")) {
-      fTOFTracker = aliRun->GetDetector("TOF")->CreateTracker();
+    AliReconstructor* tofReconstructor = (AliReconstructor*) 
+      fReconstructors.FindObject("AliTOFReconstructor");
+    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");
@@ -686,6 +753,8 @@ void AliReconstruction::CleanUp(TFile* file)
 {
 // delete trackers and the run loader and close and delete the file
 
+  fReconstructors.Delete();
+
   delete fITSVertexer;
   fITSVertexer = NULL;
   delete fITSTracker;
index ba4ded0736c29fae3dd4acac51c67c18f1e39422..bf97f1a722993ba12b2c517e9cd4a6805a397a14 100644 (file)
@@ -7,6 +7,7 @@
 
 #include <TNamed.h>
 #include <TString.h>
+#include <TObjArray.h>
 
 class AliRunLoader;
 class AliLoader;
@@ -27,8 +28,8 @@ public:
 
   void           SetGAliceFile(const char* fileName);
 
-  void           SetRunReconstruction(const char* detectors) {
-    fRunReconstruction = detectors;};
+  void           SetRunLocalReconstruction(const char* detectors) {
+    fRunLocalReconstruction = detectors;};
   void           SetRunVertexFinder(Bool_t run) {fRunVertexFinder = run;};
   void           SetRunTracking(Bool_t run) {fRunTracking = run;};
   void           SetFillESD(const char* detectors) {fFillESD = detectors;};
@@ -41,7 +42,7 @@ public:
   virtual Bool_t Run();
 
 private:
-  Bool_t         RunReconstruction(const TString& detectors);
+  Bool_t         RunLocalReconstruction(const TString& detectors);
   Bool_t         RunVertexFinder(AliESD*& esd);
   Bool_t         RunTracking(AliESD*& esd);
   Bool_t         FillESD(AliESD*& esd, const TString& detectors);
@@ -54,7 +55,7 @@ private:
   Bool_t         ReadESD(AliESD*& esd, const char* recStep) const;
   void           WriteESD(AliESD* esd, const char* recStep) const;
 
-  TString        fRunReconstruction;  // run the reconstr. for these detectors
+  TString        fRunLocalReconstruction; // run the local reconstruction for these detectors
   Bool_t         fRunVertexFinder;    // run the vertex finder
   Bool_t         fRunTracking;        // run the barrel tracking
   TString        fFillESD;            // fill ESD for these detectors
@@ -73,6 +74,8 @@ private:
   AliLoader*     fTOFLoader;          //! loader for TOF
   AliTracker*    fTOFTracker;         //! tracker for TOF
 
+  TObjArray      fReconstructors;     //! array of reconstructor objects
+
   ClassDef(AliReconstruction, 1)      // class for running the reconstruction
 };
 
diff --git a/STEER/AliReconstructor.cxx b/STEER/AliReconstructor.cxx
new file mode 100644 (file)
index 0000000..eafd610
--- /dev/null
@@ -0,0 +1,39 @@
+/**************************************************************************
+ * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ *                                                                        *
+ * Author: The ALICE Off-line Project.                                    *
+ * Contributors are mentioned in the code where appropriate.              *
+ *                                                                        *
+ * Permission to use, copy, modify and distribute this software and its   *
+ * documentation strictly for non-commercial purposes is hereby granted   *
+ * without fee, provided that the above copyright notice appears in all   *
+ * copies and that both the copyright notice and this permission notice   *
+ * appear in the supporting documentation. The authors make no claims     *
+ * about the suitability of this software for any purpose. It is          *
+ * provided "as is" without express or implied warranty.                  *
+ **************************************************************************/
+
+/* $Id$ */
+
+///////////////////////////////////////////////////////////////////////////////
+//                                                                           //
+// base class for reconstruction algorithms                                  //
+//                                                                           //
+// Derived classes should implement the virtual methods                      //
+// - Reconstruct : to perform the local reconstruction for all events        //
+// - FillESD     : to fill the ESD for the current event                     //
+//                                                                           //
+// The reconstructor classes for the barrel detectors should in addition     //
+// implement the method                                                      //
+// - CreateTracker : to create a tracker object for the barrel detector      //
+//                                                                           //
+// The ITS reconstructor should in addition implement the method             //
+// - CreateVertexer : to create an object for the vertex finding             //
+//                                                                           //
+///////////////////////////////////////////////////////////////////////////////
+
+
+#include "AliReconstructor.h"
+
+
+ClassImp(AliReconstructor)
diff --git a/STEER/AliReconstructor.h b/STEER/AliReconstructor.h
new file mode 100644 (file)
index 0000000..3524ddc
--- /dev/null
@@ -0,0 +1,28 @@
+#ifndef ALIRECONSTRUCTOR_H
+#define ALIRECONSTRUCTOR_H
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+/* $Id$ */
+
+#include <TObject.h>
+
+class AliRunLoader;
+class AliVertexer;
+class AliTracker;
+class AliESD;
+
+
+class AliReconstructor: public TObject {
+public:
+  virtual void         Reconstruct(AliRunLoader* runLoader) const = 0;
+  virtual AliVertexer* CreateVertexer(AliRunLoader* /*runLoader*/) const 
+    {return NULL;}
+  virtual AliTracker*  CreateTracker(AliRunLoader* /*runLoader*/) const 
+    {return NULL;}
+  virtual void         FillESD(AliRunLoader* runLoader, AliESD* esd) const = 0;
+
+  ClassDef(AliReconstructor, 0)   // base class for reconstruction algorithms
+};
+
+#endif
index 184821bdba505b7081298315e4077844c08f5192..64668640838fee45f0d591ab17ba4cbe7935b15f 100644 (file)
@@ -82,6 +82,7 @@
 #pragma link C++ class  AliReconstruction+;
 #pragma link C++ class  AliVertexGenFile+;
 #pragma link C++ class  AliVertexer+;
+#pragma link C++ class  AliReconstructor+;
 #endif
 
 
index 1dc42039e7d482ae1428b91168e58f048cfe7d30..4b257b6de3697d860dc9bf5cd4d279a2af5c26ce 100644 (file)
@@ -18,7 +18,9 @@ AliTrackMap.cxx AliTrackMapper.cxx AliCollisionGeometry.cxx \
 AliMemoryWatcher.cxx AliBarrelTrack.cxx \
 AliESDtrack.cxx AliESDCaloTrack.cxx AliESDMuonTrack.cxx AliESDv0.cxx AliESDcascade.cxx AliESDVertex.cxx AliESDpid.cxx \
 AliVertexer.cxx \
-AliMC.cxx AliSimulation.cxx AliReconstruction.cxx AliVertexGenFile.cxx 
+AliMC.cxx AliSimulation.cxx AliReconstruction.cxx AliVertexGenFile.cxx \
+AliReconstructor.cxx
+
 HDRS:= $(SRCS:.cxx=.h) 
 
 DHDR= STEERLinkDef.h
index e57a7f2f9d4da9018262a30f0b381027350d8faa..0215f0cca7ae0ad238f020001c630491818017be 100644 (file)
@@ -336,64 +336,6 @@ void AliTPC::Clusters2Tracks() const
  }
 
 
-//_____________________________________________________________________________
-void AliTPC::Reconstruct() const
-{
-// reconstruct clusters
-
-  AliLoader* loader = GetLoader();
-  loader->LoadRecPoints("recreate");
-  loader->LoadDigits("read");
-
-  AliTPCclustererMI clusterer(fTPCParam);
-  AliRunLoader* runLoader = loader->GetRunLoader();
-  Int_t nEvents = runLoader->GetNumberOfEvents();
-
-  for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
-    runLoader->GetEvent(iEvent);
-
-    TTree* treeClusters = loader->TreeR();
-    if (!treeClusters) {
-      loader->MakeTree("R");
-      treeClusters = loader->TreeR();
-    }
-    TTree* treeDigits = loader->TreeD();
-    if (!treeDigits) {
-      Error("Reconstruct", "Can't get digits tree !");
-      return;
-    }
-
-//    clusterer.Digits2Clusters(treeDigits, treeClusters);
-    clusterer.SetInput(treeDigits);
-    clusterer.SetOutput(treeClusters);
-    clusterer.Digits2Clusters();
-         
-    loader->WriteRecPoints("OVERWRITE");
-  }
-
-  loader->UnloadRecPoints();
-  loader->UnloadDigits();
-}
-
-//_____________________________________________________________________________
-AliTracker* AliTPC::CreateTracker() const
-{
-// create a TPC tracker
-
-  return new AliTPCtrackerMI(fTPCParam);
-}
-
-//_____________________________________________________________________________
-void AliTPC::FillESD(AliESD* esd) const
-{
-// make PID
-
-  Double_t parTPC[] = {47., 0.10, 10.};
-  AliTPCpidESD tpcPID(parTPC);
-  tpcPID.MakePID(esd);
-}
-
-
 //_____________________________________________________________________________
 void AliTPC::CreateMaterials()
 {
index cd8541dda0aba4c8414d5b049684578d8d776135..a3ffeaf983abce084b136ad928d18ece5f74c761 100644 (file)
@@ -66,10 +66,6 @@ public:
   virtual void  Digits2Clusters(Int_t eventnumber=0) const;
   virtual void  Clusters2Tracks() const;
 
-  virtual void  Reconstruct() const;
-  virtual AliTracker* CreateTracker() const;
-  virtual void  FillESD(AliESD* esd) const;
-
   Int_t         GetNsectors() const  {return fNsectors;}
   virtual void  MakeBranch(Option_t *opt=" ");
   virtual void  ResetDigits();
diff --git a/TPC/AliTPCReconstructor.cxx b/TPC/AliTPCReconstructor.cxx
new file mode 100644 (file)
index 0000000..b0923e8
--- /dev/null
@@ -0,0 +1,122 @@
+/**************************************************************************
+ * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ *                                                                        *
+ * Author: The ALICE Off-line Project.                                    *
+ * Contributors are mentioned in the code where appropriate.              *
+ *                                                                        *
+ * Permission to use, copy, modify and distribute this software and its   *
+ * documentation strictly for non-commercial purposes is hereby granted   *
+ * without fee, provided that the above copyright notice appears in all   *
+ * copies and that both the copyright notice and this permission notice   *
+ * appear in the supporting documentation. The authors make no claims     *
+ * about the suitability of this software for any purpose. It is          *
+ * provided "as is" without express or implied warranty.                  *
+ **************************************************************************/
+
+/* $Id$ */
+
+///////////////////////////////////////////////////////////////////////////////
+//                                                                           //
+// class for TPC reconstruction                                              //
+//                                                                           //
+///////////////////////////////////////////////////////////////////////////////
+
+
+#include "AliTPCReconstructor.h"
+#include "AliRunLoader.h"
+#include "AliRun.h"
+#include "AliTPC.h"
+#include "AliTPCclustererMI.h"
+#include "AliTPCtrackerMI.h"
+#include "AliTPCpidESD.h"
+
+
+ClassImp(AliTPCReconstructor)
+
+
+//_____________________________________________________________________________
+void AliTPCReconstructor::Reconstruct(AliRunLoader* runLoader) const
+{
+// reconstruct clusters
+
+  AliLoader* loader = runLoader->GetLoader("TPCLoader");
+  if (!loader) {
+    Error("Reconstruct", "TPC loader not found");
+    return;
+  }
+  loader->LoadRecPoints("recreate");
+  loader->LoadDigits("read");
+
+  AliTPCParam* param = GetTPCParam(runLoader);
+  if (!param) return;
+  AliTPCclustererMI clusterer(param);
+  Int_t nEvents = runLoader->GetNumberOfEvents();
+
+  for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
+    runLoader->GetEvent(iEvent);
+
+    TTree* treeClusters = loader->TreeR();
+    if (!treeClusters) {
+      loader->MakeTree("R");
+      treeClusters = loader->TreeR();
+    }
+    TTree* treeDigits = loader->TreeD();
+    if (!treeDigits) {
+      Error("Reconstruct", "Can't get digits tree !");
+      return;
+    }
+
+    clusterer.SetInput(treeDigits);
+    clusterer.SetOutput(treeClusters);
+    clusterer.Digits2Clusters();
+         
+    loader->WriteRecPoints("OVERWRITE");
+  }
+
+  loader->UnloadRecPoints();
+  loader->UnloadDigits();
+}
+
+//_____________________________________________________________________________
+AliTracker* AliTPCReconstructor::CreateTracker(AliRunLoader* runLoader) const
+{
+// create a TPC tracker
+
+  AliTPCParam* param = GetTPCParam(runLoader);
+  if (!param) return NULL;
+  return new AliTPCtrackerMI(param);
+}
+
+//_____________________________________________________________________________
+void AliTPCReconstructor::FillESD(AliRunLoader* /*runLoader*/, 
+                                 AliESD* esd) const
+{
+// make PID
+
+  Double_t parTPC[] = {47., 0.10, 10.};
+  AliTPCpidESD tpcPID(parTPC);
+  tpcPID.MakePID(esd);
+}
+
+
+//_____________________________________________________________________________
+AliTPCParam* AliTPCReconstructor::GetTPCParam(AliRunLoader* runLoader) const
+{
+// get the TPC parameters
+
+  if (!runLoader->GetAliRun()) runLoader->LoadgAlice();
+  if (!runLoader->GetAliRun()) {
+    Error("GetTPCParam", "couldn't get AliRun object");
+    return NULL;
+  }
+  AliTPC* tpc = (AliTPC*) runLoader->GetAliRun()->GetDetector("TPC");
+  if (!tpc) {
+    Error("GetTPCParam", "couldn't get TPC detector");
+    return NULL;
+  }
+  if (!tpc->GetParam()) {
+    Error("GetTPCParam", "no TPC parameters available");
+    return NULL;
+  }
+  return tpc->GetParam();
+}
diff --git a/TPC/AliTPCReconstructor.h b/TPC/AliTPCReconstructor.h
new file mode 100644 (file)
index 0000000..1831d75
--- /dev/null
@@ -0,0 +1,25 @@
+#ifndef ALITPCRECONSTRUCTOR_H
+#define ALITPCRECONSTRUCTOR_H
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+/* $Id$ */
+
+#include "AliReconstructor.h"
+
+class AliTPCParam;
+
+
+class AliTPCReconstructor: public AliReconstructor {
+public:
+  virtual void         Reconstruct(AliRunLoader* runLoader) const;
+  virtual AliTracker*  CreateTracker(AliRunLoader* runLoader) const;
+  virtual void         FillESD(AliRunLoader* runLoader, AliESD* esd) const;
+
+private:
+  AliTPCParam*         GetTPCParam(AliRunLoader* runLoader) const;
+
+  ClassDef(AliTPCReconstructor, 0)   // class for the TPC reconstruction
+};
+
+#endif
index c4c875a83a03d814a7d9ad3defe1166b77e21df3..18e61c29eb0a2e178fc34aadf7ca3eabafa3490d 100644 (file)
@@ -82,6 +82,8 @@
 #pragma link C++ class  AliTPCtrackPid+;
 #pragma link C++ class  AliHelix+;
 
+#pragma link C++ class  AliTPCReconstructor+;
+
 
 
 #endif
index 880dd761a8934805016bd22bcef0985e2a4aa2dc..bd7b408e7d02d42574bcfb66a77eaae9cb906e53 100644 (file)
@@ -15,7 +15,8 @@ SRCS:=  AliTPCLoader.cxx\
        AliTPCclustererMI.cxx AliTPCclusterMI.cxx  AliTPCtrackerMI.cxx \
        AliTPCpolyTrack.cxx \
        AliTPCBuffer.cxx AliTPCDDLRawData.cxx \
-       AliTPCtrackPid.cxx AliTPCpidESD.cxx AliHelix.cxx
+       AliTPCtrackPid.cxx AliTPCpidESD.cxx AliHelix.cxx \
+       AliTPCReconstructor.cxx
 
 HDRS:= $(SRCS:.cxx=.h)