]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALBiasAPD.cxx
Update on trigger code (Rachid)
[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..");
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       // assume that this info is already swapped and done for this basis?
82       if (swapSides) {
83         // C side, oriented differently than A side: swap is requested
84         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
85         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
86       }
87
88       t->SetElecId(iCol, iRow, iElecId);
89       t->SetDAC(iCol, iRow, iDAC);
90       t->SetVoltage(iCol, iRow, voltage);
91     }
92
93   } // i, SuperModule
94
95   inputFile.close();
96
97   return;
98 }
99
100 //____________________________________________________________________________
101 void AliEMCALBiasAPD::WriteTextBiasAPDInfo(const TString &txtFileName,
102                                            Bool_t swapSides)
103 {
104   // write data to txt file. ; coordinates given on SuperModule basis
105
106   std::ofstream outputFile(txtFileName.Data());
107   if (!outputFile) {
108     printf("AliEMCALBiasAPD::WriteBiasAPDInfo - Cannot open the APD output file %s\n", txtFileName.Data());
109     return;
110   }
111
112   Int_t iCol = 0;
113   Int_t iRow = 0;
114   Int_t iElecId = 0;
115   Int_t iDAC = 0;
116   Float_t voltage = 0;
117
118   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
119
120   for (Int_t i = 0; i < fNSuperModule; i++) {
121     AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[i];
122     outputFile << t->GetSuperModuleNum() << endl;
123
124     for (Int_t j=0; j<nAPDPerSM; j++) {
125       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
126       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
127
128       iElecId = t->GetElecId(iCol, iRow);
129       iDAC = t->GetDAC(iCol, iRow);
130       voltage = t->GetVoltage(iCol, iRow);
131
132       if (swapSides) {
133         // C side, oriented differently than A side: swap is requested
134         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
135         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
136       }
137
138       outputFile << iCol << " " << iRow << " " 
139                  << iElecId << " " << iDAC << " "
140                  << voltage << endl;
141     }
142
143   } // i, SuperModule
144
145   outputFile.close();
146
147   return;
148 }
149
150 //____________________________________________________________________________
151 void AliEMCALBiasAPD::ReadRootBiasAPDInfo(const TString &rootFileName,
152                                           Bool_t swapSides)
153 {
154   //Read data from root file. ; coordinates given on SuperModule basis
155   TFile inputFile(rootFileName, "read");  
156
157   TTree *tree = (TTree*) inputFile.Get("tree");
158
159   ReadTreeBiasAPDInfo(tree, swapSides);
160
161   inputFile.Close();
162
163   return;
164 }
165
166 //____________________________________________________________________________
167 void AliEMCALBiasAPD::ReadTreeBiasAPDInfo(TTree *tree,
168                                           Bool_t swapSides)
169 {
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;
173
174   Int_t iSM = 0; // SuperModule index
175   Int_t iCol = 0;
176   Int_t iRow = 0;
177   // list of values to be read
178   Int_t iElecId = 0;
179   Int_t iDAC = 0;
180   Float_t voltage = 0;     
181   // end - all values
182
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);
190
191   for (int ient=0; ient<tree->GetEntries(); ient++) {
192     tree->GetEntry(ient);
193
194     // assume the index SuperModules come in order: i=iSM
195     AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[iSM];
196     t->SetSuperModuleNum(iSM);
197
198     // assume that this info is already swapped and done for this basis?
199     if (swapSides) {
200       // C side, oriented differently than A side: swap is requested
201       iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
202       iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
203     }
204
205     t->SetElecId(iCol, iRow, iElecId);
206     t->SetDAC(iCol, iRow, iDAC);
207     t->SetVoltage(iCol, iRow, voltage);
208
209   } // 
210
211   return;
212 }
213
214 //____________________________________________________________________________
215 void AliEMCALBiasAPD::WriteRootBiasAPDInfo(const TString &rootFileName,
216                                            Bool_t swapSides)
217 {
218   // write data to root file. ; coordinates given on SuperModule basis
219   TFile destFile(rootFileName, "recreate");  
220   if (destFile.IsZombie()) {
221     return;
222   }  
223   destFile.cd();
224
225   TTree *tree = new TTree("tree","");
226
227   // variables for filling the TTree
228   Int_t iSM = 0; // SuperModule index
229   Int_t iCol = 0;
230   Int_t iRow = 0;
231   Int_t iElecId = 0;
232   Int_t iDAC = 0;
233   Float_t voltage = 0;
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");
241
242   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
243
244   for (iSM = 0; iSM < fNSuperModule; iSM++) {
245     AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[iSM];
246
247     for (Int_t j=0; j<nAPDPerSM; j++) {
248       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
249       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
250
251       iElecId = t->GetElecId(iCol, iRow);
252       iDAC = t->GetDAC(iCol, iRow);
253       voltage = t->GetVoltage(iCol, iRow);
254
255       if (swapSides) {
256         // C side, oriented differently than A side: swap is requested
257         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
258         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
259       }
260
261       tree->Fill();
262     }
263
264   } // i, SuperModule
265
266   tree->Write();
267   destFile.Close();
268
269   return;
270 }
271
272 //____________________________________________________________________________
273 AliEMCALBiasAPD::~AliEMCALBiasAPD()
274 {
275   fSuperModuleData.Delete();
276 }
277
278 //____________________________________________________________________________
279 AliEMCALSuperModuleBiasAPD * AliEMCALBiasAPD::GetSuperModuleBiasAPDNum(Int_t supModIndex)const
280 {
281   for (int i=0; i<fNSuperModule; i++) {
282     AliEMCALSuperModuleBiasAPD * t = (AliEMCALSuperModuleBiasAPD*) fSuperModuleData[i];
283     if (t->GetSuperModuleNum() == supModIndex) {
284       return t;
285     }
286   }
287
288   // if we arrived here, then nothing was found.. just return a NULL pointer 
289   return NULL;
290 }
291