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 ---
50 // --- AliRoot header files ---
51 #include "AliPHOSCalibData.h"
52 #include "AliPHOSRawDecoderv1.h"
53 #include "AliPHOSPulseGenerator.h"
55 ClassImp(AliPHOSRawDecoderv1)
57 //-----------------------------------------------------------------------------
58 AliPHOSRawDecoderv1::AliPHOSRawDecoderv1():AliPHOSRawDecoder(),
59 fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
61 //Default constructor.
64 //-----------------------------------------------------------------------------
65 AliPHOSRawDecoderv1::AliPHOSRawDecoderv1(AliRawReader* rawReader, AliAltroMapping **mapping):
66 AliPHOSRawDecoder(rawReader,mapping),
67 fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
69 //Construct a decoder object.
70 //Is is user responsibility to provide next raw event
71 //using AliRawReader::NextEvent().
74 gMinuit = new TMinuit(100);
75 fSampleParamsHigh =new TArrayD(7) ;
76 fSampleParamsHigh->AddAt(2.174,0) ;
77 fSampleParamsHigh->AddAt(0.106,1) ;
78 fSampleParamsHigh->AddAt(0.173,2) ;
79 fSampleParamsHigh->AddAt(0.06106,3) ;
80 //last two parameters are pedestal and overflow
81 fSampleParamsLow=new TArrayD(7) ;
82 fSampleParamsLow->AddAt(2.456,0) ;
83 fSampleParamsLow->AddAt(0.137,1) ;
84 fSampleParamsLow->AddAt(2.276,2) ;
85 fSampleParamsLow->AddAt(0.08246,3) ;
86 fToFit = new TList() ;
89 //-----------------------------------------------------------------------------
90 AliPHOSRawDecoderv1::~AliPHOSRawDecoderv1()
94 delete fSampleParamsLow ;
97 if(fSampleParamsHigh){
98 delete fSampleParamsHigh ;
107 //-----------------------------------------------------------------------------
108 AliPHOSRawDecoderv1::AliPHOSRawDecoderv1(const AliPHOSRawDecoderv1 &phosDecoder ):
109 AliPHOSRawDecoder(phosDecoder),
110 fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
113 fToFit = new TList() ;
114 fSampleParamsLow =new TArrayD(*(phosDecoder.fSampleParamsLow)) ;
115 fSampleParamsHigh=new TArrayD(*(phosDecoder.fSampleParamsHigh)) ;
118 //-----------------------------------------------------------------------------
119 AliPHOSRawDecoderv1& AliPHOSRawDecoderv1::operator = (const AliPHOSRawDecoderv1 &phosDecoder)
121 //Assignment operator.
123 // if(this != &phosDecoder) {
125 fToFit = new TList() ;
126 if(fSampleParamsLow){
127 fSampleParamsLow = phosDecoder.fSampleParamsLow ;
128 fSampleParamsHigh= phosDecoder.fSampleParamsHigh ;
131 fSampleParamsLow =new TArrayD(*(phosDecoder.fSampleParamsLow)) ;
132 fSampleParamsHigh=new TArrayD(*(phosDecoder.fSampleParamsHigh)) ;
137 //-----------------------------------------------------------------------------
138 Bool_t AliPHOSRawDecoderv1::NextDigit()
140 //Extract an energy deposited in the crystal,
141 //crystal' position (module,column,row),
142 //time and gain (high or low).
143 //First collects sample, then evaluates it and if it has
144 //reasonable shape, fits it with Gamma2 function and extracts
147 //Debug=====================
148 // TCanvas * c = 0; //(TCanvas*)gROOT->FindObjectAny("CSample") ;
150 // c = new TCanvas("CSample","CSample") ;
152 // TH1D * h = 0 ; //(TH1D*)gROOT->FindObjectAny("hSample") ;
154 // h=new TH1D("hSample","",200,0.,200.) ;
156 // TF1 * fff = 0 ; //(TF1*)gROOT->FindObjectAny("fff") ;
158 // 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.) ;
159 //End debug===========
161 AliCaloRawStream* in = fCaloStream;
163 Int_t iBin = fSamples->GetSize() ;
169 Float_t baseLine = 1.0;
170 const Float_t nPreSamples = 10;
172 const Float_t sampleMaxHG=102.332 ; //maximal height of HG sample with given parameterization
173 const Float_t sampleMaxLG=277.196 ; //maximal height of HG sample with given parameterization
174 const Float_t maxEtoFit=5 ; //fit only samples above this energy, accept all samples (with good aRMS) below it
178 while ( in->Next() ) {
181 tLength = in->GetTimeLength();
182 if(tLength!=fSamples->GetSize()) {
185 fSamples = new TArrayI(tLength);
186 fTimes = new TArrayI(tLength);
187 iBin= fSamples->GetSize() ;
194 // Fit the full sample
195 if((in->IsNewHWAddress() && iBin != fSamples->GetSize()) //new HW address
196 ||(iBin<=0)) { //or new signal in same address
198 //First remember new sample
199 fNewLowGainFlag = in->IsLowGain();
200 fNewModule = in->GetModule()+1;
201 fNewRow = in->GetRow() +1;
202 fNewColumn = in->GetColumn()+1;
203 fNewAmp = in->GetSignal() ;
204 fNewTime = in->GetTime() ;
206 //now handle already collected
207 Double_t pedestal =0. ;
211 pedestal = (Double_t)(pedMean/nPed);
212 fPedestalRMS=pedRMS/nPed-pedestal*pedestal ;
213 if(fPedestalRMS>0.) fPedestalRMS=TMath::Sqrt(fPedestalRMS) ;
219 //take pedestals from DB
220 pedestal = fAmpOffset ;
222 Float_t truePed = fCalibData->GetADCpedestalEmc(fModule, fColumn, fRow) ;
223 Int_t altroSettings = fCalibData->GetAltroOffsetEmc(fModule, fColumn, fRow) ;
224 pedestal += truePed - altroSettings ;
227 // printf("AliPHOSRawDecoderv1::NextDigit(): Can not read data from OCDB \n") ;
231 //calculate time and energy
238 for(Int_t i=iBin; i<fSamples->GetSize(); i++){
239 if(fSamples->At(i)>pedestal){
240 Double_t de=fSamples->At(i)-pedestal ;
246 if(de>2 && tStart==0)
248 if(maxAmp<fSamples->At(i)){
250 maxAmp=fSamples->At(i) ;
254 if(maxBin==fSamples->GetSize()-1){//bad "rising" sample
260 fEnergy=Double_t(maxAmp)-pedestal ;
261 fOverflow =0 ; //look for plato on the top of sample
262 if(fEnergy>500 && //this is not fluctuation of soft sample
263 maxBin<fSamples->GetSize()-1 && fSamples->At(maxBin+1)==maxAmp){ //and there is a plato
269 aRMS=aRMS/wts-aMean*aMean;
272 //do not take too small energies
273 if(fEnergy < baseLine)
276 //do not test quality of too soft samples
277 if(fEnergy<maxEtoFit){
278 fTime=fTimes->At(tStart);
279 if(aRMS<2.) //sigle peak
287 //Debug:=====Draw sample
288 //if(fEnergy>pedestal+10.){
289 //if(fLowGainFlag && fEnergy>2){
291 // if(!fLowGainFlag && fRow==32 && fColumn==18){
292 // TCanvas *c = new TCanvas("CSample","CSample") ;
296 // printf("fEnergy=%f, aRMS=%f \n",fEnergy,aRMS) ;
299 //======================
301 //IF sample has reasonable mean and RMS, try to fit it with gamma2
303 gMinuit->mncler(); // Reset Minuit's list of paramters
304 gMinuit->SetPrintLevel(-1) ; // No Printout
305 gMinuit->SetFCN(AliPHOSRawDecoderv1::UnfoldingChiSquare) ;
306 // To set the address of the minimization function
308 fToFit->Clear("nodelete") ;
309 Double_t b,bmin,bmax ;
311 fSampleParamsLow->AddAt(pedestal,4) ;
313 fSampleParamsLow->AddAt(double(maxAmp),5) ;
315 fSampleParamsLow->AddAt(double(1023),5) ;
316 fSampleParamsLow->AddAt(double(iBin),6) ;
317 fToFit->AddFirst((TObject*)fSampleParamsLow) ;
318 b=fSampleParamsLow->At(2) ;
323 fSampleParamsHigh->AddAt(pedestal,4) ;
325 fSampleParamsHigh->AddAt(double(maxAmp),5) ;
327 fSampleParamsHigh->AddAt(double(1023),5);
328 fSampleParamsHigh->AddAt(double(iBin),6);
329 fToFit->AddFirst((TObject*)fSampleParamsHigh) ;
330 b=fSampleParamsHigh->At(2) ;
334 fToFit->AddLast((TObject*)fSamples) ;
335 fToFit->AddLast((TObject*)fTimes) ;
337 gMinuit->SetObjectFit((TObject*)fToFit) ; // To tranfer pointer to UnfoldingChiSquare
339 gMinuit->mnparm(0, "t0", 1.*tStart, 0.01, -500., 500., ierflg) ;
341 // AliWarning(Form("Unable to set initial value for fit procedure : t0=%e\n",1.*tStart) ) ;
345 return kTRUE ; //will scan further
349 amp0=fEnergy/sampleMaxLG;
351 amp0=fEnergy/sampleMaxHG;
353 gMinuit->mnparm(1, "Energy", amp0 , 0.01*amp0, 0, 0, ierflg) ;
355 // AliWarning(Form("Unable to set initial value for fit procedure : E=%e\n", amp0)) ;
359 return kTRUE ; //will scan further
362 gMinuit->mnparm(2, "p2", b, 0.01*b, bmin, bmax, ierflg) ;
364 // AliWarning(Form("Unable to set initial value for fit procedure : E=%e\n", amp0)) ;
368 return kTRUE ; //will scan further
372 Double_t p0 = 0.0001 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
376 gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ; // force TMinuit to reduce function calls
377 gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ; // force TMinuit to use my gradient
378 // gMinuit->SetMaxIterations(100);
379 gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ; // No Warnings
381 gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ; // minimize
385 gMinuit->GetParameter(0,t0, t0err) ;
386 gMinuit->GetParameter(1,efit, err) ;
388 Double_t bfit, berr ;
389 gMinuit->GetParameter(2,bfit,berr) ;
391 //Calculate total energy
392 //this isparameterization of depetendence of pulse height on parameter b
394 efit*=99.54910 + 78.65038*bfit ;
396 efit*=80.33109+128.6433*bfit ;
398 if(efit<0. || efit > 10000.){
399 //set energy to previously found max
400 // fEnergy=0 ; //bad sample
406 //evaluate fit quality
407 Double_t fmin,fedm,errdef ;
408 Int_t npari,nparx,istat;
409 gMinuit->mnstat(fmin,fedm,errdef,npari,nparx,istat) ;
410 fQuality=fmin/(fSamples->GetSize()-iBin) ;
411 //compare quality with some parameterization
413 fQuality/=2.+0.002*fEnergy ;
416 fQuality/=0.75+0.0025*fEnergy ;
419 //Debug================
420 // Double_t n,alpha,beta ;
423 // n=fSampleParamsLow->At(0) ;
424 // alpha=fSampleParamsLow->At(1) ;
425 // beta=fSampleParamsLow->At(3) ;
426 // en=efit/(99.54910 + 78.65038*bfit) ;
429 // n=fSampleParamsHigh->At(0) ;
430 // alpha=fSampleParamsHigh->At(1) ;
431 // beta=fSampleParamsHigh->At(3) ;
432 // en=efit/(80.33109+128.6433*bfit) ;
435 //// if( fQuality > 1 && fEnergy > 20. && !fOverflow){
436 //// if(!fLowGainFlag && fRow==32 && fColumn==18){
438 // printf("Col=%d, row=%d, qual=%f, E=%f, t0=%f, b=%f\n",fColumn,fRow,fQuality,efit,t0,bfit) ;
439 // printf(" Energy = %f \n",fEnergy) ;
440 // TCanvas * c = new TCanvas("samp") ;
444 // fff->SetParameters(pedestal,en,t0,n,alpha,bfit,beta) ;
447 // fff->SetParameters(pedestal,en,t0,n,alpha,bfit,beta) ;
449 ////// for(Int_t i=1;i<=h->GetNbinsX(); i++){
450 //// Double_t x=h->GetBinCenter(i) ;
451 //// h->SetBinContent(i,h->GetBinContent(i)-fff->Eval(x)) ;
453 //// h->SetMinimum(-15.) ;
454 //// h->SetMaximum(15.) ;
456 // fff->Draw("same") ;
460 //====================
463 fTime=t0-4.024*bfit ; //-10.402*bfit+4.669*bfit*bfit ; //Correction for 70 samples
464 // fTime=t0+2.8*bfit ; //-10.402*bfit+4.669*bfit*bfit ; //Correction for 50 samples
469 fLowGainFlag = in->IsLowGain();
470 fModule = in->GetModule()+1;
471 fRow = in->GetRow() +1;
472 fColumn = in->GetColumn()+1;
474 //add previouly taken if coincides
475 if(fLowGainFlag==fNewLowGainFlag && fModule==fNewModule &&
476 fRow==fNewRow && fColumn==fNewColumn){
478 if(fPedSubtract && fNewTime < nPreSamples) {
479 pedMean += in->GetSignal();
480 pedRMS += in->GetSignal()*in->GetSignal() ;
483 fSamples->AddAt(fNewAmp,iBin);
484 fTimes->AddAt(fNewTime,iBin);
486 //Mark that we already take it
490 // Fill array with samples
492 if(fPedSubtract && (in->GetTime() < nPreSamples)) {
493 pedMean += in->GetSignal();
494 pedRMS += in->GetSignal()*in->GetSignal() ;
497 fSamples->AddAt(in->GetSignal()-10,iBin);
498 fTimes->AddAt(in->GetTime(),iBin);
500 //Debug==============
501 // h->SetBinContent(in->GetTime(),in->GetSignal()) ;
502 //EndDebug==============
508 //_____________________________________________________________________________
509 void AliPHOSRawDecoderv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
511 // Number of parameters, Gradient, Chi squared, parameters, what to do
513 TList * toFit= (TList*)gMinuit->GetObjectFit() ;
514 TArrayD * params=(TArrayD*)toFit->At(0) ;
515 TArrayI * samples = (TArrayI*)toFit->At(1) ;
516 TArrayI * times = (TArrayI*)toFit->At(2) ;
520 for(Int_t iparam = 0 ; iparam < 3 ; iparam++)
521 Grad[iparam] = 0 ; // Will evaluate gradient
526 Double_t n=params->At(0) ;
527 Double_t alpha=params->At(1) ;
528 Double_t beta=params->At(3) ;
529 Double_t ped=params->At(4) ;
531 Double_t overflow=params->At(5) ;
532 Int_t iBin = (Int_t) params->At(6) ;
533 Int_t nSamples=TMath::Min(iBin+70,samples->GetSize()) ; //Here we set number of points to fit (70)
534 // iBin - first non-zero sample
535 Int_t tStep=times->At(iBin+1)-times->At(iBin) ;
536 Double_t ddt=times->At(iBin)-t0-tStep ;
537 Double_t exp1=TMath::Exp(-alpha*ddt) ;
538 Double_t exp2=TMath::Exp(-beta*ddt) ;
539 Double_t dexp1=TMath::Exp(-alpha*tStep) ;
540 Double_t dexp2=TMath::Exp(-beta*tStep) ;
541 for(Int_t i = iBin; i<nSamples ; i++) {
542 Double_t dt=double(times->At(i))-t0 ;
543 Double_t fsample = double(samples->At(i)) ;
544 if(fsample>=overflow)
554 Double_t dtn=TMath::Power(dt,n) ;
555 Double_t dtnE=dtn*exp1 ;
556 Double_t dt2E=dt*dt*exp2 ;
557 Double_t fit=ped+en*(dtnE + b*dt2E) ;
558 // if(fit>=overflow){
559 // diff=fsample-overflow ;
560 // fret += diff*diff ;
561 // //zero gradient here
564 diff = fsample - fit ;
566 if(iflag == 2){ // calculate gradient
567 Grad[0] += en*diff*(dtnE*(n/dt-alpha)+b*dt2E*(2./dt-beta)) ; //derivative over t0
568 Grad[1] -= diff*(dtnE+b*dt2E) ;
569 Grad[2] -= en*diff*dt2E ;
574 for(Int_t iparam = 0 ; iparam < 3 ; iparam++)
577 //-----------------------------------------------------------------------------
578 Double_t AliPHOSRawDecoderv1::Gamma2(Double_t dt,Double_t en,Double_t b,TArrayD * params){ //Function for fitting samples
580 //dt-time after start
582 //function parameters
584 Double_t ped=params->At(4) ;
586 return ped ; //pedestal
588 return ped+en*(TMath::Power(dt,params->At(0))*TMath::Exp(-dt*params->At(1))+b*dt*dt*TMath::Exp(-dt*params->At(3))) ;