]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSRawDecoderv1.cxx
All time calibration moved to AliPHOSRawDigiProducer
[u/mrichter/AliRoot.git] / PHOS / AliPHOSRawDecoderv1.cxx
CommitLineData
77ea1c6f 1/**************************************************************************
2 * Copyright(c) 2007, 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
16/* $Id$ */
17
18// This class decodes the stream of ALTRO samples to extract
19// the PHOS "digits" of current event. Uses fitting procedure
20// to separate reasonable samples
21//
22// Typical use case:
23// AliRawReader* rf = new AliRawReaderDate("2006run2211.raw");
24// AliPHOSRawDecoder dc(rf);
25// while (rf->NextEvent()) {
26// dc.SetOldRCUFormat(kTRUE);
27// dc.SubtractPedestals(kTRUE);
28// while ( dc.NextDigit() ) {
29// Int_t module = dc.GetModule();
30// Int_t column = dc.GetColumn();
31// Int_t row = dc.GetRow();
32// Double_t amplitude = dc.GetEnergy();
33// Double_t time = dc.GetTime();
34// Bool_t IsLowGain = dc.IsLowGain();
35// ..........
36// }
37// }
38
39// Author: Dmitri Peressounko
40
41// --- ROOT system ---
11f2ec15 42#include "TList.h"
77ea1c6f 43#include "TMath.h"
44#include "TMinuit.h"
45
11f2ec15 46#include "TCanvas.h"
47#include "TH1.h"
48#include "TH2.h"
49#include "TF1.h"
50#include "TROOT.h"
51
77ea1c6f 52// --- AliRoot header files ---
11f2ec15 53//#include "AliLog.h"
77ea1c6f 54#include "AliPHOSRawDecoderv1.h"
55#include "AliPHOSPulseGenerator.h"
56
57
58ClassImp(AliPHOSRawDecoderv1)
59
60//-----------------------------------------------------------------------------
11f2ec15 61 AliPHOSRawDecoderv1::AliPHOSRawDecoderv1():AliPHOSRawDecoder(),
62fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
77ea1c6f 63{
64 //Default constructor.
65}
66
67//-----------------------------------------------------------------------------
68AliPHOSRawDecoderv1::AliPHOSRawDecoderv1(AliRawReader* rawReader, AliAltroMapping **mapping):
11f2ec15 69 AliPHOSRawDecoder(rawReader,mapping),
70 fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
77ea1c6f 71{
72 //Construct a decoder object.
73 //Is is user responsibility to provide next raw event
74 //using AliRawReader::NextEvent().
75
76 if(!gMinuit)
77 gMinuit = new TMinuit(100);
70d93620 78 fSampleParamsHigh =new TArrayD(7) ;
79 fSampleParamsHigh->AddAt(2.174,0) ;
80 fSampleParamsHigh->AddAt(0.106,1) ;
81 fSampleParamsHigh->AddAt(0.173,2) ;
82 fSampleParamsHigh->AddAt(0.06106,3) ;
83 //last two parameters are pedestal and overflow
84 fSampleParamsLow=new TArrayD(7) ;
85 fSampleParamsLow->AddAt(2.456,0) ;
86 fSampleParamsLow->AddAt(0.137,1) ;
87 fSampleParamsLow->AddAt(2.276,2) ;
88 fSampleParamsLow->AddAt(0.08246,3) ;
11f2ec15 89 fToFit = new TList() ;
77ea1c6f 90}
91
92//-----------------------------------------------------------------------------
93AliPHOSRawDecoderv1::~AliPHOSRawDecoderv1()
94{
95 //Destructor.
11f2ec15 96 if(fSampleParamsLow){
97 delete fSampleParamsLow ;
98 fSampleParamsLow=0 ;
99 }
100 if(fSampleParamsHigh){
101 delete fSampleParamsHigh ;
102 fSampleParamsHigh=0;
103 }
104 if(fToFit){
105 delete fToFit ;
106 fToFit=0 ;
107 }
77ea1c6f 108}
109
110//-----------------------------------------------------------------------------
111AliPHOSRawDecoderv1::AliPHOSRawDecoderv1(const AliPHOSRawDecoderv1 &phosDecoder ):
11f2ec15 112 AliPHOSRawDecoder(phosDecoder),
113 fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
77ea1c6f 114{
115 //Copy constructor.
11f2ec15 116 fToFit = new TList() ;
117 fSampleParamsLow =new TArrayD(*(phosDecoder.fSampleParamsLow)) ;
118 fSampleParamsHigh=new TArrayD(*(phosDecoder.fSampleParamsHigh)) ;
77ea1c6f 119}
120
121//-----------------------------------------------------------------------------
11f2ec15 122AliPHOSRawDecoderv1& AliPHOSRawDecoderv1::operator = (const AliPHOSRawDecoderv1 &phosDecoder)
77ea1c6f 123{
124 //Assignment operator.
125
11f2ec15 126 // if(this != &phosDecoder) {
77ea1c6f 127 // }
11f2ec15 128 fToFit = new TList() ;
129 if(fSampleParamsLow){
130 fSampleParamsLow = phosDecoder.fSampleParamsLow ;
131 fSampleParamsHigh= phosDecoder.fSampleParamsHigh ;
132 }
133 else{
134 fSampleParamsLow =new TArrayD(*(phosDecoder.fSampleParamsLow)) ;
135 fSampleParamsHigh=new TArrayD(*(phosDecoder.fSampleParamsHigh)) ;
136 }
77ea1c6f 137 return *this;
138}
139
140//-----------------------------------------------------------------------------
141Bool_t AliPHOSRawDecoderv1::NextDigit()
142{
143 //Extract an energy deposited in the crystal,
144 //crystal' position (module,column,row),
145 //time and gain (high or low).
146 //First collects sample, then evaluates it and if it has
147 //reasonable shape, fits it with Gamma2 function and extracts
148 //energy and time.
11f2ec15 149
70d93620 150//Debug=====================
39569d36 151// TCanvas * c = 0; //(TCanvas*)gROOT->FindObjectAny("CSample") ;
70d93620 152// if(!c)
153// c = new TCanvas("CSample","CSample") ;
154//
39569d36 155// TH1D * h = 0 ; //(TH1D*)gROOT->FindObjectAny("hSample") ;
70d93620 156// if(!h)
157// h=new TH1D("hSample","",200,0.,200.) ;
158//
39569d36 159// TF1 * fff = 0 ; //(TF1*)gROOT->FindObjectAny("fff") ;
70d93620 160// if(!fff)
161// fff = new TF1("fff","[0]+[1]*((abs(x-[2]))^[3]*exp(-(x-[2])*[4])+[5]*(x-[2])*(x-[2])*exp(-(x-[2])*[6]))",0.,1000.) ;
162//End debug===========
77ea1c6f 163
164 AliCaloRawStream* in = fCaloStream;
165
70d93620 166 Int_t iBin = fSamples->GetSize() ;
77ea1c6f 167 Int_t tLength = 0;
168 fEnergy = -111;
169 Float_t pedMean = 0;
170 Int_t nPed = 0;
171 Float_t baseLine = 1.0;
172 const Float_t nPreSamples = 10;
70d93620 173 fQuality= 999. ;
174 const Float_t sampleMaxHG=102.332 ; //maximal height of HG sample with given parameterization
175 const Float_t sampleMaxLG=277.196 ; //maximal height of HG sample with given parameterization
39569d36 176 const Float_t maxEtoFit=5 ; //fit only samples above this energy, accept all samples (with good aRMS) below it
77ea1c6f 177
178 while ( in->Next() ) {
179
180 if(!tLength) {
181 tLength = in->GetTimeLength();
182 if(tLength!=fSamples->GetSize()) {
183 delete fSamples ;
70d93620 184 delete fTimes ;
77ea1c6f 185 fSamples = new TArrayI(tLength);
70d93620 186 fTimes = new TArrayI(tLength);
187 iBin= fSamples->GetSize() ;
77ea1c6f 188 }
189 else{
70d93620 190 fSamples->Reset() ;
77ea1c6f 191 }
192 }
193
194 // Fit the full sample
70d93620 195 if(in->IsNewHWAddress() && iBin != fSamples->GetSize()) {
39569d36 196
197 //First remember new sample
198 fNewLowGainFlag = in->IsLowGain();
199 fNewModule = in->GetModule()+1;
200 fNewRow = in->GetRow() +1;
201 fNewColumn = in->GetColumn()+1;
202 fNewAmp = in->GetSignal() ;
203 fNewTime=in->GetTime() ;
204 //new handle already collected
77ea1c6f 205 Double_t pedestal =0. ;
206 if(fPedSubtract){
207 if (nPed > 0)
208 pedestal = (Double_t)(pedMean/nPed);
209 else
210 return kFALSE;
211 }
70d93620 212
213 //calculate time and energy
214 Int_t maxBin=0 ;
215 Int_t maxAmp=0 ;
216 Double_t aMean=0. ;
217 Double_t aRMS=0. ;
39569d36 218 Double_t wts=0 ;
70d93620 219 Int_t tStart = 0 ;
70d93620 220 for(Int_t i=iBin; i<fSamples->GetSize(); i++){
221 if(fSamples->At(i)>0){
222 Double_t de=fSamples->At(i)-pedestal ;
39569d36 223 if(de>1.){
224 aMean+=de*i ;
225 aRMS+=de*i*i ;
226 wts+=de;
227 }
70d93620 228 if(de>2 && tStart==0)
229 tStart=i ;
230 if(maxAmp<fSamples->At(i)){
231 maxBin=i ;
232 maxAmp=fSamples->At(i) ;
233 }
234 }
235 }
39569d36 236 if(maxBin==fSamples->GetSize()-1){//bad "rising" sample
70d93620 237 fEnergy=0. ;
238 fTime=-999.;
239 fQuality= 999. ;
240 return kTRUE ;
241 }
242 fEnergy=Double_t(maxAmp)-pedestal ;
243 fOverflow =0 ; //look for plato on the top of sample
244 if(fEnergy>500 && //this is not fluctuation of soft sample
245 maxBin<fSamples->GetSize()-1 && fSamples->At(maxBin+1)==maxAmp){ //and there is a plato
246 fOverflow = kTRUE ;
77ea1c6f 247 }
70d93620 248
39569d36 249 if(wts>0){
250 aMean/=wts;
251 aRMS=aRMS/wts-aMean*aMean;
77ea1c6f 252 }
39569d36 253
254 //do not take too small energies
255 if(fEnergy < baseLine)
256 fEnergy = 0;
257
258 //do not test quality of too soft samples
259 if(fEnergy<maxEtoFit){
260 fTime=fTimes->At(tStart);
261 if(aRMS<2.) //sigle peak
262 fQuality=999. ;
263 else
264 fQuality= 0. ;
265 return kTRUE ;
266 }
267
77ea1c6f 268
70d93620 269//Debug:=====Draw sample
270//if(fEnergy>pedestal+10.){
39569d36 271//if(fLowGainFlag && fEnergy>2){
272// if(!c)
273// c = new TCanvas("CSample","CSample") ;
70d93620 274// c->cd() ;
275// h->Draw() ;
276// c->Update() ;
39569d36 277// printf("fEnergy=%f, aRMS=%f \n",fEnergy,aRMS) ;
278//getchar() ;
70d93620 279//}
280//======================
281
77ea1c6f 282 //IF sample has reasonable mean and RMS, try to fit it with gamma2
77ea1c6f 283
284 gMinuit->mncler(); // Reset Minuit's list of paramters
285 gMinuit->SetPrintLevel(-1) ; // No Printout
286 gMinuit->SetFCN(AliPHOSRawDecoderv1::UnfoldingChiSquare) ;
287 // To set the address of the minimization function
11f2ec15 288
289 fToFit->Clear() ;
290 if(fLowGainFlag){
291 fSampleParamsLow->AddAt(pedestal,4) ;
70d93620 292 if(fOverflow)
293 fSampleParamsLow->AddAt(double(maxAmp),5) ;
294 else
295 fSampleParamsLow->AddAt(double(1023),5) ;
296 fSampleParamsLow->AddAt(double(iBin),6) ;
11f2ec15 297 fToFit->AddFirst((TObject*)fSampleParamsLow) ;
298 }
299 else{
300 fSampleParamsHigh->AddAt(pedestal,4) ;
70d93620 301 if(fOverflow)
302 fSampleParamsHigh->AddAt(double(maxAmp),5) ;
303 else
304 fSampleParamsHigh->AddAt(double(1023),5);
305 fSampleParamsHigh->AddAt(double(iBin),6);
11f2ec15 306 fToFit->AddFirst((TObject*)fSampleParamsHigh) ;
307 }
308 fToFit->AddLast((TObject*)fSamples) ;
70d93620 309 fToFit->AddLast((TObject*)fTimes) ;
11f2ec15 310
311 gMinuit->SetObjectFit((TObject*)fToFit) ; // To tranfer pointer to UnfoldingChiSquare
77ea1c6f 312 Int_t ierflg ;
11f2ec15 313 gMinuit->mnparm(0, "t0", 1.*tStart, 0.1, 0, 0, ierflg) ;
77ea1c6f 314 if(ierflg != 0){
11f2ec15 315// AliWarning(Form("Unable to set initial value for fit procedure : t0=%e\n",1.*tStart) ) ;
77ea1c6f 316 fEnergy=0. ;
317 fTime=-999. ;
70d93620 318 fQuality=999 ;
77ea1c6f 319 return kTRUE ; //will scan further
320 }
70d93620 321 Double_t amp0;
322 if(fLowGainFlag)
323 amp0=fEnergy/sampleMaxLG;
324 else
325 amp0=fEnergy/sampleMaxHG;
11f2ec15 326
327 gMinuit->mnparm(1, "Energy", amp0 , 0.001*amp0, 0, 0, ierflg) ;
77ea1c6f 328 if(ierflg != 0){
11f2ec15 329// AliWarning(Form("Unable to set initial value for fit procedure : E=%e\n", amp0)) ;
77ea1c6f 330 fEnergy=0. ;
331 fTime=-999. ;
70d93620 332 fQuality=999 ;
77ea1c6f 333 return kTRUE ; //will scan further
334 }
335
11f2ec15 336 Double_t p0 = 0.0001 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
77ea1c6f 337 // depends on it.
338 Double_t p1 = 1.0 ;
339 Double_t p2 = 0.0 ;
340 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
341 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
342 // gMinuit->SetMaxIterations(100);
343 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
344
345 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
70d93620 346
77ea1c6f 347
11f2ec15 348 Double_t err,t0err ;
349 Double_t t0,efit ;
350 gMinuit->GetParameter(0,t0, t0err) ;
351 gMinuit->GetParameter(1,efit, err) ;
11f2ec15 352
70d93620 353 //Calculate total energy
354 if(fLowGainFlag)
355 efit*=sampleMaxLG;
356 else
357 efit*=sampleMaxHG;
358
359 if(efit<0. || efit > 10000.){
360 fEnergy=0 ; //bad sample
361 fTime=-999.;
362 fQuality=999 ;
363 return kTRUE;
364 }
365
366 //evaluate fit quality
367 Double_t fmin,fedm,errdef ;
368 Int_t npari,nparx,istat;
369 gMinuit->mnstat(fmin,fedm,errdef,npari,nparx,istat) ;
370 fQuality=fmin/(fSamples->GetSize()-iBin) ;
371 //compare quality with some parameterization
372 fQuality/=5.*TMath::Exp(0.0025*efit) ;
373
374//Debug================
375// Double_t n,alpha,b,beta ;
376// if(fLowGainFlag){
377// n=fSampleParamsLow->At(0) ;
378// alpha=fSampleParamsLow->At(1) ;
379// b=fSampleParamsLow->At(2) ;
380// beta=fSampleParamsLow->At(3) ;
381// }
382// else{
383// n=fSampleParamsHigh->At(0) ;
384// alpha=fSampleParamsHigh->At(1) ;
39569d36 385// b=fSampleParamsHigh->At(2) ;
70d93620 386// beta=fSampleParamsHigh->At(3) ;
387// }
388//
389// if( fQuality > 5.*TMath::Exp(0.0025*efit)){
390// if( fQuality > 1.){
39569d36 391// if(fLowGainFlag){
70d93620 392// printf("Col=%d, row=%d, qual=%f, E=%f \n",fColumn,fRow,fQuality,efit) ;
11f2ec15 393// c->cd() ;
394// h->Draw() ;
395// if(fLowGainFlag){
70d93620 396// fff->SetParameters(pedestal,efit/sampleMaxLG,t0,n,alpha,b,beta) ;
11f2ec15 397// }
398// else{
70d93620 399// fff->SetParameters(pedestal,efit/sampleMaxHG,t0,n,alpha,b,beta) ;
11f2ec15 400// }
401// fff->Draw("same") ;
402// c->Update();
70d93620 403// getchar() ;
404// }
405//====================
11f2ec15 406
39569d36 407 fEnergy=efit ;
408 fTime=t0 ;
11f2ec15 409
77ea1c6f 410
411 return kTRUE;
412 }
413
414 fLowGainFlag = in->IsLowGain();
77ea1c6f 415 fModule = in->GetModule()+1;
416 fRow = in->GetRow() +1;
417 fColumn = in->GetColumn()+1;
39569d36 418
419 //add previouly taken if coincides
420 if(fLowGainFlag==fNewLowGainFlag && fModule==fNewModule &&
421 fRow==fNewRow && fColumn==fNewColumn){
422 iBin--;
423 if(fPedSubtract && fNewTime < nPreSamples) {
424 pedMean += in->GetSignal();
425 nPed++;
426 }
427 fSamples->AddAt(fNewAmp,iBin);
428 fTimes->AddAt(fNewTime,iBin);
429
430 //Mark that we already take it
431 fNewModule=-1 ;
432 }
77ea1c6f 433
77ea1c6f 434 // Fill array with samples
70d93620 435 iBin--;
436 if(fPedSubtract && (in->GetTime() < nPreSamples)) {
77ea1c6f 437 pedMean += in->GetSignal();
438 nPed++;
439 }
70d93620 440 fSamples->AddAt(in->GetSignal(),iBin);
441 fTimes->AddAt(in->GetTime(),iBin);
442
443//Debug==============
444// h->SetBinContent(iBin,in->GetSignal()) ;
77ea1c6f 445
446 } // in.Next()
447
448 return kFALSE;
449}
450//_____________________________________________________________________________
451void AliPHOSRawDecoderv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
452{
453 // Calculates the Chi square for the samples minimization
454 // Number of parameters, Gradient, Chi squared, parameters, what to do
455
11f2ec15 456 TList * toFit= (TList*)gMinuit->GetObjectFit() ;
457 TArrayD * params=(TArrayD*)toFit->At(0) ;
458 TArrayI * samples = (TArrayI*)toFit->At(1) ;
70d93620 459 TArrayI * times = (TArrayI*)toFit->At(2) ;
77ea1c6f 460
461 fret = 0. ;
462 if(iflag == 2)
11f2ec15 463 for(Int_t iparam = 0 ; iparam < 2 ; iparam++)
77ea1c6f 464 Grad[iparam] = 0 ; // Will evaluate gradient
465
11f2ec15 466 Double_t t0=x[0] ;
467 Double_t en=x[1] ;
70d93620 468 Double_t n=params->At(0) ;
11f2ec15 469 Double_t alpha=params->At(1) ;
470 Double_t b=params->At(2) ;
471 Double_t beta=params->At(3) ;
70d93620 472 Int_t overflow=(Int_t)params->At(5) ;
473 Int_t iBin = (Int_t) params->At(6) ;
474 Int_t nSamples=TMath::Min(iBin+70,samples->GetSize()) ; //Here we set number of points to fit (70)
77ea1c6f 475
70d93620 476 for(Int_t i = iBin ; i < nSamples ; i++) {
477 Int_t sample = samples->At(i) ;
478 if(sample==0 || sample==overflow) //zero or overflow - scip point
77ea1c6f 479 continue ;
70d93620 480 Double_t dt=1.*times->At(i)-t0 ;
481 Double_t diff=float(sample)-Gamma2(dt,en,params) ;
77ea1c6f 482 if(iflag == 2){ // calculate gradient
77ea1c6f 483 if(dt>=0.){
70d93620 484 Grad[0] += 2.*en*diff*(TMath::Power(dt,n-1.)*(n-alpha*dt)*TMath::Exp(-alpha*dt)+
485 b*dt*(2.-beta*dt)*TMath::Exp(-beta*dt)) ; //derivative over t0
486 Grad[1] += -2.*diff*(TMath::Power(dt,n)*TMath::Exp(-alpha*dt)+
487 b*dt*dt*TMath::Exp(-dt*beta)) ;
77ea1c6f 488 }
489 }
490 fret += diff*diff ;
491 }
77ea1c6f 492
493}
494//-----------------------------------------------------------------------------
11f2ec15 495Double_t AliPHOSRawDecoderv1::Gamma2(Double_t dt,Double_t en,TArrayD * params){ //Function for fitting samples
77ea1c6f 496 //parameters:
11f2ec15 497 //dt-time after start
77ea1c6f 498 //en-amplutude
11f2ec15 499 //function parameters
500
501 Double_t ped=params->At(4) ;
77ea1c6f 502 if(dt<0.)
11f2ec15 503 return ped ; //pedestal
77ea1c6f 504 else
70d93620 505 return ped+en*(TMath::Power(dt,params->At(0))*TMath::Exp(-dt*params->At(1))+params->At(2)*dt*dt*TMath::Exp(-dt*params->At(3))) ;
77ea1c6f 506}
507