]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALLoader.cxx
fblanco and dsilverm - add proper pedestal, and pedestal rms, calculation
[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
62 //____________________________________________________________________________ 
63 AliEMCALLoader::AliEMCALLoader()
64   : fDebug(0),
65     fHits(0),
66     fDigits(0),
67     fSDigits(0),
68     fRecPoints(0)
69 {
70   //Default constructor for EMCAL Loader Class
71
72   fDebug = 0;
73   fHits = new TClonesArray("AliEMCALHit");
74   fDigits = new TClonesArray("AliEMCALDigit");
75   fSDigits = new TClonesArray("AliEMCALDigit");
76   fRecPoints = new TObjArray();
77 }
78
79 //____________________________________________________________________________ 
80 AliEMCALLoader::AliEMCALLoader(const Char_t *detname,const Char_t *eventfoldername)
81   : AliLoader(detname,eventfoldername),
82     fDebug(0),
83     fHits(0),
84     fDigits(0),
85     fSDigits(0),
86     fRecPoints(0)
87 {
88   //Specific constructor for EMCAL Loader class
89
90   fDebug=0;
91   fHits = new TClonesArray("AliEMCALHit");
92   fDigits = new TClonesArray("AliEMCALDigit");
93   fSDigits = new TClonesArray("AliEMCALDigit");
94   fRecPoints = new TObjArray();
95 }
96
97 //____________________________________________________________________________
98 AliEMCALLoader::AliEMCALLoader(const Char_t *name, TFolder *topfolder)
99   : AliLoader(name,topfolder),
100     fDebug(0),
101     fHits(0),
102     fDigits(0),
103     fSDigits(0),
104     fRecPoints(0)
105 {
106   //Specific constructor for EMCAL Loader class
107
108   fDebug=0;
109   fHits = new TClonesArray("AliEMCALHit");
110   fDigits = new TClonesArray("AliEMCALDigit");
111   fSDigits = new TClonesArray("AliEMCALDigit");
112   fRecPoints = new TObjArray();
113 }
114
115 //____________________________________________________________________________ 
116 AliEMCALLoader::~AliEMCALLoader()
117 {
118   // Disconnect trees and remove arrays
119   if (TreeH())
120     TreeH()->SetBranchAddress(fDetectorName,0);
121   if (TreeD())
122     TreeD()->SetBranchAddress(fDetectorName,0);
123   if (TreeS())
124     TreeS()->SetBranchAddress(fDetectorName,0);
125   if (TreeR())
126     TreeR()->SetBranchAddress(fgkECARecPointsBranchName,0);
127   if (fHits) {
128     fHits->Delete();
129     delete fHits;
130   }
131   delete fDigits;
132   delete fSDigits;
133   delete fRecPoints;
134 }
135
136 //____________________________________________________________________________ 
137 AliEMCALCalibData* AliEMCALLoader::CalibData()
138
139   // Check if the instance of AliEMCALCalibData exists, if not, create it if 
140   // the OCDB is available, and finally return it.
141
142   if(!fgCalibData && (AliCDBManager::Instance()->IsDefaultStorageSet()))
143     {
144       AliCDBEntry *entry = (AliCDBEntry*) 
145         AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
146       if (entry) fgCalibData =  (AliEMCALCalibData*) entry->GetObject();
147     }
148   
149   if(!fgCalibData)
150     AliFatal("Calibration parameters not found in CDB!");
151   
152   return fgCalibData;
153   
154 }
155
156
157 //____________________________________________________________________________ 
158 Int_t AliEMCALLoader::CalibrateRaw(Double_t energy, Int_t module, 
159                                    Int_t column, Int_t row)
160 {
161   // Convert energy into digitized amplitude for a cell relId
162   // It is a user responsilibity to open CDB and set
163   // AliEMCALCalibData object by the following operators:
164   // 
165   // AliCDBLocal *loc = new AliCDBLocal("deCalibDB");
166   // AliEMCALCalibData* clb = (AliEMCALCalibData*)AliCDBStorage::Instance()
167   //    ->Get(path_to_calibdata,run_number);
168   // AliEMCALGetter* gime = AliEMCALGetter::Instance("galice.root");
169   // gime->SetCalibData(clb);
170
171   if (CalibData() == 0)
172     Warning("CalibrateRaw","Calibration DB is not initiated!");
173
174   Float_t gainFactor = 0.00305;//0.0015; // width of one ADC channel in GeV
175   Float_t pedestal   = 0.009;//0.005;  // pedestals
176
177   if(CalibData()) {
178     gainFactor = CalibData()->GetADCchannel (module,column,row);
179     pedestal   = CalibData()->GetADCpedestal(module,column,row);
180   }
181   
182   Int_t   amp = static_cast<Int_t>( (energy - pedestal) / gainFactor + 0.5 ) ; 
183   return amp;
184 }
185
186 //____________________________________________________________________________ 
187 Int_t AliEMCALLoader::GetEvent() 
188 {
189   //Method to load all of the data
190   //members of the EMCAL for a given
191   //event from the Trees
192
193   AliLoader::GetEvent();  // First call AliLoader to do all the groundwork
194   
195   // Now connect and fill TClonesArray
196
197   // Hits
198    TTree *treeH = TreeH();
199    
200    if (treeH) {
201      Int_t nEnt = treeH->GetEntries();  // TreeH has array of hits for every primary
202      fHits->Clear();
203      Int_t index = 0;
204      TClonesArray *tempArr = 0x0;
205      TBranch * branchH = treeH->GetBranch(fDetectorName);
206      branchH->SetAddress(&tempArr);
207      for (Int_t iEnt = 0; iEnt < nEnt; iEnt++) {
208        treeH->GetEntry(iEnt);
209        Int_t nHit = tempArr->GetEntriesFast();
210        for (Int_t iHit = 0; iHit < nHit; iHit++) {
211          new ((*fHits)[index]) AliEMCALHit(*((AliEMCALHit*)tempArr->At(iHit)));
212          index++;
213        }
214      }
215      branchH->ResetAddress();
216      if (tempArr) {
217        tempArr->Delete();
218        delete tempArr;
219      }
220    }
221   
222    // SDigits
223    TTree *treeS = TreeS();
224    if (treeS) {
225      TBranch * branchS = treeS->GetBranch(fDetectorName);
226      branchS->ResetAddress();
227      if (fSDigits) {
228        fSDigits->Clear();
229      }
230      branchS->SetAddress(&fSDigits);
231      treeS->GetEvent(0);
232    }
233    
234    // Digits
235    TTree *treeD = TreeD();
236    if (treeD) {
237      TBranch * branchD = treeD->GetBranch(fDetectorName);
238      branchD->ResetAddress();
239      if (fDigits) {
240        fDigits->Clear();
241      }
242      branchD->SetAddress(&fDigits);
243      treeD->GetEvent(0);
244    }
245
246    // RecPoints
247    TTree *treeR = TreeR();
248    if (treeR) {
249      TBranch * branchR = treeR->GetBranch(fgkECARecPointsBranchName);
250      branchR->ResetAddress();
251      if (fRecPoints) {
252        fRecPoints->Clear();
253      }
254      branchR->SetAddress(&fRecPoints);
255      treeR->GetEvent(0);
256    }
257    
258    return 0;
259 }