1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
27 #include <TDirectory.h>
31 #include "AliRawReader.h"
32 #include "AliRawReaderRoot.h"
33 #include "AliRawReaderDate.h"
34 #include "AliTPCRawStream.h"
35 #include "AliTPCCalROC.h"
36 #include "AliTPCROC.h"
37 #include "AliMathBase.h"
38 #include "TTreeStream.h"
39 #include "AliTPCRawStreamFast.h"
45 #include "AliTPCCalibPedestal.h"
48 ///////////////////////////////////////////////////////////////////////////////////////
49 // Implementation of the TPC pedestal and noise calibration
51 // Origin: Jens Wiechula, Marian Ivanov J.Wiechula@gsi.de, Marian.Ivanov@cern.ch
54 // *************************************************************************************
55 // * Class Description *
56 // *************************************************************************************
60 // Raw pedestal data is processed by calling one of the ProcessEvent(...) functions
61 // (see below). These in the end call the Update(...) function, where the data is filled
64 // For each ROC one TH2F histo (ROC channel vs. ADC channel) is created when
65 // it is filled for the first time (GetHistoPedestal(ROC,kTRUE)). All histos are stored in the
66 // TObjArray fHistoPedestalArray.
68 // For a fast filling of the histogram the corresponding bin number of the channel and ADC channel
69 // is computed by hand and the histogram array is accessed directly via its pointer.
70 // ATTENTION: Doing so the the entry counter of the histogram is not increased
71 // this means that e.g. the colz draw option gives an empty plot unless
72 // calling 'histo->SetEntries(1)' before drawing.
74 // After accumulating the desired statistics the Analyse() function has to be called.
75 // Whithin this function the pedestal and noise values are calculated for each pad, using
76 // the fast gaus fit function AliMathBase::FitGaus(...), and the calibration
77 // storage classes (AliTPCCalROC) are filled for each ROC.
78 // The calibration information is stored in the TObjArrays fCalRocArrayPedestal and fCalRocArrayRMS;
82 // User interface for filling data:
83 // --------------------------------
85 // To Fill information one of the following functions can be used:
87 // Bool_t ProcessEvent(eventHeaderStruct *event);
88 // - process Date event
89 // - use AliTPCRawReaderDate and call ProcessEvent(AliRawReader *rawReader)
91 // Bool_t ProcessEvent(AliRawReader *rawReader);
92 // - process AliRawReader event
93 // - use AliTPCRawStream to loop over data and call ProcessEvent(AliTPCRawStream *rawStream)
95 // Bool_t ProcessEvent(AliTPCRawStream *rawStream);
96 // - process event from AliTPCRawStream
97 // - call Update function for signal filling
99 // Int_t Update(const Int_t isector, const Int_t iRow, const Int_t
100 // iPad, const Int_t iTimeBin, const Float_t signal);
101 // - directly fill signal information (sector, row, pad, time bin, pad)
102 // to the reference histograms
104 // It is also possible to merge two independently taken calibrations using the function
106 // void Merge(AliTPCCalibPedestal *ped)
107 // - copy histograms in 'ped' if the do not exist in this instance
108 // - Add histograms in 'ped' to the histograms in this instance if the allready exist
109 // - After merging call Analyse again!
113 // -- example: filling data using root raw data:
114 // void fillPedestal(Char_t *filename)
116 // rawReader = new AliRawReaderRoot(fileName);
117 // if ( !rawReader ) return;
118 // AliTPCCalibPedestal *calib = new AliTPCCalibPedestal;
119 // while (rawReader->NextEvent()){
120 // calib->ProcessEvent(rawReader);
123 // calib->DumpToFile("PedestalData.root");
129 // What kind of information is stored and how to retrieve them:
130 // ------------------------------------------------------------
132 // - Accessing the 'Reference Histograms' (pedestal distribution histograms):
134 // TH2F *GetHistoPedestal(Int_t sector);
136 // - Accessing the calibration storage objects:
138 // AliTPCCalROC *GetCalRocPedestal(Int_t sector); - for the pedestal values, mean from gaus fit
139 // AliTPCCalROC *GetCalRocSigma(Int_t sector); - for the Noise values, sigma from guas fit
140 // AliTPCCalROC *GetCalRocMean(Int_t sector); - for the pedestal values, truncated mean
141 // AliTPCCalROC *GetCalRocRMS(Int_t sector); - for the Noise values, rms from truncated mean
143 // example for visualisation:
144 // if the file "PedestalData.root" was created using the above example one could do the following:
146 // TFile filePedestal("PedestalData.root")
147 // AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)filePedestal->Get("AliTPCCalibPedestal");
148 // ped->GetCalRocPedestal(0)->Draw("colz");
149 // ped->GetCalRocRMS(0)->Draw("colz");
151 // or use the AliTPCCalPad functionality:
152 // AliTPCCalPad padPedestal(ped->GetCalPadPedestal());
153 // AliTPCCalPad padNoise(ped->GetCalPadRMS());
154 // padPedestal->MakeHisto2D()->Draw("colz"); //Draw A-Side Pedestal Information
155 // padNoise->MakeHisto2D()->Draw("colz"); //Draw A-Side Noise Information
158 example: fill pedestal with gausschen noise
159 AliTPCCalibPedestal ped;
163 TCanvas* c1 = new TCanvas;
166 ped.GetHistoPedestal(0)->SetEntries(1); //needed in order for colz to work, reason: fast filling does not increase the entries counter
167 ped.GetHistoPedestal(0)->Draw("colz");
169 ped.GetHistoPedestal(36)->SetEntries(1); //needed in order for colz to work, reason: fast filling does not increase the entries counter
170 ped.GetHistoPedestal(36)->Draw("colz");
171 TCanvas* c2 = new TCanvas;
174 ped.GetCalRocPedestal(0)->Draw("colz");
176 ped.GetCalRocRMS(0)->Draw("colz");
178 ped.GetCalRocPedestal(36)->Draw("colz");
180 ped.GetCalRocRMS(36)->Draw("colz");
183 // Time dependent pedestals:
185 // If wished there is the possibility to calculate for each channel and time bin
186 // the mean pedestal [pedestals(t)]. This is done by
188 // 1) setting SetTimeAnalysis(kTRUE),
189 // 2) processing the data by looping over the events using ProcessEvent(..)
190 // 3) calling the Analyse() and AnalyseTime(nevents) functions (providing nevents)
191 // 4) getting the pedestals(t) using TArrayF **timePed = calibPedestal.GetTimePedestals();
192 // 5) looking at values using timePed[row][pad].At(timebin)
194 // This functionality is intended to be used on an LDC bu the detector algorithm
195 // (TPCPEDESTALda) to generate a data set used for configuration of the pattern
196 // memory for baseline subtraction in the ALTROs. Later the information should also
197 // be stored as reference data.
201 ClassImp(AliTPCCalibPedestal)
203 AliTPCCalibPedestal::AliTPCCalibPedestal() :
204 AliTPCCalibRawBase(),
209 fTimeAnalysis(kFALSE),
210 fCalRocArrayPedestal(72),
211 fCalRocArraySigma(72),
212 fHistoPedestalArray(72),
214 fCalRocArrayMean(72),
218 // default constructor
220 SetNameTitle("AliTPCCalibPedestal","AliTPCCalibPedestal");
226 //_____________________________________________________________________
227 AliTPCCalibPedestal::AliTPCCalibPedestal(const AliTPCCalibPedestal &ped) :
228 AliTPCCalibRawBase(ped),
229 fAdcMin(ped.GetAdcMin()),
230 fAdcMax(ped.GetAdcMax()),
231 fAnaMeanDown(ped.fAnaMeanDown),
232 fAnaMeanUp(ped.fAnaMeanUp),
233 fTimeAnalysis(ped.fTimeAnalysis),
234 fCalRocArrayPedestal(72),
235 fCalRocArraySigma(72),
236 fHistoPedestalArray(72),
237 fTimeSignal(ped.fTimeSignal),
238 fCalRocArrayMean(72),
244 for (Int_t iSec = 0; iSec < 72; ++iSec){
245 const AliTPCCalROC *calPed = (AliTPCCalROC*)ped.fCalRocArrayPedestal.UncheckedAt(iSec);
246 const AliTPCCalROC *calRMS = (AliTPCCalROC*)ped.fCalRocArrayRMS.UncheckedAt(iSec);
247 const TH2F *hPed = (TH2F*)ped.fHistoPedestalArray.UncheckedAt(iSec);
249 if ( calPed != 0x0 ) fCalRocArrayPedestal.AddAt(new AliTPCCalROC(*calPed), iSec);
250 if ( calRMS != 0x0 ) fCalRocArrayRMS.AddAt(new AliTPCCalROC(*calRMS), iSec);
253 TH2F *hNew = new TH2F(*hPed);
254 hNew->SetDirectory(0);
255 fHistoPedestalArray.AddAt(hNew,iSec);
259 AliTPCCalibPedestal::AliTPCCalibPedestal(const TMap *config):
260 AliTPCCalibRawBase(),
265 fTimeAnalysis(kFALSE),
266 fCalRocArrayPedestal(72),
267 fCalRocArraySigma(72),
268 fHistoPedestalArray(72),
270 fCalRocArrayMean(72),
274 // This constructor uses a TMap for setting some parametes
276 SetNameTitle("AliTPCCalibPedestal","AliTPCCalibPedestal");
279 if (config->GetValue("FirstTimeBin")) fFirstTimeBin = ((TObjString*)config->GetValue("FirstTimeBin"))->GetString().Atoi();
280 if (config->GetValue("LastTimeBin")) fLastTimeBin = ((TObjString*)config->GetValue("LastTimeBin"))->GetString().Atoi();
281 if (config->GetValue("AdcMin")) fAdcMin = ((TObjString*)config->GetValue("AdcMin"))->GetString().Atoi();
282 if (config->GetValue("AdcMax")) fAdcMax = ((TObjString*)config->GetValue("AdcMax"))->GetString().Atoi();
283 if (config->GetValue("TimeAnalysis")) SetTimeAnalysis(((TObjString*)config->GetValue("TimeAnalysis"))->GetString().Atoi());
287 //_____________________________________________________________________
288 AliTPCCalibPedestal& AliTPCCalibPedestal::operator = (const AliTPCCalibPedestal &source)
291 // assignment operator
293 if (&source == this) return *this;
294 new (this) AliTPCCalibPedestal(source);
300 //_____________________________________________________________________
301 AliTPCCalibPedestal::~AliTPCCalibPedestal()
307 fCalRocArrayPedestal.Delete();
308 fCalRocArrayRMS.Delete();
309 fCalRocArraySigma.Delete();
310 fHistoPedestalArray.Delete();
313 for (Int_t i = 0; i < 159; i++) {
314 delete [] fTimeSignal[i];
317 delete [] fTimeSignal;
321 // do not delete fMapping, because we do not own it.
326 //_____________________________________________________________________
327 void AliTPCCalibPedestal::SetTimeAnalysis(Bool_t time)
330 // Use time dependent analysis: Pedestals are analysed as a function
331 // of the drift time. There is one mean value generated for each time
332 // bin and each channel. It can be used as reference data and for
333 // configuration of the ALTRO pattern memory for baseline subtraction.
335 // ATTENTION: Use only on LDC in TPCPEDESTALda! On a LDC we get data
336 // only from one sector. For the full TPC we would need a lot of
337 // memory (36*159*140*1024*4bytes = 3.3GB)!
340 fTimeAnalysis = time;
342 if ( !fTimeAnalysis ) return;
344 // prepare array for one sector (159*140*1024*4bytes = 92MB):
345 fTimeSignal = new TArrayF*[159];
346 for (Int_t i = 0; i < 159; i++) { // padrows
347 fTimeSignal[i] = new TArrayF[140];
348 for (Int_t j = 0; j < 140; j++) { // pads per row
349 fTimeSignal[i][j].Set(1024);
350 for (Int_t k = 0; k < 1024; k++) { // time bins per pad
351 fTimeSignal[i][j].AddAt(0., k);
358 //_____________________________________________________________________
359 Int_t AliTPCCalibPedestal::Update(const Int_t icsector,
362 const Int_t icTimeBin,
363 const Float_t csignal)
366 // Signal filling method
368 if (icRow<0) return 0;
369 if (icPad<0) return 0;
370 if (icTimeBin<0) return 0;
372 // Time dependent pedestals
373 if ( fTimeAnalysis ) {
374 if ( icsector < 36 ) // IROC
375 fTimeSignal[icRow][icPad].AddAt(fTimeSignal[icRow][icPad].At(icTimeBin)+csignal, icTimeBin);
377 fTimeSignal[icRow+63][icPad].AddAt(fTimeSignal[icRow+63][icPad].At(icTimeBin)+csignal, icTimeBin);
379 //return if we are out of the specified time bin or adc range
380 if ( (icTimeBin>fLastTimeBin) || (icTimeBin<fFirstTimeBin) ) return 0;
381 if ( ((Int_t)csignal>fAdcMax) || ((Int_t)csignal<fAdcMin) ) return 0;
383 Int_t iChannel = fROC->GetRowIndexes(icsector)[icRow]+icPad; // global pad position in sector
385 // fast filling method
386 // Attention: the entry counter of the histogram is not increased
387 // this means that e.g. the colz draw option gives an empty plot
388 Int_t bin = (iChannel+1)*(fAdcMax-fAdcMin+2)+((Int_t)csignal-fAdcMin+1);
390 GetHistoPedestal(icsector,kTRUE)->GetArray()[bin]++;
396 //_____________________________________________________________________
397 Bool_t AliTPCCalibPedestal::TestEvent()
401 // fill one oroc and one iroc with random gaus
406 for (UInt_t iSec=0; iSec<72; ++iSec){
407 if (iSec%36>0) continue;
408 for (UInt_t iRow=0; iRow < fROC->GetNRows(iSec); ++iRow){
409 for (UInt_t iPad=0; iPad < fROC->GetNPads(iSec,iRow); ++iPad){
410 for (UInt_t iTimeBin=0; iTimeBin<1024; ++iTimeBin){
411 Float_t signal=(Int_t)(iRow+3+gRandom->Gaus(0,.7));
412 if ( signal>0 )Update(iSec,iRow,iPad,iTimeBin,signal);
421 //_____________________________________________________________________
422 TH2F* AliTPCCalibPedestal::GetHisto(Int_t sector, TObjArray *arr,
423 Int_t nbinsY, Float_t ymin, Float_t ymax,
424 const Char_t *type, Bool_t force)
427 // return pointer to Q histogram
428 // if force is true create a new histogram if it doesn't exist allready
430 if ( !force || arr->UncheckedAt(sector) )
431 return (TH2F*)arr->UncheckedAt(sector);
433 // if we are forced and histogram doesn't yes exist create it
434 Char_t name[255], title[255];
436 sprintf(name,"hCalib%s%.2d",type,sector);
437 sprintf(title,"%s calibration histogram sector %.2d;ADC channel;Channel (pad)",type,sector);
439 // new histogram with Q calib information. One value for each pad!
440 TH2F* hist = new TH2F(name,title,
442 fROC->GetNChannels(sector),0,fROC->GetNChannels(sector)
444 hist->SetDirectory(0);
445 arr->AddAt(hist,sector);
450 //_____________________________________________________________________
451 TH2F* AliTPCCalibPedestal::GetHistoPedestal(Int_t sector, Bool_t force)
454 // return pointer to T0 histogram
455 // if force is true create a new histogram if it doesn't exist allready
457 TObjArray *arr = &fHistoPedestalArray;
458 return GetHisto(sector, arr, fAdcMax-fAdcMin, fAdcMin, fAdcMax, "Pedestal", force);
462 //_____________________________________________________________________
463 AliTPCCalROC* AliTPCCalibPedestal::GetCalRoc(Int_t sector, TObjArray* arr, Bool_t force)
466 // return pointer to ROC Calibration
467 // if force is true create a new histogram if it doesn't exist allready
469 if ( !force || arr->UncheckedAt(sector) )
470 return (AliTPCCalROC*)arr->UncheckedAt(sector);
472 // if we are forced and the histogram doesn't yet exist create it
474 // new AliTPCCalROC for T0 information. One value for each pad!
475 AliTPCCalROC *croc = new AliTPCCalROC(sector);
476 arr->AddAt(croc,sector);
481 //_____________________________________________________________________
482 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocPedestal(Int_t sector, Bool_t force)
485 // return pointer to ROC with Pedestal data
486 // if force is true create a new histogram if it doesn't exist allready
488 TObjArray *arr = &fCalRocArrayPedestal;
489 return GetCalRoc(sector, arr, force);
493 //_____________________________________________________________________
494 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocSigma(Int_t sector, Bool_t force)
497 // return pointer to ROC with signal witdth in sigma
498 // if force is true create a new histogram if it doesn't exist allready
500 TObjArray *arr = &fCalRocArraySigma;
501 return GetCalRoc(sector, arr, force);
503 //_____________________________________________________________________
504 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocMean(Int_t sector, Bool_t force)
507 // return pointer to ROC with signal mean information
508 // if force is true create a new histogram if it doesn't exist allready
510 TObjArray *arr = &fCalRocArrayMean;
511 return GetCalRoc(sector, arr, force);
514 //_____________________________________________________________________
515 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocRMS(Int_t sector, Bool_t force)
518 // return pointer to signal width ROC Calibration
519 // if force is true create a new histogram if it doesn't exist allready
521 TObjArray *arr = &fCalRocArrayRMS;
522 return GetCalRoc(sector, arr, force);
526 //_____________________________________________________________________
527 void AliTPCCalibPedestal::Merge(AliTPCCalibPedestal *ped)
530 // Merge reference histograms of sig to the current AliTPCCalibSignal
534 for (Int_t iSec=0; iSec<72; ++iSec){
535 TH2F *hRefPedMerge = ped->GetHistoPedestal(iSec);
538 TDirectory *dir = hRefPedMerge->GetDirectory(); hRefPedMerge->SetDirectory(0);
539 TH2F *hRefPed = GetHistoPedestal(iSec);
540 if ( hRefPed ) hRefPed->Add(hRefPedMerge);
542 TH2F *hist = new TH2F(*hRefPedMerge);
543 hist->SetDirectory(0);
544 fHistoPedestalArray.AddAt(hist, iSec);
546 hRefPedMerge->SetDirectory(dir);
556 //_____________________________________________________________________
557 void AliTPCCalibPedestal::Analyse()
560 // Calculate calibration constants
563 Int_t nbinsAdc = fAdcMax-fAdcMin;
568 TH1F *hChannel=new TH1F("hChannel","hChannel",nbinsAdc,fAdcMin,fAdcMax);
572 for (Int_t iSec=0; iSec<72; ++iSec){
573 TH2F *hP = GetHistoPedestal(iSec);
576 AliTPCCalROC *rocPedestal = GetCalRocPedestal(iSec,kTRUE);
577 AliTPCCalROC *rocSigma = GetCalRocSigma(iSec,kTRUE);
578 AliTPCCalROC *rocMean = GetCalRocMean(iSec,kTRUE);
579 AliTPCCalROC *rocRMS = GetCalRocRMS(iSec,kTRUE);
581 array_hP = hP->GetArray();
582 UInt_t nChannels = fROC->GetNChannels(iSec);
584 for (UInt_t iChannel=0; iChannel<nChannels; ++iChannel){
585 Int_t offset = (nbinsAdc+2)*(iChannel+1)+1;
586 //calculate mean and sigma using a gaus fit
588 AliMathBase::FitGaus(array_hP+offset,nbinsAdc,fAdcMin,fAdcMax,¶m,&dummy);
589 // if the fitting failed set noise and pedestal to 0
590 // is now done in AliMathBase::FitGaus !
591 // if ( ret == -4 ) {
595 if ( param[1]<fAdcMin || param[1]>fAdcMax ){
599 rocPedestal->SetValue(iChannel,param[1]);
600 rocSigma->SetValue(iChannel,param[2]);
601 //calculate mean and RMS using a truncated means
602 hChannel->Set(nbinsAdc+2,array_hP+offset-1);
603 hChannel->SetEntries(param[3]);
606 if ( param[3]>0 ) AliMathBase::TruncatedMean(hChannel,¶m,fAnaMeanDown,fAnaMeanUp);
607 rocMean->SetValue(iChannel,param[1]);
608 rocRMS->SetValue(iChannel,param[2]);
615 //_____________________________________________________________________
616 void AliTPCCalibPedestal::AnalyseTime(Int_t nevents)
619 // Calculate for each channel and time bin the mean pedestal. This
620 // is used on LDC by TPCPEDESTALda to generate data used for configuration
621 // of the pattern memory for baseline subtraction in the ALTROs.
624 if ( nevents <= 0 ) return;
625 if ( fTimeAnalysis ) {
626 for (Int_t i = 0; i < 159; i++) { // padrows
627 for (Int_t j = 0; j < 140; j++) { // pads per row
628 for (Int_t k = 0; k < 1024; k++) { // time bins per pad
629 fTimeSignal[i][j].AddAt(fTimeSignal[i][j].At(k)/(Float_t)nevents, k);