1 /**************************************************************************
2 * Copyright(c) 2007, 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 **************************************************************************/
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
23 // AliRawReader* rf = new AliRawReaderDate("2006run2211.raw");
24 // AliPHOSRawDecoder dc(rf);
25 // while (rf->NextEvent()) {
26 // dc.SubtractPedestals(kTRUE);
27 // while ( dc.NextDigit() ) {
28 // Int_t module = dc.GetModule();
29 // Int_t column = dc.GetColumn();
30 // Int_t row = dc.GetRow();
31 // Double_t amplitude = dc.GetEnergy();
32 // Double_t time = dc.GetTime();
33 // Bool_t IsLowGain = dc.IsLowGain();
38 // Author: Dmitri Peressounko
40 // --- ROOT system ---
51 // --- AliRoot header files ---
53 #include "AliPHOSRawDecoderv1.h"
54 #include "AliPHOSPulseGenerator.h"
57 ClassImp(AliPHOSRawDecoderv1)
59 //-----------------------------------------------------------------------------
60 AliPHOSRawDecoderv1::AliPHOSRawDecoderv1():AliPHOSRawDecoder(),
61 fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
63 //Default constructor.
66 //-----------------------------------------------------------------------------
67 AliPHOSRawDecoderv1::AliPHOSRawDecoderv1(AliRawReader* rawReader, AliAltroMapping **mapping):
68 AliPHOSRawDecoder(rawReader,mapping),
69 fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
71 //Construct a decoder object.
72 //Is is user responsibility to provide next raw event
73 //using AliRawReader::NextEvent().
76 gMinuit = new TMinuit(100);
77 fSampleParamsHigh =new TArrayD(7) ;
78 fSampleParamsHigh->AddAt(2.174,0) ;
79 fSampleParamsHigh->AddAt(0.106,1) ;
80 fSampleParamsHigh->AddAt(0.173,2) ;
81 fSampleParamsHigh->AddAt(0.06106,3) ;
82 //last two parameters are pedestal and overflow
83 fSampleParamsLow=new TArrayD(7) ;
84 fSampleParamsLow->AddAt(2.456,0) ;
85 fSampleParamsLow->AddAt(0.137,1) ;
86 fSampleParamsLow->AddAt(2.276,2) ;
87 fSampleParamsLow->AddAt(0.08246,3) ;
88 fToFit = new TList() ;
91 //-----------------------------------------------------------------------------
92 AliPHOSRawDecoderv1::~AliPHOSRawDecoderv1()
96 delete fSampleParamsLow ;
99 if(fSampleParamsHigh){
100 delete fSampleParamsHigh ;
109 //-----------------------------------------------------------------------------
110 AliPHOSRawDecoderv1::AliPHOSRawDecoderv1(const AliPHOSRawDecoderv1 &phosDecoder ):
111 AliPHOSRawDecoder(phosDecoder),
112 fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
115 fToFit = new TList() ;
116 fSampleParamsLow =new TArrayD(*(phosDecoder.fSampleParamsLow)) ;
117 fSampleParamsHigh=new TArrayD(*(phosDecoder.fSampleParamsHigh)) ;
120 //-----------------------------------------------------------------------------
121 AliPHOSRawDecoderv1& AliPHOSRawDecoderv1::operator = (const AliPHOSRawDecoderv1 &phosDecoder)
123 //Assignment operator.
125 // if(this != &phosDecoder) {
127 fToFit = new TList() ;
128 if(fSampleParamsLow){
129 fSampleParamsLow = phosDecoder.fSampleParamsLow ;
130 fSampleParamsHigh= phosDecoder.fSampleParamsHigh ;
133 fSampleParamsLow =new TArrayD(*(phosDecoder.fSampleParamsLow)) ;
134 fSampleParamsHigh=new TArrayD(*(phosDecoder.fSampleParamsHigh)) ;
139 //-----------------------------------------------------------------------------
140 Bool_t AliPHOSRawDecoderv1::NextDigit()
142 //Extract an energy deposited in the crystal,
143 //crystal' position (module,column,row),
144 //time and gain (high or low).
145 //First collects sample, then evaluates it and if it has
146 //reasonable shape, fits it with Gamma2 function and extracts
149 //Debug=====================
150 // TCanvas * c = 0; //(TCanvas*)gROOT->FindObjectAny("CSample") ;
152 // c = new TCanvas("CSample","CSample") ;
154 // TH1D * h = 0 ; //(TH1D*)gROOT->FindObjectAny("hSample") ;
156 // h=new TH1D("hSample","",200,0.,200.) ;
158 // TF1 * fff = 0 ; //(TF1*)gROOT->FindObjectAny("fff") ;
160 // 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.) ;
161 //End debug===========
163 AliCaloRawStream* in = fCaloStream;
165 Int_t iBin = fSamples->GetSize() ;
170 Float_t baseLine = 1.0;
171 const Float_t nPreSamples = 10;
173 const Float_t sampleMaxHG=102.332 ; //maximal height of HG sample with given parameterization
174 const Float_t sampleMaxLG=277.196 ; //maximal height of HG sample with given parameterization
175 const Float_t maxEtoFit=5 ; //fit only samples above this energy, accept all samples (with good aRMS) below it
177 while ( in->Next() ) {
180 tLength = in->GetTimeLength();
181 if(tLength!=fSamples->GetSize()) {
184 fSamples = new TArrayI(tLength);
185 fTimes = new TArrayI(tLength);
186 iBin= fSamples->GetSize() ;
193 // Fit the full sample
194 if((in->IsNewHWAddress() && iBin != fSamples->GetSize()) //new HW address
195 ||(iBin<=0)) { //or new signal in same address
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() ;
205 //new handle already collected
206 Double_t pedestal =0. ;
209 pedestal = (Double_t)(pedMean/nPed);
214 //calculate time and energy
221 for(Int_t i=iBin; i<fSamples->GetSize(); i++){
222 if(fSamples->At(i)>0){
223 Double_t de=fSamples->At(i)-pedestal ;
229 if(de>2 && tStart==0)
231 if(maxAmp<fSamples->At(i)){
233 maxAmp=fSamples->At(i) ;
237 if(maxBin==fSamples->GetSize()-1){//bad "rising" sample
243 fEnergy=Double_t(maxAmp)-pedestal ;
244 fOverflow =0 ; //look for plato on the top of sample
245 if(fEnergy>500 && //this is not fluctuation of soft sample
246 maxBin<fSamples->GetSize()-1 && fSamples->At(maxBin+1)==maxAmp){ //and there is a plato
252 aRMS=aRMS/wts-aMean*aMean;
255 //do not take too small energies
256 if(fEnergy < baseLine)
259 //do not test quality of too soft samples
260 if(fEnergy<maxEtoFit){
261 fTime=fTimes->At(tStart);
262 if(aRMS<2.) //sigle peak
270 //Debug:=====Draw sample
271 //if(fEnergy>pedestal+10.){
272 //if(fLowGainFlag && fEnergy>2){
274 // if(!fLowGainFlag && fRow==32 && fColumn==18){
275 // TCanvas *c = new TCanvas("CSample","CSample") ;
279 // printf("fEnergy=%f, aRMS=%f \n",fEnergy,aRMS) ;
282 //======================
284 //IF sample has reasonable mean and RMS, try to fit it with gamma2
286 gMinuit->mncler(); // Reset Minuit's list of paramters
287 gMinuit->SetPrintLevel(-1) ; // No Printout
288 gMinuit->SetFCN(AliPHOSRawDecoderv1::UnfoldingChiSquare) ;
289 // To set the address of the minimization function
292 Double_t b,bmin,bmax ;
294 fSampleParamsLow->AddAt(pedestal,4) ;
296 fSampleParamsLow->AddAt(double(maxAmp),5) ;
298 fSampleParamsLow->AddAt(double(1023),5) ;
299 fSampleParamsLow->AddAt(double(iBin),6) ;
300 fToFit->AddFirst((TObject*)fSampleParamsLow) ;
301 b=fSampleParamsLow->At(2) ;
306 fSampleParamsHigh->AddAt(pedestal,4) ;
308 fSampleParamsHigh->AddAt(double(maxAmp),5) ;
310 fSampleParamsHigh->AddAt(double(1023),5);
311 fSampleParamsHigh->AddAt(double(iBin),6);
312 fToFit->AddFirst((TObject*)fSampleParamsHigh) ;
313 b=fSampleParamsHigh->At(2) ;
317 fToFit->AddLast((TObject*)fSamples) ;
318 fToFit->AddLast((TObject*)fTimes) ;
320 gMinuit->SetObjectFit((TObject*)fToFit) ; // To tranfer pointer to UnfoldingChiSquare
322 gMinuit->mnparm(0, "t0", 1.*tStart, 0.01, -500., 500., ierflg) ;
324 // AliWarning(Form("Unable to set initial value for fit procedure : t0=%e\n",1.*tStart) ) ;
328 return kTRUE ; //will scan further
332 amp0=fEnergy/sampleMaxLG;
334 amp0=fEnergy/sampleMaxHG;
336 gMinuit->mnparm(1, "Energy", amp0 , 0.01*amp0, 0, 0, ierflg) ;
338 // AliWarning(Form("Unable to set initial value for fit procedure : E=%e\n", amp0)) ;
342 return kTRUE ; //will scan further
345 gMinuit->mnparm(2, "p2", b, 0.01*b, bmin, bmax, ierflg) ;
347 // AliWarning(Form("Unable to set initial value for fit procedure : E=%e\n", amp0)) ;
351 return kTRUE ; //will scan further
355 Double_t p0 = 0.0001 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
359 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
360 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
361 // gMinuit->SetMaxIterations(100);
362 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
364 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
368 gMinuit->GetParameter(0,t0, t0err) ;
369 gMinuit->GetParameter(1,efit, err) ;
371 Double_t bfit, berr ;
372 gMinuit->GetParameter(2,bfit,berr) ;
374 //Calculate total energy
375 //this isparameterization of depetendence of pulse height on parameter b
377 efit*=99.54910 + 78.65038*bfit ;
379 efit*=80.33109+128.6433*bfit ;
381 if(efit<0. || efit > 10000.){
382 //set energy to previously found max
383 // fEnergy=0 ; //bad sample
389 //evaluate fit quality
390 Double_t fmin,fedm,errdef ;
391 Int_t npari,nparx,istat;
392 gMinuit->mnstat(fmin,fedm,errdef,npari,nparx,istat) ;
393 fQuality=fmin/(fSamples->GetSize()-iBin) ;
394 //compare quality with some parameterization
396 fQuality/=2.+0.002*fEnergy ;
399 fQuality/=0.75+0.0025*fEnergy ;
402 //Debug================
403 // Double_t n,alpha,beta ;
406 // n=fSampleParamsLow->At(0) ;
407 // alpha=fSampleParamsLow->At(1) ;
408 // beta=fSampleParamsLow->At(3) ;
409 // en=efit/(99.54910 + 78.65038*bfit) ;
412 // n=fSampleParamsHigh->At(0) ;
413 // alpha=fSampleParamsHigh->At(1) ;
414 // beta=fSampleParamsHigh->At(3) ;
415 // en=efit/(80.33109+128.6433*bfit) ;
418 //// if( fQuality > 1 && fEnergy > 20. && !fOverflow){
419 //// if(!fLowGainFlag && fRow==32 && fColumn==18){
421 // printf("Col=%d, row=%d, qual=%f, E=%f, t0=%f, b=%f\n",fColumn,fRow,fQuality,efit,t0,bfit) ;
422 // printf(" Energy = %f \n",fEnergy) ;
423 // TCanvas * c = new TCanvas("samp") ;
427 // fff->SetParameters(pedestal,en,t0,n,alpha,bfit,beta) ;
430 // fff->SetParameters(pedestal,en,t0,n,alpha,bfit,beta) ;
432 ////// for(Int_t i=1;i<=h->GetNbinsX(); i++){
433 //// Double_t x=h->GetBinCenter(i) ;
434 //// h->SetBinContent(i,h->GetBinContent(i)-fff->Eval(x)) ;
436 //// h->SetMinimum(-15.) ;
437 //// h->SetMaximum(15.) ;
439 // fff->Draw("same") ;
443 //====================
446 fTime=t0-4.024*bfit ; //-10.402*bfit+4.669*bfit*bfit ; //Correction for 70 samples
447 // fTime=t0+2.8*bfit ; //-10.402*bfit+4.669*bfit*bfit ; //Correction for 50 samples
453 fLowGainFlag = in->IsLowGain();
454 fModule = in->GetModule()+1;
455 fRow = in->GetRow() +1;
456 fColumn = in->GetColumn()+1;
458 //add previouly taken if coincides
459 if(fLowGainFlag==fNewLowGainFlag && fModule==fNewModule &&
460 fRow==fNewRow && fColumn==fNewColumn){
462 if(fPedSubtract && fNewTime < nPreSamples) {
463 pedMean += in->GetSignal();
466 fSamples->AddAt(fNewAmp,iBin);
467 fTimes->AddAt(fNewTime,iBin);
469 //Mark that we already take it
473 // Fill array with samples
475 if(fPedSubtract && (in->GetTime() < nPreSamples)) {
476 pedMean += in->GetSignal();
479 fSamples->AddAt(in->GetSignal(),iBin);
480 fTimes->AddAt(in->GetTime(),iBin);
482 //Debug==============
483 // h->SetBinContent(in->GetTime(),in->GetSignal()) ;
484 //EndDebug==============
490 //_____________________________________________________________________________
491 void AliPHOSRawDecoderv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
493 // Number of parameters, Gradient, Chi squared, parameters, what to do
495 TList * toFit= (TList*)gMinuit->GetObjectFit() ;
496 TArrayD * params=(TArrayD*)toFit->At(0) ;
497 TArrayI * samples = (TArrayI*)toFit->At(1) ;
498 TArrayI * times = (TArrayI*)toFit->At(2) ;
502 for(Int_t iparam = 0 ; iparam < 3 ; iparam++)
503 Grad[iparam] = 0 ; // Will evaluate gradient
508 Double_t n=params->At(0) ;
509 Double_t alpha=params->At(1) ;
510 Double_t beta=params->At(3) ;
511 Double_t ped=params->At(4) ;
513 Double_t overflow=params->At(5) ;
514 Int_t iBin = (Int_t) params->At(6) ;
515 Int_t nSamples=TMath::Min(iBin+70,samples->GetSize()) ; //Here we set number of points to fit (70)
516 // iBin - first non-zero sample
517 Int_t tStep=times->At(iBin+1)-times->At(iBin) ;
518 Double_t ddt=times->At(iBin)-t0-tStep ;
519 Double_t exp1=TMath::Exp(-alpha*ddt) ;
520 Double_t exp2=TMath::Exp(-beta*ddt) ;
521 Double_t dexp1=TMath::Exp(-alpha*tStep) ;
522 Double_t dexp2=TMath::Exp(-beta*tStep) ;
523 for(Int_t i = iBin; i<nSamples ; i++) {
524 Double_t dt=double(times->At(i))-t0 ;
525 Double_t fsample = double(samples->At(i)) ;
526 if(fsample>=overflow)
536 Double_t dtn=TMath::Power(dt,n) ;
537 Double_t dtnE=dtn*exp1 ;
538 Double_t dt2E=dt*dt*exp2 ;
539 Double_t fit=ped+en*(dtnE + b*dt2E) ;
540 // if(fit>=overflow){
541 // diff=fsample-overflow ;
542 // fret += diff*diff ;
543 // //zero gradient here
546 diff = fsample - fit ;
548 if(iflag == 2){ // calculate gradient
549 Grad[0] += en*diff*(dtnE*(n/dt-alpha)+b*dt2E*(2./dt-beta)) ; //derivative over t0
550 Grad[1] -= diff*(dtnE+b*dt2E) ;
551 Grad[2] -= en*diff*dt2E ;
556 for(Int_t iparam = 0 ; iparam < 3 ; iparam++)
559 //-----------------------------------------------------------------------------
560 Double_t AliPHOSRawDecoderv1::Gamma2(Double_t dt,Double_t en,Double_t b,TArrayD * params){ //Function for fitting samples
562 //dt-time after start
564 //function parameters
566 Double_t ped=params->At(4) ;
568 return ped ; //pedestal
570 return ped+en*(TMath::Power(dt,params->At(0))*TMath::Exp(-dt*params->At(1))+b*dt*dt*TMath::Exp(-dt*params->At(3))) ;