]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDInput.h
New detector numbering scheme (common for DAQ/HLT/Offline). All the subdetectors...
[u/mrichter/AliRoot.git] / FMD / AliFMDInput.h
1 #ifndef AliFMDInput_H
2 #define AliFMDInput_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights
4  * reserved. 
5  *
6  * See cxx source for full Copyright notice                               
7  */
8 //___________________________________________________________________
9 //
10 // The classes defined here, are utility classes for reading in data
11 // for the FMD.  They are  put in a seperate library to not polute the
12 // normal libraries.  The classes are intended to be used as base
13 // classes for customized class that do some sort of analysis on the
14 // various types of data produced by the FMD. 
15 /** @file    AliFMDInput.h
16     @author  Christian Holm Christensen <cholm@nbi.dk>
17     @date    Mon Mar 27 12:42:40 2006
18     @brief   FMD utility classes for reading FMD data
19 */
20 //___________________________________________________________________
21 /** @defgroup FMD_util Utility classes. 
22
23     The classes defined here, are utility classes for reading in data
24     for the FMD.  They are put in a seperate library to not polute the
25     normal libraries.  The classes are intended to be used as base
26     classes for customized class that do some sort of analysis on the
27     various types of data produced by the FMD.
28 */
29 #include <TObject.h>
30 #ifndef ROOT_TString
31 # include <TString.h>
32 #endif
33 class AliRunLoader;
34 class AliLoader;
35 class AliStack;
36 class AliRun;
37 class AliRawReader;
38 class AliFMD;
39 class AliFMDHit;
40 class AliFMDDigit;
41 class AliFMDSDigit;
42 class AliFMDRecPoint;
43 class AliESD;
44 class AliESDFMD;
45 class TString;
46 class TClonesArray;
47 class TTree;
48 class TGeoManager;
49 class TParticle;
50 class TChain;
51
52 //___________________________________________________________________
53 /** @class AliFMDInput 
54     @brief Base class for reading in various FMD data. 
55     The class loops over all found events.  For each event the
56     specified data is read in.   The class then loops over all
57     elements of the read data, and process these with user defined
58     code. 
59     @code 
60     struct DigitInput : public AliFMDInput 
61     { 
62       DigitInput() 
63       { 
64          // Load digits 
65          AddLoad(kDigits);
66          // Make a histogram 
67          fHist = new TH1F("adc", "ADC spectra", 1024, -.5, 1023.5);
68       }
69       // Process one digit. 
70       Bool_t ProcessDigit(AliFMDDigit* d)
71       { 
72          fHist->Fill(d->Counts());
73          return kTRUE;
74       }
75       // After processing all events, display spectrum
76       Bool_t Finish()
77       {
78         fHist->Draw();
79       }
80       TH1F* fHist;
81     };
82
83     void AdcSpectrum()
84     { 
85        DigitInput di;
86        di.Run();
87     }
88     @endcode
89     This class allows for writing small scripts, that can be compiled
90     with AcLIC, to do all sorts of tests, quick prototyping, and so
91     on.  It has proven to be quiet useful.   One can load more than
92     one type of data in one derived class, to for example to make
93     comparisons between hits and reconstructed points.   See also the
94     various scripts in @c FMD/scripts. 
95     @ingroup FMD_util
96  */
97 class AliFMDInput : public TObject
98 {
99 public:
100   /** The kinds of data that can be read in. */
101   enum ETrees {
102     kHits       = 1,  // Hits
103     kKinematics,      // Kinematics (from sim)
104     kDigits,          // Digits
105     kSDigits,         // Summable digits 
106     kHeader,          // Header information 
107     kRecPoints,       // Reconstructed points
108     kESD,             // Load ESD's
109     kRaw,             // Read raw data 
110     kGeometry         // Not really a tree 
111   };
112   /** CTOR  */
113   AliFMDInput();
114   /** CTOR
115       @param gAliceFile galice file  */
116   AliFMDInput(const char* gAliceFile);
117   /** DTOR */
118   virtual ~AliFMDInput() {}
119
120   /** Add a data type to load 
121       @param tree Data to load */
122   virtual void   AddLoad(ETrees tree)     { SETBIT(fTreeMask, tree); }
123   /** Remove a data type to load 
124       @param tree Data to @e not load */
125   virtual void   RemoveLoad(ETrees tree)  { CLRBIT(fTreeMask, tree); }
126   /** @return # of available events */
127   virtual Int_t  NEvents() const;
128
129   /** Initialize the class.  If a user class overloads this member
130       function, then this @e must be explicitly called
131       @return @c false on error */
132   virtual Bool_t Init();
133   /** Callled at the beginning of each event. If a user class
134       overloads this member  function, then this @e must be explicitly
135       called. 
136       @param event Event number
137       @return @c false on error */
138   virtual Bool_t Begin(Int_t event);
139   /** Process one event.  This loops over all the loaded data.   Users
140       can overload this member function, but then it's @e strongly
141       recommended to explicitly call this classes version. 
142       @return @c false on error  */
143   virtual Bool_t Event();
144   /** Called at the end of each event. 
145       @return @c false on error  */
146   virtual Bool_t End();
147   /** Called at the end of the run.
148       @return  @c false on error */
149   virtual Bool_t Finish() { return kTRUE; }
150   /** Run a full job. 
151       @return @c false on error  */
152   virtual Bool_t Run();
153
154   /** Loop over all hits, and call ProcessHit with that hit, and
155       optionally the corresponding kinematics track. 
156       @return @c false on error  */
157   virtual Bool_t ProcessHits();
158   /** Loop over all digits, and call ProcessDigit for each digit.
159       @return @c false on error  */
160   virtual Bool_t ProcessDigits();
161   /** Loop over all summable digits, and call ProcessSDigit for each
162       digit. 
163       @return @c false on error  */
164   virtual Bool_t ProcessSDigits();
165   /** Loop over all digits read from raw data files, and call
166       ProcessRawDigit for each digit. 
167       @return @c false on error  */
168   virtual Bool_t ProcessRawDigits();
169   /** Loop over all reconstructed points, and call ProcessRecPoint for
170       each reconstructed point. 
171       @return @c false on error  */
172   virtual Bool_t ProcessRecPoints();
173
174   /** Process one hit, and optionally it's corresponding kinematics
175       track.  Users should over this to process each hit. 
176       @return  @c false on error   */
177   virtual Bool_t ProcessHit(AliFMDHit*, TParticle*)  { return kTRUE; }
178   /** Process one digit.  Users should over this to process each digit. 
179       @return  @c false on error   */
180   virtual Bool_t ProcessDigit(AliFMDDigit*)          { return kTRUE; }
181   /** Process one summable digit.  Users should over this to process
182       each summable digit.  
183       @return  @c false on error   */
184   virtual Bool_t ProcessSDigit(AliFMDSDigit*)        { return kTRUE; }
185   /** Process one digit from raw data files.  Users should over this
186       to process each raw digit.  
187       @return  @c false on error   */
188   virtual Bool_t ProcessRawDigit(AliFMDDigit*)       { return kTRUE; }
189   /** Process one reconstructed point.  Users should over this to
190       process each reconstructed point.  
191       @return  @c false on error   */
192   virtual Bool_t ProcessRecPoint(AliFMDRecPoint*)    { return kTRUE; }
193   /** Process ESD data for the FMD.  Users should overload this to
194       deal with ESD data. 
195       @return  @c false on error  */
196   virtual Bool_t ProcessESD(AliESDFMD*)              { return kTRUE; }
197   
198 protected:
199   /** Copy ctor 
200       @param o Object to copy from  */
201   AliFMDInput(const AliFMDInput& o) : TObject(o) {}
202   /** Assignement operator 
203       @return  REference to this */
204   AliFMDInput& operator=(const AliFMDInput&) { return *this; }
205
206   TString       fGAliceFile; // File name of gAlice file
207   AliRunLoader* fLoader;     // Loader of FMD data 
208   AliRun*       fRun;        // Run information
209   AliStack*     fStack;      // Stack of particles 
210   AliLoader*    fFMDLoader;  // Loader of FMD data 
211   AliRawReader* fReader;     // Raw data reader 
212   AliFMD*       fFMD;        // FMD object
213   AliESD*       fMainESD;    // ESD Object
214   AliESDFMD*    fESD;        // FMD ESD data  
215   TTree*        fTreeE;      // Header tree 
216   TTree*        fTreeH;      // Hits tree
217   TTree*        fTreeD;      // Digit tree 
218   TTree*        fTreeS;      // SDigit tree 
219   TTree*        fTreeR;      // RecPoint tree
220   TTree*        fTreeA;      // Raw data tree
221   TChain*       fChainE;     // Chain of ESD's
222   TClonesArray* fArrayE;     // Event info array
223   TClonesArray* fArrayH;     // Hit info array
224   TClonesArray* fArrayD;     // Digit info array
225   TClonesArray* fArrayS;     // SDigit info array
226   TClonesArray* fArrayR;     // Rec points info array
227   TClonesArray* fArrayA;     // Raw data (digits) info array
228   TGeoManager*  fGeoManager; // Geometry manager
229   Int_t         fTreeMask;   // Which tree's to load
230   Bool_t        fIsInit;     // Have we been initialized 
231   ClassDef(AliFMDInput,0)  //Hits for detector FMD
232 };
233
234
235 #endif
236 //____________________________________________________________________
237 //
238 // Local Variables:
239 //   mode: C++
240 // End:
241 //
242 // EOF
243 //