]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALCalibTimeDepCorrection.cxx
Correction of SA track rejection
[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..");
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       // assume that this info is already swapped and done for this basis?
135       if (swapSides) {
136         // C side, oriented differently than A side: swap is requested
137         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
138         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
139       }
140
141       // set size of TArray
142       t->GetCorrection(iCol,iRow)->Set(nCorr);
143       for (Int_t k=0; k<nCorr; k++) {
144         inputFile >> correction;
145         // add to TArray
146         t->GetCorrection(iCol,iRow)->AddAt(correction, k);
147       }
148     }
149
150   } // i, SuperModule
151
152   inputFile.close();
153
154   return;
155 }
156
157 //____________________________________________________________________________
158 void AliEMCALCalibTimeDepCorrection::WriteTextInfo(const TString &txtFileName,
159                                                    Bool_t swapSides)
160 {
161   // write data to txt file. ; coordinates given on SuperModule basis
162
163   std::ofstream outputFile(txtFileName.Data());
164   if (!outputFile) {
165     printf("AliEMCALCalibTimeDepCorrection::WriteTextInfo - Cannot open the APD output file %s\n", txtFileName.Data());
166     return;
167   }
168
169   Int_t iCol = 0;
170   Int_t iRow = 0;
171   Int_t nCorr = 0;
172
173   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
174
175   // write header info 
176   outputFile << fStartTime << " " << fNTimeBins << " " << fTimeBinSize << endl;
177
178   for (Int_t i = 0; i < fNSuperModule; i++) {
179     AliEMCALSuperModuleCalibTimeDepCorrection * t = (AliEMCALSuperModuleCalibTimeDepCorrection*) fSuperModuleData[i];
180     outputFile << t->GetSuperModuleNum() << endl;
181
182     for (Int_t j=0; j<nAPDPerSM; j++) {
183       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
184       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
185
186       nCorr = t->GetCorrection(iCol,iRow)->GetSize();
187
188       if (swapSides) {
189         // C side, oriented differently than A side: swap is requested
190         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
191         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
192       }
193
194       outputFile << iCol << " " << iRow << " " << nCorr << endl;
195       for (Int_t k=0; k<nCorr; k++) {
196         outputFile << t->GetCorrection(iCol,iRow)->At(k) << " ";
197       }
198       outputFile << endl;
199
200     }
201
202   } // i, SuperModule
203
204   outputFile.close();
205
206   return;
207 }
208
209 //____________________________________________________________________________
210 void AliEMCALCalibTimeDepCorrection::ReadRootInfo(const TString &rootFileName,
211                                                   Bool_t swapSides)
212 {
213   //Read data from root file. ; coordinates given on SuperModule basis
214   TFile inputFile(rootFileName, "read");  
215
216   TTree *treeGlob = (TTree*) inputFile.Get("treeGlob");
217   TTree *treeCorr = (TTree*) inputFile.Get("treeCorr");
218
219   ReadTreeInfo(treeGlob, treeCorr, swapSides);
220
221   inputFile.Close();
222
223   return;
224 }
225
226 //____________________________________________________________________________
227 void AliEMCALCalibTimeDepCorrection::ReadTreeInfo(TTree *treeGlob, TTree *treeCorr,
228                                                   Bool_t swapSides)
229 {
230  // how many SuperModule's worth of entries / APDs do we have?
231   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
232   fNSuperModule = treeCorr->GetEntries() / nAPDPerSM;
233
234   // global variables : only one entry
235   treeGlob->SetBranchAddress("fStartTime", &fStartTime);
236   treeGlob->SetBranchAddress("fNTimeBins", &fNTimeBins);
237   treeGlob->SetBranchAddress("fTimeBinSize", &fTimeBinSize);
238   treeGlob->GetEntry(0); 
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 = (AliEMCALSuperModuleCalibTimeDepCorrection*) fSuperModuleData[iSM];
262     t->SetSuperModuleNum(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->GetCorrection(iCol,iRow)->Set(nCorr);
274     for (Int_t k=0; k<nCorr; k++) {
275       // add to TArray
276       t->GetCorrection(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 = (AliEMCALSuperModuleCalibTimeDepCorrection*) 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->GetCorrection(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->GetCorrection(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   fSuperModuleData.Delete();
361 }
362
363 //____________________________________________________________________________
364 AliEMCALSuperModuleCalibTimeDepCorrection * AliEMCALCalibTimeDepCorrection::GetSuperModuleCalibTimeDepCorrectionNum(Int_t supModIndex)const
365 {
366   for (int i=0; i<fNSuperModule; i++) {
367     AliEMCALSuperModuleCalibTimeDepCorrection * t = (AliEMCALSuperModuleCalibTimeDepCorrection*) fSuperModuleData[i];
368     if (t->GetSuperModuleNum() == supModIndex) {
369       return t;
370     }
371   }
372
373   // if we arrived here, then nothing was found.. just return a NULL pointer 
374   return NULL;
375 }
376