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..\n");
76 t->SetSuperModuleNum(iSM);
78 for (Int_t j=0; j<nAPDPerSM; j++) {
79 inputFile >> iCol >> iRow >> iElecId >> iDAC >> voltage;
81 // check that input values are not out bounds
82 if (iCol<0 || iCol>(AliEMCALGeoParams::fgkEMCALCols-1) ||
83 iRow<0 || iRow>(AliEMCALGeoParams::fgkEMCALRows-1) ) {
84 printf("AliEMCALBiasAPD::ReadBiasAPDInfo - Error while reading input file; j %d iCol %d iRow %d\n", j, iCol, iRow);
88 // assume that this info is already swapped and done for this basis?
90 // C side, oriented differently than A side: swap is requested
91 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
92 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
95 t->SetElecId(iCol, iRow, iElecId);
96 t->SetDAC(iCol, iRow, iDAC);
97 t->SetVoltage(iCol, iRow, voltage);
107 //____________________________________________________________________________
108 void AliEMCALBiasAPD::WriteTextBiasAPDInfo(const TString &txtFileName,
111 // write data to txt file. ; coordinates given on SuperModule basis
113 std::ofstream outputFile(txtFileName.Data());
115 printf("AliEMCALBiasAPD::WriteBiasAPDInfo - Cannot open the APD output file %s\n", txtFileName.Data());
125 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
127 for (Int_t i = 0; i < fNSuperModule; i++) {
128 AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[i];
129 outputFile << t->GetSuperModuleNum() << endl;
131 for (Int_t j=0; j<nAPDPerSM; j++) {
132 iCol = j / AliEMCALGeoParams::fgkEMCALRows;
133 iRow = j % AliEMCALGeoParams::fgkEMCALRows;
135 iElecId = t->GetElecId(iCol, iRow);
136 iDAC = t->GetDAC(iCol, iRow);
137 voltage = t->GetVoltage(iCol, iRow);
140 // C side, oriented differently than A side: swap is requested
141 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
142 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
145 outputFile << iCol << " " << iRow << " "
146 << iElecId << " " << iDAC << " "
157 //____________________________________________________________________________
158 void AliEMCALBiasAPD::ReadRootBiasAPDInfo(const TString &rootFileName,
161 //Read data from root file. ; coordinates given on SuperModule basis
162 TFile inputFile(rootFileName, "read");
164 TTree *tree = (TTree*) inputFile.Get("tree");
166 ReadTreeBiasAPDInfo(tree, swapSides);
173 //____________________________________________________________________________
174 void AliEMCALBiasAPD::ReadTreeBiasAPDInfo(TTree *tree,
177 // how many SuperModule's worth of entries / APDs do we have?
178 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
179 fNSuperModule = tree->GetEntries() / nAPDPerSM;
181 Int_t iSM = 0; // SuperModule index
184 // list of values to be read
190 // declare the branches
191 tree->SetBranchAddress("iSM", &iSM);
192 tree->SetBranchAddress("iCol", &iCol);
193 tree->SetBranchAddress("iRow", &iRow);
194 tree->SetBranchAddress("iElecId", &iElecId);
195 tree->SetBranchAddress("iDAC", &iDAC);
196 tree->SetBranchAddress("voltage", &voltage);
198 for (int ient=0; ient<tree->GetEntries(); ient++) {
199 tree->GetEntry(ient);
201 // assume the index SuperModules come in order: i=iSM
202 AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[iSM];
203 t->SetSuperModuleNum(iSM);
205 // assume that this info is already swapped and done for this basis?
207 // C side, oriented differently than A side: swap is requested
208 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
209 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
212 t->SetElecId(iCol, iRow, iElecId);
213 t->SetDAC(iCol, iRow, iDAC);
214 t->SetVoltage(iCol, iRow, voltage);
221 //____________________________________________________________________________
222 void AliEMCALBiasAPD::WriteRootBiasAPDInfo(const TString &rootFileName,
225 // write data to root file. ; coordinates given on SuperModule basis
226 TFile destFile(rootFileName, "recreate");
227 if (destFile.IsZombie()) {
232 TTree *tree = new TTree("tree","");
234 // variables for filling the TTree
235 Int_t iSM = 0; // SuperModule index
241 // declare the branches
242 tree->Branch("iSM", &iSM, "iSM/I");
243 tree->Branch("iCol", &iCol, "iCol/I");
244 tree->Branch("iRow", &iRow, "iRow/I");
245 tree->Branch("iElecId", &iElecId, "iElecId/I");
246 tree->Branch("iDAC", &iDAC, "iDAC/I");
247 tree->Branch("voltage", &voltage, "voltage/F");
249 Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
251 for (iSM = 0; iSM < fNSuperModule; iSM++) {
252 AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[iSM];
254 for (Int_t j=0; j<nAPDPerSM; j++) {
255 iCol = j / AliEMCALGeoParams::fgkEMCALRows;
256 iRow = j % AliEMCALGeoParams::fgkEMCALRows;
258 iElecId = t->GetElecId(iCol, iRow);
259 iDAC = t->GetDAC(iCol, iRow);
260 voltage = t->GetVoltage(iCol, iRow);
263 // C side, oriented differently than A side: swap is requested
264 iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
265 iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
279 //____________________________________________________________________________
280 AliEMCALBiasAPD::~AliEMCALBiasAPD()
282 fSuperModuleData.Delete();
285 //____________________________________________________________________________
286 AliEMCALSuperModuleBiasAPD * AliEMCALBiasAPD::GetSuperModuleBiasAPDNum(Int_t supModIndex)const
287 { // getter via index
288 for (int i=0; i<fNSuperModule; i++) {
289 AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[i];
290 if (t->GetSuperModuleNum() == supModIndex) {
295 // if we arrived here, then nothing was found.. just return a NULL pointer