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