]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALCalibReference.cxx
Correction of SA track rejection
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALCalibReference.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 basis for reference calibrations
19 //
20
21 #include <fstream>
22 #include <TString.h>
23 #include <TFile.h>
24 #include <TTree.h>
25
26 #include "AliEMCALCalibReference.h"
27
28 using namespace std;
29
30 ClassImp(AliEMCALCalibReference)
31
32 //____________________________________________________________________________
33 AliEMCALCalibReference::AliEMCALCalibReference(const int nSM) : 
34   fNSuperModule(nSM),
35   fSuperModuleData()
36 {
37   //Default constructor.
38   for (int i=0; i<fNSuperModule; i++) {
39     fSuperModuleData.Add(new AliEMCALSuperModuleCalibReference(i));
40   }
41   fSuperModuleData.Compress(); // compress the TObjArray
42 }
43
44 //____________________________________________________________________________
45 void AliEMCALCalibReference::ReadTextCalibReferenceInfo(Int_t nSM, const TString &txtFileName,
46                                             Bool_t swapSides)
47 {
48   //Read data from txt file. ; coordinates given on SuperModule basis
49
50   std::ifstream inputFile(txtFileName.Data());
51   if (!inputFile) {
52     printf("AliEMCALCalibReference::ReadCalibReferenceInfo - Cannot open the APD info file %s\n", txtFileName.Data());
53     return;
54   }
55
56   fNSuperModule = nSM;
57
58   Int_t iSM = 0; // SuperModule index
59   Int_t iCol = 0;
60   Int_t iRow = 0;
61   Int_t id = 0;
62
63   // list of values to be read
64   // first: overall values for the whole SuperModule
65   Int_t iReferenceTime; 
66   // second: additional info for LED Reference and SM temperature
67   Float_t rLEDRefAmp;
68   Float_t rLEDRefAmpRMS;
69   Int_t iLEDRefHighLow;
70   Float_t temperature;
71   Float_t temperatureRMS;
72   // third: info for each tower
73   Int_t iHighLow; // 
74   Float_t rLEDAmp; // low gain eq. amplitude
75   Float_t rLEDAmpRMS; //
76   // end - all values
77
78   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
79
80   for (Int_t i = 0; i < fNSuperModule; i++) {
81     AliEMCALSuperModuleCalibReference * t = (AliEMCALSuperModuleCalibReference*) fSuperModuleData[i];
82     if (!inputFile) {
83       printf("AliEMCALCalibReference::ReadCalibReferenceInfo - Error while reading input file; likely EOF..");
84       return;
85     }
86     inputFile >> iSM;
87     t->SetSuperModuleNum(iSM);
88
89     // first: overall values for the whole SuperModule
90     inputFile >> iReferenceTime;
91     t->SetReferenceTime(iReferenceTime);
92
93     // second: additional info for LED Reference and SM temperature
94     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALLEDRefs; j++) {
95       inputFile >> id >> iLEDRefHighLow >> rLEDRefAmp >> rLEDRefAmpRMS;
96       t->SetLEDRefHighLow(id, iLEDRefHighLow);
97       t->SetLEDRefAmp(id, rLEDRefAmp);
98       t->SetLEDRefAmpRMS(id, rLEDRefAmpRMS);
99     }
100     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALTempSensors; j++) {
101       inputFile >> id >> temperature >> temperatureRMS;
102       t->SetTemperature(id, temperature);
103       t->SetTemperatureRMS(id, temperatureRMS);
104     }
105
106     // third: info for each tower
107     for (Int_t j=0; j<nAPDPerSM; j++) {
108       inputFile >> iCol >> iRow 
109                 >> iHighLow >> rLEDAmp >> rLEDAmpRMS;
110
111       // assume that this info is already swapped and done for this basis?
112       if (swapSides) {
113         // C side, oriented differently than A side: swap is requested
114         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
115         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
116       }
117
118       AliEMCALCalibReferenceVal * v = t->GetAPDVal(iCol, iRow);
119
120       v->SetHighLow(iHighLow);
121       v->SetLEDAmp(rLEDAmp);
122       v->SetLEDAmpRMS(rLEDAmpRMS);
123     }
124
125   } // i, SuperModule
126
127   inputFile.close();
128
129   return;
130 }
131
132 //____________________________________________________________________________
133 void AliEMCALCalibReference::WriteTextCalibReferenceInfo(const TString &txtFileName,
134                                              Bool_t swapSides)
135 {
136   // write data to txt file. ; coordinates given on SuperModule basis
137
138   std::ofstream outputFile(txtFileName.Data());
139   if (!outputFile) {
140     printf("AliEMCALCalibReference::WriteCalibReferenceInfo - Cannot open the APD output file %s\n", txtFileName.Data());
141     return;
142   }
143
144   Int_t iCol = 0;
145   Int_t iRow = 0;
146
147   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
148
149   for (Int_t i = 0; i < fNSuperModule; i++) {
150     AliEMCALSuperModuleCalibReference * t = (AliEMCALSuperModuleCalibReference*) fSuperModuleData[i];
151
152     // first: overall values for the whole SuperModule
153     outputFile << t->GetSuperModuleNum() << endl;
154     outputFile << t->GetReferenceTime() << endl;
155
156     // second: additional info for LED Reference and SM temperature
157     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALLEDRefs; j++) {
158       outputFile << j << " " << t->GetLEDRefHighLow(j) 
159                  << " " << t->GetLEDRefAmp(j) << " " << t->GetLEDRefAmpRMS(j) 
160                  << endl;
161     }
162     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALTempSensors; j++) {
163       outputFile << j << " " << t->GetTemperature(j) << " " << t->GetTemperatureRMS(j) << endl;
164     }
165
166     // third: info for each tower
167     for (Int_t j=0; j<nAPDPerSM; j++) {
168       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
169       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
170
171       AliEMCALCalibReferenceVal * v = t->GetAPDVal(iCol, iRow);
172
173       if (swapSides) {
174         // C side, oriented differently than A side: swap is requested
175         iCol = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
176         iRow = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
177       }
178
179       outputFile << iCol << " " << iRow 
180                  << " " << v->GetHighLow() 
181                  << " " << v->GetLEDAmp() 
182                  << " " << v->GetLEDAmpRMS() << endl;
183     }
184
185   } // i, SuperModule
186
187   outputFile.close();
188
189   return;
190 }
191
192 //____________________________________________________________________________
193 void AliEMCALCalibReference::ReadRootCalibReferenceInfo(const TString &rootFileName,
194                                             Bool_t swapSides)
195 {
196   //Read data from root file. ; coordinates given on SuperModule basis
197   TFile inputFile(rootFileName, "read");  
198
199   TTree *tree = (TTree*) inputFile.Get("tree");
200
201   ReadTreeCalibReferenceInfo(tree, swapSides);
202
203   inputFile.Close();
204
205   return;
206 }
207
208 //____________________________________________________________________________
209 void AliEMCALCalibReference::ReadTreeCalibReferenceInfo(TTree *tree,
210                                             Bool_t swapSides)
211 {
212   // how many SuperModule's worth of info do we have?
213   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
214   fNSuperModule = tree->GetEntries();
215
216   Int_t iSM = 0; // SuperModule index
217   // list of values to be read
218   // first: overall values for the whole SuperModule
219   Int_t iReferenceTime= {0}; 
220   // second: additional info for LED Reference and SM temperature
221   Float_t rLEDRefAmp[AliEMCALGeoParams::fgkEMCALLEDRefs]= {0};
222   Float_t rLEDRefAmpRMS[AliEMCALGeoParams::fgkEMCALLEDRefs]= {0};
223   Int_t iLEDRefHighLow[AliEMCALGeoParams::fgkEMCALLEDRefs]= {0};
224   Float_t temperature[AliEMCALGeoParams::fgkEMCALTempSensors]= {0};
225   Float_t temperatureRMS[AliEMCALGeoParams::fgkEMCALTempSensors]= {0};
226   // third: info for each tower
227   Int_t iHighLow[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
228   Float_t rLEDAmp[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
229   Float_t rLEDAmpRMS[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
230   // end - all values
231
232   // just to make the initializations of the arrays are done correctly, let's use memset
233   memset(rLEDRefAmp, 0, sizeof(rLEDRefAmp)); 
234   memset(rLEDRefAmpRMS, 0, sizeof(rLEDRefAmpRMS)); 
235   memset(iLEDRefHighLow, 0, sizeof(iLEDRefHighLow)); 
236   memset(temperature, 0, sizeof(temperature)); 
237   memset(temperatureRMS, 0, sizeof(temperatureRMS)); 
238   memset(iHighLow, 0, sizeof(iHighLow)); 
239   memset(rLEDAmp, 0, sizeof(rLEDAmp)); 
240   memset(rLEDAmpRMS, 0, sizeof(rLEDAmpRMS)); 
241
242   // declare the branches
243   tree->SetBranchAddress("iSM", &iSM);
244   tree->SetBranchAddress("ReferenceTime", &iReferenceTime);
245   //
246   tree->SetBranchAddress("LEDRefAmp", rLEDRefAmp);
247   tree->SetBranchAddress("LEDRefAmpRMS", rLEDRefAmpRMS);
248   tree->SetBranchAddress("LEDRefHighLow", iLEDRefHighLow);
249   tree->SetBranchAddress("Temperature", temperature);
250   tree->SetBranchAddress("TemperatureRMS", temperatureRMS);
251   //
252   tree->SetBranchAddress("HighLow", iHighLow);
253   tree->SetBranchAddress("LEDAmp", rLEDAmp);
254   tree->SetBranchAddress("LEDAmpRMS", rLEDAmpRMS);
255
256   // indices for looping over the towers
257   Int_t iCol = 0;
258   Int_t iRow = 0;
259
260   for (int ient=0; ient<tree->GetEntries(); ient++) {
261     tree->GetEntry(ient);
262
263     // assume the index SuperModules come in order: i=iSM
264     AliEMCALSuperModuleCalibReference * t = (AliEMCALSuperModuleCalibReference*) fSuperModuleData[iSM];
265
266     t->SetSuperModuleNum(iSM);
267     // first, overall values
268     t->SetReferenceTime(iReferenceTime);
269
270     // second: additional info for LED references and SM temperatures
271     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALLEDRefs; j++) {
272       t->SetLEDRefAmp(j, rLEDRefAmp[j]);
273       t->SetLEDRefAmpRMS(j, rLEDRefAmpRMS[j]);
274       t->SetLEDRefHighLow(j, iLEDRefHighLow[j]);
275     }
276     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALTempSensors; j++) {
277       t->SetTemperature(j, temperature[j]);
278       t->SetTemperatureRMS(j, temperatureRMS[j]);
279     }
280
281     // third: info for each tower
282     for (Int_t j=0; j<nAPDPerSM; j++) {
283       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
284       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
285
286       // help variables: possibly modified or swapped indices
287       int iColMod = iCol;
288       int iRowMod = iRow;
289       // assume that this info is already swapped and done for this basis?
290       if (swapSides) {
291         // C side, oriented differently than A side: swap is requested
292         iColMod = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
293         iRowMod = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
294       }
295
296       AliEMCALCalibReferenceVal * v = t->GetAPDVal(iColMod, iRowMod);
297
298       v->SetHighLow(iHighLow[iCol][iRow]);
299       v->SetLEDAmp(rLEDAmp[iCol][iRow]);
300       v->SetLEDAmpRMS(rLEDAmpRMS[iCol][iRow]);
301     }
302
303   } // loop over entries
304
305   return;
306 }
307
308 //____________________________________________________________________________
309 void AliEMCALCalibReference::WriteRootCalibReferenceInfo(const TString &rootFileName,
310                                              Bool_t swapSides)
311 {
312   // write data to root file. ; coordinates given on SuperModule basis
313   TFile destFile(rootFileName, "recreate");  
314   if (destFile.IsZombie()) {
315     return;
316   }  
317   destFile.cd();
318
319   TTree *tree = new TTree("tree","");
320
321   // variables for filling the TTree
322   Int_t iSM = 0; // SuperModule index
323   // list of values to be written
324   // first: overall values for the whole SuperModule
325   Int_t iReferenceTime = 0; 
326   // second: additional info for LED Reference and SM temperature
327   Float_t rLEDRefAmp[AliEMCALGeoParams::fgkEMCALLEDRefs] = {0};
328   Float_t rLEDRefAmpRMS[AliEMCALGeoParams::fgkEMCALLEDRefs]= {0};
329   Int_t iLEDRefHighLow[AliEMCALGeoParams::fgkEMCALLEDRefs]= {0};
330   Float_t temperature[AliEMCALGeoParams::fgkEMCALTempSensors]= {0};
331   Float_t temperatureRMS[AliEMCALGeoParams::fgkEMCALTempSensors]= {0};
332   // third: info for each tower
333   Int_t iHighLow[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
334   Float_t rLEDAmp[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
335   Float_t rLEDAmpRMS[AliEMCALGeoParams::fgkEMCALCols][AliEMCALGeoParams::fgkEMCALRows]; 
336   // end - all values
337
338   // just to make the initializations of the arrays are done correctly, let's use memset
339   memset(rLEDRefAmp, 0, sizeof(rLEDRefAmp)); 
340   memset(rLEDRefAmpRMS, 0, sizeof(rLEDRefAmpRMS)); 
341   memset(iLEDRefHighLow, 0, sizeof(iLEDRefHighLow)); 
342   memset(temperature, 0, sizeof(temperature)); 
343   memset(temperatureRMS, 0, sizeof(temperatureRMS)); 
344   memset(iHighLow, 0, sizeof(iHighLow)); 
345   memset(rLEDAmp, 0, sizeof(rLEDAmp)); 
346   memset(rLEDAmpRMS, 0, sizeof(rLEDAmpRMS)); 
347
348   Int_t nAPDPerSM = AliEMCALGeoParams::fgkEMCALCols * AliEMCALGeoParams::fgkEMCALRows;
349   // for looping over towers
350   Int_t iCol = 0;
351   Int_t iRow = 0;
352
353   // declare the branches
354   // first
355   tree->Branch("iSM", &iSM, "iSM/I");
356   tree->Branch("ReferenceTime", &iReferenceTime, "ReferenceTime/I");
357   // second  
358   tree->Branch( "LEDRefAmp", &rLEDRefAmp, Form("LEDRefAmp[%d]/F", AliEMCALGeoParams::fgkEMCALLEDRefs) );
359   tree->Branch( "LEDRefAmpRMS", &rLEDRefAmpRMS, Form("LEDRefAmpRMS[%d]/F", AliEMCALGeoParams::fgkEMCALLEDRefs) );
360   tree->Branch( "LEDRefHighLow", &iLEDRefHighLow, Form("LEDRefHighLow[%d]/I", AliEMCALGeoParams::fgkEMCALLEDRefs) );
361   tree->Branch( "Temperature", &temperature, Form("Temperature[%d]/F", AliEMCALGeoParams::fgkEMCALTempSensors) );
362   tree->Branch( "TemperatureRMS", &temperatureRMS, Form("TemperatureRMS[%d]/F", AliEMCALGeoParams::fgkEMCALTempSensors) );
363   // third: info for each tower; see if a 2D array works OK or if we'll have to use 1D arrays instead 
364   tree->Branch( "HighLow", &iHighLow, Form("HighLow[%d][%d]/I", AliEMCALGeoParams::fgkEMCALCols, AliEMCALGeoParams::fgkEMCALRows) );
365   tree->Branch( "LEDAmp", &rLEDAmp, Form("LEDAmp[%d][%d]/F", AliEMCALGeoParams::fgkEMCALCols, AliEMCALGeoParams::fgkEMCALRows) );
366   tree->Branch( "LEDAmpRMS", &rLEDAmpRMS, Form("LEDAmpRMS[%d][%d]/F", AliEMCALGeoParams::fgkEMCALCols, AliEMCALGeoParams::fgkEMCALRows) );
367
368   for (iSM = 0; iSM < fNSuperModule; iSM++) {
369     AliEMCALSuperModuleCalibReference * t = (AliEMCALSuperModuleCalibReference*) fSuperModuleData[iSM];
370
371     iSM = t->GetSuperModuleNum();
372     // first, overall values
373     iReferenceTime = t->GetReferenceTime();
374
375     // second: additional info for LED references and SM temperatures
376     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALLEDRefs; j++) {
377       rLEDRefAmp[j] = t->GetLEDRefAmp(j);
378       rLEDRefAmpRMS[j] = t->GetLEDRefAmpRMS(j);
379       iLEDRefHighLow[j] = t->GetLEDRefHighLow(j);
380     }
381     for (Int_t j=0; j<AliEMCALGeoParams::fgkEMCALTempSensors; j++) {
382       temperature[j] = t->GetTemperature(j);
383       temperatureRMS[j] = t->GetTemperatureRMS(j);
384     }
385
386     // third: info for each tower
387     for (Int_t j=0; j<nAPDPerSM; j++) {
388       iCol = j / AliEMCALGeoParams::fgkEMCALRows;
389       iRow = j % AliEMCALGeoParams::fgkEMCALRows;
390
391       // help variables: possibly modified or swapped indices
392       int iColMod = iCol;
393       int iRowMod = iRow;
394       // assume that this info is already swapped and done for this basis?
395       if (swapSides) {
396         // C side, oriented differently than A side: swap is requested
397         iColMod = AliEMCALGeoParams::fgkEMCALCols-1 - iCol;
398         iRowMod = AliEMCALGeoParams::fgkEMCALRows-1 - iRow;
399       }
400
401       AliEMCALCalibReferenceVal * v = t->GetAPDVal(iCol, iRow);
402
403       iHighLow[iColMod][iRowMod] = v->GetHighLow();
404       rLEDAmp[iColMod][iRowMod] = v->GetLEDAmp();
405       rLEDAmpRMS[iColMod][iRowMod] = v->GetLEDAmpRMS();
406     }
407
408     tree->Fill();
409   } // i, SuperModule
410
411   tree->Write();
412   destFile.Close();
413
414   return;
415 }
416
417 //____________________________________________________________________________
418 AliEMCALCalibReference::~AliEMCALCalibReference()
419 {
420   fSuperModuleData.Delete();
421 }
422
423 //____________________________________________________________________________
424 AliEMCALSuperModuleCalibReference * AliEMCALCalibReference::GetSuperModuleCalibReferenceNum(Int_t supModIndex)const
425 {
426   for (int i=0; i<fNSuperModule; i++) {
427     AliEMCALSuperModuleCalibReference * t = (AliEMCALSuperModuleCalibReference*) fSuperModuleData[i];
428     if (t->GetSuperModuleNum() == supModIndex) {
429       return t;
430     }
431   }
432
433   // if we arrived here, then nothing was found.. just return a NULL pointer 
434   return NULL;
435 }
436