]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSRawFitterv1.cxx
Update in the T0 configuration @P2.
[u/mrichter/AliRoot.git] / PHOS / AliPHOSRawFitterv1.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 extracts the signal parameters (energy, time, quality)
19 // from ALTRO samples. Energy is in ADC counts, time is in time bin units.
20 // A fitting algorithm evaluates the energy and the time from Minuit minimization
21 // 
22 // Typical use case:
23 //     AliPHOSRawFitter *fitter=new AliPHOSRawFitter();
24 //     fitter->SetChannelGeo(module,cellX,cellZ,caloFlag);
25 //     fitter->SetCalibData(fgCalibData) ;
26 //     fitter->Eval(sig,sigStart,sigLength);
27 //     Double_t amplitude = fitter.GetEnergy();
28 //     Double_t time      = fitter.GetTime();
29 //     Bool_t   isLowGain = fitter.GetCaloFlag()==0;
30
31 // Author: Dmitri Peressounko (Oct.2008)
32 // Modified: Yuri Kharlov (Jul.2009)
33
34 // --- ROOT system ---
35 #include "TArrayD.h"
36 #include "TArrayI.h"
37 #include "TList.h"
38 #include "TMath.h"
39 #include "TMinuit.h"
40
41 // --- AliRoot header files ---
42 #include "AliLog.h"
43 #include "AliPHOSCalibData.h"
44 #include "AliPHOSRawFitterv1.h"
45 #include "AliPHOSPulseGenerator.h"
46
47 ClassImp(AliPHOSRawFitterv1)
48
49 //-----------------------------------------------------------------------------
50 AliPHOSRawFitterv1::AliPHOSRawFitterv1():
51   AliPHOSRawFitterv0(),
52   fSampleParamsLow(0x0),
53   fSampleParamsHigh(0x0),
54   fToFit(0x0)
55 {
56   //Default constructor.
57   if(!gMinuit) 
58     gMinuit = new TMinuit(100);
59   fSampleParamsHigh =new TArrayD(7) ;
60   fSampleParamsHigh->AddAt(2.174,0) ;
61   fSampleParamsHigh->AddAt(0.106,1) ;
62   fSampleParamsHigh->AddAt(0.173,2) ;
63   fSampleParamsHigh->AddAt(0.06106,3) ;
64   //last two parameters are pedestal and overflow
65   fSampleParamsLow=new TArrayD(7) ;
66   fSampleParamsLow->AddAt(2.456,0) ;
67   fSampleParamsLow->AddAt(0.137,1) ;
68   fSampleParamsLow->AddAt(2.276,2) ;
69   fSampleParamsLow->AddAt(0.08246,3) ;
70   fToFit = new TList() ;
71 }
72
73 //-----------------------------------------------------------------------------
74 AliPHOSRawFitterv1::~AliPHOSRawFitterv1()
75 {
76   //Destructor.
77   if(fSampleParamsLow){
78     delete fSampleParamsLow ; 
79     fSampleParamsLow=0 ;
80   }
81   if(fSampleParamsHigh){
82     delete fSampleParamsHigh ;
83     fSampleParamsHigh=0;
84   }
85   if(fToFit){
86     delete fToFit ;
87     fToFit=0 ;
88   }
89 }
90
91 //-----------------------------------------------------------------------------
92 AliPHOSRawFitterv1::AliPHOSRawFitterv1(const AliPHOSRawFitterv1 &phosFitter ):
93   AliPHOSRawFitterv0(phosFitter), 
94   fSampleParamsLow(0x0),
95   fSampleParamsHigh(0x0),
96   fToFit(0x0)
97 {
98   //Copy constructor.
99   fToFit = new TList() ;
100   fSampleParamsLow =new TArrayD(*(phosFitter.fSampleParamsLow)) ;
101   fSampleParamsHigh=new TArrayD(*(phosFitter.fSampleParamsHigh)) ;
102 }
103
104 //-----------------------------------------------------------------------------
105 AliPHOSRawFitterv1& AliPHOSRawFitterv1::operator = (const AliPHOSRawFitterv1 &phosFitter)
106 {
107   //Assignment operator.
108
109   fToFit = new TList() ;
110   if(fSampleParamsLow){
111     fSampleParamsLow = phosFitter.fSampleParamsLow ;
112     fSampleParamsHigh= phosFitter.fSampleParamsHigh ;
113   }
114   else{
115     fSampleParamsLow =new TArrayD(*(phosFitter.fSampleParamsLow)) ; 
116     fSampleParamsHigh=new TArrayD(*(phosFitter.fSampleParamsHigh)) ;
117   }
118   return *this;
119 }
120
121 //-----------------------------------------------------------------------------
122 Bool_t AliPHOSRawFitterv1::Eval(const UShort_t *signal, Int_t sigStart, Int_t sigLength)
123 {
124   //Extract an energy deposited in the crystal,
125   //crystal' position (module,column,row),
126   //time and gain (high or low).
127   //First collects sample, then evaluates it and if it has
128   //reasonable shape, fits it with Gamma2 function and extracts 
129   //energy and time.
130
131   if (fCaloFlag == 2 || fNBunches > 1) {
132     fQuality = 1000;
133     return kTRUE;
134   }
135
136   Float_t pedMean = 0;
137   Float_t pedRMS  = 0;
138   Int_t   nPed    = 0;
139   const Float_t kBaseLine   = 1.0;
140   const Int_t   kPreSamples = 10;
141   
142   TArrayI *fSamples = new TArrayI(sigLength); // array of sample amplitudes
143   TArrayI *fTimes   = new TArrayI(sigLength); // array of sample time stamps
144   for (Int_t i=0; i<sigLength; i++) {
145     if (i<kPreSamples) {
146       nPed++;
147       pedMean += signal[i];
148       pedRMS  += signal[i]*signal[i] ;
149     }
150     fSamples->AddAt(signal[i],sigLength-i-1);
151     fTimes  ->AddAt(i ,i);
152   }
153
154   fEnergy = -111;
155   fQuality= 999. ;
156   const Float_t sampleMaxHG=102.332 ;  //maximal height of HG sample with given parameterization
157   const Float_t sampleMaxLG=277.196 ;  //maximal height of LG sample with given parameterization
158   const Float_t maxEtoFit=5 ; //fit only samples above this energy, accept all samples (with good aRMS) below it
159   Double_t pedestal = 0;
160
161   if (fPedSubtract) {
162     if (nPed > 0) {
163       fPedestalRMS=(pedRMS - pedMean*pedMean/nPed)/nPed ;
164       if(fPedestalRMS > 0.) 
165         fPedestalRMS = TMath::Sqrt(fPedestalRMS) ;
166       fEnergy -= (Double_t)(pedMean/nPed); // pedestal subtraction
167     }
168     else
169       return kFALSE;
170   }
171   else {
172     //take pedestals from DB
173     pedestal = (Double_t) fAmpOffset ;
174     if (fCalibData) {
175       Float_t truePed       = fCalibData->GetADCpedestalEmc(fModule, fCellZ, fCellX) ;
176       Int_t   altroSettings = fCalibData->GetAltroOffsetEmc(fModule, fCellZ, fCellX) ;
177       pedestal += truePed - altroSettings ;
178     }
179     else{
180       AliWarning(Form("Can not read data from OCDB")) ;
181     }
182     fEnergy-=pedestal ;
183   }
184
185   if (fEnergy < kBaseLine) fEnergy = 0;
186   //Evaluate time
187   Int_t iStart = 0;
188   while(iStart<sigLength && fSamples->At(iStart)-pedestal <kBaseLine) iStart++ ;
189   fTime = sigStart-sigLength+iStart; 
190   
191   //calculate time and energy
192   Int_t    maxBin=0 ;
193   Int_t    maxAmp=0 ;
194   Double_t aMean =0. ;
195   Double_t aRMS  =0. ;
196   Double_t wts   =0 ;
197   Int_t    tStart=0 ;
198
199   for (Int_t i=0; i<sigLength; i++){
200     if(signal[i] > pedestal){
201       Double_t de = signal[i] - pedestal ;
202       if(de > 1.) {
203         aMean += de*i ;
204         aRMS  += de*i*i ;
205         wts   += de; 
206       }
207       if(de > 2 && tStart==0) 
208         tStart = i ;
209       if(maxAmp < signal[i]){
210         maxBin = i ;
211         maxAmp = signal[i] ;
212       }
213     }
214   }
215
216   if (maxBin==sigLength-1){//bad "rising" sample
217     fEnergy =    0. ;
218     fTime   = -999. ;
219     fQuality=  999. ;
220     return kTRUE ;
221   }
222
223   fEnergy=Double_t(maxAmp)-pedestal ;
224   fOverflow =0 ;  //look for plato on the top of sample
225   if (fEnergy>500 &&  //this is not fluctuation of soft sample
226      maxBin<sigLength-1 && fSamples->At(maxBin+1)==maxAmp){ //and there is a plato
227     fOverflow = kTRUE ;
228   }
229   
230   if (wts > 0) {
231     aMean /= wts; 
232     aRMS   = aRMS/wts - aMean*aMean;
233   }
234
235   //do not take too small energies
236   if (fEnergy < kBaseLine) 
237     fEnergy = 0;
238   
239   //do not test quality of too soft samples
240   if (fEnergy < maxEtoFit){
241     fTime = tStart;
242     if (aRMS < 2.) //sigle peak
243       fQuality = 999. ;
244     else
245       fQuality =   0. ;
246     return kTRUE ;
247   }
248       
249   // if sample has reasonable mean and RMS, try to fit it with gamma2
250   
251   gMinuit->mncler();                     // Reset Minuit's list of paramters
252   gMinuit->SetPrintLevel(-1) ;           // No Printout
253   gMinuit->SetFCN(AliPHOSRawFitterv1::UnfoldingChiSquare) ;  
254   // To set the address of the minimization function 
255   
256   fToFit->Clear("nodelete") ;
257   Double_t b=0,bmin=0,bmax=0 ;
258   if      (fCaloFlag == 0){ // Low gain
259     fSampleParamsLow->AddAt(pedestal,4) ;
260     if (fOverflow)
261       fSampleParamsLow->AddAt(double(maxAmp),5) ;
262     else
263       fSampleParamsLow->AddAt(double(1023),5) ;
264     fSampleParamsLow->AddAt(double(iStart),6) ;
265     fToFit->AddFirst((TObject*)fSampleParamsLow) ; 
266     b=fSampleParamsLow->At(2) ;
267     bmin=0.5 ;
268     bmax=10. ;
269   }
270   else if (fCaloFlag == 1){ // High gain
271     fSampleParamsHigh->AddAt(pedestal,4) ;
272     if (fOverflow)
273       fSampleParamsHigh->AddAt(double(maxAmp),5) ;
274     else
275       fSampleParamsHigh->AddAt(double(1023),5);
276     fSampleParamsHigh->AddAt(double(iStart),6);
277     fToFit->AddFirst((TObject*)fSampleParamsHigh) ; 
278     b=fSampleParamsHigh->At(2) ;
279     bmin=0.05 ;
280     bmax=0.4 ;
281   }
282   fToFit->AddLast((TObject*)fSamples) ;
283   fToFit->AddLast((TObject*)fTimes) ;
284   
285   gMinuit->SetObjectFit((TObject*)fToFit) ;         // To tranfer pointer to UnfoldingChiSquare
286   Int_t ierflg ;
287   gMinuit->mnparm(0, "t0",  1.*tStart, 0.01, -500., 500., ierflg) ;
288   if(ierflg != 0){
289     //    AliWarning(Form("Unable to set initial value for fit procedure : t0=%e\n",1.*tStart) ) ;
290     fEnergy =   0. ;
291     fTime   =-999. ;
292     fQuality= 999. ;
293     return kTRUE ; //will scan further
294   }
295   Double_t amp0=0; 
296   if      (fCaloFlag == 0) // Low gain
297     amp0 = fEnergy/sampleMaxLG;
298   else if (fCaloFlag == 1) // High gain
299     amp0 = fEnergy/sampleMaxHG;
300   
301   gMinuit->mnparm(1, "Energy", amp0 , 0.01*amp0, 0, 0, ierflg) ;
302   if(ierflg != 0){
303     //    AliWarning(Form("Unable to set initial value for fit procedure : E=%e\n", amp0)) ;
304     fEnergy =   0. ;
305     fTime   =-999. ;
306     fQuality= 999. ;
307     return kTRUE ; //will scan further
308   }
309   
310   gMinuit->mnparm(2, "p2", b, 0.01*b, bmin, bmax, ierflg) ;
311   if(ierflg != 0){                                         
312     //        AliWarning(Form("Unable to set initial value for fit procedure : E=%e\n", amp0)) ;  
313     fEnergy =   0. ;
314     fTime   =-999. ;
315     fQuality= 999. ;
316     return kTRUE ; //will scan further  
317   }             
318   
319   Double_t p0 = 0.0001 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
320   //  depends on it. 
321   Double_t p1 = 1.0 ;
322   Double_t p2 = 0.0 ;
323   gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls  
324   gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient  
325   //    gMinuit->SetMaxIterations(100);
326   gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
327   
328   gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize 
329   
330   Double_t err,t0err ;
331   Double_t t0,efit ;
332   gMinuit->GetParameter(0,t0, t0err) ;
333   gMinuit->GetParameter(1,efit, err) ;
334   
335   Double_t bfit, berr ;
336   gMinuit->GetParameter(2,bfit,berr) ;
337   
338   //Calculate total energy
339   //this is parameterization of dependence of pulse height on parameter b
340   if(fCaloFlag == 0) // Low gain
341     efit *= 99.54910 + 78.65038*bfit ;
342   else if(fCaloFlag == 1) // High gain
343     efit *= 80.33109 + 128.6433*bfit ;
344   
345   if(efit < 0. || efit > 10000.){
346     //set energy to previously found max
347     fTime   =-999.;
348     fQuality= 999 ;
349     return kTRUE;
350   }                                                                             
351   
352   //evaluate fit quality
353   Double_t fmin,fedm,errdef ;
354   Int_t npari,nparx,istat;
355   gMinuit->mnstat(fmin,fedm,errdef,npari,nparx,istat) ;
356   fQuality = fmin/sigLength ;
357   //compare quality with some parameterization
358   if      (fCaloFlag == 0) // Low gain
359     fQuality /= 2.00 + 0.0020*fEnergy ;
360   else if (fCaloFlag == 1) // High gain
361     fQuality /= 0.75 + 0.0025*fEnergy ;
362   
363   fEnergy = efit ;
364   fTime  += t0 - 4.024*bfit ; //-10.402*bfit+4.669*bfit*bfit ; //Correction for 70 samples
365 //  fTime  += sigStart;
366   
367   delete fSamples ;
368   delete fTimes ;
369   return kTRUE;
370 }
371 //_____________________________________________________________________________
372 void AliPHOSRawFitterv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
373 {
374   // Number of parameters, Gradient, Chi squared, parameters, what to do
375
376   TList * toFit= (TList*)gMinuit->GetObjectFit() ;
377   TArrayD * params=(TArrayD*)toFit->At(0) ; 
378   TArrayI * samples = (TArrayI*)toFit->At(1) ;
379   TArrayI * times = (TArrayI*)toFit->At(2) ;
380
381   fret = 0. ;     
382   if(iflag == 2)
383     for(Int_t iparam = 0 ; iparam < 3 ; iparam++)    
384       Grad[iparam] = 0 ; // Will evaluate gradient
385   
386   Double_t t0=x[0] ;
387   Double_t en=x[1] ;
388   Double_t b=x[2] ;
389   Double_t n=params->At(0) ;
390   Double_t alpha=params->At(1) ;
391   Double_t beta=params->At(3) ;
392   Double_t ped=params->At(4) ;
393
394   Double_t overflow=params->At(5) ;
395   Int_t iBin = (Int_t) params->At(6) ;
396   Int_t nSamples=TMath::Min(iBin+70,samples->GetSize()) ; //Here we set number of points to fit (70)
397   // iBin - first non-zero sample 
398   Int_t tStep=times->At(iBin+1)-times->At(iBin) ;
399   Double_t ddt=times->At(iBin)-t0-tStep ;
400   Double_t exp1=TMath::Exp(-alpha*ddt) ;
401   Double_t exp2=TMath::Exp(-beta*ddt) ;
402   Double_t dexp1=TMath::Exp(-alpha*tStep) ;
403   Double_t dexp2=TMath::Exp(-beta*tStep) ;
404   for(Int_t i = iBin; i<nSamples ; i++) {
405     Double_t dt=double(times->At(i))-t0 ;
406     Double_t fsample = double(samples->At(i)) ;
407     if(fsample>=overflow)
408       continue ;
409     Double_t diff ;
410     exp1*=dexp1 ;
411     exp2*=dexp2 ;
412     if(dt<=0.){
413       diff=fsample - ped ; 
414       fret += diff*diff ;
415       continue ;
416     }
417     Double_t dtn=TMath::Power(dt,n) ;
418     Double_t dtnE=dtn*exp1 ;
419     Double_t dt2E=dt*dt*exp2 ;
420     Double_t fit=ped+en*(dtnE + b*dt2E) ;
421     diff = fsample - fit ;
422     fret += diff*diff ;
423     if(iflag == 2){  // calculate gradient
424       Grad[0] += en*diff*(dtnE*(n/dt-alpha)+b*dt2E*(2./dt-beta))  ; //derivative over t0
425       Grad[1] -= diff*(dtnE+b*dt2E) ;
426       Grad[2] -= en*diff*dt2E ;
427     }
428   }
429   if(iflag == 2)
430     for(Int_t iparam = 0 ; iparam < 3 ; iparam++)    
431       Grad[iparam] *= 2. ; 
432 }
433 //-----------------------------------------------------------------------------
434 Double_t AliPHOSRawFitterv1::Gamma2(Double_t dt,Double_t en,Double_t b,TArrayD * params){  //Function for fitting samples
435   //parameters:
436   //dt-time after start
437   //en-amplutude
438   //function parameters
439   
440   Double_t ped=params->At(4) ;
441   if(dt<0.)
442     return ped ; //pedestal
443   else
444     return ped+en*(TMath::Power(dt,params->At(0))*TMath::Exp(-dt*params->At(1))+b*dt*dt*TMath::Exp(-dt*params->At(3))) ;
445 }