]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALCalibReference.cxx
Coverity corrections
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALCalibReference.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 reference calibrations
19 //
20
21 #include <fstream>
22 #include <TString.h>
23 #include <TFile.h>
24 #include <TTree.h>
25
26 #include "AliEMCALCalibReference.h"
27
28 using namespace std;
29
30 ClassImp(AliEMCALCalibReference)
31
32 //____________________________________________________________________________
33 AliEMCALCalibReference::AliEMCALCalibReference(const int nSM) : 
34   fNSuperModule(nSM),
35   fSuperModuleData()
36 {
37   //Default constructor.
38   for (int i=0; i<fNSuperModule; i++) {
39     fSuperModuleData.Add(new AliEMCALSuperModuleCalibReference(i));
40   }
41   fSuperModuleData.Compress(); // compress the TObjArray
42   fSuperModuleData.SetOwner(kTRUE); 
43 }
44
45 //____________________________________________________________________________
46 void AliEMCALCalibReference::ReadTextCalibReferenceInfo(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("AliEMCALCalibReference::ReadCalibReferenceInfo - 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   Int_t id = 0;
63
64   // list of values to be read
65   // first: overall values for the whole SuperModule
66   Int_t iReferenceTime = 0; 
67   // second: additional info for LED Reference and SM temperature
68   Float_t rLEDRefAmp = 0;
69   Float_t rLEDRefAmpRMS = 0;
70   Int_t iLEDRefHighLow = 0;
71   Float_t temperature = 0;
72   Float_t temperatureRMS = 0;
73   // third: info for each tower
74   Int_t iHighLow = 0; // 
75   Float_t rLEDAmp = 0; // low gain eq. amplitude
76   Float_t rLEDAmpRMS = 0; //
77   // end - all values
78
79   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
80
81   for (Int_t i = 0; i < fNSuperModule; i++) {
82     AliEMCALSuperModuleCalibReference * t = (AliEMCALSuperModuleCalibReference*) fSuperModuleData[i];
83     if (!inputFile) {
84       printf("AliEMCALCalibReference::ReadCalibReferenceInfo - Error while reading input file; likely EOF..");
85       return;
86     }
87     inputFile >> iSM;
88     t->SetSuperModuleNum(iSM);
89
90     // first: overall values for the whole SuperModule
91     inputFile >> iReferenceTime;
92     t->SetReferenceTime(iReferenceTime);
93
94     // second: additional info for LED Reference and SM temperature
95     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALLEDRefs; j++) {
96       inputFile >> id >> iLEDRefHighLow >> rLEDRefAmp >> rLEDRefAmpRMS;
97       t->SetLEDRefHighLow(id, iLEDRefHighLow);
98       t->SetLEDRefAmp(id, rLEDRefAmp);
99       t->SetLEDRefAmpRMS(id, rLEDRefAmpRMS);
100     }
101     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALTempSensors; j++) {
102       inputFile >> id >> temperature >> temperatureRMS;
103       t->SetTemperature(id, temperature);
104       t->SetTemperatureRMS(id, temperatureRMS);
105     }
106
107     // third: info for each tower
108     for (Int_t j=0; j<nAPDPerSM; j++) {
109       inputFile >> iCol >> iRow 
110                 >> iHighLow >> rLEDAmp >> rLEDAmpRMS;
111
112       // assume that this info is already swapped and done for this basis?
113       if (swapSides) {
114         // C side, oriented differently than A side: swap is requested
115         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
116         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
117       }
118
119       AliEMCALCalibReferenceVal * v = t->GetAPDVal(iCol, iRow);
120
121       v->SetHighLow(iHighLow);
122       v->SetLEDAmp(rLEDAmp);
123       v->SetLEDAmpRMS(rLEDAmpRMS);
124     }
125
126   } // i, SuperModule
127
128   inputFile.close();
129
130   return;
131 }
132
133 //____________________________________________________________________________
134 void AliEMCALCalibReference::WriteTextCalibReferenceInfo(const TString &txtFileName,
135                                              Bool_t swapSides)
136 {
137   // write data to txt file. ; coordinates given on SuperModule basis
138
139   std::ofstream outputFile(txtFileName.Data());
140   if (!outputFile) {
141     printf("AliEMCALCalibReference::WriteCalibReferenceInfo - Cannot open the APD output file %s\n", txtFileName.Data());
142     return;
143   }
144
145   Int_t iCol = 0;
146   Int_t iRow = 0;
147
148   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
149
150   for (Int_t i = 0; i < fNSuperModule; i++) {
151     AliEMCALSuperModuleCalibReference * t = (AliEMCALSuperModuleCalibReference*) fSuperModuleData[i];
152
153     // first: overall values for the whole SuperModule
154     outputFile << t->GetSuperModuleNum() << endl;
155     outputFile << t->GetReferenceTime() << endl;
156
157     // second: additional info for LED Reference and SM temperature
158     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALLEDRefs; j++) {
159       outputFile << j << " " << t->GetLEDRefHighLow(j) 
160                  << " " << t->GetLEDRefAmp(j) << " " << t->GetLEDRefAmpRMS(j) 
161                  << endl;
162     }
163     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALTempSensors; j++) {
164       outputFile << j << " " << t->GetTemperature(j) << " " << t->GetTemperatureRMS(j) << endl;
165     }
166
167     // third: info for each tower
168     for (Int_t j=0; j<nAPDPerSM; j++) {
169       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
170       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
171
172       AliEMCALCalibReferenceVal * v = t->GetAPDVal(iCol, iRow);
173
174       if (swapSides) {
175         // C side, oriented differently than A side: swap is requested
176         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
177         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
178       }
179
180       outputFile << iCol << " " << iRow 
181                  << " " << v->GetHighLow() 
182                  << " " << v->GetLEDAmp() 
183                  << " " << v->GetLEDAmpRMS() << endl;
184     }
185
186   } // i, SuperModule
187
188   outputFile.close();
189
190   return;
191 }
192
193 //____________________________________________________________________________
194 void AliEMCALCalibReference::ReadRootCalibReferenceInfo(const TString &rootFileName,
195                                             Bool_t swapSides)
196 {
197   //Read data from root file. ; coordinates given on SuperModule basis
198   TFile inputFile(rootFileName, "read");  
199
200   TTree *tree = (TTree*) inputFile.Get("tree");
201
202   ReadTreeCalibReferenceInfo(tree, swapSides);
203
204   inputFile.Close();
205
206   return;
207 }
208
209 //____________________________________________________________________________
210 void AliEMCALCalibReference::ReadTreeCalibReferenceInfo(TTree *tree,
211                                             Bool_t swapSides)
212 {
213   // how many SuperModule's worth of info do we have?
214   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
215   fNSuperModule = tree->GetEntries();
216
217   Int_t iSM = 0; // SuperModule index
218   // list of values to be read
219   // first: overall values for the whole SuperModule
220   Int_t iReferenceTime= 0; 
221   // second: additional info for LED Reference and SM temperature
222   Float_t rLEDRefAmp[AliEMCALGeoParams::fgkEMCALLEDRefs]= {0};
223   Float_t rLEDRefAmpRMS[AliEMCALGeoParams::fgkEMCALLEDRefs]= {0};
224   Int_t iLEDRefHighLow[AliEMCALGeoParams::fgkEMCALLEDRefs]= {0};
225   Float_t temperature[AliEMCALGeoParams::fgkEMCALTempSensors]= {0};
226   Float_t temperatureRMS[AliEMCALGeoParams::fgkEMCALTempSensors]= {0};
227   // third: info for each tower
228   Int_t iHighLow[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
229   Float_t rLEDAmp[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
230   Float_t rLEDAmpRMS[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
231   // end - all values
232
233   // just to make the initializations of the arrays are done correctly, let's use memset
234   memset(rLEDRefAmp, 0, sizeof(rLEDRefAmp)); 
235   memset(rLEDRefAmpRMS, 0, sizeof(rLEDRefAmpRMS)); 
236   memset(iLEDRefHighLow, 0, sizeof(iLEDRefHighLow)); 
237   memset(temperature, 0, sizeof(temperature)); 
238   memset(temperatureRMS, 0, sizeof(temperatureRMS)); 
239   memset(iHighLow, 0, sizeof(iHighLow)); 
240   memset(rLEDAmp, 0, sizeof(rLEDAmp)); 
241   memset(rLEDAmpRMS, 0, sizeof(rLEDAmpRMS)); 
242
243   // declare the branches
244   tree->SetBranchAddress("iSM", &iSM);
245   tree->SetBranchAddress("ReferenceTime", &iReferenceTime);
246   //
247   tree->SetBranchAddress("LEDRefAmp", rLEDRefAmp);
248   tree->SetBranchAddress("LEDRefAmpRMS", rLEDRefAmpRMS);
249   tree->SetBranchAddress("LEDRefHighLow", iLEDRefHighLow);
250   tree->SetBranchAddress("Temperature", temperature);
251   tree->SetBranchAddress("TemperatureRMS", temperatureRMS);
252   //
253   tree->SetBranchAddress("HighLow", iHighLow);
254   tree->SetBranchAddress("LEDAmp", rLEDAmp);
255   tree->SetBranchAddress("LEDAmpRMS", rLEDAmpRMS);
256
257   // indices for looping over the towers
258   Int_t iCol = 0;
259   Int_t iRow = 0;
260
261   for (int ient=0; ient<tree->GetEntries(); ient++) {
262     tree->GetEntry(ient);
263
264     // assume the index SuperModules come in order: i=iSM
265     AliEMCALSuperModuleCalibReference * t = (AliEMCALSuperModuleCalibReference*) fSuperModuleData[iSM];
266
267     t->SetSuperModuleNum(iSM);
268     // first, overall values
269     t->SetReferenceTime(iReferenceTime);
270
271     // second: additional info for LED references and SM temperatures
272     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALLEDRefs; j++) {
273       t->SetLEDRefAmp(j, rLEDRefAmp[j]);
274       t->SetLEDRefAmpRMS(j, rLEDRefAmpRMS[j]);
275       t->SetLEDRefHighLow(j, iLEDRefHighLow[j]);
276     }
277     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALTempSensors; j++) {
278       t->SetTemperature(j, temperature[j]);
279       t->SetTemperatureRMS(j, temperatureRMS[j]);
280     }
281
282     // third: info for each tower
283     for (Int_t j=0; j<nAPDPerSM; j++) {
284       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
285       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
286
287       // help variables: possibly modified or swapped indices
288       int iColMod = iCol;
289       int iRowMod = iRow;
290       // assume that this info is already swapped and done for this basis?
291       if (swapSides) {
292         // C side, oriented differently than A side: swap is requested
293         iColMod = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
294         iRowMod = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
295       }
296
297       AliEMCALCalibReferenceVal * v = t->GetAPDVal(iColMod, iRowMod);
298
299       v->SetHighLow(iHighLow[iCol][iRow]);
300       v->SetLEDAmp(rLEDAmp[iCol][iRow]);
301       v->SetLEDAmpRMS(rLEDAmpRMS[iCol][iRow]);
302     }
303
304   } // loop over entries
305
306   return;
307 }
308
309 //____________________________________________________________________________
310 void AliEMCALCalibReference::WriteRootCalibReferenceInfo(const TString &rootFileName,
311                                              Bool_t swapSides)
312 {
313   // write data to root file. ; coordinates given on SuperModule basis
314   TFile destFile(rootFileName, "recreate");  
315   if (destFile.IsZombie()) {
316     return;
317   }  
318   destFile.cd();
319
320   TTree *tree = new TTree("tree","");
321
322   // variables for filling the TTree
323   Int_t iSM = 0; // SuperModule index
324   // list of values to be written
325   // first: overall values for the whole SuperModule
326   Int_t iReferenceTime = 0; 
327   // second: additional info for LED Reference and SM temperature
328   Float_t rLEDRefAmp[AliEMCALGeoParams::fgkEMCALLEDRefs] = {0};
329   Float_t rLEDRefAmpRMS[AliEMCALGeoParams::fgkEMCALLEDRefs]= {0};
330   Int_t iLEDRefHighLow[AliEMCALGeoParams::fgkEMCALLEDRefs]= {0};
331   Float_t temperature[AliEMCALGeoParams::fgkEMCALTempSensors]= {0};
332   Float_t temperatureRMS[AliEMCALGeoParams::fgkEMCALTempSensors]= {0};
333   // third: info for each tower
334   Int_t iHighLow[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
335   Float_t rLEDAmp[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
336   Float_t rLEDAmpRMS[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
337   // end - all values
338
339   // just to make the initializations of the arrays are done correctly, let's use memset
340   memset(rLEDRefAmp, 0, sizeof(rLEDRefAmp)); 
341   memset(rLEDRefAmpRMS, 0, sizeof(rLEDRefAmpRMS)); 
342   memset(iLEDRefHighLow, 0, sizeof(iLEDRefHighLow)); 
343   memset(temperature, 0, sizeof(temperature)); 
344   memset(temperatureRMS, 0, sizeof(temperatureRMS)); 
345   memset(iHighLow, 0, sizeof(iHighLow)); 
346   memset(rLEDAmp, 0, sizeof(rLEDAmp)); 
347   memset(rLEDAmpRMS, 0, sizeof(rLEDAmpRMS)); 
348
349   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
350   // for looping over towers
351   Int_t iCol = 0;
352   Int_t iRow = 0;
353
354   // declare the branches
355   // first
356   tree->Branch("iSM", &iSM, "iSM/I");
357   tree->Branch("ReferenceTime", &iReferenceTime, "ReferenceTime/I");
358   // second  
359   tree->Branch( "LEDRefAmp", &rLEDRefAmp, Form("LEDRefAmp[%d]/F", AliEMCALGeoParams::fgkEMCALLEDRefs) );
360   tree->Branch( "LEDRefAmpRMS", &rLEDRefAmpRMS, Form("LEDRefAmpRMS[%d]/F", AliEMCALGeoParams::fgkEMCALLEDRefs) );
361   tree->Branch( "LEDRefHighLow", &iLEDRefHighLow, Form("LEDRefHighLow[%d]/I", AliEMCALGeoParams::fgkEMCALLEDRefs) );
362   tree->Branch( "Temperature", &temperature, Form("Temperature[%d]/F", AliEMCALGeoParams::fgkEMCALTempSensors) );
363   tree->Branch( "TemperatureRMS", &temperatureRMS, Form("TemperatureRMS[%d]/F", AliEMCALGeoParams::fgkEMCALTempSensors) );
364   // third: info for each tower; see if a 2D array works OK or if we'll have to use 1D arrays instead 
365   tree->Branch( "HighLow", &iHighLow, Form("HighLow[%d][%d]/I", AliEMCALGeoParams::fgkEMCALCols, AliEMCALGeoParams::fgkEMCALRows) );
366   tree->Branch( "LEDAmp", &rLEDAmp, Form("LEDAmp[%d][%d]/F", AliEMCALGeoParams::fgkEMCALCols, AliEMCALGeoParams::fgkEMCALRows) );
367   tree->Branch( "LEDAmpRMS", &rLEDAmpRMS, Form("LEDAmpRMS[%d][%d]/F", AliEMCALGeoParams::fgkEMCALCols, AliEMCALGeoParams::fgkEMCALRows) );
368
369   for (iSM = 0; iSM < fNSuperModule; iSM++) {
370     AliEMCALSuperModuleCalibReference * t = (AliEMCALSuperModuleCalibReference*) fSuperModuleData[iSM];
371
372     iSM = t->GetSuperModuleNum();
373     // first, overall values
374     iReferenceTime = t->GetReferenceTime();
375
376     // second: additional info for LED references and SM temperatures
377     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALLEDRefs; j++) {
378       rLEDRefAmp[j] = t->GetLEDRefAmp(j);
379       rLEDRefAmpRMS[j] = t->GetLEDRefAmpRMS(j);
380       iLEDRefHighLow[j] = t->GetLEDRefHighLow(j);
381     }
382     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALTempSensors; j++) {
383       temperature[j] = t->GetTemperature(j);
384       temperatureRMS[j] = t->GetTemperatureRMS(j);
385     }
386
387     // third: info for each tower
388     for (Int_t j=0; j<nAPDPerSM; j++) {
389       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
390       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
391
392       // help variables: possibly modified or swapped indices
393       int iColMod = iCol;
394       int iRowMod = iRow;
395       // assume that this info is already swapped and done for this basis?
396       if (swapSides) {
397         // C side, oriented differently than A side: swap is requested
398         iColMod = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
399         iRowMod = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
400       }
401
402       AliEMCALCalibReferenceVal * v = t->GetAPDVal(iCol, iRow);
403
404       iHighLow[iColMod][iRowMod] = v->GetHighLow();
405       rLEDAmp[iColMod][iRowMod] = v->GetLEDAmp();
406       rLEDAmpRMS[iColMod][iRowMod] = v->GetLEDAmpRMS();
407     }
408
409     tree->Fill();
410   } // i, SuperModule
411
412   tree->Write();
413   destFile.Close();
414
415   return;
416 }
417
418 //____________________________________________________________________________
419 AliEMCALCalibReference::~AliEMCALCalibReference()
420 {
421   fSuperModuleData.Delete();
422 }
423
424 //____________________________________________________________________________
425 AliEMCALSuperModuleCalibReference * AliEMCALCalibReference::GetSuperModuleCalibReferenceNum(Int_t supModIndex)const
426 {
427   for (int i=0; i<fNSuperModule; i++) {
428     AliEMCALSuperModuleCalibReference * t = (AliEMCALSuperModuleCalibReference*) fSuperModuleData[i];
429     if (t->GetSuperModuleNum() == supModIndex) {
430       return t;
431     }
432   }
433
434   // if we arrived here, then nothing was found.. just return a NULL pointer 
435   return NULL;
436 }
437