8#include "Math/SpecFuncMathMore.h"
27bool
TGRSIFunctions::CheckParameterErrors(const TFitResultPtr& fitres, std::
string opt)
36 if(fitres.Get() ==
nullptr) { return true; }
37 std::transform(opt.begin(), opt.end(), opt.begin(), [](
unsigned char c) { return std::tolower(c); });
38 bool quiet = (opt.find(
'q') != std::string::npos);
39 bool verbose = (opt.find(
'v') != std::string::npos);
40 if(quiet && verbose) {
41 std::cout <<
"Don't know how to be quiet and verbose at once (" << opt <<
"), going to be verbose!" << std::endl;
47 auto covarianceMatrix = fitres->GetCovarianceMatrix();
49 if(covarianceMatrix.GetNrows() !=
static_cast<Int_t
>(fitres->NPar()) || covarianceMatrix.GetNcols() !=
static_cast<Int_t
>(fitres->NPar())) { return false; }
50 for(
unsigned int i = 0; i < fitres->NPar(); ++i) {
51 if(std::fabs(fitres->ParError(i) - TMath::Sqrt(covarianceMatrix[i][i])) > 0.1) {
53 std::cout <<
RED <<
"Parameter " << i <<
" = " << fitres->GetParameterName(i) <<
" is at or close to the limit, error is " << fitres->ParError(i) <<
", and square root of diagonal is " << TMath::Sqrt(covarianceMatrix[i][i]) <<
RESET_COLOR << std::endl;
59 std::cout <<
GREEN <<
"Parameter " << i <<
" = " << fitres->GetParameterName(i) <<
" is not close to the limit, error is " << fitres->ParError(i) <<
", and square root of diagonal is " << TMath::Sqrt(covarianceMatrix[i][i]) <<
RESET_COLOR << std::endl;
70 Double_t x = time[0] - par[0];
71 Double_t e = exp(-x / par[1]);
76 s += par[6] * (1 - exp(-x / par[2])) * e;
77 s += par[7] * (1 - exp(-x / par[3])) * e;
78 s += par[8] * (1 - exp(-x / par[4])) * e;
87 for(Int_t i = 0; i <= order; i++) {
88 result += par[i] * TMath::Power(x[0] - par[order + 1], i);
104 Double_t height = par[0];
105 Double_t centroid = par[1];
106 Double_t sigma = par[2];
107 Double_t step = par[5];
109 return TMath::Abs(step) * height / 100. * TMath::Erfc((x - centroid) / (TMath::Sqrt(2.) * sigma));
114 return StepFunction(dim, par) + PolyBg(dim, &par[6], 2);
120 return Gaus(dim, par) + SkewedGaus(dim, par);
126 double result = Gaus(dim, par) + SkewedGaus(dim, par) + StepFunction(dim, par) + PolyBg(dim, &par[6], 2);
127 if(std::isfinite(result)) {
return result; }
134 int nPeaks =
static_cast<int>(par[0] + 0.5);
135 double result = PolyBg(dim, &par[1], 2);
136 for(
int i = 0; i < nPeaks; i++) {
137 std::array<Double_t, 6> tmpPar;
138 tmpPar[0] = par[6 * i + 5];
139 tmpPar[1] = par[6 * i + 6];
140 tmpPar[2] = par[6 * i + 7];
141 tmpPar[3] = par[6 * i + 8];
142 tmpPar[4] = par[6 * i + 9];
143 tmpPar[5] = par[6 * i + 10];
144 result += PhotoPeak(dim, tmpPar.data()) + StepFunction(dim, tmpPar.data());
161 Double_t height = par[0];
162 Double_t centroid = par[1];
163 Double_t sigma = par[2];
164 Double_t relHeight = par[4];
166 return height * (1. - relHeight / 100.) * TMath::Gaus(x, centroid, sigma);
182 Double_t height = par[0];
183 Double_t centroid = par[1];
184 Double_t sigma = par[2];
185 Double_t beta = par[3];
186 Double_t relHeight = par[4];
191 return relHeight * height / 100. * (TMath::Exp((x - centroid) / beta)) *
192 (TMath::Erfc(((x - centroid) / (TMath::Sqrt(2.) * sigma)) + sigma / (TMath::Sqrt(2.) * beta)));
201 double result = par[1] + dim[0] * par[2];
202 for(
int i = 0; i < par[0]; i++) {
203 std::array<Double_t, 5> tmpPar;
204 tmpPar[0] = par[5 * i + 3];
205 tmpPar[1] = par[5 * i + 4];
206 tmpPar[2] = par[5 * i + 5];
207 tmpPar[3] = par[5 * i + 6];
208 tmpPar[4] = par[5 * i + 7];
209 result += SkewedGaus(dim, tmpPar.data());
225 return par[0] * (TMath::Exp((x[0] - par[1]) / par[3])) *
226 (TMath::Erfc(((x[0] - par[1]) / (TMath::Sqrt(2.) * par[2])) + par[2] / (TMath::Sqrt(2.) * par[3])));
235 double result = par[1] + dim[0] * par[2];
236 for(
int i = 0; i < static_cast<int>(par[0] + 0.5); i++) {
237 std::array<Double_t, 4> tmpPar;
238 tmpPar[0] = par[4 * i + 3];
239 tmpPar[1] = par[4 * i + 4];
240 tmpPar[2] = par[4 * i + 5];
241 tmpPar[3] = par[4 * i + 6];
256 for(
int i = 0; i < 10; i++) {
257 dy = 5 * pars[3] / 10.;
259 y = x[0] - 2.5 * pars[3] + dy * i;
263 for(
int n = 0; n < static_cast<int>(pars[0] + 0.5); n++) {
265 spec += pars[3 * n + 4] * TMath::Landau(-y, -pars[3 * n + 5], pars[3 * n + 6]) /
266 TMath::Landau(0, 0, 100);
269 gaus = TMath::Gaus(-x[0], -y, pars[3]) /
270 sqrt(2 * TMath::Pi() * pars[3] * pars[3]);
271 conv += gaus * spec * dy;
285 for(
int i = 0; i < 50; i++) {
286 dy = 8 * pars[3] / 50.;
287 y = x[0] - 4 * pars[3] + dy * i;
289 spec = pars[1] + pars[2] * y;
290 for(
int n = 0; n < static_cast<int>(pars[0] + 0.5); n++) {
291 spec += pars[3 * n + 4] * TMath::Landau(-y, -pars[3 * n + 5], pars[3 * n + 6]) / TMath::Landau(0, 0, 100);
294 gaus = TMath::Gaus(-x[0], -y, pars[3]) / sqrt(2 * TMath::Pi() * pars[3] * pars[3]);
295 conv += gaus * spec * dy;
309 double result = par[1] + dim[0] * par[2];
310 for(
int i = 0; i < static_cast<int>(par[0] + 0.5); i++) {
311 amp = par[3 * i + 3];
312 mean = par[3 * i + 4];
313 sigma = par[3 * i + 5];
314 result += amp * TMath::Gaus(dim[0], mean, sigma);
333 if(par.size() < (nChain * 3)) {
334 std::cout <<
"not enough parameters passed to function" << std::endl;
338 Double_t totalActivity = 0.;
342 for(UInt_t n = 0; n < nChain; n++) {
344 Double_t firstterm = 1.;
346 for(UInt_t j = 0; j < n - 1; j++) {
347 firstterm *= par[1 + 3 * j];
349 Double_t secondterm = 0.;
350 for(UInt_t i = 0; i < n; i++) {
352 for(UInt_t j = i; j < n; j++) {
354 for(UInt_t p = i; p < n; p++) {
356 denom *= par[1 + 3 * p] - par[1 + 3 * j];
359 sum += par[0 + 3 * i] / par[1 + 3 * i] * TMath::Exp(-par[1 + 3 * j] * dim[0]) / denom;
363 par[2 + 3 * n] = par[1 + 3 * n] * firstterm * secondterm;
364 totalActivity += par[2 + 3 * n];
366 return totalActivity;
375 return dim[0] / (1. - dim[0] * deadtime / (binWidth * 1000000.));
384 return function / (1. + function * deadtime / (binWidth * 1000000.));
390 Double_t val = p[0] * (1 + p[1] * ::ROOT::Math::legendre(2, x[0]) + p[2] * ::ROOT::Math::legendre(4, x[0]));
393 std::cout <<
"Mathmore feature of ROOT is missing, " << __PRETTY_FUNCTION__ <<
" will always return 1!" << std::endl;
401 for(
int i = 0; i < 9; ++i) {
402 sum += par[i] * TMath::Power(TMath::Log(dim[0]), i);
405 return TMath::Exp(sum);
418 return TMath::Sqrt(TMath::Pi()) * par[0] * par[3] / 2 * TMath::Exp(par[3] / 2 * (2 * par[1] + par[3] * TMath::Power(par[2], 2) - 2 * x[0])) * TMath::Erfc((par[1] + par[3] * TMath::Power(par[2], 2) - x[0]) / (TMath::Sqrt(2) * par[2])) + par[4];
432 return TMath::Sqrt(TMath::Pi()) * par[0] * par[3] / 2 * TMath::Exp(par[3] / 2 * (2 * par[1] + par[3] * TMath::Power(par[2], 2) - 2 * x[0])) * TMath::Erfc((par[1] + par[3] * TMath::Power(par[2], 2) - x[0]) / (TMath::Sqrt(2) * par[2])) + TMath::Sqrt(TMath::Pi()) * par[4] * par[5] / 2 * TMath::Exp(par[5] / 2 * (2 * par[1] + par[5] * TMath::Power(par[2], 2) - 2 * x[0])) * TMath::Erfc((par[1] + par[5] * TMath::Power(par[2], 2) - x[0]) / (TMath::Sqrt(2) * par[2]));
441 if(2 * j1 != floor(2 * j1) ||
442 2 * j2 != floor(2 * j2) ||
443 2 * j != floor(2 * j) ||
444 2 * m1 != floor(2 * m1) ||
445 2 * m2 != floor(2 * m2) ||
446 2 * m != floor(2 * m)) {
452 if((j1 - m1) != floor(j1 - m1)) {
455 if((j2 - m2) != floor(j2 - m2)) {
458 if(j - m != floor(j - m)) {
461 if(j > (j1 + j2) || j < abs(j1 - j2)) {
474 double term1 = TMath::Sqrt((((2 * j + 1) / TMath::Factorial(
static_cast<Int_t
>(j1 + j2 + j + 1))) * TMath::Factorial(
static_cast<Int_t
>(j2 + j - j1)) * TMath::Factorial(
static_cast<Int_t
>(j + j1 - j2)) * TMath::Factorial(
static_cast<Int_t
>(j1 + j2 - j)) * TMath::Factorial(
static_cast<Int_t
>(j1 + m1)) * TMath::Factorial(
static_cast<Int_t
>(j1 - m1)) * TMath::Factorial(
static_cast<Int_t
>(j2 + m2)) * TMath::Factorial(
static_cast<Int_t
>(j2 - m2)) * TMath::Factorial(
static_cast<Int_t
>(j + m)) * TMath::Factorial(
static_cast<Int_t
>(j - m))));
476 for(
int k = 0; k <= 99; k++) {
477 if((j1 + j2 - j - k < 0) || (j1 - m1 - k < 0) || (j2 + m2 - k < 0)) {
481 if((j - j1 - m2 + k < 0) || (j - j2 + m1 + k < 0)) {
483 const auto a1 =
static_cast<Int_t
>(j - j1 - m2);
484 const auto a2 =
static_cast<Int_t
>(j - j2 + m1);
485 k = TMath::Max(-TMath::Min(a1, a2) - 1, k);
487 double term = TMath::Factorial(
static_cast<Int_t
>(j1 + j2 - j - k)) * TMath::Factorial(
static_cast<Int_t
>(j - j1 - m2 + k)) * TMath::Factorial(
static_cast<Int_t
>(j - j2 + m1 + k)) * TMath::Factorial(
static_cast<Int_t
>(j1 - m1 - k)) * TMath::Factorial(
static_cast<Int_t
>(j2 + m2 - k)) * TMath::Factorial(k);
504 return TMath::Power(-1.,
static_cast<int>(a + b + d + c)) * ::ROOT::Math::wigner_6j(
static_cast<int>(2 * a),
static_cast<int>(2 * b),
static_cast<int>(2 * e),
static_cast<int>(2 * d),
static_cast<int>(2 * c),
static_cast<int>(2 * f));
506 std::cout <<
"Mathmore feature of ROOT is missing, " << __PRETTY_FUNCTION__ <<
" will always return 1!" << std::endl;
521 return TMath::Power((-1), (jf - ji - 1)) * TMath::Power((2 * L1 + 1) * (2 * L2 + 1) * (2 * ji + 1), (1. / 2.)) * cg * w;
528 double f1 =
F(k, ji, L1, L1, jf);
529 double f2 =
F(k, ji, L1, L2, jf);
530 double f3 =
F(k, ji, L2, L2, jf);
532 return (1. / (1. + TMath::Power(delta, 2))) * (f1 + 2 * delta * f2 + delta * delta * f3);
537 double f1 =
F(k, jf, L1, L1, ji);
538 double f2 =
F(k, jf, L1, L2, ji);
539 double f3 =
F(k, jf, L2, L2, ji);
540 return (1. / (1. + TMath::Power(delta, 2))) * (f1 + TMath::Power((-1), L1 + L2) * 2 * delta * f2 + delta * delta * f3);
545 return B(2, j2, j1, l1a, l1b, delta1) *
A(2, j3, j2, l2a, l2b, delta2);
550 return B(4, j2, j1, l1a, l1b, delta1) *
A(4, j3, j2, l2a, l2b, delta2);
NamespaceImp(GRootFunctions) Double_t GRootFunctions
Double_t MultiPhotoPeakBG(Double_t *dim, Double_t *par)
Double_t MultiSkewedGausWithBG(Double_t *dim, Double_t *par)
Double_t MultiSkewedGausWithBG2(Double_t *dim, Double_t *par)
double F(double k, double jf, double L1, double L2, double ji)
Double_t MultiGausWithBG(Double_t *dim, Double_t *par)
Double_t StepFunction(Double_t *dim, Double_t *par)
double CalculateA2(double j1, double j2, double j3, double l1a, double l1b, double l2a, double l2b, double delta1, double delta2)
Double_t PhotoPeak(Double_t *dim, Double_t *par)
Double_t LegendrePolynomial(Double_t *x, Double_t *p)
Double_t CsIFitFunction(Double_t *time, Double_t *par)
Double_t Bateman(std::vector< Double_t > &dim, std::vector< Double_t > &par, UInt_t nChain=1, Double_t SecondsPerBin=1.0)
Double_t LanGausHighRes(Double_t *x, Double_t *pars)
double CalculateA4(double j1, double j2, double j3, double l1a, double l1b, double l2a, double l2b, double delta1, double delta2)
double RacahW(double a, double b, double c, double d, double e, double f)
Double_t PhotoPeakBG(Double_t *dim, Double_t *par)
Double_t ConvolutedDecay(Double_t *x, Double_t *par)
Double_t DeadTimeAffect(Double_t function, Double_t deadtime, Double_t binWidth=1.0)
Double_t PolyBg(Double_t *x, Double_t *par, Int_t order)
double B(double k, double ji, double jf, double L1, double L2, double delta)
Double_t SkewedGaus2(Double_t *x, Double_t *par)
Double_t DeadTimeCorrect(Double_t *dim, Double_t deadtime, Double_t binWidth=1.0)
double A(double k, double ji, double jf, double L1, double L2, double delta)
Double_t StepBG(Double_t *dim, Double_t *par)
double ClebschGordan(double j1, double m1, double j2, double m2, double j, double m)
Double_t Gaus(Double_t *dim, Double_t *par)
Double_t PhotoEfficiency(Double_t *dim, Double_t *par)
Double_t SkewedGaus(Double_t *dim, Double_t *par)
Double_t ConvolutedDecay2(Double_t *x, Double_t *par)
Double_t LanGaus(Double_t *x, Double_t *pars)