]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FMD/AliFMDInput.h
Fixed compiler warnings
[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 */
02a27b50 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.
c2fc1258 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*/
a1f80595 20//___________________________________________________________________
9f662337 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*/
8f6ee336 29#include <TObject.h>
a1f80595 30#ifndef ROOT_TString
31# include <TString.h>
32#endif
7003fbd0 33#ifndef ROOT_TArrayF
34# include <TArrayF.h>
35#endif
8f6ee336 36class AliRunLoader;
37class AliLoader;
38class AliStack;
39class AliRun;
d760ea03 40class AliRawReader;
8f6ee336 41class AliFMD;
42class AliFMDHit;
bf000c32 43class AliFMDDigit;
44class AliFMDSDigit;
45class AliFMDRecPoint;
df137876 46class AliESDEvent;
bf000c32 47class AliESDFMD;
9b98d361 48class AliHeader;
8f6ee336 49class TString;
50class TClonesArray;
51class TTree;
52class TGeoManager;
53class TParticle;
bf000c32 54class TChain;
a1f80595 55
56//___________________________________________________________________
9f662337 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
62 code.
63 @code
64 struct DigitInput : public AliFMDInput
65 {
66 DigitInput()
67 {
68 // Load digits
69 AddLoad(kDigits);
70 // Make a histogram
71 fHist = new TH1F("adc", "ADC spectra", 1024, -.5, 1023.5);
72 }
73 // Process one digit.
74 Bool_t ProcessDigit(AliFMDDigit* d)
75 {
76 fHist->Fill(d->Counts());
77 return kTRUE;
78 }
79 // After processing all events, display spectrum
80 Bool_t Finish()
81 {
82 fHist->Draw();
83 }
84 TH1F* fHist;
85 };
86
87 void AdcSpectrum()
88 {
89 DigitInput di;
90 di.Run();
91 }
92 @endcode
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.
99 @ingroup FMD_util
100 */
a1f80595 101class AliFMDInput : public TObject
102{
103public:
9f662337 104 /** The kinds of data that can be read in. */
a1f80595 105 enum ETrees {
8f6ee336 106 kHits = 1, // Hits
107 kKinematics, // Kinematics (from sim)
108 kDigits, // Digits
109 kSDigits, // Summable digits
110 kHeader, // Header information
111 kRecPoints, // Reconstructed points
bf000c32 112 kESD, // Load ESD's
d760ea03 113 kRaw, // Read raw data
f95a63c4 114 kGeometry, // Not really a tree
115 kTracks // Hits and tracs - for BG study
a1f80595 116 };
9f662337 117 /** CTOR */
a1f80595 118 AliFMDInput();
9f662337 119 /** CTOR
120 @param gAliceFile galice file */
a1f80595 121 AliFMDInput(const char* gAliceFile);
9f662337 122 /** DTOR */
a1f80595 123 virtual ~AliFMDInput() {}
124
9f662337 125 /** Add a data type to load
126 @param tree Data to load */
a1f80595 127 virtual void AddLoad(ETrees tree) { SETBIT(fTreeMask, tree); }
9f662337 128 /** Remove a data type to load
129 @param tree Data to @e not load */
a1f80595 130 virtual void RemoveLoad(ETrees tree) { CLRBIT(fTreeMask, tree); }
9f662337 131 /** @return # of available events */
a1f80595 132 virtual Int_t NEvents() const;
133
9f662337 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 */
a1f80595 137 virtual Bool_t Init();
9f662337 138 /** Callled at the beginning of each event. If a user class
139 overloads this member function, then this @e must be explicitly
140 called.
141 @param event Event number
142 @return @c false on error */
a1f80595 143 virtual Bool_t Begin(Int_t event);
9f662337 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 */
bf000c32 148 virtual Bool_t Event();
9f662337 149 /** Called at the end of each event.
150 @return @c false on error */
a1f80595 151 virtual Bool_t End();
9f662337 152 /** Called at the end of the run.
153 @return @c false on error */
a1f80595 154 virtual Bool_t Finish() { return kTRUE; }
9f662337 155 /** Run a full job.
156 @return @c false on error */
a1f80595 157 virtual Bool_t Run();
bf000c32 158
9f662337 159 /** Loop over all hits, and call ProcessHit with that hit, and
160 optionally the corresponding kinematics track.
161 @return @c false on error */
bf000c32 162 virtual Bool_t ProcessHits();
f95a63c4 163 /** Loop over all tracks, and call ProcessTrack with each hit for
164 that track
165 @return @c false on error */
166 virtual Bool_t ProcessTracks();
9f662337 167 /** Loop over all digits, and call ProcessDigit for each digit.
168 @return @c false on error */
bf000c32 169 virtual Bool_t ProcessDigits();
9f662337 170 /** Loop over all summable digits, and call ProcessSDigit for each
171 digit.
172 @return @c false on error */
bf000c32 173 virtual Bool_t ProcessSDigits();
9f662337 174 /** Loop over all digits read from raw data files, and call
175 ProcessRawDigit for each digit.
176 @return @c false on error */
d760ea03 177 virtual Bool_t ProcessRawDigits();
9f662337 178 /** Loop over all reconstructed points, and call ProcessRecPoint for
179 each reconstructed point.
180 @return @c false on error */
bf000c32 181 virtual Bool_t ProcessRecPoints();
a9579262 182 /** Loop over all ESD data, and call ProcessESD for each entry.
183 @return @c false on error */
184 virtual Bool_t ProcessESDs();
bf000c32 185
9f662337 186 /** Process one hit, and optionally it's corresponding kinematics
187 track. Users should over this to process each hit.
a9579262 188 @param h Hit
189 @param p Associated track
9f662337 190 @return @c false on error */
a9579262 191 virtual Bool_t ProcessHit(AliFMDHit* h, TParticle* p);
f95a63c4 192 /** Process one hit per track. Users should over this to process
193 each hit.
194 @param i Track number
195 @param p Track
196 @param h Associated Hit
197 @return @c false on error */
198 virtual Bool_t ProcessTrack(Int_t i, TParticle* p, AliFMDHit* h);
a9579262 199 /** Process one digit. Users should over this to process each
200 digit.
201 @param digit Digit
9f662337 202 @return @c false on error */
a9579262 203 virtual Bool_t ProcessDigit(AliFMDDigit* digit);
9f662337 204 /** Process one summable digit. Users should over this to process
205 each summable digit.
a9579262 206 @param sdigit Summable digit
9f662337 207 @return @c false on error */
a9579262 208 virtual Bool_t ProcessSDigit(AliFMDSDigit* sdigit);
9f662337 209 /** Process one digit from raw data files. Users should over this
210 to process each raw digit.
a9579262 211 @param digit Raw digit
9f662337 212 @return @c false on error */
a9579262 213 virtual Bool_t ProcessRawDigit(AliFMDDigit* digit);
9f662337 214 /** Process one reconstructed point. Users should over this to
215 process each reconstructed point.
a9579262 216 @param point Reconstructed point
9f662337 217 @return @c false on error */
a9579262 218 virtual Bool_t ProcessRecPoint(AliFMDRecPoint* point);
9f662337 219 /** Process ESD data for the FMD. Users should overload this to
220 deal with ESD data.
a9579262 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
9f662337 227 @return @c false on error */
a9579262 228 virtual Bool_t ProcessESD(UShort_t, Char_t, UShort_t, UShort_t,
229 Float_t, Float_t);
69893a66 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);
039842fe 236
237 /** Set the raw data input
238 @param file File name - if empty, assume simulated raw. */
239 void SetRawFile(const char* file) { if (file) fRawFile = file; }
240
a1f80595 241protected:
02a27b50 242 /** Copy ctor
243 @param o Object to copy from */
b5ee4425 244 AliFMDInput(const AliFMDInput& o)
245 : TObject(o),
246 fGAliceFile(""),
247 fLoader(0),
248 fRun(0),
249 fStack(0),
250 fFMDLoader(0),
251 fReader(0),
252 fFMD(0),
b5ee4425 253 fESD(0),
df137876 254 fESDEvent(0),
b5ee4425 255 fTreeE(0),
256 fTreeH(0),
257 fTreeD(0),
258 fTreeS(0),
259 fTreeR(0),
260 fTreeA(0),
261 fChainE(0),
262 fArrayE(0),
263 fArrayH(0),
264 fArrayD(0),
265 fArrayS(0),
266 fArrayR(0),
267 fArrayA(0),
17e542eb 268 fHeader(0),
b5ee4425 269 fGeoManager(0),
270 fTreeMask(0),
6ce810fc 271 fRawFile(""),
17e542eb 272 fIsInit(kFALSE),
273 fEventCount(0)
b5ee4425 274 {}
02a27b50 275 /** Assignement operator
276 @return REference to this */
277 AliFMDInput& operator=(const AliFMDInput&) { return *this; }
278
a1f80595 279 TString fGAliceFile; // File name of gAlice file
280 AliRunLoader* fLoader; // Loader of FMD data
281 AliRun* fRun; // Run information
282 AliStack* fStack; // Stack of particles
283 AliLoader* fFMDLoader; // Loader of FMD data
d760ea03 284 AliRawReader* fReader; // Raw data reader
a1f80595 285 AliFMD* fFMD; // FMD object
bf000c32 286 AliESDFMD* fESD; // FMD ESD data
df137876 287 AliESDEvent* fESDEvent; // ESD Event object.
a1f80595 288 TTree* fTreeE; // Header tree
289 TTree* fTreeH; // Hits tree
290 TTree* fTreeD; // Digit tree
291 TTree* fTreeS; // SDigit tree
292 TTree* fTreeR; // RecPoint tree
d760ea03 293 TTree* fTreeA; // Raw data tree
bf000c32 294 TChain* fChainE; // Chain of ESD's
a1f80595 295 TClonesArray* fArrayE; // Event info array
296 TClonesArray* fArrayH; // Hit info array
297 TClonesArray* fArrayD; // Digit info array
298 TClonesArray* fArrayS; // SDigit info array
d760ea03 299 TClonesArray* fArrayR; // Rec points info array
300 TClonesArray* fArrayA; // Raw data (digits) info array
9b98d361 301 AliHeader* fHeader; // Header
8f6ee336 302 TGeoManager* fGeoManager; // Geometry manager
a1f80595 303 Int_t fTreeMask; // Which tree's to load
039842fe 304 TString fRawFile; // Raw input file
02a27b50 305 Bool_t fIsInit; // Have we been initialized
69893a66 306 Int_t fEventCount; // Event counter
a1f80595 307 ClassDef(AliFMDInput,0) //Hits for detector FMD
308};
309
a9579262 310inline Bool_t AliFMDInput::ProcessHit(AliFMDHit*,TParticle*) { return kTRUE; }
f95a63c4 311inline Bool_t AliFMDInput::ProcessTrack(Int_t,TParticle*,
312 AliFMDHit*) { return kTRUE; }
a9579262 313inline Bool_t AliFMDInput::ProcessDigit(AliFMDDigit*) { return kTRUE; }
314inline Bool_t AliFMDInput::ProcessSDigit(AliFMDSDigit*) { return kTRUE; }
315inline Bool_t AliFMDInput::ProcessRawDigit(AliFMDDigit*) { return kTRUE; }
316inline Bool_t AliFMDInput::ProcessRecPoint(AliFMDRecPoint*) { return kTRUE; }
317inline Bool_t AliFMDInput::ProcessESD(UShort_t,Char_t,UShort_t,UShort_t,
318 Float_t,Float_t) { return kTRUE; }
319
a1f80595 320
a1f80595 321#endif
322//____________________________________________________________________
323//
324// Local Variables:
325// mode: C++
326// End:
327//
328// EOF
329//