]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSRawDecoderv1.cxx
3-par fitting for better quality estimate
[u/mrichter/AliRoot.git] / PHOS / AliPHOSRawDecoderv1.cxx
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 ---
42 #include "TList.h"
43 #include "TMath.h"
44 #include "TMinuit.h"
45
46 #include "TCanvas.h"
47 #include "TH1.h"
48 #include "TH2.h"
49 #include "TF1.h"
50 #include "TROOT.h"
51
52 // --- AliRoot header files ---
53 //#include "AliLog.h"
54 #include "AliPHOSRawDecoderv1.h"
55 #include "AliPHOSPulseGenerator.h"
56
57
58 ClassImp(AliPHOSRawDecoderv1)
59
60 //-----------------------------------------------------------------------------
61   AliPHOSRawDecoderv1::AliPHOSRawDecoderv1():AliPHOSRawDecoder(),
62 fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
63 {
64   //Default constructor.
65 }
66
67 //-----------------------------------------------------------------------------
68 AliPHOSRawDecoderv1::AliPHOSRawDecoderv1(AliRawReader* rawReader,  AliAltroMapping **mapping):
69   AliPHOSRawDecoder(rawReader,mapping),
70   fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
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);
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) ;
89   fToFit = new TList() ;
90 }
91
92 //-----------------------------------------------------------------------------
93 AliPHOSRawDecoderv1::~AliPHOSRawDecoderv1()
94 {
95   //Destructor.
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   }
108 }
109
110 //-----------------------------------------------------------------------------
111 AliPHOSRawDecoderv1::AliPHOSRawDecoderv1(const AliPHOSRawDecoderv1 &phosDecoder ):
112   AliPHOSRawDecoder(phosDecoder), 
113   fSampleParamsLow(0x0),fSampleParamsHigh(0x0),fToFit(0x0)
114 {
115   //Copy constructor.
116   fToFit = new TList() ;
117   fSampleParamsLow =new TArrayD(*(phosDecoder.fSampleParamsLow)) ;
118   fSampleParamsHigh=new TArrayD(*(phosDecoder.fSampleParamsHigh)) ;
119 }
120
121 //-----------------------------------------------------------------------------
122 AliPHOSRawDecoderv1& AliPHOSRawDecoderv1::operator = (const AliPHOSRawDecoderv1 &phosDecoder)
123 {
124   //Assignment operator.
125
126   //  if(this != &phosDecoder) {
127   //  }
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   }
137   return *this;
138 }
139
140 //-----------------------------------------------------------------------------
141 Bool_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.
149
150 //Debug=====================
151 //  TCanvas * c = 0; //(TCanvas*)gROOT->FindObjectAny("CSample") ;
152 //  if(!c)
153 //    c = new TCanvas("CSample","CSample") ;
154 // 
155 //  TH1D * h = 0 ; //(TH1D*)gROOT->FindObjectAny("hSample") ;
156 //  if(!h)
157 //    h=new TH1D("hSample","",200,0.,200.) ;
158 // 
159 //  TF1 * fff = 0 ; //(TF1*)gROOT->FindObjectAny("fff") ;
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===========
163   
164   AliCaloRawStream* in = fCaloStream;
165   
166   Int_t    iBin     = fSamples->GetSize() ;
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;
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
176   const Float_t maxEtoFit=5 ; //fit only samples above this energy, accept all samples (with good aRMS) below it
177   
178   while ( in->Next() ) { 
179     
180     if(!tLength) {
181       tLength = in->GetTimeLength();
182       if(tLength!=fSamples->GetSize()) {
183         delete fSamples ;
184         delete fTimes ;
185         fSamples = new TArrayI(tLength);
186         fTimes = new TArrayI(tLength);
187         iBin= fSamples->GetSize() ;
188       }
189       else{
190         fSamples->Reset() ;
191       }
192     }
193     
194     // Fit the full sample
195     if(in->IsNewHWAddress() && iBin != fSamples->GetSize()) {
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 
205       Double_t pedestal =0. ;
206       if(fPedSubtract){ 
207         if (nPed > 0)
208           pedestal = (Double_t)(pedMean/nPed); 
209         else
210           return kFALSE;
211       }
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. ;                                                                                                                   
218       Double_t wts=0 ;                                                                                                                       
219       Int_t tStart = 0 ;                                                                                                                   
220       for(Int_t i=iBin; i<fSamples->GetSize(); i++){
221         if(fSamples->At(i)>0){                                                                                                             
222           Double_t de=fSamples->At(i)-pedestal ;                                                                                           
223           if(de>1.){
224             aMean+=de*i ;                                                                                                                      
225             aRMS+=de*i*i ;                                                                                                                    
226             wts+=de; 
227           }                                                                                                                         
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       }
236       if(maxBin==fSamples->GetSize()-1){//bad "rising" sample
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 ;
247       }
248       
249       if(wts>0){
250         aMean/=wts; 
251         aRMS=aRMS/wts-aMean*aMean;
252       }
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
268       
269 //Debug:=====Draw sample
270 //if(fEnergy>pedestal+10.){
271 //if(fLowGainFlag && fEnergy>2){
272 //  if(!c)
273 //    if(!fLowGainFlag && fRow==32 && fColumn==18){
274 //    TCanvas *c = new TCanvas("CSample","CSample") ;
275 //    c->cd() ;
276 //    h->Draw() ;
277 //    c->Update() ;
278 // printf("fEnergy=%f, aRMS=%f \n",fEnergy,aRMS) ;   
279 //getchar() ;
280 //}
281 //======================
282
283       //IF sample has reasonable mean and RMS, try to fit it with gamma2
284         
285         gMinuit->mncler();                     // Reset Minuit's list of paramters
286         gMinuit->SetPrintLevel(-1) ;           // No Printout
287         gMinuit->SetFCN(AliPHOSRawDecoderv1::UnfoldingChiSquare) ;  
288         // To set the address of the minimization function 
289         
290        fToFit->Clear() ;
291        Double_t b,bmin,bmax ;
292        if(fLowGainFlag){
293          fSampleParamsLow->AddAt(pedestal,4) ;
294          if(fOverflow)
295            fSampleParamsLow->AddAt(double(maxAmp),5) ;
296          else
297            fSampleParamsLow->AddAt(double(1023),5) ;
298          fSampleParamsLow->AddAt(double(iBin),6) ;
299          fToFit->AddFirst((TObject*)fSampleParamsLow) ; 
300          b=fSampleParamsLow->At(2) ;
301          bmin=0.5 ;
302          bmax=10. ;
303        }
304        else{
305          fSampleParamsHigh->AddAt(pedestal,4) ;
306          if(fOverflow)
307            fSampleParamsHigh->AddAt(double(maxAmp),5) ;
308          else
309            fSampleParamsHigh->AddAt(double(1023),5);
310          fSampleParamsHigh->AddAt(double(iBin),6);
311          fToFit->AddFirst((TObject*)fSampleParamsHigh) ; 
312          b=fSampleParamsHigh->At(2) ;
313          bmin=0.05 ;
314          bmax=0.4 ;
315         }
316         fToFit->AddLast((TObject*)fSamples) ;
317         fToFit->AddLast((TObject*)fTimes) ;
318
319         gMinuit->SetObjectFit((TObject*)fToFit) ;         // To tranfer pointer to UnfoldingChiSquare
320         Int_t ierflg ;
321         gMinuit->mnparm(0, "t0",  1.*tStart, 0.01, 0, 0, ierflg) ;
322         if(ierflg != 0){
323 //        AliWarning(Form("Unable to set initial value for fit procedure : t0=%e\n",1.*tStart) ) ;
324           fEnergy=0. ;
325           fTime=-999. ;
326           fQuality=999 ;
327           return kTRUE ; //will scan further
328         }
329         Double_t amp0; 
330         if(fLowGainFlag)
331           amp0=fEnergy/sampleMaxLG;
332         else
333           amp0=fEnergy/sampleMaxHG;
334
335         gMinuit->mnparm(1, "Energy", amp0 , 0.01*amp0, 0, 0, ierflg) ;
336         if(ierflg != 0){
337 //        AliWarning(Form("Unable to set initial value for fit procedure : E=%e\n", amp0)) ;
338           fEnergy=0. ;
339           fTime=-999. ;
340           fQuality=999 ;
341           return kTRUE ; //will scan further
342         }
343
344         gMinuit->mnparm(2, "p2", b, 0.01*b, bmin, bmax, ierflg) ;
345         if(ierflg != 0){                                         
346 //        AliWarning(Form("Unable to set initial value for fit procedure : E=%e\n", amp0)) ;  
347           fEnergy=0. ;           
348           fTime=-999. ;         
349           fQuality=999 ;       
350           return kTRUE ; //will scan further  
351         }             
352  
353         
354         Double_t p0 = 0.0001 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
355         //  depends on it. 
356         Double_t p1 = 1.0 ;
357         Double_t p2 = 0.0 ;
358         gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls  
359         gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient  
360         //      gMinuit->SetMaxIterations(100);
361         gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
362         
363         gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize 
364
365         
366         Double_t err,t0err ;
367         Double_t t0,efit ;
368         gMinuit->GetParameter(0,t0, t0err) ;    
369         gMinuit->GetParameter(1,efit, err) ;    
370
371         Double_t bfit, berr ;
372         gMinuit->GetParameter(2,bfit,berr) ;
373
374         //Calculate total energy
375         //this isparameterization of depetendence of pulse height on parameter b
376         if(fLowGainFlag)
377           efit*=99.54910 + 78.65038*bfit ;
378         else
379           efit*=80.33109+128.6433*bfit ;
380
381         if(efit<0. || efit > 10000.){                                                                          
382           fEnergy=0 ; //bad sample                                                    
383           fTime=-999.;                                                                
384           fQuality=999 ;                                                              
385           return kTRUE;
386         }                                                                             
387  
388         //evaluate fit quality
389         Double_t fmin,fedm,errdef ;
390         Int_t npari,nparx,istat;
391         gMinuit->mnstat(fmin,fedm,errdef,npari,nparx,istat) ;
392         fQuality=fmin/(fSamples->GetSize()-iBin) ;
393         //compare quality with some parameterization
394         if(fLowGainFlag){
395           fQuality/=2.+0.002*fEnergy ;
396         }
397         else{
398           fQuality/=0.75+0.0025*fEnergy ;
399         }
400
401 //Debug================
402 //        Double_t n,alpha,beta ;
403 //        if(fLowGainFlag){
404 //          n=fSampleParamsLow->At(0) ;
405 //          alpha=fSampleParamsLow->At(1) ;
406 //          beta=fSampleParamsLow->At(3) ;
407 //        }
408 //        else{
409 //          n=fSampleParamsHigh->At(0) ;
410 //          alpha=fSampleParamsHigh->At(1) ;
411 //          beta=fSampleParamsHigh->At(3) ;
412 //        }
413 //
414 //    if( fQuality > 1 && fEnergy > 20. && !fOverflow){
415 //    if(!fLowGainFlag && fRow==32 && fColumn==18){
416 //    printf("Col=%d, row=%d, qual=%f, E=%f, t0=%f, b=%f\n",fColumn,fRow,fQuality,efit,t0,bfit) ;
417 //    TCanvas * c = new TCanvas("samp") ;
418 //    c->cd() ;
419 //    h->Draw() ;
420 //    if(fLowGainFlag){
421 //      fff->SetParameters(pedestal,efit/sampleMaxLG,t0,n,alpha,bfit,beta) ;
422 //    }
423 //    else{
424 //     fff->SetParameters(pedestal,efit/sampleMaxHG,t0,n,alpha,bfit,beta) ;
425 //    }
426 ////    for(Int_t i=1;i<=h->GetNbinsX(); i++){
427 //       Double_t x=h->GetBinCenter(i) ;
428 //       h->SetBinContent(i,h->GetBinContent(i)-fff->Eval(x)) ;
429 //    }
430 //    h->SetMinimum(-15.) ;
431 //    h->SetMaximum(15.) ;
432 //    h->Draw() ;
433 //    fff->Draw("same") ;
434 //    c->Update();
435 //    getchar() ;
436 //    }
437 //====================
438
439       fEnergy=efit ;
440       fTime=t0 ;
441
442       
443       return kTRUE;
444     }
445     
446     fLowGainFlag = in->IsLowGain();
447     fModule = in->GetModule()+1;
448     fRow    = in->GetRow()   +1;
449     fColumn = in->GetColumn()+1;
450
451     //add previouly taken if coincides
452     if(fLowGainFlag==fNewLowGainFlag && fModule==fNewModule &&
453        fRow==fNewRow && fColumn==fNewColumn){
454        iBin--;                                                                                                                                
455        if(fPedSubtract && fNewTime < nPreSamples) {                                                                                    
456          pedMean += in->GetSignal();                                                                                                          
457          nPed++;                                                                                                                              
458        }                                                                                                                                      
459        fSamples->AddAt(fNewAmp,iBin);                                                                                                 
460        fTimes->AddAt(fNewTime,iBin);                                                                                                     
461     
462        //Mark that we already take it
463        fNewModule=-1 ;
464     }
465     
466     // Fill array with samples
467     iBin--;
468     if(fPedSubtract && (in->GetTime() < nPreSamples)) {
469       pedMean += in->GetSignal();
470       nPed++;
471     }
472     fSamples->AddAt(in->GetSignal(),iBin);
473     fTimes->AddAt(in->GetTime(),iBin);
474  
475 //Debug==============
476 //    h->SetBinContent(iBin,in->GetSignal()) ;
477     
478   } // in.Next()
479   
480   return kFALSE;
481 }
482 //_____________________________________________________________________________
483 void AliPHOSRawDecoderv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
484 {
485   // Number of parameters, Gradient, Chi squared, parameters, what to do
486
487   TList * toFit= (TList*)gMinuit->GetObjectFit() ;
488   TArrayD * params=(TArrayD*)toFit->At(0) ; 
489   TArrayI * samples = (TArrayI*)toFit->At(1) ;
490   TArrayI * times = (TArrayI*)toFit->At(2) ;
491
492   fret = 0. ;     
493   if(iflag == 2)
494     for(Int_t iparam = 0 ; iparam < 3 ; iparam++)    
495       Grad[iparam] = 0 ; // Will evaluate gradient
496   
497   Double_t t0=x[0] ;
498   Double_t en=x[1] ;
499   Double_t b=x[2] ;
500   Double_t n=params->At(0) ;
501   Double_t alpha=params->At(1) ;
502   Double_t beta=params->At(3) ;
503   Double_t ped=params->At(4) ;
504
505   Double_t overflow=params->At(5) ;
506   Int_t iBin = (Int_t) params->At(6) ;
507   Int_t nSamples=TMath::Min(iBin+70,samples->GetSize()) ; //Here we set number of points to fit (70)
508   
509   Int_t tStep=times->At(iBin+1)-times->At(iBin) ;
510   Double_t ddt=TMath::Ceil(t0)-t0-tStep ;
511   Double_t exp1=TMath::Exp(-alpha*ddt) ;
512   Double_t exp2=TMath::Exp(-beta*ddt) ;
513   Double_t dexp1=TMath::Exp(-alpha*tStep) ;
514   Double_t dexp2=TMath::Exp(-beta*tStep) ;
515   for(Int_t i = iBin; i<nSamples ; i++) {
516     Double_t fsample = double(samples->At(i)) ;
517 //    if(fsample>=overflow)
518 //      continue ;
519     Double_t dt=double(times->At(i))-t0 ;
520     Double_t diff ;
521     if(dt<=0.){
522       diff=fsample - ped ; 
523       fret += diff*diff ;
524       continue ;
525     }
526     exp1*=dexp1 ;
527     exp2*=dexp2 ;
528     Double_t dtn=TMath::Power(dt,n) ;
529     Double_t dtnE=dtn*exp1 ;
530     Double_t dt2E=dt*dt*exp2 ;
531     Double_t fit=ped+en*(dtnE + b*dt2E) ;
532     if(fit>=overflow){
533       diff=fsample-overflow ;
534       fret += diff*diff ;
535       //zero gradient here
536     }
537     else{
538       diff = fsample - fit ;
539       fret += diff*diff ;
540       if(iflag == 2){  // calculate gradient
541         Grad[0] += en*diff*(dtnE*(n/dt-alpha)+b*dt2E*(2./dt-beta))  ; //derivative over t0
542         Grad[1] -= diff*(dtnE+b*dt2E) ;
543         Grad[2] -= en*diff*dt2E ;
544       }
545     }
546   }
547   if(iflag == 2)
548     for(Int_t iparam = 0 ; iparam < 3 ; iparam++)    
549       Grad[iparam] *= 2. ; 
550 }
551 //-----------------------------------------------------------------------------
552 Double_t AliPHOSRawDecoderv1::Gamma2(Double_t dt,Double_t en,Double_t b,TArrayD * params){  //Function for fitting samples
553   //parameters:
554   //dt-time after start
555   //en-amplutude
556   //function parameters
557   
558   Double_t ped=params->At(4) ;
559   if(dt<0.)
560     return ped ; //pedestal
561   else
562     return ped+en*(TMath::Power(dt,params->At(0))*TMath::Exp(-dt*params->At(1))+b*dt*dt*TMath::Exp(-dt*params->At(3))) ;
563 }
564