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