]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALBiasAPD.cxx
New CE DA ... (Jens)
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALBiasAPD.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id: $ */
17
18 // Objects of this class contain info on APD bias settings/voltages
19 //
20
21 #include <fstream>
22 #include <TString.h>
23 #include <TFile.h>
24 #include <TTree.h>
25
26 #include "AliEMCALBiasAPD.h"
27
28 using namespace std;
29
30 ClassImp(AliEMCALBiasAPD)
31
32 //____________________________________________________________________________
33 AliEMCALBiasAPD::AliEMCALBiasAPD(const int nSM) : 
34   fNSuperModule(nSM), // make space for everyone 
35   fSuperModuleData()
36 {
37   //Default constructor.
38   for (int i=0; i<fNSuperModule; i++) {
39     fSuperModuleData.Add(new AliEMCALSuperModuleBiasAPD(i));
40   }
41   fSuperModuleData.Compress(); // compress the TObjArray
42   fSuperModuleData.SetOwner(kTRUE); 
43 }
44
45 //____________________________________________________________________________
46 void AliEMCALBiasAPD::ReadTextBiasAPDInfo(Int_t nSM, const TString &txtFileName,
47                                           Bool_t swapSides)
48 {
49   //Read data from txt file. ; coordinates given on SuperModule basis
50
51   std::ifstream inputFile(txtFileName.Data());
52   if (!inputFile) {
53     printf("AliEMCALBiasAPD::ReadBiasAPDInfo - Cannot open the APD info file %s\n", txtFileName.Data());
54     return;
55   }
56
57   fNSuperModule = nSM;
58
59   Int_t iSM = 0; // SuperModule index
60   Int_t iCol = 0;
61   Int_t iRow = 0;
62   Int_t iElecId = 0;
63   Int_t iDAC = 0;
64   Float_t voltage = 0;
65
66   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
67
68   for (Int_t i = 0; i < fNSuperModule; i++) {
69     AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[i];
70
71     if (!inputFile) {
72       printf("AliEMCALBiasAPD::ReadBiasAPDInfo - Error while reading input file; likely EOF..\n");
73       return;
74     }
75     inputFile >> iSM;
76     t->SetSuperModuleNum(iSM);
77
78     for (Int_t j=0; j<nAPDPerSM; j++) {
79       inputFile >> iCol >> iRow >> iElecId >> iDAC >> voltage;
80
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);
85       return;
86       }
87
88       // assume that this info is already swapped and done for this basis?
89       if (swapSides) {
90         // C side, oriented differently than A side: swap is requested
91         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
92         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
93       }
94
95       t->SetElecId(iCol, iRow, iElecId);
96       t->SetDAC(iCol, iRow, iDAC);
97       t->SetVoltage(iCol, iRow, voltage);
98     }
99
100   } // i, SuperModule
101
102   inputFile.close();
103
104   return;
105 }
106
107 //____________________________________________________________________________
108 void AliEMCALBiasAPD::WriteTextBiasAPDInfo(const TString &txtFileName,
109                                            Bool_t swapSides)
110 {
111   // write data to txt file. ; coordinates given on SuperModule basis
112
113   std::ofstream outputFile(txtFileName.Data());
114   if (!outputFile) {
115     printf("AliEMCALBiasAPD::WriteBiasAPDInfo - Cannot open the APD output file %s\n", txtFileName.Data());
116     return;
117   }
118
119   Int_t iCol = 0;
120   Int_t iRow = 0;
121   Int_t iElecId = 0;
122   Int_t iDAC = 0;
123   Float_t voltage = 0;
124
125   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
126
127   for (Int_t i = 0; i < fNSuperModule; i++) {
128     AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[i];
129     outputFile << t->GetSuperModuleNum() << endl;
130
131     for (Int_t j=0; j<nAPDPerSM; j++) {
132       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
133       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
134
135       iElecId = t->GetElecId(iCol, iRow);
136       iDAC = t->GetDAC(iCol, iRow);
137       voltage = t->GetVoltage(iCol, iRow);
138
139       if (swapSides) {
140         // C side, oriented differently than A side: swap is requested
141         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
142         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
143       }
144
145       outputFile << iCol << " " << iRow << " " 
146                  << iElecId << " " << iDAC << " "
147                  << voltage << endl;
148     }
149
150   } // i, SuperModule
151
152   outputFile.close();
153
154   return;
155 }
156
157 //____________________________________________________________________________
158 void AliEMCALBiasAPD::ReadRootBiasAPDInfo(const TString &rootFileName,
159                                           Bool_t swapSides)
160 {
161   //Read data from root file. ; coordinates given on SuperModule basis
162   TFile inputFile(rootFileName, "read");  
163
164   TTree *tree = (TTree*) inputFile.Get("tree");
165
166   ReadTreeBiasAPDInfo(tree, swapSides);
167
168   inputFile.Close();
169
170   return;
171 }
172
173 //____________________________________________________________________________
174 void AliEMCALBiasAPD::ReadTreeBiasAPDInfo(TTree *tree,
175                                           Bool_t swapSides)
176 {
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;
180
181   Int_t iSM = 0; // SuperModule index
182   Int_t iCol = 0;
183   Int_t iRow = 0;
184   // list of values to be read
185   Int_t iElecId = 0;
186   Int_t iDAC = 0;
187   Float_t voltage = 0;     
188   // end - all values
189
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);
197
198   for (int ient=0; ient<tree->GetEntries(); ient++) {
199     tree->GetEntry(ient);
200
201     // assume the index SuperModules come in order: i=iSM
202     AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[iSM];
203     t->SetSuperModuleNum(iSM);
204
205     // assume that this info is already swapped and done for this basis?
206     if (swapSides) {
207       // C side, oriented differently than A side: swap is requested
208       iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
209       iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
210     }
211
212     t->SetElecId(iCol, iRow, iElecId);
213     t->SetDAC(iCol, iRow, iDAC);
214     t->SetVoltage(iCol, iRow, voltage);
215
216   } // 
217
218   return;
219 }
220
221 //____________________________________________________________________________
222 void AliEMCALBiasAPD::WriteRootBiasAPDInfo(const TString &rootFileName,
223                                            Bool_t swapSides)
224 {
225   // write data to root file. ; coordinates given on SuperModule basis
226   TFile destFile(rootFileName, "recreate");  
227   if (destFile.IsZombie()) {
228     return;
229   }  
230   destFile.cd();
231
232   TTree *tree = new TTree("tree","");
233
234   // variables for filling the TTree
235   Int_t iSM = 0; // SuperModule index
236   Int_t iCol = 0;
237   Int_t iRow = 0;
238   Int_t iElecId = 0;
239   Int_t iDAC = 0;
240   Float_t voltage = 0;
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");
248
249   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
250
251   for (iSM = 0; iSM < fNSuperModule; iSM++) {
252     AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[iSM];
253
254     for (Int_t j=0; j<nAPDPerSM; j++) {
255       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
256       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
257
258       iElecId = t->GetElecId(iCol, iRow);
259       iDAC = t->GetDAC(iCol, iRow);
260       voltage = t->GetVoltage(iCol, iRow);
261
262       if (swapSides) {
263         // C side, oriented differently than A side: swap is requested
264         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
265         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
266       }
267
268       tree->Fill();
269     }
270
271   } // i, SuperModule
272
273   tree->Write();
274   destFile.Close();
275
276   return;
277 }
278
279 //____________________________________________________________________________
280 AliEMCALBiasAPD::~AliEMCALBiasAPD()
281 {
282   fSuperModuleData.Delete();
283 }
284
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) {
291       return t;
292     }
293   }
294
295   // if we arrived here, then nothing was found.. just return a NULL pointer 
296   return NULL;
297 }
298