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
26 #include "AliEMCALBiasAPD.h"
30 ClassImp(AliEMCALBiasAPD)
32 //____________________________________________________________________________
33 AliEMCALBiasAPD::AliEMCALBiasAPD(const int nSM) :
34 fNSuperModule(nSM), // make space for everyone
37 //Default constructor.
38 for (int i=0; i<fNSuperModule; i++) {
39 fSuperModuleData.Add(new AliEMCALSuperModuleBiasAPD(i));
41 fSuperModuleData.Compress(); // compress the TObjArray
44 //____________________________________________________________________________
45 void AliEMCALBiasAPD::ReadTextBiasAPDInfo(Int_t nSM, const TString &txtFileName,
48 //Read data from txt file. ; coordinates given on SuperModule basis
50 std::ifstream inputFile(txtFileName.Data());
52 printf("AliEMCALBiasAPD::ReadBiasAPDInfo - Cannot open the APD info file %s\n", txtFileName.Data());
58 Int_t iSM = 0; // SuperModule index
65 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
67 for (Int_t i = 0; i < fNSuperModule; i++) {
68 AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[i];
71 printf("AliEMCALBiasAPD::ReadBiasAPDInfo - Error while reading input file; likely EOF..");
75 t->SetSuperModuleNum(iSM);
77 for (Int_t j=0; j<nAPDPerSM; j++) {
78 inputFile >> iCol >> iRow >> iElecId >> iDAC >> voltage;
80 // assume that this info is already swapped and done for this basis?
82 // C side, oriented differently than A side: swap is requested
83 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
84 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
87 t->SetElecId(iCol, iRow, iElecId);
88 t->SetDAC(iCol, iRow, iDAC);
89 t->SetVoltage(iCol, iRow, voltage);
99 //____________________________________________________________________________
100 void AliEMCALBiasAPD::WriteTextBiasAPDInfo(const TString &txtFileName,
103 // write data to txt file. ; coordinates given on SuperModule basis
105 std::ofstream outputFile(txtFileName.Data());
107 printf("AliEMCALBiasAPD::WriteBiasAPDInfo - Cannot open the APD output file %s\n", txtFileName.Data());
117 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
119 for (Int_t i = 0; i < fNSuperModule; i++) {
120 AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[i];
121 outputFile << t->GetSuperModuleNum() << endl;
123 for (Int_t j=0; j<nAPDPerSM; j++) {
124 iCol = j / AliEMCALGeoParams::fgkEMCALRows;
125 iRow = j % AliEMCALGeoParams::fgkEMCALRows;
127 iElecId = t->GetElecId(iCol, iRow);
128 iDAC = t->GetDAC(iCol, iRow);
129 voltage = t->GetVoltage(iCol, iRow);
132 // C side, oriented differently than A side: swap is requested
133 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
134 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
137 outputFile << iCol << " " << iRow << " "
138 << iElecId << " " << iDAC << " "
149 //____________________________________________________________________________
150 void AliEMCALBiasAPD::ReadRootBiasAPDInfo(const TString &rootFileName,
153 //Read data from root file. ; coordinates given on SuperModule basis
154 TFile inputFile(rootFileName, "read");
156 TTree *tree = (TTree*) inputFile.Get("tree");
158 ReadTreeBiasAPDInfo(tree, swapSides);
165 //____________________________________________________________________________
166 void AliEMCALBiasAPD::ReadTreeBiasAPDInfo(TTree *tree,
169 // how many SuperModule's worth of entries / APDs do we have?
170 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
171 fNSuperModule = tree->GetEntries() / nAPDPerSM;
173 Int_t iSM = 0; // SuperModule index
176 // list of values to be read
182 // declare the branches
183 tree->SetBranchAddress("iSM", &iSM);
184 tree->SetBranchAddress("iCol", &iCol);
185 tree->SetBranchAddress("iRow", &iRow);
186 tree->SetBranchAddress("iElecId", &iElecId);
187 tree->SetBranchAddress("iDAC", &iDAC);
188 tree->SetBranchAddress("voltage", &voltage);
190 for (int ient=0; ient<tree->GetEntries(); ient++) {
191 tree->GetEntry(ient);
193 // assume the index SuperModules come in order: i=iSM
194 AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[iSM];
195 t->SetSuperModuleNum(iSM);
197 // assume that this info is already swapped and done for this basis?
199 // C side, oriented differently than A side: swap is requested
200 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
201 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
204 t->SetElecId(iCol, iRow, iElecId);
205 t->SetDAC(iCol, iRow, iDAC);
206 t->SetVoltage(iCol, iRow, voltage);
213 //____________________________________________________________________________
214 void AliEMCALBiasAPD::WriteRootBiasAPDInfo(const TString &rootFileName,
217 // write data to root file. ; coordinates given on SuperModule basis
218 TFile destFile(rootFileName, "recreate");
219 if (destFile.IsZombie()) {
224 TTree *tree = new TTree("tree","");
226 // variables for filling the TTree
227 Int_t iSM = 0; // SuperModule index
233 // declare the branches
234 tree->Branch("iSM", &iSM, "iSM/I");
235 tree->Branch("iCol", &iCol, "iCol/I");
236 tree->Branch("iRow", &iRow, "iRow/I");
237 tree->Branch("iElecId", &iElecId, "iElecId/I");
238 tree->Branch("iDAC", &iDAC, "iDAC/I");
239 tree->Branch("voltage", &voltage, "voltage/F");
241 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
243 for (iSM = 0; iSM < fNSuperModule; iSM++) {
244 AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[iSM];
246 for (Int_t j=0; j<nAPDPerSM; j++) {
247 iCol = j / AliEMCALGeoParams::fgkEMCALRows;
248 iRow = j % AliEMCALGeoParams::fgkEMCALRows;
250 iElecId = t->GetElecId(iCol, iRow);
251 iDAC = t->GetDAC(iCol, iRow);
252 voltage = t->GetVoltage(iCol, iRow);
255 // C side, oriented differently than A side: swap is requested
256 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
257 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
271 //____________________________________________________________________________
272 AliEMCALBiasAPD::~AliEMCALBiasAPD()
274 fSuperModuleData.Delete();
277 //____________________________________________________________________________
278 AliEMCALSuperModuleBiasAPD * AliEMCALBiasAPD::GetSuperModuleBiasAPDNum(Int_t supModIndex)const
280 for (int i=0; i<fNSuperModule; i++) {
281 AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[i];
282 if (t->GetSuperModuleNum() == supModIndex) {
287 // if we arrived here, then nothing was found.. just return a NULL pointer