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()),
210 fCalRocArrayPedestal(72),
212 fHistoPedestalArray(72),
216 // default constructor
221 //_____________________________________________________________________
222 AliTPCCalibPedestal::AliTPCCalibPedestal(const AliTPCCalibPedestal &ped) : /*FOLD00*/
224 fFirstTimeBin(ped.GetFirstTimeBin()),
225 fLastTimeBin(ped.GetLastTimeBin()),
226 fAdcMin(ped.GetAdcMin()),
227 fAdcMax(ped.GetAdcMax()),
228 fOldRCUformat(ped.fOldRCUformat),
229 fTimeAnalysis(ped.fTimeAnalysis),
230 fROC(AliTPCROC::Instance()),
231 fCalRocArrayPedestal(72),
233 fHistoPedestalArray(72),
234 fTimeSignal(ped.fTimeSignal)
239 for (Int_t iSec = 0; iSec < 72; ++iSec){
240 const AliTPCCalROC *calPed = (AliTPCCalROC*)ped.fCalRocArrayPedestal.UncheckedAt(iSec);
241 const AliTPCCalROC *calRMS = (AliTPCCalROC*)ped.fCalRocArrayRMS.UncheckedAt(iSec);
242 const TH2F *hPed = (TH2F*)ped.fHistoPedestalArray.UncheckedAt(iSec);
244 if ( calPed != 0x0 ) fCalRocArrayPedestal.AddAt(new AliTPCCalROC(*calPed), iSec);
245 if ( calRMS != 0x0 ) fCalRocArrayRMS.AddAt(new AliTPCCalROC(*calRMS), iSec);
248 TH2F *hNew = new TH2F(*hPed);
249 hNew->SetDirectory(0);
250 fHistoPedestalArray.AddAt(hNew,iSec);
256 //_____________________________________________________________________
257 AliTPCCalibPedestal& AliTPCCalibPedestal::operator = (const AliTPCCalibPedestal &source)
260 // assignment operator
262 if (&source == this) return *this;
263 new (this) AliTPCCalibPedestal(source);
269 //_____________________________________________________________________
270 AliTPCCalibPedestal::~AliTPCCalibPedestal() /*FOLD00*/
276 fCalRocArrayPedestal.Delete();
277 fCalRocArrayRMS.Delete();
278 fHistoPedestalArray.Delete();
281 for (Int_t i = 0; i < 159; i++) {
282 delete [] fTimeSignal[i];
285 delete [] fTimeSignal;
291 //_____________________________________________________________________
292 void AliTPCCalibPedestal::SetTimeAnalysis(Bool_t time)
295 // Use time dependent analysis: Pedestals are analysed as a function
296 // of the drift time. There is one mean value generated for each time
297 // bin and each channel. It can be used as reference data and for
298 // configuration of the ALTRO pattern memory for baseline subtraction.
300 // ATTENTION: Use only on LDC in TPCPEDESTALda! On a LDC we get data
301 // only from one sector. For the full TPC we would need a lot of
302 // memory (36*159*140*1024*4bytes = 3.3GB)!
305 fTimeAnalysis = time;
307 if ( !fTimeAnalysis ) return;
309 // prepare array for one sector (159*140*1024*4bytes = 92MB):
310 fTimeSignal = new TArrayF*[159];
311 for (Int_t i = 0; i < 159; i++) { // padrows
312 fTimeSignal[i] = new TArrayF[140];
313 for (Int_t j = 0; j < 140; j++) { // pads per row
314 fTimeSignal[i][j].Set(1024);
315 for (Int_t k = 0; k < 1024; k++) { // time bins per pad
316 fTimeSignal[i][j].AddAt(0., k);
323 //_____________________________________________________________________
324 Int_t AliTPCCalibPedestal::Update(const Int_t icsector, /*FOLD00*/
327 const Int_t icTimeBin,
328 const Float_t csignal)
331 // Signal filling method
334 // Time dependent pedestals
335 if ( fTimeAnalysis ) {
336 if ( icsector < 36 ) // IROC
337 fTimeSignal[icRow][icPad].AddAt(fTimeSignal[icRow][icPad].At(icTimeBin)+csignal, icTimeBin);
339 fTimeSignal[icRow+63][icPad].AddAt(fTimeSignal[icRow+63][icPad].At(icTimeBin)+csignal, icTimeBin);
341 //return if we are out of the specified time bin or adc range
342 if ( (icTimeBin>fLastTimeBin) || (icTimeBin<fFirstTimeBin) ) return 0;
343 if ( ((Int_t)csignal>fAdcMax) || ((Int_t)csignal<fAdcMin) ) return 0;
345 Int_t iChannel = fROC->GetRowIndexes(icsector)[icRow]+icPad; // global pad position in sector
347 // fast filling methode.
348 // Attention: the entry counter of the histogram is not increased
349 // this means that e.g. the colz draw option gives an empty plot
350 Int_t bin = (iChannel+1)*(fAdcMax-fAdcMin+2)+((Int_t)csignal-fAdcMin+1);
352 GetHistoPedestal(icsector,kTRUE)->GetArray()[bin]++;
358 //_____________________________________________________________________
359 Bool_t AliTPCCalibPedestal::ProcessEvent(AliTPCRawStream *rawStream)
362 // Event Processing loop - AliTPCRawStream
365 rawStream->SetOldRCUFormat(fOldRCUformat);
367 Bool_t withInput = kFALSE;
369 while (rawStream->Next()) {
371 Int_t iSector = rawStream->GetSector(); // current ROC
372 Int_t iRow = rawStream->GetRow(); // current row
373 Int_t iPad = rawStream->GetPad(); // current pad
374 Int_t iTimeBin = rawStream->GetTime(); // current time bin
375 Float_t signal = rawStream->GetSignal(); // current ADC signal
377 Update(iSector,iRow,iPad,iTimeBin,signal);
385 //_____________________________________________________________________
386 Bool_t AliTPCCalibPedestal::ProcessEvent(AliRawReader *rawReader)
389 // Event processing loop - AliRawReader
392 AliTPCRawStream rawStream(rawReader);
393 rawReader->Select("TPC");
394 return ProcessEvent(&rawStream);
398 //_____________________________________________________________________
399 Bool_t AliTPCCalibPedestal::ProcessEvent(eventHeaderStruct *event)
402 // process date event
405 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
406 Bool_t result=ProcessEvent(rawReader);
412 //_____________________________________________________________________
413 Bool_t AliTPCCalibPedestal::TestEvent() /*FOLD00*/
417 // fill one oroc and one iroc with random gaus
422 for (UInt_t iSec=0; iSec<72; ++iSec){
423 if (iSec%36>0) continue;
424 for (UInt_t iRow=0; iRow < fROC->GetNRows(iSec); ++iRow){
425 for (UInt_t iPad=0; iPad < fROC->GetNPads(iSec,iRow); ++iPad){
426 for (UInt_t iTimeBin=0; iTimeBin<1024; ++iTimeBin){
427 Float_t signal=(Int_t)(iRow+3+gRandom->Gaus(0,.7));
428 if ( signal>0 )Update(iSec,iRow,iPad,iTimeBin,signal);
437 //_____________________________________________________________________
438 TH2F* AliTPCCalibPedestal::GetHisto(Int_t sector, TObjArray *arr, /*FOLD00*/
439 Int_t nbinsY, Float_t ymin, Float_t ymax,
440 Char_t *type, Bool_t force)
443 // return pointer to Q histogram
444 // if force is true create a new histogram if it doesn't exist allready
446 if ( !force || arr->UncheckedAt(sector) )
447 return (TH2F*)arr->UncheckedAt(sector);
449 // if we are forced and histogram doesn't yes exist create it
450 Char_t name[255], title[255];
452 sprintf(name,"hCalib%s%.2d",type,sector);
453 sprintf(title,"%s calibration histogram sector %.2d;ADC channel;Channel (pad)",type,sector);
455 // new histogram with Q calib information. One value for each pad!
456 TH2F* hist = new TH2F(name,title,
458 fROC->GetNChannels(sector),0,fROC->GetNChannels(sector)
460 hist->SetDirectory(0);
461 arr->AddAt(hist,sector);
466 //_____________________________________________________________________
467 TH2F* AliTPCCalibPedestal::GetHistoPedestal(Int_t sector, Bool_t force) /*FOLD00*/
470 // return pointer to T0 histogram
471 // if force is true create a new histogram if it doesn't exist allready
473 TObjArray *arr = &fHistoPedestalArray;
474 return GetHisto(sector, arr, fAdcMax-fAdcMin, fAdcMin, fAdcMax, "Pedestal", force);
478 //_____________________________________________________________________
479 AliTPCCalROC* AliTPCCalibPedestal::GetCalRoc(Int_t sector, TObjArray* arr, Bool_t force) /*FOLD00*/
482 // return pointer to ROC Calibration
483 // if force is true create a new histogram if it doesn't exist allready
485 if ( !force || arr->UncheckedAt(sector) )
486 return (AliTPCCalROC*)arr->UncheckedAt(sector);
488 // if we are forced and the histogram doesn't yet exist create it
490 // new AliTPCCalROC for T0 information. One value for each pad!
491 AliTPCCalROC *croc = new AliTPCCalROC(sector);
492 arr->AddAt(croc,sector);
497 //_____________________________________________________________________
498 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocPedestal(Int_t sector, Bool_t force) /*FOLD00*/
501 // return pointer to Carge ROC Calibration
502 // if force is true create a new histogram if it doesn't exist allready
504 TObjArray *arr = &fCalRocArrayPedestal;
505 return GetCalRoc(sector, arr, force);
509 //_____________________________________________________________________
510 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocRMS(Int_t sector, Bool_t force) /*FOLD00*/
513 // return pointer to signal width ROC Calibration
514 // if force is true create a new histogram if it doesn't exist allready
516 TObjArray *arr = &fCalRocArrayRMS;
517 return GetCalRoc(sector, arr, force);
521 //_____________________________________________________________________
522 void AliTPCCalibPedestal::Merge(AliTPCCalibPedestal *ped)
525 // Merge reference histograms of sig to the current AliTPCCalibSignal
529 for (Int_t iSec=0; iSec<72; ++iSec){
530 TH2F *hRefPedMerge = ped->GetHistoPedestal(iSec);
533 TDirectory *dir = hRefPedMerge->GetDirectory(); hRefPedMerge->SetDirectory(0);
534 TH2F *hRefPed = GetHistoPedestal(iSec);
535 if ( hRefPed ) hRefPed->Add(hRefPedMerge);
537 TH2F *hist = new TH2F(*hRefPedMerge);
538 hist->SetDirectory(0);
539 fHistoPedestalArray.AddAt(hist, iSec);
541 hRefPedMerge->SetDirectory(dir);
551 //_____________________________________________________________________
552 void AliTPCCalibPedestal::Analyse() /*FOLD00*/
555 // Calculate calibration constants
558 Int_t nbinsAdc = fAdcMax-fAdcMin;
565 for (Int_t iSec=0; iSec<72; ++iSec){
566 TH2F *hP = GetHistoPedestal(iSec);
569 AliTPCCalROC *rocPedestal = GetCalRocPedestal(iSec,kTRUE);
570 AliTPCCalROC *rocRMS = GetCalRocRMS(iSec,kTRUE);
572 array_hP = hP->GetArray();
573 UInt_t nChannels = fROC->GetNChannels(iSec);
575 for (UInt_t iChannel=0; iChannel<nChannels; ++iChannel){
576 Int_t offset = (nbinsAdc+2)*(iChannel+1)+1;
577 Double_t ret = AliMathBase::FitGaus(array_hP+offset,nbinsAdc,fAdcMin,fAdcMax,¶m,&dummy);
578 // if the fitting failed set noise and pedestal to 0
583 rocPedestal->SetValue(iChannel,param[1]);
584 rocRMS->SetValue(iChannel,param[2]);
590 //_____________________________________________________________________
591 void AliTPCCalibPedestal::AnalyseTime(Int_t nevents)
594 // Calculate for each channel and time bin the mean pedestal. This
595 // is used on LDC by TPCPEDESTALda to generate data used for configuration
596 // of the pattern memory for baseline subtraction in the ALTROs.
599 if ( nevents <= 0 ) return;
600 if ( fTimeAnalysis ) {
601 for (Int_t i = 0; i < 159; i++) { // padrows
602 for (Int_t j = 0; j < 140; j++) { // pads per row
603 for (Int_t k = 0; k < 1024; k++) { // time bins per pad
604 fTimeSignal[i][j].AddAt(fTimeSignal[i][j].At(k)/(Float_t)nevents, k);
612 //_____________________________________________________________________
613 void AliTPCCalibPedestal::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append) /*FOLD00*/
616 // Write class to file
627 TDirectory *backup = gDirectory;
628 TFile f(filename,option.Data());
630 if ( !sDir.IsNull() ){
631 f.mkdir(sDir.Data());
637 if ( backup ) backup->cd();