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>
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"
43 #include "AliTPCCalibPedestal.h"
46 ///////////////////////////////////////////////////////////////////////////////////////
47 // Implementation of the TPC pedestal and noise calibration
49 // Origin: Jens Wiechula, Marian Ivanov J.Wiechula@gsi.de, Marian.Ivanov@cern.ch
52 // *************************************************************************************
53 // * Class Description *
54 // *************************************************************************************
58 // Raw pedestal data is processed by calling one of the ProcessEvent(...) functions
59 // (see below). These in the end call the Update(...) function, where the data is filled
62 // For each ROC one TH2F histo (ROC channel vs. ADC channel) is created when
63 // it is filled for the first time (GetHistoPedestal(ROC,kTRUE)). All histos are stored in the
64 // TObjArray fHistoPedestalArray.
66 // For a fast filling of the histogram the corresponding bin number of the channel and ADC channel
67 // is computed by hand and the histogram array is accessed directly via its pointer.
68 // ATTENTION: Doing so the the entry counter of the histogram is not increased
69 // this means that e.g. the colz draw option gives an empty plot unless
70 // calling 'histo->SetEntries(1)' before drawing.
72 // After accumulating the desired statistics the Analyse() function has to be called.
73 // Whithin this function the pedestal and noise values are calculated for each pad, using
74 // the fast gaus fit function AliMathBase::FitGaus(...), and the calibration
75 // storage classes (AliTPCCalROC) are filled for each ROC.
76 // The calibration information is stored in the TObjArrays fCalRocArrayPedestal and fCalRocArrayRMS;
80 // User interface for filling data:
81 // --------------------------------
83 // To Fill information one of the following functions can be used:
85 // Bool_t ProcessEvent(eventHeaderStruct *event);
86 // - process Date event
87 // - use AliTPCRawReaderDate and call ProcessEvent(AliRawReader *rawReader)
89 // Bool_t ProcessEvent(AliRawReader *rawReader);
90 // - process AliRawReader event
91 // - use AliTPCRawStream to loop over data and call ProcessEvent(AliTPCRawStream *rawStream)
93 // Bool_t ProcessEvent(AliTPCRawStream *rawStream);
94 // - process event from AliTPCRawStream
95 // - call Update function for signal filling
97 // Int_t Update(const Int_t isector, const Int_t iRow, const Int_t
98 // iPad, const Int_t iTimeBin, const Float_t signal);
99 // - directly fill signal information (sector, row, pad, time bin, pad)
100 // to the reference histograms
102 // It is also possible to merge two independently taken calibrations using the function
104 // void Merge(AliTPCCalibPedestal *ped)
105 // - copy histograms in 'ped' if the do not exist in this instance
106 // - Add histograms in 'ped' to the histograms in this instance if the allready exist
107 // - After merging call Analyse again!
111 // -- example: filling data using root raw data:
112 // void fillPedestal(Char_t *filename)
114 // rawReader = new AliRawReaderRoot(fileName);
115 // if ( !rawReader ) return;
116 // AliTPCCalibPedestal *calib = new AliTPCCalibPedestal;
117 // while (rawReader->NextEvent()){
118 // calib->ProcessEvent(rawReader);
121 // calib->DumpToFile("PedestalData.root");
127 // What kind of information is stored and how to retrieve them:
128 // ------------------------------------------------------------
130 // - Accessing the 'Reference Histograms' (pedestal distribution histograms):
132 // TH2F *GetHistoPedestal(Int_t sector);
134 // - Accessing the calibration storage objects:
136 // AliTPCCalROC *GetCalRocPedestal(Int_t sector); - for the pedestal values
137 // AliTPCCalROC *GetCalRocNoise(Int_t sector); - for the Noise values
139 // example for visualisation:
140 // if the file "PedestalData.root" was created using the above example one could do the following:
142 // TFile filePedestal("PedestalData.root")
143 // AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)filePedestal->Get("AliTPCCalibPedestal");
144 // ped->GetCalRocPedestal(0)->Draw("colz");
145 // ped->GetCalRocRMS(0)->Draw("colz");
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
155 example: fill pedestal with gausschen noise
156 AliTPCCalibPedestal ped;
160 TCanvas* c1 = new TCanvas;
163 ped.GetHistoPedestal(0)->SetEntries(1); //needed in order for colz to work, reason: fast filling does not increase the entries counter
164 ped.GetHistoPedestal(0)->Draw("colz");
166 ped.GetHistoPedestal(36)->SetEntries(1); //needed in order for colz to work, reason: fast filling does not increase the entries counter
167 ped.GetHistoPedestal(36)->Draw("colz");
168 TCanvas* c2 = new TCanvas;
171 ped.GetCalRocPedestal(0)->Draw("colz");
173 ped.GetCalRocRMS(0)->Draw("colz");
175 ped.GetCalRocPedestal(36)->Draw("colz");
177 ped.GetCalRocRMS(36)->Draw("colz");
181 // Time dependent pedestals:
183 // If wished there is the possibility to calculate for each channel and time bin
184 // the mean pedestal [pedestals(t)]. This is done by
186 // 1) setting SetTimeAnalysis(kTRUE),
187 // 2) processing the data by looping over the events using ProcessEvent(..)
188 // 3) calling the Analyse() and AnalyseTime(nevents) functions (providing nevents)
189 // 4) getting the pedestals(t) using TArrayF **timePed = calibPedestal.GetTimePedestals();
190 // 5) looking at values using timePed[row][pad].At(timebin)
192 // This functionality is intended to be used on an LDC bu the detector algorithm
193 // (TPCPEDESTALda) to generate a data set used for configuration of the pattern
194 // memory for baseline subtraction in the ALTROs. Later the information should also
195 // be stored as reference data.
199 ClassImp(AliTPCCalibPedestal)
201 AliTPCCalibPedestal::AliTPCCalibPedestal() : /*FOLD00*/
207 fOldRCUformat(kTRUE),
208 fTimeAnalysis(kFALSE),
209 fROC(AliTPCROC::Instance()),
211 fCalRocArrayPedestal(72),
213 fHistoPedestalArray(72),
217 // default constructor
222 //_____________________________________________________________________
223 AliTPCCalibPedestal::AliTPCCalibPedestal(const AliTPCCalibPedestal &ped) : /*FOLD00*/
225 fFirstTimeBin(ped.GetFirstTimeBin()),
226 fLastTimeBin(ped.GetLastTimeBin()),
227 fAdcMin(ped.GetAdcMin()),
228 fAdcMax(ped.GetAdcMax()),
229 fOldRCUformat(ped.fOldRCUformat),
230 fTimeAnalysis(ped.fTimeAnalysis),
231 fROC(AliTPCROC::Instance()),
233 fCalRocArrayPedestal(72),
235 fHistoPedestalArray(72),
236 fTimeSignal(ped.fTimeSignal)
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);
246 if ( calPed != 0x0 ) fCalRocArrayPedestal.AddAt(new AliTPCCalROC(*calPed), iSec);
247 if ( calRMS != 0x0 ) fCalRocArrayRMS.AddAt(new AliTPCCalROC(*calRMS), iSec);
250 TH2F *hNew = new TH2F(*hPed);
251 hNew->SetDirectory(0);
252 fHistoPedestalArray.AddAt(hNew,iSec);
258 //_____________________________________________________________________
259 AliTPCCalibPedestal& AliTPCCalibPedestal::operator = (const AliTPCCalibPedestal &source)
262 // assignment operator
264 if (&source == this) return *this;
265 new (this) AliTPCCalibPedestal(source);
271 //_____________________________________________________________________
272 AliTPCCalibPedestal::~AliTPCCalibPedestal() /*FOLD00*/
278 fCalRocArrayPedestal.Delete();
279 fCalRocArrayRMS.Delete();
280 fHistoPedestalArray.Delete();
283 for (Int_t i = 0; i < 159; i++) {
284 delete [] fTimeSignal[i];
287 delete [] fTimeSignal;
291 // do not delete fMapping, because we do not own it.
296 //_____________________________________________________________________
297 void AliTPCCalibPedestal::SetTimeAnalysis(Bool_t time)
300 // Use time dependent analysis: Pedestals are analysed as a function
301 // of the drift time. There is one mean value generated for each time
302 // bin and each channel. It can be used as reference data and for
303 // configuration of the ALTRO pattern memory for baseline subtraction.
305 // ATTENTION: Use only on LDC in TPCPEDESTALda! On a LDC we get data
306 // only from one sector. For the full TPC we would need a lot of
307 // memory (36*159*140*1024*4bytes = 3.3GB)!
310 fTimeAnalysis = time;
312 if ( !fTimeAnalysis ) return;
314 // prepare array for one sector (159*140*1024*4bytes = 92MB):
315 fTimeSignal = new TArrayF*[159];
316 for (Int_t i = 0; i < 159; i++) { // padrows
317 fTimeSignal[i] = new TArrayF[140];
318 for (Int_t j = 0; j < 140; j++) { // pads per row
319 fTimeSignal[i][j].Set(1024);
320 for (Int_t k = 0; k < 1024; k++) { // time bins per pad
321 fTimeSignal[i][j].AddAt(0., k);
328 //_____________________________________________________________________
329 Int_t AliTPCCalibPedestal::Update(const Int_t icsector, /*FOLD00*/
332 const Int_t icTimeBin,
333 const Float_t csignal)
336 // Signal filling method
339 // Time dependent pedestals
340 if ( fTimeAnalysis ) {
341 if ( icsector < 36 ) // IROC
342 fTimeSignal[icRow][icPad].AddAt(fTimeSignal[icRow][icPad].At(icTimeBin)+csignal, icTimeBin);
344 fTimeSignal[icRow+63][icPad].AddAt(fTimeSignal[icRow+63][icPad].At(icTimeBin)+csignal, icTimeBin);
346 //return if we are out of the specified time bin or adc range
347 if ( (icTimeBin>fLastTimeBin) || (icTimeBin<fFirstTimeBin) ) return 0;
348 if ( ((Int_t)csignal>fAdcMax) || ((Int_t)csignal<fAdcMin) ) return 0;
350 Int_t iChannel = fROC->GetRowIndexes(icsector)[icRow]+icPad; // global pad position in sector
352 // fast filling method
353 // Attention: the entry counter of the histogram is not increased
354 // this means that e.g. the colz draw option gives an empty plot
355 Int_t bin = (iChannel+1)*(fAdcMax-fAdcMin+2)+((Int_t)csignal-fAdcMin+1);
357 GetHistoPedestal(icsector,kTRUE)->GetArray()[bin]++;
363 //_____________________________________________________________________
364 Bool_t AliTPCCalibPedestal::ProcessEvent(AliTPCRawStream *rawStream)
367 // Event Processing loop - AliTPCRawStream
370 rawStream->SetOldRCUFormat(fOldRCUformat);
372 Bool_t withInput = kFALSE;
374 while (rawStream->Next()) {
376 Int_t iSector = rawStream->GetSector(); // current ROC
377 Int_t iRow = rawStream->GetRow(); // current row
378 Int_t iPad = rawStream->GetPad(); // current pad
379 Int_t iTimeBin = rawStream->GetTime(); // current time bin
380 Float_t signal = rawStream->GetSignal(); // current ADC signal
382 Update(iSector,iRow,iPad,iTimeBin,signal);
390 //_____________________________________________________________________
391 Bool_t AliTPCCalibPedestal::ProcessEvent(AliRawReader *rawReader)
394 // Event processing loop - AliRawReader
397 // if fMapping is NULL the rawstream will crate its own mapping
398 AliTPCRawStream rawStream(rawReader, (AliAltroMapping**)fMapping);
399 rawReader->Select("TPC");
400 return ProcessEvent(&rawStream);
404 //_____________________________________________________________________
405 Bool_t AliTPCCalibPedestal::ProcessEvent(eventHeaderStruct *event)
408 // process date event
411 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
412 Bool_t result=ProcessEvent(rawReader);
418 //_____________________________________________________________________
419 Bool_t AliTPCCalibPedestal::TestEvent() /*FOLD00*/
423 // fill one oroc and one iroc with random gaus
428 for (UInt_t iSec=0; iSec<72; ++iSec){
429 if (iSec%36>0) continue;
430 for (UInt_t iRow=0; iRow < fROC->GetNRows(iSec); ++iRow){
431 for (UInt_t iPad=0; iPad < fROC->GetNPads(iSec,iRow); ++iPad){
432 for (UInt_t iTimeBin=0; iTimeBin<1024; ++iTimeBin){
433 Float_t signal=(Int_t)(iRow+3+gRandom->Gaus(0,.7));
434 if ( signal>0 )Update(iSec,iRow,iPad,iTimeBin,signal);
443 //_____________________________________________________________________
444 TH2F* AliTPCCalibPedestal::GetHisto(Int_t sector, TObjArray *arr, /*FOLD00*/
445 Int_t nbinsY, Float_t ymin, Float_t ymax,
446 Char_t *type, Bool_t force)
449 // return pointer to Q histogram
450 // if force is true create a new histogram if it doesn't exist allready
452 if ( !force || arr->UncheckedAt(sector) )
453 return (TH2F*)arr->UncheckedAt(sector);
455 // if we are forced and histogram doesn't yes exist create it
456 Char_t name[255], title[255];
458 sprintf(name,"hCalib%s%.2d",type,sector);
459 sprintf(title,"%s calibration histogram sector %.2d;ADC channel;Channel (pad)",type,sector);
461 // new histogram with Q calib information. One value for each pad!
462 TH2F* hist = new TH2F(name,title,
464 fROC->GetNChannels(sector),0,fROC->GetNChannels(sector)
466 hist->SetDirectory(0);
467 arr->AddAt(hist,sector);
472 //_____________________________________________________________________
473 TH2F* AliTPCCalibPedestal::GetHistoPedestal(Int_t sector, Bool_t force) /*FOLD00*/
476 // return pointer to T0 histogram
477 // if force is true create a new histogram if it doesn't exist allready
479 TObjArray *arr = &fHistoPedestalArray;
480 return GetHisto(sector, arr, fAdcMax-fAdcMin, fAdcMin, fAdcMax, "Pedestal", force);
484 //_____________________________________________________________________
485 AliTPCCalROC* AliTPCCalibPedestal::GetCalRoc(Int_t sector, TObjArray* arr, Bool_t force) /*FOLD00*/
488 // return pointer to ROC Calibration
489 // if force is true create a new histogram if it doesn't exist allready
491 if ( !force || arr->UncheckedAt(sector) )
492 return (AliTPCCalROC*)arr->UncheckedAt(sector);
494 // if we are forced and the histogram doesn't yet exist create it
496 // new AliTPCCalROC for T0 information. One value for each pad!
497 AliTPCCalROC *croc = new AliTPCCalROC(sector);
498 arr->AddAt(croc,sector);
503 //_____________________________________________________________________
504 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocPedestal(Int_t sector, Bool_t force) /*FOLD00*/
507 // return pointer to Carge ROC Calibration
508 // if force is true create a new histogram if it doesn't exist allready
510 TObjArray *arr = &fCalRocArrayPedestal;
511 return GetCalRoc(sector, arr, force);
515 //_____________________________________________________________________
516 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocRMS(Int_t sector, Bool_t force) /*FOLD00*/
519 // return pointer to signal width ROC Calibration
520 // if force is true create a new histogram if it doesn't exist allready
522 TObjArray *arr = &fCalRocArrayRMS;
523 return GetCalRoc(sector, arr, force);
527 //_____________________________________________________________________
528 void AliTPCCalibPedestal::Merge(AliTPCCalibPedestal *ped)
531 // Merge reference histograms of sig to the current AliTPCCalibSignal
535 for (Int_t iSec=0; iSec<72; ++iSec){
536 TH2F *hRefPedMerge = ped->GetHistoPedestal(iSec);
539 TDirectory *dir = hRefPedMerge->GetDirectory(); hRefPedMerge->SetDirectory(0);
540 TH2F *hRefPed = GetHistoPedestal(iSec);
541 if ( hRefPed ) hRefPed->Add(hRefPedMerge);
543 TH2F *hist = new TH2F(*hRefPedMerge);
544 hist->SetDirectory(0);
545 fHistoPedestalArray.AddAt(hist, iSec);
547 hRefPedMerge->SetDirectory(dir);
557 //_____________________________________________________________________
558 void AliTPCCalibPedestal::Analyse() /*FOLD00*/
561 // Calculate calibration constants
564 Int_t nbinsAdc = fAdcMax-fAdcMin;
571 for (Int_t iSec=0; iSec<72; ++iSec){
572 TH2F *hP = GetHistoPedestal(iSec);
575 AliTPCCalROC *rocPedestal = GetCalRocPedestal(iSec,kTRUE);
576 AliTPCCalROC *rocRMS = GetCalRocRMS(iSec,kTRUE);
578 array_hP = hP->GetArray();
579 UInt_t nChannels = fROC->GetNChannels(iSec);
581 for (UInt_t iChannel=0; iChannel<nChannels; ++iChannel){
582 Int_t offset = (nbinsAdc+2)*(iChannel+1)+1;
583 Double_t ret = AliMathBase::FitGaus(array_hP+offset,nbinsAdc,fAdcMin,fAdcMax,¶m,&dummy);
584 // if the fitting failed set noise and pedestal to 0
589 rocPedestal->SetValue(iChannel,param[1]);
590 rocRMS->SetValue(iChannel,param[2]);
596 //_____________________________________________________________________
597 void AliTPCCalibPedestal::AnalyseTime(Int_t nevents)
600 // Calculate for each channel and time bin the mean pedestal. This
601 // is used on LDC by TPCPEDESTALda to generate data used for configuration
602 // of the pattern memory for baseline subtraction in the ALTROs.
605 if ( nevents <= 0 ) return;
606 if ( fTimeAnalysis ) {
607 for (Int_t i = 0; i < 159; i++) { // padrows
608 for (Int_t j = 0; j < 140; j++) { // pads per row
609 for (Int_t k = 0; k < 1024; k++) { // time bins per pad
610 fTimeSignal[i][j].AddAt(fTimeSignal[i][j].At(k)/(Float_t)nevents, k);
618 //_____________________________________________________________________
619 void AliTPCCalibPedestal::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append) /*FOLD00*/
622 // Write class to file
633 TDirectory *backup = gDirectory;
634 TFile f(filename,option.Data());
636 if ( !sDir.IsNull() ){
637 f.mkdir(sDir.Data());
643 if ( backup ) backup->cd();