]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliReconstruction.h
d402f9bbc7e5e358e84f0ea9accbd3f1b27ab609
[u/mrichter/AliRoot.git] / STEER / AliReconstruction.h
1 #ifndef ALIRECONSTRUCTION_H
2 #define ALIRECONSTRUCTION_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 /* $Id$ */
7
8 ///////////////////////////////////////////////////////////////////////////////
9 //                                                                           //
10 // class for running the reconstruction                                      //
11 // Clusters and tracks are created for all detectors and all events by       //
12 // typing:                                                                   //
13 //                                                                           //
14 //   AliReconstruction rec;                                                  //
15 //   rec.Run();                                                              //
16 //                                                                           //
17 ///////////////////////////////////////////////////////////////////////////////
18
19
20 #include <TNamed.h>
21 #include <TString.h>
22 #include <TObjArray.h>
23 #include "AliReconstructor.h"
24 #include "AliDetector.h"
25
26 class AliRunLoader;
27 class AliLoader;
28 class AliVertexer;
29 class AliTracker;
30 class AliESD;
31 class TFile;
32
33
34 class AliReconstruction: public TNamed {
35 public:
36   AliReconstruction(const char* gAliceFilename = "galice.root",
37                     const char* name = "AliReconstruction", 
38                     const char* title = "reconstruction");
39   AliReconstruction(const AliReconstruction& rec);
40   AliReconstruction& operator = (const AliReconstruction& rec);
41   virtual ~AliReconstruction();
42
43   void           SetGAliceFile(const char* fileName);
44   void           SetOption(const char* detector, const char* option);
45
46   void           SetRunLocalReconstruction(const char* detectors) {
47     fRunLocalReconstruction = detectors;};
48   void           SetRunVertexFinder(Bool_t run) {fRunVertexFinder = run;};
49   void           SetRunTracking(Bool_t run) {fRunTracking = run;};
50   void           SetFillESD(const char* detectors) {fFillESD = detectors;};
51
52   void           SetStopOnError(Bool_t stopOnError) 
53     {fStopOnError = stopOnError;}
54   void           SetCheckPointLevel(Int_t checkPointLevel)
55     {fCheckPointLevel = checkPointLevel;}
56
57   virtual Bool_t Run();
58
59 private:
60   class AliDummyReconstructor: public AliReconstructor {
61   public:
62     AliDummyReconstructor(AliDetector* detector) {fDetector = detector;};
63     virtual ~AliDummyReconstructor() {};
64
65     virtual void         Reconstruct(AliRunLoader* /*runLoader*/) const
66       {fDetector->Reconstruct();};
67     virtual AliVertexer* CreateVertexer(AliRunLoader* /*runLoader*/) const 
68       {return fDetector->CreateVertexer();}
69     virtual AliTracker*  CreateTracker(AliRunLoader* /*runLoader*/) const 
70       {return fDetector->CreateTracker();}
71     virtual void         FillESD(AliRunLoader* /*runLoader*/, AliESD* esd) const
72       {fDetector->FillESD(esd);};
73
74     virtual const char*  GetDetectorName() const
75       {return fDetector->GetName();};
76   private:
77     AliDummyReconstructor(const AliDummyReconstructor &drc):
78       AliReconstructor(drc)
79       {Fatal("copy ctor","Not implemented\n");}
80     AliDummyReconstructor & operator=(const AliDummyReconstructor &)
81       {Fatal("= operator","Not implemented\n"); return *this;}
82     AliDetector*         fDetector;   // detector object
83   };
84
85   Bool_t         RunLocalReconstruction(const TString& detectors);
86   Bool_t         RunVertexFinder(AliESD*& esd);
87   Bool_t         RunTracking(AliESD*& esd);
88   Bool_t         FillESD(AliESD*& esd, const TString& detectors);
89
90   Bool_t         IsSelected(TString detName, TString& detectors) const;
91   AliReconstructor* GetReconstructor(const char* detName) const;
92   Bool_t         CreateVertexer();
93   Bool_t         CreateTrackers();
94   void           CleanUp(TFile* file = NULL);
95
96   Bool_t         ReadESD(AliESD*& esd, const char* recStep) const;
97   void           WriteESD(AliESD* esd, const char* recStep) const;
98
99   TString        fRunLocalReconstruction; // run the local reconstruction for these detectors
100   Bool_t         fRunVertexFinder;    // run the vertex finder
101   Bool_t         fRunTracking;        // run the barrel tracking
102   TString        fFillESD;            // fill ESD for these detectors
103   TString        fGAliceFileName;     // name of the galice file
104   Bool_t         fStopOnError;        // stop or continue on errors
105   Int_t          fCheckPointLevel;    // level of ESD check points
106
107   AliRunLoader*  fRunLoader;          //! current run loader object
108   AliLoader*     fITSLoader;          //! loader for ITS
109   AliVertexer*   fITSVertexer;        //! vertexer for ITS
110   AliTracker*    fITSTracker;         //! tracker for ITS
111   AliLoader*     fTPCLoader;          //! loader for TPC
112   AliTracker*    fTPCTracker;         //! tracker for TPC
113   AliLoader*     fTRDLoader;          //! loader for TRD
114   AliTracker*    fTRDTracker;         //! tracker for TRD
115   AliLoader*     fTOFLoader;          //! loader for TOF
116   AliTracker*    fTOFTracker;         //! tracker for TOF
117
118   static const Int_t fgkNDetectors = 15;   //! number of detectors
119   static const char* fgkDetectorName[fgkNDetectors]; //! names of detectors
120   TObjArray      fReconstructors;     //! array of reconstructor objects
121   TObjArray      fOptions;            // options for reconstructor objects
122
123   ClassDef(AliReconstruction, 1)      // class for running the reconstruction
124 };
125
126 #endif