]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/Base/AliTPCCalibPedestal.cxx
Fixes for wrong use of const causing PW.CAST_TO_QUALIFIED_TYPE defect in Coverity
[u/mrichter/AliRoot.git] / TPC / Base / AliTPCCalibPedestal.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
17 /* $Id$ */
18
19
20 //Root includes
21 #include <TH1F.h>
22 #include <TH2F.h>
23 #include <TString.h>
24 #include <TMath.h>
25 #include <TF1.h>
26 #include <TRandom.h>
27 #include <TDirectory.h>
28 #include <TFile.h>
29 #include <TMap.h>
30 //AliRoot includes
31 #include "AliRawReader.h"
32 #include "AliRawReaderRoot.h"
33 #include "AliRawReaderDate.h"
34 #include "AliTPCCalROC.h"
35 #include "AliTPCROC.h"
36 #include "AliMathBase.h"
37 #include "TTreeStream.h"
38
39 //date
40 #include "event.h"
41
42 //header file
43 #include "AliTPCCalibPedestal.h"
44
45
46 ///////////////////////////////////////////////////////////////////////////////////////
47 //          Implementation of the TPC pedestal and noise calibration
48 //
49 //   Origin: Jens Wiechula, Marian Ivanov   J.Wiechula@gsi.de, Marian.Ivanov@cern.ch
50 // 
51 // 
52 // *************************************************************************************
53 // *                                Class Description                                  *
54 // *************************************************************************************
55 //
56 // Working principle:
57 // ------------------
58 // Raw pedestal data is processed by calling one of the ProcessEvent(...) functions
59 // (see below). These in the end call the Update(...) function, where the data is filled
60 // into histograms.
61 //
62 // For each ROC one TH2F histo (ROC channel vs. ADC channel) is created when
63 // it is filled for the first time (GetHistoPedestal(ROC,kTRUE)). All histos are stored in the
64 // TObjArray fHistoPedestalArray.
65 //
66 // For a fast filling of the histogram the corresponding bin number of the channel and ADC channel
67 // is computed by hand and the histogram array is accessed directly via its pointer.
68 // ATTENTION: Doing so the the entry counter of the histogram is not increased
69 //            this means that e.g. the colz draw option gives an empty plot unless
70 //          calling 'histo->SetEntries(1)' before drawing.
71 //
72 // After accumulating the desired statistics the Analyse() function has to be called.
73 // Whithin this function the pedestal and noise values are calculated for each pad, using
74 // the fast gaus fit function  AliMathBase::FitGaus(...), and the calibration
75 // storage classes (AliTPCCalROC) are filled for each ROC.
76 // The calibration information is stored in the TObjArrays fCalRocArrayPedestal and fCalRocArrayRMS;
77 //
78 //
79 //
80 // User interface for filling data:
81 // --------------------------------
82 //
83 // To Fill information one of the following functions can be used:
84 //
85 // Bool_t ProcessEvent(eventHeaderStruct *event);
86 //   - process Date event
87 //   - use AliTPCRawReaderDate and call ProcessEvent(AliRawReader *rawReader)
88 //
89 // Bool_t ProcessEvent(AliRawReader *rawReader);
90 //  - process AliRawReader event
91 //   - use AliTPCRawStreamV3 to loop over data and call ProcessEvent(AliTPCRawStreamV3 *rawStream)
92 //
93 // Bool_t ProcessEvent(AliTPCRawStreamV3 *rawStream);
94 //   - process event from AliTPCRawStreamV3
95 //   - call Update function for signal filling
96 //
97 // Int_t Update(const Int_t isector, const Int_t iRow, const Int_t
98 //              iPad,  const Int_t iTimeBin, const Float_t signal);
99 //   - directly  fill signal information (sector, row, pad, time bin, pad)
100 //     to the reference histograms
101 //
102 // It is also possible to merge two independently taken calibrations using the function
103 //
104 // void Merge(AliTPCCalibPedestal *ped)
105 //   - copy histograms in 'ped' if the do not exist in this instance
106 //   - Add histograms in 'ped' to the histograms in this instance if the allready exist
107 //   - After merging call Analyse again!
108 //
109 //
110 //
111 // -- example: filling data using root raw data:
112 // void fillPedestal(Char_t *filename)
113 // {
114 //    rawReader = new AliRawReaderRoot(fileName);
115 //    if ( !rawReader ) return;
116 //    AliTPCCalibPedestal *calib = new AliTPCCalibPedestal;
117 //    while (rawReader->NextEvent()){
118 //      calib->ProcessEvent(rawReader);
119 //    }
120 //    calib->Analyse();
121 //    calib->DumpToFile("PedestalData.root");
122 //    delete rawReader;
123 //    delete calib;
124 // }
125 //
126 //
127 // What kind of information is stored and how to retrieve them:
128 // ------------------------------------------------------------
129 //
130 // - Accessing the 'Reference Histograms' (pedestal distribution histograms):
131 //
132 // TH2F *GetHistoPedestal(Int_t sector);
133 //
134 // - Accessing the calibration storage objects:
135 //
136 // AliTPCCalROC *GetCalRocPedestal(Int_t sector);  - for the pedestal values, mean from gaus fit
137 // AliTPCCalROC *GetCalRocSigma(Int_t sector);     - for the Noise values, sigma from guas fit
138 // AliTPCCalROC *GetCalRocMean(Int_t sector);  - for the pedestal values, truncated mean
139 // AliTPCCalROC *GetCalRocRMS(Int_t sector);     - for the Noise values, rms from truncated mean
140 //
141 // example for visualisation:
142 // if the file "PedestalData.root" was created using the above example one could do the following:
143 //
144 // TFile filePedestal("PedestalData.root")
145 // AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)filePedestal->Get("AliTPCCalibPedestal");
146 // ped->GetCalRocPedestal(0)->Draw("colz");
147 // ped->GetCalRocRMS(0)->Draw("colz");
148 //
149 // or use the AliTPCCalPad functionality:
150 // AliTPCCalPad padPedestal(ped->GetCalPadPedestal());
151 // AliTPCCalPad padNoise(ped->GetCalPadRMS());
152 // padPedestal->MakeHisto2D()->Draw("colz");  //Draw A-Side Pedestal Information
153 // padNoise->MakeHisto2D()->Draw("colz");  //Draw A-Side Noise Information
154 //
155 /*
156  example: fill pedestal with gausschen noise
157  AliTPCCalibPedestal ped;
158  ped.TestEvent();
159  ped.Analyse();
160  //Draw output;
161  TCanvas* c1 = new TCanvas;
162  c1->Divide(1,2);
163  c1->cd(1);
164  ped.GetHistoPedestal(0)->SetEntries(1); //needed in order for colz to work, reason: fast filling does not increase the entries counter
165  ped.GetHistoPedestal(0)->Draw("colz");
166  c1->cd(2);
167  ped.GetHistoPedestal(36)->SetEntries(1); //needed in order for colz to work, reason: fast filling does not increase the entries counter
168  ped.GetHistoPedestal(36)->Draw("colz");
169  TCanvas* c2 = new TCanvas;
170  c2->Divide(2,2);
171  c2->cd(1);
172  ped.GetCalRocPedestal(0)->Draw("colz");
173  c2->cd(2);
174  ped.GetCalRocRMS(0)->Draw("colz");
175  c2->cd(3);
176  ped.GetCalRocPedestal(36)->Draw("colz");
177  c2->cd(4);
178  ped.GetCalRocRMS(36)->Draw("colz");
179 */
180 //
181 // Time dependent pedestals:
182 //
183 // If wished there is the possibility to calculate for each channel and time bin
184 // the mean pedestal [pedestals(t)]. This is done by
185 //
186 // 1) setting SetTimeAnalysis(kTRUE),
187 // 2) processing the data by looping over the events using ProcessEvent(..)
188 // 3) calling the Analyse() and AnalyseTime(nevents) functions (providing nevents)
189 // 4) getting the pedestals(t) using   TArrayF **timePed = calibPedestal.GetTimePedestals();
190 // 5) looking at values using   timePed[row][pad].At(timebin)
191 //
192 // This functionality is intended to be used on an LDC bu the detector algorithm
193 // (TPCPEDESTALda) to generate a data set used for configuration of the pattern
194 // memory for baseline subtraction in the ALTROs. Later the information should also
195 // be stored as reference data.
196 //
197
198
199 ClassImp(AliTPCCalibPedestal)
200
201 AliTPCCalibPedestal::AliTPCCalibPedestal() : 
202   AliTPCCalibRawBase(),
203   fAdcMin(1),
204   fAdcMax(100),
205   fAnaMeanDown(0.),
206   fAnaMeanUp(1.),
207   fTimeAnalysis(kFALSE),
208   fCalRocArrayPedestal(72),
209   fCalRocArraySigma(72),
210   fHistoPedestalArray(72),
211   fTimeSignal(NULL),
212   fCalRocArrayMean(72),
213   fCalRocArrayRMS(72)
214 {
215   //
216   // default constructor
217   //
218   SetNameTitle("AliTPCCalibPedestal","AliTPCCalibPedestal");
219   fFirstTimeBin=60;
220   fLastTimeBin=1000;
221 }
222
223
224 //_____________________________________________________________________
225 AliTPCCalibPedestal::AliTPCCalibPedestal(const AliTPCCalibPedestal &ped) : 
226   AliTPCCalibRawBase(ped),
227   fAdcMin(ped.GetAdcMin()),
228   fAdcMax(ped.GetAdcMax()),
229   fAnaMeanDown(ped.fAnaMeanDown),
230   fAnaMeanUp(ped.fAnaMeanUp),
231   fTimeAnalysis(ped.fTimeAnalysis),
232   fCalRocArrayPedestal(72),
233   fCalRocArraySigma(72),
234   fHistoPedestalArray(72),
235   fTimeSignal(ped.fTimeSignal),
236   fCalRocArrayMean(72),
237   fCalRocArrayRMS(72)
238 {
239   //
240   // copy constructor
241   //
242   for (Int_t iSec = 0; iSec < 72; ++iSec){
243     const AliTPCCalROC *calPed = (AliTPCCalROC*)ped.fCalRocArrayPedestal.UncheckedAt(iSec);
244     const AliTPCCalROC *calRMS = (AliTPCCalROC*)ped.fCalRocArrayRMS.UncheckedAt(iSec);
245     const TH2F         *hPed   = (TH2F*)ped.fHistoPedestalArray.UncheckedAt(iSec);
246     
247     if ( calPed != 0x0 ) fCalRocArrayPedestal.AddAt(new AliTPCCalROC(*calPed), iSec);
248     if ( calRMS != 0x0 ) fCalRocArrayRMS.AddAt(new AliTPCCalROC(*calRMS), iSec);
249     
250     if ( hPed != 0x0 ){
251       TH2F *hNew = new TH2F(*hPed);
252       hNew->SetDirectory(0);
253       fHistoPedestalArray.AddAt(hNew,iSec);
254     }
255   }
256 }
257 AliTPCCalibPedestal::AliTPCCalibPedestal(const TMap *config): 
258   AliTPCCalibRawBase(),
259   fAdcMin(1),
260   fAdcMax(100),
261   fAnaMeanDown(0.),
262   fAnaMeanUp(1.),
263   fTimeAnalysis(kFALSE),
264   fCalRocArrayPedestal(72),
265   fCalRocArraySigma(72),
266   fHistoPedestalArray(72),
267   fTimeSignal(NULL),
268   fCalRocArrayMean(72),
269   fCalRocArrayRMS(72)  
270 {
271  //
272  // This constructor uses a TMap for setting some parametes
273  //
274   SetNameTitle("AliTPCCalibPedestal","AliTPCCalibPedestal");
275   fFirstTimeBin=60;
276   fLastTimeBin=1000;
277   if (config->GetValue("FirstTimeBin")) fFirstTimeBin = ((TObjString*)config->GetValue("FirstTimeBin"))->GetString().Atoi();
278   if (config->GetValue("LastTimeBin"))  fLastTimeBin = ((TObjString*)config->GetValue("LastTimeBin"))->GetString().Atoi();
279   if (config->GetValue("AdcMin"))       fAdcMin = ((TObjString*)config->GetValue("AdcMin"))->GetString().Atoi();
280   if (config->GetValue("AdcMax"))       fAdcMax = ((TObjString*)config->GetValue("AdcMax"))->GetString().Atoi();
281   if (config->GetValue("TimeAnalysis")) SetTimeAnalysis(((TObjString*)config->GetValue("TimeAnalysis"))->GetString().Atoi());
282
283
284
285 //_____________________________________________________________________
286 AliTPCCalibPedestal& AliTPCCalibPedestal::operator = (const  AliTPCCalibPedestal &source)
287 {
288   //
289   // assignment operator
290   //
291   if (&source == this) return *this;
292   new (this) AliTPCCalibPedestal(source);
293
294   return *this;
295 }
296
297
298 //_____________________________________________________________________
299 AliTPCCalibPedestal::~AliTPCCalibPedestal() 
300 {
301   //
302   // destructor
303   //
304
305   fCalRocArrayPedestal.Delete();
306   fCalRocArrayRMS.Delete();
307   fCalRocArraySigma.Delete();
308   fHistoPedestalArray.Delete();
309
310   if ( fTimeSignal ) {
311     for (Int_t i = 0; i < 159; i++) {
312       delete [] fTimeSignal[i];
313       fTimeSignal[i] = 0;
314     }
315     delete [] fTimeSignal;
316     fTimeSignal = 0;
317   }
318
319   // do not delete fMapping, because we do not own it.
320
321 }
322
323
324 //_____________________________________________________________________
325 void AliTPCCalibPedestal::SetTimeAnalysis(Bool_t time)
326 {
327   //
328   // Use time dependent analysis: Pedestals are analysed as a function
329   // of the drift time. There is one mean value generated for each time
330   // bin and each channel. It can be used as reference data and for
331   // configuration of the ALTRO pattern memory for baseline subtraction.
332   //
333   // ATTENTION: Use only on LDC in TPCPEDESTALda! On a LDC we get data
334   // only from one sector. For the full TPC we would need a lot of
335   // memory (36*159*140*1024*4bytes = 3.3GB)!
336   //
337
338   fTimeAnalysis = time;
339
340   if ( !fTimeAnalysis ) return;
341
342   // prepare array for one sector (159*140*1024*4bytes = 92MB):
343   fTimeSignal = new TArrayF*[159];
344   for (Int_t i = 0; i < 159; i++) {  // padrows
345     fTimeSignal[i] = new TArrayF[140];
346     for (Int_t j = 0; j < 140; j++) {  // pads per row
347       fTimeSignal[i][j].Set(1024);
348       for (Int_t k = 0; k < 1024; k++) {  // time bins per pad
349         fTimeSignal[i][j].AddAt(0., k);
350       }
351     }      
352   }
353 }
354
355
356 //_____________________________________________________________________
357 Int_t AliTPCCalibPedestal::Update(const Int_t icsector, 
358                                   const Int_t icRow,
359                                   const Int_t icPad,
360                                   const Int_t icTimeBin,
361                                   const Float_t csignal)
362 {
363   //
364   // Signal filling method
365   //
366   if (icRow<0) return 0;
367   if (icPad<0) return 0;
368   if (icTimeBin<0) return 0;
369  
370   // Time dependent pedestals
371   if ( fTimeAnalysis ) {
372     if ( icsector < 36 ) // IROC
373       fTimeSignal[icRow][icPad].AddAt(fTimeSignal[icRow][icPad].At(icTimeBin)+csignal, icTimeBin);
374     else 
375       fTimeSignal[icRow+63][icPad].AddAt(fTimeSignal[icRow+63][icPad].At(icTimeBin)+csignal, icTimeBin);
376   }
377   //return if we are out of the specified time bin or adc range
378   if ( (icTimeBin>fLastTimeBin) || (icTimeBin<fFirstTimeBin) ) return 0;
379   if ( ((Int_t)csignal>fAdcMax) || ((Int_t)csignal<fAdcMin)  ) return 0;
380
381   Int_t iChannel  = fROC->GetRowIndexes(icsector)[icRow]+icPad; //  global pad position in sector
382
383   // fast filling method
384   // Attention: the entry counter of the histogram is not increased
385   //            this means that e.g. the colz draw option gives an empty plot
386   Int_t bin = (iChannel+1)*(fAdcMax-fAdcMin+2)+((Int_t)csignal-fAdcMin+1);
387
388   GetHistoPedestal(icsector,kTRUE)->GetArray()[bin]++;
389
390   return 0;
391 }
392
393
394 //_____________________________________________________________________
395 Bool_t AliTPCCalibPedestal::TestEvent() 
396 {
397   //
398   //  Test event loop
399   // fill one oroc and one iroc with random gaus
400   //
401
402   gRandom->SetSeed(0);
403
404   for (UInt_t iSec=0; iSec<72; ++iSec){
405     if (iSec%36>0) continue;
406     for (UInt_t iRow=0; iRow < fROC->GetNRows(iSec); ++iRow){
407       for (UInt_t iPad=0; iPad < fROC->GetNPads(iSec,iRow); ++iPad){
408         for (UInt_t iTimeBin=0; iTimeBin<1024; ++iTimeBin){
409           Float_t signal=(Int_t)(iRow+3+gRandom->Gaus(0,.7));
410           if ( signal>0 )Update(iSec,iRow,iPad,iTimeBin,signal);
411         }
412       }
413     }
414   }
415   return kTRUE;
416 }
417
418
419 //_____________________________________________________________________
420 TH2F* AliTPCCalibPedestal::GetHisto(Int_t sector, TObjArray *arr, 
421                                     Int_t nbinsY, Float_t ymin, Float_t ymax,
422                                     const Char_t *type, Bool_t force)
423 {
424     //
425     // return pointer to Q histogram
426     // if force is true create a new histogram if it doesn't exist allready
427     //
428     if ( !force || arr->UncheckedAt(sector) )
429       return (TH2F*)arr->UncheckedAt(sector);
430
431     // if we are forced and histogram doesn't yes exist create it
432     // new histogram with Q calib information. One value for each pad!
433     TH2F* hist = new TH2F(Form("hCalib%s%.2d",type,sector),
434                           Form("%s calibration histogram sector %.2d;ADC channel;Channel (pad)",type,sector),
435                           nbinsY, ymin, ymax,
436                           fROC->GetNChannels(sector),0,fROC->GetNChannels(sector)
437                          );
438     hist->SetDirectory(0);
439     arr->AddAt(hist,sector);
440     return hist;
441 }
442
443
444 //_____________________________________________________________________
445 TH2F* AliTPCCalibPedestal::GetHistoPedestal(Int_t sector, Bool_t force) 
446 {
447     //
448     // return pointer to T0 histogram
449     // if force is true create a new histogram if it doesn't exist allready
450     //
451     TObjArray *arr = &fHistoPedestalArray;
452     return GetHisto(sector, arr, fAdcMax-fAdcMin, fAdcMin, fAdcMax, "Pedestal", force);
453 }
454
455
456 //_____________________________________________________________________
457 AliTPCCalROC* AliTPCCalibPedestal::GetCalRoc(Int_t sector, TObjArray* arr, Bool_t force) 
458 {
459     //
460     // return pointer to ROC Calibration
461     // if force is true create a new histogram if it doesn't exist allready
462     //
463     if ( !force || arr->UncheckedAt(sector) )
464         return (AliTPCCalROC*)arr->UncheckedAt(sector);
465
466     // if we are forced and the histogram doesn't yet exist create it
467
468     // new AliTPCCalROC for T0 information. One value for each pad!
469     AliTPCCalROC *croc = new AliTPCCalROC(sector);
470     arr->AddAt(croc,sector);
471     return croc;
472 }
473
474
475 //_____________________________________________________________________
476 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocPedestal(Int_t sector, Bool_t force) 
477 {
478     //
479     // return pointer to ROC with Pedestal data
480     // if force is true create a new histogram if it doesn't exist allready
481     //
482     TObjArray *arr = &fCalRocArrayPedestal;
483     return GetCalRoc(sector, arr, force);
484 }
485
486
487 //_____________________________________________________________________
488 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocSigma(Int_t sector, Bool_t force) 
489 {
490     //
491     // return pointer to  ROC with signal witdth in sigma
492     // if force is true create a new histogram if it doesn't exist allready
493     //
494     TObjArray *arr = &fCalRocArraySigma;
495     return GetCalRoc(sector, arr, force);
496 }
497 //_____________________________________________________________________
498 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocMean(Int_t sector, Bool_t force)
499 {
500   //
501     // return pointer to ROC with signal mean information
502     // if force is true create a new histogram if it doesn't exist allready
503   //
504   TObjArray *arr = &fCalRocArrayMean;
505   return GetCalRoc(sector, arr, force);
506 }
507
508 //_____________________________________________________________________
509 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocRMS(Int_t sector, Bool_t force) 
510 {
511   //
512     // return pointer to signal width ROC Calibration
513     // if force is true create a new histogram if it doesn't exist allready
514   //
515   TObjArray *arr = &fCalRocArrayRMS;
516   return GetCalRoc(sector, arr, force);
517 }
518
519
520 //_____________________________________________________________________
521 void AliTPCCalibPedestal::Merge(AliTPCCalibPedestal * const ped)
522 {
523   //
524   //  Merge reference histograms of sig to the current AliTPCCalibPedestal
525   //
526   MergeBase(ped);
527   // merge histograms
528   for (Int_t iSec=0; iSec<72; ++iSec){
529     TH2F *hRefPedMerge   = ped->GetHistoPedestal(iSec);
530     
531     if ( hRefPedMerge ){
532       TDirectory *dir = hRefPedMerge->GetDirectory(); hRefPedMerge->SetDirectory(0);
533       TH2F *hRefPed   = GetHistoPedestal(iSec);
534       if ( hRefPed ) hRefPed->Add(hRefPedMerge);
535       else {
536         TH2F *hist = new TH2F(*hRefPedMerge);
537         hist->SetDirectory(0);
538         fHistoPedestalArray.AddAt(hist, iSec);
539       }
540       hRefPedMerge->SetDirectory(dir);
541     }
542   }
543   
544   // merge array
545   // ...
546   
547 }
548
549 //_____________________________________________________________________
550 Long64_t AliTPCCalibPedestal::Merge(TCollection * const list)
551 {
552   //
553   // Merge all objects of this type in list
554   //
555   
556   Long64_t nmerged=1;
557   
558   TIter next(list);
559   AliTPCCalibPedestal *ce=0;
560   TObject *o=0;
561   
562   while ( (o=next()) ){
563     ce=dynamic_cast<AliTPCCalibPedestal*>(o);
564     if (ce){
565       Merge(ce);
566       ++nmerged;
567     }
568   }
569   
570   return nmerged;
571 }
572
573 //_____________________________________________________________________
574 void AliTPCCalibPedestal::Analyse() 
575 {
576   //
577   //  Calculate calibration constants
578   //
579
580   Int_t nbinsAdc = fAdcMax-fAdcMin;
581
582   TVectorD param(4);
583   TMatrixD dummy(3,3);
584
585   TH1F *hChannel=new TH1F("hChannel","hChannel",nbinsAdc,fAdcMin,fAdcMax);
586   
587   Float_t *arrayhP=0;  
588
589   for (Int_t iSec=0; iSec<72; ++iSec){
590     TH2F *hP = GetHistoPedestal(iSec);
591     if ( !hP ) continue;
592
593     AliTPCCalROC *rocPedestal = GetCalRocPedestal(iSec,kTRUE);
594     AliTPCCalROC *rocSigma    = GetCalRocSigma(iSec,kTRUE);
595     AliTPCCalROC *rocMean     = GetCalRocMean(iSec,kTRUE);
596     AliTPCCalROC *rocRMS      = GetCalRocRMS(iSec,kTRUE);
597
598     arrayhP = hP->GetArray();
599     UInt_t nChannels = fROC->GetNChannels(iSec);
600
601     for (UInt_t iChannel=0; iChannel<nChannels; ++iChannel){
602       Int_t offset = (nbinsAdc+2)*(iChannel+1)+1;
603       //calculate mean and sigma using a gaus fit
604       //Double_t ret =
605       AliMathBase::FitGaus(arrayhP+offset,nbinsAdc,fAdcMin,fAdcMax,&param,&dummy);
606       // if the fitting failed set noise and pedestal to 0
607       // is now done in AliMathBase::FitGaus !
608 //       if ( ret == -4 ) {
609 //      param[1]=0;
610 //      param[2]=0;
611 //       }
612       if ( param[1]<fAdcMin || param[1]>fAdcMax ){
613         param[1]=0;
614         param[2]=0;
615       }
616       rocPedestal->SetValue(iChannel,param[1]);
617       rocSigma->SetValue(iChannel,param[2]);
618       //calculate mean and RMS using a truncated means
619       hChannel->Set(nbinsAdc+2,arrayhP+offset-1);
620       hChannel->SetEntries(param[3]);
621       param[1]=0;
622       param[2]=0;
623       if ( param[3]>0 ) AliMathBase::TruncatedMean(hChannel,&param,fAnaMeanDown,fAnaMeanUp);
624       rocMean->SetValue(iChannel,param[1]);
625       rocRMS->SetValue(iChannel,param[2]);
626     }
627   }
628   delete hChannel;
629 }
630
631
632 //_____________________________________________________________________
633 void AliTPCCalibPedestal::AnalyseTime(Int_t nevents)
634 {
635   //
636   // Calculate for each channel and time bin the mean pedestal. This
637   // is used on LDC by TPCPEDESTALda to generate data used for configuration
638   // of the pattern memory for baseline subtraction in the ALTROs.
639   //
640
641   if ( nevents <= 0 ) return;
642   if ( fTimeAnalysis ) {
643     for (Int_t i = 0; i < 159; i++) {  // padrows
644       for (Int_t j = 0; j < 140; j++) {  // pads per row
645         for (Int_t k = 0; k < 1024; k++) {  // time bins per pad
646           fTimeSignal[i][j].AddAt(fTimeSignal[i][j].At(k)/(Float_t)nevents, k);
647         }
648       }
649     }
650   }
651 }