]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliReconstruction.h
increased flexibility of tracking, removal of dummy reconstructors, creation of recon...
[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 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           SetInput(const char* input) {fInput = input;};
45   void           SetOption(const char* detector, const char* option);
46
47   void           SetRunLocalReconstruction(const char* detectors) {
48     fRunLocalReconstruction = detectors;};
49   void           SetRunVertexFinder(Bool_t run) {fRunVertexFinder = run;};
50   void           SetRunTracking(const char* detectors) {
51     fRunTracking = detectors;};
52   void           SetFillESD(const char* detectors) {fFillESD = detectors;};
53   void           SetRunReconstruction(const char* detectors) {
54     SetRunLocalReconstruction(detectors); 
55     SetRunTracking(detectors);
56     SetFillESD(detectors);};
57
58   void           SetStopOnError(Bool_t stopOnError) 
59     {fStopOnError = stopOnError;}
60   void           SetCheckPointLevel(Int_t checkPointLevel)
61     {fCheckPointLevel = checkPointLevel;}
62
63   virtual Bool_t Run(const char* input = NULL);
64
65 private:
66   Bool_t         RunLocalReconstruction(const TString& detectors);
67   Bool_t         RunVertexFinder(AliESD*& esd);
68   Bool_t         RunTracking(AliESD*& esd);
69   Bool_t         FillESD(AliESD*& esd, const TString& detectors);
70
71   Bool_t         IsSelected(TString detName, TString& detectors) const;
72   AliReconstructor* GetReconstructor(Int_t iDet);
73   Bool_t         CreateVertexer();
74   Bool_t         CreateTrackers(const TString& detectors);
75   void           CleanUp(TFile* file = NULL);
76
77   Bool_t         ReadESD(AliESD*& esd, const char* recStep) const;
78   void           WriteESD(AliESD* esd, const char* recStep) const;
79
80   TString        fRunLocalReconstruction; // run the local reconstruction for these detectors
81   Bool_t         fRunVertexFinder;    // run the vertex finder
82   TString        fRunTracking;        // run the tracking for these detectors
83   TString        fFillESD;            // fill ESD for these detectors
84   TString        fGAliceFileName;     // name of the galice file
85   TString        fInput;              // name of input file or directory
86   Bool_t         fStopOnError;        // stop or continue on errors
87   Int_t          fCheckPointLevel;    // level of ESD check points
88   TObjArray      fOptions;            // options for reconstructor objects
89
90   AliRunLoader*  fRunLoader;          //! current run loader object
91   AliRawReader*  fRawReader;          //! current raw data reader
92
93   static const Int_t fgkNDetectors = 15;   //! number of detectors
94   static const char* fgkDetectorName[fgkNDetectors]; //! names of detectors
95   AliReconstructor*  fReconstructor[fgkNDetectors];  //! array of reconstructor objects
96   AliLoader*     fLoader[fgkNDetectors];   //! detector loaders
97   AliVertexer*   fVertexer;                //! vertexer for ITS
98   AliTracker*    fTracker[fgkNDetectors];  //! trackers
99
100   ClassDef(AliReconstruction, 3)      // class for running the reconstruction
101 };
102
103 #endif