]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/PHOS/AliHLTPHOSBaselineAnalyzer.cxx
Coding conventions and removal of obsolete files
[u/mrichter/AliRoot.git] / HLT / PHOS / AliHLTPHOSBaselineAnalyzer.cxx
1  /**************************************************************************
2  * This file is property of and copyright by the ALICE HLT Project        * 
3  * All rights reserved.                                                   *
4  *                                                                        *
5  * Primary Authors: Oystein Djuvsland                                     *
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
17 /**
18  * Measures the baselines and calculates relevant values
19  *
20  * @file   AliHLTPHOSBaselineAnalyzer.cxx
21  * @author Oystein Djuvsland
22  * @date
23  * @brief  Baseline analyzer for PHOS HLT
24  */
25
26 #include "AliHLTPHOSBaselineAnalyzer.h"
27 #include "AliHLTPHOSBase.h"
28 #include "AliHLTPHOSBaseline.h"
29 #include "AliHLTPHOSValidCellDataStruct.h"
30 #include "AliHLTPHOSSanityInspector.h"
31 #include "AliHLTPHOSConstants.h"
32
33 #include "AliHLTPHOSRcuCellEnergyDataStruct.h"
34 #include "TTree.h"
35 #include "TClonesArray.h"
36 #include "TFile.h"
37 #include "TH1F.h"
38 #include "TH2F.h"
39
40 using namespace PhosHLTConst;
41
42 ClassImp(AliHLTPHOSBaselineAnalyzer); 
43
44 AliHLTPHOSBaselineAnalyzer::AliHLTPHOSBaselineAnalyzer() : 
45   AliHLTPHOSBase(),
46   fSampleStart(5),
47   fMaxCrazyDifference(0),       
48   fMaxSignal(0),    
49   fTreePtr(0),
50   fSanityInspector(0)      
51 {  
52   //see headerfile for documentation
53   fSanityInspector = new AliHLTPHOSSanityInspector();
54   
55   char histName[128];
56   char histTitle[128];
57   for(UInt_t x = 0; x < N_XCOLUMNS_MOD; x++)
58   {
59     for(UInt_t z = 0; z < N_ZROWS_MOD; z++)
60     {
61       for(UInt_t gain = 0; gain < N_GAINS; gain++)
62       { 
63         sprintf(histName, "sample_value_channel_x_col_%d_z_row_%d_gain_%d", x, z, gain);
64         sprintf(histTitle, "Distribution of Sample Values for Channel X: %d - Z: %d - Gain: %d", x, z, gain);
65         fChannelHistogramsPtr[x][z][gain] = new TH1F(histName, histTitle, 1024, 0, 1023);
66         sprintf(histName, "fixed_sample_value_channel_x_col_%d_z_row_%d_gain_%d", x, z, gain);
67         sprintf(histTitle, "Distribution of Corrected Sample Values for Channel X: %d - Z: %d - Gain: %d", x, z, gain);
68         fFixedChannelHistogramsPtr[x][z][gain] = new TH1F(histName, histTitle, 1024, 0, 1023);
69       }
70     }
71   }
72   
73   fRMSHistogramPtr = new TH1F("RMSHist", "RMS Values for All Channels", 1023, 0, 1023);
74   fRMSMapHighGainHistogramPtr = new TH2F("RMSHGMapHist", "Map of RMS Values for High Gain Channels", 64, 0, 63, 56, 0, 55);
75   fRMSMapLowGainHistogramPtr = new TH2F("RMSLGMapHist", "Map of RMS Values for Low Gain Channels", 64, 0, 63, 56, 0, 55);
76   
77   fFixedRMSHistogramPtr = new TH1F("fixedRMSHist", "Corrected RMS Values for All Channels", 1023, 0, 1023);
78   fFixedRMSMapHighGainHistogramPtr = new TH2F("fixedRMSHGMapHist", "Map of Corrected RMS Values for High Gain Channels", 64, 0, 63, 56, 0, 55);
79   fFixedRMSMapLowGainHistogramPtr = new TH2F("fixedRMSLGMapHist", "Map of Corrected RMS Values for Low Gain Channels", 64, 0, 63, 56, 0, 55);
80   
81   ResetBaselines();
82   ResetAccumulatedBaselines();
83
84
85 AliHLTPHOSBaselineAnalyzer::~AliHLTPHOSBaselineAnalyzer() 
86
87   //see headerfile for documentation
88
89
90
91 void
92 AliHLTPHOSBaselineAnalyzer::CalculateRcuBaselines(AliHLTPHOSRcuCellEnergyDataStruct* rcuData)
93 {
94   //see headerfile for documentation
95   Int_t xOff = rcuData->fRcuX * N_XCOLUMNS_RCU;
96   Int_t zOff = rcuData->fRcuZ * N_ZROWS_RCU;
97   Float_t baseline = 0;
98
99   AliHLTPHOSValidCellDataStruct *data = rcuData->fValidData;
100
101   for(Int_t i = 0; i < rcuData->fCnt; i++)
102     {
103       baseline = CalculateChannelBaseline(&(data[i]), xOff, zOff);
104       if(!(baseline < 0))
105       {
106         CalculateAccumulatedChannelBaseline(data[i].fX + xOff, data[i].fZ + zOff, data[i].fGain, baseline);
107       }
108     }
109 }
110   
111
112
113 Float_t
114 AliHLTPHOSBaselineAnalyzer::CalculateChannelBaseline(AliHLTPHOSValidCellDataStruct *cellData, Int_t xOff, Int_t zOff) 
115
116   //see headerfile for documentation
117   Int_t crazyness = 0;
118   Int_t total = 0;
119   Int_t *data = cellData->fData;
120   for(Int_t i = fSampleStart; i < fNSamples; i++)
121   { 
122     /*
123     if(cellData->fGain == 0)
124     {
125       fChannelLowGainHistogramsPtr[cellData->fX + xOff][cellData->fZ + zOff]->Fill(data[i]);
126     }
127     else
128     {
129       fChannelHighGainHistogramsPtr[cellData->fX + xOff][cellData->fZ + zOff]->Fill(data[i]);
130     }
131     */
132     fChannelHistogramsPtr[cellData->fX + xOff][cellData->fZ + zOff][cellData->fGain]->Fill(data[i]);
133   }
134   
135   if(cellData->fCrazyness == 0)
136   {
137        crazyness = fSanityInspector->CheckAndHealInsanity(data, fNSamples);
138   }
139   if(crazyness < 0)
140     return crazyness;
141   
142   for(Int_t j = fSampleStart; j < fNSamples; j++)
143     {
144       if(data[j] > fMaxSignal)
145         return -data[j];
146     
147       fFixedChannelHistogramsPtr[cellData->fX + xOff][cellData->fZ + zOff][cellData->fGain]->Fill(data[j]);
148       total += data[j];
149     }
150   fBaselines[cellData->fX + xOff][cellData->fZ + zOff][cellData->fGain] = (float)total / (fNSamples - fSampleStart);
151     
152   return (float)total/(fNSamples - fSampleStart);
153 }
154
155 void 
156 AliHLTPHOSBaselineAnalyzer::CalculateAccumulatedChannelBaseline(Int_t x, Int_t z, Int_t gain, Float_t baseline)
157 {
158   //see headerfile for documentation
159   Float_t oldBaseline = fAccumulatedBaselines[x][z][gain][0];
160   Float_t nEntries = fAccumulatedBaselines[x][z][gain][1];
161   
162   Float_t total = nEntries * oldBaseline + baseline;
163   nEntries++;
164   fAccumulatedBaselines[x][z][gain][0] = total / nEntries;
165   fAccumulatedBaselines[x][z][gain][1] = nEntries;
166 }
167
168 void
169 AliHLTPHOSBaselineAnalyzer::CalculateChannelsBaselineRMS()
170 {
171
172   //see headerfile for documentation
173   for(UInt_t x = 0; x < N_XCOLUMNS_MOD; x++)
174   {
175     for(UInt_t z = 0; z < N_ZROWS_MOD; z++)
176     {
177       for(UInt_t gain = 0; gain < N_GAINS; gain++)
178       {
179         fRMSHistogramPtr->Fill(fChannelHistogramsPtr[x][z][gain]->GetRMS());  
180         if(gain == 1)
181           fRMSMapHighGainHistogramPtr->SetBinContent(x, z, fChannelHistogramsPtr[x][z][gain]->GetRMS());
182         else
183           fRMSMapLowGainHistogramPtr->SetBinContent(x, z, fChannelHistogramsPtr[x][z][gain]->GetRMS());
184                    
185         fFixedRMSHistogramPtr->Fill(fFixedChannelHistogramsPtr[x][z][gain]->GetRMS());  
186         if(gain == 1)
187           fFixedRMSMapHighGainHistogramPtr->SetBinContent(x, z, fFixedChannelHistogramsPtr[x][z][gain]->GetRMS());
188         else
189           fFixedRMSMapLowGainHistogramPtr->SetBinContent(x, z, fFixedChannelHistogramsPtr[x][z][gain]->GetRMS());
190       }
191     }
192   }
193 }    
194  
195
196 void 
197 AliHLTPHOSBaselineAnalyzer::SetRootObjects(TTree *tree, TClonesArray *baselineArray)
198 {
199   //see headerfile for documentation
200   fTreePtr = tree;
201   fBaselineArrayPtr = baselineArray;
202   tree->Branch("Baselines", &fBaselineArrayPtr);
203 }
204   
205
206 void
207 AliHLTPHOSBaselineAnalyzer::FillTree()
208 {
209   //see headerfile for documentation
210   AliHLTPHOSBaseline *baseline;
211   Int_t n = 0;
212
213   for(UInt_t x = 0; x < N_XCOLUMNS_MOD; x++)
214     {
215       for(UInt_t z = 0; z < N_ZROWS_MOD; z++)
216         {
217           for(UInt_t gain = 0; gain < N_GAINS; gain++)
218             {
219               baseline = (AliHLTPHOSBaseline*)fBaselineArrayPtr->New(n);
220               baseline->SetX(x);
221               baseline->SetZ(z);
222               baseline->SetGain(gain);
223               baseline->SetBaseline(fAccumulatedBaselines[x][z][gain][0]);
224               baseline->SetEntries((Int_t)fAccumulatedBaselines[x][z][gain][1]);
225               n++;
226             }
227         }
228     }
229   fTreePtr->Fill();
230   fBaselineArrayPtr->Clear();
231 }
232
233 void
234 AliHLTPHOSBaselineAnalyzer::WriteAccumulatedBaselines(const Char_t* filename)
235 {
236   //see headerfile for documentation
237   TFile *baselineFile = TFile::Open(filename, "recreate");
238   fTreePtr->Write();
239   baselineFile->Close();
240 }
241
242 void 
243 AliHLTPHOSBaselineAnalyzer::WriteChannelHistograms(const Char_t* filename)
244 {
245   //see headerfile for documentation
246   TFile *file = TFile::Open(filename, "recreate");
247   
248   for(UInt_t x = 0; x < N_XCOLUMNS_MOD; x++)
249   {
250     for(UInt_t z = 0; z < N_ZROWS_MOD; z++)
251     {
252       for(UInt_t gain = 0; gain < N_GAINS; gain++)
253       {
254         fChannelHistogramsPtr[x][z][gain]->Write();
255       }
256     }
257   }
258   file->Close();
259 }
260
261
262 void
263 AliHLTPHOSBaselineAnalyzer::WriteRMSHistogram(const Char_t* filename)
264 {
265   //see headerfile for documentation
266   TFile *file = TFile::Open(filename, "recreate");
267   fRMSHistogramPtr->Write();
268   fRMSMapHighGainHistogramPtr->Write();
269   fRMSMapLowGainHistogramPtr->Write();
270
271   fFixedRMSHistogramPtr->Write();
272   fFixedRMSMapHighGainHistogramPtr->Write();
273   fFixedRMSMapLowGainHistogramPtr->Write();
274   
275   file->Close();
276 }
277
278 void 
279 AliHLTPHOSBaselineAnalyzer::ResetBaselines()
280 {
281   //see headerfile for documentation
282    for(UInt_t x = 0; x < N_XCOLUMNS_MOD; x++)
283     {
284       for(UInt_t z = 0; z < N_ZROWS_MOD; z++)
285         {
286           for(UInt_t gain = 0; gain < N_GAINS; gain++)
287             {
288               fBaselines[x][z][gain] = 0;    
289             }
290         }
291     }
292 }
293
294 void 
295 AliHLTPHOSBaselineAnalyzer::ResetChannelCount()
296 {
297   //see headerfile for documentation
298   fChannelCount = 0;
299 }
300     
301 void 
302 AliHLTPHOSBaselineAnalyzer::ResetAccumulatedBaselines()
303 {
304    for(UInt_t x = 0; x < N_XCOLUMNS_MOD; x++)
305     {
306       for(UInt_t z = 0; z < N_ZROWS_MOD; z++)
307         {
308           for(UInt_t gain = 0; gain < N_GAINS; gain++)
309             {
310              fAccumulatedBaselines[x][z][gain][0] = 0;
311              fAccumulatedBaselines[x][z][gain][1] = 0;
312             }
313         }
314     }
315 }
316
317 void 
318 AliHLTPHOSBaselineAnalyzer::SetMaxCrazyDifference(Int_t diff)
319 {
320   //see headerfile for documentation
321   fMaxCrazyDifference = diff;  
322   fSanityInspector->SetMaxDifference(diff);
323 }
324
325
326
327 /*   
328 void 
329 AliHLTPHOSBaselineAnalyzer::SetChannelsHistograms(TH1F *channelLowGainHistArray[N_XCOLUMNS_MOD][N_ZROWS_MOD], TH1F *channelHighGainHistArray[N_XCOLUMNS_MOD][N_ZROWS_MOD])
330
331   for(UInt_t x = 0; x < N_XCOLUMNS_MOD; x++)
332   {
333     for(UInt_t z = 0; z < N_ZROWS_MOD; z++)
334     {
335       for(UInt_t gain = 0; gain < N_GAINS; gain++)
336       { 
337         fChannelLowGainHistogramsPtr[x][z][gain] = channelLowGainHistArray[x][z][gain];
338         fChannelHighGainHistogramsPtr[x][z][gain] = channelHighGainHistArray[x][z][gain];
339       }
340     }
341   }
342 }
343 */