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