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