68#define GSL_FN_EVAL(F,x) (*((F)->function))(x,(F)->params)
95 double epsabs,
double epsrel,
size_t limit,
98 double *
result,
double *abserr);
103 double epsabs,
double epsrel,
size_t limit,
105 double *
result,
double * abserr) ;
109 double epsabs,
double epsrel,
size_t limit,
111 double *
result,
double *abserr) ;
116 double epsabs,
double epsrel,
size_t limit,
118 double *
result,
double *abserr) ;
123 double epsabs,
double epsrel,
size_t limit,
125 double *
result,
double *abserr) ;
135namespace RooFit_internal {
138 static void registerIntegrator()
145struct Roo_reg_AGKInteg1D {
146 Roo_reg_AGKInteg1D() { Roo_internal_AGKInteg1D::registerIntegrator(); }
158 RooRealVar maxSeg(
"maxSeg",
"maximum number of segments", 100);
159 RooCategory method(
"method",
"Integration method for each segment");
170 return std::make_unique<RooAdaptiveGaussKronrodIntegrator1D>(function, config);
173 fact.
registerPlugin(
"RooAdaptiveGaussKronrodIntegrator1D", creator, {maxSeg, method},
179 oocoutI(
nullptr,Integration) <<
"RooAdaptiveGaussKronrodIntegrator1D has been registered " << std::endl;
188 :
RooAbsIntegrator(function), _useIntegrandLimits(true), _epsAbs(config.epsRel()), _epsRel(config.epsAbs())
206 _useIntegrandLimits(false),
207 _epsAbs(config.epsRel()),
208 _epsRel(config.epsAbs()),
255 oocoutE(
nullptr, Integration) <<
"RooAdaptiveGaussKronrodIntegrator1D::setLimits: cannot override integrand's limits" << endl;
282 if (!infLo && !infHi) {
284 }
else if (infLo && infHi) {
286 }
else if (infLo && !infHi) {
382#define GSL_EBADTOL 13
384#define GSL_ERROR(a,b) oocoutE(nullptr,Integration) << "RooAdaptiveGaussKronrodIntegrator1D::integral() ERROR: " << a << endl ; return b ;
385#define GSL_DBL_MIN 2.2250738585072014e-308
386#define GSL_DBL_MAX 1.7976931348623157e+308
387#define GSL_DBL_EPSILON 2.2204460492503131e-16
390#define GSL_EMAXITER 3
393#define GSL_EDIVERGE 6
396#define GSL_ERROR_VAL(reason, gsl_errno, value) return value ;
398#define GSL_MAX(a,b) ((a) > (b) ? (a) : (b))
414#define GSL_COERCE_DBL(x) (gsl_coerce_double(x))
427 double *
result,
double *abserr,
428 double *defabs,
double *resabs);
431 double *
result,
double *abserr,
432 double *resabs,
double *resasc);
435 double *
result,
double *abserr,
436 double *resabs,
double *resasc);
439 double *
result,
double *abserr,
440 double *resabs,
double *resasc);
443 double *
result,
double *abserr,
444 double *resabs,
double *resasc);
447 double *
result,
double *abserr,
448 double *resabs,
double *resasc);
451 double *
result,
double *abserr,
452 double *resabs,
double *resasc);
455 double *cheb12,
double *cheb24);
473 const double wg[],
const double wgk[],
474 double fv1[],
double fv2[],
476 double *
result,
double * abserr,
477 double * resabs,
double * resasc);
482 double epsabs,
double epsrel,
size_t limit,
485 double *
result,
double *abserr);
496 workspace->
nrmax = 0;
500 workspace->
rlist[0] = 0.0;
501 workspace->
elist[0] = 0.0;
502 workspace->
order[0] = 0;
503 workspace->
level[0] = 0;
511 double result,
double error);
515 double result,
double error)
519 workspace->
elist[0] = error;
529 const size_t last = workspace->
size - 1;
530 const size_t limit = workspace->
limit;
532 double * elist = workspace->
elist;
533 size_t * order = workspace->
order;
541 size_t i_nrmax = workspace->
nrmax;
542 size_t i_maxerr = order[i_nrmax] ;
550 workspace->
i = i_maxerr ;
554 errmax = elist[i_maxerr] ;
561 while (i_nrmax > 0 && errmax > elist[order[i_nrmax - 1]])
563 order[i_nrmax] = order[i_nrmax - 1] ;
571 if(last < (limit/2 + 2))
577 top = limit - last + 1;
588 while (i < top && errmax < elist[order[i]])
590 order[i-1] = order[i] ;
594 order[i-1] = i_maxerr ;
598 errmin = elist[last] ;
602 while (k > i - 2 && errmin >= elist[order[k]])
604 order[k+1] = order[k] ;
612 i_maxerr = order[i_nrmax] ;
614 workspace->
i = i_maxerr ;
615 workspace->
nrmax = i_nrmax ;
623 double a1,
double b1,
double area1,
double error1,
624 double a2,
double b2,
double area2,
double error2);
628 double *
a,
double *
b,
double *
r,
double *
e);
634 double a1,
double b1,
double area1,
double error1,
635 double a2,
double b2,
double area2,
double error2)
637 double * alist = workspace->
alist ;
638 double * blist = workspace->
blist ;
639 double * rlist = workspace->
rlist ;
640 double * elist = workspace->
elist ;
641 size_t * level = workspace->
level ;
643 const size_t i_max = workspace->
i ;
644 const size_t i_new = workspace->
size ;
646 const size_t new_level = workspace->
level[i_max] + 1;
653 rlist[i_max] = area2;
654 elist[i_max] = error2;
655 level[i_max] = new_level;
659 rlist[i_new] = area1;
660 elist[i_new] = error1;
661 level[i_new] = new_level;
666 rlist[i_max] = area1;
667 elist[i_max] = error1;
668 level[i_max] = new_level;
672 rlist[i_new] = area2;
673 elist[i_new] = error2;
674 level[i_new] = new_level;
689 double *
a,
double *
b,
double *
r,
double *
e)
691 const size_t i = workspace->
i;
692 double * alist = workspace->
alist;
693 double * blist = workspace->
blist;
694 double * rlist = workspace->
rlist;
695 double * elist = workspace->
elist;
709 const double *
const rlist = workspace->
rlist ;
710 const size_t n = workspace->
size;
713 double result_sum = 0;
715 for (k = 0; k <
n; k++)
717 result_sum += rlist[k];
732 double tmp = (1 + 100 *
e) * (std::abs(a2) + 1000 * u);
734 int status = std::abs(a1) <= tmp && std::abs(b2) <= tmp;
742 const double a,
const double b,
743 const double epsabs,
const double epsrel,
746 double *
result,
double * abserr,
752 double epsabs,
double epsrel,
size_t limit,
755 double *
result,
double * abserr)
791 status =
qag (
f,
a,
b, epsabs, epsrel, limit,
801 const double a,
const double b,
802 const double epsabs,
const double epsrel,
805 double *
result,
double *abserr,
815 size_t iteration = 0;
816 int roundoff_type1 = 0;
817 int roundoff_type2 = 0;
829 if (limit > workspace->
limit)
834 if (epsabs <= 0 && (epsrel < 50 *
GSL_DBL_EPSILON || epsrel < 0.5e-28))
836 GSL_ERROR (
"tolerance cannot be acheieved with given epsabs and epsrel",
842 q (
f,
a,
b, &result0, &abserr0, &resabs0, &resasc0);
848 tolerance =
GSL_MAX_DBL (epsabs, epsrel * std::abs(result0));
854 if (abserr0 <= round_off && abserr0 > tolerance)
859 GSL_ERROR (
"cannot reach tolerance because of roundoff error "
862 else if ((abserr0 <= tolerance && abserr0 != resasc0) || abserr0 == 0.0)
905 retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
908 b1 = 0.5 * (a_i + b_i);
912 q (
f, a1, b1, &area1, &error1, &resabs1, &resasc1);
913 q (
f, a2, b2, &area2, &error2, &resabs2, &resasc2);
915 area12 = area1 + area2;
916 error12 = error1 + error2;
918 errsum += (error12 - e_i);
919 area += area12 - r_i;
921 if (resasc1 != error1 && resasc2 != error2)
923 double delta = r_i - area12;
925 if (std::abs(delta) <= 1.0e-5 * std::abs(area12) && error12 >= 0.99 * e_i)
929 if (iteration >= 10 && error12 > e_i)
935 tolerance =
GSL_MAX_DBL (epsabs, epsrel * std::abs(area));
937 if (errsum > tolerance)
939 if (roundoff_type1 >= 6 || roundoff_type2 >= 20)
953 update (workspace, a1, b1, area1, error1, a2, b2, area2, error2);
955 retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
960 while (iteration < limit && !error_type && errsum > tolerance);
965 if (errsum <= tolerance)
969 else if (error_type == 2)
971 GSL_ERROR (
"roundoff error prevents tolerance from being achieved",
974 else if (error_type == 3)
976 GSL_ERROR (
"bad integrand behavior found in the integration interval",
979 else if (iteration == limit)
989static double rescale_error (
double err,
const double result_abs,
const double result_asc) ;
994 err = std::abs(err) ;
996 if (result_asc != 0 && err != 0)
998 double scale = std::pow((200 * err / result_asc), 1.5) ;
1002 err = result_asc * scale ;
1025 const double xgk[],
const double wg[],
const double wgk[],
1026 double fv1[],
double fv2[],
1028 double *
result,
double *abserr,
1029 double *resabs,
double *resasc)
1032 const double center = 0.5 * (
a +
b);
1033 const double half_length = 0.5 * (
b -
a);
1034 const double abs_half_length = std::abs(half_length);
1037 double result_gauss = 0;
1038 double result_kronrod = f_center * wgk[
n - 1];
1040 double result_abs = std::abs(result_kronrod);
1041 double result_asc = 0;
1049 result_gauss = f_center * wg[
n / 2 - 1];
1052 for (j = 0; j < (
n - 1) / 2; j++)
1054 const int jtw = j * 2 + 1;
1055 const double abscissa = half_length * xgk[jtw];
1056 const double fval1 =
GSL_FN_EVAL (
f, center - abscissa);
1057 const double fval2 =
GSL_FN_EVAL (
f, center + abscissa);
1058 const double fsum = fval1 + fval2;
1061 result_gauss += wg[j] * fsum;
1062 result_kronrod += wgk[jtw] * fsum;
1063 result_abs += wgk[jtw] * (std::abs(fval1) + std::abs(fval2));
1066 for (j = 0; j <
n / 2; j++)
1069 const double abscissa = half_length * xgk[jtwm1];
1070 const double fval1 =
GSL_FN_EVAL (
f, center - abscissa);
1071 const double fval2 =
GSL_FN_EVAL (
f, center + abscissa);
1074 result_kronrod += wgk[jtwm1] * (fval1 + fval2);
1075 result_abs += wgk[jtwm1] * (std::abs(fval1) + std::abs(fval2));
1078 mean = result_kronrod * 0.5;
1080 result_asc = wgk[
n - 1] * std::abs(f_center - mean);
1082 for (j = 0; j <
n - 1; j++)
1084 result_asc += wgk[j] * (std::abs(fv1[j] - mean) + std::abs(fv2[j] - mean));
1089 err = (result_kronrod - result_gauss) * half_length;
1091 result_kronrod *= half_length;
1092 result_abs *= abs_half_length;
1093 result_asc *= abs_half_length;
1095 *
result = result_kronrod;
1096 *resabs = result_abs;
1097 *resasc = result_asc;
1108 0.991455371120812639206854697526329,
1109 0.949107912342758524526189684047851,
1110 0.864864423359769072789712788640926,
1111 0.741531185599394439863864773280788,
1112 0.586087235467691130294144838258730,
1113 0.405845151377397166906606412076961,
1114 0.207784955007898467600689403773245,
1115 0.000000000000000000000000000000000
1123 0.129484966168869693270611432679082,
1124 0.279705391489276667901467771423780,
1125 0.381830050505118944950369775488975,
1126 0.417959183673469387755102040816327
1131 0.022935322010529224963732008058970,
1132 0.063092092629978553290700663189204,
1133 0.104790010322250183839876322541518,
1134 0.140653259715525918745189590510238,
1135 0.169004726639267902826583426598550,
1136 0.190350578064785409913256402421014,
1137 0.204432940075298892414161999234649,
1138 0.209482141084727828012999174891714
1143 double *
result,
double *abserr,
1144 double *resabs,
double *resasc)
1149 gsl_integration_qk (8,
xgkA,
wgA,
wgkA, fv1, fv2,
f,
a,
b,
result, abserr, resabs, resasc);
1158 0.995657163025808080735527280689003,
1159 0.973906528517171720077964012084452,
1160 0.930157491355708226001207180059508,
1161 0.865063366688984510732096688423493,
1162 0.780817726586416897063717578345042,
1163 0.679409568299024406234327365114874,
1164 0.562757134668604683339000099272694,
1165 0.433395394129247190799265943165784,
1166 0.294392862701460198131126603103866,
1167 0.148874338981631210884826001129720,
1168 0.000000000000000000000000000000000
1176 0.066671344308688137593568809893332,
1177 0.149451349150580593145776339657697,
1178 0.219086362515982043995534934228163,
1179 0.269266719309996355091226921569469,
1180 0.295524224714752870173892994651338
1185 0.011694638867371874278064396062192,
1186 0.032558162307964727478818972459390,
1187 0.054755896574351996031381300244580,
1188 0.075039674810919952767043140916190,
1189 0.093125454583697605535065465083366,
1190 0.109387158802297641899210590325805,
1191 0.123491976262065851077958109831074,
1192 0.134709217311473325928054001771707,
1193 0.142775938577060080797094273138717,
1194 0.147739104901338491374841515972068,
1195 0.149445554002916905664936468389821
1201 double *
result,
double *abserr,
1202 double *resabs,
double *resasc)
1207 gsl_integration_qk (11,
xgkB,
wgB,
wgkB, fv1, fv2,
f,
a,
b,
result, abserr, resabs, resasc);
1216 0.998002298693397060285172840152271,
1217 0.987992518020485428489565718586613,
1218 0.967739075679139134257347978784337,
1219 0.937273392400705904307758947710209,
1220 0.897264532344081900882509656454496,
1221 0.848206583410427216200648320774217,
1222 0.790418501442465932967649294817947,
1223 0.724417731360170047416186054613938,
1224 0.650996741297416970533735895313275,
1225 0.570972172608538847537226737253911,
1226 0.485081863640239680693655740232351,
1227 0.394151347077563369897207370981045,
1228 0.299180007153168812166780024266389,
1229 0.201194093997434522300628303394596,
1230 0.101142066918717499027074231447392,
1231 0.000000000000000000000000000000000
1239 0.030753241996117268354628393577204,
1240 0.070366047488108124709267416450667,
1241 0.107159220467171935011869546685869,
1242 0.139570677926154314447804794511028,
1243 0.166269205816993933553200860481209,
1244 0.186161000015562211026800561866423,
1245 0.198431485327111576456118326443839,
1246 0.202578241925561272880620199967519
1251 0.005377479872923348987792051430128,
1252 0.015007947329316122538374763075807,
1253 0.025460847326715320186874001019653,
1254 0.035346360791375846222037948478360,
1255 0.044589751324764876608227299373280,
1256 0.053481524690928087265343147239430,
1257 0.062009567800670640285139230960803,
1258 0.069854121318728258709520077099147,
1259 0.076849680757720378894432777482659,
1260 0.083080502823133021038289247286104,
1261 0.088564443056211770647275443693774,
1262 0.093126598170825321225486872747346,
1263 0.096642726983623678505179907627589,
1264 0.099173598721791959332393173484603,
1265 0.100769845523875595044946662617570,
1266 0.101330007014791549017374792767493
1271 double *
result,
double *abserr,
1272 double *resabs,
double *resasc)
1277 gsl_integration_qk (16,
xgkC,
wgC,
wgkC, fv1, fv2,
f,
a,
b,
result, abserr, resabs, resasc);
1286 0.998859031588277663838315576545863,
1287 0.993128599185094924786122388471320,
1288 0.981507877450250259193342994720217,
1289 0.963971927277913791267666131197277,
1290 0.940822633831754753519982722212443,
1291 0.912234428251325905867752441203298,
1292 0.878276811252281976077442995113078,
1293 0.839116971822218823394529061701521,
1294 0.795041428837551198350638833272788,
1295 0.746331906460150792614305070355642,
1296 0.693237656334751384805490711845932,
1297 0.636053680726515025452836696226286,
1298 0.575140446819710315342946036586425,
1299 0.510867001950827098004364050955251,
1300 0.443593175238725103199992213492640,
1301 0.373706088715419560672548177024927,
1302 0.301627868114913004320555356858592,
1303 0.227785851141645078080496195368575,
1304 0.152605465240922675505220241022678,
1305 0.076526521133497333754640409398838,
1306 0.000000000000000000000000000000000
1314 0.017614007139152118311861962351853,
1315 0.040601429800386941331039952274932,
1316 0.062672048334109063569506535187042,
1317 0.083276741576704748724758143222046,
1318 0.101930119817240435036750135480350,
1319 0.118194531961518417312377377711382,
1320 0.131688638449176626898494499748163,
1321 0.142096109318382051329298325067165,
1322 0.149172986472603746787828737001969,
1323 0.152753387130725850698084331955098
1328 0.003073583718520531501218293246031,
1329 0.008600269855642942198661787950102,
1330 0.014626169256971252983787960308868,
1331 0.020388373461266523598010231432755,
1332 0.025882133604951158834505067096153,
1333 0.031287306777032798958543119323801,
1334 0.036600169758200798030557240707211,
1335 0.041668873327973686263788305936895,
1336 0.046434821867497674720231880926108,
1337 0.050944573923728691932707670050345,
1338 0.055195105348285994744832372419777,
1339 0.059111400880639572374967220648594,
1340 0.062653237554781168025870122174255,
1341 0.065834597133618422111563556969398,
1342 0.068648672928521619345623411885368,
1343 0.071054423553444068305790361723210,
1344 0.073030690332786667495189417658913,
1345 0.074582875400499188986581418362488,
1346 0.075704497684556674659542775376617,
1347 0.076377867672080736705502835038061,
1348 0.076600711917999656445049901530102
1353 double *
result,
double *abserr,
1354 double *resabs,
double *resasc)
1359 gsl_integration_qk (21,
xgkD,
wgD,
wgkD, fv1, fv2,
f,
a,
b,
result, abserr, resabs, resasc);
1368 0.999262104992609834193457486540341,
1369 0.995556969790498097908784946893902,
1370 0.988035794534077247637331014577406,
1371 0.976663921459517511498315386479594,
1372 0.961614986425842512418130033660167,
1373 0.942974571228974339414011169658471,
1374 0.920747115281701561746346084546331,
1375 0.894991997878275368851042006782805,
1376 0.865847065293275595448996969588340,
1377 0.833442628760834001421021108693570,
1378 0.797873797998500059410410904994307,
1379 0.759259263037357630577282865204361,
1380 0.717766406813084388186654079773298,
1381 0.673566368473468364485120633247622,
1382 0.626810099010317412788122681624518,
1383 0.577662930241222967723689841612654,
1384 0.526325284334719182599623778158010,
1385 0.473002731445714960522182115009192,
1386 0.417885382193037748851814394594572,
1387 0.361172305809387837735821730127641,
1388 0.303089538931107830167478909980339,
1389 0.243866883720988432045190362797452,
1390 0.183718939421048892015969888759528,
1391 0.122864692610710396387359818808037,
1392 0.061544483005685078886546392366797,
1393 0.000000000000000000000000000000000
1401 0.011393798501026287947902964113235,
1402 0.026354986615032137261901815295299,
1403 0.040939156701306312655623487711646,
1404 0.054904695975835191925936891540473,
1405 0.068038333812356917207187185656708,
1406 0.080140700335001018013234959669111,
1407 0.091028261982963649811497220702892,
1408 0.100535949067050644202206890392686,
1409 0.108519624474263653116093957050117,
1410 0.114858259145711648339325545869556,
1411 0.119455763535784772228178126512901,
1412 0.122242442990310041688959518945852,
1413 0.123176053726715451203902873079050
1418 0.001987383892330315926507851882843,
1419 0.005561932135356713758040236901066,
1420 0.009473973386174151607207710523655,
1421 0.013236229195571674813656405846976,
1422 0.016847817709128298231516667536336,
1423 0.020435371145882835456568292235939,
1424 0.024009945606953216220092489164881,
1425 0.027475317587851737802948455517811,
1426 0.030792300167387488891109020215229,
1427 0.034002130274329337836748795229551,
1428 0.037116271483415543560330625367620,
1429 0.040083825504032382074839284467076,
1430 0.042872845020170049476895792439495,
1431 0.045502913049921788909870584752660,
1432 0.047982537138836713906392255756915,
1433 0.050277679080715671963325259433440,
1434 0.052362885806407475864366712137873,
1435 0.054251129888545490144543370459876,
1436 0.055950811220412317308240686382747,
1437 0.057437116361567832853582693939506,
1438 0.058689680022394207961974175856788,
1439 0.059720340324174059979099291932562,
1440 0.060539455376045862945360267517565,
1441 0.061128509717053048305859030416293,
1442 0.061471189871425316661544131965264,
1443 0.061580818067832935078759824240066
1450 double *
result,
double *abserr,
1451 double *resabs,
double *resasc)
1456 gsl_integration_qk (26,
xgkE,
wgE,
wgkE, fv1, fv2,
f,
a,
b,
result, abserr, resabs, resasc);
1465 0.999484410050490637571325895705811,
1466 0.996893484074649540271630050918695,
1467 0.991630996870404594858628366109486,
1468 0.983668123279747209970032581605663,
1469 0.973116322501126268374693868423707,
1470 0.960021864968307512216871025581798,
1471 0.944374444748559979415831324037439,
1472 0.926200047429274325879324277080474,
1473 0.905573307699907798546522558925958,
1474 0.882560535792052681543116462530226,
1475 0.857205233546061098958658510658944,
1476 0.829565762382768397442898119732502,
1477 0.799727835821839083013668942322683,
1478 0.767777432104826194917977340974503,
1479 0.733790062453226804726171131369528,
1480 0.697850494793315796932292388026640,
1481 0.660061064126626961370053668149271,
1482 0.620526182989242861140477556431189,
1483 0.579345235826361691756024932172540,
1484 0.536624148142019899264169793311073,
1485 0.492480467861778574993693061207709,
1486 0.447033769538089176780609900322854,
1487 0.400401254830394392535476211542661,
1488 0.352704725530878113471037207089374,
1489 0.304073202273625077372677107199257,
1490 0.254636926167889846439805129817805,
1491 0.204525116682309891438957671002025,
1492 0.153869913608583546963794672743256,
1493 0.102806937966737030147096751318001,
1494 0.051471842555317695833025213166723,
1495 0.000000000000000000000000000000000
1503 0.007968192496166605615465883474674,
1504 0.018466468311090959142302131912047,
1505 0.028784707883323369349719179611292,
1506 0.038799192569627049596801936446348,
1507 0.048402672830594052902938140422808,
1508 0.057493156217619066481721689402056,
1509 0.065974229882180495128128515115962,
1510 0.073755974737705206268243850022191,
1511 0.080755895229420215354694938460530,
1512 0.086899787201082979802387530715126,
1513 0.092122522237786128717632707087619,
1514 0.096368737174644259639468626351810,
1515 0.099593420586795267062780282103569,
1516 0.101762389748405504596428952168554,
1517 0.102852652893558840341285636705415
1522 0.001389013698677007624551591226760,
1523 0.003890461127099884051267201844516,
1524 0.006630703915931292173319826369750,
1525 0.009273279659517763428441146892024,
1526 0.011823015253496341742232898853251,
1527 0.014369729507045804812451432443580,
1528 0.016920889189053272627572289420322,
1529 0.019414141193942381173408951050128,
1530 0.021828035821609192297167485738339,
1531 0.024191162078080601365686370725232,
1532 0.026509954882333101610601709335075,
1533 0.028754048765041292843978785354334,
1534 0.030907257562387762472884252943092,
1535 0.032981447057483726031814191016854,
1536 0.034979338028060024137499670731468,
1537 0.036882364651821229223911065617136,
1538 0.038678945624727592950348651532281,
1539 0.040374538951535959111995279752468,
1540 0.041969810215164246147147541285970,
1541 0.043452539701356069316831728117073,
1542 0.044814800133162663192355551616723,
1543 0.046059238271006988116271735559374,
1544 0.047185546569299153945261478181099,
1545 0.048185861757087129140779492298305,
1546 0.049055434555029778887528165367238,
1547 0.049795683427074206357811569379942,
1548 0.050405921402782346840893085653585,
1549 0.050881795898749606492297473049805,
1550 0.051221547849258772170656282604944,
1551 0.051426128537459025933862879215781,
1552 0.051494729429451567558340433647099
1557 double *
result,
double *abserr,
1558 double *resabs,
double *resasc)
1563 gsl_integration_qk (31,
xgkF,
wgF,
wgkF, fv1, fv2,
f,
a,
b,
result, abserr, resabs, resasc);
1573 GSL_ERROR_VAL (
"workspace length n must be positive integer",
1581 GSL_ERROR_VAL (
"failed to allocate space for workspace struct",
1585 w->alist =
reinterpret_cast<double *
>(
malloc (
n *
sizeof (
double)));
1587 if (
w->alist ==
nullptr)
1595 w->blist =
reinterpret_cast<double *
>(
malloc (
n *
sizeof (
double)));
1597 if (
w->blist ==
nullptr)
1606 w->rlist =
reinterpret_cast<double *
>(
malloc (
n *
sizeof (
double)));
1608 if (
w->rlist ==
nullptr)
1619 w->elist =
reinterpret_cast<double *
>(
malloc (
n *
sizeof (
double)));
1621 if (
w->elist ==
nullptr)
1632 w->order =
reinterpret_cast<size_t *
>(
malloc (
n *
sizeof (
size_t)));
1634 if (
w->order ==
nullptr)
1646 w->level =
reinterpret_cast<size_t *
>(
malloc (
n *
sizeof (
size_t)));
1648 if (
w->level ==
nullptr)
1663 w->maximum_level = 0 ;
1689 workspace->
nrmax = 0;
1690 workspace->
i = workspace->
order[0] ;
1706 int id = workspace->
nrmax;
1709 const size_t * level = workspace->
level;
1710 const size_t * order = workspace->
order;
1712 size_t limit = workspace->
limit ;
1713 size_t last = workspace->
size - 1 ;
1715 if (last > (1 + limit / 2))
1717 jupbnd = limit + 1 - last;
1724 for (k =
id; k <= jupbnd; k++)
1726 size_t i_max = order[workspace->
nrmax];
1728 workspace->
i = i_max ;
1744 size_t i = workspace->
i ;
1745 const size_t * level = workspace->
level;
1808 double *epstab = table->
rlist2;
1809 double *res3la = table->
res3la;
1810 const size_t n = table->
n - 1;
1812 const double current = epstab[
n];
1817 const size_t newelm =
n / 2;
1818 const size_t n_orig =
n;
1822 const size_t nres_orig = table->
nres;
1834 epstab[
n + 2] = epstab[
n];
1837 for (i = 0; i < newelm; i++)
1839 double res = epstab[
n - 2 * i + 2];
1840 double e0 = epstab[
n - 2 * i - 2];
1841 double e1 = epstab[
n - 2 * i - 1];
1844 double e1abs = std::abs(e1);
1845 double delta2 = e2 - e1;
1846 double err2 = std::abs(delta2);
1848 double delta3 = e1 - e0;
1849 double err3 = std::abs(delta3);
1858 if (err2 <= tol2 && err3 <= tol3)
1864 absolute = err2 + err3;
1870 e3 = epstab[
n - 2 * i];
1871 epstab[
n - 2 * i] = e1;
1873 err1 = std::abs(delta1);
1879 if (err1 <= tol1 || err2 <= tol2 || err3 <= tol3)
1885 ss = (1 / delta1 + 1 / delta2) - 1 / delta3;
1891 if (std::abs(ss * e1) <= 0.0001)
1901 epstab[
n - 2 * i] = res;
1904 const double error = err2 + std::abs(res - e2) + err3;
1906 if (error <= *abserr)
1917 const size_t limexp = 50 - 1;
1919 if (n_final == limexp)
1921 n_final = 2 * (limexp / 2);
1925 if (n_orig % 2 == 1)
1927 for (i = 0; i <= newelm; i++)
1929 epstab[1 + i * 2] = epstab[i * 2 + 3];
1934 for (i = 0; i <= newelm; i++)
1936 epstab[i * 2] = epstab[i * 2 + 2];
1940 if (n_orig != n_final)
1942 for (i = 0; i <= n_final; i++)
1944 epstab[i] = epstab[n_orig - n_final + i];
1948 table->
n = n_final + 1;
1952 res3la[nres_orig] = *
result;
1957 *abserr = (std::abs(*
result - res3la[2]) + std::abs(*
result - res3la[1])
1958 + std::abs(*
result - res3la[0]));
1960 res3la[0] = res3la[1];
1961 res3la[1] = res3la[2];
1971 table->
nres = nres_orig + 1;
1995 b,
const double epsabs,
const double epsrel,
const size_t limit,
2002 double epsabs,
double epsrel,
size_t limit,
2004 double *
result,
double * abserr)
2006 int status =
qags (
f,
a,
b, epsabs, epsrel, limit,
2020static double i_transform (
double t,
void *params);
2024 double epsabs,
double epsrel,
size_t limit,
2026 double *
result,
double *abserr)
2035 status =
qags (&f_transform, 0.0, 1.0,
2036 epsabs, epsrel, limit,
2048 double x = (1 - t) / t;
2068 double epsabs,
double epsrel,
size_t limit,
2070 double *
result,
double *abserr)
2077 transform_params.
b =
b ;
2078 transform_params.
f =
f ;
2081 f_transform.
params = &transform_params;
2083 status =
qags (&f_transform, 0.0, 1.0,
2084 epsabs, epsrel, limit,
2098 double x =
b - (1 - t) / t;
2117 double epsabs,
double epsrel,
size_t limit,
2119 double *
result,
double *abserr)
2126 transform_params.
a =
a ;
2127 transform_params.
f =
f ;
2130 f_transform.
params = &transform_params;
2132 status =
qags (&f_transform, 0.0, 1.0,
2133 epsabs, epsrel, limit,
2147 double x =
a + (1 - t) / t;
2156 const double a,
const double b,
2157 const double epsabs,
const double epsrel,
2160 double *
result,
double *abserr,
2174 double error_over_large_intervals = 0;
2179 int roundoff_type1 = 0;
2180 int roundoff_type2 = 0;
2181 int roundoff_type3 = 0;
2183 int error_type2 = 0;
2185 size_t iteration = 0;
2187 int positive_integrand = 0;
2188 int extrapolate = 0;
2189 int disallow_extrapolation = 0;
2200 if (limit > workspace->
limit)
2207 if (epsabs <= 0 && (epsrel < 50 *
GSL_DBL_EPSILON || epsrel < 0.5e-28))
2209 GSL_ERROR (
"tolerance cannot be achieved with given epsabs and epsrel",
2215 q (
f,
a,
b, &result0, &abserr0, &resabs0, &resasc0);
2219 tolerance =
GSL_MAX_DBL (epsabs, epsrel * std::abs(result0));
2221 if (abserr0 <= 100 * GSL_DBL_EPSILON * resabs0 && abserr0 > tolerance)
2226 GSL_ERROR (
"cannot reach tolerance because of roundoff error"
2229 else if ((abserr0 <= tolerance && abserr0 != resasc0) || abserr0 == 0.0)
2236 else if (limit == 1)
2261 size_t current_level;
2284 retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
2286 current_level = workspace->
level[workspace->
i] + 1;
2289 b1 = 0.5 * (a_i + b_i);
2295 q (
f, a1, b1, &area1, &error1, &resabs1, &resasc1);
2296 q (
f, a2, b2, &area2, &error2, &resabs2, &resasc2);
2298 area12 = area1 + area2;
2299 error12 = error1 + error2;
2309 errsum = errsum + error12 - e_i;
2310 area = area + area12 - r_i;
2312 tolerance =
GSL_MAX_DBL (epsabs, epsrel * std::abs(area));
2314 if (resasc1 != error1 && resasc2 != error2)
2316 double delta = r_i - area12;
2318 if (std::abs(delta) <= 1.0e-5 * std::abs(area12) && error12 >= 0.99 * e_i)
2329 if (iteration > 10 && error12 > e_i)
2337 if (roundoff_type1 + roundoff_type2 >= 10 || roundoff_type3 >= 20)
2342 if (roundoff_type2 >= 5)
2357 update (workspace, a1, b1, area1, error1, a2, b2, area2, error2);
2359 if (errsum <= tolerance)
2361 goto compute_result;
2369 if (iteration >= limit - 1)
2377 error_over_large_intervals = errsum;
2383 if (disallow_extrapolation)
2388 error_over_large_intervals += -last_e_i;
2390 if (current_level < workspace->maximum_level)
2392 error_over_large_intervals += error12;
2404 workspace->
nrmax = 1;
2407 if (!error_type2 && error_over_large_intervals > ertest)
2417 qelg (&table, &reseps, &abseps);
2421 if (ktmin > 5 && err_ext < 0.001 * errsum)
2426 if (abseps < err_ext)
2431 correc = error_over_large_intervals;
2432 ertest =
GSL_MAX_DBL (epsabs, epsrel * std::abs(reseps));
2433 if (err_ext <= ertest)
2441 disallow_extrapolation = 1;
2444 if (error_type == 5)
2453 error_over_large_intervals = errsum;
2456 while (iteration < limit);
2462 goto compute_result;
2464 if (error_type || error_type2)
2474 if (res_ext != 0.0 && area != 0.0)
2476 if (err_ext / std::abs(res_ext) > errsum / std::abs(area))
2477 goto compute_result;
2479 else if (err_ext > errsum)
2481 goto compute_result;
2483 else if (area == 0.0)
2492 double max_area =
GSL_MAX_DBL (std::abs(res_ext), std::abs(area));
2494 if (!positive_integrand && max_area < 0.01 * resabs0)
2499 double ratio = res_ext / area;
2501 if (ratio < 0.01 || ratio > 100.0 || errsum > std::abs(area))
2519 if (error_type == 0)
2523 else if (error_type == 1)
2527 else if (error_type == 2)
2529 GSL_ERROR (
"cannot reach tolerance because of roundoff error",
2532 else if (error_type == 3)
2534 GSL_ERROR (
"bad integrand behavior found in the integration interval",
2537 else if (error_type == 4)
2539 GSL_ERROR (
"roundoff error detected in the extrapolation table",
2542 else if (error_type == 5)
2544 GSL_ERROR (
"integral is divergent, or slowly convergent",
static void update(gsl_integration_workspace *workspace, double a1, double b1, double area1, double error1, double a2, double b2, double area2, double error2)
void gsl_integration_qk61(const gsl_function *f, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
gsl_integration_workspace * gsl_integration_workspace_alloc(const size_t n)
static int increase_nrmax(gsl_integration_workspace *workspace)
static double rescale_error(double err, const double result_abs, const double result_asc)
void gsl_integration_qk51(const gsl_function *f, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
static void retrieve(const gsl_integration_workspace *workspace, double *a, double *b, double *r, double *e)
static void set_initial_result(gsl_integration_workspace *workspace, double result, double error)
static const double wgD[11]
static void qelg(struct extrapolation_table *table, double *result, double *abserr)
static double il_transform(double t, void *params)
void gsl_integration_qk(const int n, const double xgk[], const double wg[], const double wgk[], double fv1[], double fv2[], const gsl_function *f, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
static void initialise_table(struct extrapolation_table *table)
static const double wgC[8]
static int test_positivity(double result, double resabs)
static const double wgkE[26]
static const double wgkB[11]
static int large_interval(gsl_integration_workspace *workspace)
static const double wgkF[31]
static void qpsrt(gsl_integration_workspace *workspace)
#define GSL_ERROR_VAL(reason, gsl_errno, value)
double GSL_MAX_DBL(double a, double b)
static const double wgB[5]
double RooAdaptiveGaussKronrodIntegrator1D_GSL_GlueFunction(double x, void *data)
Glue function interacing to GSL code.
static int qag(const gsl_function *f, const double a, const double b, const double epsabs, const double epsrel, const size_t limit, gsl_integration_workspace *workspace, double *result, double *abserr, gsl_integration_rule *q)
static const double xgkE[26]
int gsl_integration_qagil(gsl_function *f, double b, double epsabs, double epsrel, size_t limit, gsl_integration_workspace *workspace, double *result, double *abserr)
static Roo_reg_AGKInteg1D instance
int gsl_integration_qag(const gsl_function *f, double a, double b, double epsabs, double epsrel, size_t limit, int key, gsl_integration_workspace *workspace, double *result, double *abserr)
static void reset_nrmax(gsl_integration_workspace *workspace)
#define GSL_COERCE_DBL(x)
static const double wgkC[16]
static const double xgkA[8]
static const double wgF[15]
void gsl_integration_qk21(const gsl_function *f, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
int gsl_integration_qagi(gsl_function *f, double epsabs, double epsrel, size_t limit, gsl_integration_workspace *workspace, double *result, double *abserr)
void gsl_integration_rule(const gsl_function *f, double a, double b, double *result, double *abserr, double *defabs, double *resabs)
static int qags(const gsl_function *f, const double a, const double b, const double epsabs, const double epsrel, const size_t limit, gsl_integration_workspace *workspace, double *result, double *abserr, gsl_integration_rule *q)
static double sum_results(const gsl_integration_workspace *workspace)
void gsl_integration_qk15(const gsl_function *f, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
static int subinterval_too_small(double a1, double a2, double b2)
int gsl_integration_qagiu(gsl_function *f, double a, double epsabs, double epsrel, size_t limit, gsl_integration_workspace *workspace, double *result, double *abserr)
static const double xgkD[21]
static const double wgA[4]
void gsl_integration_qk41(const gsl_function *f, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
void gsl_integration_qcheb(gsl_function *f, double a, double b, double *cheb12, double *cheb24)
static const double wgkD[21]
void gsl_integration_workspace_free(gsl_integration_workspace *w)
int gsl_integration_qags(const gsl_function *f, double a, double b, double epsabs, double epsrel, size_t limit, gsl_integration_workspace *workspace, double *result, double *abserr)
#define GSL_FN_EVAL(F, x)
static const double xgkF[31]
static const double wgE[13]
static double i_transform(double t, void *params)
double gsl_coerce_double(const double x)
void gsl_integration_qk31(const gsl_function *f, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
static double iu_transform(double t, void *params)
static const double xgkC[16]
static void append_table(struct extrapolation_table *table, double y)
static const double xgkB[11]
static const double wgkA[8]
static void initialise(gsl_integration_workspace *workspace, double a, double b)
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
double getRealValue(const char *name, double defVal=0.0, bool verbose=false) const
Get value of a RooAbsReal stored in set with given name.
Int_t getCatIndex(const char *name, Int_t defVal=0, bool verbose=false) const
Get index value of a RooAbsCategory stored in set with given name.
Abstract interface for evaluating a real-valued function of one real variable and performing numerica...
virtual double getMaxLimit(UInt_t dimension) const =0
virtual double getMinLimit(UInt_t dimension) const =0
UInt_t getDimension() const
Abstract interface for integrators of real-valued functions that implement the RooAbsFunc interface.
bool isValid() const
Is integrator in valid state.
double integrand(const double x[]) const
Return value of integrand at given observable values.
const RooAbsFunc * _function
Pointer to function binding of integrand.
const RooAbsFunc * integrand() const
Return integrand function binding.
bool _valid
Is integrator in valid state?
Implements the Gauss-Kronrod integration algorithm.
double integral(const double *yvec=nullptr) override
Calculate and return integral at at given parameter values.
RooAdaptiveGaussKronrodIntegrator1D(const RooAbsFunc &function, const RooNumIntConfig &config)
Constructor taking a function binding and a configuration object.
bool initialize()
Initialize integrator allocate buffers and setup GSL workspace.
double _epsAbs
Current coordinate.
bool setLimits(double *xmin, double *xmax) override
Change our integration limits.
bool checkLimits() const override
Check that our integration range is finite and otherwise return false.
friend double RooAdaptiveGaussKronrodIntegrator1D_GSL_GlueFunction(double x, void *data)
Glue function interacing to GSL code.
~RooAdaptiveGaussKronrodIntegrator1D() override
Destructor.
double _xmax
Lower integration bound.
static void registerIntegrator(RooNumIntFactory &fact)
Register this class with RooNumIntConfig as a possible choice of numeric integrator for one-dimension...
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Object to represent discrete states.
bool setIndex(Int_t index, bool printError=true) override
Set value by specifying the index code of the desired state.
bool defineType(const std::string &label)
Define a state with given name.
Holds the configuration parameters of the various numeric integrators used by RooRealIntegral.
const RooArgSet & getConfigSection(const char *name) const
Retrieve configuration information specific to integrator with given name.
Factory to instantiate numeric integrators from a given function binding and a given configuration.
static RooNumIntFactory & instance()
Static method returning reference to singleton instance of factory.
bool registerPlugin(std::string const &name, Creator const &creator, const RooArgSet &defConfig, bool canIntegrate1D, bool canIntegrate2D, bool canIntegrateND, bool canIntegrateOpenEnded, const char *depName="")
Method accepting registration of a prototype numeric integrator along with a RooArgSet of its default...
static constexpr int isInfinite(double x)
Return true if x is infinite by RooNumber internal specification.
Variable that can be changed from the outside.
double(* function)(double x, void *params)