]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDInput.h
fb27fc65b64d6b534d39cf7ae946f74095399c27
[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 AliESDEvent;
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     kTracks           // Hits and tracs - for BG study  
112   };
113   /** CTOR  */
114   AliFMDInput();
115   /** CTOR
116       @param gAliceFile galice file  */
117   AliFMDInput(const char* gAliceFile);
118   /** DTOR */
119   virtual ~AliFMDInput() {}
120
121   /** Add a data type to load 
122       @param tree Data to load */
123   virtual void   AddLoad(ETrees tree)     { SETBIT(fTreeMask, tree); }
124   /** Remove a data type to load 
125       @param tree Data to @e not load */
126   virtual void   RemoveLoad(ETrees tree)  { CLRBIT(fTreeMask, tree); }
127   /** @return # of available events */
128   virtual Int_t  NEvents() const;
129
130   /** Initialize the class.  If a user class overloads this member
131       function, then this @e must be explicitly called
132       @return @c false on error */
133   virtual Bool_t Init();
134   /** Callled at the beginning of each event. If a user class
135       overloads this member  function, then this @e must be explicitly
136       called. 
137       @param event Event number
138       @return @c false on error */
139   virtual Bool_t Begin(Int_t event);
140   /** Process one event.  This loops over all the loaded data.   Users
141       can overload this member function, but then it's @e strongly
142       recommended to explicitly call this classes version. 
143       @return @c false on error  */
144   virtual Bool_t Event();
145   /** Called at the end of each event. 
146       @return @c false on error  */
147   virtual Bool_t End();
148   /** Called at the end of the run.
149       @return  @c false on error */
150   virtual Bool_t Finish() { return kTRUE; }
151   /** Run a full job. 
152       @return @c false on error  */
153   virtual Bool_t Run();
154
155   /** Loop over all hits, and call ProcessHit with that hit, and
156       optionally the corresponding kinematics track. 
157       @return @c false on error  */
158   virtual Bool_t ProcessHits();
159   /** Loop over all tracks, and call ProcessTrack with each hit for
160       that track
161       @return @c false on error  */
162   virtual Bool_t ProcessTracks();
163   /** Loop over all digits, and call ProcessDigit for each digit.
164       @return @c false on error  */
165   virtual Bool_t ProcessDigits();
166   /** Loop over all summable digits, and call ProcessSDigit for each
167       digit. 
168       @return @c false on error  */
169   virtual Bool_t ProcessSDigits();
170   /** Loop over all digits read from raw data files, and call
171       ProcessRawDigit for each digit. 
172       @return @c false on error  */
173   virtual Bool_t ProcessRawDigits();
174   /** Loop over all reconstructed points, and call ProcessRecPoint for
175       each reconstructed point. 
176       @return @c false on error  */
177   virtual Bool_t ProcessRecPoints();
178   /** Loop over all ESD data, and call ProcessESD for each entry.
179       @return  @c false on error  */
180   virtual Bool_t ProcessESDs();
181
182   /** Process one hit, and optionally it's corresponding kinematics
183       track.  Users should over this to process each hit. 
184       @param h Hit 
185       @param p Associated track
186       @return  @c false on error   */
187   virtual Bool_t ProcessHit(AliFMDHit* h, TParticle* p);
188   /** Process one hit per track. Users should over this to process
189       each hit. 
190       @param i Track number 
191       @param p Track  
192       @param h Associated Hit
193       @return  @c false on error   */
194   virtual Bool_t ProcessTrack(Int_t i, TParticle* p, AliFMDHit* h);
195   /** Process one digit.  Users should over this to process each
196       digit. 
197       @param digit Digit
198       @return  @c false on error   */
199   virtual Bool_t ProcessDigit(AliFMDDigit* digit);
200   /** Process one summable digit.  Users should over this to process
201       each summable digit.  
202       @param sdigit Summable digit
203       @return  @c false on error   */
204   virtual Bool_t ProcessSDigit(AliFMDSDigit* sdigit);
205   /** Process one digit from raw data files.  Users should over this
206       to process each raw digit.  
207       @param digit Raw digit
208       @return  @c false on error   */
209   virtual Bool_t ProcessRawDigit(AliFMDDigit* digit);
210   /** Process one reconstructed point.  Users should over this to
211       process each reconstructed point.  
212       @param point Reconstructed point 
213       @return  @c false on error   */
214   virtual Bool_t ProcessRecPoint(AliFMDRecPoint* point);
215   /** Process ESD data for the FMD.  Users should overload this to
216       deal with ESD data. 
217       @param d    Detector number (1-3)
218       @param r    Ring identifier ('I' or 'O')
219       @param s    Sector number (0-19, or 0-39)
220       @param t    Strip number (0-511, or 0-255)
221       @param eta  Psuedo-rapidity 
222       @param mult Psuedo-multiplicity 
223       @return  @c false on error  */
224   virtual Bool_t ProcessESD(UShort_t, Char_t, UShort_t, UShort_t, 
225                             Float_t, Float_t);
226   /** Service function to make a logarithmic axis. 
227       @param n    Number of bins 
228       @param min  Minimum of axis 
229       @param max  Maximum of axis. 
230       @return An array with the bin boundaries. */
231   static TArrayF MakeLogScale(Int_t n, Double_t min, Double_t max);
232 protected:
233   /** Copy ctor 
234       @param o Object to copy from  */
235   AliFMDInput(const AliFMDInput& o) 
236     : TObject(o),
237       fGAliceFile(""),
238       fLoader(0),
239       fRun(0),
240       fStack(0),
241       fFMDLoader(0),
242       fReader(0),
243       fFMD(0),
244       fESD(0),
245       fESDEvent(0),
246       fTreeE(0),
247       fTreeH(0),
248       fTreeD(0),
249       fTreeS(0),
250       fTreeR(0),
251       fTreeA(0),
252       fChainE(0),
253       fArrayE(0),
254       fArrayH(0),
255       fArrayD(0),
256       fArrayS(0),
257       fArrayR(0),
258       fArrayA(0),
259       fGeoManager(0),
260       fTreeMask(0),
261       fIsInit(kFALSE)
262   {}
263   /** Assignement operator 
264       @return  REference to this */
265   AliFMDInput& operator=(const AliFMDInput&) { return *this; }
266
267   TString       fGAliceFile; // File name of gAlice file
268   AliRunLoader* fLoader;     // Loader of FMD data 
269   AliRun*       fRun;        // Run information
270   AliStack*     fStack;      // Stack of particles 
271   AliLoader*    fFMDLoader;  // Loader of FMD data 
272   AliRawReader* fReader;     // Raw data reader 
273   AliFMD*       fFMD;        // FMD object
274   AliESDFMD*    fESD;        // FMD ESD data  
275   AliESDEvent*  fESDEvent;   // ESD Event object. 
276   TTree*        fTreeE;      // Header tree 
277   TTree*        fTreeH;      // Hits tree
278   TTree*        fTreeD;      // Digit tree 
279   TTree*        fTreeS;      // SDigit tree 
280   TTree*        fTreeR;      // RecPoint tree
281   TTree*        fTreeA;      // Raw data tree
282   TChain*       fChainE;     // Chain of ESD's
283   TClonesArray* fArrayE;     // Event info array
284   TClonesArray* fArrayH;     // Hit info array
285   TClonesArray* fArrayD;     // Digit info array
286   TClonesArray* fArrayS;     // SDigit info array
287   TClonesArray* fArrayR;     // Rec points info array
288   TClonesArray* fArrayA;     // Raw data (digits) info array
289   TGeoManager*  fGeoManager; // Geometry manager
290   Int_t         fTreeMask;   // Which tree's to load
291   Bool_t        fIsInit;     // Have we been initialized 
292   Int_t         fEventCount; // Event counter 
293   ClassDef(AliFMDInput,0)  //Hits for detector FMD
294 };
295
296 inline Bool_t AliFMDInput::ProcessHit(AliFMDHit*,TParticle*) { return kTRUE; }
297 inline Bool_t AliFMDInput::ProcessTrack(Int_t,TParticle*,
298                                         AliFMDHit*) { return kTRUE; }
299 inline Bool_t AliFMDInput::ProcessDigit(AliFMDDigit*) { return kTRUE; }
300 inline Bool_t AliFMDInput::ProcessSDigit(AliFMDSDigit*) { return kTRUE; }
301 inline Bool_t AliFMDInput::ProcessRawDigit(AliFMDDigit*) { return kTRUE; }
302 inline Bool_t AliFMDInput::ProcessRecPoint(AliFMDRecPoint*) { return kTRUE; }
303 inline Bool_t AliFMDInput::ProcessESD(UShort_t,Char_t,UShort_t,UShort_t,
304                                       Float_t,Float_t) { return kTRUE; }
305
306
307 #endif
308 //____________________________________________________________________
309 //
310 // Local Variables:
311 //   mode: C++
312 // End:
313 //
314 // EOF
315 //