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