1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 // Objects of this class contain basis for absolute calibrations
26 #include "AliEMCALCalibAbs.h"
30 ClassImp(AliEMCALCalibAbs)
32 //____________________________________________________________________________
33 AliEMCALCalibAbs::AliEMCALCalibAbs(const int nSM) :
37 //Default constructor.
38 for (int i=0; i<fNSuperModule; i++) {
39 fSuperModuleData.Add(new AliEMCALSuperModuleCalibAbs(i));
41 fSuperModuleData.Compress(); // compress the TObjArray
42 fSuperModuleData.SetOwner(kTRUE);
45 //____________________________________________________________________________
46 void AliEMCALCalibAbs::ReadTextCalibAbsInfo(Int_t nSM, const TString &txtFileName,
49 //Read data from txt file. ; coordinates given on SuperModule basis
51 std::ifstream inputFile(txtFileName.Data());
53 printf("AliEMCALCalibAbs::ReadCalibAbsInfo - Cannot open the APD info file %s\n", txtFileName.Data());
59 Int_t iSM = 0; // SuperModule index
63 // list of values to be read
64 // first: overall values for the whole SuperModule
65 Int_t iCalibMethod = 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
72 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
74 for (Int_t i = 0; i < fNSuperModule; i++) {
75 AliEMCALSuperModuleCalibAbs * t = (AliEMCALSuperModuleCalibAbs*) fSuperModuleData[i];
77 printf("AliEMCALCalibAbs::ReadCalibAbsInfo - Error while reading input file; likely EOF..");
81 t->SetSuperModuleNum(iSM);
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);
89 // third: info for each tower
90 for (Int_t j=0; j<nAPDPerSM; j++) {
91 inputFile >> iCol >> iRow
94 // assume that this info is already swapped and done for this basis?
96 // C side, oriented differently than A side: swap is requested
97 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
98 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
101 t->SetRelativeCalib(iCol, iRow, relativeCalib);
111 //____________________________________________________________________________
112 void AliEMCALCalibAbs::WriteTextCalibAbsInfo(const TString &txtFileName,
115 // write data to txt file. ; coordinates given on SuperModule basis
117 std::ofstream outputFile(txtFileName.Data());
119 printf("AliEMCALCalibAbs::WriteCalibAbsInfo - Cannot open the APD output file %s\n", txtFileName.Data());
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];
131 // first: overall values for the whole SuperModule
132 outputFile << t->GetSuperModuleNum() << endl;
133 outputFile << t->GetCalibMethod() << " "
134 << t->GetCalibPass() << " "
135 << t->GetAbsoluteCalib() << endl;
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;
142 relativeCalib = t->GetRelativeCalib(iCol, iRow);
145 // C side, oriented differently than A side: swap is requested
146 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
147 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
150 outputFile << iCol << " " << iRow
151 << " " << relativeCalib << endl;
161 //____________________________________________________________________________
162 void AliEMCALCalibAbs::ReadRootCalibAbsInfo(const TString &rootFileName,
165 //Read data from root file. ; coordinates given on SuperModule basis
166 TFile inputFile(rootFileName, "read");
168 TTree *tree = (TTree*) inputFile.Get("tree");
170 ReadTreeCalibAbsInfo(tree, swapSides);
177 //____________________________________________________________________________
178 void AliEMCALCalibAbs::ReadTreeCalibAbsInfo(TTree *tree,
181 // how many SuperModule's worth of info do we have?
182 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
183 fNSuperModule = tree->GetEntries();
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];
195 // just to make the initializations of the arrays are done correctly, let's use memset
196 memset(relativeCalib, 0, sizeof(relativeCalib));
198 // declare the branches
199 tree->SetBranchAddress("iSM", &iSM);
200 tree->SetBranchAddress("CalibMethod", &iCalibMethod);
201 tree->SetBranchAddress("CalibPass", &iCalibPass);
202 tree->SetBranchAddress("AbsoluteCalib", &absoluteCalib);
204 tree->SetBranchAddress("RelativeCalib", relativeCalib);
206 // indices for looping over the towers
210 for (int ient=0; ient<tree->GetEntries(); ient++) {
211 tree->GetEntry(ient);
213 // assume the index SuperModules come in order: i=iSM
214 AliEMCALSuperModuleCalibAbs * t = (AliEMCALSuperModuleCalibAbs*) fSuperModuleData[iSM];
216 t->SetSuperModuleNum(iSM);
217 // first, overall values
218 t->SetCalibMethod(iCalibMethod);
219 t->SetCalibPass(iCalibPass);
220 t->SetAbsoluteCalib(absoluteCalib);
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;
227 // help variables: possibly modified or swapped indices
230 // assume that this info is already swapped and done for this basis?
232 // C side, oriented differently than A side: swap is requested
233 iColMod = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
234 iRowMod = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
237 t->SetRelativeCalib(iColMod, iRowMod, relativeCalib[iCol][iRow]);
240 } // loop over entries
245 //____________________________________________________________________________
246 void AliEMCALCalibAbs::WriteRootCalibAbsInfo(const TString &rootFileName,
249 // write data to root file. ; coordinates given on SuperModule basis
250 TFile destFile(rootFileName, "recreate");
251 if (destFile.IsZombie()) {
256 TTree *tree = new TTree("tree","");
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];
269 // just to make the initializations of the arrays are done correctly, let's use memset
270 memset(relativeCalib, 0, sizeof(relativeCalib));
272 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
273 // for looping over towers
277 // declare the branches
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) );
286 for (iSM = 0; iSM < fNSuperModule; iSM++) {
287 AliEMCALSuperModuleCalibAbs * t = (AliEMCALSuperModuleCalibAbs*) fSuperModuleData[iSM];
289 iSM = t->GetSuperModuleNum();
290 // first, overall values
291 iCalibMethod = t->GetCalibMethod();
292 iCalibPass = t->GetCalibPass();
293 absoluteCalib = t->GetAbsoluteCalib();
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;
300 // help variables: possibly modified or swapped indices
303 // assume that this info is already swapped and done for this basis?
305 // C side, oriented differently than A side: swap is requested
306 iColMod = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
307 iRowMod = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
310 relativeCalib[iColMod][iRowMod] = t->GetRelativeCalib(iCol, iRow);
322 //____________________________________________________________________________
323 AliEMCALCalibAbs::~AliEMCALCalibAbs()
325 fSuperModuleData.Delete();
328 //____________________________________________________________________________
329 AliEMCALSuperModuleCalibAbs * AliEMCALCalibAbs::GetSuperModuleCalibAbsNum(Int_t supModIndex)const
331 for (int i=0; i<fNSuperModule; i++) {
332 AliEMCALSuperModuleCalibAbs * t = (AliEMCALSuperModuleCalibAbs*) fSuperModuleData[i];
333 if (t->GetSuperModuleNum() == supModIndex) {
338 // if we arrived here, then nothing was found.. just return a NULL pointer