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