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