GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TPeakFitter.cxx
Go to the documentation of this file.
1#include "TPeakFitter.h"
2#include "Math/Minimizer.h"
3#include "TH1.h"
4#include "TClass.h"
5
6#include "Globals.h"
7
8TPeakFitter::TPeakFitter(const Double_t& rangeLow, const Double_t& rangeHigh)
9 : fRangeLow(rangeLow), fRangeHigh(rangeHigh)
10{
11 fBGToFit = new TF1("fbg", this, &TPeakFitter::DefaultBackgroundFunction, fRangeLow, fRangeHigh, 4, "TPeakFitter", "DefaultBackgroundFunction");
12 fBGToFit->FixParameter(3, 0);
13 fBGToFit->SetLineColor(static_cast<Color_t>(kRed + fColorIndex));
14}
15
16void TPeakFitter::Print(Option_t* opt) const
17{
18 /// Print information from the fit, opt is passed along to each individual peaks Print function.
19 if(fTotalFitFunction != nullptr) {
20 double xMin = 0.;
21 double xMax = 0.;
22 fTotalFitFunction->GetRange(xMin, xMax);
23 std::cout << "Fitting range: " << xMin << " to " << xMax << std::endl;
24 } else {
25 std::cout << "Range: " << fRangeLow << " to " << fRangeHigh << std::endl;
26 }
27 Int_t counter = 0;
28 if(!fPeaksToFit.empty()) {
29 std::cout << "Peaks: " << std::endl;
30 for(auto* peak : fPeaksToFit) {
31 std::cout << "Peak #" << counter++ << " - ";
32 peak->Print(opt);
33 }
34 std::cout << "Chi2/Ndf = " << fPeaksToFit.front()->GetChi2() << "/" << fPeaksToFit.front()->GetNDF() << " = " << fPeaksToFit.front()->GetReducedChi2() << std::endl;
35 } else {
36 std::cout << "Could not find peak to print" << std::endl;
37 }
38 if(fBGToFit != nullptr) {
39 std::cout << "BG: " << std::endl;
40 fBGToFit->Print(opt);
41 } else {
42 std::cout << "Could not find a BG to print" << std::endl;
43 }
44}
45
47{
48 /// Print the range of the fit and the parameters of each peak on a single line.
49 std::cout << "Range: " << fRangeLow << " to " << fRangeHigh << " - ";
50 int i = 0;
51 for(auto* peak : fPeaksToFit) {
52 std::cout << peak->IsA()->GetName() << " #" << i++ << " ";
53 peak->PrintParameters();
54 }
55 std::cout << std::endl;
56}
57
59{
60 Int_t n_par = 0;
61 for(auto* peak : fPeaksToFit) {
62 n_par += peak->GetNParameters();
63 }
64 if(fBGToFit != nullptr) {
65 n_par += fBGToFit->GetNpar();
66 }
67
68 return n_par;
69}
70
71void TPeakFitter::SetRange(const Double_t& low, const Double_t& high)
72{
73 fRangeLow = low;
74 fRangeHigh = high;
75}
76
77TFitResultPtr TPeakFitter::Fit(TH1* fit_hist, Option_t* opt)
78{
79 /// Fit the histogram. Recognized options are "q" for a quiet fit, "retryfit" to retry a fit without parameter limits
80 /// if one of the parameters got close to its limit. All options (apart from "retryfit") are passed on to the actual
81 /// call to TH1::Fit.
82 if(fPeaksToFit.empty()) {
83 std::cout << "No peaks provided!" << std::endl;
84 return {};
85 }
86
87 TString options = opt;
88 options.ToLower();
89 bool quiet = options.Contains("q");
90 bool verbose = options.Contains("v");
91 if(quiet && verbose) {
92 std::cout << "Don't know how to be quiet and verbose at once (" << opt << "), going to be verbose!" << std::endl;
93 quiet = false;
94 }
95 bool retryFit = options.Contains("retryfit");
96 options.ReplaceAll("retryfit", "");
97
98 // store original range
99 int firstBin = fit_hist->GetXaxis()->GetFirst();
100 int lastBin = fit_hist->GetXaxis()->GetLast();
101 // set range to fit range of peak fitter
102 fit_hist->GetXaxis()->SetRangeUser(fRangeLow, fRangeHigh);
103
104 ROOT::Math::MinimizerOptions::SetDefaultMinimizer("Minuit2", "Combination");
105 TVirtualFitter::SetMaxIterations(100000);
106 TVirtualFitter::SetPrecision(1e-4);
107 if(fTotalFitFunction == nullptr) {
108 fTotalFitFunction = new TF1("total_fit", this, &TPeakFitter::FitFunction, fRangeLow, fRangeHigh, GetNParameters(), "TPeakFitter", "FitFunction");
109 }
110 fTotalFitFunction->SetLineColor(static_cast<Color_t>(kMagenta + fColorIndex));
112 fBGToFit->SetRange(fRangeLow, fRangeHigh);
113 // We need to initialize all of the parameters based on the peak parameters
114 if(fit_hist != nullptr && (!fInitFlag || fit_hist != fLastHistFit)) {
115 if(verbose) { std::cout << "Initializing Fit...." << std::endl; }
117 InitializeParameters(fit_hist);
118 fInitFlag = true;
119 fLastHistFit = fit_hist;
120 }
122 if(verbose) { PrintParameters(); }
123 // Do a first fit to get closer on the parameters
124 fit_hist->Fit(fTotalFitFunction, "RNQ0");
125 // Now try the fit with the options provided
126 fit_hist->Fit(fTotalFitFunction, Form("RNQO%s", options.Data()));
127 // Now do a final fit with integrated bins
128 TFitResultPtr fit_res = fit_hist->Fit(fTotalFitFunction, Form("SRI%s", options.Data()));
129
130 // check the parameter errors
131 if(!TGRSIFunctions::CheckParameterErrors(fit_res, options.Data())) {
132 if(retryFit) {
133 // fit again with all parameters released
134 if(!quiet) { std::cout << GREEN << "Re-fitting with released parameters (without any limits)" << RESET_COLOR << std::endl; }
135 for(int i = 0; i < fTotalFitFunction->GetNpar(); ++i) {
136 if(i == 1) { continue; } // skipping centroid, which should always be parameter 1
137 fTotalFitFunction->ReleaseParameter(i);
138 }
139 fit_res = fit_hist->Fit(fTotalFitFunction, Form("SRI%s", options.Data()));
140 TGRSIFunctions::CheckParameterErrors(fit_res, options.Data());
141 } else {
142 // re-try using minos instead of minuit (only thing from opt that might play a role it 'q' or 'v')
143 if(!quiet) { std::cout << YELLOW << "Re-fitting with \"E\" option to get better error estimation using Minos technique." << RESET_COLOR << std::endl; }
144 fit_res = fit_hist->Fit(fTotalFitFunction, Form("SRIE%s", options.Data()));
145 }
146 }
147
148 if(fit_res.Get() != nullptr) {
149 fTotalFitFunction->SetFitResult(*fit_res); // Not sure if this is needed
150
151 // Once we do the fit, we want to update all of the Peak parameters.
152 UpdatePeakParameters(fit_res, fit_hist);
153 } else if(!quiet) {
154 std::cout << RED << "Failed to get fit result, be aware that the results and especially the error estimates might be wrong!" << RESET_COLOR << std::endl;
155 }
156 fit_hist->Draw("hist");
157 fTotalFitFunction->Draw("same");
158 fPeaksToFit.front()->DrawBackground("same");
159
160 // always print the result of the fit even if not verbose
161 if(!quiet) {
162 std::cout << "****************" << std::endl
163 << "Summary of Fit: " << std::endl;
164 Print();
165 }
166 // restore old range
167 fit_hist->GetXaxis()->SetRange(firstBin, lastBin);
168 // reset the default back to minuit in case we switched to minos in between
169 ROOT::Math::MinimizerOptions::SetDefaultMinimizer("Minuit2", "Combination");
170
171 return fit_res;
172}
173
174void TPeakFitter::UpdatePeakParameters(const TFitResultPtr& fit_res, TH1* fit_hist)
175{
176 // We enter this function once we have performed a fit and need to tell the peaks
177 // what their new parameters should be.
178 Int_t param_counter = 0;
179
180 TF1* global_bg = new TF1("global_bg", this, &TPeakFitter::BackgroundFunction, fRangeLow, fRangeHigh, fTotalFitFunction->GetNpar(), "TPeakFitter", "BackgroundFunction");
181 global_bg->SetParameters(fTotalFitFunction->GetParameters());
182
183 // Start by looping through all of the peaks in the fitter.
184 for(auto* p_it : fPeaksToFit) {
185 TF1* peak_func = p_it->GetFitFunction();
186 // We need to make a clone of the main fit function, and set all of the non-peak parameters to 0.
187 // This will allow us to do the proper integrals. We have to do the same with the covariance matrix
188 // We get the covariance matrix that we are using to update the individual peaks
189 TF1* total_function_copy = new TF1(*fTotalFitFunction);
190 TMatrixDSym covariance_matrix = fit_res->GetCovarianceMatrix();
191 bool goodCovarianceMatrix = true;
192 if(covariance_matrix.GetNrows() < peak_func->GetNpar() || covariance_matrix.GetNcols() < peak_func->GetNpar()) {
193 goodCovarianceMatrix = false;
194 }
195 // Now we need to remove all of the parts of the covariance peak that has nothing to do with this particular peak
196 // This would be the diagonals of the background, and all of the other peaks.
197 Int_t param_to_zero_counter = 0;
198 std::vector<Int_t> param_to_zero_list;
199 for(auto* other_p_it : fPeaksToFit) {
200 if(other_p_it != p_it) {
201 for(int i = 0; i < other_p_it->GetFitFunction()->GetNpar(); ++i) {
202 param_to_zero_list.push_back(param_to_zero_counter);
203 ++param_to_zero_counter;
204 }
205 } else {
206 for(int i = 0; i < peak_func->GetNpar(); ++i) {
207 if(p_it->IsBackgroundParameter(i)) {
208 param_to_zero_list.push_back(param_to_zero_counter);
209 }
210 ++param_to_zero_counter;
211 }
212 Double_t low_range = 0.;
213 Double_t high_range = 0.;
214 fTotalFitFunction->GetRange(low_range, high_range);
215 peak_func->SetRange(low_range, high_range);
216 }
217 }
218 // We have now taken care to zero all of the non-peak parameters in the peak list, so we want to remove the total background contribution
219 for(int i = param_to_zero_counter; i < fTotalFitFunction->GetNpar(); ++i) {
220 param_to_zero_list.push_back(i);
221 }
222
223 // zero other non-peak portions of matrix
224 for(auto i : param_to_zero_list) {
225 for(auto j : param_to_zero_list) {
226 if(goodCovarianceMatrix) { covariance_matrix(i, j) = 0.0; }
227 }
228 total_function_copy->SetParameter(i, 0.0);
229 }
230 if(peak_func != nullptr) {
231 for(int i = 0; i < peak_func->GetNpar(); ++i) {
232 peak_func->SetParameter(i, fTotalFitFunction->GetParameter(param_counter));
233 peak_func->SetParError(i, fTotalFitFunction->GetParError(param_counter));
234 Double_t low = 0.;
235 Double_t high = 0.;
236 fTotalFitFunction->GetParLimits(param_counter, low, high);
237 peak_func->SetParLimits(i, low, high);
238 ++param_counter;
239 // Lets do some integrals meow.
240 }
241 p_it->SetArea(total_function_copy->Integral(p_it->Centroid() - p_it->Width() * 5., p_it->Centroid() + p_it->Width() * 5., 1e-8) / fit_hist->GetBinWidth(1));
242 if(goodCovarianceMatrix) {
243 p_it->SetAreaErr(total_function_copy->IntegralError(p_it->Centroid() - p_it->Width() * 5., p_it->Centroid() + p_it->Width() * 5., total_function_copy->GetParameters(), covariance_matrix.GetMatrixArray(), 1E-5) / fit_hist->GetBinWidth(1));
244 } else {
245 std::cout << "Not setting area error because we don't have a good covariance matrix!" << std::endl;
246 }
247 }
248 total_function_copy->Delete();
249 }
250 if(fBGToFit != nullptr) {
251 for(int i = 0; i < fBGToFit->GetNpar(); ++i) {
252 fBGToFit->SetParameter(i, fTotalFitFunction->GetParameter(param_counter));
253 fBGToFit->SetParError(i, fTotalFitFunction->GetParError(param_counter));
254 ++param_counter;
255 }
256 }
257 // We now have a copy of the global background. Copy and send to every peak
258 for(auto* p_it : fPeaksToFit) {
259 p_it->SetGlobalBackground(new TF1(*global_bg));
260 p_it->SetChi2(fTotalFitFunction->GetChisquare());
261 p_it->SetNDF(fTotalFitFunction->GetNDF());
262 }
263 global_bg->Delete();
264}
265
267{
268 for(auto* p_it : fPeaksToFit) {
269 p_it->InitializeParameters(fit_hist, fRangeLow, fRangeHigh);
270 }
271}
272
274{
275 /// This functions gets the parameters and their limits from the peak functions and
276 /// sets them for the total fit function
277 Int_t param_counter = 0;
278 Int_t peak_counter = 0;
279 for(auto* p_it : fPeaksToFit) {
280 TF1* peak_func = p_it->GetFitFunction();
281 if(peak_func != nullptr) {
282 for(int i = 0; i < peak_func->GetNpar(); ++i) {
283 fTotalFitFunction->SetParName(param_counter, Form("%s_%i", peak_func->GetParName(i), peak_counter));
284 fTotalFitFunction->SetParameter(param_counter, peak_func->GetParameter(i));
285 fTotalFitFunction->SetParError(param_counter, peak_func->GetParError(i));
286 Double_t limit_low = 0.;
287 Double_t limit_high = 0.;
288 peak_func->GetParLimits(i, limit_low, limit_high);
289 fTotalFitFunction->SetParLimits(param_counter, limit_low, limit_high);
290 ++param_counter;
291 }
292 ++peak_counter;
293 }
294 }
295 if(fBGToFit != nullptr) {
296 for(int i = 0; i < fBGToFit->GetNpar(); ++i) {
297 fTotalFitFunction->SetParName(param_counter, fBGToFit->GetParName(i));
298 fTotalFitFunction->SetParameter(param_counter, fBGToFit->GetParameter(i));
299 fTotalFitFunction->SetParError(param_counter, fBGToFit->GetParError(i));
300 Double_t limit_low = 0.;
301 Double_t limit_high = 0.;
302 fBGToFit->GetParLimits(i, limit_low, limit_high);
303 fTotalFitFunction->SetParLimits(param_counter, limit_low, limit_high);
304 ++param_counter;
305 }
306 }
307}
308
309Double_t TPeakFitter::FitFunction(Double_t* dim, Double_t* par)
310{
311 // I want to use the EvalPar command here in order to ge the individual peaks
312 Double_t sum = 0;
313 Int_t params_so_far = 0;
314 for(auto* p_it : fPeaksToFit) {
315 TF1* peak_func = p_it->GetFitFunction();
316 sum += peak_func->EvalPar(dim, &par[params_so_far]);
317 params_so_far += peak_func->GetNpar();
318 }
319 sum += fBGToFit->EvalPar(dim, &par[params_so_far]);
320
321 return sum;
322}
323
324Double_t TPeakFitter::BackgroundFunction(Double_t* dim, Double_t* par)
325{
326 Double_t sum = 0;
327 Int_t params_so_far = 0;
328 for(auto* p_it : fPeaksToFit) {
329 TF1* peak_func = p_it->GetBackgroundFunction();
330 sum += peak_func->EvalPar(dim, &par[params_so_far]);
331 params_so_far += peak_func->GetNpar();
332 }
333 sum += fBGToFit->EvalPar(dim, &par[params_so_far]);
334
335 return sum;
336}
337
339{
340 fBGToFit->SetParName(0, "A");
341 fBGToFit->SetParName(1, "B");
342 fBGToFit->SetParName(2, "C");
343 fBGToFit->SetParName(3, "bg_offset");
344
345 Double_t lowLimit = 0.;
346 Double_t highLimit = 0.;
347 fBGToFit->GetParLimits(0, lowLimit, highLimit);
348
349 double value = fBGToFit->GetParameter(0);
350 if(value == 0 && lowLimit == 0 && highLimit == 0) {
351 fBGToFit->SetParLimits(0, 0.0, fit_hist->GetBinContent(fit_hist->FindBin(fRangeHigh)) * 100.);
352 fBGToFit->SetParameter("A", fit_hist->GetBinContent(fit_hist->FindBin(fRangeHigh)));
353 }
354 fBGToFit->GetParLimits(1, lowLimit, highLimit);
355
356 value = fBGToFit->GetParameter(1);
357 if(value == 0 && lowLimit == 0 && highLimit == 0) {
358 fBGToFit->SetParameter("B", fit_hist->GetBinWidth(1) * ((fit_hist->GetBinContent(fit_hist->FindBin(fRangeLow)) - fit_hist->GetBinContent(fit_hist->FindBin(fRangeHigh))) / (fRangeLow - fRangeHigh)));
359 }
360 fBGToFit->GetParLimits(2, lowLimit, highLimit);
361
362 value = fBGToFit->GetParameter(2);
363 if(value == 0 && lowLimit == 0 && highLimit == 0) {
364 fBGToFit->SetParameter("C", 0.0000);
365 fBGToFit->FixParameter(2, 0.00);
366 }
367 fBGToFit->GetParLimits(3, lowLimit, highLimit);
368
369 value = fBGToFit->GetParameter(3);
370 if(value == 0 && lowLimit == 0 && highLimit == 0) {
371 fBGToFit->SetParLimits(3, fRangeLow, fRangeHigh);
372 if(!fPeaksToFit.empty()) {
373 fBGToFit->SetParameter("bg_offset", fPeaksToFit.front()->Centroid());
374 } else {
375 fBGToFit->SetParameter("bg_offset", (fRangeHigh + fRangeLow) / 2.);
376 }
377 }
378}
379
380Double_t TPeakFitter::DefaultBackgroundFunction(Double_t* dim, Double_t* par) // NOLINT(readability-non-const-parameter)
381{
382 Double_t x = dim[0];
383 Double_t A = par[0];
384 Double_t B = par[1];
385 Double_t C = par[2];
386 Double_t bg_offset = par[3];
387
388 return A + B * (x - bg_offset) + C * TMath::Power(x - bg_offset, 2.);
389}
390
391void TPeakFitter::DrawPeaks(Option_t*) const
392{
393 // First we are going to draw the background
394 TF1* bg_to_draw = new TF1;
395 fTotalFitFunction->Copy(*bg_to_draw);
396 bg_to_draw->SetLineColor(static_cast<Color_t>(kRed + fColorIndex));
397
398 for(auto* p_it : fPeaksToFit) {
399 TF1* peak_func = p_it->GetFitFunction();
400 TF1* total_function_copy = new TF1;
401 fTotalFitFunction->Copy(*total_function_copy);
402 // Now we need to remove all of the parts of the covariance peak that has nothing to do with this particular peak
403 // This would be the diagonals of the background, and all of the other peaks.
404 Int_t param_to_zero_counter = 0;
405 for(auto* other_p_it : fPeaksToFit) {
406 if(other_p_it != p_it) {
407 for(int i = 0; i < other_p_it->GetFitFunction()->GetNpar(); ++i) {
408 total_function_copy->SetParameter(param_to_zero_counter, 0.0);
409 ++param_to_zero_counter;
410 }
411 } else {
412 for(int i = 0; i < peak_func->GetNpar(); ++i) {
413 if(!p_it->IsBackgroundParameter(i)) {
414 bg_to_draw->SetParameter(param_to_zero_counter, 0.0);
415 }
416 ++param_to_zero_counter;
417 }
418 }
419 }
420 total_function_copy->Draw("same");
421 }
422 bg_to_draw->Draw("same");
423}
#define RED
Definition Globals.h:9
#define YELLOW
Definition Globals.h:7
#define GREEN
Definition Globals.h:8
#define RESET_COLOR
Definition Globals.h:5
void UpdatePeakParameters(const TFitResultPtr &fit_res, TH1 *fit_hist)
TF1 * fTotalFitFunction
Definition TPeakFitter.h:91
void PrintParameters() const
void SetRange(const Double_t &low, const Double_t &high)
void Print(Option_t *opt="") const override
Int_t GetNParameters() const
void InitializeBackgroundParameters(TH1 *fit_hist)
Double_t FitFunction(Double_t *dim, Double_t *par)
Double_t fRangeLow
Definition TPeakFitter.h:93
void InitializeParameters(TH1 *fit_hist)
int fColorIndex
this index is added to the colors kRed for the total function and kMagenta for the individual peaks
MultiplePeak_t fPeaksToFit
Definition TPeakFitter.h:88
void DrawPeaks(Option_t *="") const
Double_t BackgroundFunction(Double_t *dim, Double_t *par)
Double_t fRangeHigh
Definition TPeakFitter.h:94
TFitResultPtr Fit(TH1 *fit_hist, Option_t *opt="")
Double_t DefaultBackgroundFunction(Double_t *dim, Double_t *par)
TH1 * fLastHistFit
TF1 * fBGToFit
Definition TPeakFitter.h:89
void UpdateFitterParameters()
bool CheckParameterErrors(const TFitResultPtr &fitres, std::string opt="")