3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights
6 * See cxx source for full Copyright notice
8 //___________________________________________________________________
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
20 //___________________________________________________________________
21 /** @defgroup FMD_util Utility classes.
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.
56 //___________________________________________________________________
57 /** @class AliFMDInput
58 @brief Base class for reading in various FMD data.
59 The class loops over all found events. For each event the
60 specified data is read in. The class then loops over all
61 elements of the read data, and process these with user defined
64 struct DigitInput : public AliFMDInput
71 fHist = new TH1F("adc", "ADC spectra", 1024, -.5, 1023.5);
74 Bool_t ProcessDigit(AliFMDDigit* d)
76 fHist->Fill(d->Counts());
79 // After processing all events, display spectrum
93 This class allows for writing small scripts, that can be compiled
94 with AcLIC, to do all sorts of tests, quick prototyping, and so
95 on. It has proven to be quiet useful. One can load more than
96 one type of data in one derived class, to for example to make
97 comparisons between hits and reconstructed points. See also the
98 various scripts in @c FMD/scripts.
101 class AliFMDInput : public TObject
104 /** The kinds of data that can be read in. */
107 kKinematics, // Kinematics (from sim)
109 kSDigits, // Summable digits
110 kHeader, // Header information
111 kRecPoints, // Reconstructed points
113 kRaw, // Read raw data
114 kGeometry, // Not really a tree
115 kTracks // Hits and tracs - for BG study
120 @param gAliceFile galice file */
121 AliFMDInput(const char* gAliceFile);
123 virtual ~AliFMDInput() {}
125 /** Add a data type to load
126 @param tree Data to load */
127 virtual void AddLoad(ETrees tree) { SETBIT(fTreeMask, tree); }
128 /** Remove a data type to load
129 @param tree Data to @e not load */
130 virtual void RemoveLoad(ETrees tree) { CLRBIT(fTreeMask, tree); }
131 /** @return # of available events */
132 virtual Int_t NEvents() const;
134 /** Initialize the class. If a user class overloads this member
135 function, then this @e must be explicitly called
136 @return @c false on error */
137 virtual Bool_t Init();
138 /** Callled at the beginning of each event. If a user class
139 overloads this member function, then this @e must be explicitly
141 @param event Event number
142 @return @c false on error */
143 virtual Bool_t Begin(Int_t event);
144 /** Process one event. This loops over all the loaded data. Users
145 can overload this member function, but then it's @e strongly
146 recommended to explicitly call this classes version.
147 @return @c false on error */
148 virtual Bool_t Event();
149 /** Called at the end of each event.
150 @return @c false on error */
151 virtual Bool_t End();
152 /** Called at the end of the run.
153 @return @c false on error */
154 virtual Bool_t Finish() { return kTRUE; }
156 @return @c false on error */
157 virtual Bool_t Run();
159 /** Loop over all hits, and call ProcessHit with that hit, and
160 optionally the corresponding kinematics track.
161 @return @c false on error */
162 virtual Bool_t ProcessHits();
163 /** Loop over all tracks, and call ProcessTrack with each hit for
165 @return @c false on error */
166 virtual Bool_t ProcessTracks();
167 /** Loop over all digits, and call ProcessDigit for each digit.
168 @return @c false on error */
169 virtual Bool_t ProcessDigits();
170 /** Loop over all summable digits, and call ProcessSDigit for each
172 @return @c false on error */
173 virtual Bool_t ProcessSDigits();
174 /** Loop over all digits read from raw data files, and call
175 ProcessRawDigit for each digit.
176 @return @c false on error */
177 virtual Bool_t ProcessRawDigits();
178 /** Loop over all reconstructed points, and call ProcessRecPoint for
179 each reconstructed point.
180 @return @c false on error */
181 virtual Bool_t ProcessRecPoints();
182 /** Loop over all ESD data, and call ProcessESD for each entry.
183 @return @c false on error */
184 virtual Bool_t ProcessESDs();
186 /** Process one hit, and optionally it's corresponding kinematics
187 track. Users should over this to process each hit.
189 @param p Associated track
190 @return @c false on error */
191 virtual Bool_t ProcessHit(AliFMDHit* h, TParticle* p);
192 /** Process one hit per track. Users should over this to process
194 @param i Track number
196 @param h Associated Hit
197 @return @c false on error */
198 virtual Bool_t ProcessTrack(Int_t i, TParticle* p, AliFMDHit* h);
199 /** Process one digit. Users should over this to process each
202 @return @c false on error */
203 virtual Bool_t ProcessDigit(AliFMDDigit* digit);
204 /** Process one summable digit. Users should over this to process
206 @param sdigit Summable digit
207 @return @c false on error */
208 virtual Bool_t ProcessSDigit(AliFMDSDigit* sdigit);
209 /** Process one digit from raw data files. Users should over this
210 to process each raw digit.
211 @param digit Raw digit
212 @return @c false on error */
213 virtual Bool_t ProcessRawDigit(AliFMDDigit* digit);
214 /** Process one reconstructed point. Users should over this to
215 process each reconstructed point.
216 @param point Reconstructed point
217 @return @c false on error */
218 virtual Bool_t ProcessRecPoint(AliFMDRecPoint* point);
219 /** Process ESD data for the FMD. Users should overload this to
221 @param d Detector number (1-3)
222 @param r Ring identifier ('I' or 'O')
223 @param s Sector number (0-19, or 0-39)
224 @param t Strip number (0-511, or 0-255)
225 @param eta Psuedo-rapidity
226 @param mult Psuedo-multiplicity
227 @return @c false on error */
228 virtual Bool_t ProcessESD(UShort_t, Char_t, UShort_t, UShort_t,
230 /** Service function to make a logarithmic axis.
231 @param n Number of bins
232 @param min Minimum of axis
233 @param max Maximum of axis.
234 @return An array with the bin boundaries. */
235 static TArrayF MakeLogScale(Int_t n, Double_t min, Double_t max);
238 @param o Object to copy from */
239 AliFMDInput(const AliFMDInput& o)
267 /** Assignement operator
268 @return REference to this */
269 AliFMDInput& operator=(const AliFMDInput&) { return *this; }
271 TString fGAliceFile; // File name of gAlice file
272 AliRunLoader* fLoader; // Loader of FMD data
273 AliRun* fRun; // Run information
274 AliStack* fStack; // Stack of particles
275 AliLoader* fFMDLoader; // Loader of FMD data
276 AliRawReader* fReader; // Raw data reader
277 AliFMD* fFMD; // FMD object
278 AliESDFMD* fESD; // FMD ESD data
279 AliESDEvent* fESDEvent; // ESD Event object.
280 TTree* fTreeE; // Header tree
281 TTree* fTreeH; // Hits tree
282 TTree* fTreeD; // Digit tree
283 TTree* fTreeS; // SDigit tree
284 TTree* fTreeR; // RecPoint tree
285 TTree* fTreeA; // Raw data tree
286 TChain* fChainE; // Chain of ESD's
287 TClonesArray* fArrayE; // Event info array
288 TClonesArray* fArrayH; // Hit info array
289 TClonesArray* fArrayD; // Digit info array
290 TClonesArray* fArrayS; // SDigit info array
291 TClonesArray* fArrayR; // Rec points info array
292 TClonesArray* fArrayA; // Raw data (digits) info array
293 AliHeader* fHeader; // Header
294 TGeoManager* fGeoManager; // Geometry manager
295 Int_t fTreeMask; // Which tree's to load
296 Bool_t fIsInit; // Have we been initialized
297 Int_t fEventCount; // Event counter
298 ClassDef(AliFMDInput,0) //Hits for detector FMD
301 inline Bool_t AliFMDInput::ProcessHit(AliFMDHit*,TParticle*) { return kTRUE; }
302 inline Bool_t AliFMDInput::ProcessTrack(Int_t,TParticle*,
303 AliFMDHit*) { return kTRUE; }
304 inline Bool_t AliFMDInput::ProcessDigit(AliFMDDigit*) { return kTRUE; }
305 inline Bool_t AliFMDInput::ProcessSDigit(AliFMDSDigit*) { return kTRUE; }
306 inline Bool_t AliFMDInput::ProcessRawDigit(AliFMDDigit*) { return kTRUE; }
307 inline Bool_t AliFMDInput::ProcessRecPoint(AliFMDRecPoint*) { return kTRUE; }
308 inline Bool_t AliFMDInput::ProcessESD(UShort_t,Char_t,UShort_t,UShort_t,
309 Float_t,Float_t) { return kTRUE; }
313 //____________________________________________________________________