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