]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALCalibAbs.cxx
COVERITY reports many FORWARD_NULL, corrected; AliEMCALv2: initialization of fHits...
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALCalibAbs.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 basis for absolute calibrations
19 //
20
21 #include <fstream>
22 #include <TString.h>
23 #include <TFile.h>
24 #include <TTree.h>
25
26 #include "AliEMCALCalibAbs.h"
27
28 using namespace std;
29
30 ClassImp(AliEMCALCalibAbs)
31
32 //____________________________________________________________________________
33 AliEMCALCalibAbs::AliEMCALCalibAbs(const int nSM) : 
34   fNSuperModule(nSM),
35   fSuperModuleData()
36 {
37   //Default constructor.
38   for (int i=0; i<fNSuperModule; i++) {
39     fSuperModuleData.Add(new AliEMCALSuperModuleCalibAbs(i));
40   }
41   fSuperModuleData.Compress(); // compress the TObjArray
42   fSuperModuleData.SetOwner(kTRUE); 
43 }
44
45 //____________________________________________________________________________
46 void AliEMCALCalibAbs::ReadTextCalibAbsInfo(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("AliEMCALCalibAbs::ReadCalibAbsInfo - 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   // first: overall values for the whole SuperModule
65   Int_t iCalibMethod = 0; 
66   Int_t iCalibPass = 0; 
67   Float_t absoluteCalib = 0; 
68   // third: info for each tower
69   Float_t relativeCalib = 0; // (ADC>GeV relative gain/conversion), value around 1
70   // end - all values
71
72   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
73
74   for (Int_t i = 0; i < fNSuperModule; i++) {
75     AliEMCALSuperModuleCalibAbs * t = (AliEMCALSuperModuleCalibAbs*) fSuperModuleData[i];
76     if (!inputFile) {
77       printf("AliEMCALCalibAbs::ReadCalibAbsInfo - Error while reading input file; likely EOF..");
78       return;
79     }
80     inputFile >> iSM;
81     t->SetSuperModuleNum(iSM);
82
83     // first: overall values for the whole SuperModule
84     inputFile >> iCalibMethod >> iCalibPass >> absoluteCalib;
85     t->SetCalibMethod(iCalibMethod);
86     t->SetCalibPass(iCalibPass);
87     t->SetAbsoluteCalib(absoluteCalib);
88
89     // third: info for each tower
90     for (Int_t j=0; j<nAPDPerSM; j++) {
91       inputFile >> iCol >> iRow 
92                 >> relativeCalib;
93
94       // assume that this info is already swapped and done for this basis?
95       if (swapSides) {
96         // C side, oriented differently than A side: swap is requested
97         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
98         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
99       }
100
101       t->SetRelativeCalib(iCol, iRow, relativeCalib);
102     }
103
104   } // i, SuperModule
105
106   inputFile.close();
107
108   return;
109 }
110
111 //____________________________________________________________________________
112 void AliEMCALCalibAbs::WriteTextCalibAbsInfo(const TString &txtFileName,
113                                              Bool_t swapSides)
114 {
115   // write data to txt file. ; coordinates given on SuperModule basis
116
117   std::ofstream outputFile(txtFileName.Data());
118   if (!outputFile) {
119     printf("AliEMCALCalibAbs::WriteCalibAbsInfo - Cannot open the APD output file %s\n", txtFileName.Data());
120     return;
121   }
122
123   Int_t iCol = 0;
124   Int_t iRow = 0;
125
126   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
127   Float_t relativeCalib = 0;
128   for (Int_t i = 0; i < fNSuperModule; i++) {
129     AliEMCALSuperModuleCalibAbs * t = (AliEMCALSuperModuleCalibAbs*) fSuperModuleData[i];
130
131     // first: overall values for the whole SuperModule
132     outputFile << t->GetSuperModuleNum() << endl;
133     outputFile << t->GetCalibMethod() << " " 
134                << t->GetCalibPass() << " " 
135                << t->GetAbsoluteCalib() << endl;
136
137     // third: info for each tower
138     for (Int_t j=0; j<nAPDPerSM; j++) {
139       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
140       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
141
142       relativeCalib = t->GetRelativeCalib(iCol, iRow);
143
144       if (swapSides) {
145         // C side, oriented differently than A side: swap is requested
146         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
147         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
148       }
149
150       outputFile << iCol << " " << iRow 
151                  << " " << relativeCalib << endl;
152     }
153
154   } // i, SuperModule
155
156   outputFile.close();
157
158   return;
159 }
160
161 //____________________________________________________________________________
162 void AliEMCALCalibAbs::ReadRootCalibAbsInfo(const TString &rootFileName,
163                                             Bool_t swapSides)
164 {
165   //Read data from root file. ; coordinates given on SuperModule basis
166   TFile inputFile(rootFileName, "read");  
167
168   TTree *tree = (TTree*) inputFile.Get("tree");
169
170   ReadTreeCalibAbsInfo(tree, swapSides);
171
172   inputFile.Close();
173
174   return;
175 }
176
177 //____________________________________________________________________________
178 void AliEMCALCalibAbs::ReadTreeCalibAbsInfo(TTree *tree,
179                                             Bool_t swapSides)
180 {
181   // how many SuperModule's worth of info do we have?
182   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
183   fNSuperModule = tree->GetEntries();
184
185   Int_t iSM = 0; // SuperModule index
186   // list of values to be read
187   // first: overall values for the whole SuperModule
188   Int_t iCalibMethod = 0; 
189   Int_t iCalibPass = 0; 
190   Float_t absoluteCalib = 0; 
191   // third: info for each tower
192   Float_t relativeCalib[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
193   // end - all values
194
195   // just to make the initializations of the arrays are done correctly, let's use memset
196   memset(relativeCalib, 0, sizeof(relativeCalib)); 
197
198   // declare the branches
199   tree->SetBranchAddress("iSM", &iSM);
200   tree->SetBranchAddress("CalibMethod", &iCalibMethod);
201   tree->SetBranchAddress("CalibPass", &iCalibPass);
202   tree->SetBranchAddress("AbsoluteCalib", &absoluteCalib);
203   //
204   tree->SetBranchAddress("RelativeCalib", relativeCalib);
205
206   // indices for looping over the towers
207   Int_t iCol = 0;
208   Int_t iRow = 0;
209
210   for (int ient=0; ient<tree->GetEntries(); ient++) {
211     tree->GetEntry(ient);
212
213     // assume the index SuperModules come in order: i=iSM
214     AliEMCALSuperModuleCalibAbs * t = (AliEMCALSuperModuleCalibAbs*) fSuperModuleData[iSM];
215
216     t->SetSuperModuleNum(iSM);
217     // first, overall values
218     t->SetCalibMethod(iCalibMethod);
219     t->SetCalibPass(iCalibPass);
220     t->SetAbsoluteCalib(absoluteCalib);
221
222     // third: info for each tower
223     for (Int_t j=0; j<nAPDPerSM; j++) {
224       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
225       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
226
227       // help variables: possibly modified or swapped indices
228       int iColMod = iCol;
229       int iRowMod = iRow;
230       // assume that this info is already swapped and done for this basis?
231       if (swapSides) {
232         // C side, oriented differently than A side: swap is requested
233         iColMod = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
234         iRowMod = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
235       }
236
237       t->SetRelativeCalib(iColMod, iRowMod, relativeCalib[iCol][iRow]);
238     }
239
240   } // loop over entries
241
242   return;
243 }
244
245 //____________________________________________________________________________
246 void AliEMCALCalibAbs::WriteRootCalibAbsInfo(const TString &rootFileName,
247                                              Bool_t swapSides)
248 {
249   // write data to root file. ; coordinates given on SuperModule basis
250   TFile destFile(rootFileName, "recreate");  
251   if (destFile.IsZombie()) {
252     return;
253   }  
254   destFile.cd();
255
256   TTree *tree = new TTree("tree","");
257
258   // variables for filling the TTree
259   Int_t iSM = 0; // SuperModule index
260   // list of values to be written
261   // first: overall values for the whole SuperModule
262   Int_t iCalibMethod = 0; 
263   Int_t iCalibPass = 0; 
264   Float_t absoluteCalib = 0; 
265   // third: info for each tower
266   Float_t relativeCalib[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
267   // end - all values
268
269   // just to make the initializations of the arrays are done correctly, let's use memset
270   memset(relativeCalib, 0, sizeof(relativeCalib)); 
271
272   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
273   // for looping over towers
274   Int_t iCol = 0;
275   Int_t iRow = 0;
276
277   // declare the branches
278   // first
279   tree->Branch("iSM", &iSM, "iSM/I");
280   tree->Branch("CalibMethod", &iCalibMethod, "CalibMethod/I");
281   tree->Branch("CalibPass", &iCalibPass, "CalibPass/I");
282   tree->Branch("AbsoluteCalib", &absoluteCalib, "AbsoluteCalib/F");
283   // third: info for each tower; see if a 2D array works OK or if we'll have to use 1D arrays instead 
284   tree->Branch( "RelativeCalib", &relativeCalib, Form("RelativeCalib[%d][%d]/F", AliEMCALGeoParams::fgkEMCALCols, AliEMCALGeoParams::fgkEMCALRows) );
285
286   for (iSM = 0; iSM < fNSuperModule; iSM++) {
287     AliEMCALSuperModuleCalibAbs * t = (AliEMCALSuperModuleCalibAbs*) fSuperModuleData[iSM];
288
289     iSM = t->GetSuperModuleNum();
290     // first, overall values
291     iCalibMethod = t->GetCalibMethod();
292     iCalibPass = t->GetCalibPass();
293     absoluteCalib = t->GetAbsoluteCalib();
294
295     // third: info for each tower
296     for (Int_t j=0; j<nAPDPerSM; j++) {
297       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
298       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
299
300       // help variables: possibly modified or swapped indices
301       int iColMod = iCol;
302       int iRowMod = iRow;
303       // assume that this info is already swapped and done for this basis?
304       if (swapSides) {
305         // C side, oriented differently than A side: swap is requested
306         iColMod = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
307         iRowMod = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
308       }
309
310       relativeCalib[iColMod][iRowMod] = t->GetRelativeCalib(iCol, iRow);
311     }
312
313     tree->Fill();
314   } // i, SuperModule
315
316   tree->Write();
317   destFile.Close();
318
319   return;
320 }
321
322 //____________________________________________________________________________
323 AliEMCALCalibAbs::~AliEMCALCalibAbs()
324 {
325   fSuperModuleData.Delete();
326 }
327
328 //____________________________________________________________________________
329 AliEMCALSuperModuleCalibAbs * AliEMCALCalibAbs::GetSuperModuleCalibAbsNum(Int_t supModIndex)const
330 {
331   for (int i=0; i<fNSuperModule; i++) {
332     AliEMCALSuperModuleCalibAbs * t = (AliEMCALSuperModuleCalibAbs*) fSuperModuleData[i];
333     if (t->GetSuperModuleNum() == supModIndex) {
334       return t;
335     }
336   }
337
338   // if we arrived here, then nothing was found.. just return a NULL pointer 
339   return NULL;
340 }
341