]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EMCAL/AliEMCALCalibTempCoeff.cxx
Add new functionalities (Laurent)
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALCalibTempCoeff.cxx
CommitLineData
3a282863 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// Objects of this class contain temperature-dependence coefficients
19//
20
21#include <fstream>
22#include <TString.h>
23#include <TFile.h>
24#include <TTree.h>
25
26#include "AliEMCALCalibTempCoeff.h"
27
28using namespace std;
29
30ClassImp(AliEMCALCalibTempCoeff)
31
32//____________________________________________________________________________
33AliEMCALCalibTempCoeff::AliEMCALCalibTempCoeff(const int nSM) :
34 fNSuperModule(nSM),
35 fSuperModuleData()
36{
37 //Default constructor.
38 for (int i=0; i<fNSuperModule; i++) {
39 fSuperModuleData.Add(new AliEMCALSuperModuleCalibTempCoeff(i));
40 }
41 fSuperModuleData.Compress(); // compress the TObjArray
42 fSuperModuleData.SetOwner(kTRUE);
43}
44
45//____________________________________________________________________________
46void AliEMCALCalibTempCoeff::ReadTextCalibTempCoeffInfo(Int_t nSM, const TString &txtFileName,
47 Bool_t swapSides)
48{
49 //Read data from txt file. ; coordinates given on SuperModule basis
50
51 std::ifstream inputFile(txtFileName.Data());
52 if (!inputFile) {
53 printf("AliEMCALCalibTempCoeff::ReadCalibTempCoeffInfo - Cannot open the APD info file %s\n", txtFileName.Data());
54 return;
55 }
56
57 fNSuperModule = nSM;
58
59 Int_t iSM = 0; // SuperModule index
60 Int_t iCol = 0;
61 Int_t iRow = 0;
62
63 // list of values to be read
64 Int_t iSrc = 0;
65 Float_t tempCoeff = 0;
66 // end - all values
67
68 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
69
70 for (Int_t i = 0; i < fNSuperModule; i++) {
71 AliEMCALSuperModuleCalibTempCoeff * t = (AliEMCALSuperModuleCalibTempCoeff*) fSuperModuleData[i];
72 if (!inputFile) {
73 printf("AliEMCALCalibTempCoeff::ReadCalibTempCoeffInfo - Error while reading input file; likely EOF..\n");
74 return;
75 }
76 inputFile >> iSM;
77 t->SetSuperModuleNum(iSM);
78
79 // info for each tower
80 for (Int_t j=0; j<nAPDPerSM; j++) {
81 inputFile >> iCol >> iRow >> tempCoeff >> iSrc;
82
83 // check that input values are not out bounds
84 if (iCol<0 || iCol>(AliEMCALGeoParams::fgkEMCALCols-1) ||
85 iRow<0 || iRow>(AliEMCALGeoParams::fgkEMCALRows-1) ) {
86 printf("AliEMCALCalibTempCoeff::ReadCalibTempCoeffInfo - Error while reading input file; j %d iCol %d iRow %d\n", j, iCol, iRow);
87 return;
88 }
89
90 // assume that this info is already swapped and done for this basis?
91 if (swapSides) {
92 // C side, oriented differently than A side: swap is requested
93 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
94 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
95 }
96
97 t->SetTC(iCol, iRow, tempCoeff);
98 t->SetSrc(iCol, iRow, iSrc);
99 }
100
101 } // i, SuperModule
102
103 inputFile.close();
104
105 return;
106}
107
108//____________________________________________________________________________
109void AliEMCALCalibTempCoeff::WriteTextCalibTempCoeffInfo(const TString &txtFileName,
110 Bool_t swapSides)
111{
112 // write data to txt file. ; coordinates given on SuperModule basis
113
114 std::ofstream outputFile(txtFileName.Data());
115 if (!outputFile) {
116 printf("AliEMCALCalibTempCoeff::WriteCalibTempCoeffInfo - Cannot open the APD output file %s\n", txtFileName.Data());
117 return;
118 }
119
120 Int_t iCol = 0;
121 Int_t iRow = 0;
122
123 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
124 Float_t tempCoeff = 0;
125 Int_t iSrc = 0;
126
127 for (Int_t i = 0; i < fNSuperModule; i++) {
128 AliEMCALSuperModuleCalibTempCoeff * t = (AliEMCALSuperModuleCalibTempCoeff*) fSuperModuleData[i];
129
60f4d62b 130 outputFile << t->GetSuperModuleNum() << endl;
131
3a282863 132 // info for each tower
133 for (Int_t j=0; j<nAPDPerSM; j++) {
134 iCol = j / AliEMCALGeoParams::fgkEMCALRows;
135 iRow = j % AliEMCALGeoParams::fgkEMCALRows;
136
137 tempCoeff = t->GetTC(iCol, iRow);
138 iSrc = t->GetSrc(iCol, iRow);
139
140 if (swapSides) {
141 // C side, oriented differently than A side: swap is requested
142 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
143 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
144 }
145
146 outputFile << iCol << " " << iRow
147 << " " << tempCoeff
148 << " " << iSrc << endl;
149 }
150
151 } // i, SuperModule
152
153 outputFile.close();
154
155 return;
156}
157
158//____________________________________________________________________________
159void AliEMCALCalibTempCoeff::ReadRootCalibTempCoeffInfo(const TString &rootFileName,
160 Bool_t swapSides)
161{
162 //Read data from root file. ; coordinates given on SuperModule basis
163 TFile inputFile(rootFileName, "read");
164
165 TTree *tree = (TTree*) inputFile.Get("tree");
166
167 ReadTreeCalibTempCoeffInfo(tree, swapSides);
168
169 inputFile.Close();
170
171 return;
172}
173
174//____________________________________________________________________________
175void AliEMCALCalibTempCoeff::ReadTreeCalibTempCoeffInfo(TTree *tree,
176 Bool_t swapSides)
177{
178 // how many SuperModule's worth of info do we have?
179 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
180 fNSuperModule = tree->GetEntries();
181
182 Int_t iSM = 0; // SuperModule index
183 // list of values to be read
184 // info for each tower
185 Float_t tempCoeff[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows];
186 Int_t iSrc[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows];
187 // end - all values
188
189 // just to make the initializations of the arrays are done correctly, let's use memset
190 memset(tempCoeff, 0, sizeof(tempCoeff));
191 memset(iSrc, 0, sizeof(iSrc));
192
193 // declare the branches
194 tree->SetBranchAddress("iSM", &iSM);
195 tree->SetBranchAddress("TempCoeff", tempCoeff);
196 tree->SetBranchAddress("Src", iSrc);
197
198 // indices for looping over the towers
199 Int_t iCol = 0;
200 Int_t iRow = 0;
201
202 for (int ient=0; ient<tree->GetEntries(); ient++) {
203 tree->GetEntry(ient);
204
205 // assume the index SuperModules come in order: i=iSM
206 AliEMCALSuperModuleCalibTempCoeff * t = (AliEMCALSuperModuleCalibTempCoeff*) fSuperModuleData[iSM];
207
208 t->SetSuperModuleNum(iSM);
209
210 // third: info for each tower
211 for (Int_t j=0; j<nAPDPerSM; j++) {
212 iCol = j / AliEMCALGeoParams::fgkEMCALRows;
213 iRow = j % AliEMCALGeoParams::fgkEMCALRows;
214
215 // help variables: possibly modified or swapped indices
216 int iColMod = iCol;
217 int iRowMod = iRow;
218 // assume that this info is already swapped and done for this basis?
219 if (swapSides) {
220 // C side, oriented differently than A side: swap is requested
221 iColMod = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
222 iRowMod = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
223 }
224
225 t->SetTC(iColMod, iRowMod, tempCoeff[iCol][iRow]);
226 t->SetSrc(iColMod, iRowMod, iSrc[iCol][iRow]);
227 }
228
229 } // loop over entries
230
231 return;
232}
233
234//____________________________________________________________________________
235void AliEMCALCalibTempCoeff::WriteRootCalibTempCoeffInfo(const TString &rootFileName,
236 Bool_t swapSides)
237{
238 // write data to root file. ; coordinates given on SuperModule basis
239 TFile destFile(rootFileName, "recreate");
240 if (destFile.IsZombie()) {
241 return;
242 }
243 destFile.cd();
244
245 TTree *tree = new TTree("tree","");
246
247 // variables for filling the TTree
248 Int_t iSM = 0; // SuperModule index
249 // list of values to be written
250
251 Float_t tempCoeff[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows];
252 Int_t iSrc[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; // end - all values
253
254 // just to make the initializations of the arrays are done correctly, let's use memset
255 memset(tempCoeff, 0, sizeof(tempCoeff));
256 memset(iSrc, 0, sizeof(iSrc));
257
258 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
259 // for looping over towers
260 Int_t iCol = 0;
261 Int_t iRow = 0;
262
263 // declare the branches
264 // first
265 tree->Branch("iSM", &iSM, "iSM/I");
266 // info for each tower; see if a 2D array works OK or if we'll have to use 1D arrays instead
267 tree->Branch( "TempCoeff", &tempCoeff, Form("TempCoeff[%d][%d]/F", AliEMCALGeoParams::fgkEMCALCols, AliEMCALGeoParams::fgkEMCALRows) );
268 tree->Branch( "Src", &iSrc, Form("Src[%d][%d]/I", AliEMCALGeoParams::fgkEMCALCols, AliEMCALGeoParams::fgkEMCALRows) );
269
270 for (iSM = 0; iSM < fNSuperModule; iSM++) {
271 AliEMCALSuperModuleCalibTempCoeff * t = (AliEMCALSuperModuleCalibTempCoeff*) fSuperModuleData[iSM];
272
273 iSM = t->GetSuperModuleNum();
274
275 // info for each tower
276 for (Int_t j=0; j<nAPDPerSM; j++) {
277 iCol = j / AliEMCALGeoParams::fgkEMCALRows;
278 iRow = j % AliEMCALGeoParams::fgkEMCALRows;
279
280 // help variables: possibly modified or swapped indices
281 int iColMod = iCol;
282 int iRowMod = iRow;
283 // assume that this info is already swapped and done for this basis?
284 if (swapSides) {
285 // C side, oriented differently than A side: swap is requested
286 iColMod = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
287 iRowMod = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
288 }
289
290 tempCoeff[iColMod][iRowMod] = t->GetTC(iCol, iRow);
291 iSrc[iColMod][iRowMod] = t->GetSrc(iCol, iRow);
292 }
293
294 tree->Fill();
295 } // i, SuperModule
296
297 tree->Write();
298 destFile.Close();
299
300 return;
301}
302
303//____________________________________________________________________________
304AliEMCALCalibTempCoeff::~AliEMCALCalibTempCoeff()
305{
306 fSuperModuleData.Delete();
307}
308
309//____________________________________________________________________________
310AliEMCALSuperModuleCalibTempCoeff * AliEMCALCalibTempCoeff::GetSuperModuleCalibTempCoeffNum(Int_t supModIndex)const
afae9650 311{ // getter via index
3a282863 312 for (int i=0; i<fNSuperModule; i++) {
313 AliEMCALSuperModuleCalibTempCoeff * t = (AliEMCALSuperModuleCalibTempCoeff*) fSuperModuleData[i];
314 if (t->GetSuperModuleNum() == supModIndex) {
315 return t;
316 }
317 }
318
319 // if we arrived here, then nothing was found.. just return a NULL pointer
320 return NULL;
321}
322