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 info on APD bias settings/voltages
25 #include "AliEMCALCalibTimeDepCorrection.h"
29 ClassImp(AliEMCALCalibTimeDepCorrection)
31 //____________________________________________________________________________
32 AliEMCALCalibTimeDepCorrection::AliEMCALCalibTimeDepCorrection() :
36 //Default constructor.
39 //____________________________________________________________________________
40 void AliEMCALCalibTimeDepCorrection::InitCorrection(Int_t nSM, Int_t nBins, Float_t val=1.0)
42 // This methods assumes that you are using SuperModules 0..nSM-1
44 if (fSuperModuleData) delete [] fSuperModuleData;
45 fSuperModuleData = new AliEMCALSuperModuleCalibTimeDepCorrection[fNSuperModule];
47 Int_t iSM = 0; // SuperModule index
51 Float_t Correction = val;
53 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
55 for (Int_t i = 0; i < fNSuperModule; i++) {
56 AliEMCALSuperModuleCalibTimeDepCorrection &t = fSuperModuleData[i];
57 t.fSuperModuleNum = iSM; // assume SMs are coming in order
59 for (Int_t j=0; j<nAPDPerSM; j++) {
62 t.fCorrection[iCol][iRow].Set(nCorr);
63 for (Int_t k=0; k<nCorr; k++) {
65 t.fCorrection[iCol][iRow].AddAt(Correction, k); // AddAt = SetAt..
74 //____________________________________________________________________________
75 void AliEMCALCalibTimeDepCorrection::SetCorrection(Int_t supModIndex, Int_t iCol, Int_t iRow, Int_t iBin, Float_t val)
76 { // if you call for non-existing data, there may be a crash..
77 fSuperModuleData[supModIndex].fCorrection[iCol][iRow].AddAt(val, iBin); // AddAt = SetAt..
81 //____________________________________________________________________________
82 Float_t AliEMCALCalibTimeDepCorrection::GetCorrection(Int_t supModIndex, Int_t iCol, Int_t iRow, Int_t iBin)const
83 { // if you call for non-existing data, there may be a crash..
84 return fSuperModuleData[supModIndex].fCorrection[iCol][iRow].At(iBin);
87 //____________________________________________________________________________
88 void AliEMCALCalibTimeDepCorrection::ReadCalibTimeDepCorrectionInfo(Int_t nSM, const TString &txtFileName,
91 //Read data from txt file. ; coordinates given on SuperModule basis
93 std::ifstream inputFile(txtFileName.Data());
95 printf("AliEMCALCalibTimeDepCorrection::ReadCalibTimeDepCorrectionInfo - Cannot open the APD info file %s\n", txtFileName.Data());
100 if (fSuperModuleData) delete [] fSuperModuleData;
101 fSuperModuleData = new AliEMCALSuperModuleCalibTimeDepCorrection[fNSuperModule];
103 Int_t iSM = 0; // SuperModule index
107 Float_t Correction = 0;
109 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
111 for (Int_t i = 0; i < fNSuperModule; i++) {
112 AliEMCALSuperModuleCalibTimeDepCorrection &t = fSuperModuleData[i];
114 printf("AliEMCALCalibTimeDepCorrection::ReadCalibTimeDepCorrectionInfo - Error while reading input file; likely EOF..");
118 t.fSuperModuleNum = iSM;
120 for (Int_t j=0; j<nAPDPerSM; j++) {
121 inputFile >> iCol >> iRow >> nCorr;
123 // assume that this info is already swapped and done for this basis?
125 // C side, oriented differently than A side: swap is requested
126 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
127 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
130 // set size of TArray
131 t.fCorrection[iCol][iRow].Set(nCorr);
132 for (Int_t k=0; k<nCorr; k++) {
133 inputFile >> Correction;
135 t.fCorrection[iCol][iRow].AddAt(Correction, k);
146 //____________________________________________________________________________
147 void AliEMCALCalibTimeDepCorrection::WriteCalibTimeDepCorrectionInfo(const TString &txtFileName,
150 // write data to txt file. ; coordinates given on SuperModule basis
152 std::ofstream outputFile(txtFileName.Data());
154 printf("AliEMCALCalibTimeDepCorrection::WriteCalibTimeDepCorrectionInfo - Cannot open the APD output file %s\n", txtFileName.Data());
162 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
164 for (Int_t i = 0; i < fNSuperModule; i++) {
165 AliEMCALSuperModuleCalibTimeDepCorrection &t = fSuperModuleData[i];
166 outputFile << t.fSuperModuleNum << endl;
168 for (Int_t j=0; j<nAPDPerSM; j++) {
169 iCol = j / AliEMCALGeoParams::fgkEMCALRows;
170 iRow = j % AliEMCALGeoParams::fgkEMCALRows;
172 nCorr = t.fCorrection[iCol][iRow].GetSize();
175 // C side, oriented differently than A side: swap is requested
176 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
177 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
180 outputFile << iCol << " " << iRow << " " << nCorr << endl;
181 for (Int_t k=0; k<nCorr; k++) {
182 outputFile << t.fCorrection[iCol][iRow].At(k) << " ";
195 //____________________________________________________________________________
196 AliEMCALCalibTimeDepCorrection::~AliEMCALCalibTimeDepCorrection()
198 delete [] fSuperModuleData;
201 //____________________________________________________________________________
202 AliEMCALCalibTimeDepCorrection::AliEMCALSuperModuleCalibTimeDepCorrection AliEMCALCalibTimeDepCorrection::GetSuperModuleCalibTimeDepCorrectionId(Int_t supModIndex)const
204 AliEMCALSuperModuleCalibTimeDepCorrection t; // just to maybe prevent a crash, but we are returning something not-initialized so maybe not better really..
205 if (!fSuperModuleData)
208 return fSuperModuleData[supModIndex];
211 //____________________________________________________________________________
212 AliEMCALCalibTimeDepCorrection::AliEMCALSuperModuleCalibTimeDepCorrection AliEMCALCalibTimeDepCorrection::GetSuperModuleCalibTimeDepCorrectionNum(Int_t supModIndex)const
214 AliEMCALSuperModuleCalibTimeDepCorrection t; // just to maybe prevent a crash, but we are returning something not-initialized so maybe not better really..
215 if (!fSuperModuleData)
218 for (int i=0; i<fNSuperModule; i++) {
219 if (fSuperModuleData[i].fSuperModuleNum == supModIndex) {
220 return fSuperModuleData[i];