]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCCalibPedestal.cxx
Adding analysis task involving for TPC calibration
[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 //AliRoot includes
30 #include "AliRawReader.h"
31 #include "AliRawReaderRoot.h"
32 #include "AliRawReaderDate.h"
33 #include "AliTPCRawStream.h"
34 #include "AliTPCCalROC.h"
35 #include "AliTPCROC.h"
36 #include "AliMathBase.h"
37 #include "TTreeStream.h"
38 #include "AliTPCRawStreamFast.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
138 // AliTPCCalROC *GetCalRocNoise(Int_t sector);     - for the Noise values
139 //
140 // example for visualisation:
141 // if the file "PedestalData.root" was created using the above example one could do the following:
142 //
143 // TFile filePedestal("PedestalData.root")
144 // AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)filePedestal->Get("AliTPCCalibPedestal");
145 // ped->GetCalRocPedestal(0)->Draw("colz");
146 // ped->GetCalRocRMS(0)->Draw("colz");
147 //
148 // or use the AliTPCCalPad functionality:
149 // AliTPCCalPad padPedestal(ped->GetCalPadPedestal());
150 // AliTPCCalPad padNoise(ped->GetCalPadRMS());
151 // padPedestal->MakeHisto2D()->Draw("colz");  //Draw A-Side Pedestal Information
152 // padNoise->MakeHisto2D()->Draw("colz");  //Draw A-Side Noise Information
153 //
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 //
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() : /*FOLD00*/
203   TObject(),
204   fFirstTimeBin(60),
205   fLastTimeBin(1000),
206   fAdcMin(1),
207   fAdcMax(100),
208   fOldRCUformat(kTRUE),
209   fTimeAnalysis(kFALSE),
210   fROC(AliTPCROC::Instance()),
211   fMapping(NULL),
212   fCalRocArrayPedestal(72),
213   fCalRocArrayRMS(72),
214   fHistoPedestalArray(72),
215   fTimeSignal(NULL)
216 {
217   //
218   // default constructor
219   //
220 }
221
222
223 //_____________________________________________________________________
224 AliTPCCalibPedestal::AliTPCCalibPedestal(const AliTPCCalibPedestal &ped) : /*FOLD00*/
225   TObject(ped),
226   fFirstTimeBin(ped.GetFirstTimeBin()),
227   fLastTimeBin(ped.GetLastTimeBin()),
228   fAdcMin(ped.GetAdcMin()),
229   fAdcMax(ped.GetAdcMax()),
230   fOldRCUformat(ped.fOldRCUformat),
231   fTimeAnalysis(ped.fTimeAnalysis),
232   fROC(AliTPCROC::Instance()),
233   fMapping(NULL),
234   fCalRocArrayPedestal(72),
235   fCalRocArrayRMS(72),
236   fHistoPedestalArray(72),
237   fTimeSignal(ped.fTimeSignal)
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
258
259 //_____________________________________________________________________
260 AliTPCCalibPedestal& AliTPCCalibPedestal::operator = (const  AliTPCCalibPedestal &source)
261 {
262   //
263   // assignment operator
264   //
265   if (&source == this) return *this;
266   new (this) AliTPCCalibPedestal(source);
267
268   return *this;
269 }
270
271
272 //_____________________________________________________________________
273 AliTPCCalibPedestal::~AliTPCCalibPedestal() /*FOLD00*/
274 {
275   //
276   // destructor
277   //
278
279   fCalRocArrayPedestal.Delete();
280   fCalRocArrayRMS.Delete();
281   fHistoPedestalArray.Delete();
282
283   if ( fTimeSignal ) {
284     for (Int_t i = 0; i < 159; i++) {
285       delete [] fTimeSignal[i];
286       fTimeSignal[i] = 0;
287     }
288     delete [] fTimeSignal;
289     fTimeSignal = 0;
290   }
291
292   // do not delete fMapping, because we do not own it.
293
294 }
295
296
297 //_____________________________________________________________________
298 void AliTPCCalibPedestal::SetTimeAnalysis(Bool_t time)
299 {
300   //
301   // Use time dependent analysis: Pedestals are analysed as a function
302   // of the drift time. There is one mean value generated for each time
303   // bin and each channel. It can be used as reference data and for
304   // configuration of the ALTRO pattern memory for baseline subtraction.
305   //
306   // ATTENTION: Use only on LDC in TPCPEDESTALda! On a LDC we get data
307   // only from one sector. For the full TPC we would need a lot of
308   // memory (36*159*140*1024*4bytes = 3.3GB)!
309   //
310
311   fTimeAnalysis = time;
312
313   if ( !fTimeAnalysis ) return;
314
315   // prepare array for one sector (159*140*1024*4bytes = 92MB):
316   fTimeSignal = new TArrayF*[159];
317   for (Int_t i = 0; i < 159; i++) {  // padrows
318     fTimeSignal[i] = new TArrayF[140];
319     for (Int_t j = 0; j < 140; j++) {  // pads per row
320       fTimeSignal[i][j].Set(1024);
321       for (Int_t k = 0; k < 1024; k++) {  // time bins per pad
322         fTimeSignal[i][j].AddAt(0., k);
323       }
324     }      
325   }
326 }
327
328
329 //_____________________________________________________________________
330 Int_t AliTPCCalibPedestal::Update(const Int_t icsector, /*FOLD00*/
331                                   const Int_t icRow,
332                                   const Int_t icPad,
333                                   const Int_t icTimeBin,
334                                   const Float_t csignal)
335 {
336   //
337   // Signal filling method
338   //
339
340   // Time dependent pedestals
341   if ( fTimeAnalysis ) {
342     if ( icsector < 36 ) // IROC
343       fTimeSignal[icRow][icPad].AddAt(fTimeSignal[icRow][icPad].At(icTimeBin)+csignal, icTimeBin);
344     else 
345       fTimeSignal[icRow+63][icPad].AddAt(fTimeSignal[icRow+63][icPad].At(icTimeBin)+csignal, icTimeBin);
346   }
347   //return if we are out of the specified time bin or adc range
348   if ( (icTimeBin>fLastTimeBin) || (icTimeBin<fFirstTimeBin) ) return 0;
349   if ( ((Int_t)csignal>fAdcMax) || ((Int_t)csignal<fAdcMin)  ) return 0;
350
351   Int_t iChannel  = fROC->GetRowIndexes(icsector)[icRow]+icPad; //  global pad position in sector
352
353   // fast filling method
354   // Attention: the entry counter of the histogram is not increased
355   //            this means that e.g. the colz draw option gives an empty plot
356   Int_t bin = (iChannel+1)*(fAdcMax-fAdcMin+2)+((Int_t)csignal-fAdcMin+1);
357
358   GetHistoPedestal(icsector,kTRUE)->GetArray()[bin]++;
359
360   return 0;
361 }
362
363
364 //_____________________________________________________________________
365 Bool_t AliTPCCalibPedestal::ProcessEventFast(AliTPCRawStreamFast *rawStreamFast)
366 {
367   //
368   // Event Processing loop - AliTPCRawStream
369   //
370   Bool_t withInput = kFALSE;
371
372   while ( rawStreamFast->NextDDL() ){
373       while ( rawStreamFast->NextChannel() ){
374           Int_t isector  = rawStreamFast->GetSector();                       //  current sector
375           Int_t iRow     = rawStreamFast->GetRow();                          //  current row
376           Int_t iPad     = rawStreamFast->GetPad();                          //  current pad
377           Int_t startTbin = (Int_t)rawStreamFast->GetStartTimeBin();
378           Int_t endTbin = (Int_t)rawStreamFast->GetEndTimeBin();
379
380           while ( rawStreamFast->NextBunch() ){
381               for (Int_t iTimeBin = startTbin; iTimeBin < endTbin; iTimeBin++){
382                   Float_t signal=(Float_t)rawStreamFast->GetSignals()[iTimeBin-startTbin];
383                   Update(isector,iRow,iPad,iTimeBin+1,signal);
384                   withInput = kTRUE;
385               }
386           }
387       }
388   }
389
390   return withInput;
391 }
392 //_____________________________________________________________________
393 Bool_t AliTPCCalibPedestal::ProcessEventFast(AliRawReader *rawReader)
394 {
395   //
396   //  Event processing loop - AliRawReader
397   //
398   AliTPCRawStreamFast *rawStreamFast = new AliTPCRawStreamFast(rawReader, (AliAltroMapping**)fMapping);
399   Bool_t res=ProcessEventFast(rawStreamFast);
400   delete rawStreamFast;
401   return res;
402 }
403
404 //_____________________________________________________________________
405 Bool_t AliTPCCalibPedestal::ProcessEvent(AliTPCRawStream *rawStream)
406 {
407   //
408   // Event Processing loop - AliTPCRawStream
409   //
410
411   rawStream->SetOldRCUFormat(fOldRCUformat);
412
413   Bool_t withInput = kFALSE;
414
415   while (rawStream->Next()) {
416
417     Int_t iSector  = rawStream->GetSector();      //  current ROC
418     Int_t iRow     = rawStream->GetRow();         //  current row
419     Int_t iPad     = rawStream->GetPad();         //  current pad
420     Int_t iTimeBin = rawStream->GetTime();        //  current time bin
421     Float_t signal = rawStream->GetSignal();      //  current ADC signal
422     
423     Update(iSector,iRow,iPad,iTimeBin,signal);
424     withInput = kTRUE;
425   }
426
427   return withInput;
428 }
429
430
431 //_____________________________________________________________________
432 Bool_t AliTPCCalibPedestal::ProcessEvent(AliRawReader *rawReader)
433 {
434   //
435   //  Event processing loop - AliRawReader
436   //
437
438   // if fMapping is NULL the rawstream will crate its own mapping
439   AliTPCRawStream rawStream(rawReader, (AliAltroMapping**)fMapping);
440   rawReader->Select("TPC");
441   return ProcessEvent(&rawStream);
442 }
443
444
445 //_____________________________________________________________________
446 Bool_t AliTPCCalibPedestal::ProcessEvent(eventHeaderStruct *event)
447 {
448   //
449   //  process date event
450   //
451
452   AliRawReader *rawReader = new AliRawReaderDate((void*)event);
453   Bool_t result=ProcessEvent(rawReader);
454   delete rawReader;
455   return result;
456 }
457
458
459 //_____________________________________________________________________
460 Bool_t AliTPCCalibPedestal::TestEvent() /*FOLD00*/
461 {
462   //
463   //  Test event loop
464   // fill one oroc and one iroc with random gaus
465   //
466
467     gRandom->SetSeed(0);
468
469     for (UInt_t iSec=0; iSec<72; ++iSec){
470         if (iSec%36>0) continue;
471         for (UInt_t iRow=0; iRow < fROC->GetNRows(iSec); ++iRow){
472             for (UInt_t iPad=0; iPad < fROC->GetNPads(iSec,iRow); ++iPad){
473                 for (UInt_t iTimeBin=0; iTimeBin<1024; ++iTimeBin){
474                     Float_t signal=(Int_t)(iRow+3+gRandom->Gaus(0,.7));
475                     if ( signal>0 )Update(iSec,iRow,iPad,iTimeBin,signal);
476                 }
477             }
478         }
479     }
480     return kTRUE;
481 }
482
483
484 //_____________________________________________________________________
485 TH2F* AliTPCCalibPedestal::GetHisto(Int_t sector, TObjArray *arr, /*FOLD00*/
486                                     Int_t nbinsY, Float_t ymin, Float_t ymax,
487                                     Char_t *type, Bool_t force)
488 {
489     //
490     // return pointer to Q histogram
491     // if force is true create a new histogram if it doesn't exist allready
492     //
493     if ( !force || arr->UncheckedAt(sector) )
494         return (TH2F*)arr->UncheckedAt(sector);
495
496     // if we are forced and histogram doesn't yes exist create it
497     Char_t name[255], title[255];
498
499     sprintf(name,"hCalib%s%.2d",type,sector);
500     sprintf(title,"%s calibration histogram sector %.2d;ADC channel;Channel (pad)",type,sector);
501
502     // new histogram with Q calib information. One value for each pad!
503     TH2F* hist = new TH2F(name,title,
504                           nbinsY, ymin, ymax,
505                           fROC->GetNChannels(sector),0,fROC->GetNChannels(sector)
506                          );
507     hist->SetDirectory(0);
508     arr->AddAt(hist,sector);
509     return hist;
510 }
511
512
513 //_____________________________________________________________________
514 TH2F* AliTPCCalibPedestal::GetHistoPedestal(Int_t sector, Bool_t force) /*FOLD00*/
515 {
516     //
517     // return pointer to T0 histogram
518     // if force is true create a new histogram if it doesn't exist allready
519     //
520     TObjArray *arr = &fHistoPedestalArray;
521     return GetHisto(sector, arr, fAdcMax-fAdcMin, fAdcMin, fAdcMax, "Pedestal", force);
522 }
523
524
525 //_____________________________________________________________________
526 AliTPCCalROC* AliTPCCalibPedestal::GetCalRoc(Int_t sector, TObjArray* arr, Bool_t force) /*FOLD00*/
527 {
528     //
529     // return pointer to ROC Calibration
530     // if force is true create a new histogram if it doesn't exist allready
531     //
532     if ( !force || arr->UncheckedAt(sector) )
533         return (AliTPCCalROC*)arr->UncheckedAt(sector);
534
535     // if we are forced and the histogram doesn't yet exist create it
536
537     // new AliTPCCalROC for T0 information. One value for each pad!
538     AliTPCCalROC *croc = new AliTPCCalROC(sector);
539     arr->AddAt(croc,sector);
540     return croc;
541 }
542
543
544 //_____________________________________________________________________
545 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocPedestal(Int_t sector, Bool_t force) /*FOLD00*/
546 {
547     //
548     // return pointer to Carge ROC Calibration
549     // if force is true create a new histogram if it doesn't exist allready
550     //
551     TObjArray *arr = &fCalRocArrayPedestal;
552     return GetCalRoc(sector, arr, force);
553 }
554
555
556 //_____________________________________________________________________
557 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocRMS(Int_t sector, Bool_t force) /*FOLD00*/
558 {
559     //
560     // return pointer to signal width ROC Calibration
561     // if force is true create a new histogram if it doesn't exist allready
562     //
563     TObjArray *arr = &fCalRocArrayRMS;
564     return GetCalRoc(sector, arr, force);
565 }
566
567
568 //_____________________________________________________________________
569 void AliTPCCalibPedestal::Merge(AliTPCCalibPedestal *ped)
570 {
571   //
572   //  Merge reference histograms of sig to the current AliTPCCalibSignal
573   //
574
575   // merge histograms
576   for (Int_t iSec=0; iSec<72; ++iSec){
577     TH2F *hRefPedMerge   = ped->GetHistoPedestal(iSec);
578
579     if ( hRefPedMerge ){
580       TDirectory *dir = hRefPedMerge->GetDirectory(); hRefPedMerge->SetDirectory(0);
581       TH2F *hRefPed   = GetHistoPedestal(iSec);
582       if ( hRefPed ) hRefPed->Add(hRefPedMerge);
583       else {
584         TH2F *hist = new TH2F(*hRefPedMerge);
585         hist->SetDirectory(0);
586         fHistoPedestalArray.AddAt(hist, iSec);
587       }
588       hRefPedMerge->SetDirectory(dir);
589     }
590   }
591
592   // merge array
593   // ...
594
595 }
596
597
598 //_____________________________________________________________________
599 void AliTPCCalibPedestal::Analyse() /*FOLD00*/
600 {
601   //
602   //  Calculate calibration constants
603   //
604
605   Int_t nbinsAdc = fAdcMax-fAdcMin;
606
607   TVectorD param(3);
608   TMatrixD dummy(3,3);
609
610   Float_t *array_hP=0;  
611
612   for (Int_t iSec=0; iSec<72; ++iSec){
613     TH2F *hP = GetHistoPedestal(iSec);
614     if ( !hP ) continue;
615
616     AliTPCCalROC *rocPedestal = GetCalRocPedestal(iSec,kTRUE);
617     AliTPCCalROC *rocRMS      = GetCalRocRMS(iSec,kTRUE);
618
619     array_hP = hP->GetArray();
620     UInt_t nChannels = fROC->GetNChannels(iSec);
621
622     for (UInt_t iChannel=0; iChannel<nChannels; ++iChannel){
623       Int_t offset = (nbinsAdc+2)*(iChannel+1)+1;
624       Double_t ret = AliMathBase::FitGaus(array_hP+offset,nbinsAdc,fAdcMin,fAdcMax,&param,&dummy);
625       // if the fitting failed set noise and pedestal to 0
626       if ( ret == -4 ) {
627         param[1]=0;
628         param[2]=0;
629       }
630       rocPedestal->SetValue(iChannel,param[1]);
631       rocRMS->SetValue(iChannel,param[2]);
632     }
633   }
634 }
635
636
637 //_____________________________________________________________________
638 void AliTPCCalibPedestal::AnalyseTime(Int_t nevents)
639 {
640   //
641   // Calculate for each channel and time bin the mean pedestal. This
642   // is used on LDC by TPCPEDESTALda to generate data used for configuration
643   // of the pattern memory for baseline subtraction in the ALTROs.
644   //
645
646   if ( nevents <= 0 ) return;
647   if ( fTimeAnalysis ) {
648     for (Int_t i = 0; i < 159; i++) {  // padrows
649       for (Int_t j = 0; j < 140; j++) {  // pads per row
650         for (Int_t k = 0; k < 1024; k++) {  // time bins per pad
651           fTimeSignal[i][j].AddAt(fTimeSignal[i][j].At(k)/(Float_t)nevents, k);
652         }
653       }
654     }
655   }
656 }
657
658
659 //_____________________________________________________________________
660 void AliTPCCalibPedestal::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append) /*FOLD00*/
661 {
662   //
663   //  Write class to file
664   //
665
666   TString sDir(dir);
667   TString option;
668
669   if ( append )
670     option = "update";
671   else
672     option = "recreate";
673
674   TDirectory *backup = gDirectory;
675   TFile f(filename,option.Data());
676   f.cd();
677   if ( !sDir.IsNull() ){
678     f.mkdir(sDir.Data());
679     f.cd(sDir);
680   }
681   this->Write();
682   f.Close();
683
684   if ( backup ) backup->cd();
685 }