GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TPeak.cxx
Go to the documentation of this file.
1#include "TPeak.h"
2#include "TGraph.h"
3
4#include "Math/Minimizer.h"
5
6Bool_t TPeak::fLogLikelihoodFlag = true;
7TPeak* TPeak::fLastFit = nullptr;
8
9TPeak::TPeak(Double_t cent, Double_t xlow, Double_t xhigh, TF1* background)
10 : TGRSIFit("photopeakbg", TGRSIFunctions::PhotoPeakBG, xlow, xhigh, 10), fResiduals(new TGraph)
11{
12 Clear();
13 Bool_t outOfRangeFlag = false;
14
15 if(cent > xhigh) {
16 std::cout << "centroid is higher than range" << std::endl;
17 outOfRangeFlag = true;
18 } else if(cent < xlow) {
19 std::cout << "centroid is lower than range" << std::endl;
20 outOfRangeFlag = true;
21 }
22
23 // This fixes things if your user is like me and screws up a lot.
24 if(outOfRangeFlag) {
25 if(xlow > cent) {
26 std::swap(xlow, cent);
27 }
28 if(xlow > xhigh) {
29 std::swap(xlow, xhigh);
30 }
31 if(cent > xhigh) {
32 std::swap(cent, xhigh);
33 }
34 std::cout << "Something about your range was wrong. Assuming:" << std::endl;
35 std::cout << "centroid: " << cent << " \t range: " << xlow << " to " << xhigh << std::endl;
36 }
37
38 SetRange(xlow, xhigh);
39 // We also need to make initial guesses at parameters
40 // We need nice ways to access parameters etc.
41 // Need to make a TMultipeak-like thing (does a helper class come into play then?)
42
43 // Set the fit function to be a radware style photo peak.
44 // This function might be unnecessary. Will revist this later. rd.
45 SetName(Form("Chan%d_%d_to_%d", static_cast<Int_t>(cent), static_cast<Int_t>(xlow),
46 static_cast<Int_t>(xhigh))); // Gives a default name to the peak
47
48 // We need to set parameter names now.
49 InitNames();
50 SetParameter("centroid", cent);
51
52 // Check to see if background is good.
53 if(background != nullptr) {
54 fBackground = background;
55 fOwnBgFlag = false;
56 } else {
57 fBackground = new TF1(
58 Form("background%d_%d_to_%d", static_cast<Int_t>(cent), static_cast<Int_t>(xlow), static_cast<Int_t>(xhigh)),
59 TGRSIFunctions::StepBG, xlow, xhigh, 10);
61 fOwnBgFlag = true;
62 }
63
64 fBackground->SetNpx(1000);
65 fBackground->SetLineStyle(2);
66 fBackground->SetLineColor(kBlack);
67}
68
69TPeak::TPeak() : TGRSIFit("photopeakbg", TGRSIFunctions::PhotoPeakBG, 0, 1000, 10),
70 fBackground(new TF1("background", TGRSIFunctions::StepBG, 0, 1000, 10)), fResiduals(new TGraph)
71{
72 InitNames();
73 fBackground->SetNpx(1000);
74 fBackground->SetLineStyle(2);
75 fBackground->SetLineColor(kBlack);
77}
78
80{
81 if(fOwnBgFlag) {
82 delete fBackground;
83 }
84 delete fResiduals;
85}
86
87TPeak::TPeak(const TPeak& copy) : TGRSIFit(copy)
88{
89 copy.Copy(*this);
90}
91
93{
94 SetParName(0, "Height");
95 SetParName(1, "centroid");
96 SetParName(2, "sigma");
97 SetParName(3, "beta");
98 SetParName(4, "R");
99 SetParName(5, "step");
100 SetParName(6, "A");
101 SetParName(7, "B");
102 SetParName(8, "C");
103 SetParName(9, "bg_offset");
104}
105
106void TPeak::Copy(TObject& obj) const
107{
108 TGRSIFit::Copy(obj);
109
110 if(static_cast<TPeak&>(obj).fBackground == nullptr) {
111 static_cast<TPeak&>(obj).fBackground = new TF1(*fBackground);
112 }
113 if(static_cast<TPeak&>(obj).fResiduals == nullptr) {
114 static_cast<TPeak&>(obj).fResiduals = new TGraph(*fResiduals);
115 }
116
117 static_cast<TPeak&>(obj).fArea = fArea;
118 static_cast<TPeak&>(obj).fDArea = fDArea;
119
120 static_cast<TPeak&>(obj).fChi2 = fChi2;
121 static_cast<TPeak&>(obj).fNdf = fNdf;
122
123 *(static_cast<TPeak&>(obj).fBackground) = *fBackground;
124 // We are making a direct copy of the background so the ownership is that of the copy
125 static_cast<TPeak&>(obj).fOwnBgFlag = true;
126 *(static_cast<TPeak&>(obj).fResiduals) = *fResiduals;
127
128 static_cast<TPeak&>(obj).SetHist(GetHist());
129}
130
131Bool_t TPeak::InitParams(TH1* fitHist)
132{
133 // Makes initial guesses at parameters for the fit. Uses the histogram to
134 if(fitHist == nullptr) {
135 return false;
136 }
137 Double_t xlow = 0.;
138 Double_t xhigh = 0.;
139 GetRange(xlow, xhigh);
140 Int_t bin = fitHist->FindBin(GetParameter("centroid"));
141 Int_t binlow = fitHist->GetXaxis()->FindBin(xlow);
142 Int_t binhigh = fitHist->GetXaxis()->FindBin(xhigh);
143 SetParLimits(0, 0, fitHist->GetMaximum());
144 Double_t low = 0.;
145 Double_t high = 0.;
146 GetParLimits(1, low, high);
147 if(low == high && low == 0.) {
148 SetParLimits(1, xlow, xhigh);
149 }
150 GetParLimits(2, low, high);
151 if(low == high && low == 0.) {
152 SetParLimits(2, 0.5, (xhigh - xlow)); // sigma should be less than the window width - JKS
153 }
154 SetParLimits(3, 0.000001, 10);
155 SetParLimits(4, 0.000001, 100); // this is a percentage. no reason for it to go to 500% - JKS
156 // Step size is allow to vary to anything. If it goes below 0, the code will fix it to 0
157 SetParLimits(5, 0.0, 1.0E2);
158 SetParLimits(6, 0.0, fitHist->GetBinContent(bin) * 100.);
159 SetParLimits(9, xlow, xhigh);
160
161 if((fitHist == nullptr) && (GetHist() != nullptr)) {
162 fitHist = GetHist();
163 }
164
165 if(fitHist == nullptr) {
166 std::cout << "No histogram is associated yet, no initial guesses made" << std::endl;
167 return false;
168 }
169 // Make initial guesses
170 // Actually set the parameters in the photopeak function
171 // Fixing has to come after setting
172 // Might have to include bin widths eventually
173 // The centroid should already be set by this point in the ctor
174 SetParameter("Height", fitHist->GetBinContent(bin));
175 SetParameter("centroid", GetParameter("centroid"));
176 SetParameter("sigma", TMath::Sqrt(9.0 + 4. * GetParameter("centroid") / 1000.) / 2.35);
177 SetParameter("beta", GetParameter("sigma") / 2.0);
178 SetParameter("R", 0.001);
179 SetParameter("step", 1.0);
180 SetParameter("A", fitHist->GetBinContent(binhigh));
181 SetParameter("B", (fitHist->GetBinContent(binlow) - fitHist->GetBinContent(binhigh)) / (xlow - xhigh));
182 SetParameter("C", 0.0000);
183 SetParameter("bg_offset", GetParameter("centroid"));
184 FixParameter(8, 0.00);
185 FixParameter(3, GetParameter("beta"));
186 FixParameter(4, 0.00);
188 return true;
189}
190
191Bool_t TPeak::Fit(TH1* fitHist, Option_t* opt)
192{
193 TString options = opt;
194 options.ToLower();
195 bool quiet = options.Contains("q");
196 bool verbose = options.Contains("v");
197 if(quiet && verbose) {
198 std::cout << "Don't know how to be quiet and verbose at once (" << opt << "), going to be verbose!" << std::endl;
199 quiet = false;
200 }
201 bool retryFit = options.Contains("retryfit");
202 options.ReplaceAll("retryfit", "");
203 if(!verbose && !quiet) { options.Append("q"); }
204
205 if(fitHist == nullptr && GetHist() == nullptr) {
206 std::cout << "No hist passed, trying something... ";
207 fitHist = fHistogram;
208 }
209 if(fitHist == nullptr) {
210 std::cout << "No histogram associated with Peak" << std::endl;
211 return false;
212 }
213 if(!IsInitialized()) {
214 InitParams(fitHist);
215 }
216 ROOT::Math::MinimizerOptions::SetDefaultMinimizer("Minuit2", "Combination");
217 TVirtualFitter::SetMaxIterations(100000);
218 TVirtualFitter::SetPrecision(1e-3);
219
220 SetHist(fitHist);
221
222 // Now that it is initialized, let's fit it.
223 // Just in case the range changed, we should reset the centroid and bg energy limits
224 // check first if the parameter has been fixed!
225 std::vector<double> lowerLimit(GetNpar());
226 std::vector<double> upperLimit(GetNpar());
227 for(int i = 0; i < GetNpar(); ++i) {
228 GetParLimits(i, lowerLimit[i], upperLimit[i]);
229 if(i < 2) { // height, and centroid
230 FixParameter(i, GetParameter(i));
231 }
232 }
233
234 TFitResultPtr fitres;
235 // Log likelihood is the proper fitting technique UNLESS the data is a result of an addition or subtraction.
237 fitres = fitHist->Fit(this, Form("%sRLSN", options.Data())); // The RS needs to always be there
238 } else {
239 fitres = fitHist->Fit(this, Form("%sRSN", options.Data())); // The RS needs to always be there
240 }
241
242 // Check fit exited successfully before continuing
243 if(static_cast<int>(fitres) == -1) { return false; }
244
245 for(int i = 0; i < GetNpar(); ++i) {
246 SetParLimits(i, lowerLimit[i], upperLimit[i]);
247 }
248
249 // Log likelihood is the proper fitting technique UNLESS the data is a result of an addition or subtraction.
251 fitres = fitHist->Fit(this, Form("%sRLS", options.Data())); // The RS needs to always be there
252 } else {
253 fitres = fitHist->Fit(this, Form("%sRS", options.Data())); // The RS needs to always be there
254 }
255
256 // Check fit exited successfully before continuing
257 if(static_cast<int>(fitres) == -1) { return false; }
258
259 // After performing this fit I want to put something here that takes the fit result (good,bad,etc)
260 // for printing out. RD
261
262 if(fitres->ParError(2) != fitres->ParError(2)) { // Check to see if nan
263 if(fitres->Parameter(3) < 1) {
264 InitParams(fitHist);
265 FixParameter(4, 0);
266 FixParameter(3, 1);
267 if(verbose) { std::cout << "Beta may have broken the fit, retrying with R=0" << std::endl; }
268 // Leaving the log-likelihood argument out so users are not constrained to just using that. - JKS
269 fitHist->GetListOfFunctions()->Last()->Delete();
271 fitres = fitHist->Fit(this, Form("%sRLS", options.Data())); // The RS needs to always be there
272 } else {
273 fitres = fitHist->Fit(this, Form("%sRS", options.Data()));
274 }
275 }
276 }
277
278 // check parameter errors
279 if(!TGRSIFunctions::CheckParameterErrors(fitres, options.Data())) {
280 if(retryFit) {
281 // fit again with all parameters released
282 if(!quiet) { std::cout << GREEN << "Re-fitting with released parameters (without any limits)" << RESET_COLOR << std::endl; }
283 for(int i = 0; i < GetNpar(); ++i) {
284 ReleaseParameter(i);
285 }
287 fitres = fitHist->Fit(this, Form("%sRLS", options.Data())); // The RS needs to always be there
288 } else {
289 fitres = fitHist->Fit(this, Form("%sRS", options.Data()));
290 }
291 } else {
292 // re-try using minos instead of minuit
293 if(!quiet) { std::cout << YELLOW << "Re-fitting with \"E\" option to get better error estimation using Minos technique." << RESET_COLOR << std::endl; }
295 fitres = fitHist->Fit(this, Form("%sERLS", options.Data())); // The RS needs to always be there
296 } else {
297 fitres = fitHist->Fit(this, Form("%sERS", options.Data()));
298 }
299 }
300 }
302
303 Double_t binWidth = fitHist->GetBinWidth(1);
304 Double_t width = GetParameter("sigma");
305 if(verbose) {
306 std::cout << "Chi^2/NDF = " << fitres->Chi2() / fitres->Ndf() << std::endl;
307 }
308 fChi2 = fitres->Chi2();
309 fNdf = fitres->Ndf();
310 Double_t xlow = 0.;
311 Double_t xhigh = 0.;
312 GetRange(xlow, xhigh);
313 Double_t int_low = xlow - 10. * width; // making the integration bounds a bit smaller, but still large enough. -JKS
314 Double_t int_high = xhigh + 10. * width;
315
316 // Make a function that does not include the background
317 // Integrate the background.
318 // TPeak* tmppeak = new TPeak(*this);
319 auto* tmppeak = new TPeak;
320 Copy(*tmppeak);
321 tmppeak->SetParameter("step", 0.0);
322 tmppeak->SetParameter("A", 0.0);
323 tmppeak->SetParameter("B", 0.0);
324 tmppeak->SetParameter("C", 0.0);
325 tmppeak->SetParameter("bg_offset", 0.0);
326 tmppeak->SetRange(int_low, int_high); // This will help get the true area of the gaussian 200 ~ infinity in a gaus
327 tmppeak->SetName("tmppeak");
328
329 // SOMETHING IS WRONG WITH THESE UNCERTAINTIES
330 // This is where we will do integrals and stuff.
331 fArea = (tmppeak->Integral(int_low, int_high)) / binWidth;
332 // Set the background values in the covariance matrix to 0, while keeping their covariance errors
333 TMatrixDSym CovMat = fitres->GetCovarianceMatrix();
334 CovMat(5, 5) = 0.0;
335 CovMat(6, 6) = 0.0;
336 CovMat(7, 7) = 0.0;
337 CovMat(8, 8) = 0.0;
338 CovMat(9, 9) = 0.0;
339 fDArea = (tmppeak->IntegralError(int_low, int_high, tmppeak->GetParameters(), CovMat.GetMatrixArray())) / binWidth;
340
341 if(verbose) {
342 std::cout << "Integral: " << fArea << " +/- " << fDArea << std::endl;
343 }
344 // Set the background for drawing later
345 fBackground->SetParameters(GetParameters());
346 // To DO: put a flag in signalling that the errors are not to be trusted if we have a bad cov matrix
347 Copy(*fitHist->GetListOfFunctions()->Last());
348
349 // always print result of the fit even if not verbose
350 if(!quiet) { Print("+"); }
351 delete tmppeak;
352 fLastFit = this;
353 return true;
354}
355
356void TPeak::Clear(Option_t*)
357{
358 // Clear the TPeak including functions and histogram, does not
359 // currently clear inherited members such as name.
360 // want to make a clear that doesn't clear everything
361 fArea = 0.0;
362 fDArea = 0.0;
363 fChi2 = 0.0;
364 fNdf = 0.0;
366 // Do deep clean stuff maybe? require an option?
367}
368
369void TPeak::Print(Option_t* opt) const
370{
371 // Prints TPeak properties. To see More properties use the option "+"
372 std::cout << "Name: " << GetName() << std::endl;
373 std::cout << "Centroid: " << GetParameter("centroid") << " +/- " << GetParError(GetParNumber("centroid")) << std::endl;
374 std::cout << "Area: " << fArea << " +/- " << fDArea << std::endl;
375 std::cout << "FWHM: " << GetParameter("sigma") * 2.3548 << " +/- " << GetParError(GetParNumber("sigma")) * 2.3548 << std::endl;
376 std::cout << "Chi^2/NDF: " << fChi2 / fNdf << std::endl;
377 if(strchr(opt, '+') != nullptr) {
378 TF1::Print();
379 TGRSIFit::Print(opt); // Polymorphise this a bit better
380 }
381}
382
383void TPeak::DrawBackground(Option_t* opt) const
384{
385 fBackground->Draw(opt);
386}
387
389{
390 if(GetHist() == nullptr) {
391 std::cout << "No hist set" << std::endl;
392 return;
393 }
394 if(fChi2 < 0.000000001) {
395 std::cout << "No fit performed" << std::endl;
396 return;
397 }
398
399 Double_t xlow = 0.;
400 Double_t xhigh = 0.;
401 GetRange(xlow, xhigh);
402 Int_t nbins = GetHist()->GetXaxis()->GetNbins();
403 auto* res = new Double_t[nbins];
404 auto* bin = new Double_t[nbins];
405 Int_t points = 0;
406 fResiduals->Clear();
407
408 for(int i = 1; i <= nbins; i++) {
409 if(GetHist()->GetBinCenter(i) <= xlow || GetHist()->GetBinCenter(i) >= xhigh) {
410 continue;
411 }
412 res[points] = (GetHist()->GetBinContent(i) - Eval(GetHist()->GetBinCenter(i))) +
413 GetParameter("Height") / 2; /// GetHist()->GetBinError(i));// + GetParameter("Height") + 10.;
414 bin[points] = GetHist()->GetBinCenter(i);
415 fResiduals->SetPoint(i, bin[i], res[i]);
416
417 points++;
418 }
419
420 fResiduals->Draw();
421
422 delete[] res;
423 delete[] bin;
424}
425
427{
428 if(GetHist() == nullptr) {
429 std::cout << "No hist set" << std::endl;
430 return false;
431 }
432 if(fChi2 < 0.000000001) {
433 std::cout << "No fit performed" << std::endl;
434 return false;
435 }
436 return true;
437}
438
440{
441 if(!GoodStatus()) { return 0.; }
442
443 Double_t width = GetParameter("sigma");
444 Double_t xlow = 0.;
445 Double_t xhigh = 0.;
446 GetRange(xlow, xhigh);
447 Double_t int_low = xlow - 10. * width; // making the integration bounds a bit smaller, but still large enough. -JKS
448 Double_t int_high = xhigh + 10. * width;
449 return GetIntegralArea(int_low, int_high);
450}
451
452Double_t TPeak::GetIntegralArea(Double_t int_low, Double_t int_high)
453{
454 if(!GoodStatus()) { return 0.; }
455
456 // pull appropriate properties from peak and histogram
457 TH1* hist = GetHist();
458
459 // use those properties to integrate the histogram
460 Int_t binlow = hist->FindBin(int_low);
461 Int_t binhigh = hist->FindBin(int_high);
462 Double_t binWidth = hist->GetBinWidth(binlow);
463 Double_t hist_integral = hist->Integral(binlow, binhigh);
464 Double_t xlow = hist->GetXaxis()->GetBinLowEdge(binlow);
465 Double_t xhigh = hist->GetXaxis()->GetBinUpEdge(binhigh);
466
467 // ... and then integrate the background
468 Double_t bg_area = (Background()->Integral(xlow, xhigh)) / binWidth;
469
470 // calculate the peak area and error
471 Double_t peakarea = hist_integral - bg_area;
472
473 return peakarea;
474}
475
476Double_t TPeak::GetIntegralAreaErr(Double_t int_low, Double_t int_high)
477{
478 if(!GoodStatus()) { return 0.; }
479
480 // pull appropriate properties from peak and histogram
481 TH1* hist = GetHist();
482
483 // use those properties to integrate the histogram
484 Int_t binlow = hist->FindBin(int_low);
485 Int_t binhigh = hist->FindBin(int_high);
486 Double_t binWidth = hist->GetBinWidth(binlow);
487 Double_t hist_integral = hist->Integral(binlow, binhigh);
488 Double_t xlow = hist->GetXaxis()->GetBinLowEdge(binlow);
489 Double_t xhigh = hist->GetXaxis()->GetBinUpEdge(binhigh);
490
491 // ... and then integrate the background
492 Double_t bg_area = (Background()->Integral(xlow, xhigh)) / binWidth;
493
494 // calculate the peak error
495 Double_t peakerr = sqrt(hist_integral + bg_area);
496
497 return peakerr;
498}
499
501{
502 if(!GoodStatus()) { return 0.; }
503
504 Double_t width = GetParameter("sigma");
505 Double_t xlow = 0.;
506 Double_t xhigh = 0.;
507 GetRange(xlow, xhigh);
508 Double_t int_low = xlow - 10. * width; // making the integration bounds a bit smaller, but still large enough. -JKS
509 Double_t int_high = xhigh + 10. * width;
510 return GetIntegralAreaErr(int_low, int_high);
511}
512
513void TPeak::CheckArea(Double_t int_low, Double_t int_high)
514{
515 if(!GoodStatus()) { return; }
516
517 // calculate the peak area and error
518 Double_t peakarea = GetIntegralArea(int_low, int_high);
519 Double_t peakerr = GetIntegralAreaErr(int_low, int_high);
520
521 // now print properties
522 std::cout << "TPeak integral: " << fArea << " +/- " << fDArea << std::endl;
523 std::cout << "Histogram - BG integral: " << peakarea << " +/- " << peakerr << std::endl;
524 if(std::abs(peakarea - fArea) < (fDArea + peakerr)) {
525 std::cout << DGREEN << "Areas are consistent." << RESET_COLOR << std::endl;
526 } else if(std::abs(peakarea - fArea) < 2 * (fDArea + peakerr)) {
527 std::cout << DYELLOW << "Areas are consistent within 2 sigma." << RESET_COLOR << std::endl;
528 } else {
529 std::cout << DRED << "Areas are inconsistent." << RESET_COLOR << std::endl;
530 }
531}
532
534{
535 if(!GoodStatus()) { return; }
536
537 // calculate the peak area and error
538 Double_t peakarea = GetIntegralArea();
539 Double_t peakerr = GetIntegralAreaErr();
540
541 // now print properties
542 std::cout << "TPeak integral: " << fArea << " +/- " << fDArea << std::endl;
543 std::cout << "Histogram - BG integral: " << peakarea << " +/- " << peakerr << std::endl;
544 if(std::abs(peakarea - fArea) < (fDArea + peakerr)) {
545 std::cout << DGREEN << "Areas are consistent." << RESET_COLOR << std::endl;
546 } else if(std::abs(peakarea - fArea) < 2 * (fDArea + peakerr)) {
547 std::cout << DYELLOW << "Areas are consistent within 2 sigma." << RESET_COLOR << std::endl;
548 } else {
549 std::cout << DRED << "Areas are inconsistent." << RESET_COLOR << std::endl;
550 }
551}
#define DYELLOW
Definition Globals.h:16
#define DRED
Definition Globals.h:18
#define DGREEN
Definition Globals.h:17
#define YELLOW
Definition Globals.h:7
#define GREEN
Definition Globals.h:8
#define RESET_COLOR
Definition Globals.h:5
TH1D * hist
Definition UserFillObj.h:3
Bool_t AddToGlobalList(Bool_t yes=kTRUE) override
Definition TGRSIFit.cxx:62
virtual void SetHist(TH1 *hist)
Definition TGRSIFit.h:53
void Clear(Option_t *opt="") override
Definition TGRSIFit.cxx:36
Bool_t IsInitialized() const
Definition TGRSIFit.h:67
void SetInitialized(Bool_t flag=true)
Definition TGRSIFit.h:68
void Copy(TObject &obj) const override
Definition TGRSIFit.cxx:20
virtual TH1 * GetHist() const
Definition TGRSIFit.h:57
void Print(Option_t *opt="") const override
Definition TGRSIFit.cxx:27
Definition TPeak.h:28
Bool_t Fit(TH1 *fitHist, Option_t *opt="")
Definition TPeak.cxx:191
Double_t GetIntegralArea()
Definition TPeak.cxx:439
~TPeak()
Definition TPeak.cxx:79
TF1 * Background() const
Definition TPeak.h:96
void InitNames()
Definition TPeak.cxx:92
Double_t GetIntegralAreaErr()
Definition TPeak.cxx:500
Double_t fDArea
Definition TPeak.h:119
void DrawResiduals()
Definition TPeak.cxx:388
void CheckArea()
Definition TPeak.cxx:533
Double_t fChi2
Definition TPeak.h:120
Bool_t fOwnBgFlag
Definition TPeak.h:122
TGraph * fResiduals
Definition TPeak.h:128
void Print(Option_t *opt="") const override
Definition TPeak.cxx:369
TF1 * fBackground
Definition TPeak.h:127
void DrawBackground(Option_t *opt="SAME") const
Definition TPeak.cxx:383
bool GoodStatus()
Definition TPeak.cxx:426
void Copy(TObject &obj) const override
Definition TPeak.cxx:106
Double_t fNdf
Definition TPeak.h:121
Double_t fArea
Definition TPeak.h:118
Bool_t InitParams(TH1 *fitHist=nullptr) override
Definition TPeak.cxx:131
static TPeak * fLastFit
!
Definition TPeak.h:125
void Clear(Option_t *opt="") override
Definition TPeak.cxx:356
static bool fLogLikelihoodFlag
!
Definition TPeak.h:124
static Bool_t GetLogLikelihoodFlag()
Definition TPeak.h:103
TPeak()
Definition TPeak.cxx:69
bool CheckParameterErrors(const TFitResultPtr &fitres, std::string opt="")
Double_t StepBG(Double_t *dim, Double_t *par)