]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSRawFitterv1.cxx
Another histos for lumi
[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();
379c5c09 24// fitter->SetChannelGeo(module,cellX,cellZ,caloFlag);
25// fitter->SetCalibData(fgCalibData) ;
1dfadc32 26// fitter->Eval(sig,sigStart,sigLength);
379c5c09 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 ---
33373c87 35#include "TArrayI.h"
379c5c09 36#include "TList.h"
37#include "TMath.h"
38#include "TMinuit.h"
379c5c09 39
40// --- AliRoot header files ---
41#include "AliLog.h"
42#include "AliPHOSCalibData.h"
43#include "AliPHOSRawFitterv1.h"
44#include "AliPHOSPulseGenerator.h"
45
46ClassImp(AliPHOSRawFitterv1)
47
48//-----------------------------------------------------------------------------
49AliPHOSRawFitterv1::AliPHOSRawFitterv1():
50 AliPHOSRawFitterv0(),
51 fSampleParamsLow(0x0),
52 fSampleParamsHigh(0x0),
53 fToFit(0x0)
54{
55 //Default constructor.
56 if(!gMinuit)
57 gMinuit = new TMinuit(100);
58 fSampleParamsHigh =new TArrayD(7) ;
59 fSampleParamsHigh->AddAt(2.174,0) ;
60 fSampleParamsHigh->AddAt(0.106,1) ;
61 fSampleParamsHigh->AddAt(0.173,2) ;
62 fSampleParamsHigh->AddAt(0.06106,3) ;
63 //last two parameters are pedestal and overflow
64 fSampleParamsLow=new TArrayD(7) ;
65 fSampleParamsLow->AddAt(2.456,0) ;
66 fSampleParamsLow->AddAt(0.137,1) ;
67 fSampleParamsLow->AddAt(2.276,2) ;
68 fSampleParamsLow->AddAt(0.08246,3) ;
69 fToFit = new TList() ;
70}
71
72//-----------------------------------------------------------------------------
73AliPHOSRawFitterv1::~AliPHOSRawFitterv1()
74{
75 //Destructor.
7dac04a3 76 //Destructor
379c5c09 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//-----------------------------------------------------------------------------
92AliPHOSRawFitterv1::AliPHOSRawFitterv1(const AliPHOSRawFitterv1 &phosFitter ):
7dac04a3 93 AliPHOSRawFitterv0(phosFitter),
379c5c09 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//-----------------------------------------------------------------------------
105AliPHOSRawFitterv1& AliPHOSRawFitterv1::operator = (const AliPHOSRawFitterv1 &phosFitter)
106{
107 //Assignment operator.
7dac04a3 108 if(this != &phosFitter) {
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 }
379c5c09 118 }
119 return *this;
120}
121
122//-----------------------------------------------------------------------------
92236b27 123Bool_t AliPHOSRawFitterv1::Eval(const UShort_t *signal, Int_t sigStart, Int_t sigLength)
379c5c09 124{
125 //Extract an energy deposited in the crystal,
126 //crystal' position (module,column,row),
127 //time and gain (high or low).
128 //First collects sample, then evaluates it and if it has
129 //reasonable shape, fits it with Gamma2 function and extracts
130 //energy and time.
131
92236b27 132 if (fCaloFlag == 2 || fNBunches > 1) {
379c5c09 133 fQuality = 1000;
134 return kTRUE;
135 }
136
137 Float_t pedMean = 0;
138 Float_t pedRMS = 0;
139 Int_t nPed = 0;
140 const Float_t kBaseLine = 1.0;
141 const Int_t kPreSamples = 10;
142
92236b27 143 TArrayI *fSamples = new TArrayI(sigLength); // array of sample amplitudes
144 TArrayI *fTimes = new TArrayI(sigLength); // array of sample time stamps
145 for (Int_t i=0; i<sigLength; i++) {
379c5c09 146 if (i<kPreSamples) {
147 nPed++;
92236b27 148 pedMean += signal[i];
149 pedRMS += signal[i]*signal[i] ;
379c5c09 150 }
33373c87 151 fSamples->AddAt(signal[i],sigLength-i-1);
152 fTimes ->AddAt(i ,i);
379c5c09 153 }
154
379c5c09 155 fEnergy = -111;
156 fQuality= 999. ;
157 const Float_t sampleMaxHG=102.332 ; //maximal height of HG sample with given parameterization
158 const Float_t sampleMaxLG=277.196 ; //maximal height of LG sample with given parameterization
159 const Float_t maxEtoFit=5 ; //fit only samples above this energy, accept all samples (with good aRMS) below it
160 Double_t pedestal = 0;
161
162 if (fPedSubtract) {
163 if (nPed > 0) {
164 fPedestalRMS=(pedRMS - pedMean*pedMean/nPed)/nPed ;
165 if(fPedestalRMS > 0.)
166 fPedestalRMS = TMath::Sqrt(fPedestalRMS) ;
167 fEnergy -= (Double_t)(pedMean/nPed); // pedestal subtraction
168 }
169 else
170 return kFALSE;
171 }
172 else {
173 //take pedestals from DB
174 pedestal = (Double_t) fAmpOffset ;
175 if (fCalibData) {
176 Float_t truePed = fCalibData->GetADCpedestalEmc(fModule, fCellZ, fCellX) ;
177 Int_t altroSettings = fCalibData->GetAltroOffsetEmc(fModule, fCellZ, fCellX) ;
178 pedestal += truePed - altroSettings ;
179 }
180 else{
181 AliWarning(Form("Can not read data from OCDB")) ;
182 }
183 fEnergy-=pedestal ;
184 }
185
186 if (fEnergy < kBaseLine) fEnergy = 0;
33373c87 187 //Evaluate time
188 Int_t iStart = 0;
189 while(iStart<sigLength && fSamples->At(iStart)-pedestal <kBaseLine) iStart++ ;
190 fTime = sigStart-sigLength+iStart;
379c5c09 191
192 //calculate time and energy
92236b27 193 Int_t maxBin=0 ;
194 Int_t maxAmp=0 ;
195 Double_t aMean =0. ;
196 Double_t aRMS =0. ;
197 Double_t wts =0 ;
198 Int_t tStart=0 ;
199
200 for (Int_t i=0; i<sigLength; i++){
201 if(signal[i] > pedestal){
202 Double_t de = signal[i] - pedestal ;
203 if(de > 1.) {
204 aMean += de*i ;
205 aRMS += de*i*i ;
206 wts += de;
379c5c09 207 }
92236b27 208 if(de > 2 && tStart==0)
209 tStart = i ;
210 if(maxAmp < signal[i]){
211 maxBin = i ;
212 maxAmp = signal[i] ;
379c5c09 213 }
214 }
215 }
92236b27 216
217 if (maxBin==sigLength-1){//bad "rising" sample
218 fEnergy = 0. ;
219 fTime = -999. ;
220 fQuality= 999. ;
379c5c09 221 return kTRUE ;
222 }
92236b27 223
379c5c09 224 fEnergy=Double_t(maxAmp)-pedestal ;
225 fOverflow =0 ; //look for plato on the top of sample
92236b27 226 if (fEnergy>500 && //this is not fluctuation of soft sample
227 maxBin<sigLength-1 && fSamples->At(maxBin+1)==maxAmp){ //and there is a plato
379c5c09 228 fOverflow = kTRUE ;
229 }
92236b27 230
231 if (wts > 0) {
232 aMean /= wts;
233 aRMS = aRMS/wts - aMean*aMean;
379c5c09 234 }
235
236 //do not take too small energies
92236b27 237 if (fEnergy < kBaseLine)
379c5c09 238 fEnergy = 0;
239
240 //do not test quality of too soft samples
92236b27 241 if (fEnergy < maxEtoFit){
242 fTime = tStart;
243 if (aRMS < 2.) //sigle peak
244 fQuality = 999. ;
379c5c09 245 else
92236b27 246 fQuality = 0. ;
379c5c09 247 return kTRUE ;
248 }
249
92236b27 250 // if sample has reasonable mean and RMS, try to fit it with gamma2
379c5c09 251
252 gMinuit->mncler(); // Reset Minuit's list of paramters
253 gMinuit->SetPrintLevel(-1) ; // No Printout
254 gMinuit->SetFCN(AliPHOSRawFitterv1::UnfoldingChiSquare) ;
7dac04a3 255
379c5c09 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) ;
33373c87 266 fSampleParamsLow->AddAt(double(iStart),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);
33373c87 278 fSampleParamsHigh->AddAt(double(iStart),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) ) ;
7dac04a3 292 fEnergy = 0. ;
92236b27 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
de0cfd46 332 Double_t err=0.,t0err=0. ;
333 Double_t t0=0.,efit=0. ;
92236b27 334 gMinuit->GetParameter(0,t0, t0err) ;
335 gMinuit->GetParameter(1,efit, err) ;
379c5c09 336
de0cfd46 337 Double_t bfit=0., berr=0. ;
379c5c09 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
de0cfd46 355 Double_t fmin=0.,fedm=0.,errdef=0. ;
379c5c09 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 ;
33373c87 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}
7dac04a3 373//-----------------------------------------------------------------------------
374Double_t AliPHOSRawFitterv1::Gamma2(Double_t dt,Double_t en,Double_t b,TArrayD * params){ //Function for fitting samples
375 //parameters:
376 //dt-time after start
377 //en-amplutude
378 //function parameters
379
380 Double_t ped=params->At(4) ;
381 if(dt<0.)
382 return ped ; //pedestal
383 else
384 return ped+en*(TMath::Power(dt,params->At(0))*TMath::Exp(-dt*params->At(1))+b*dt*dt*TMath::Exp(-dt*params->At(3))) ;
385}
379c5c09 386//_____________________________________________________________________________
387void AliPHOSRawFitterv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
388{
389 // Number of parameters, Gradient, Chi squared, parameters, what to do
390
391 TList * toFit= (TList*)gMinuit->GetObjectFit() ;
392 TArrayD * params=(TArrayD*)toFit->At(0) ;
393 TArrayI * samples = (TArrayI*)toFit->At(1) ;
394 TArrayI * times = (TArrayI*)toFit->At(2) ;
395
396 fret = 0. ;
397 if(iflag == 2)
398 for(Int_t iparam = 0 ; iparam < 3 ; iparam++)
399 Grad[iparam] = 0 ; // Will evaluate gradient
400
401 Double_t t0=x[0] ;
402 Double_t en=x[1] ;
403 Double_t b=x[2] ;
404 Double_t n=params->At(0) ;
405 Double_t alpha=params->At(1) ;
406 Double_t beta=params->At(3) ;
7dac04a3 407 // Double_t ped=params->At(4) ;
379c5c09 408
409 Double_t overflow=params->At(5) ;
410 Int_t iBin = (Int_t) params->At(6) ;
411 Int_t nSamples=TMath::Min(iBin+70,samples->GetSize()) ; //Here we set number of points to fit (70)
412 // iBin - first non-zero sample
413 Int_t tStep=times->At(iBin+1)-times->At(iBin) ;
414 Double_t ddt=times->At(iBin)-t0-tStep ;
415 Double_t exp1=TMath::Exp(-alpha*ddt) ;
416 Double_t exp2=TMath::Exp(-beta*ddt) ;
417 Double_t dexp1=TMath::Exp(-alpha*tStep) ;
418 Double_t dexp2=TMath::Exp(-beta*tStep) ;
419 for(Int_t i = iBin; i<nSamples ; i++) {
420 Double_t dt=double(times->At(i))-t0 ;
421 Double_t fsample = double(samples->At(i)) ;
7dac04a3 422 Double_t diff=0. ;
379c5c09 423 exp1*=dexp1 ;
424 exp2*=dexp2 ;
de0cfd46 425// if(fsample>=overflow)
426// continue ;
379c5c09 427 if(dt<=0.){
7dac04a3 428 diff=fsample ;
379c5c09 429 fret += diff*diff ;
430 continue ;
431 }
7dac04a3 432
379c5c09 433 Double_t dtn=TMath::Power(dt,n) ;
434 Double_t dtnE=dtn*exp1 ;
435 Double_t dt2E=dt*dt*exp2 ;
7dac04a3 436 Double_t fit=en*(dtnE + b*dt2E) ;
de0cfd46 437 if(fsample>=overflow && fit>=overflow)
438 continue ;
439
379c5c09 440 diff = fsample - fit ;
441 fret += diff*diff ;
442 if(iflag == 2){ // calculate gradient
443 Grad[0] += en*diff*(dtnE*(n/dt-alpha)+b*dt2E*(2./dt-beta)) ; //derivative over t0
444 Grad[1] -= diff*(dtnE+b*dt2E) ;
445 Grad[2] -= en*diff*dt2E ;
446 }
7dac04a3 447
379c5c09 448 }
449 if(iflag == 2)
450 for(Int_t iparam = 0 ; iparam < 3 ; iparam++)
451 Grad[iparam] *= 2. ;
452}