]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliReconstruction.h
Removing run loader
[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
24 class AliReconstructor;
25 class AliRunLoader;
26 class AliRawReader;
27 class AliLoader;
28 class AliTracker;
29 class AliVertexer;
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           SetInput(const char* input) {fInput = input;};
45   void           SetEventRange(Int_t firstEvent = 0, Int_t lastEvent = -1) 
46     {fFirstEvent = firstEvent; fLastEvent = lastEvent;};
47   void           SetOption(const char* detector, const char* option);
48
49   void           SetRunLocalReconstruction(const char* detectors) {
50     fRunLocalReconstruction = detectors;};
51   void           SetRunVertexFinder(Bool_t run) {fRunVertexFinder = run;};
52   void           SetRunTracking(const char* detectors) {
53     fRunTracking = detectors;};
54   void           SetFillESD(const char* detectors) {fFillESD = detectors;};
55   void           SetRunReconstruction(const char* detectors) {
56     SetRunLocalReconstruction(detectors); 
57     SetRunTracking(detectors);
58     SetFillESD(detectors);};
59
60    void SetUniformFieldTracking(){fUniformField=kTRUE;} 
61    void SetNonuniformFieldTracking(){fUniformField=kFALSE;} 
62
63   void           SetStopOnError(Bool_t stopOnError) 
64     {fStopOnError = stopOnError;}
65   void           SetCheckPointLevel(Int_t checkPointLevel)
66     {fCheckPointLevel = checkPointLevel;}
67
68   virtual Bool_t Run(const char* input, 
69                      Int_t firstEvent, Int_t lastEvent = -1);
70   Bool_t         Run(const char* input = NULL)
71     {return Run(input, fFirstEvent, fLastEvent);};
72   Bool_t         Run(Int_t firstEvent, Int_t lastEvent = -1)
73     {return Run(NULL, firstEvent, lastEvent);};
74
75 private:
76   Bool_t         RunLocalReconstruction(const TString& detectors);
77   Bool_t         RunLocalEventReconstruction(const TString& detectors);
78   Bool_t         RunVertexFinder(AliESD*& esd);
79   Bool_t         RunHLTTracking(AliESD*& esd);
80   Bool_t         RunTracking(AliESD*& esd);
81   Bool_t         FillESD(AliESD*& esd, const TString& detectors);
82
83   Bool_t         IsSelected(TString detName, TString& detectors) const;
84   Bool_t         InitRunLoader();
85   AliReconstructor* GetReconstructor(Int_t iDet);
86   Bool_t         CreateVertexer();
87   Bool_t         CreateTrackers(const TString& detectors);
88   void           CleanUp(TFile* file = NULL, TFile* fileOld = NULL);
89
90   Bool_t         ReadESD(AliESD*& esd, const char* recStep) const;
91   void           WriteESD(AliESD* esd, const char* recStep) const;
92
93   TString        fRunLocalReconstruction; // run the local reconstruction for these detectors
94   Bool_t         fUniformField;       // uniform field tracking flag
95   Bool_t         fRunVertexFinder;    // run the vertex finder
96   Bool_t         fRunHLTTracking;     // run the HLT tracking
97   TString        fRunTracking;        // run the tracking for these detectors
98   TString        fFillESD;            // fill ESD for these detectors
99   TString        fGAliceFileName;     // name of the galice file
100   TString        fInput;              // name of input file or directory
101   Int_t          fFirstEvent;         // index of first event to be reconstr.
102   Int_t          fLastEvent;          // index of last event to be reconstr.
103   Bool_t         fStopOnError;        // stop or continue on errors
104   Int_t          fCheckPointLevel;    // level of ESD check points
105   TObjArray      fOptions;            // options for reconstructor objects
106
107   AliRunLoader*  fRunLoader;          //! current run loader object
108   AliRawReader*  fRawReader;          //! current raw data reader
109
110   static const Int_t fgkNDetectors = 15;   //! number of detectors
111   static const char* fgkDetectorName[fgkNDetectors]; //! names of detectors
112   AliReconstructor*  fReconstructor[fgkNDetectors];  //! array of reconstructor objects
113   AliLoader*     fLoader[fgkNDetectors];   //! detector loaders
114   AliVertexer*   fVertexer;                //! vertexer for ITS
115   AliTracker*    fTracker[fgkNDetectors];  //! trackers
116
117   ClassDef(AliReconstruction, 5)      // class for running the reconstruction
118 };
119
120 #endif