]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALLoader.h
Adding MUON HLT code to the repository.
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALLoader.h
1 #ifndef ALIEMCALLOADER_H
2 #define ALIEMCALLOADER_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 /* $Id$ */
7
8 //_________________________________________________________________________
9 //  A singleton that returns various objects 
10 //  Should be used on the analysis stage to avoid confusing between different
11 //  branches of reconstruction tree: e.g. reading RecPoints and TS made from 
12 //  another set of RecPoints.
13 // 
14 //  The objects are retrived from folders.  
15 //*-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (RRC KI & SUBATECH)
16 //    
17
18
19 // --- ROOT system ---
20 #include "TClonesArray.h"
21 #include "TFolder.h"  
22 #include "TTree.h"
23 class TString ;
24 class TParticle ;
25 class TTask ;
26
27 // --- Standard library ---
28
29 // --- AliRoot header files ---
30
31 #include "AliRun.h"
32 #include "AliRunLoader.h"
33 #include "AliEMCALDigit.h"
34 #include "AliEMCALRecPoint.h"
35 #include "AliEMCALTrackSegment.h"
36 #include "AliEMCALClusterizer.h"
37 #include "AliEMCALTrackSegmentMaker.h" 
38 #include "AliEMCALPID.h"
39 class AliLoader ;
40 class AliEMCAL ; 
41 class AliEMCALHit ;
42 class AliEMCALRecParticle ; 
43 class AliEMCALGeometry ;
44 class AliEMCALDigitizer ; 
45 class AliEMCALSDigitizer ; 
46 class AliEMCALCalibrationDB ;
47
48
49 //
50
51 class AliEMCALLoader : public AliLoader {
52   
53  public:
54
55   AliEMCALLoader();
56   AliEMCALLoader(const AliEMCALLoader & obj):AliLoader(obj){}
57   AliEMCALLoader(const Char_t *detname,const Char_t *eventfoldername); 
58   
59   virtual ~AliEMCALLoader() ; 
60
61   // assignement operator requested by coding convention, but not needed
62   const AliEMCALLoader & operator = (const AliEMCALLoader & ) {return *this;}
63
64   Int_t   GetEvent();//extends the method on EMCAL RecPart posting
65   Int_t   SetEvent();//extends the method on EMCAL RecPart posting
66   
67   Bool_t  BranchExists(const TString& recName);
68   Int_t   LoadHits(Option_t* opt=""); //reads  from disk and sends them to folder; array as well as tree
69   Int_t   LoadSDigits(Option_t* opt="");
70   Int_t   LoadDigits(Option_t* opt=""); //reads Digits from disk and sends them to folder; array as well as tree
71   Int_t   LoadRecPoints(Option_t* opt=""); //reads RecPoints from disk and sends them to folder; array as well as tree
72   Int_t   LoadTracks(Option_t* opt="");  //reads Tracks from disk and sends them to folder; array as well as tree
73   Int_t   LoadRecParticles(Option_t* opt="");
74  
75   void    UnloadRecParticles();
76   void    UnloadTracks();
77   
78   Int_t   PostHits();  //Posts the 
79   Int_t   PostSDigits();
80   Int_t   PostDigits();
81   Int_t   PostRecPoints();
82   Int_t   PostTracks();
83   Int_t   PostRecParticles();
84   
85   void    CleanFolders();//cleans all the stuff loaded by this detector + calls AliLoader::Clean
86
87   void    CleanHits();
88   void    CleanSDigits();
89   void    CleanDigits();
90   void    CleanRecPoints();
91   void    CleanTracks();
92   void    CleanRecParticles();
93
94 //up to now it is only here -> no definition about global/incremental tracking/PID
95  
96 //   Int_t   WriteRecParticles(Option_t* opt="");//writes the reconstructed particles
97 //   Int_t   WritePID(Option_t* opt="");//writes the task for PID to file
98 //   Bool_t  PostPID  (AliEMCALPID * pid) const {return kTRUE;}
99
100 /*******************************************************************/
101 /*******************************************************************/
102 /*******************************************************************/
103
104   TObject** HitsRef(){return GetDetectorDataRef(Hits());}
105   TObject** SDigitsRef(){return GetDetectorDataRef(SDigits());}
106   TObject** DigitsRef(){return GetDetectorDataRef(Digits());}
107   TObject** ECARecPointsRef(){return GetDetectorDataRef(ECARecPoints());}
108   TObject** TracksRef(){return GetDetectorDataRef(TrackSegments());}
109   TObject** RecParticlesRef(){return GetDetectorDataRef(RecParticles());}
110
111   void   Track(Int_t itrack) ;
112
113   static AliEMCALGeometry* GetEMCALGeometry();
114   static AliEMCALLoader* GetEMCALLoader(const  char* eventfoldername);
115
116   //Method to be used when digitizing under AliRunDigitizer, who opens all files etc.
117   Int_t  EventNumber()       { return (Int_t) GetRunLoader()->GetEventNumber();}
118   Int_t  MaxEvent()          { return (Int_t) GetRunLoader()->TreeE()->GetEntries();}
119
120   const AliEMCAL *         EMCAL();
121   const AliEMCALGeometry  *EMCALGeometry() ; 
122
123   /*********************************************/
124   /************    TClonesArrays     ***********/
125   /*********************************************/
126   /****   H i t s  ****/
127   TClonesArray*  Hits(void);
128   const AliEMCALHit*    Hit(Int_t index);
129   void MakeHitsArray();
130   /****   S D i g i t s  ****/ 
131   TClonesArray*  SDigits();
132   const AliEMCALDigit*  SDigit(Int_t index);
133   void MakeSDigitsArray();
134   /****  D i g i t s  ****/
135   TClonesArray*   Digits();
136   const AliEMCALDigit *  Digit(Int_t index);
137   void MakeDigitsArray();
138   /****  R e c P o i n t s  ****/
139   TObjArray * ECARecPoints();
140   const AliEMCALRecPoint * ECARecPoint(Int_t index) ;
141   void MakeRecPointsArray();
142   /****   T r a c k S e g m e n t s ****/
143   TClonesArray * TrackSegments();
144   const AliEMCALTrackSegment * TrackSegment(Int_t index);
145   void MakeTrackSegmentsArray();
146   /****  R e c P a r t ic l e s   ****/
147   TClonesArray * RecParticles() ;
148   const AliEMCALRecParticle * RecParticle(Int_t index);
149   void MakeRecParticlesArray();
150
151   /*********************************************/
152   /************    T A S K S      **************/
153   /*********************************************/
154   // 
155   //  AliEMCALSDigitizer*  EMCALSDigitizer(TString name = AliConfig::GetDefaultEventFolderName());
156   //AliEMCALDigitizer*   EMCALDigitizer()  { return  dynamic_cast<AliEMCALDigitizer*>(Digitizer()) ;}
157
158   AliEMCALClusterizer* Clusterizer () const {return dynamic_cast<AliEMCALClusterizer*>(Reconstructioner()) ;}
159   Int_t PostClusterizer(TTask* clust) const {return PostReconstructioner(clust);}
160   Int_t LoadClusterizer(Option_t * opt="") const {return LoadReconstructioner(opt);}
161   Int_t WriteClusterizer(Option_t * opt="") const {return WriteReconstructioner(opt);}
162
163   AliEMCALPID * PID() const {return dynamic_cast<AliEMCALPID*>(PIDTask()) ;}
164   Int_t PostPID(TTask* pid) const {return PostPIDTask(pid);}
165   Int_t LoadPID(Option_t * opt="") const {return LoadPIDTask(opt);}
166   Int_t WritePID(Option_t * opt="") const {return WritePIDTask(opt);}
167
168
169   AliEMCALTrackSegmentMaker * TrackSegmentMaker () const  { return dynamic_cast<AliEMCALTrackSegmentMaker *>(Tracker()) ;}
170   Int_t PostTrackSegmentMaker(TTask* segmaker) const {return PostTracker(segmaker);}
171   Int_t LoadTrackSegmentMaker(Option_t * opt="") const {return LoadTracker(opt);}
172   Int_t WriteTrackSegmentMaker(Option_t * opt="") const  {return WriteTracker(opt);}
173
174   
175   void   SetDebug(Int_t level) {fDebug = level;} // Set debug level
176   void   SetBranchTitle(const TString& btitle);
177   
178   AliEMCALCalibrationDB * CalibrationDB(){return  fcdb; }
179   //void ReadCalibrationDB(const char * name, const char * filename);
180   
181
182   static TString HitsName() { return fgkHitsName ; }   //Name for TClonesArray with hits from one event
183   static TString SDigitsName() { return fgkSDigitsName ;} //Name for TClonesArray 
184   static TString DigitsName() { return fgkDigitsName ;} //Name for TClonesArray 
185   static TString ECARecPointsName() { return fgkECARecPointsName ;} //Name for TClonesArray y 
186   static TString TracksName() { return fgkTracksName ;} //Name for TClonesArray 
187   static TString RecParticlesName() { return fgkRecParticlesName ;} //Name for TClonesArray
188   static TString ECARecPointsBranchName() { return fgkECARecPointsBranchName ;} //Name for branch
189   static TString TrackSegmentsBranchName() { return fgkTrackSegmentsBranchName ;} //Name for branch
190   static TString RecParticlesBranchName() { return fgkRecParticlesBranchName ;} //Name for branch
191
192 protected:
193   TString fBranchTitle;            //Title of the branch
194   Bool_t  fRecParticlesLoaded;     //Flag signing if Reconstructed Particles are loaded
195   Bool_t  fTracksLoaded;           //Flag signing if Tracks are loaded
196   TString fRecParticlesFileOption; //Loading Option for Reconstructed Particles
197   AliEMCALCalibrationDB * fcdb ;       //!
198
199 private:
200
201   Int_t ReadHits();
202   Int_t ReadDigits();
203   Int_t ReadSDigits();
204   Int_t ReadRecPoints();
205   Int_t ReadTracks();
206   Int_t ReadRecParticles();
207
208   static const TString fgkHitsName;//Name for TClonesArray with hits from one event
209   static const TString fgkSDigitsName;//Name for TClonesArray 
210   static const TString fgkDigitsName;//Name for TClonesArray 
211   static const TString fgkECARecPointsName;//Name for TClonesArray 
212   static const TString fgkTracksName;//Name for TClonesArray 
213   static const TString fgkRecParticlesName;//Name for TClonesArray
214
215   static const TString fgkECARecPointsBranchName;//Name for branch
216   static const TString fgkTrackSegmentsBranchName;//Name for branch
217   static const TString fgkRecParticlesBranchName;//Name for branch
218   Int_t  fDebug ;             // Debug level
219  
220
221   
222   ClassDef(AliEMCALLoader,4)  // Algorithm class that provides methods to retrieve objects from a list knowing the index 
223
224 };
225
226 /******************************************************************************/
227 /****************    I N L I N E S     ****************************************/
228 /******************************************************************************/
229
230 inline TClonesArray* AliEMCALLoader::Hits()  
231 {
232  return (TClonesArray*)GetDetectorData(fgkHitsName);
233 }
234 /******************************************************************************/
235
236 inline const AliEMCALHit* AliEMCALLoader::Hit(Int_t index)  
237 {
238   const TClonesArray* tcarr = Hits();
239   if (tcarr)
240     return (const AliEMCALHit*) tcarr->At(index);
241   return 0x0; 
242 }
243 /******************************************************************************/
244
245 inline TClonesArray* AliEMCALLoader::SDigits()
246 {
247    return dynamic_cast<TClonesArray*>(GetDetectorData(fgkSDigitsName));
248 }
249 /******************************************************************************/
250
251 inline const AliEMCALDigit*  AliEMCALLoader::SDigit(Int_t index)
252 {
253   const TClonesArray* tcarr = SDigits();
254   if (tcarr)
255     return (const AliEMCALDigit*) tcarr->At(index);
256   return 0x0; 
257 }
258 /******************************************************************************/
259
260 inline TClonesArray* AliEMCALLoader::Digits()
261 {
262  return dynamic_cast<TClonesArray*>(GetDetectorData(fgkDigitsName));
263 }
264 /******************************************************************************/
265
266 inline const AliEMCALDigit*  AliEMCALLoader::Digit(Int_t index)
267 {
268   const TClonesArray* tcarr = Digits();
269   if (tcarr)
270     return (const AliEMCALDigit*) tcarr->At(index);
271   return 0x0; 
272 }
273
274 /******************************************************************************/
275
276 inline TObjArray * AliEMCALLoader::ECARecPoints()
277 {
278  return dynamic_cast<TObjArray*>(GetDetectorData(fgkECARecPointsName));
279 }
280
281 /******************************************************************************/
282
283 inline const AliEMCALRecPoint * AliEMCALLoader::ECARecPoint(Int_t index)
284 {
285   TObjArray* tcarr = ECARecPoints();
286   if (tcarr)
287     return dynamic_cast<const AliEMCALRecPoint*>(tcarr->At(index));
288   return 0x0; 
289 }
290
291 /******************************************************************************/
292
293 inline TClonesArray * AliEMCALLoader::TrackSegments()
294 {
295  return dynamic_cast<TClonesArray*>(GetDetectorData(fgkTracksName));
296 }
297 /******************************************************************************/
298
299 inline const AliEMCALTrackSegment * AliEMCALLoader::TrackSegment(Int_t index)
300 {
301   const TClonesArray* tcarr = TrackSegments();
302   if (tcarr)
303     return (const AliEMCALTrackSegment*) tcarr->At(index);
304   return 0x0; 
305 }
306 /******************************************************************************/
307
308 inline TClonesArray * AliEMCALLoader::RecParticles() 
309 {
310  return dynamic_cast<TClonesArray*>(GetDetectorData(fgkRecParticlesName)); 
311 }
312 /******************************************************************************/
313
314 inline const AliEMCALRecParticle* AliEMCALLoader::RecParticle(Int_t index)
315 {
316   TClonesArray* tcarr = RecParticles();
317   if (tcarr)
318     return (const AliEMCALRecParticle*) tcarr->At(index);
319   return 0x0;  
320 }
321
322 #endif // AliEMCALLOADER_H