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