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.
52 //___________________________________________________________________
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
60 struct DigitInput : public AliFMDInput
67 fHist = new TH1F("adc", "ADC spectra", 1024, -.5, 1023.5);
70 Bool_t ProcessDigit(AliFMDDigit* d)
72 fHist->Fill(d->Counts());
75 // After processing all events, display spectrum
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.
97 class AliFMDInput : public TObject
100 /** The kinds of data that can be read in. */
103 kKinematics, // Kinematics (from sim)
105 kSDigits, // Summable digits
106 kHeader, // Header information
107 kRecPoints, // Reconstructed points
109 kRaw, // Read raw data
110 kGeometry // Not really a tree
115 @param gAliceFile galice file */
116 AliFMDInput(const char* gAliceFile);
118 virtual ~AliFMDInput() {}
120 /** Add a data type to load
121 @param tree Data to load */
122 virtual void AddLoad(ETrees tree) { SETBIT(fTreeMask, tree); }
123 /** Remove a data type to load
124 @param tree Data to @e not load */
125 virtual void RemoveLoad(ETrees tree) { CLRBIT(fTreeMask, tree); }
126 /** @return # of available events */
127 virtual Int_t NEvents() const;
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 */
132 virtual Bool_t Init();
133 /** Callled at the beginning of each event. If a user class
134 overloads this member function, then this @e must be explicitly
136 @param event Event number
137 @return @c false on error */
138 virtual Bool_t Begin(Int_t event);
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 */
143 virtual Bool_t Event();
144 /** Called at the end of each event.
145 @return @c false on error */
146 virtual Bool_t End();
147 /** Called at the end of the run.
148 @return @c false on error */
149 virtual Bool_t Finish() { return kTRUE; }
151 @return @c false on error */
152 virtual Bool_t Run();
154 /** Loop over all hits, and call ProcessHit with that hit, and
155 optionally the corresponding kinematics track.
156 @return @c false on error */
157 virtual Bool_t ProcessHits();
158 /** Loop over all digits, and call ProcessDigit for each digit.
159 @return @c false on error */
160 virtual Bool_t ProcessDigits();
161 /** Loop over all summable digits, and call ProcessSDigit for each
163 @return @c false on error */
164 virtual Bool_t ProcessSDigits();
165 /** Loop over all digits read from raw data files, and call
166 ProcessRawDigit for each digit.
167 @return @c false on error */
168 virtual Bool_t ProcessRawDigits();
169 /** Loop over all reconstructed points, and call ProcessRecPoint for
170 each reconstructed point.
171 @return @c false on error */
172 virtual Bool_t ProcessRecPoints();
173 /** Loop over all ESD data, and call ProcessESD for each entry.
174 @return @c false on error */
175 virtual Bool_t ProcessESDs();
177 /** Process one hit, and optionally it's corresponding kinematics
178 track. Users should over this to process each hit.
180 @param p Associated track
181 @return @c false on error */
182 virtual Bool_t ProcessHit(AliFMDHit* h, TParticle* p);
183 /** Process one digit. Users should over this to process each
186 @return @c false on error */
187 virtual Bool_t ProcessDigit(AliFMDDigit* digit);
188 /** Process one summable digit. Users should over this to process
190 @param sdigit Summable digit
191 @return @c false on error */
192 virtual Bool_t ProcessSDigit(AliFMDSDigit* sdigit);
193 /** Process one digit from raw data files. Users should over this
194 to process each raw digit.
195 @param digit Raw digit
196 @return @c false on error */
197 virtual Bool_t ProcessRawDigit(AliFMDDigit* digit);
198 /** Process one reconstructed point. Users should over this to
199 process each reconstructed point.
200 @param point Reconstructed point
201 @return @c false on error */
202 virtual Bool_t ProcessRecPoint(AliFMDRecPoint* point);
203 /** Process ESD data for the FMD. Users should overload this to
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
211 @return @c false on error */
212 virtual Bool_t ProcessESD(UShort_t, Char_t, UShort_t, UShort_t,
217 @param o Object to copy from */
218 AliFMDInput(const AliFMDInput& o)
246 /** Assignement operator
247 @return REference to this */
248 AliFMDInput& operator=(const AliFMDInput&) { return *this; }
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
255 AliRawReader* fReader; // Raw data reader
256 AliFMD* fFMD; // FMD object
257 AliESD* fMainESD; // ESD Object
258 AliESDFMD* fESD; // FMD ESD data
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
264 TTree* fTreeA; // Raw data tree
265 TChain* fChainE; // Chain of ESD's
266 TClonesArray* fArrayE; // Event info array
267 TClonesArray* fArrayH; // Hit info array
268 TClonesArray* fArrayD; // Digit info array
269 TClonesArray* fArrayS; // SDigit info array
270 TClonesArray* fArrayR; // Rec points info array
271 TClonesArray* fArrayA; // Raw data (digits) info array
272 TGeoManager* fGeoManager; // Geometry manager
273 Int_t fTreeMask; // Which tree's to load
274 Bool_t fIsInit; // Have we been initialized
275 ClassDef(AliFMDInput,0) //Hits for detector FMD
278 inline Bool_t AliFMDInput::ProcessHit(AliFMDHit*,TParticle*) { return kTRUE; }
279 inline Bool_t AliFMDInput::ProcessDigit(AliFMDDigit*) { return kTRUE; }
280 inline Bool_t AliFMDInput::ProcessSDigit(AliFMDSDigit*) { return kTRUE; }
281 inline Bool_t AliFMDInput::ProcessRawDigit(AliFMDDigit*) { return kTRUE; }
282 inline Bool_t AliFMDInput::ProcessRecPoint(AliFMDRecPoint*) { return kTRUE; }
283 inline Bool_t AliFMDInput::ProcessESD(UShort_t,Char_t,UShort_t,UShort_t,
284 Float_t,Float_t) { return kTRUE; }
288 //____________________________________________________________________