]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FMD/AliFMDInput.h
Revert unwanted changes concerning PID in HLT
[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*/
42f1b2f5 29#include <TNamed.h>
a1f80595 30#ifndef ROOT_TString
31# include <TString.h>
32#endif
7003fbd0 33#ifndef ROOT_TArrayF
34# include <TArrayF.h>
35#endif
08d168d9 36class AliTrackReference;
8f6ee336 37class AliRunLoader;
38class AliLoader;
39class AliStack;
40class AliRun;
d760ea03 41class AliRawReader;
5cf05dbb 42class AliFMDRawReader;
8f6ee336 43class AliFMD;
44class AliFMDHit;
bf000c32 45class AliFMDDigit;
46class AliFMDSDigit;
47class AliFMDRecPoint;
df137876 48class AliESDEvent;
bf000c32 49class AliESDFMD;
9b98d361 50class AliHeader;
8f6ee336 51class TString;
52class TClonesArray;
53class TTree;
54class TGeoManager;
55class TParticle;
bf000c32 56class TChain;
a1f80595 57
58//___________________________________________________________________
9f662337 59/** @class AliFMDInput
60 @brief Base class for reading in various FMD data.
61 The class loops over all found events. For each event the
62 specified data is read in. The class then loops over all
63 elements of the read data, and process these with user defined
64 code.
65 @code
66 struct DigitInput : public AliFMDInput
67 {
68 DigitInput()
69 {
70 // Load digits
71 AddLoad(kDigits);
72 // Make a histogram
73 fHist = new TH1F("adc", "ADC spectra", 1024, -.5, 1023.5);
74 }
75 // Process one digit.
76 Bool_t ProcessDigit(AliFMDDigit* d)
77 {
78 fHist->Fill(d->Counts());
79 return kTRUE;
80 }
81 // After processing all events, display spectrum
82 Bool_t Finish()
83 {
84 fHist->Draw();
85 }
86 TH1F* fHist;
87 };
88
89 void AdcSpectrum()
90 {
91 DigitInput di;
92 di.Run();
93 }
94 @endcode
95 This class allows for writing small scripts, that can be compiled
96 with AcLIC, to do all sorts of tests, quick prototyping, and so
97 on. It has proven to be quiet useful. One can load more than
98 one type of data in one derived class, to for example to make
99 comparisons between hits and reconstructed points. See also the
100 various scripts in @c FMD/scripts.
101 @ingroup FMD_util
102 */
42f1b2f5 103class AliFMDInput : public TNamed
a1f80595 104{
105public:
9f662337 106 /** The kinds of data that can be read in. */
a1f80595 107 enum ETrees {
8f6ee336 108 kHits = 1, // Hits
109 kKinematics, // Kinematics (from sim)
110 kDigits, // Digits
111 kSDigits, // Summable digits
112 kHeader, // Header information
113 kRecPoints, // Reconstructed points
bf000c32 114 kESD, // Load ESD's
d760ea03 115 kRaw, // Read raw data
f95a63c4 116 kGeometry, // Not really a tree
08d168d9 117 kTracks, // Hits and tracs - for BG study
e064ab4a 118 kTrackRefs, // Track references - also for BG study
119 kRawCalib // Read raws and calibrate them
a1f80595 120 };
9f662337 121 /** CTOR */
a1f80595 122 AliFMDInput();
9f662337 123 /** CTOR
124 @param gAliceFile galice file */
a1f80595 125 AliFMDInput(const char* gAliceFile);
9f662337 126 /** DTOR */
a1f80595 127 virtual ~AliFMDInput() {}
128
9f662337 129 /** Add a data type to load
130 @param tree Data to load */
a1f80595 131 virtual void AddLoad(ETrees tree) { SETBIT(fTreeMask, tree); }
9f662337 132 /** Remove a data type to load
133 @param tree Data to @e not load */
a1f80595 134 virtual void RemoveLoad(ETrees tree) { CLRBIT(fTreeMask, tree); }
9f662337 135 /** @return # of available events */
a1f80595 136 virtual Int_t NEvents() const;
137
9f662337 138 /** Initialize the class. If a user class overloads this member
139 function, then this @e must be explicitly called
140 @return @c false on error */
a1f80595 141 virtual Bool_t Init();
9f662337 142 /** Callled at the beginning of each event. If a user class
143 overloads this member function, then this @e must be explicitly
144 called.
145 @param event Event number
146 @return @c false on error */
a1f80595 147 virtual Bool_t Begin(Int_t event);
9f662337 148 /** Process one event. This loops over all the loaded data. Users
149 can overload this member function, but then it's @e strongly
150 recommended to explicitly call this classes version.
151 @return @c false on error */
bf000c32 152 virtual Bool_t Event();
9f662337 153 /** Called at the end of each event.
154 @return @c false on error */
a1f80595 155 virtual Bool_t End();
9f662337 156 /** Called at the end of the run.
157 @return @c false on error */
a1f80595 158 virtual Bool_t Finish() { return kTRUE; }
9f662337 159 /** Run a full job.
160 @return @c false on error */
a1f80595 161 virtual Bool_t Run();
bf000c32 162
9f662337 163 /** Loop over all hits, and call ProcessHit with that hit, and
164 optionally the corresponding kinematics track.
165 @return @c false on error */
bf000c32 166 virtual Bool_t ProcessHits();
08d168d9 167 /** Loop over all track refs, and call ProcessTrackRef with that hit, and
168 optionally the corresponding kinematics track.
169 @return @c false on error */
170 virtual Bool_t ProcessTrackRefs();
f95a63c4 171 /** Loop over all tracks, and call ProcessTrack with each hit for
172 that track
173 @return @c false on error */
174 virtual Bool_t ProcessTracks();
9f662337 175 /** Loop over all digits, and call ProcessDigit for each digit.
176 @return @c false on error */
bf000c32 177 virtual Bool_t ProcessDigits();
9f662337 178 /** Loop over all summable digits, and call ProcessSDigit for each
179 digit.
180 @return @c false on error */
bf000c32 181 virtual Bool_t ProcessSDigits();
9f662337 182 /** Loop over all digits read from raw data files, and call
183 ProcessRawDigit for each digit.
184 @return @c false on error */
d760ea03 185 virtual Bool_t ProcessRawDigits();
e064ab4a 186 /** Loop over all digits read from raw data files, and call
187 ProcessRawDigit for each digit.
188 @return @c false on error */
189 virtual Bool_t ProcessRawCalibDigits();
9f662337 190 /** Loop over all reconstructed points, and call ProcessRecPoint for
191 each reconstructed point.
192 @return @c false on error */
bf000c32 193 virtual Bool_t ProcessRecPoints();
a9579262 194 /** Loop over all ESD data, and call ProcessESD for each entry.
195 @return @c false on error */
196 virtual Bool_t ProcessESDs();
bf000c32 197
9f662337 198 /** Process one hit, and optionally it's corresponding kinematics
199 track. Users should over this to process each hit.
a9579262 200 @param h Hit
201 @param p Associated track
9f662337 202 @return @c false on error */
a9579262 203 virtual Bool_t ProcessHit(AliFMDHit* h, TParticle* p);
08d168d9 204 /** Process one track reference, and optionally it's corresponding kinematics
205 track. Users should overload this to process each track reference.
206 @param trackRef Track Reference
207 @param track Associated track
208 @return @c false on error */
209 virtual Bool_t ProcessTrackRef(AliTrackReference* trackRef, TParticle* track);
f95a63c4 210 /** Process one hit per track. Users should over this to process
211 each hit.
212 @param i Track number
213 @param p Track
214 @param h Associated Hit
215 @return @c false on error */
216 virtual Bool_t ProcessTrack(Int_t i, TParticle* p, AliFMDHit* h);
a9579262 217 /** Process one digit. Users should over this to process each
218 digit.
219 @param digit Digit
9f662337 220 @return @c false on error */
a9579262 221 virtual Bool_t ProcessDigit(AliFMDDigit* digit);
9f662337 222 /** Process one summable digit. Users should over this to process
223 each summable digit.
a9579262 224 @param sdigit Summable digit
9f662337 225 @return @c false on error */
a9579262 226 virtual Bool_t ProcessSDigit(AliFMDSDigit* sdigit);
9f662337 227 /** Process one digit from raw data files. Users should over this
228 to process each raw digit.
a9579262 229 @param digit Raw digit
9f662337 230 @return @c false on error */
a9579262 231 virtual Bool_t ProcessRawDigit(AliFMDDigit* digit);
e064ab4a 232 /** Process one digit from raw data files. Users should over this
233 to process each raw digit.
234 @param digit Raw digit
235 @return @c false on error */
236 virtual Bool_t ProcessRawCalibDigit(AliFMDDigit* digit);
9f662337 237 /** Process one reconstructed point. Users should over this to
238 process each reconstructed point.
a9579262 239 @param point Reconstructed point
9f662337 240 @return @c false on error */
a9579262 241 virtual Bool_t ProcessRecPoint(AliFMDRecPoint* point);
9f662337 242 /** Process ESD data for the FMD. Users should overload this to
243 deal with ESD data.
a9579262 244 @param d Detector number (1-3)
245 @param r Ring identifier ('I' or 'O')
246 @param s Sector number (0-19, or 0-39)
247 @param t Strip number (0-511, or 0-255)
248 @param eta Psuedo-rapidity
249 @param mult Psuedo-multiplicity
9f662337 250 @return @c false on error */
a9579262 251 virtual Bool_t ProcessESD(UShort_t, Char_t, UShort_t, UShort_t,
252 Float_t, Float_t);
69893a66 253 /** Service function to make a logarithmic axis.
254 @param n Number of bins
255 @param min Minimum of axis
256 @param max Maximum of axis.
257 @return An array with the bin boundaries. */
258 static TArrayF MakeLogScale(Int_t n, Double_t min, Double_t max);
039842fe 259
260 /** Set the raw data input
261 @param file File name - if empty, assume simulated raw. */
262 void SetRawFile(const char* file) { if (file) fRawFile = file; }
263
a1f80595 264protected:
02a27b50 265 /** Copy ctor
266 @param o Object to copy from */
b5ee4425 267 AliFMDInput(const AliFMDInput& o)
42f1b2f5 268 : TNamed(o),
b5ee4425 269 fGAliceFile(""),
270 fLoader(0),
271 fRun(0),
272 fStack(0),
273 fFMDLoader(0),
274 fReader(0),
5cf05dbb 275 fFMDReader(0),
b5ee4425 276 fFMD(0),
b5ee4425 277 fESD(0),
df137876 278 fESDEvent(0),
b5ee4425 279 fTreeE(0),
280 fTreeH(0),
08d168d9 281 fTreeTR(0),
b5ee4425 282 fTreeD(0),
283 fTreeS(0),
284 fTreeR(0),
285 fTreeA(0),
286 fChainE(0),
287 fArrayE(0),
288 fArrayH(0),
08d168d9 289 fArrayTR(0),
b5ee4425 290 fArrayD(0),
291 fArrayS(0),
292 fArrayR(0),
293 fArrayA(0),
17e542eb 294 fHeader(0),
b5ee4425 295 fGeoManager(0),
296 fTreeMask(0),
6ce810fc 297 fRawFile(""),
17e542eb 298 fIsInit(kFALSE),
299 fEventCount(0)
b5ee4425 300 {}
02a27b50 301 /** Assignement operator
302 @return REference to this */
303 AliFMDInput& operator=(const AliFMDInput&) { return *this; }
304
5cf05dbb 305 TString fGAliceFile; // File name of gAlice file
306 AliRunLoader* fLoader; // Loader of FMD data
307 AliRun* fRun; // Run information
308 AliStack* fStack; // Stack of particles
309 AliLoader* fFMDLoader; // Loader of FMD data
310 AliRawReader* fReader; // Raw data reader
311 AliFMDRawReader* fFMDReader; // FMD raw reader
312 AliFMD* fFMD; // FMD object
313 AliESDFMD* fESD; // FMD ESD data
314 AliESDEvent* fESDEvent; // ESD Event object.
315 TTree* fTreeE; // Header tree
316 TTree* fTreeH; // Hits tree
08d168d9 317 TTree* fTreeTR; // Track Reference tree
5cf05dbb 318 TTree* fTreeD; // Digit tree
319 TTree* fTreeS; // SDigit tree
320 TTree* fTreeR; // RecPoint tree
321 TTree* fTreeA; // Raw data tree
322 TChain* fChainE; // Chain of ESD's
323 TClonesArray* fArrayE; // Event info array
324 TClonesArray* fArrayH; // Hit info array
08d168d9 325 TClonesArray* fArrayTR; // Hit info array
5cf05dbb 326 TClonesArray* fArrayD; // Digit info array
327 TClonesArray* fArrayS; // SDigit info array
328 TClonesArray* fArrayR; // Rec points info array
329 TClonesArray* fArrayA; // Raw data (digits) info array
330 AliHeader* fHeader; // Header
331 TGeoManager* fGeoManager; // Geometry manager
332 Int_t fTreeMask; // Which tree's to load
333 TString fRawFile; // Raw input file
334 Bool_t fIsInit; // Have we been initialized
335 Int_t fEventCount; // Event counter
a1f80595 336 ClassDef(AliFMDInput,0) //Hits for detector FMD
337};
338
a9579262 339inline Bool_t AliFMDInput::ProcessHit(AliFMDHit*,TParticle*) { return kTRUE; }
faf80567 340inline Bool_t AliFMDInput::ProcessTrackRef(AliTrackReference*,
341 TParticle*) { return kTRUE; }
f95a63c4 342inline Bool_t AliFMDInput::ProcessTrack(Int_t,TParticle*,
343 AliFMDHit*) { return kTRUE; }
a9579262 344inline Bool_t AliFMDInput::ProcessDigit(AliFMDDigit*) { return kTRUE; }
345inline Bool_t AliFMDInput::ProcessSDigit(AliFMDSDigit*) { return kTRUE; }
346inline Bool_t AliFMDInput::ProcessRawDigit(AliFMDDigit*) { return kTRUE; }
e064ab4a 347inline Bool_t AliFMDInput::ProcessRawCalibDigit(AliFMDDigit*) { return kTRUE; }
a9579262 348inline Bool_t AliFMDInput::ProcessRecPoint(AliFMDRecPoint*) { return kTRUE; }
349inline Bool_t AliFMDInput::ProcessESD(UShort_t,Char_t,UShort_t,UShort_t,
350 Float_t,Float_t) { return kTRUE; }
351
a1f80595 352
a1f80595 353#endif
354//____________________________________________________________________
355//
356// Local Variables:
357// mode: C++
358// End:
359//
360// EOF
361//