]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALLoader.cxx
Completely Updated (Mario Sitta)
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALLoader.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 //  A singleton. This class should be used in the analysis stage to get 
20 //  reconstructed objects: Digits, RecPoints, TrackSegments and RecParticles,
21 //  instead of directly reading them from galice.root file. 
22 //                  
23 //  MvL Feb 2006:
24 //  The AliEMCALLoader now holds the TClonesArray and TObjArray for reading
25 //  Hits, Dgits, SDigits and RecPoints. Filling is managed in the GetEvent()
26 //  method.
27 //
28 //  Creation/writing of files is managed by the relevant parts of the 
29 //  reconstruction software (AliEMCALDigitiser etx)
30 //
31 //*-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (RRC KI & SUBATECH)
32 //*--         Completely redesigned by Dmitri Peressounko March 2001  
33 //
34 //*-- YS June 2001 : renamed the original AliEMCALIndexToObject and make
35 //*--         systematic usage of TFolders without changing the interface
36 // 
37 //*-- Marco van Leeuwen, Jan 2006: complete revision to simplify reading
38 //*--         and fit better in general ALICE scheme
39 //
40 //////////////////////////////////////////////////////////////////////////////
41
42 // --- ROOT system ---
43 #include "TMath.h"
44 #include "TTree.h"
45
46 // --- Standard library ---
47
48 // --- AliRoot header files ---
49 #include "AliEMCALLoader.h"
50 #include "AliLog.h"
51 #include "AliCDBLocal.h"
52 #include "AliCDBStorage.h"
53 #include "AliCDBManager.h"
54 #include "AliCDBEntry.h"
55 #include "AliEMCALHit.h"
56
57 ClassImp(AliEMCALLoader)
58   
59 const TString AliEMCALLoader::fgkECARecPointsBranchName("EMCALECARP");//Name for branch with ECA Reconstructed Points
60 AliEMCALCalibData* AliEMCALLoader::fgCalibData = 0; //calibation data
61 AliEMCALRecParam * AliEMCALLoader::fgRecParam  = 0; //reconstruction parameters
62
63 //____________________________________________________________________________ 
64 AliEMCALLoader::AliEMCALLoader()
65   : fDebug(0),
66     fHits(0),
67     fDigits(0),
68     fSDigits(0),
69     fRecPoints(0)
70 {
71   //Default constructor for EMCAL Loader Class
72
73   fDebug = 0;
74   fHits = new TClonesArray("AliEMCALHit");
75   fDigits = new TClonesArray("AliEMCALDigit");
76   fSDigits = new TClonesArray("AliEMCALDigit");
77   fRecPoints = new TObjArray();
78 }
79
80 //____________________________________________________________________________ 
81 AliEMCALLoader::AliEMCALLoader(const Char_t *detname,const Char_t *eventfoldername)
82   : AliLoader(detname,eventfoldername),
83     fDebug(0),
84     fHits(0),
85     fDigits(0),
86     fSDigits(0),
87     fRecPoints(0)
88 {
89   //Specific constructor for EMCAL Loader class
90
91   fDebug=0;
92   fHits = new TClonesArray("AliEMCALHit");
93   fDigits = new TClonesArray("AliEMCALDigit");
94   fSDigits = new TClonesArray("AliEMCALDigit");
95   fRecPoints = new TObjArray();
96 }
97
98 //____________________________________________________________________________
99 AliEMCALLoader::AliEMCALLoader(const Char_t *name, TFolder *topfolder)
100   : AliLoader(name,topfolder),
101     fDebug(0),
102     fHits(0),
103     fDigits(0),
104     fSDigits(0),
105     fRecPoints(0)
106 {
107   //Specific constructor for EMCAL Loader class
108
109   fDebug=0;
110   fHits = new TClonesArray("AliEMCALHit");
111   fDigits = new TClonesArray("AliEMCALDigit");
112   fSDigits = new TClonesArray("AliEMCALDigit");
113   fRecPoints = new TObjArray();
114 }
115
116 //____________________________________________________________________________
117 AliEMCALLoader::AliEMCALLoader(const AliEMCALLoader & obj)
118   : AliLoader(obj),
119     fDebug(obj.fDebug),
120     fHits(obj.fHits),
121     fDigits(obj.fDigits),
122     fSDigits(obj.fSDigits),
123     fRecPoints(obj.fRecPoints)
124 {
125   //copy ctor
126 }
127
128 //____________________________________________________________________________ 
129 AliEMCALLoader::~AliEMCALLoader()
130 {
131   // Disconnect trees and remove arrays
132   if (TreeH())
133     TreeH()->SetBranchAddress(fDetectorName,0);
134   if (TreeD())
135     TreeD()->SetBranchAddress(fDetectorName,0);
136   if (TreeS())
137     TreeS()->SetBranchAddress(fDetectorName,0);
138   if (TreeR())
139     TreeR()->SetBranchAddress(fgkECARecPointsBranchName,0);
140   if (fHits) {
141     fHits->Delete();
142     delete fHits;
143   }
144   delete fDigits;
145   delete fSDigits;
146   delete fRecPoints;
147 }
148
149 //____________________________________________________________________________ 
150 AliEMCALCalibData* AliEMCALLoader::CalibData()
151
152   // Check if the instance of AliEMCALCalibData exists, if not, create it if 
153   // the OCDB is available, and finally return it.
154
155   if(!fgCalibData && (AliCDBManager::Instance()->IsDefaultStorageSet()))
156     {
157       AliCDBEntry *entry = (AliCDBEntry*) 
158         AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
159       if (entry) fgCalibData =  (AliEMCALCalibData*) entry->GetObject();
160     }
161   
162   if(!fgCalibData)
163     AliFatal("Calibration parameters not found in CDB!");
164   
165   return fgCalibData;
166   
167 }
168
169 //____________________________________________________________________________ 
170 AliEMCALRecParam* AliEMCALLoader::RecParam()
171
172   // Check if the instance of AliEMCALRecParam exists, if not, create it if 
173   // the OCDB is available, and finally return it.
174
175   if(!fgRecParam && (AliCDBManager::Instance()->IsDefaultStorageSet())) {
176     AliCDBEntry *entry = (AliCDBEntry*) 
177       AliCDBManager::Instance()->Get("EMCAL/Config/RecParam/");
178     if (entry) fgRecParam =  (AliEMCALRecParam*) entry->GetObject();
179   }
180   
181   if(!fgRecParam)
182     AliWarning("Recostruction parameters not found in CDB!");
183   
184   return fgRecParam;
185   
186 }
187
188 //____________________________________________________________________________ 
189 Int_t AliEMCALLoader::CalibrateRaw(Double_t energy, Int_t module, 
190                                    Int_t column, Int_t row)
191 {
192   // Convert energy into digitized amplitude for a cell relId
193   // It is a user responsilibity to open CDB and set
194   // AliEMCALCalibData object by the following operators:
195   // 
196   // AliCDBLocal *loc = new AliCDBLocal("deCalibDB");
197   // AliEMCALCalibData* clb = (AliEMCALCalibData*)AliCDBStorage::Instance()
198   //    ->Get(path_to_calibdata,run_number);
199   // AliEMCALGetter* gime = AliEMCALGetter::Instance("galice.root");
200   // gime->SetCalibData(clb);
201
202   if (CalibData() == 0)
203     Warning("CalibrateRaw","Calibration DB is not initiated!");
204
205   Float_t gainFactor = 0.00305;//0.0015; // width of one ADC channel in GeV
206   Float_t pedestal   = 0.009;//0.005;  // pedestals
207
208   if(CalibData()) {
209     gainFactor = CalibData()->GetADCchannel (module,column,row);
210     pedestal   = CalibData()->GetADCpedestal(module,column,row);
211   }
212   
213   Int_t   amp = static_cast<Int_t>( (energy - pedestal) / gainFactor + 0.5 ) ; 
214   return amp;
215 }
216
217 //____________________________________________________________________________ 
218 Int_t AliEMCALLoader::GetEvent() 
219 {
220   //Method to load all of the data
221   //members of the EMCAL for a given
222   //event from the Trees
223
224   AliLoader::GetEvent();  // First call AliLoader to do all the groundwork
225   
226   // Now connect and fill TClonesArray
227
228   // Hits
229    TTree *treeH = TreeH();
230    
231    if (treeH) {
232      Int_t nEnt = treeH->GetEntries();  // TreeH has array of hits for every primary
233      fHits->Clear();
234      Int_t index = 0;
235      TClonesArray *tempArr = 0x0;
236      TBranch * branchH = treeH->GetBranch(fDetectorName);
237      branchH->SetAddress(&tempArr);
238      for (Int_t iEnt = 0; iEnt < nEnt; iEnt++) {
239        treeH->GetEntry(iEnt);
240        Int_t nHit = tempArr->GetEntriesFast();
241        for (Int_t iHit = 0; iHit < nHit; iHit++) {
242          new ((*fHits)[index]) AliEMCALHit(*((AliEMCALHit*)tempArr->At(iHit)));
243          index++;
244        }
245      }
246      branchH->ResetAddress();
247      if (tempArr) {
248        tempArr->Delete();
249        delete tempArr;
250      }
251    }
252   
253    // SDigits
254    TTree *treeS = TreeS();
255    if (treeS) {
256      TBranch * branchS = treeS->GetBranch(fDetectorName);
257      branchS->ResetAddress();
258      if (fSDigits) {
259        fSDigits->Delete();
260        delete fSDigits;
261        fSDigits = 0x0;
262      }
263      branchS->SetAddress(&fSDigits);
264      treeS->GetEvent(0);
265    }
266    
267    // Digits
268    TTree *treeD = TreeD();
269    if (treeD) {
270      TBranch * branchD = treeD->GetBranch(fDetectorName);
271      branchD->ResetAddress();
272      if (fDigits) {
273        fDigits->Delete();
274        delete fDigits;
275        fDigits = 0x0;
276      }
277      branchD->SetAddress(&fDigits);
278      treeD->GetEvent(0);
279    }
280
281    // RecPoints
282    TTree *treeR = TreeR();
283    if (treeR) {
284      TBranch * branchR = treeR->GetBranch(fgkECARecPointsBranchName);
285      branchR->ResetAddress();
286      if (fRecPoints) {
287        fRecPoints->Delete();
288        delete fRecPoints;
289        fRecPoints = 0x0;
290      }
291      branchR->SetAddress(&fRecPoints);
292      treeR->GetEvent(0);
293    }
294    
295    return 0;
296 }