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 */
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.
14 // The objects are retrived from folders.
15 //*-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (RRC KI & SUBATECH)
19 // --- ROOT system ---
20 #include "TClonesArray.h"
27 // --- Standard library ---
29 // --- AliRoot header files ---
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"
43 class AliEMCALRecParticle ;
44 class AliEMCALGeometry ;
45 class AliEMCALDigitizer ;
46 class AliEMCALSDigitizer ;
47 class AliEMCALCalibrationDB ;
52 class AliEMCALLoader : public AliLoader {
57 AliEMCALLoader(const AliEMCALLoader & obj):AliLoader(obj){}
58 AliEMCALLoader(const Char_t *detname,const Char_t *eventfoldername);
60 virtual ~AliEMCALLoader() ;
62 // assignement operator requested by coding convention, but not needed
63 const AliEMCALLoader & operator = (const AliEMCALLoader & ) {return *this;}
65 Int_t GetEvent();//extends the method on EMCAL RecPart posting
66 Int_t SetEvent();//extends the method on EMCAL RecPart posting
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="");
76 void UnloadRecParticles();
79 Int_t PostHits(); //Posts the
82 Int_t PostRecPoints();
84 Int_t PostRecParticles();
86 void CleanFolders();//cleans all the stuff loaded by this detector + calls AliLoader::Clean
91 void CleanRecPoints();
93 void CleanRecParticles();
95 //up to now it is only here -> no definition about global/incremental tracking/PID
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 // Bool_t PostQA (void) const ; //it was empty anyway
102 /*******************************************************************/
103 /*******************************************************************/
104 /*******************************************************************/
106 TObject** HitsRef(){return GetDetectorDataRef(Hits());}
107 TObject** SDigitsRef(){return GetDetectorDataRef(SDigits());}
108 TObject** DigitsRef(){return GetDetectorDataRef(Digits());}
109 TObject** PRERecPointsRef(){return GetDetectorDataRef(PRERecPoints());}
110 TObject** ECARecPointsRef(){return GetDetectorDataRef(ECARecPoints());}
111 TObject** HCARecPointsRef(){return GetDetectorDataRef(HCARecPoints());}
112 TObject** TracksRef(){return GetDetectorDataRef(TrackSegments());}
113 TObject** RecParticlesRef(){return GetDetectorDataRef(RecParticles());}
114 TObject** AlarmsRef(){return GetDetectorDataRef(Alarms());}
115 void Track(Int_t itrack) ;
117 static AliEMCALGeometry* GetEMCALGeometry();
118 static AliEMCALLoader* GetEMCALLoader(const char* eventfoldername);
120 //Method to be used when digitizing under AliRunDigitizer, who opens all files etc.
121 Int_t EventNumber() { return (Int_t) GetRunLoader()->GetEventNumber();}
122 Int_t MaxEvent() { return (Int_t) GetRunLoader()->TreeE()->GetEntries();}
124 const AliEMCAL * EMCAL();
125 const AliEMCALGeometry *EMCALGeometry() ;
127 // TFolder * Alarms() const { return (TFolder*)(ReturnO("Alarms", 0)); }
128 TObjArray * Alarms();
130 /*********************************************/
131 /************ TClonesArrays ***********/
132 /*********************************************/
134 TClonesArray* Hits(void);
135 const AliEMCALHit* Hit(Int_t index);
136 void MakeHitsArray();
137 /**** S D i g i t s ****/
138 TClonesArray* SDigits();
139 const AliEMCALDigit* SDigit(Int_t index);
140 void MakeSDigitsArray();
141 /**** D i g i t s ****/
142 TClonesArray* Digits();
143 const AliEMCALDigit * Digit(Int_t index);
144 void MakeDigitsArray();
145 /**** R e c P o i n t s ****/
146 TObjArray * PRERecPoints();
147 TObjArray * ECARecPoints();
148 TObjArray * HCARecPoints();
149 const AliEMCALRecPoint * PRERecPoint(Int_t index) ;
150 const AliEMCALTowerRecPoint * ECARecPoint(Int_t index) ;
151 const AliEMCALTowerRecPoint * HCARecPoint(Int_t index) ;
152 void MakeRecPointsArray();
153 /**** T r a c k S e g m e n t s ****/
154 TClonesArray * TrackSegments();
155 const AliEMCALTrackSegment * TrackSegment(Int_t index);
156 void MakeTrackSegmentsArray();
157 /**** R e c P a r t ic l e s ****/
158 TClonesArray * RecParticles() ;
159 const AliEMCALRecParticle * RecParticle(Int_t index);
160 void MakeRecParticlesArray();
162 /*********************************************/
163 /************ T A S K S **************/
164 /*********************************************/
166 // AliEMCALSDigitizer* EMCALSDigitizer(TString name = AliConfig::fgkDefaultEventFolderName);
167 //AliEMCALDigitizer* EMCALDigitizer() { return dynamic_cast<AliEMCALDigitizer*>(Digitizer()) ;}
169 AliEMCALClusterizer* Clusterizer () {return dynamic_cast<AliEMCALClusterizer*>(Reconstructioner()) ;}
170 Int_t PostClusterizer(TTask* clust){return PostReconstructioner(clust);}
171 Int_t LoadClusterizer(Option_t * opt="") {return LoadReconstructioner(opt);}
172 Int_t WriteClusterizer(Option_t * opt="") {return WriteReconstructioner(opt);}
174 AliEMCALPID * PID (){return dynamic_cast<AliEMCALPID*>(PIDTask()) ;}
175 Int_t PostPID(TTask* pid){return PostPIDTask(pid);}
176 Int_t LoadPID(Option_t * opt="") {return LoadPIDTask(opt);}
177 Int_t WritePID(Option_t * opt="") {return WritePIDTask(opt);}
180 AliEMCALTrackSegmentMaker * TrackSegmentMaker () { return dynamic_cast<AliEMCALTrackSegmentMaker *>(Tracker()) ;}
181 Int_t PostTrackSegmentMaker(TTask* segmaker){return PostTracker(segmaker);}
182 Int_t LoadTrackSegmentMaker(Option_t * opt="") {return LoadTracker(opt);}
183 Int_t WriteTrackSegmentMaker(Option_t * opt="") {return WriteTracker(opt);}
186 void SetDebug(Int_t level) {fDebug = level;} // Set debug level
187 void SetBranchTitle(const TString& btitle);
189 AliEMCALCalibrationDB * CalibrationDB(){return fcdb; }
190 //void ReadCalibrationDB(const char * name, const char * filename);
193 static TString HitsName() { return fgkHitsName ; } //Name for TClonesArray with hits from one event
194 static TString SDigitsName() { return fgkSDigitsName ;} //Name for TClonesArray
195 static TString DigitsName() { return fgkDigitsName ;} //Name for TClonesArray
196 static TString PreRecPointsName() { return fgkPRERecPointsName ;} //Name for TClonesArray
197 static TString ECARecPointsName() { return fgkECARecPointsName ;} //Name for TClonesArray
198 static TString HCARecPointsName() { return fgkHCARecPointsName ;} //Name for TClonesArray
199 static TString TracksName() { return fgkTracksName ;} //Name for TClonesArray
200 static TString RecParticlesName() { return fgkRecParticlesName ;} //Name for TClonesArray
201 static TString PRERecPointsBranchName() { return fgkPRERecPointsBranchName ;} //Name for branch
202 static TString ECARecPointsBranchName() { return fgkECARecPointsBranchName ;} //Name for branch
203 static TString HCARecPointsBranchName() { return fgkHCARecPointsBranchName ;} //Name for branch
204 static TString TrackSegmentsBranchName() { return fgkTrackSegmentsBranchName ;} //Name for branch
205 static TString RecParticlesBranchName() { return fgkRecParticlesBranchName ;} //Name for branch
208 TString fBranchTitle; //Title of the branch
209 Bool_t fRecParticlesLoaded; //Flag signing if Reconstructed Particles are loaded
210 Bool_t fTracksLoaded; //Flag signing if Tracks are loaded
211 TString fRecParticlesFileOption; //Loading Option for Reconstructed Particles
212 AliEMCALCalibrationDB * fcdb ; //!
219 Int_t ReadRecPoints();
221 Int_t ReadRecParticles();
225 static const TString fgkHitsName;//Name for TClonesArray with hits from one event
226 static const TString fgkSDigitsName;//Name for TClonesArray
227 static const TString fgkDigitsName;//Name for TClonesArray
228 static const TString fgkPRERecPointsName;//Name for TClonesArray
229 static const TString fgkECARecPointsName;//Name for TClonesArray
230 static const TString fgkHCARecPointsName;//Name for TClonesArray
231 static const TString fgkTracksName;//Name for TClonesArray
232 static const TString fgkRecParticlesName;//Name for TClonesArray
234 static const TString fgkPRERecPointsBranchName;//Name for branch
235 static const TString fgkECARecPointsBranchName;//Name for branch
236 static const TString fgkHCARecPointsBranchName;//Name for branch
237 static const TString fgkTrackSegmentsBranchName;//Name for branch
238 static const TString fgkRecParticlesBranchName;//Name for branch
239 Int_t fDebug ; // Debug level
243 ClassDef(AliEMCALLoader,3) // Algorithm class that provides methods to retrieve objects from a list knowing the index
247 /******************************************************************************/
248 /**************** I N L I N E S ****************************************/
249 /******************************************************************************/
251 inline TClonesArray* AliEMCALLoader::Hits()
253 return (TClonesArray*)GetDetectorData(fgkHitsName);
255 /******************************************************************************/
257 inline const AliEMCALHit* AliEMCALLoader::Hit(Int_t index)
259 const TClonesArray* tcarr = Hits();
261 return (const AliEMCALHit*) tcarr->At(index);
264 /******************************************************************************/
266 inline TClonesArray* AliEMCALLoader::SDigits()
268 return dynamic_cast<TClonesArray*>(GetDetectorData(fgkSDigitsName));
270 /******************************************************************************/
272 inline const AliEMCALDigit* AliEMCALLoader::SDigit(Int_t index)
274 const TClonesArray* tcarr = SDigits();
276 return (const AliEMCALDigit*) tcarr->At(index);
279 /******************************************************************************/
281 inline TClonesArray* AliEMCALLoader::Digits()
283 return dynamic_cast<TClonesArray*>(GetDetectorData(fgkDigitsName));
285 /******************************************************************************/
287 inline const AliEMCALDigit* AliEMCALLoader::Digit(Int_t index)
289 const TClonesArray* tcarr = Digits();
291 return (const AliEMCALDigit*) tcarr->At(index);
295 /******************************************************************************/
297 inline TObjArray * AliEMCALLoader::PRERecPoints()
299 return dynamic_cast<TObjArray*>(GetDetectorData(fgkPRERecPointsName));
301 /******************************************************************************/
303 inline const AliEMCALRecPoint * AliEMCALLoader::PRERecPoint(Int_t index)
305 TObjArray* tcarr = PRERecPoints();
307 return dynamic_cast<const AliEMCALRecPoint*>(tcarr->At(index));
311 /******************************************************************************/
313 inline TObjArray * AliEMCALLoader::ECARecPoints()
315 return dynamic_cast<TObjArray*>(GetDetectorData(fgkECARecPointsName));
318 /******************************************************************************/
320 inline const AliEMCALTowerRecPoint * AliEMCALLoader::ECARecPoint(Int_t index)
322 TObjArray* tcarr = ECARecPoints();
324 return dynamic_cast<const AliEMCALTowerRecPoint*>(tcarr->At(index));
328 /******************************************************************************/
330 inline TObjArray * AliEMCALLoader::HCARecPoints()
332 return dynamic_cast<TObjArray*>(GetDetectorData(fgkHCARecPointsName));
335 /******************************************************************************/
337 inline const AliEMCALTowerRecPoint * AliEMCALLoader::HCARecPoint(Int_t index)
339 TObjArray* tcarr = HCARecPoints();
341 return dynamic_cast<const AliEMCALTowerRecPoint*>(tcarr->At(index));
345 /******************************************************************************/
347 inline TClonesArray * AliEMCALLoader::TrackSegments()
349 return dynamic_cast<TClonesArray*>(GetDetectorData(fgkTracksName));
351 /******************************************************************************/
353 inline const AliEMCALTrackSegment * AliEMCALLoader::TrackSegment(Int_t index)
355 const TClonesArray* tcarr = TrackSegments();
357 return (const AliEMCALTrackSegment*) tcarr->At(index);
360 /******************************************************************************/
362 inline TClonesArray * AliEMCALLoader::RecParticles()
364 return dynamic_cast<TClonesArray*>(GetDetectorData(fgkRecParticlesName));
366 /******************************************************************************/
368 inline const AliEMCALRecParticle* AliEMCALLoader::RecParticle(Int_t index)
370 TClonesArray* tcarr = RecParticles();
372 return (const AliEMCALRecParticle*) tcarr->At(index);
375 /******************************************************************************/
376 inline TObjArray * AliEMCALLoader::Alarms()
377 { return (TObjArray*)(GetQAFolder()->FindObject(fDetectorName));}
379 #endif // AliEMCALLOADER_H