]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDInput.h
New EMCAL cosmic trigger as defined by Federico A.
[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     kRawCalib,        // Read raws and calibrate them
120     kUser
121   };
122   /** CTOR  */
123   AliFMDInput();
124   /** CTOR
125       @param gAliceFile galice file  */
126   AliFMDInput(const char* gAliceFile);
127   /** DTOR */
128   virtual ~AliFMDInput() {}
129
130   /** Add a data type to load 
131       @param tree Data to load */
132   virtual void   AddLoad(ETrees tree)     { SETBIT(fTreeMask, tree); }
133   /** Remove a data type to load 
134       @param tree Data to @e not load */
135   virtual void   RemoveLoad(ETrees tree)  { CLRBIT(fTreeMask, tree); }
136   /** @return # of available events */
137   virtual Int_t  NEvents() const;
138
139   /** Initialize the class.  If a user class overloads this member
140       function, then this @e must be explicitly called
141       @return @c false on error */
142   virtual Bool_t Init();
143   /** Callled at the beginning of each event. If a user class
144       overloads this member  function, then this @e must be explicitly
145       called. 
146       @param event Event number
147       @return @c false on error */
148   virtual Bool_t Begin(Int_t event);
149   /** Process one event.  This loops over all the loaded data.   Users
150       can overload this member function, but then it's @e strongly
151       recommended to explicitly call this classes version. 
152       @return @c false on error  */
153   virtual Bool_t Event();
154   /** Called at the end of each event. 
155       @return @c false on error  */
156   virtual Bool_t End();
157   /** Called at the end of the run.
158       @return  @c false on error */
159   virtual Bool_t Finish() { return kTRUE; }
160   /** Run a full job. 
161       @return @c false on error  */
162   virtual Bool_t Run();
163
164   /** Loop over all hits, and call ProcessHit with that hit, and
165       optionally the corresponding kinematics track. 
166       @return @c false on error  */
167   virtual Bool_t ProcessHits();
168   /** Loop over all track refs, and call ProcessTrackRef with that hit, and
169       optionally the corresponding kinematics track. 
170       @return @c false on error  */
171   virtual Bool_t ProcessTrackRefs();
172   /** Loop over all tracks, and call ProcessTrack with each hit for
173       that track
174       @return @c false on error  */
175   virtual Bool_t ProcessTracks();
176   /** Loop over all digits, and call ProcessDigit for each digit.
177       @return @c false on error  */
178   virtual Bool_t ProcessDigits();
179   /** Loop over all summable digits, and call ProcessSDigit for each
180       digit. 
181       @return @c false on error  */
182   virtual Bool_t ProcessSDigits();
183   /** Loop over all digits read from raw data files, and call
184       ProcessRawDigit for each digit. 
185       @return @c false on error  */
186   virtual Bool_t ProcessRawDigits();
187   /** Loop over all digits read from raw data files, and call
188       ProcessRawDigit for each digit. 
189       @return @c false on error  */
190   virtual Bool_t ProcessRawCalibDigits();
191   /** Loop over all reconstructed points, and call ProcessRecPoint for
192       each reconstructed point. 
193       @return @c false on error  */
194   virtual Bool_t ProcessRecPoints();
195   /** Loop over all ESD data, and call ProcessESD for each entry.
196       @return  @c false on error  */
197   virtual Bool_t ProcessESDs();
198   /** Loop over all strips and ask user routine to supply the data.
199       @return  @c false on error  */
200   virtual Bool_t ProcessUsers();
201
202   /** Process one hit, and optionally it's corresponding kinematics
203       track.  Users should over this to process each hit. 
204       @param h Hit 
205       @param p Associated track
206       @return  @c false on error   */
207   virtual Bool_t ProcessHit(AliFMDHit* h, TParticle* p);
208   /** Process one track reference, and optionally it's corresponding kinematics
209       track.  Users should overload this to process each track reference. 
210       @param trackRef Track Reference 
211       @param track Associated track
212       @return  @c false on error   */
213   virtual Bool_t ProcessTrackRef(AliTrackReference* trackRef, TParticle* track);
214   /** Process one hit per track. Users should over this to process
215       each hit. 
216       @param i Track number 
217       @param p Track  
218       @param h Associated Hit
219       @return  @c false on error   */
220   virtual Bool_t ProcessTrack(Int_t i, TParticle* p, AliFMDHit* h);
221   /** Process one digit.  Users should over this to process each
222       digit. 
223       @param digit Digit
224       @return  @c false on error   */
225   virtual Bool_t ProcessDigit(AliFMDDigit* digit);
226   /** Process one summable digit.  Users should over this to process
227       each summable digit.  
228       @param sdigit Summable digit
229       @return  @c false on error   */
230   virtual Bool_t ProcessSDigit(AliFMDSDigit* sdigit);
231   /** Process one digit from raw data files.  Users should over this
232       to process each raw digit.  
233       @param digit Raw digit
234       @return  @c false on error   */
235   virtual Bool_t ProcessRawDigit(AliFMDDigit* digit);
236   /** Process one digit from raw data files.  Users should over this
237       to process each raw digit.  
238       @param digit Raw digit
239       @return  @c false on error   */
240   virtual Bool_t ProcessRawCalibDigit(AliFMDDigit* digit);
241   /** Process one reconstructed point.  Users should over this to
242       process each reconstructed point.  
243       @param point Reconstructed point 
244       @return  @c false on error   */
245   virtual Bool_t ProcessRecPoint(AliFMDRecPoint* point);
246   /** Process ESD data for the FMD.  Users should overload this to
247       deal with ESD data. 
248       @param d    Detector number (1-3)
249       @param r    Ring identifier ('I' or 'O')
250       @param s    Sector number (0-19, or 0-39)
251       @param t    Strip number (0-511, or 0-255)
252       @param eta  Psuedo-rapidity 
253       @param mult Psuedo-multiplicity 
254       @return  @c false on error  */
255   virtual Bool_t ProcessESD(UShort_t d, Char_t r, UShort_t s, UShort_t t, 
256                             Float_t eta, Float_t mult);
257   /** Process User data for the FMD.  Users should overload this to
258       deal with ESD data. 
259       @param d    Detector number (1-3)
260       @param r    Ring identifier ('I' or 'O')
261       @param s    Sector number (0-19, or 0-39)
262       @param t    Strip number (0-511, or 0-255)
263       @param v    Value
264       @return  @c false on error  */
265   virtual Bool_t ProcessUser(UShort_t d, Char_t r, UShort_t s, UShort_t t, 
266                              Float_t v);
267   /** Service function to make a logarithmic axis. 
268       @param n    Number of bins 
269       @param min  Minimum of axis 
270       @param max  Maximum of axis. 
271       @return An array with the bin boundaries. */
272   static TArrayF MakeLogScale(Int_t n, Double_t min, Double_t max);
273
274   /** Set the raw data input 
275       @param file File name - if empty, assume simulated raw. */
276   void SetRawFile(const char* file) { if (file) fRawFile = file; }
277      
278 protected:
279   /** Copy ctor 
280       @param o Object to copy from  */
281   AliFMDInput(const AliFMDInput& o) 
282     : TNamed(o),
283       fGAliceFile(""),
284       fLoader(0),
285       fRun(0),
286       fStack(0),
287       fFMDLoader(0),
288       fReader(0),
289       fFMDReader(0),
290       fFMD(0),
291       fESD(0),
292       fESDEvent(0),
293       fTreeE(0),
294       fTreeH(0),
295       fTreeTR(0),
296       fTreeD(0),
297       fTreeS(0),
298       fTreeR(0),
299       fTreeA(0),
300       fChainE(0),
301       fArrayE(0),
302       fArrayH(0),
303       fArrayTR(0),
304       fArrayD(0),
305       fArrayS(0),
306       fArrayR(0),
307       fArrayA(0),
308       fHeader(0),
309       fGeoManager(0),
310       fTreeMask(0),
311       fRawFile(""),      
312       fIsInit(kFALSE),
313       fEventCount(0)
314   {}
315   /** Assignement operator 
316       @return  REference to this */
317   AliFMDInput& operator=(const AliFMDInput&) { return *this; }
318   /** 
319    * Get user supplued data
320    * 
321    * @param d Detector
322    * @param r Ring
323    * @param s Sector
324    * @param t Strip
325    * 
326    * @return Value
327    */
328   virtual Float_t GetSignal(UShort_t d, Char_t r, UShort_t s, UShort_t t);
329
330   TString          fGAliceFile; // File name of gAlice file
331   AliRunLoader*    fLoader;     // Loader of FMD data 
332   AliRun*          fRun;        // Run information
333   AliStack*        fStack;      // Stack of particles 
334   AliLoader*       fFMDLoader;  // Loader of FMD data 
335   AliRawReader*    fReader;     // Raw data reader 
336   AliFMDRawReader* fFMDReader;  // FMD raw reader
337   AliFMD*          fFMD;        // FMD object
338   AliESDFMD*       fESD;        // FMD ESD data  
339   AliESDEvent*     fESDEvent;   // ESD Event object. 
340   TTree*           fTreeE;      // Header tree 
341   TTree*           fTreeH;      // Hits tree
342   TTree*           fTreeTR;     // Track Reference tree
343   TTree*           fTreeD;      // Digit tree 
344   TTree*           fTreeS;      // SDigit tree 
345   TTree*           fTreeR;      // RecPoint tree
346   TTree*           fTreeA;      // Raw data tree
347   TChain*          fChainE;     // Chain of ESD's
348   TClonesArray*    fArrayE;     // Event info array
349   TClonesArray*    fArrayH;     // Hit info array
350   TClonesArray*    fArrayTR;    // Hit info array
351   TClonesArray*    fArrayD;     // Digit info array
352   TClonesArray*    fArrayS;     // SDigit info array
353   TClonesArray*    fArrayR;     // Rec points info array
354   TClonesArray*    fArrayA;     // Raw data (digits) info array
355   AliHeader*       fHeader;     // Header 
356   TGeoManager*     fGeoManager; // Geometry manager
357   Int_t            fTreeMask;   // Which tree's to load
358   TString          fRawFile;    // Raw input file
359   Bool_t           fIsInit;     // Have we been initialized 
360   Int_t            fEventCount; // Event counter 
361   ClassDef(AliFMDInput,0)  //Hits for detector FMD
362 };
363
364 inline Bool_t AliFMDInput::ProcessHit(AliFMDHit*,TParticle*) { return kTRUE; }
365 inline Bool_t AliFMDInput::ProcessTrackRef(AliTrackReference*, 
366                                            TParticle*) { return kTRUE; }
367 inline Bool_t AliFMDInput::ProcessTrack(Int_t,TParticle*,
368                                         AliFMDHit*) { return kTRUE; }
369 inline Bool_t AliFMDInput::ProcessDigit(AliFMDDigit*) { return kTRUE; }
370 inline Bool_t AliFMDInput::ProcessSDigit(AliFMDSDigit*) { return kTRUE; }
371 inline Bool_t AliFMDInput::ProcessRawDigit(AliFMDDigit*) { return kTRUE; }
372 inline Bool_t AliFMDInput::ProcessRawCalibDigit(AliFMDDigit*) { return kTRUE; }
373 inline Bool_t AliFMDInput::ProcessRecPoint(AliFMDRecPoint*) { return kTRUE; }
374 inline Bool_t AliFMDInput::ProcessESD(UShort_t,Char_t,UShort_t,UShort_t,
375                                       Float_t,Float_t) { return kTRUE; }
376 inline Bool_t AliFMDInput::ProcessUser(UShort_t,Char_t,UShort_t,UShort_t,
377                                        Float_t) { return kTRUE; }
378 inline Float_t AliFMDInput::GetSignal(UShort_t, Char_t, UShort_t, UShort_t) { 
379   return 0.; }
380
381
382 #endif
383 //____________________________________________________________________
384 //
385 // Local Variables:
386 //   mode: C++
387 // End:
388 //
389 // EOF
390 //