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