]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALCalibTempCoeff.cxx
Update the mult corr histograms
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALCalibTempCoeff.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 // 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
28 using namespace std;
29
30 ClassImp(AliEMCALCalibTempCoeff)
31
32 //____________________________________________________________________________
33 AliEMCALCalibTempCoeff::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 //____________________________________________________________________________
46 void 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 //____________________________________________________________________________
109 void 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
130     // info for each tower
131     for (Int_t j=0; j<nAPDPerSM; j++) {
132       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
133       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
134
135       tempCoeff = t->GetTC(iCol, iRow);
136       iSrc = t->GetSrc(iCol, iRow);
137
138       if (swapSides) {
139         // C side, oriented differently than A side: swap is requested
140         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
141         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
142       }
143
144       outputFile << iCol << " " << iRow 
145                  << " " << tempCoeff 
146                  << " " << iSrc << endl;
147     }
148
149   } // i, SuperModule
150
151   outputFile.close();
152
153   return;
154 }
155
156 //____________________________________________________________________________
157 void AliEMCALCalibTempCoeff::ReadRootCalibTempCoeffInfo(const TString &rootFileName,
158                                             Bool_t swapSides)
159 {
160   //Read data from root file. ; coordinates given on SuperModule basis
161   TFile inputFile(rootFileName, "read");  
162
163   TTree *tree = (TTree*) inputFile.Get("tree");
164
165   ReadTreeCalibTempCoeffInfo(tree, swapSides);
166
167   inputFile.Close();
168
169   return;
170 }
171
172 //____________________________________________________________________________
173 void AliEMCALCalibTempCoeff::ReadTreeCalibTempCoeffInfo(TTree *tree,
174                                             Bool_t swapSides)
175 {
176   // how many SuperModule's worth of info do we have?
177   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
178   fNSuperModule = tree->GetEntries();
179
180   Int_t iSM = 0; // SuperModule index
181   // list of values to be read
182   // info for each tower
183   Float_t tempCoeff[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
184   Int_t iSrc[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
185   // end - all values
186
187   // just to make the initializations of the arrays are done correctly, let's use memset
188   memset(tempCoeff, 0, sizeof(tempCoeff)); 
189   memset(iSrc, 0, sizeof(iSrc)); 
190
191   // declare the branches
192   tree->SetBranchAddress("iSM", &iSM);
193   tree->SetBranchAddress("TempCoeff", tempCoeff);
194   tree->SetBranchAddress("Src", iSrc);
195
196   // indices for looping over the towers
197   Int_t iCol = 0;
198   Int_t iRow = 0;
199
200   for (int ient=0; ient<tree->GetEntries(); ient++) {
201     tree->GetEntry(ient);
202
203     // assume the index SuperModules come in order: i=iSM
204     AliEMCALSuperModuleCalibTempCoeff * t = (AliEMCALSuperModuleCalibTempCoeff*) fSuperModuleData[iSM];
205
206     t->SetSuperModuleNum(iSM);
207
208     // third: info for each tower
209     for (Int_t j=0; j<nAPDPerSM; j++) {
210       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
211       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
212
213       // help variables: possibly modified or swapped indices
214       int iColMod = iCol;
215       int iRowMod = iRow;
216       // assume that this info is already swapped and done for this basis?
217       if (swapSides) {
218         // C side, oriented differently than A side: swap is requested
219         iColMod = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
220         iRowMod = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
221       }
222
223       t->SetTC(iColMod, iRowMod, tempCoeff[iCol][iRow]);
224       t->SetSrc(iColMod, iRowMod, iSrc[iCol][iRow]);
225     }
226
227   } // loop over entries
228
229   return;
230 }
231
232 //____________________________________________________________________________
233 void AliEMCALCalibTempCoeff::WriteRootCalibTempCoeffInfo(const TString &rootFileName,
234                                              Bool_t swapSides)
235 {
236   // write data to root file. ; coordinates given on SuperModule basis
237   TFile destFile(rootFileName, "recreate");  
238   if (destFile.IsZombie()) {
239     return;
240   }  
241   destFile.cd();
242
243   TTree *tree = new TTree("tree","");
244
245   // variables for filling the TTree
246   Int_t iSM = 0; // SuperModule index
247   // list of values to be written
248
249   Float_t tempCoeff[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
250   Int_t iSrc[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows];   // end - all values
251
252   // just to make the initializations of the arrays are done correctly, let's use memset
253   memset(tempCoeff, 0, sizeof(tempCoeff)); 
254   memset(iSrc, 0, sizeof(iSrc)); 
255
256   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
257   // for looping over towers
258   Int_t iCol = 0;
259   Int_t iRow = 0;
260
261   // declare the branches
262   // first
263   tree->Branch("iSM", &iSM, "iSM/I");
264   // info for each tower; see if a 2D array works OK or if we'll have to use 1D arrays instead 
265   tree->Branch( "TempCoeff", &tempCoeff, Form("TempCoeff[%d][%d]/F", AliEMCALGeoParams::fgkEMCALCols, AliEMCALGeoParams::fgkEMCALRows) );
266   tree->Branch( "Src", &iSrc, Form("Src[%d][%d]/I", AliEMCALGeoParams::fgkEMCALCols, AliEMCALGeoParams::fgkEMCALRows) );
267
268   for (iSM = 0; iSM < fNSuperModule; iSM++) {
269     AliEMCALSuperModuleCalibTempCoeff * t = (AliEMCALSuperModuleCalibTempCoeff*) fSuperModuleData[iSM];
270
271     iSM = t->GetSuperModuleNum();
272
273     // info for each tower
274     for (Int_t j=0; j<nAPDPerSM; j++) {
275       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
276       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
277
278       // help variables: possibly modified or swapped indices
279       int iColMod = iCol;
280       int iRowMod = iRow;
281       // assume that this info is already swapped and done for this basis?
282       if (swapSides) {
283         // C side, oriented differently than A side: swap is requested
284         iColMod = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
285         iRowMod = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
286       }
287
288       tempCoeff[iColMod][iRowMod] = t->GetTC(iCol, iRow);
289       iSrc[iColMod][iRowMod] = t->GetSrc(iCol, iRow);
290     }
291
292     tree->Fill();
293   } // i, SuperModule
294
295   tree->Write();
296   destFile.Close();
297
298   return;
299 }
300
301 //____________________________________________________________________________
302 AliEMCALCalibTempCoeff::~AliEMCALCalibTempCoeff()
303 {
304   fSuperModuleData.Delete();
305 }
306
307 //____________________________________________________________________________
308 AliEMCALSuperModuleCalibTempCoeff * AliEMCALCalibTempCoeff::GetSuperModuleCalibTempCoeffNum(Int_t supModIndex)const
309 {
310   for (int i=0; i<fNSuperModule; i++) {
311     AliEMCALSuperModuleCalibTempCoeff * t = (AliEMCALSuperModuleCalibTempCoeff*) fSuperModuleData[i];
312     if (t->GetSuperModuleNum() == supModIndex) {
313       return t;
314     }
315   }
316
317   // if we arrived here, then nothing was found.. just return a NULL pointer 
318   return NULL;
319 }
320