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