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
42 fSuperModuleData.SetOwner(kTRUE);
45 //____________________________________________________________________________
46 void AliEMCALBiasAPD::ReadTextBiasAPDInfo(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("AliEMCALBiasAPD::ReadBiasAPDInfo - Cannot open the APD info file %s\n", txtFileName.Data());
59 Int_t iSM = 0; // SuperModule index
66 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
68 for (Int_t i = 0; i < fNSuperModule; i++) {
69 AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[i];
72 printf("AliEMCALBiasAPD::ReadBiasAPDInfo - Error while reading input file; likely EOF..");
76 t->SetSuperModuleNum(iSM);
78 for (Int_t j=0; j<nAPDPerSM; j++) {
79 inputFile >> iCol >> iRow >> iElecId >> iDAC >> voltage;
81 // assume that this info is already swapped and done for this basis?
83 // C side, oriented differently than A side: swap is requested
84 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
85 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
88 t->SetElecId(iCol, iRow, iElecId);
89 t->SetDAC(iCol, iRow, iDAC);
90 t->SetVoltage(iCol, iRow, voltage);
100 //____________________________________________________________________________
101 void AliEMCALBiasAPD::WriteTextBiasAPDInfo(const TString &txtFileName,
104 // write data to txt file. ; coordinates given on SuperModule basis
106 std::ofstream outputFile(txtFileName.Data());
108 printf("AliEMCALBiasAPD::WriteBiasAPDInfo - Cannot open the APD output file %s\n", txtFileName.Data());
118 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
120 for (Int_t i = 0; i < fNSuperModule; i++) {
121 AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[i];
122 outputFile << t->GetSuperModuleNum() << endl;
124 for (Int_t j=0; j<nAPDPerSM; j++) {
125 iCol = j / AliEMCALGeoParams::fgkEMCALRows;
126 iRow = j % AliEMCALGeoParams::fgkEMCALRows;
128 iElecId = t->GetElecId(iCol, iRow);
129 iDAC = t->GetDAC(iCol, iRow);
130 voltage = t->GetVoltage(iCol, iRow);
133 // C side, oriented differently than A side: swap is requested
134 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
135 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
138 outputFile << iCol << " " << iRow << " "
139 << iElecId << " " << iDAC << " "
150 //____________________________________________________________________________
151 void AliEMCALBiasAPD::ReadRootBiasAPDInfo(const TString &rootFileName,
154 //Read data from root file. ; coordinates given on SuperModule basis
155 TFile inputFile(rootFileName, "read");
157 TTree *tree = (TTree*) inputFile.Get("tree");
159 ReadTreeBiasAPDInfo(tree, swapSides);
166 //____________________________________________________________________________
167 void AliEMCALBiasAPD::ReadTreeBiasAPDInfo(TTree *tree,
170 // how many SuperModule's worth of entries / APDs do we have?
171 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
172 fNSuperModule = tree->GetEntries() / nAPDPerSM;
174 Int_t iSM = 0; // SuperModule index
177 // list of values to be read
183 // declare the branches
184 tree->SetBranchAddress("iSM", &iSM);
185 tree->SetBranchAddress("iCol", &iCol);
186 tree->SetBranchAddress("iRow", &iRow);
187 tree->SetBranchAddress("iElecId", &iElecId);
188 tree->SetBranchAddress("iDAC", &iDAC);
189 tree->SetBranchAddress("voltage", &voltage);
191 for (int ient=0; ient<tree->GetEntries(); ient++) {
192 tree->GetEntry(ient);
194 // assume the index SuperModules come in order: i=iSM
195 AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[iSM];
196 t->SetSuperModuleNum(iSM);
198 // assume that this info is already swapped and done for this basis?
200 // C side, oriented differently than A side: swap is requested
201 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
202 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
205 t->SetElecId(iCol, iRow, iElecId);
206 t->SetDAC(iCol, iRow, iDAC);
207 t->SetVoltage(iCol, iRow, voltage);
214 //____________________________________________________________________________
215 void AliEMCALBiasAPD::WriteRootBiasAPDInfo(const TString &rootFileName,
218 // write data to root file. ; coordinates given on SuperModule basis
219 TFile destFile(rootFileName, "recreate");
220 if (destFile.IsZombie()) {
225 TTree *tree = new TTree("tree","");
227 // variables for filling the TTree
228 Int_t iSM = 0; // SuperModule index
234 // declare the branches
235 tree->Branch("iSM", &iSM, "iSM/I");
236 tree->Branch("iCol", &iCol, "iCol/I");
237 tree->Branch("iRow", &iRow, "iRow/I");
238 tree->Branch("iElecId", &iElecId, "iElecId/I");
239 tree->Branch("iDAC", &iDAC, "iDAC/I");
240 tree->Branch("voltage", &voltage, "voltage/F");
242 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
244 for (iSM = 0; iSM < fNSuperModule; iSM++) {
245 AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[iSM];
247 for (Int_t j=0; j<nAPDPerSM; j++) {
248 iCol = j / AliEMCALGeoParams::fgkEMCALRows;
249 iRow = j % AliEMCALGeoParams::fgkEMCALRows;
251 iElecId = t->GetElecId(iCol, iRow);
252 iDAC = t->GetDAC(iCol, iRow);
253 voltage = t->GetVoltage(iCol, iRow);
256 // C side, oriented differently than A side: swap is requested
257 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
258 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
272 //____________________________________________________________________________
273 AliEMCALBiasAPD::~AliEMCALBiasAPD()
275 fSuperModuleData.Delete();
278 //____________________________________________________________________________
279 AliEMCALSuperModuleBiasAPD * AliEMCALBiasAPD::GetSuperModuleBiasAPDNum(Int_t supModIndex)const
281 for (int i=0; i<fNSuperModule; i++) {
282 AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[i];
283 if (t->GetSuperModuleNum() == supModIndex) {
288 // if we arrived here, then nothing was found.. just return a NULL pointer