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