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