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