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