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