]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALCalibTimeDepCorrection.cxx
a macro with TPC,ITS, ITSSPD vertex, global vertex and v0's
[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=1.0)
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     t->SetSuperModuleNum(iSM); // assume SMs are coming in order
64
65     for (Int_t j=0; j<nAPDPerSM; j++) {
66
67       // set size of TArray
68       t->GetCorrection(iCol,iRow)->Set(nCorr);
69       for (Int_t k=0; k<nCorr; k++) {
70         // add to TArray
71         t->GetCorrection(iCol,iRow)->AddAt(correction, k); // AddAt = SetAt..
72       }
73     }
74
75   } // i, SuperModule
76
77   return;
78 }
79
80 //____________________________________________________________________________
81 void AliEMCALCalibTimeDepCorrection::SetCorrection(Int_t supModIndex, Int_t iCol, Int_t iRow, Int_t iBin, Float_t val)
82 { // if you call for non-existing data, there may be a crash..
83   ((AliEMCALSuperModuleCalibTimeDepCorrection*)fSuperModuleData[supModIndex])->GetCorrection(iCol,iRow)->AddAt(val, iBin); // AddAt = SetAt..
84  return; 
85 }
86
87 //____________________________________________________________________________
88 Float_t AliEMCALCalibTimeDepCorrection::GetCorrection(Int_t supModIndex, Int_t iCol, Int_t iRow, Int_t iBin) const
89 { // if you call for non-existing data, there may be a crash..
90   return ((AliEMCALSuperModuleCalibTimeDepCorrection*)fSuperModuleData[supModIndex])->GetCorrection(iCol,iRow)->At(iBin);
91 }
92
93 //____________________________________________________________________________
94 void AliEMCALCalibTimeDepCorrection::ReadTextInfo(Int_t nSM, const TString &txtFileName,
95                                                   Bool_t swapSides)
96 {
97   //Read data from txt file. ; coordinates given on SuperModule basis
98
99   std::ifstream inputFile(txtFileName.Data());
100   if (!inputFile) {
101     printf("AliEMCALCalibTimeDepCorrection::ReadTextInfo - Cannot open the APD info file %s\n", txtFileName.Data());
102     return;
103   }
104
105   fNSuperModule = nSM;
106
107   Int_t iSM = 0; // SuperModule index
108   Int_t iCol = 0;
109   Int_t iRow = 0;
110   Int_t nCorr = 0;
111   Float_t correction = 0;
112
113   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
114
115   // read header info 
116   inputFile >> fStartTime >> fNTimeBins >> fTimeBinSize;
117
118   for (Int_t i = 0; i < fNSuperModule; i++) {
119     AliEMCALSuperModuleCalibTimeDepCorrection * t = (AliEMCALSuperModuleCalibTimeDepCorrection*) fSuperModuleData[i];
120
121     if (!inputFile) {
122       printf("AliEMCALCalibTimeDepCorrection::ReadTextInfo - Error while reading input file; likely EOF..");
123       return;
124     }
125     inputFile >> iSM;
126     t->SetSuperModuleNum(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->GetCorrection(iCol,iRow)->Set(nCorr);
140       for (Int_t k=0; k<nCorr; k++) {
141         inputFile >> correction;
142         // add to TArray
143         t->GetCorrection(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 = (AliEMCALSuperModuleCalibTimeDepCorrection*) fSuperModuleData[i];
177     outputFile << t->GetSuperModuleNum() << 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->GetCorrection(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->GetCorrection(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   Int_t iSM = 0; // SuperModule index
238   Int_t iCol = 0;
239   Int_t iRow = 0;
240   // list of values to be read
241   Int_t nCorr = 0;
242   Float_t correction[fgkMaxTimeBins] = {0};
243   // make sure it's really initialized correctly
244   memset(correction, 0, sizeof(correction)); // better safe than sorry
245   // end - all values
246
247   // declare the branches
248   treeCorr->SetBranchAddress("iSM", &iSM);
249   treeCorr->SetBranchAddress("iCol", &iCol);
250   treeCorr->SetBranchAddress("iRow", &iRow);
251   treeCorr->SetBranchAddress("nCorr", &nCorr);
252   treeCorr->SetBranchAddress("correction", correction);
253
254   for (int ient=0; ient<treeCorr->GetEntries(); ient++) {
255     treeCorr->GetEntry(ient);
256
257     // assume the index SuperModules come in order: i=iSM
258     AliEMCALSuperModuleCalibTimeDepCorrection * t = (AliEMCALSuperModuleCalibTimeDepCorrection*) fSuperModuleData[iSM];
259     t->SetSuperModuleNum(iSM);
260
261     // assume that this info is already swapped and done for this basis?
262     if (swapSides) {
263       // C side, oriented differently than A side: swap is requested
264       iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
265       iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
266     }
267
268
269     // set size of TArray
270     t->GetCorrection(iCol,iRow)->Set(nCorr);
271     for (Int_t k=0; k<nCorr; k++) {
272       // add to TArray
273       t->GetCorrection(iCol,iRow)->AddAt(correction[k], k);
274     }
275
276   } // entry
277
278   return;
279 }
280
281 //____________________________________________________________________________
282 void AliEMCALCalibTimeDepCorrection::WriteRootInfo(const TString &rootFileName,
283                                                    Bool_t swapSides)
284 {
285   // write data to root file. ; coordinates given on SuperModule basis
286   TFile destFile(rootFileName, "recreate");  
287   if (destFile.IsZombie()) {
288     return;
289   }  
290   destFile.cd();
291
292   TTree *treeGlob = new TTree("treeGlob","");
293   TTree *treeCorr = new TTree("treeCorr","");
294
295   // global part only has one entry
296   treeGlob->Branch("fStartTime", &fStartTime, "fStartTime/I"); // really unsigned int..
297   treeGlob->Branch("fNTimeBins", &fNTimeBins, "fNTimeBins/I");
298   treeGlob->Branch("fTimeBinSize", &fTimeBinSize, "fTimeBinSize/I");
299   treeGlob->Fill();
300
301   // variables for filling the TTree
302   Int_t iSM = 0; // SuperModule index
303   Int_t iCol = 0;
304   Int_t iRow = 0;
305   Int_t nCorr = 0;
306   Float_t correction[fgkMaxTimeBins] = {0};
307   // make sure it's really initialized correctly
308   memset(correction, 0, sizeof(correction)); // better safe than sorry
309
310   // declare the branches
311   treeCorr->Branch("iSM", &iSM, "iSM/I");
312   treeCorr->Branch("iCol", &iCol, "iCol/I");
313   treeCorr->Branch("iRow", &iRow, "iRow/I");
314   treeCorr->Branch("nCorr", &nCorr, "nCorr/I");
315   treeCorr->Branch("correction", &correction, "correction[nCorr]/F");
316
317   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
318
319   for (iSM = 0; iSM < fNSuperModule; iSM++) {
320     AliEMCALSuperModuleCalibTimeDepCorrection * t = (AliEMCALSuperModuleCalibTimeDepCorrection*) fSuperModuleData[iSM];
321
322     for (Int_t j=0; j<nAPDPerSM; j++) {
323       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
324       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
325
326       nCorr = t->GetCorrection(iCol,iRow)->GetSize();
327       if (nCorr > fgkMaxTimeBins) {
328         printf("AliEMCALCalibTimeDepCorrection::WriteRootInfo - too many correction/timebins %d kept\n", nCorr);
329         return;
330       }
331
332       if (swapSides) {
333         // C side, oriented differently than A side: swap is requested
334         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
335         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
336       }
337
338       for (Int_t k=0; k<nCorr; k++) {
339         correction[k] = t->GetCorrection(iCol,iRow)->At(k);
340       }
341
342       treeCorr->Fill();
343     }
344
345   } // i, SuperModule
346
347   treeGlob->Write();
348   treeCorr->Write();
349   destFile.Close();
350
351   return;
352 }
353
354 //____________________________________________________________________________
355 AliEMCALCalibTimeDepCorrection::~AliEMCALCalibTimeDepCorrection()
356 {
357   fSuperModuleData.Delete();
358 }
359
360 //____________________________________________________________________________
361 AliEMCALSuperModuleCalibTimeDepCorrection * AliEMCALCalibTimeDepCorrection::GetSuperModuleCalibTimeDepCorrectionNum(Int_t supModIndex)const
362 {
363   for (int i=0; i<fNSuperModule; i++) {
364     AliEMCALSuperModuleCalibTimeDepCorrection * t = (AliEMCALSuperModuleCalibTimeDepCorrection*) fSuperModuleData[i];
365     if (t->GetSuperModuleNum() == supModIndex) {
366       return t;
367     }
368   }
369
370   // if we arrived here, then nothing was found.. just return a NULL pointer 
371   return NULL;
372 }
373