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