3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights
6 * See cxx source for full Copyright notice
8 //___________________________________________________________________
9 /** @defgroup FMD_util Utility classes.
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.
40 //___________________________________________________________________
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
48 struct DigitInput : public AliFMDInput
55 fHist = new TH1F("adc", "ADC spectra", 1024, -.5, 1023.5);
58 Bool_t ProcessDigit(AliFMDDigit* d)
60 fHist->Fill(d->Counts());
63 // After processing all events, display spectrum
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.
85 class AliFMDInput : public TObject
88 /** The kinds of data that can be read in. */
91 kKinematics, // Kinematics (from sim)
93 kSDigits, // Summable digits
94 kHeader, // Header information
95 kRecPoints, // Reconstructed points
97 kRaw, // Read raw data
98 kGeometry // Not really a tree
103 @param gAliceFile galice file */
104 AliFMDInput(const char* gAliceFile);
106 virtual ~AliFMDInput() {}
108 /** Add a data type to load
109 @param tree Data to load */
110 virtual void AddLoad(ETrees tree) { SETBIT(fTreeMask, tree); }
111 /** Remove a data type to load
112 @param tree Data to @e not load */
113 virtual void RemoveLoad(ETrees tree) { CLRBIT(fTreeMask, tree); }
114 /** @return # of available events */
115 virtual Int_t NEvents() const;
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 */
120 virtual Bool_t Init();
121 /** Callled at the beginning of each event. If a user class
122 overloads this member function, then this @e must be explicitly
124 @param event Event number
125 @return @c false on error */
126 virtual Bool_t Begin(Int_t event);
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 */
131 virtual Bool_t Event();
132 /** Called at the end of each event.
133 @return @c false on error */
134 virtual Bool_t End();
135 /** Called at the end of the run.
136 @return @c false on error */
137 virtual Bool_t Finish() { return kTRUE; }
139 @return @c false on error */
140 virtual Bool_t Run();
142 /** Loop over all hits, and call ProcessHit with that hit, and
143 optionally the corresponding kinematics track.
144 @return @c false on error */
145 virtual Bool_t ProcessHits();
146 /** Loop over all digits, and call ProcessDigit for each digit.
147 @return @c false on error */
148 virtual Bool_t ProcessDigits();
149 /** Loop over all summable digits, and call ProcessSDigit for each
151 @return @c false on error */
152 virtual Bool_t ProcessSDigits();
153 /** Loop over all digits read from raw data files, and call
154 ProcessRawDigit for each digit.
155 @return @c false on error */
156 virtual Bool_t ProcessRawDigits();
157 /** Loop over all reconstructed points, and call ProcessRecPoint for
158 each reconstructed point.
159 @return @c false on error */
160 virtual Bool_t ProcessRecPoints();
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 */
165 virtual Bool_t ProcessHit(AliFMDHit*, TParticle*) { return kTRUE; }
166 /** Process one digit. Users should over this to process each digit.
167 @return @c false on error */
168 virtual Bool_t ProcessDigit(AliFMDDigit*) { return kTRUE; }
169 /** Process one summable digit. Users should over this to process
171 @return @c false on error */
172 virtual Bool_t ProcessSDigit(AliFMDSDigit*) { return kTRUE; }
173 /** Process one digit from raw data files. Users should over this
174 to process each raw digit.
175 @return @c false on error */
176 virtual Bool_t ProcessRawDigit(AliFMDDigit*) { return kTRUE; }
177 /** Process one reconstructed point. Users should over this to
178 process each reconstructed point.
179 @return @c false on error */
180 virtual Bool_t ProcessRecPoint(AliFMDRecPoint*) { return kTRUE; }
181 /** Process ESD data for the FMD. Users should overload this to
183 @return @c false on error */
184 virtual Bool_t ProcessESD(AliESDFMD*) { return kTRUE; }
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
192 AliRawReader* fReader; // Raw data reader
193 AliFMD* fFMD; // FMD object
194 AliESD* fMainESD; // ESD Object
195 AliESDFMD* fESD; // FMD ESD data
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
201 TTree* fTreeA; // Raw data tree
202 TChain* fChainE; // Chain of ESD's
203 TClonesArray* fArrayE; // Event info array
204 TClonesArray* fArrayH; // Hit info array
205 TClonesArray* fArrayD; // Digit info array
206 TClonesArray* fArrayS; // SDigit info array
207 TClonesArray* fArrayR; // Rec points info array
208 TClonesArray* fArrayA; // Raw data (digits) info array
209 TGeoManager* fGeoManager; // Geometry manager
210 Int_t fTreeMask; // Which tree's to load
212 ClassDef(AliFMDInput,0) //Hits for detector FMD
216 //____________________________________________________________________
218 class AliFMDInputHits : public AliFMDInput
221 AliFMDInputHits(const char* file="galice.root")
222 : AliFMDInput(file) { AddLoad(kHits); }
223 ClassDef(AliFMDInputHits, 0);
226 //____________________________________________________________________
228 class AliFMDInputDigits : public AliFMDInput
231 AliFMDInputDigits(const char* file="galice.root")
232 : AliFMDInput(file) { AddLoad(kDigits); }
233 ClassDef(AliFMDInputDigits, 0);
236 //____________________________________________________________________
238 class AliFMDInputSDigits : public AliFMDInput
241 AliFMDInputSDigits(const char* file="galice.root")
242 : AliFMDInput(file) { AddLoad(kSDigits); }
243 ClassDef(AliFMDInputSDigits, 0);
246 //____________________________________________________________________
247 class AliFMDInputRaw : public AliFMDInput
250 AliFMDInputRaw(const char* file="galice.root")
251 : AliFMDInput(file) { AddLoad(kRaw); }
252 ClassDef(AliFMDInputRaw, 0);
255 //____________________________________________________________________
256 class AliFMDRecPoint;
257 class AliFMDInputRecPoints : public AliFMDInput
260 AliFMDInputRecPoints(const char* file="galice.root")
261 : AliFMDInput(file) { AddLoad(kRecPoints); }
262 ClassDef(AliFMDInputRecPoints, 0);
266 //____________________________________________________________________