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 **************************************************************************/
17 //-------------------------------------------------------
18 // Implementation of the TPC pedestal and noise calibration
20 // Origin: Jens Wiechula, Marian Ivanov J.Wiechula@gsi.de, Marian.Ivanov@cern.ch
23 //-------------------------------------------------------
31 #include <TObjArray.h>
43 #include <TDirectory.h>
49 #include "AliRawReader.h"
50 #include "AliRawReaderRoot.h"
51 #include "AliTPCRawStream.h"
52 #include "AliTPCCalROC.h"
53 #include "AliTPCCalPad.h"
54 #include "AliTPCROC.h"
55 #include "AliTPCCalibPedestal.h"
56 #include "AliMathBase.h"
58 #include "TTreeStream.h"
62 ClassImp(AliTPCCalibPedestal) /*FOLD00*/
70 AliTPCCalibPedestal::AliTPCCalibPedestal() : /*FOLD00*/
76 fROC(AliTPCROC::Instance()),
77 fCalRocArrayPedestal(72),
79 fHistoPedestalArray(72),
84 // AliTPCSignal default constructor
87 TDirectory *backup = gDirectory;
88 fDebugStreamer = new TTreeSRedirector("deb2.root");
89 if ( backup ) backup->cd(); //we don't want to be cd'd to the debug streamer
92 //_____________________________________________________________________
93 AliTPCCalibPedestal::~AliTPCCalibPedestal() /*FOLD00*/
98 if ( fDebugStreamer ) delete fDebugStreamer;
105 //_____________________________________________________________________
106 Int_t AliTPCCalibPedestal::Update(const Int_t icsector, /*FOLD00*/
109 const Int_t icTimeBin,
110 const Float_t csignal)
113 // Signal filling methode
115 if ( (icTimeBin>fLastTimeBin) || (icTimeBin<fFirstTimeBin) ) return 0;
117 Int_t iChannel = fROC->GetRowIndexes(icsector)[icRow]+icPad; // global pad position in sector
119 // dirty and fast filling methode. No boundary checking!!!
120 // Attention: the entry counter of the histogram is not increased
121 // this means that e.g. the colz draw option gives an empty plot
122 Int_t bin = (iChannel+1)*(fAdcMax-fAdcMin+2)+((Int_t)csignal-fAdcMin+1);
124 GetHistoPedestal(icsector,kTRUE)->GetArray()[bin]++;
128 //_____________________________________________________________________
129 Bool_t AliTPCCalibPedestal::ProcessEvent(AliRawReader *rawReader) /*FOLD00*/
132 // simple event processing loop
136 AliTPCRawStream input(rawReader);
138 rawReader->Select("TPC");
140 input.SetOldRCUFormat(1);
141 printf("start event processing\n");
144 Bool_t withInput = kFALSE;
146 while (input.Next()) {
148 Int_t isector = input.GetSector(); // current sector
149 Int_t iRow = input.GetRow(); // current row
150 Int_t iPad = input.GetPad(); // current pad
151 Int_t iTimeBin = input.GetTime(); // current time bin
152 Float_t signal = input.GetSignal(); // current ADC signal
154 Update(isector,iRow,iPad,iTimeBin,signal);
158 printf("end event processing\n");
160 fDebugStreamer->GetFile()->Write();
163 //_____________________________________________________________________
164 Bool_t AliTPCCalibPedestal::TestEvent() /*FOLD00*/
172 for (UInt_t iSec=0; iSec<72; iSec++){
173 if (iSec%36>0) continue;
174 for (UInt_t iRow=0; iRow < fROC->GetNRows(iSec); iRow++){
175 for (UInt_t iPad=0; iPad < fROC->GetNPads(iSec,iRow); iPad++){
176 for (UInt_t iTimeBin=0; iTimeBin<1024; iTimeBin++){
177 Float_t signal=(Int_t)(iRow+5+gRandom->Gaus(0,.7));
178 if ( signal>0 )Update(iSec,iRow,iPad,iTimeBin,signal);
185 //_____________________________________________________________________
186 TH2S* AliTPCCalibPedestal::GetHisto(Int_t sector, TObjArray *arr, /*FOLD00*/
187 Int_t nbinsY, Float_t ymin, Float_t ymax,
188 Char_t *type, Bool_t force)
191 // return pointer to Q histogram
192 // if force is true create a new histogram if it doesn't exist allready
194 if ( !force || arr->UncheckedAt(sector) )
195 return (TH2S*)arr->UncheckedAt(sector);
197 // if we are forced and histogram doesn't yes exist create it
198 Char_t name[255], title[255];
200 sprintf(name,"hCalib%s%.2d",type,sector);
201 sprintf(title,"%s calibration histogram sector %.2d",type,sector);
203 // new histogram with Q calib information. One value for each pad!
204 TH2S* hist = new TH2S(name,title,
206 fROC->GetNChannels(sector),0,fROC->GetNChannels(sector)
208 hist->SetDirectory(0);
209 arr->AddAt(hist,sector);
212 //_____________________________________________________________________
213 TH2S* AliTPCCalibPedestal::GetHistoPedestal(Int_t sector, Bool_t force) /*FOLD00*/
216 // return pointer to T0 histogram
217 // if force is true create a new histogram if it doesn't exist allready
219 TObjArray *arr = &fHistoPedestalArray;
220 return GetHisto(sector, arr, fAdcMax-fAdcMin, fAdcMin, fAdcMax, "Pedestal", force);
222 //_____________________________________________________________________
223 AliTPCCalROC* AliTPCCalibPedestal::GetCalRoc(Int_t sector, TObjArray* arr, Bool_t force) /*FOLD00*/
226 // return pointer to ROC Calibration
227 // if force is true create a new histogram if it doesn't exist allready
229 if ( !force || arr->UncheckedAt(sector) )
230 return (AliTPCCalROC*)arr->UncheckedAt(sector);
232 // if we are forced and histogram doesn't yes exist create it
234 // new AliTPCCalROC for T0 information. One value for each pad!
235 AliTPCCalROC *croc = new AliTPCCalROC(sector);
237 for ( UInt_t iChannel = 0; iChannel<croc->GetNchannels(); iChannel++){
238 croc->SetValue(iChannel, 0);
240 arr->AddAt(croc,sector);
243 //_____________________________________________________________________
244 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocPedestal(Int_t sector, Bool_t force) /*FOLD00*/
247 // return pointer to Carge ROC Calibration
248 // if force is true create a new histogram if it doesn't exist allready
250 TObjArray *arr = &fCalRocArrayPedestal;
251 return GetCalRoc(sector, arr, force);
253 //_____________________________________________________________________
254 AliTPCCalROC* AliTPCCalibPedestal::GetCalRocRMS(Int_t sector, Bool_t force) /*FOLD00*/
257 // return pointer to signal width ROC Calibration
258 // if force is true create a new histogram if it doesn't exist allready
260 TObjArray *arr = &fCalRocArrayRMS;
261 return GetCalRoc(sector, arr, force);
263 //_____________________________________________________________________
264 void AliTPCCalibPedestal::Analyse() /*FOLD00*/
267 // Calculate calibration constants
270 Int_t nbinsAdc = fAdcMax-fAdcMin;
272 TH1F *py = new TH1F("htemp_py","htemp_py", nbinsAdc, fAdcMin, fAdcMax);
273 TF1 *gaus = new TF1("fit","gaus");
279 array_py = py->GetArray();
281 for (Int_t iSec=0; iSec<72; iSec++){
282 TH2S *hP = GetHistoPedestal(iSec);
285 AliTPCCalROC *rocPedestal = GetCalRocPedestal(iSec,kTRUE);
286 AliTPCCalROC *rocRMS = GetCalRocRMS(iSec,kTRUE);
288 array_hP = hP->GetArray();
289 UInt_t nChannels = fROC->GetNChannels(iSec);
291 for (UInt_t iChannel=0; iChannel<nChannels; iChannel++){
293 // set bin content of py in a dirty but fast way
294 for (Int_t iAdc=0;iAdc<nbinsAdc;iAdc++)
295 array_py[iAdc+1] = (Float_t)array_hP[(iChannel+1)*(nbinsAdc+2)+(iAdc+1)];
297 gaus->SetParameters(0,0);
300 rocPedestal->SetValue(iChannel,gaus->GetParameter(1));
301 rocRMS->SetValue(iChannel,gaus->GetParameter(2));
303 //AliMathBase::FitGaus(nbinsAdc,array_hP+(iChannel+1)*(nbinsAdc+2),¶m)
304 //rocPedestal->SetValue(iChannel,param[1]);
305 //rocRMS->SetValue(iChannel,param[2]);
310 delete fDebugStreamer;
311 fDebugStreamer = 0x0;
313 //_____________________________________________________________________
314 void AliTPCCalibPedestal::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append) /*FOLD00*/
317 // Write class to file
320 TDirectory *backup = gDirectory;
329 TFile f(filename,option.Data());
330 if ( !sDir.IsNull() ){
331 f.mkdir(sDir.Data());
334 gDirectory->WriteTObject(this);
337 if ( backup ) backup->cd();