]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDInput.h
AliEveEventManager
[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 <TNamed.h>
30 #ifndef ROOT_TString
31 # include <TString.h>
32 #endif
33 #ifndef ROOT_TArrayF
34 # include <TArrayF.h>
35 #endif
36 class AliRunLoader;
37 class AliLoader;
38 class AliStack;
39 class AliRun;
40 class AliRawReader;
41 class AliFMDRawReader;
42 class AliFMD;
43 class AliFMDHit;
44 class AliFMDDigit;
45 class AliFMDSDigit;
46 class AliFMDRecPoint;
47 class AliESDEvent;
48 class AliESDFMD;
49 class AliHeader;
50 class TString;
51 class TClonesArray;
52 class TTree;
53 class TGeoManager;
54 class TParticle;
55 class TChain;
56
57 //___________________________________________________________________
58 /** @class AliFMDInput 
59     @brief Base class for reading in various FMD data. 
60     The class loops over all found events.  For each event the
61     specified data is read in.   The class then loops over all
62     elements of the read data, and process these with user defined
63     code. 
64     @code 
65     struct DigitInput : public AliFMDInput 
66     { 
67       DigitInput() 
68       { 
69          // Load digits 
70          AddLoad(kDigits);
71          // Make a histogram 
72          fHist = new TH1F("adc", "ADC spectra", 1024, -.5, 1023.5);
73       }
74       // Process one digit. 
75       Bool_t ProcessDigit(AliFMDDigit* d)
76       { 
77          fHist->Fill(d->Counts());
78          return kTRUE;
79       }
80       // After processing all events, display spectrum
81       Bool_t Finish()
82       {
83         fHist->Draw();
84       }
85       TH1F* fHist;
86     };
87
88     void AdcSpectrum()
89     { 
90        DigitInput di;
91        di.Run();
92     }
93     @endcode
94     This class allows for writing small scripts, that can be compiled
95     with AcLIC, to do all sorts of tests, quick prototyping, and so
96     on.  It has proven to be quiet useful.   One can load more than
97     one type of data in one derived class, to for example to make
98     comparisons between hits and reconstructed points.   See also the
99     various scripts in @c FMD/scripts. 
100     @ingroup FMD_util
101  */
102 class AliFMDInput : public TNamed
103 {
104 public:
105   /** The kinds of data that can be read in. */
106   enum ETrees {
107     kHits       = 1,  // Hits
108     kKinematics,      // Kinematics (from sim)
109     kDigits,          // Digits
110     kSDigits,         // Summable digits 
111     kHeader,          // Header information 
112     kRecPoints,       // Reconstructed points
113     kESD,             // Load ESD's
114     kRaw,             // Read raw data 
115     kGeometry,        // Not really a tree 
116     kTracks           // Hits and tracs - for BG study  
117   };
118   /** CTOR  */
119   AliFMDInput();
120   /** CTOR
121       @param gAliceFile galice file  */
122   AliFMDInput(const char* gAliceFile);
123   /** DTOR */
124   virtual ~AliFMDInput() {}
125
126   /** Add a data type to load 
127       @param tree Data to load */
128   virtual void   AddLoad(ETrees tree)     { SETBIT(fTreeMask, tree); }
129   /** Remove a data type to load 
130       @param tree Data to @e not load */
131   virtual void   RemoveLoad(ETrees tree)  { CLRBIT(fTreeMask, tree); }
132   /** @return # of available events */
133   virtual Int_t  NEvents() const;
134
135   /** Initialize the class.  If a user class overloads this member
136       function, then this @e must be explicitly called
137       @return @c false on error */
138   virtual Bool_t Init();
139   /** Callled at the beginning of each event. If a user class
140       overloads this member  function, then this @e must be explicitly
141       called. 
142       @param event Event number
143       @return @c false on error */
144   virtual Bool_t Begin(Int_t event);
145   /** Process one event.  This loops over all the loaded data.   Users
146       can overload this member function, but then it's @e strongly
147       recommended to explicitly call this classes version. 
148       @return @c false on error  */
149   virtual Bool_t Event();
150   /** Called at the end of each event. 
151       @return @c false on error  */
152   virtual Bool_t End();
153   /** Called at the end of the run.
154       @return  @c false on error */
155   virtual Bool_t Finish() { return kTRUE; }
156   /** Run a full job. 
157       @return @c false on error  */
158   virtual Bool_t Run();
159
160   /** Loop over all hits, and call ProcessHit with that hit, and
161       optionally the corresponding kinematics track. 
162       @return @c false on error  */
163   virtual Bool_t ProcessHits();
164   /** Loop over all tracks, and call ProcessTrack with each hit for
165       that track
166       @return @c false on error  */
167   virtual Bool_t ProcessTracks();
168   /** Loop over all digits, and call ProcessDigit for each digit.
169       @return @c false on error  */
170   virtual Bool_t ProcessDigits();
171   /** Loop over all summable digits, and call ProcessSDigit for each
172       digit. 
173       @return @c false on error  */
174   virtual Bool_t ProcessSDigits();
175   /** Loop over all digits read from raw data files, and call
176       ProcessRawDigit for each digit. 
177       @return @c false on error  */
178   virtual Bool_t ProcessRawDigits();
179   /** Loop over all reconstructed points, and call ProcessRecPoint for
180       each reconstructed point. 
181       @return @c false on error  */
182   virtual Bool_t ProcessRecPoints();
183   /** Loop over all ESD data, and call ProcessESD for each entry.
184       @return  @c false on error  */
185   virtual Bool_t ProcessESDs();
186
187   /** Process one hit, and optionally it's corresponding kinematics
188       track.  Users should over this to process each hit. 
189       @param h Hit 
190       @param p Associated track
191       @return  @c false on error   */
192   virtual Bool_t ProcessHit(AliFMDHit* h, TParticle* p);
193   /** Process one hit per track. Users should over this to process
194       each hit. 
195       @param i Track number 
196       @param p Track  
197       @param h Associated Hit
198       @return  @c false on error   */
199   virtual Bool_t ProcessTrack(Int_t i, TParticle* p, AliFMDHit* h);
200   /** Process one digit.  Users should over this to process each
201       digit. 
202       @param digit Digit
203       @return  @c false on error   */
204   virtual Bool_t ProcessDigit(AliFMDDigit* digit);
205   /** Process one summable digit.  Users should over this to process
206       each summable digit.  
207       @param sdigit Summable digit
208       @return  @c false on error   */
209   virtual Bool_t ProcessSDigit(AliFMDSDigit* sdigit);
210   /** Process one digit from raw data files.  Users should over this
211       to process each raw digit.  
212       @param digit Raw digit
213       @return  @c false on error   */
214   virtual Bool_t ProcessRawDigit(AliFMDDigit* digit);
215   /** Process one reconstructed point.  Users should over this to
216       process each reconstructed point.  
217       @param point Reconstructed point 
218       @return  @c false on error   */
219   virtual Bool_t ProcessRecPoint(AliFMDRecPoint* point);
220   /** Process ESD data for the FMD.  Users should overload this to
221       deal with ESD data. 
222       @param d    Detector number (1-3)
223       @param r    Ring identifier ('I' or 'O')
224       @param s    Sector number (0-19, or 0-39)
225       @param t    Strip number (0-511, or 0-255)
226       @param eta  Psuedo-rapidity 
227       @param mult Psuedo-multiplicity 
228       @return  @c false on error  */
229   virtual Bool_t ProcessESD(UShort_t, Char_t, UShort_t, UShort_t, 
230                             Float_t, Float_t);
231   /** Service function to make a logarithmic axis. 
232       @param n    Number of bins 
233       @param min  Minimum of axis 
234       @param max  Maximum of axis. 
235       @return An array with the bin boundaries. */
236   static TArrayF MakeLogScale(Int_t n, Double_t min, Double_t max);
237
238   /** Set the raw data input 
239       @param file File name - if empty, assume simulated raw. */
240   void SetRawFile(const char* file) { if (file) fRawFile = file; }
241      
242 protected:
243   /** Copy ctor 
244       @param o Object to copy from  */
245   AliFMDInput(const AliFMDInput& o) 
246     : TNamed(o),
247       fGAliceFile(""),
248       fLoader(0),
249       fRun(0),
250       fStack(0),
251       fFMDLoader(0),
252       fReader(0),
253       fFMDReader(0),
254       fFMD(0),
255       fESD(0),
256       fESDEvent(0),
257       fTreeE(0),
258       fTreeH(0),
259       fTreeD(0),
260       fTreeS(0),
261       fTreeR(0),
262       fTreeA(0),
263       fChainE(0),
264       fArrayE(0),
265       fArrayH(0),
266       fArrayD(0),
267       fArrayS(0),
268       fArrayR(0),
269       fArrayA(0),
270       fHeader(0),
271       fGeoManager(0),
272       fTreeMask(0),
273       fRawFile(""),      
274       fIsInit(kFALSE),
275       fEventCount(0)
276   {}
277   /** Assignement operator 
278       @return  REference to this */
279   AliFMDInput& operator=(const AliFMDInput&) { return *this; }
280
281   TString          fGAliceFile; // File name of gAlice file
282   AliRunLoader*    fLoader;     // Loader of FMD data 
283   AliRun*          fRun;        // Run information
284   AliStack*        fStack;      // Stack of particles 
285   AliLoader*       fFMDLoader;  // Loader of FMD data 
286   AliRawReader*    fReader;     // Raw data reader 
287   AliFMDRawReader* fFMDReader;  // FMD raw reader
288   AliFMD*          fFMD;        // FMD object
289   AliESDFMD*       fESD;        // FMD ESD data  
290   AliESDEvent*     fESDEvent;   // ESD Event object. 
291   TTree*           fTreeE;      // Header tree 
292   TTree*           fTreeH;      // Hits tree
293   TTree*           fTreeD;      // Digit tree 
294   TTree*           fTreeS;      // SDigit tree 
295   TTree*           fTreeR;      // RecPoint tree
296   TTree*           fTreeA;      // Raw data tree
297   TChain*          fChainE;     // Chain of ESD's
298   TClonesArray*    fArrayE;     // Event info array
299   TClonesArray*    fArrayH;     // Hit info array
300   TClonesArray*    fArrayD;     // Digit info array
301   TClonesArray*    fArrayS;     // SDigit info array
302   TClonesArray*    fArrayR;     // Rec points info array
303   TClonesArray*    fArrayA;     // Raw data (digits) info array
304   AliHeader*       fHeader;     // Header 
305   TGeoManager*     fGeoManager; // Geometry manager
306   Int_t            fTreeMask;   // Which tree's to load
307   TString          fRawFile;    // Raw input file
308   Bool_t           fIsInit;     // Have we been initialized 
309   Int_t            fEventCount; // Event counter 
310   ClassDef(AliFMDInput,0)  //Hits for detector FMD
311 };
312
313 inline Bool_t AliFMDInput::ProcessHit(AliFMDHit*,TParticle*) { return kTRUE; }
314 inline Bool_t AliFMDInput::ProcessTrack(Int_t,TParticle*,
315                                         AliFMDHit*) { return kTRUE; }
316 inline Bool_t AliFMDInput::ProcessDigit(AliFMDDigit*) { return kTRUE; }
317 inline Bool_t AliFMDInput::ProcessSDigit(AliFMDSDigit*) { return kTRUE; }
318 inline Bool_t AliFMDInput::ProcessRawDigit(AliFMDDigit*) { return kTRUE; }
319 inline Bool_t AliFMDInput::ProcessRecPoint(AliFMDRecPoint*) { return kTRUE; }
320 inline Bool_t AliFMDInput::ProcessESD(UShort_t,Char_t,UShort_t,UShort_t,
321                                       Float_t,Float_t) { return kTRUE; }
322
323
324 #endif
325 //____________________________________________________________________
326 //
327 // Local Variables:
328 //   mode: C++
329 // End:
330 //
331 // EOF
332 //