]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDInput.h
Using detector quality flag (taken from ALICE logbook) to decide whether to rpodcue...
[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 class TSystemDirectory;
58
59 //___________________________________________________________________
60 /** @class AliFMDInput 
61     @brief Base class for reading in various FMD data. 
62     The class loops over all found events.  For each event the
63     specified data is read in.   The class then loops over all
64     elements of the read data, and process these with user defined
65     code. 
66     @code 
67     struct DigitInput : public AliFMDInput 
68     { 
69       DigitInput() 
70       { 
71          // Load digits 
72          AddLoad(kDigits);
73          // Make a histogram 
74          fHist = new TH1F("adc", "ADC spectra", 1024, -.5, 1023.5);
75       }
76       // Process one digit. 
77       Bool_t ProcessDigit(AliFMDDigit* d)
78       { 
79          fHist->Fill(d->Counts());
80          return kTRUE;
81       }
82       // After processing all events, display spectrum
83       Bool_t Finish()
84       {
85         fHist->Draw();
86       }
87       TH1F* fHist;
88     };
89
90     void AdcSpectrum()
91     { 
92        DigitInput di;
93        di.Run();
94     }
95     @endcode
96     This class allows for writing small scripts, that can be compiled
97     with AcLIC, to do all sorts of tests, quick prototyping, and so
98     on.  It has proven to be quiet useful.   One can load more than
99     one type of data in one derived class, to for example to make
100     comparisons between hits and reconstructed points.   See also the
101     various scripts in @c FMD/scripts. 
102     @ingroup FMD_util
103  */
104 class AliFMDInput : public TNamed
105 {
106 public:
107   /** The kinds of data that can be read in. */
108   enum ETrees {
109     kHits       = 1,  // Hits
110     kKinematics,      // Kinematics (from sim)
111     kDigits,          // Digits
112     kSDigits,         // Summable digits 
113     kHeader,          // Header information 
114     kRecPoints,       // Reconstructed points
115     kESD,             // Load ESD's
116     kRaw,             // Read raw data 
117     kGeometry,        // Not really a tree 
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   /** @return true if passed tree is loaded */
139   virtual Bool_t IsLoaded(ETrees tree)const { return TESTBIT(fTreeMask, tree); }
140   /** 
141    * Set the trees to load.  
142    * 
143    * @param mask Bit mask of trees to load.  Should be constructed
144    * like for example 
145    *
146    * @code 
147    * UInt_t mask = ((1 << AliFMDInput::kHits) | 
148    *                (1 << AliFMDInput::kDigits) | 
149    *                (1 << AliFMDInput::kSDigits));
150    * @endcode
151    */
152   virtual void SetLoads(UInt_t mask);
153   /** 
154    * Set the trees to load. 
155    * 
156    * @param mask A comma or space separated list of trees to load.
157    * The case is not important, and a short from of the tree name can
158    * be used.  
159    */
160   virtual void SetLoads(const char* mask);
161   /** 
162    * Get a string that represents the loaded trees 
163    * 
164    * @param dataOnly If true, then only show data 
165    * 
166    * @return String representing the loaded trees. 
167    */
168   virtual const char* LoadedString(Bool_t dataOnly=false) const;
169
170   /** Initialize the class.  If a user class overloads this member
171       function, then this @e must be explicitly called
172       @return @c false on error */
173   virtual Bool_t Init();
174   /** Callled at the beginning of each event. If a user class
175       overloads this member  function, then this @e must be explicitly
176       called. 
177       @param event Event number
178       @return @c false on error */
179   virtual Bool_t Begin(Int_t event);
180   /** Process one event.  This loops over all the loaded data.   Users
181       can overload this member function, but then it's @e strongly
182       recommended to explicitly call this classes version. 
183       @return @c false on error  */
184   virtual Bool_t Event();
185   /** Called at the end of each event. 
186       @return @c false on error  */
187   virtual Bool_t End();
188   /** Called at the end of the run.
189       @return  @c false on error */
190   virtual Bool_t Finish() { return kTRUE; }
191   /** Run a full job. 
192       @return @c false on error  */
193   virtual Bool_t Run(UInt_t maxEvents=0);
194
195   /** Loop over all hits, and call ProcessHit with that hit, and
196       optionally the corresponding kinematics track. 
197       @return @c false on error  */
198   virtual Bool_t ProcessHits();
199   /** Loop over all track refs, and call ProcessTrackRef with that hit, and
200       optionally the corresponding kinematics track. 
201       @return @c false on error  */
202   virtual Bool_t ProcessTrackRefs();
203   /** Loop over all tracks, and call ProcessTrack with each hit for
204       that track
205       @return @c false on error  */
206   virtual Bool_t ProcessTracks();
207   /** Loop over all tracks, and call ProcessTrack with each hit for
208       that track
209       @return @c false on error  */
210   virtual Bool_t ProcessStack();
211   /** Loop over all digits, and call ProcessDigit for each digit.
212       @return @c false on error  */
213   virtual Bool_t ProcessDigits();
214   /** Loop over all summable digits, and call ProcessSDigit for each
215       digit. 
216       @return @c false on error  */
217   virtual Bool_t ProcessSDigits();
218   /** Loop over all digits read from raw data files, and call
219       ProcessRawDigit for each digit. 
220       @return @c false on error  */
221   virtual Bool_t ProcessRawDigits();
222   /** Loop over all digits read from raw data files, and call
223       ProcessRawDigit for each digit. 
224       @return @c false on error  */
225   virtual Bool_t ProcessRawCalibDigits();
226   /** Loop over all reconstructed points, and call ProcessRecPoint for
227       each reconstructed point. 
228       @return @c false on error  */
229   virtual Bool_t ProcessRecPoints();
230   /** Loop over all ESD data, and call ProcessESD for each entry.
231       @return  @c false on error  */
232   virtual Bool_t ProcessESDs();
233   /** Loop over all strips and ask user routine to supply the data.
234       @return  @c false on error  */
235   virtual Bool_t ProcessUsers();
236
237   /** Process one hit, and optionally it's corresponding kinematics
238       track.  Users should over this to process each hit. 
239       @param h Hit 
240       @param p Associated track
241       @return  @c false on error   */
242   virtual Bool_t ProcessHit(AliFMDHit* h, TParticle* p);
243   /** Process one track reference, and optionally it's corresponding kinematics
244       track.  Users should overload this to process each track reference. 
245       @param trackRef Track Reference 
246       @param track Associated track
247       @return  @c false on error   */
248   virtual Bool_t ProcessTrackRef(AliTrackReference* trackRef, TParticle* track);
249   /** Process one hit per track. Users should over this to process
250       each hit. 
251       @param i Track number 
252       @param p Track  
253       @param h Associated Hit
254       @return  @c false on error   */
255   virtual Bool_t ProcessTrack(Int_t i, TParticle* p, AliFMDHit* h);
256   /** Process stack particle 
257       @param i Track number 
258       @param p Track  
259       @return  @c false on error   */
260   virtual Bool_t ProcessParticle(Int_t i , TParticle* p);
261   /** Process one digit.  Users should over this to process each
262       digit. 
263       @param digit Digit
264       @return  @c false on error   */
265   virtual Bool_t ProcessDigit(AliFMDDigit* digit);
266   /** Process one summable digit.  Users should over this to process
267       each summable digit.  
268       @param sdigit Summable digit
269       @return  @c false on error   */
270   virtual Bool_t ProcessSDigit(AliFMDSDigit* sdigit);
271   /** Process one digit from raw data files.  Users should over this
272       to process each raw digit.  
273       @param digit Raw digit
274       @return  @c false on error   */
275   virtual Bool_t ProcessRawDigit(AliFMDDigit* digit);
276   /** Process one digit from raw data files.  Users should over this
277       to process each raw digit.  
278       @param digit Raw digit
279       @return  @c false on error   */
280   virtual Bool_t ProcessRawCalibDigit(AliFMDDigit* digit);
281   /** Process one reconstructed point.  Users should over this to
282       process each reconstructed point.  
283       @param point Reconstructed point 
284       @return  @c false on error   */
285   virtual Bool_t ProcessRecPoint(AliFMDRecPoint* point);
286   /** Process ESD data for the FMD.  Users should overload this to
287       deal with ESD data. 
288       @param d    Detector number (1-3)
289       @param r    Ring identifier ('I' or 'O')
290       @param s    Sector number (0-19, or 0-39)
291       @param t    Strip number (0-511, or 0-255)
292       @param eta  Psuedo-rapidity 
293       @param mult Psuedo-multiplicity 
294       @return  @c false on error  */
295   virtual Bool_t ProcessESD(UShort_t d, Char_t r, UShort_t s, UShort_t t, 
296                             Float_t eta, Float_t mult);
297   /** Process User data for the FMD.  Users should overload this to
298       deal with ESD data. 
299       @param d    Detector number (1-3)
300       @param r    Ring identifier ('I' or 'O')
301       @param s    Sector number (0-19, or 0-39)
302       @param t    Strip number (0-511, or 0-255)
303       @param v    Value
304       @return  @c false on error  */
305   virtual Bool_t ProcessUser(UShort_t d, Char_t r, UShort_t s, UShort_t t, 
306                              Float_t v);
307   /** Service function to make a logarithmic axis. 
308       @param n    Number of bins 
309       @param min  Minimum of axis 
310       @param max  Maximum of axis. 
311       @return An array with the bin boundaries. */
312   static TArrayF MakeLogScale(Int_t n, Double_t min, Double_t max);
313
314   /** Set the raw data input 
315       @param file File name - if empty, assume simulated raw. */
316   void SetRawFile(const char* file) { if (file) fRawFile = file; }
317   void SetInputDir(const char* dir) { fInputDir = (dir && dir[0] != '\0') 
318       ? dir : "."; }
319   /** 
320    * Parse a string as a load option
321    * 
322    * @param what String to pass
323    * 
324    * @return Load option value, or 0 in case of error
325    */
326   static ETrees ParseLoad(const char* what);     
327 protected:
328   /** Copy ctor 
329       @param o Object to copy from  */
330   AliFMDInput(const AliFMDInput& o) 
331     : TNamed(o),
332       fGAliceFile(""),
333       fLoader(0),
334       fRun(0),
335       fStack(0),
336       fFMDLoader(0),
337       fReader(0),
338       fFMDReader(0),
339       fFMD(0),
340       fESD(0),
341       fESDEvent(0),
342       fTreeE(0),
343       fTreeH(0),
344       fTreeTR(0),
345       fTreeD(0),
346       fTreeS(0),
347       fTreeR(0),
348       fTreeA(0),
349       fChainE(0),
350       fArrayE(0),
351       fArrayH(0),
352       fArrayTR(0),
353       fArrayD(0),
354       fArrayS(0),
355       fArrayR(0),
356       fArrayA(0),
357       fHeader(0),
358       fGeoManager(0),
359       fTreeMask(0),
360       fRawFile(""),      
361       fInputDir("."),
362       fIsInit(kFALSE),
363       fEventCount(0),
364       fNEvents(-1)
365   {}
366   /** Assignement operator 
367       @return  REference to this */
368   AliFMDInput& operator=(const AliFMDInput&) { return *this; }
369   /** 
370    * Get user supplued data
371    * 
372    * @param d Detector
373    * @param r Ring
374    * @param s Sector
375    * @param t Strip
376    * 
377    * @return Value
378    */
379   virtual Float_t GetSignal(UShort_t d, Char_t r, UShort_t s, UShort_t t);
380
381   static const char* TreeName(ETrees tree, bool shortest=false);
382
383   /** 
384    * Make a chain of specified data 
385    * 
386    * @param what       What data to chain.  Possible values are 
387    *                   - ESD Event summary data (AliESD)
388    *                   - MC  Simulation data (galice)
389    * @param datadir    Data directory to scan 
390    * @param recursive  Whether to recurse into sub-directories 
391    * 
392    * @return Pointer to newly create chain, or null
393    */
394   static TChain* MakeChain(const char* what, const char* datadir, 
395                            bool recursive=false);
396   /** 
397    * Scan a directory (optionally recursive) for data files to add to
398    * the chain.  Only ROOT files, and files which name contain the
399    * passed pattern are considered.
400    * 
401    * @param dir        Directory to scan
402    * @param chain      Chain to add data to 
403    * @param pattern    Pattern that the file name must contain
404    * @param recursive  Whether to scan recursively 
405    */
406   static void ScanDirectory(TSystemDirectory* dir, 
407                             const TString& olddir, 
408                             TChain* chain, 
409                             const char* pattern, bool recursive);
410
411   TString          fGAliceFile; // File name of gAlice file
412   AliRunLoader*    fLoader;     // Loader of FMD data 
413   AliRun*          fRun;        // Run information
414   AliStack*        fStack;      // Stack of particles 
415   AliLoader*       fFMDLoader;  // Loader of FMD data 
416   AliRawReader*    fReader;     // Raw data reader 
417   AliFMDRawReader* fFMDReader;  // FMD raw reader
418   AliFMD*          fFMD;        // FMD object
419   AliESDFMD*       fESD;        // FMD ESD data  
420   AliESDEvent*     fESDEvent;   // ESD Event object. 
421   TTree*           fTreeE;      // Header tree 
422   TTree*           fTreeH;      // Hits tree
423   TTree*           fTreeTR;     // Track Reference tree
424   TTree*           fTreeD;      // Digit tree 
425   TTree*           fTreeS;      // SDigit tree 
426   TTree*           fTreeR;      // RecPoint tree
427   TTree*           fTreeA;      // Raw data tree
428   TChain*          fChainE;     // Chain of ESD's
429   TClonesArray*    fArrayE;     // Event info array
430   TClonesArray*    fArrayH;     // Hit info array
431   TClonesArray*    fArrayTR;    // Hit info array
432   TClonesArray*    fArrayD;     // Digit info array
433   TClonesArray*    fArrayS;     // SDigit info array
434   TClonesArray*    fArrayR;     // Rec points info array
435   TClonesArray*    fArrayA;     // Raw data (digits) info array
436   AliHeader*       fHeader;     // Header 
437   TGeoManager*     fGeoManager; // Geometry manager
438   Int_t            fTreeMask;   // Which tree's to load
439   TString          fRawFile;    // Raw input file
440   TString          fInputDir;   // Input directory
441   Bool_t           fIsInit;     // Have we been initialized 
442   Int_t            fEventCount; // Event counter 
443   Int_t            fNEvents;    // The maximum number of events
444   static const ETrees fgkAllLoads[kUser+1]; // List of all possible loads
445   ClassDef(AliFMDInput,0)  //Hits for detector FMD
446 };
447
448 inline Bool_t AliFMDInput::ProcessHit(AliFMDHit*,TParticle*) { return kTRUE; }
449 inline Bool_t AliFMDInput::ProcessTrackRef(AliTrackReference*, 
450                                            TParticle*) { return kTRUE; }
451 inline Bool_t AliFMDInput::ProcessTrack(Int_t,TParticle*,
452                                         AliFMDHit*) { return kTRUE; }
453 inline Bool_t AliFMDInput::ProcessParticle(Int_t,TParticle*) { return kTRUE; }
454 inline Bool_t AliFMDInput::ProcessDigit(AliFMDDigit*) { return kTRUE; }
455 inline Bool_t AliFMDInput::ProcessSDigit(AliFMDSDigit*) { return kTRUE; }
456 inline Bool_t AliFMDInput::ProcessRawDigit(AliFMDDigit*) { return kTRUE; }
457 inline Bool_t AliFMDInput::ProcessRawCalibDigit(AliFMDDigit*) { return kTRUE; }
458 inline Bool_t AliFMDInput::ProcessRecPoint(AliFMDRecPoint*) { return kTRUE; }
459 inline Bool_t AliFMDInput::ProcessESD(UShort_t,Char_t,UShort_t,UShort_t,
460                                       Float_t,Float_t) { return kTRUE; }
461 inline Bool_t AliFMDInput::ProcessUser(UShort_t,Char_t,UShort_t,UShort_t,
462                                        Float_t) { return kTRUE; }
463 inline Float_t AliFMDInput::GetSignal(UShort_t, Char_t, UShort_t, UShort_t) { 
464   return 0.; }
465
466
467 #endif
468 //____________________________________________________________________
469 //
470 // Local Variables:
471 //   mode: C++
472 // End:
473 //
474 // EOF
475 //