]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FMD/AliFMDInput.h
Additiona TOF data members (S.Arcelli)
[u/mrichter/AliRoot.git] / FMD / AliFMDInput.h
CommitLineData
a1f80595 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//___________________________________________________________________
9f662337 9/** @defgroup FMD_util Utility classes.
10
11 The classes defined here, are utility classes for reading in data
12 for the FMD. They are put in a seperate library to not polute the
13 normal libraries. The classes are intended to be used as base
14 classes for customized class that do some sort of analysis on the
15 various types of data produced by the FMD.
16*/
8f6ee336 17#include <TObject.h>
a1f80595 18#ifndef ROOT_TString
19# include <TString.h>
20#endif
8f6ee336 21class AliRunLoader;
22class AliLoader;
23class AliStack;
24class AliRun;
d760ea03 25class AliRawReader;
8f6ee336 26class AliFMD;
27class AliFMDHit;
bf000c32 28class AliFMDDigit;
29class AliFMDSDigit;
30class AliFMDRecPoint;
31class AliESD;
32class AliESDFMD;
8f6ee336 33class TString;
34class TClonesArray;
35class TTree;
36class TGeoManager;
37class TParticle;
bf000c32 38class TChain;
a1f80595 39
40//___________________________________________________________________
9f662337 41/** @class AliFMDInput
42 @brief Base class for reading in various FMD data.
43 The class loops over all found events. For each event the
44 specified data is read in. The class then loops over all
45 elements of the read data, and process these with user defined
46 code.
47 @code
48 struct DigitInput : public AliFMDInput
49 {
50 DigitInput()
51 {
52 // Load digits
53 AddLoad(kDigits);
54 // Make a histogram
55 fHist = new TH1F("adc", "ADC spectra", 1024, -.5, 1023.5);
56 }
57 // Process one digit.
58 Bool_t ProcessDigit(AliFMDDigit* d)
59 {
60 fHist->Fill(d->Counts());
61 return kTRUE;
62 }
63 // After processing all events, display spectrum
64 Bool_t Finish()
65 {
66 fHist->Draw();
67 }
68 TH1F* fHist;
69 };
70
71 void AdcSpectrum()
72 {
73 DigitInput di;
74 di.Run();
75 }
76 @endcode
77 This class allows for writing small scripts, that can be compiled
78 with AcLIC, to do all sorts of tests, quick prototyping, and so
79 on. It has proven to be quiet useful. One can load more than
80 one type of data in one derived class, to for example to make
81 comparisons between hits and reconstructed points. See also the
82 various scripts in @c FMD/scripts.
83 @ingroup FMD_util
84 */
a1f80595 85class AliFMDInput : public TObject
86{
87public:
9f662337 88 /** The kinds of data that can be read in. */
a1f80595 89 enum ETrees {
8f6ee336 90 kHits = 1, // Hits
91 kKinematics, // Kinematics (from sim)
92 kDigits, // Digits
93 kSDigits, // Summable digits
94 kHeader, // Header information
95 kRecPoints, // Reconstructed points
bf000c32 96 kESD, // Load ESD's
d760ea03 97 kRaw, // Read raw data
8f6ee336 98 kGeometry // Not really a tree
a1f80595 99 };
9f662337 100 /** CTOR */
a1f80595 101 AliFMDInput();
9f662337 102 /** CTOR
103 @param gAliceFile galice file */
a1f80595 104 AliFMDInput(const char* gAliceFile);
9f662337 105 /** DTOR */
a1f80595 106 virtual ~AliFMDInput() {}
107
9f662337 108 /** Add a data type to load
109 @param tree Data to load */
a1f80595 110 virtual void AddLoad(ETrees tree) { SETBIT(fTreeMask, tree); }
9f662337 111 /** Remove a data type to load
112 @param tree Data to @e not load */
a1f80595 113 virtual void RemoveLoad(ETrees tree) { CLRBIT(fTreeMask, tree); }
9f662337 114 /** @return # of available events */
a1f80595 115 virtual Int_t NEvents() const;
116
9f662337 117 /** Initialize the class. If a user class overloads this member
118 function, then this @e must be explicitly called
119 @return @c false on error */
a1f80595 120 virtual Bool_t Init();
9f662337 121 /** Callled at the beginning of each event. If a user class
122 overloads this member function, then this @e must be explicitly
123 called.
124 @param event Event number
125 @return @c false on error */
a1f80595 126 virtual Bool_t Begin(Int_t event);
9f662337 127 /** Process one event. This loops over all the loaded data. Users
128 can overload this member function, but then it's @e strongly
129 recommended to explicitly call this classes version.
130 @return @c false on error */
bf000c32 131 virtual Bool_t Event();
9f662337 132 /** Called at the end of each event.
133 @return @c false on error */
a1f80595 134 virtual Bool_t End();
9f662337 135 /** Called at the end of the run.
136 @return @c false on error */
a1f80595 137 virtual Bool_t Finish() { return kTRUE; }
9f662337 138 /** Run a full job.
139 @return @c false on error */
a1f80595 140 virtual Bool_t Run();
bf000c32 141
9f662337 142 /** Loop over all hits, and call ProcessHit with that hit, and
143 optionally the corresponding kinematics track.
144 @return @c false on error */
bf000c32 145 virtual Bool_t ProcessHits();
9f662337 146 /** Loop over all digits, and call ProcessDigit for each digit.
147 @return @c false on error */
bf000c32 148 virtual Bool_t ProcessDigits();
9f662337 149 /** Loop over all summable digits, and call ProcessSDigit for each
150 digit.
151 @return @c false on error */
bf000c32 152 virtual Bool_t ProcessSDigits();
9f662337 153 /** Loop over all digits read from raw data files, and call
154 ProcessRawDigit for each digit.
155 @return @c false on error */
d760ea03 156 virtual Bool_t ProcessRawDigits();
9f662337 157 /** Loop over all reconstructed points, and call ProcessRecPoint for
158 each reconstructed point.
159 @return @c false on error */
bf000c32 160 virtual Bool_t ProcessRecPoints();
161
9f662337 162 /** Process one hit, and optionally it's corresponding kinematics
163 track. Users should over this to process each hit.
164 @return @c false on error */
bf000c32 165 virtual Bool_t ProcessHit(AliFMDHit*, TParticle*) { return kTRUE; }
9f662337 166 /** Process one digit. Users should over this to process each digit.
167 @return @c false on error */
bf000c32 168 virtual Bool_t ProcessDigit(AliFMDDigit*) { return kTRUE; }
9f662337 169 /** Process one summable digit. Users should over this to process
170 each summable digit.
171 @return @c false on error */
bf000c32 172 virtual Bool_t ProcessSDigit(AliFMDSDigit*) { return kTRUE; }
9f662337 173 /** Process one digit from raw data files. Users should over this
174 to process each raw digit.
175 @return @c false on error */
d760ea03 176 virtual Bool_t ProcessRawDigit(AliFMDDigit*) { return kTRUE; }
9f662337 177 /** Process one reconstructed point. Users should over this to
178 process each reconstructed point.
179 @return @c false on error */
bf000c32 180 virtual Bool_t ProcessRecPoint(AliFMDRecPoint*) { return kTRUE; }
9f662337 181 /** Process ESD data for the FMD. Users should overload this to
182 deal with ESD data.
183 @return @c false on error */
bf000c32 184 virtual Bool_t ProcessESD(AliESDFMD*) { return kTRUE; }
185
a1f80595 186protected:
187 TString fGAliceFile; // File name of gAlice file
188 AliRunLoader* fLoader; // Loader of FMD data
189 AliRun* fRun; // Run information
190 AliStack* fStack; // Stack of particles
191 AliLoader* fFMDLoader; // Loader of FMD data
d760ea03 192 AliRawReader* fReader; // Raw data reader
a1f80595 193 AliFMD* fFMD; // FMD object
bf000c32 194 AliESD* fMainESD; // ESD Object
195 AliESDFMD* fESD; // FMD ESD data
a1f80595 196 TTree* fTreeE; // Header tree
197 TTree* fTreeH; // Hits tree
198 TTree* fTreeD; // Digit tree
199 TTree* fTreeS; // SDigit tree
200 TTree* fTreeR; // RecPoint tree
d760ea03 201 TTree* fTreeA; // Raw data tree
bf000c32 202 TChain* fChainE; // Chain of ESD's
a1f80595 203 TClonesArray* fArrayE; // Event info array
204 TClonesArray* fArrayH; // Hit info array
205 TClonesArray* fArrayD; // Digit info array
206 TClonesArray* fArrayS; // SDigit info array
d760ea03 207 TClonesArray* fArrayR; // Rec points info array
208 TClonesArray* fArrayA; // Raw data (digits) info array
8f6ee336 209 TGeoManager* fGeoManager; // Geometry manager
a1f80595 210 Int_t fTreeMask; // Which tree's to load
211 Bool_t fIsInit;
212 ClassDef(AliFMDInput,0) //Hits for detector FMD
213};
214
215
216//____________________________________________________________________
217class AliFMDHit;
218class AliFMDInputHits : public AliFMDInput
219{
220public:
221 AliFMDInputHits(const char* file="galice.root")
222 : AliFMDInput(file) { AddLoad(kHits); }
a1f80595 223 ClassDef(AliFMDInputHits, 0);
224};
225
226//____________________________________________________________________
227class AliFMDDigit;
228class AliFMDInputDigits : public AliFMDInput
229{
230public:
231 AliFMDInputDigits(const char* file="galice.root")
232 : AliFMDInput(file) { AddLoad(kDigits); }
a1f80595 233 ClassDef(AliFMDInputDigits, 0);
234};
235
236//____________________________________________________________________
237class AliFMDSDigit;
238class AliFMDInputSDigits : public AliFMDInput
239{
240public:
241 AliFMDInputSDigits(const char* file="galice.root")
242 : AliFMDInput(file) { AddLoad(kSDigits); }
a1f80595 243 ClassDef(AliFMDInputSDigits, 0);
244};
245
d760ea03 246//____________________________________________________________________
247class AliFMDInputRaw : public AliFMDInput
248{
249public:
250 AliFMDInputRaw(const char* file="galice.root")
251 : AliFMDInput(file) { AddLoad(kRaw); }
252 ClassDef(AliFMDInputRaw, 0);
253};
254
a1f80595 255//____________________________________________________________________
bf000c32 256class AliFMDRecPoint;
a1f80595 257class AliFMDInputRecPoints : public AliFMDInput
258{
259public:
260 AliFMDInputRecPoints(const char* file="galice.root")
261 : AliFMDInput(file) { AddLoad(kRecPoints); }
a1f80595 262 ClassDef(AliFMDInputRecPoints, 0);
263};
264
265#endif
266//____________________________________________________________________
267//
268// Local Variables:
269// mode: C++
270// End:
271//
272// EOF
273//