intermediate calibration class commit: struct changed to classes and added root i/o
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALCalibTimeDepCorrection.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 info on time-dependent corrections
19 //
20
21 #include <fstream>
22 #include <TString.h>
23 #include <TArrayF.h>
24 #include <TFile.h>
25 #include <TTree.h>
26
27 #include "AliEMCALCalibTimeDepCorrection.h"
28
29 using namespace std;
30
31 ClassImp(AliEMCALCalibTimeDepCorrection)
32
33 //____________________________________________________________________________
34 AliEMCALCalibTimeDepCorrection::AliEMCALCalibTimeDepCorrection() : 
35   fNSuperModule(0),
36   fSuperModuleData(0),
37   fStartTime(0),
38   fNTimeBins(0),
39   fTimeBinSize(0)
40 {
41   //Default constructor.
42 }
43
44 //____________________________________________________________________________
45 void AliEMCALCalibTimeDepCorrection::InitCorrection(Int_t nSM, Int_t nBins, Float_t val=1.0)
46 {
47   // This methods assumes that you are using SuperModules 0..nSM-1
48   fNSuperModule = nSM;
49   if (fSuperModuleData) delete [] fSuperModuleData;
50   fSuperModuleData = new AliEMCALSuperModuleCalibTimeDepCorrection[fNSuperModule];
51
52   Int_t iSM = 0; // SuperModule index
53   Int_t iCol = 0;
54   Int_t iRow = 0;
55   Int_t nCorr = nBins;
56   Float_t Correction = val;
57
58   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
59
60   for (Int_t i = 0; i < fNSuperModule; i++) {
61     AliEMCALSuperModuleCalibTimeDepCorrection &t = fSuperModuleData[i];
62     t.fSuperModuleNum = iSM; // assume SMs are coming in order
63
64     for (Int_t j=0; j<nAPDPerSM; j++) {
65
66       // set size of TArray
67       t.fCorrection[iCol][iRow].Set(nCorr);
68       for (Int_t k=0; k<nCorr; k++) {
69         // add to TArray
70         t.fCorrection[iCol][iRow].AddAt(Correction, k); // AddAt = SetAt..
71       }
72     }
73
74   } // i, SuperModule
75
76   return;
77 }
78
79 //____________________________________________________________________________
80 void AliEMCALCalibTimeDepCorrection::SetCorrection(Int_t supModIndex, Int_t iCol, Int_t iRow, Int_t iBin, Float_t val)
81 { // if you call for non-existing data, there may be a crash..
82   fSuperModuleData[supModIndex].fCorrection[iCol][iRow].AddAt(val, iBin); // AddAt = SetAt..
83  return; 
84 }
85
86 //____________________________________________________________________________
87 Float_t AliEMCALCalibTimeDepCorrection::GetCorrection(Int_t supModIndex, Int_t iCol, Int_t iRow, Int_t iBin)const
88 { // if you call for non-existing data, there may be a crash..
89   return fSuperModuleData[supModIndex].fCorrection[iCol][iRow].At(iBin);
90 }
91
92 //____________________________________________________________________________
93 void AliEMCALCalibTimeDepCorrection::ReadTextInfo(Int_t nSM, const TString &txtFileName,
94                                                   Bool_t swapSides)
95 {
96   //Read data from txt file. ; coordinates given on SuperModule basis
97
98   std::ifstream inputFile(txtFileName.Data());
99   if (!inputFile) {
100     printf("AliEMCALCalibTimeDepCorrection::ReadTextInfo - Cannot open the APD info file %s\n", txtFileName.Data());
101     return;
102   }
103
104   fNSuperModule = nSM;
105   if (fSuperModuleData) delete [] fSuperModuleData;
106   fSuperModuleData = new AliEMCALSuperModuleCalibTimeDepCorrection[fNSuperModule];
107
108   Int_t iSM = 0; // SuperModule index
109   Int_t iCol = 0;
110   Int_t iRow = 0;
111   Int_t nCorr = 0;
112   Float_t Correction = 0;
113
114   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
115
116   // read header info 
117   inputFile >> fStartTime >> fNTimeBins >> fTimeBinSize;
118
119   for (Int_t i = 0; i < fNSuperModule; i++) {
120     AliEMCALSuperModuleCalibTimeDepCorrection &t = fSuperModuleData[i];
121     if (!inputFile) {
122       printf("AliEMCALCalibTimeDepCorrection::ReadTextInfo - Error while reading input file; likely EOF..");
123       return;
124     }
125     inputFile >> iSM;
126     t.fSuperModuleNum = iSM;
127
128     for (Int_t j=0; j<nAPDPerSM; j++) {
129       inputFile >> iCol >> iRow >> nCorr;
130
131       // assume that this info is already swapped and done for this basis?
132       if (swapSides) {
133         // C side, oriented differently than A side: swap is requested
134         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
135         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
136       }
137
138       // set size of TArray
139       t.fCorrection[iCol][iRow].Set(nCorr);
140       for (Int_t k=0; k<nCorr; k++) {
141         inputFile >> Correction;
142         // add to TArray
143         t.fCorrection[iCol][iRow].AddAt(Correction, k);
144       }
145     }
146
147   } // i, SuperModule
148
149   inputFile.close();
150
151   return;
152 }
153
154 //____________________________________________________________________________
155 void AliEMCALCalibTimeDepCorrection::WriteTextInfo(const TString &txtFileName,
156                                                    Bool_t swapSides)
157 {
158   // write data to txt file. ; coordinates given on SuperModule basis
159
160   std::ofstream outputFile(txtFileName.Data());
161   if (!outputFile) {
162     printf("AliEMCALCalibTimeDepCorrection::WriteTextInfo - Cannot open the APD output file %s\n", txtFileName.Data());
163     return;
164   }
165
166   Int_t iCol = 0;
167   Int_t iRow = 0;
168   Int_t nCorr = 0;
169
170   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
171
172   // write header info 
173   outputFile << fStartTime << " " << fNTimeBins << " " << fTimeBinSize << endl;
174
175   for (Int_t i = 0; i < fNSuperModule; i++) {
176     AliEMCALSuperModuleCalibTimeDepCorrection &t = fSuperModuleData[i];
177     outputFile << t.fSuperModuleNum << endl;
178
179     for (Int_t j=0; j<nAPDPerSM; j++) {
180       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
181       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
182
183       nCorr = t.fCorrection[iCol][iRow].GetSize();
184
185       if (swapSides) {
186         // C side, oriented differently than A side: swap is requested
187         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
188         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
189       }
190
191       outputFile << iCol << " " << iRow << " " << nCorr << endl;
192       for (Int_t k=0; k<nCorr; k++) {
193         outputFile << t.fCorrection[iCol][iRow].At(k) << " ";
194       }
195       outputFile << endl;
196
197     }
198
199   } // i, SuperModule
200
201   outputFile.close();
202
203   return;
204 }
205
206 //____________________________________________________________________________
207 void AliEMCALCalibTimeDepCorrection::ReadRootInfo(const TString &rootFileName,
208                                                   Bool_t swapSides)
209 {
210   //Read data from root file. ; coordinates given on SuperModule basis
211   TFile inputFile(rootFileName, "read");  
212
213   TTree *treeGlob = (TTree*) inputFile.Get("treeGlob");
214   TTree *treeCorr = (TTree*) inputFile.Get("treeCorr");
215
216   ReadTreeInfo(treeGlob, treeCorr, swapSides);
217
218   inputFile.Close();
219
220   return;
221 }
222
223 //____________________________________________________________________________
224 void AliEMCALCalibTimeDepCorrection::ReadTreeInfo(TTree *treeGlob, TTree *treeCorr,
225                                                   Bool_t swapSides)
226 {
227  // how many SuperModule's worth of entries / APDs do we have?
228   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
229   fNSuperModule = treeCorr->GetEntries() / nAPDPerSM;
230
231   // global variables : only one entry
232   treeGlob->SetBranchAddress("fStartTime", &fStartTime);
233   treeGlob->SetBranchAddress("fNTimeBins", &fNTimeBins);
234   treeGlob->SetBranchAddress("fTimeBinSize", &fTimeBinSize);
235   treeGlob->GetEntry(0); 
236
237   if (fSuperModuleData) delete [] fSuperModuleData;
238   fSuperModuleData = new AliEMCALSuperModuleCalibTimeDepCorrection[fNSuperModule];
239
240   Int_t iSM = 0; // SuperModule index
241   Int_t iCol = 0;
242   Int_t iRow = 0;
243   // list of values to be read
244   Int_t nCorr = 0;
245   Float_t correction[fgkMaxTimeBins] = {0};
246   // make sure it's really initialized correctly
247   memset(correction, 0, sizeof(correction)); // better safe than sorry
248   // end - all values
249
250   // declare the branches
251   treeCorr->SetBranchAddress("iSM", &iSM);
252   treeCorr->SetBranchAddress("iCol", &iCol);
253   treeCorr->SetBranchAddress("iRow", &iRow);
254   treeCorr->SetBranchAddress("nCorr", &nCorr);
255   treeCorr->SetBranchAddress("correction", correction);
256
257   for (int ient=0; ient<treeCorr->GetEntries(); ient++) {
258     treeCorr->GetEntry(ient);
259
260     // assume the index SuperModules come in order: i=iSM
261     AliEMCALSuperModuleCalibTimeDepCorrection &t = fSuperModuleData[iSM];
262     t.fSuperModuleNum = iSM;
263
264     // assume that this info is already swapped and done for this basis?
265     if (swapSides) {
266       // C side, oriented differently than A side: swap is requested
267       iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
268       iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
269     }
270
271
272     // set size of TArray
273     t.fCorrection[iCol][iRow].Set(nCorr);
274     for (Int_t k=0; k<nCorr; k++) {
275       // add to TArray
276       t.fCorrection[iCol][iRow].AddAt(correction[k], k);
277     }
278
279   } // entry
280
281   return;
282 }
283
284 //____________________________________________________________________________
285 void AliEMCALCalibTimeDepCorrection::WriteRootInfo(const TString &rootFileName,
286                                                    Bool_t swapSides)
287 {
288   // write data to root file. ; coordinates given on SuperModule basis
289   TFile destFile(rootFileName, "recreate");  
290   if (destFile.IsZombie()) {
291     return;
292   }  
293   destFile.cd();
294
295   TTree *treeGlob = new TTree("treeGlob","");
296   TTree *treeCorr = new TTree("treeCorr","");
297
298   // global part only has one entry
299   treeGlob->Branch("fStartTime", &fStartTime, "fStartTime/I"); // really unsigned int..
300   treeGlob->Branch("fNTimeBins", &fNTimeBins, "fNTimeBins/I");
301   treeGlob->Branch("fTimeBinSize", &fTimeBinSize, "fTimeBinSize/I");
302   treeGlob->Fill();
303
304   // variables for filling the TTree
305   Int_t iSM = 0; // SuperModule index
306   Int_t iCol = 0;
307   Int_t iRow = 0;
308   Int_t nCorr = 0;
309   Float_t correction[fgkMaxTimeBins] = {0};
310   // make sure it's really initialized correctly
311   memset(correction, 0, sizeof(correction)); // better safe than sorry
312
313   // declare the branches
314   treeCorr->Branch("iSM", &iSM, "iSM/I");
315   treeCorr->Branch("iCol", &iCol, "iCol/I");
316   treeCorr->Branch("iRow", &iRow, "iRow/I");
317   treeCorr->Branch("nCorr", &nCorr, "nCorr/I");
318   treeCorr->Branch("correction", &correction, "correction[nCorr]/F");
319
320   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
321
322   for (iSM = 0; iSM < fNSuperModule; iSM++) {
323     AliEMCALSuperModuleCalibTimeDepCorrection &t = fSuperModuleData[iSM];
324
325     for (Int_t j=0; j<nAPDPerSM; j++) {
326       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
327       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
328
329       nCorr = t.fCorrection[iCol][iRow].GetSize();
330       if (nCorr > fgkMaxTimeBins) {
331         printf("AliEMCALCalibTimeDepCorrection::WriteRootInfo - too many correction/timebins %d kept\n", nCorr);
332         return;
333       }
334
335       if (swapSides) {
336         // C side, oriented differently than A side: swap is requested
337         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
338         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
339       }
340
341       for (Int_t k=0; k<nCorr; k++) {
342         correction[k] = t.fCorrection[iCol][iRow].At(k);
343       }
344
345       treeCorr->Fill();
346     }
347
348   } // i, SuperModule
349
350   treeGlob->Write();
351   treeCorr->Write();
352   destFile.Close();
353
354   return;
355 }
356
357 //____________________________________________________________________________
358 AliEMCALCalibTimeDepCorrection::~AliEMCALCalibTimeDepCorrection()
359 {
360   delete [] fSuperModuleData;
361 }
362
363 //____________________________________________________________________________
364 AliEMCALSuperModuleCalibTimeDepCorrection AliEMCALCalibTimeDepCorrection::GetSuperModuleCalibTimeDepCorrectionId(Int_t supModIndex)const
365 {
366   AliEMCALSuperModuleCalibTimeDepCorrection t;  // just to maybe prevent a crash, but we are returning something not-initialized so maybe not better really..
367   if (!fSuperModuleData)
368     return t;
369
370   return fSuperModuleData[supModIndex];
371 }
372
373 //____________________________________________________________________________
374 AliEMCALSuperModuleCalibTimeDepCorrection AliEMCALCalibTimeDepCorrection::GetSuperModuleCalibTimeDepCorrectionNum(Int_t supModIndex)const
375 {
376   AliEMCALSuperModuleCalibTimeDepCorrection t;  // just to maybe prevent a crash, but we are returning something not-initialized so maybe not better really..
377   if (!fSuperModuleData)
378     return t;
379
380   for (int i=0; i<fNSuperModule; i++) {
381     if (fSuperModuleData[i].fSuperModuleNum == supModIndex) {
382       return fSuperModuleData[i];
383     }
384   }
385
386   return t;
387 }
388