GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TMultiPeak.cxx
Go to the documentation of this file.
1#include "TMultiPeak.h"
2
3#include "TGRSIFunctions.h"
4
5#include "Math/MinimizerOptions.h"
6#include "TVirtualFitter.h"
7#include "TFitResult.h"
8
9#include <algorithm>
10
12
13TMultiPeak::TMultiPeak(Double_t xlow, Double_t xhigh, const std::vector<Double_t>& centroids, Option_t*)
14 : TGRSIFit("multipeakbg", this, &TMultiPeak::MultiPhotoPeakBG, xlow, xhigh, centroids.size() * 6 + 5, "TMultiPeak", "MultiPhotoPeakBG")
15{
16 std::cout << "Warning, the class TMultiPeak is deprecated (use TPeakFitter instead)!" << std::endl;
17 Clear();
18 // We make the background first so we can send it to the TPeaks.
19 fBackground = new TF1(Form("MPbackground_%d_to_%d", static_cast<Int_t>(xlow), static_cast<Int_t>(xhigh)), this, &TMultiPeak::MultiStepBG, xlow, xhigh, centroids.size() * 6 + 5, "TMuliPeak", "MultiStepBG");
20 fBackground->SetNpx(1000);
21 fBackground->SetLineStyle(2);
22 fBackground->SetLineColor(kBlack);
24
25 for(double cent : centroids) {
26 Bool_t out_of_range_flag = false;
27 if(cent > xhigh) {
28 std::cout << "centroid " << cent << " is higher than range" << std::endl;
29 out_of_range_flag = true;
30 } else if(cent < xlow) {
31 std::cout << "centroid " << cent << " is lower than range" << std::endl;
32 out_of_range_flag = true;
33 }
34 if(out_of_range_flag) {
35 std::cout << "ignoring peak at " << cent << ", make a new multi peak with the corrected energy" << std::endl;
36 } else {
37 auto* peak = new TPeak(cent, xlow, xhigh, fBackground);
38 peak->AddToGlobalList(kFALSE);
39 fPeakVec.push_back(peak);
40 }
41 }
42 SetRange(xlow, xhigh);
43
44 SetName(Form("MultiPeak_%d_to_%d", static_cast<Int_t>(xlow),
45 static_cast<Int_t>(xhigh))); // Gives a default name to the peak
46 SortPeaks(); // Defaults to sorting by TPeak::CompareEnergy
47 InitNames();
48}
49
50TMultiPeak::TMultiPeak() : TGRSIFit("multipeakbg", this, &TMultiPeak::MultiPhotoPeakBG, 0, 1000, 10, "TMultiPeak", "MultiPhotoPeakBG")
51{
52 std::cout << "Warning, the class TMultiPeak is deprecated (use TPeakFitter instead)!" << std::endl;
53 // I don't think this constructor should be used, RD.
54 InitNames();
55 fBackground = new TF1("background", this, &TMultiPeak::MultiStepBG, 1000, 10, 10, "TMultiPeak", "MultiStepBG"); // This is a weird nonsense line.
56 fBackground->SetNpx(1000);
57 fBackground->SetLineStyle(2);
58 fBackground->SetLineColor(kBlack);
60}
61
63{
64 delete fBackground;
65
66 for(auto& peak : fPeakVec) {
67 delete peak;
68 }
69}
70
71void TMultiPeak::SortPeaks(Bool_t (*SortFunction)(const TPeak*, const TPeak*))
72{
73 std::sort(fPeakVec.begin(), fPeakVec.end(), SortFunction);
74}
75
77{
78 SetParName(0, "N_Peaks");
79 SetParName(1, "A");
80 SetParName(2, "B");
81 SetParName(3, "C");
82 SetParName(4, "bg_offset");
83
84 for(int i = 0; i < static_cast<int>(fPeakVec.size()); ++i) {
85 SetParName(6 * i + 5, Form("Height_%i", i));
86 SetParName(6 * i + 6, Form("Centroid_%i", i));
87 SetParName(6 * i + 7, Form("Sigma_%i", i));
88 SetParName(6 * i + 8, Form("Beta_%i", i));
89 SetParName(6 * i + 9, Form("R_%i", i));
90 SetParName(6 * i + 10, Form("Step_%i", i));
91 }
92 FixParameter(0, static_cast<double>(fPeakVec.size()));
93}
94
96{
97 copy.Copy(*this);
98}
99
100void TMultiPeak::Copy(TObject& obj) const
101{
102 TGRSIFit::Copy(obj);
103 auto* mpobj = static_cast<TMultiPeak*>(&obj);
104 if((mpobj->fBackground) == nullptr) {
105 mpobj->fBackground = new TF1(*(fBackground));
106 } else {
107 *(mpobj->fBackground) = *fBackground;
108 }
109
111
112 // Copy all of the TPeaks.
113 for(auto* i : fPeakVec) {
114 auto* peak = new TPeak(*i);
115 peak->AddToGlobalList(kFALSE);
116 mpobj->fPeakVec.push_back(peak);
117 }
118}
119
120Bool_t TMultiPeak::InitParams(TH1* fithist)
121{
122 // Makes initial guesses at parameters for the fit. Uses the histogram to make the initial guesses
123 if((fithist == nullptr) && (GetHist() != nullptr)) {
124 fithist = GetHist();
125 }
126
127 if(fithist == nullptr) {
128 std::cout << "No histogram is associated yet, no initial guesses made" << std::endl;
129 return false;
130 }
131
132 FixParameter(0, static_cast<double>(fPeakVec.size()));
133 // This is the range for the fit.
134 Double_t xlow = 0;
135 Double_t xhigh = 0;
136 GetRange(xlow, xhigh);
137 Int_t binlow = fithist->GetXaxis()->FindBin(xlow);
138 Int_t binhigh = fithist->GetXaxis()->FindBin(xhigh);
139
140 // Initialize background
141 SetParLimits(1, 0.0, fithist->GetBinContent(binlow) * 100.0);
142 SetParameter("A", fithist->GetBinContent(binlow));
143 SetParameter("B", (fithist->GetBinContent(binlow) - fithist->GetBinContent(binhigh)) / (xlow - xhigh));
144 SetParameter("C", 0.0000);
145 SetParameter("bg_offset", (xhigh + xlow) / 2.0);
146
147 FixParameter(3, 0);
148
149 // We need to initialize parameters for every peak in the fit
150 for(int i = 0; i < static_cast<int>(fPeakVec.size()); ++i) {
151 Double_t centroid = fPeakVec.at(i)->GetCentroid();
152 Int_t bin = fithist->GetXaxis()->FindBin(centroid);
153 SetParLimits(6 * i + 5, 0, fithist->GetBinContent(bin) * 5.);
154 SetParLimits(6 * i + 6, centroid - 4, centroid + 4);
155 // SetParLimits(6*i+7,0.1,xhigh-xlow);//This will be linked to other peaks eventually.
156 SetParLimits(6 * i + 7, 0.1, 1.5); // This will be linked to other peaks eventually.
157 SetParLimits(6 * i + 8, 0.000001, 10);
158 SetParLimits(6 * i + 9, 0.000001, 100);
159 SetParLimits(6 * i + 10, 0.0, 1.0E2);
160 // Step size is allow to vary to anything. If it goes below 0, the code will fix it to 0
161
162 // Now set the actual paramter to start the fit from these points
163 SetParameter(Form("Height_%i", i), fithist->GetBinContent(bin));
164 SetParameter(Form("Centroid_%i", i), centroid);
165 // SetParameter("sigma",(xhigh-xlow)*0.5); // slightly more robust starting value for sigma -JKS
166 // SetParameter("sigma",1.0/binWidth); // slightly more robust starting value for sigma -JKS
167 SetParameter(Form("Sigma_%i", i), TMath::Sqrt(9.0 + 4. * GetParameter(Form("Centroid_%i", i)) / 1000.));
168 SetParameter(Form("Beta_%i", i), GetParameter(Form("Sigma_%i", i)) / 2.0);
169 SetParameter(Form("R_%i", i), 1.0);
170 SetParameter(Form("Step_%i", i), 1.0);
171
172 // Fix beta and R. These will be released if they are needed (or can be asked to be released).
173 FixParameter(6 * i + 8, GetParameter(Form("Beta_%i", i)));
174 FixParameter(6 * i + 9, 0.00);
175 }
176
178 return true;
179}
180
181Bool_t TMultiPeak::Fit(TH1* fithist, Option_t* opt)
182{
183 TString optstr = opt;
184 if((fithist == nullptr) && (GetHist() == nullptr)) {
185 std::cout << "No hist passed, trying something...";
186 fithist = fHistogram;
187 }
188 if(fithist == nullptr) {
189 std::cout << "No histogram associated with Peak" << std::endl;
190 return false;
191 }
192 if(!IsInitialized()) {
193 InitParams(fithist);
194 }
195
196 TVirtualFitter::SetMaxIterations(100000);
197 ROOT::Math::MinimizerOptions::SetDefaultMinimizer("Minuit2", "Combination");
198 SetHist(fithist);
199
200 TString options(opt);
201 bool print_flag = true;
202 if(options.Contains("Q")) {
203 print_flag = false;
204 }
205
206 // Now that it is initialized, let's fit it.
207
208 TFitResultPtr fitres;
209 // Log likelihood is the proper fitting technique UNLESS the data is a result of an addition or subtraction.
211 fitres = fithist->Fit(this, Form("%sLRSQ", opt)); // The RS needs to always be there
212 } else {
213 fitres = fithist->Fit(this, Form("%sRSQ", opt)); // The RS needs to always be there
214 }
215
216 std::vector<Double_t> sigma_list(static_cast<int>(GetParameter(0) + 0.5));
217 for(int i = 0; i < static_cast<int>(GetParameter(0) + 0.5); ++i) {
218 // Get Median sigma
219 sigma_list[i] = GetParameter(6 * i + 7);
220 }
221 std::sort(sigma_list.begin(), sigma_list.end(), std::greater<>());
222 Double_t median = sigma_list.at(sigma_list.size() / 2);
223
224 Double_t range_low = 0.;
225 Double_t range_high = 0.;
226 for(int i = 0; i < static_cast<int>(GetParameter(0) + 0.5); ++i) {
227 GetParLimits(6 * i + 7, range_low, range_high);
228 if(range_low != range_high) {
229 SetParLimits(6 * i + 7, median * .95, median * 1.05);
230 }
231 }
232
233 // Refit
235 fitres = fithist->Fit(this, Form("%sLRS", opt)); // The RS needs to always be there
236 } else {
237 fitres = fithist->Fit(this, Form("%sRS", opt)); // The RS needs to always be there
238 }
239
240 // After performing this fit I want to put something here that takes the fit result (good,bad,etc)
241 // for printing out. RD
242
243 // This removes the background parts of the fit form the integral error, while maintaining the covariance between the
244 // fits and the background.
245 TMatrixDSym CovMat = fitres->GetCovarianceMatrix();
246 CovMat(0, 0) = 0.0;
247 CovMat(1, 1) = 0.0;
248 CovMat(2, 2) = 0.0;
249 CovMat(3, 3) = 0.0;
250 CovMat(4, 4) = 0.0;
251
252 // This copies the parameters background but the background function doesn't have peaks
254 // We now make a copy of the covariance matrix that has completel 0 diagonals so that we can remove the other peaks
255 // form the integral error.
256 // This is done by adding back the peak of interest on the diagonal when it is integrated.
257 TMatrixDSym emptyCovMat = CovMat;
258 for(size_t i = 0; i < fPeakVec.size(); ++i) {
259 emptyCovMat(6 * i + 5, 6 * i + 5) = 0.0;
260 emptyCovMat(6 * i + 6, 6 * i + 6) = 0.0;
261 emptyCovMat(6 * i + 7, 6 * i + 7) = 0.0;
262 emptyCovMat(6 * i + 8, 6 * i + 8) = 0.0;
263 emptyCovMat(6 * i + 9, 6 * i + 9) = 0.0;
264 emptyCovMat(6 * i + 10, 6 * i + 10) = 0.0;
265 }
266
267 if(print_flag) {
268 std::cout << "Chi^2/NDF = " << fitres->Chi2() / fitres->Ndf() << std::endl;
269 }
270 // We will now set the parameters of each of the peaks based on the fits.
271 for(int i = 0; i < static_cast<int>(fPeakVec.size()); ++i) {
272 auto* tmpMp = new TMultiPeak(*this);
273 tmpMp->ClearParameters(); // We need to clear all of the parameters so that we can add the ones we want back in
274 Double_t binWidth = fithist->GetBinWidth(1); // This assumes all bins are the same width, so we can just take the first bin instead of the bin of the centroid.
275 TPeak* peak = fPeakVec.at(i);
276 TMatrixDSym tmpCovMat = emptyCovMat;
277 peak->SetParameter("Height", GetParameter(Form("Height_%i", i)));
278 peak->SetParameter("centroid", GetParameter(Form("Centroid_%i", i)));
279 peak->SetParameter("sigma", GetParameter(Form("Sigma_%i", i)));
280 peak->SetParameter("beta", GetParameter(Form("Beta_%i", i)));
281 peak->SetParameter("R", GetParameter(Form("R_%i", i)));
282 peak->SetParameter("step", GetParameter(Form("Step_%i", i)));
283 peak->SetParameter("A", 0.0);
284 peak->SetParameter("B", 0.0);
285 peak->SetParameter("C", 0.0);
286 peak->SetParameter("bg_offset", 0.0);
287 peak->SetChi2(fitres->Chi2());
288 peak->SetNdf(fitres->Ndf());
289
290 // Set the important diagonals for the integral of the covariance matrix
291 tmpCovMat(6 * i + 5, 6 * i + 5) = CovMat(6 * i + 5, 6 * i + 5);
292 tmpCovMat(6 * i + 6, 6 * i + 6) = CovMat(6 * i + 6, 6 * i + 6);
293 tmpCovMat(6 * i + 7, 6 * i + 7) = CovMat(6 * i + 7, 6 * i + 7);
294 tmpCovMat(6 * i + 8, 6 * i + 8) = CovMat(6 * i + 8, 6 * i + 8);
295 tmpCovMat(6 * i + 9, 6 * i + 9) = CovMat(6 * i + 9, 6 * i + 9);
296
297 tmpMp->SetParameter("N_Peaks", GetParameter("N_Peaks"));
298 tmpMp->SetParameter(Form("Height_%i", i), GetParameter(Form("Height_%i", i)));
299 tmpMp->SetParameter(Form("Centroid_%i", i), GetParameter(Form("Centroid_%i", i)));
300 tmpMp->SetParameter(Form("Sigma_%i", i), GetParameter(Form("Sigma_%i", i)));
301 tmpMp->SetParameter(Form("Beta_%i", i), GetParameter(Form("Beta_%i", i)));
302 tmpMp->SetParameter(Form("R_%i", i), GetParameter(Form("R_%i", i)));
303
304 Double_t width = GetParameter(Form("Sigma_%i", i));
305 Double_t xlow = 0.;
306 Double_t xhigh = 0.;
307 GetRange(xlow, xhigh);
308 Double_t int_low = xlow - 10. * width; // making the integration bounds a bit smaller, but still large enough. -JKS
309 Double_t int_high = xhigh + 10. * width;
310
311 // Make a function that does not include the background
312 // Intgrate the background.
313 tmpMp->SetRange(int_low, int_high); // This will help get the true area of the gaussian 200 ~ infinity in a gaus
314 // peak->SetName("tmppeak");
315
316 // This is where we will do integrals and stuff.
317 peak->SetArea((tmpMp->Integral(int_low, int_high)) / binWidth);
318 peak->SetAreaErr((tmpMp->IntegralError(int_low, int_high, tmpMp->GetParameters(), tmpCovMat.GetMatrixArray())) /
319 binWidth);
320 peak->SetParameter("centroid", GetParameter(Form("Centroid_%i", i)));
321 peak->SetParError(peak->GetParNumber("centroid"), GetParError(GetParNumber(Form("Centroid_%i", i))));
322 peak->SetParameter("sigma", GetParameter(GetParNumber(Form("Sigma_%i", i))));
323 peak->SetParError(peak->GetParNumber("sigma"), GetParError(GetParNumber(Form("Sigma_%i", i))));
324 if(print_flag) {
325 std::cout << "Integral: " << peak->GetArea() << " +- " << peak->GetAreaErr() << std::endl;
326 }
327 }
328
329 return true;
330}
331
332void TMultiPeak::Clear(Option_t* opt)
333{
334 TGRSIFit::Clear(opt);
335 for(auto& peak : fPeakVec) {
336 delete peak;
337 peak = nullptr;
338 }
339 fPeakVec.clear();
340}
341
342void TMultiPeak::Print(Option_t* opt) const
343{
344 /// Prints TMultiPeak properties. To see More properties use the option "+"
345 std::cout << "Name: " << GetName() << std::endl;
346 std::cout << "Number of Peaks: " << fPeakVec.size() << std::endl;
347 TF1::Print();
348 for(int i = 0; i < static_cast<int>(fPeakVec.size()); ++i) {
349 std::cout << "Peak: " << i << std::endl;
350 fPeakVec.at(i)->Print(opt);
351 std::cout << std::endl;
352 }
353}
354
355Double_t TMultiPeak::MultiPhotoPeakBG(Double_t* dim, Double_t* par)
356{
357 // Limits need to be imposed or error states may occour.
358 //
359 // General background.
360 int npeaks = static_cast<int>(par[0] + 0.5);
361 double result = TGRSIFunctions::PolyBg(dim, &par[1], 2); // polynomial background. uses par[1->4]
362 for(int i = 0; i < npeaks; ++i) { // par[0] is number of peaks
363 std::array<Double_t, 6> tmp_par;
364 tmp_par[0] = par[6 * i + 5]; // height of photopeak
365 tmp_par[1] = par[6 * i + 6]; // Peak Centroid of non skew gaus
366 tmp_par[2] = par[6 * i + 7]; // standard deviation of gaussian
367 tmp_par[3] = par[6 * i + 8]; //"skewedness" of the skewed gaussian
368 tmp_par[4] = par[6 * i + 9]; // relative height of skewed gaussian
369 tmp_par[5] = par[6 * i + 10]; // Size of step in background
370 result += TGRSIFunctions::PhotoPeak(dim, tmp_par.data()) + TGRSIFunctions::StepFunction(dim, tmp_par.data());
371 }
372 return result;
373}
374
375Double_t TMultiPeak::MultiStepBG(Double_t* dim, Double_t* par)
376{
377 // Limits need to be imposed or error states may occour.
378 //
379 // General background.
380 int npeaks = static_cast<int>(par[0] + 0.5);
381 double result = TGRSIFunctions::PolyBg(dim, &par[1], 2); // polynomial background. uses par[1->4]
382 for(int i = 0; i < npeaks; i++) { // par[0] is number of peaks
383 std::array<Double_t, 6> tmp_par;
384 tmp_par[0] = par[6 * i + 5]; // height of photopeak
385 tmp_par[1] = par[6 * i + 6]; // Peak Centroid of non skew gaus
386 tmp_par[2] = par[6 * i + 7]; // standard deviation of gaussian
387 tmp_par[3] = par[6 * i + 8]; //"skewedness" of the skewed gaussian
388 tmp_par[4] = par[6 * i + 9]; // relative height of skewed gaussian
389 tmp_par[5] = par[6 * i + 10]; // Size of step in background
390 result += TGRSIFunctions::StepFunction(dim, tmp_par.data());
391 }
392 return result;
393}
394
395Double_t TMultiPeak::SinglePeakBG(Double_t* dim, Double_t* par)
396{
397 // Limits need to be imposed or error states may occour.
398 //
399 // General background.
400
401 int npeaks = static_cast<int>(par[0] + 0.5);
402 double result = TGRSIFunctions::PolyBg(dim, &par[1], 2); // polynomial background. uses par[1->4]
403 for(int i = 0; i < npeaks; i++) { // par[0] is number of peaks
404 std::array<Double_t, 6> tmp_par;
405 tmp_par[0] = par[6 * i + 5]; // height of photopeak
406 tmp_par[1] = par[6 * i + 6]; // Peak Centroid of non skew gaus
407 tmp_par[2] = par[6 * i + 7]; // standard deviation of gaussian
408 tmp_par[3] = par[6 * i + 8]; //"skewedness" of the skewed gaussian
409 tmp_par[4] = par[6 * i + 9]; // relative height of skewed gaussian
410 tmp_par[5] = par[6 * i + 10]; // Size of step in background
411 result += TGRSIFunctions::StepFunction(dim, tmp_par.data());
412 }
413 result += TGRSIFunctions::PhotoPeak(dim, &par[npeaks * 6 + 5]);
414 return result;
415}
416
418{
419 if(idx < fPeakVec.size()) {
420 return fPeakVec.at(idx);
421 }
422 std::cout << "No matching peak at index " << idx << std::endl;
423
424 return nullptr;
425}
426
428{
429 size_t closest_idx = 0;
430 Double_t closest_so_far = 1000000.;
431 for(size_t i = 0; i < fPeakVec.size(); ++i) {
432 if(std::abs(energy - fPeakVec.at(i)->GetCentroid()) < closest_so_far) {
433 closest_so_far = std::abs(energy - fPeakVec.at(i)->GetCentroid());
434 closest_idx = i;
435 }
436 }
437
438 return GetPeak(closest_idx);
439}
440
442{
443 // Draws the individual TPeaks that make up the TMultiPeak. ROOT makes this a complicated process. The result on the
444 // histogram might have memory issues.
445 Double_t xlow = 0.;
446 Double_t xhigh = 0.;
447 GetRange(xlow, xhigh);
448 int npeaks = fPeakVec.size();
449 for(auto* peak : fPeakVec) {
450 // Should be good enough to draw between -2 and +2 fwhm
451 Double_t centroid = peak->GetCentroid();
452 Double_t range = 2. * peak->GetFWHM();
453
454 auto* sum = new TF1(Form("tmp%s", peak->GetName()), this, &TMultiPeak::SinglePeakBG, centroid - range, centroid + range,
455 fPeakVec.size() * 6 + 11, "TMultiPeak", "SinglePeakBG");
456
457 for(int j = 0; j < GetNpar(); ++j) {
458 sum->SetParameter(j, GetParameter(j));
459 }
460 for(int j = 0; j < 5; ++j) {
461 sum->SetParameter(npeaks * 6 + 5 + j, peak->GetParameter(j));
462 }
463
464 sum->SetNpx(1000);
465 sum->SetLineStyle(2);
466 sum->SetLineColor(kMagenta);
467
468 sum->SetRange(centroid - range, centroid + range);
469 sum->DrawCopy("SAME");
470 delete sum;
471 }
472}
Bool_t AddToGlobalList(Bool_t yes=kTRUE) override
Definition TGRSIFit.cxx:64
virtual void SetHist(TH1 *hist)
Definition TGRSIFit.h:42
void Clear(Option_t *opt="") override
Definition TGRSIFit.cxx:38
virtual void CopyParameters(TF1 *copy) const
Definition TGRSIFit.cxx:52
Bool_t IsInitialized() const
Definition TGRSIFit.h:56
void SetInitialized(Bool_t flag=true)
Definition TGRSIFit.h:57
void Copy(TObject &obj) const override
Definition TGRSIFit.cxx:22
virtual TH1 * GetHist() const
Definition TGRSIFit.h:46
TF1 * fBackground
Definition TMultiPeak.h:61
void InitNames()
TPeak * GetPeakClosestTo(Double_t energy)
static bool fLogLikelihoodFlag
!
Definition TMultiPeak.h:59
void Copy(TObject &obj) const override
Double_t MultiPhotoPeakBG(Double_t *dim, Double_t *par)
Double_t SinglePeakBG(Double_t *dim, Double_t *par)
Bool_t Fit(TH1 *fithist, Option_t *opt="")
TPeak * GetPeak(UInt_t idx)
std::vector< TPeak * > fPeakVec
Definition TMultiPeak.h:60
static bool GetLogLikelihoodFlag()
Definition TMultiPeak.h:49
void Print(Option_t *opt="") const override
bool InitParams(TH1 *fithist) override
void DrawPeaks()
void SortPeaks(Bool_t(*SortFunction)(const TPeak *, const TPeak *)=TPeak::CompareEnergy)
void Clear(Option_t *opt="") override
Double_t MultiStepBG(Double_t *dim, Double_t *par)
Definition TPeak.h:23
void SetNdf(Double_t Ndf)
Definition TPeak.h:87
void SetArea(Double_t area)
Definition TPeak.h:79
Double_t GetAreaErr() const
Definition TPeak.h:49
void SetAreaErr(Double_t areaErr)
Definition TPeak.h:80
void SetChi2(Double_t chi2)
Definition TPeak.h:86
Double_t GetArea() const
Definition TPeak.h:48
Double_t StepFunction(Double_t *dim, Double_t *par)
Double_t PhotoPeak(Double_t *dim, Double_t *par)
Double_t PolyBg(Double_t *x, Double_t *par, Int_t order)