72#define GSL_FN_EVAL(F,x) (*((F)->function))(x,(F)->params)
99 double epsabs,
double epsrel,
size_t limit,
102 double *
result,
double *abserr);
107 double epsabs,
double epsrel,
size_t limit,
109 double *
result,
double * abserr) ;
113 double epsabs,
double epsrel,
size_t limit,
115 double *
result,
double *abserr) ;
120 double epsabs,
double epsrel,
size_t limit,
122 double *
result,
double *abserr) ;
127 double epsabs,
double epsrel,
size_t limit,
129 double *
result,
double *abserr) ;
158 RooRealVar maxSeg(
"maxSeg",
"maximum number of segments", 100);
159 RooCategory method(
"method",
"Integration method for each segment");
169 oocoutI(
nullptr,Integration) <<
"RooAdaptiveGaussKronrodIntegrator1D has been registered " << std::endl;
179 _epsAbs(config.epsRel()),
180 _epsRel(config.epsAbs())
200 _epsAbs(config.epsRel()),
201 _epsRel(config.epsAbs()),
260 coutE(Integration) <<
"RooAdaptiveGaussKronrodIntegrator1D::setLimits: cannot override integrand's limits" << endl;
287 if (!infLo && !infHi) {
289 }
else if (infLo && infHi) {
291 }
else if (infLo && !infHi) {
309 return instance->integrand(instance->xvec(
x)) ;
386#define GSL_EBADTOL 13
388#define GSL_ERROR(a,b) oocoutE(nullptr,Integration) << "RooAdaptiveGaussKronrodIntegrator1D::integral() ERROR: " << a << endl ; return b ;
389#define GSL_DBL_MIN 2.2250738585072014e-308
390#define GSL_DBL_MAX 1.7976931348623157e+308
391#define GSL_DBL_EPSILON 2.2204460492503131e-16
394#define GSL_EMAXITER 3
397#define GSL_EDIVERGE 6
400#define GSL_ERROR_VAL(reason, gsl_errno, value) return value ;
402#define GSL_MAX(a,b) ((a) > (b) ? (a) : (b))
418#define GSL_COERCE_DBL(x) (gsl_coerce_double(x))
431 double *
result,
double *abserr,
432 double *defabs,
double *resabs);
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 *
result,
double *abserr,
456 double *resabs,
double *resasc);
459 double *cheb12,
double *cheb24);
477 const double wg[],
const double wgk[],
478 double fv1[],
double fv2[],
480 double *
result,
double * abserr,
481 double * resabs,
double * resasc);
486 double epsabs,
double epsrel,
size_t limit,
489 double *
result,
double *abserr);
500 workspace->
nrmax = 0;
504 workspace->
rlist[0] = 0.0;
505 workspace->
elist[0] = 0.0;
506 workspace->
order[0] = 0;
507 workspace->
level[0] = 0;
515 double result,
double error);
519 double result,
double error)
523 workspace->
elist[0] = error;
533 const size_t last = workspace->
size - 1;
534 const size_t limit = workspace->
limit;
536 double * elist = workspace->
elist;
537 size_t * order = workspace->
order;
543 size_t i_nrmax = workspace->
nrmax;
544 size_t i_maxerr = order[i_nrmax] ;
552 workspace->
i = i_maxerr ;
556 errmax = elist[i_maxerr] ;
563 while (i_nrmax > 0 && errmax > elist[order[i_nrmax - 1]])
565 order[i_nrmax] = order[i_nrmax - 1] ;
573 if(last < (limit/2 + 2))
579 top = limit - last + 1;
590 while (i < top && errmax < elist[order[i]])
592 order[i-1] = order[i] ;
596 order[i-1] = i_maxerr ;
600 errmin = elist[last] ;
604 while (k > i - 2 && errmin >= elist[order[k]])
606 order[k+1] = order[k] ;
614 i_maxerr = order[i_nrmax] ;
616 workspace->
i = i_maxerr ;
617 workspace->
nrmax = i_nrmax ;
625 double a1,
double b1,
double area1,
double error1,
626 double a2,
double b2,
double area2,
double error2);
630 double *
a,
double *
b,
double *
r,
double *
e);
636 double a1,
double b1,
double area1,
double error1,
637 double a2,
double b2,
double area2,
double error2)
639 double * alist = workspace->
alist ;
640 double * blist = workspace->
blist ;
641 double * rlist = workspace->
rlist ;
642 double * elist = workspace->
elist ;
643 size_t * level = workspace->
level ;
645 const size_t i_max = workspace->
i ;
646 const size_t i_new = workspace->
size ;
648 const size_t new_level = workspace->
level[i_max] + 1;
655 rlist[i_max] = area2;
656 elist[i_max] = error2;
657 level[i_max] = new_level;
661 rlist[i_new] = area1;
662 elist[i_new] = error1;
663 level[i_new] = new_level;
668 rlist[i_max] = area1;
669 elist[i_max] = error1;
670 level[i_max] = new_level;
674 rlist[i_new] = area2;
675 elist[i_new] = error2;
676 level[i_new] = new_level;
691 double *
a,
double *
b,
double *
r,
double *
e)
693 const size_t i = workspace->
i;
694 double * alist = workspace->
alist;
695 double * blist = workspace->
blist;
696 double * rlist = workspace->
rlist;
697 double * elist = workspace->
elist;
711 const double *
const rlist = workspace->
rlist ;
712 const size_t n = workspace->
size;
715 double result_sum = 0;
717 for (k = 0; k <
n; k++)
719 result_sum += rlist[k];
734 double tmp = (1 + 100 *
e) * (std::abs(a2) + 1000 * u);
736 int status = std::abs(a1) <= tmp && std::abs(b2) <= tmp;
744 const double a,
const double b,
745 const double epsabs,
const double epsrel,
748 double *
result,
double * abserr,
754 double epsabs,
double epsrel,
size_t limit,
757 double *
result,
double * abserr)
793 status =
qag (
f,
a,
b, epsabs, epsrel, limit,
803 const double a,
const double b,
804 const double epsabs,
const double epsrel,
807 double *
result,
double *abserr,
811 double result0, abserr0, resabs0, resasc0;
813 size_t iteration = 0;
814 int roundoff_type1 = 0, roundoff_type2 = 0, error_type = 0;
825 if (limit > workspace->
limit)
830 if (epsabs <= 0 && (epsrel < 50 *
GSL_DBL_EPSILON || epsrel < 0.5e-28))
832 GSL_ERROR (
"tolerance cannot be acheived with given epsabs and epsrel",
838 q (
f,
a,
b, &result0, &abserr0, &resabs0, &resasc0);
844 tolerance =
GSL_MAX_DBL (epsabs, epsrel * std::abs(result0));
850 if (abserr0 <= round_off && abserr0 > tolerance)
855 GSL_ERROR (
"cannot reach tolerance because of roundoff error "
858 else if ((abserr0 <= tolerance && abserr0 != resasc0) || abserr0 == 0.0)
880 double a1, b1, a2, b2;
881 double a_i, b_i, r_i, e_i;
882 double area1 = 0, area2 = 0, area12 = 0;
883 double error1 = 0, error2 = 0, error12 = 0;
884 double resasc1, resasc2;
885 double resabs1, resabs2;
889 retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
892 b1 = 0.5 * (a_i + b_i);
896 q (
f, a1, b1, &area1, &error1, &resabs1, &resasc1);
897 q (
f, a2, b2, &area2, &error2, &resabs2, &resasc2);
899 area12 = area1 + area2;
900 error12 = error1 + error2;
902 errsum += (error12 - e_i);
903 area += area12 - r_i;
905 if (resasc1 != error1 && resasc2 != error2)
907 double delta = r_i - area12;
909 if (std::abs(delta) <= 1.0e-5 * std::abs(area12) && error12 >= 0.99 * e_i)
913 if (iteration >= 10 && error12 > e_i)
919 tolerance =
GSL_MAX_DBL (epsabs, epsrel * std::abs(area));
921 if (errsum > tolerance)
923 if (roundoff_type1 >= 6 || roundoff_type2 >= 20)
937 update (workspace, a1, b1, area1, error1, a2, b2, area2, error2);
939 retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
944 while (iteration < limit && !error_type && errsum > tolerance);
949 if (errsum <= tolerance)
953 else if (error_type == 2)
955 GSL_ERROR (
"roundoff error prevents tolerance from being achieved",
958 else if (error_type == 3)
960 GSL_ERROR (
"bad integrand behavior found in the integration interval",
963 else if (iteration == limit)
973static double rescale_error (
double err,
const double result_abs,
const double result_asc) ;
978 err = std::abs(err) ;
980 if (result_asc != 0 && err != 0)
982 double scale =
TMath::Power((200 * err / result_asc), 1.5) ;
986 err = result_asc * scale ;
1009 const double xgk[],
const double wg[],
const double wgk[],
1010 double fv1[],
double fv2[],
1012 double *
result,
double *abserr,
1013 double *resabs,
double *resasc)
1016 const double center = 0.5 * (
a +
b);
1017 const double half_length = 0.5 * (
b -
a);
1018 const double abs_half_length = std::abs(half_length);
1021 double result_gauss = 0;
1022 double result_kronrod = f_center * wgk[
n - 1];
1024 double result_abs = std::abs(result_kronrod);
1025 double result_asc = 0;
1026 double mean = 0, err = 0;
1032 result_gauss = f_center * wg[
n / 2 - 1];
1035 for (j = 0; j < (
n - 1) / 2; j++)
1037 const int jtw = j * 2 + 1;
1038 const double abscissa = half_length * xgk[jtw];
1039 const double fval1 =
GSL_FN_EVAL (
f, center - abscissa);
1040 const double fval2 =
GSL_FN_EVAL (
f, center + abscissa);
1041 const double fsum = fval1 + fval2;
1044 result_gauss += wg[j] * fsum;
1045 result_kronrod += wgk[jtw] * fsum;
1046 result_abs += wgk[jtw] * (std::abs(fval1) + std::abs(fval2));
1049 for (j = 0; j <
n / 2; j++)
1052 const double abscissa = half_length * xgk[jtwm1];
1053 const double fval1 =
GSL_FN_EVAL (
f, center - abscissa);
1054 const double fval2 =
GSL_FN_EVAL (
f, center + abscissa);
1057 result_kronrod += wgk[jtwm1] * (fval1 + fval2);
1058 result_abs += wgk[jtwm1] * (std::abs(fval1) + std::abs(fval2));
1061 mean = result_kronrod * 0.5;
1063 result_asc = wgk[
n - 1] * std::abs(f_center - mean);
1065 for (j = 0; j <
n - 1; j++)
1067 result_asc += wgk[j] * (std::abs(fv1[j] - mean) + std::abs(fv2[j] - mean));
1072 err = (result_kronrod - result_gauss) * half_length;
1074 result_kronrod *= half_length;
1075 result_abs *= abs_half_length;
1076 result_asc *= abs_half_length;
1078 *
result = result_kronrod;
1079 *resabs = result_abs;
1080 *resasc = result_asc;
1091 0.991455371120812639206854697526329,
1092 0.949107912342758524526189684047851,
1093 0.864864423359769072789712788640926,
1094 0.741531185599394439863864773280788,
1095 0.586087235467691130294144838258730,
1096 0.405845151377397166906606412076961,
1097 0.207784955007898467600689403773245,
1098 0.000000000000000000000000000000000
1106 0.129484966168869693270611432679082,
1107 0.279705391489276667901467771423780,
1108 0.381830050505118944950369775488975,
1109 0.417959183673469387755102040816327
1114 0.022935322010529224963732008058970,
1115 0.063092092629978553290700663189204,
1116 0.104790010322250183839876322541518,
1117 0.140653259715525918745189590510238,
1118 0.169004726639267902826583426598550,
1119 0.190350578064785409913256402421014,
1120 0.204432940075298892414161999234649,
1121 0.209482141084727828012999174891714
1126 double *
result,
double *abserr,
1127 double *resabs,
double *resasc)
1129 double fv1[8], fv2[8];
1131 gsl_integration_qk (8,
xgkA,
wgA,
wgkA, fv1, fv2,
f,
a,
b,
result, abserr, resabs, resasc);
1140 0.995657163025808080735527280689003,
1141 0.973906528517171720077964012084452,
1142 0.930157491355708226001207180059508,
1143 0.865063366688984510732096688423493,
1144 0.780817726586416897063717578345042,
1145 0.679409568299024406234327365114874,
1146 0.562757134668604683339000099272694,
1147 0.433395394129247190799265943165784,
1148 0.294392862701460198131126603103866,
1149 0.148874338981631210884826001129720,
1150 0.000000000000000000000000000000000
1158 0.066671344308688137593568809893332,
1159 0.149451349150580593145776339657697,
1160 0.219086362515982043995534934228163,
1161 0.269266719309996355091226921569469,
1162 0.295524224714752870173892994651338
1167 0.011694638867371874278064396062192,
1168 0.032558162307964727478818972459390,
1169 0.054755896574351996031381300244580,
1170 0.075039674810919952767043140916190,
1171 0.093125454583697605535065465083366,
1172 0.109387158802297641899210590325805,
1173 0.123491976262065851077958109831074,
1174 0.134709217311473325928054001771707,
1175 0.142775938577060080797094273138717,
1176 0.147739104901338491374841515972068,
1177 0.149445554002916905664936468389821
1183 double *
result,
double *abserr,
1184 double *resabs,
double *resasc)
1186 double fv1[11], fv2[11];
1188 gsl_integration_qk (11,
xgkB,
wgB,
wgkB, fv1, fv2,
f,
a,
b,
result, abserr, resabs, resasc);
1197 0.998002298693397060285172840152271,
1198 0.987992518020485428489565718586613,
1199 0.967739075679139134257347978784337,
1200 0.937273392400705904307758947710209,
1201 0.897264532344081900882509656454496,
1202 0.848206583410427216200648320774217,
1203 0.790418501442465932967649294817947,
1204 0.724417731360170047416186054613938,
1205 0.650996741297416970533735895313275,
1206 0.570972172608538847537226737253911,
1207 0.485081863640239680693655740232351,
1208 0.394151347077563369897207370981045,
1209 0.299180007153168812166780024266389,
1210 0.201194093997434522300628303394596,
1211 0.101142066918717499027074231447392,
1212 0.000000000000000000000000000000000
1220 0.030753241996117268354628393577204,
1221 0.070366047488108124709267416450667,
1222 0.107159220467171935011869546685869,
1223 0.139570677926154314447804794511028,
1224 0.166269205816993933553200860481209,
1225 0.186161000015562211026800561866423,
1226 0.198431485327111576456118326443839,
1227 0.202578241925561272880620199967519
1232 0.005377479872923348987792051430128,
1233 0.015007947329316122538374763075807,
1234 0.025460847326715320186874001019653,
1235 0.035346360791375846222037948478360,
1236 0.044589751324764876608227299373280,
1237 0.053481524690928087265343147239430,
1238 0.062009567800670640285139230960803,
1239 0.069854121318728258709520077099147,
1240 0.076849680757720378894432777482659,
1241 0.083080502823133021038289247286104,
1242 0.088564443056211770647275443693774,
1243 0.093126598170825321225486872747346,
1244 0.096642726983623678505179907627589,
1245 0.099173598721791959332393173484603,
1246 0.100769845523875595044946662617570,
1247 0.101330007014791549017374792767493
1252 double *
result,
double *abserr,
1253 double *resabs,
double *resasc)
1255 double fv1[16], fv2[16];
1257 gsl_integration_qk (16,
xgkC,
wgC,
wgkC, fv1, fv2,
f,
a,
b,
result, abserr, resabs, resasc);
1266 0.998859031588277663838315576545863,
1267 0.993128599185094924786122388471320,
1268 0.981507877450250259193342994720217,
1269 0.963971927277913791267666131197277,
1270 0.940822633831754753519982722212443,
1271 0.912234428251325905867752441203298,
1272 0.878276811252281976077442995113078,
1273 0.839116971822218823394529061701521,
1274 0.795041428837551198350638833272788,
1275 0.746331906460150792614305070355642,
1276 0.693237656334751384805490711845932,
1277 0.636053680726515025452836696226286,
1278 0.575140446819710315342946036586425,
1279 0.510867001950827098004364050955251,
1280 0.443593175238725103199992213492640,
1281 0.373706088715419560672548177024927,
1282 0.301627868114913004320555356858592,
1283 0.227785851141645078080496195368575,
1284 0.152605465240922675505220241022678,
1285 0.076526521133497333754640409398838,
1286 0.000000000000000000000000000000000
1294 0.017614007139152118311861962351853,
1295 0.040601429800386941331039952274932,
1296 0.062672048334109063569506535187042,
1297 0.083276741576704748724758143222046,
1298 0.101930119817240435036750135480350,
1299 0.118194531961518417312377377711382,
1300 0.131688638449176626898494499748163,
1301 0.142096109318382051329298325067165,
1302 0.149172986472603746787828737001969,
1303 0.152753387130725850698084331955098
1308 0.003073583718520531501218293246031,
1309 0.008600269855642942198661787950102,
1310 0.014626169256971252983787960308868,
1311 0.020388373461266523598010231432755,
1312 0.025882133604951158834505067096153,
1313 0.031287306777032798958543119323801,
1314 0.036600169758200798030557240707211,
1315 0.041668873327973686263788305936895,
1316 0.046434821867497674720231880926108,
1317 0.050944573923728691932707670050345,
1318 0.055195105348285994744832372419777,
1319 0.059111400880639572374967220648594,
1320 0.062653237554781168025870122174255,
1321 0.065834597133618422111563556969398,
1322 0.068648672928521619345623411885368,
1323 0.071054423553444068305790361723210,
1324 0.073030690332786667495189417658913,
1325 0.074582875400499188986581418362488,
1326 0.075704497684556674659542775376617,
1327 0.076377867672080736705502835038061,
1328 0.076600711917999656445049901530102
1333 double *
result,
double *abserr,
1334 double *resabs,
double *resasc)
1336 double fv1[21], fv2[21];
1338 gsl_integration_qk (21,
xgkD,
wgD,
wgkD, fv1, fv2,
f,
a,
b,
result, abserr, resabs, resasc);
1347 0.999262104992609834193457486540341,
1348 0.995556969790498097908784946893902,
1349 0.988035794534077247637331014577406,
1350 0.976663921459517511498315386479594,
1351 0.961614986425842512418130033660167,
1352 0.942974571228974339414011169658471,
1353 0.920747115281701561746346084546331,
1354 0.894991997878275368851042006782805,
1355 0.865847065293275595448996969588340,
1356 0.833442628760834001421021108693570,
1357 0.797873797998500059410410904994307,
1358 0.759259263037357630577282865204361,
1359 0.717766406813084388186654079773298,
1360 0.673566368473468364485120633247622,
1361 0.626810099010317412788122681624518,
1362 0.577662930241222967723689841612654,
1363 0.526325284334719182599623778158010,
1364 0.473002731445714960522182115009192,
1365 0.417885382193037748851814394594572,
1366 0.361172305809387837735821730127641,
1367 0.303089538931107830167478909980339,
1368 0.243866883720988432045190362797452,
1369 0.183718939421048892015969888759528,
1370 0.122864692610710396387359818808037,
1371 0.061544483005685078886546392366797,
1372 0.000000000000000000000000000000000
1380 0.011393798501026287947902964113235,
1381 0.026354986615032137261901815295299,
1382 0.040939156701306312655623487711646,
1383 0.054904695975835191925936891540473,
1384 0.068038333812356917207187185656708,
1385 0.080140700335001018013234959669111,
1386 0.091028261982963649811497220702892,
1387 0.100535949067050644202206890392686,
1388 0.108519624474263653116093957050117,
1389 0.114858259145711648339325545869556,
1390 0.119455763535784772228178126512901,
1391 0.122242442990310041688959518945852,
1392 0.123176053726715451203902873079050
1397 0.001987383892330315926507851882843,
1398 0.005561932135356713758040236901066,
1399 0.009473973386174151607207710523655,
1400 0.013236229195571674813656405846976,
1401 0.016847817709128298231516667536336,
1402 0.020435371145882835456568292235939,
1403 0.024009945606953216220092489164881,
1404 0.027475317587851737802948455517811,
1405 0.030792300167387488891109020215229,
1406 0.034002130274329337836748795229551,
1407 0.037116271483415543560330625367620,
1408 0.040083825504032382074839284467076,
1409 0.042872845020170049476895792439495,
1410 0.045502913049921788909870584752660,
1411 0.047982537138836713906392255756915,
1412 0.050277679080715671963325259433440,
1413 0.052362885806407475864366712137873,
1414 0.054251129888545490144543370459876,
1415 0.055950811220412317308240686382747,
1416 0.057437116361567832853582693939506,
1417 0.058689680022394207961974175856788,
1418 0.059720340324174059979099291932562,
1419 0.060539455376045862945360267517565,
1420 0.061128509717053048305859030416293,
1421 0.061471189871425316661544131965264,
1422 0.061580818067832935078759824240066
1429 double *
result,
double *abserr,
1430 double *resabs,
double *resasc)
1432 double fv1[26], fv2[26];
1434 gsl_integration_qk (26,
xgkE,
wgE,
wgkE, fv1, fv2,
f,
a,
b,
result, abserr, resabs, resasc);
1443 0.999484410050490637571325895705811,
1444 0.996893484074649540271630050918695,
1445 0.991630996870404594858628366109486,
1446 0.983668123279747209970032581605663,
1447 0.973116322501126268374693868423707,
1448 0.960021864968307512216871025581798,
1449 0.944374444748559979415831324037439,
1450 0.926200047429274325879324277080474,
1451 0.905573307699907798546522558925958,
1452 0.882560535792052681543116462530226,
1453 0.857205233546061098958658510658944,
1454 0.829565762382768397442898119732502,
1455 0.799727835821839083013668942322683,
1456 0.767777432104826194917977340974503,
1457 0.733790062453226804726171131369528,
1458 0.697850494793315796932292388026640,
1459 0.660061064126626961370053668149271,
1460 0.620526182989242861140477556431189,
1461 0.579345235826361691756024932172540,
1462 0.536624148142019899264169793311073,
1463 0.492480467861778574993693061207709,
1464 0.447033769538089176780609900322854,
1465 0.400401254830394392535476211542661,
1466 0.352704725530878113471037207089374,
1467 0.304073202273625077372677107199257,
1468 0.254636926167889846439805129817805,
1469 0.204525116682309891438957671002025,
1470 0.153869913608583546963794672743256,
1471 0.102806937966737030147096751318001,
1472 0.051471842555317695833025213166723,
1473 0.000000000000000000000000000000000
1481 0.007968192496166605615465883474674,
1482 0.018466468311090959142302131912047,
1483 0.028784707883323369349719179611292,
1484 0.038799192569627049596801936446348,
1485 0.048402672830594052902938140422808,
1486 0.057493156217619066481721689402056,
1487 0.065974229882180495128128515115962,
1488 0.073755974737705206268243850022191,
1489 0.080755895229420215354694938460530,
1490 0.086899787201082979802387530715126,
1491 0.092122522237786128717632707087619,
1492 0.096368737174644259639468626351810,
1493 0.099593420586795267062780282103569,
1494 0.101762389748405504596428952168554,
1495 0.102852652893558840341285636705415
1500 0.001389013698677007624551591226760,
1501 0.003890461127099884051267201844516,
1502 0.006630703915931292173319826369750,
1503 0.009273279659517763428441146892024,
1504 0.011823015253496341742232898853251,
1505 0.014369729507045804812451432443580,
1506 0.016920889189053272627572289420322,
1507 0.019414141193942381173408951050128,
1508 0.021828035821609192297167485738339,
1509 0.024191162078080601365686370725232,
1510 0.026509954882333101610601709335075,
1511 0.028754048765041292843978785354334,
1512 0.030907257562387762472884252943092,
1513 0.032981447057483726031814191016854,
1514 0.034979338028060024137499670731468,
1515 0.036882364651821229223911065617136,
1516 0.038678945624727592950348651532281,
1517 0.040374538951535959111995279752468,
1518 0.041969810215164246147147541285970,
1519 0.043452539701356069316831728117073,
1520 0.044814800133162663192355551616723,
1521 0.046059238271006988116271735559374,
1522 0.047185546569299153945261478181099,
1523 0.048185861757087129140779492298305,
1524 0.049055434555029778887528165367238,
1525 0.049795683427074206357811569379942,
1526 0.050405921402782346840893085653585,
1527 0.050881795898749606492297473049805,
1528 0.051221547849258772170656282604944,
1529 0.051426128537459025933862879215781,
1530 0.051494729429451567558340433647099
1535 double *
result,
double *abserr,
1536 double *resabs,
double *resasc)
1538 double fv1[31], fv2[31];
1540 gsl_integration_qk (31,
xgkF,
wgF,
wgkF, fv1, fv2,
f,
a,
b,
result, abserr, resabs, resasc);
1550 GSL_ERROR_VAL (
"workspace length n must be positive integer",
1559 GSL_ERROR_VAL (
"failed to allocate space for workspace struct",
1563 w->alist = (
double *)
malloc (
n *
sizeof (
double));
1573 w->blist = (
double *)
malloc (
n *
sizeof (
double));
1584 w->rlist = (
double *)
malloc (
n *
sizeof (
double));
1597 w->elist = (
double *)
malloc (
n *
sizeof (
double));
1610 w->order = (
size_t *)
malloc (
n *
sizeof (
size_t));
1624 w->level = (
size_t *)
malloc (
n *
sizeof (
size_t));
1641 w->maximum_level = 0 ;
1667 workspace->
nrmax = 0;
1668 workspace->
i = workspace->
order[0] ;
1684 int id = workspace->
nrmax;
1687 const size_t * level = workspace->
level;
1688 const size_t * order = workspace->
order;
1690 size_t limit = workspace->
limit ;
1691 size_t last = workspace->
size - 1 ;
1693 if (last > (1 + limit / 2))
1695 jupbnd = limit + 1 - last;
1702 for (k =
id; k <= jupbnd; k++)
1704 size_t i_max = order[workspace->
nrmax];
1706 workspace->
i = i_max ;
1722 size_t i = workspace->
i ;
1723 const size_t * level = workspace->
level;
1786 double *epstab = table->
rlist2;
1787 double *res3la = table->
res3la;
1788 const size_t n = table->
n - 1;
1790 const double current = epstab[
n];
1795 const size_t newelm =
n / 2;
1796 const size_t n_orig =
n;
1800 const size_t nres_orig = table->
nres;
1812 epstab[
n + 2] = epstab[
n];
1815 for (i = 0; i < newelm; i++)
1817 double res = epstab[
n - 2 * i + 2];
1818 double e0 = epstab[
n - 2 * i - 2];
1819 double e1 = epstab[
n - 2 * i - 1];
1822 double e1abs = std::abs(e1);
1823 double delta2 = e2 - e1;
1824 double err2 = std::abs(delta2);
1826 double delta3 = e1 - e0;
1827 double err3 = std::abs(delta3);
1830 double e3, delta1, err1, tol1, ss;
1832 if (err2 <= tol2 && err3 <= tol3)
1838 absolute = err2 + err3;
1844 e3 = epstab[
n - 2 * i];
1845 epstab[
n - 2 * i] = e1;
1847 err1 = std::abs(delta1);
1853 if (err1 <= tol1 || err2 <= tol2 || err3 <= tol3)
1859 ss = (1 / delta1 + 1 / delta2) - 1 / delta3;
1865 if (std::abs(ss * e1) <= 0.0001)
1875 epstab[
n - 2 * i] = res;
1878 const double error = err2 + std::abs(res - e2) + err3;
1880 if (error <= *abserr)
1891 const size_t limexp = 50 - 1;
1893 if (n_final == limexp)
1895 n_final = 2 * (limexp / 2);
1899 if (n_orig % 2 == 1)
1901 for (i = 0; i <= newelm; i++)
1903 epstab[1 + i * 2] = epstab[i * 2 + 3];
1908 for (i = 0; i <= newelm; i++)
1910 epstab[i * 2] = epstab[i * 2 + 2];
1914 if (n_orig != n_final)
1916 for (i = 0; i <= n_final; i++)
1918 epstab[i] = epstab[n_orig - n_final + i];
1922 table->
n = n_final + 1;
1926 res3la[nres_orig] = *
result;
1931 *abserr = (std::abs(*
result - res3la[2]) + std::abs(*
result - res3la[1])
1932 + std::abs(*
result - res3la[0]));
1934 res3la[0] = res3la[1];
1935 res3la[1] = res3la[2];
1945 table->
nres = nres_orig + 1;
1969 b,
const double epsabs,
const double epsrel,
const size_t limit,
1976 double epsabs,
double epsrel,
size_t limit,
1978 double *
result,
double * abserr)
1980 int status =
qags (
f,
a,
b, epsabs, epsrel, limit,
1994static double i_transform (
double t,
void *params);
1998 double epsabs,
double epsrel,
size_t limit,
2000 double *
result,
double *abserr)
2009 status =
qags (&f_transform, 0.0, 1.0,
2010 epsabs, epsrel, limit,
2022 double x = (1 - t) / t;
2042 double epsabs,
double epsrel,
size_t limit,
2044 double *
result,
double *abserr)
2051 transform_params.
b =
b ;
2052 transform_params.
f =
f ;
2055 f_transform.
params = &transform_params;
2057 status =
qags (&f_transform, 0.0, 1.0,
2058 epsabs, epsrel, limit,
2072 double x =
b - (1 - t) / t;
2091 double epsabs,
double epsrel,
size_t limit,
2093 double *
result,
double *abserr)
2100 transform_params.
a =
a ;
2101 transform_params.
f =
f ;
2104 f_transform.
params = &transform_params;
2106 status =
qags (&f_transform, 0.0, 1.0,
2107 epsabs, epsrel, limit,
2121 double x =
a + (1 - t) / t;
2130 const double a,
const double b,
2131 const double epsabs,
const double epsrel,
2134 double *
result,
double *abserr,
2137 double area, errsum;
2138 double res_ext, err_ext;
2139 double result0, abserr0, resabs0, resasc0;
2143 double error_over_large_intervals = 0;
2144 double reseps = 0, abseps = 0, correc = 0;
2146 int roundoff_type1 = 0, roundoff_type2 = 0, roundoff_type3 = 0;
2147 int error_type = 0, error_type2 = 0;
2149 size_t iteration = 0;
2151 int positive_integrand = 0;
2152 int extrapolate = 0;
2153 int disallow_extrapolation = 0;
2164 if (limit > workspace->
limit)
2171 if (epsabs <= 0 && (epsrel < 50 *
GSL_DBL_EPSILON || epsrel < 0.5e-28))
2173 GSL_ERROR (
"tolerance cannot be acheived with given epsabs and epsrel",
2179 q (
f,
a,
b, &result0, &abserr0, &resabs0, &resasc0);
2183 tolerance =
GSL_MAX_DBL (epsabs, epsrel * std::abs(result0));
2185 if (abserr0 <= 100 * GSL_DBL_EPSILON * resabs0 && abserr0 > tolerance)
2190 GSL_ERROR (
"cannot reach tolerance because of roundoff error"
2193 else if ((abserr0 <= tolerance && abserr0 != resasc0) || abserr0 == 0.0)
2200 else if (limit == 1)
2225 size_t current_level;
2226 double a1, b1, a2, b2;
2227 double a_i, b_i, r_i, e_i;
2228 double area1 = 0, area2 = 0, area12 = 0;
2229 double error1 = 0, error2 = 0, error12 = 0;
2230 double resasc1, resasc2;
2231 double resabs1, resabs2;
2236 retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
2238 current_level = workspace->
level[workspace->
i] + 1;
2241 b1 = 0.5 * (a_i + b_i);
2247 q (
f, a1, b1, &area1, &error1, &resabs1, &resasc1);
2248 q (
f, a2, b2, &area2, &error2, &resabs2, &resasc2);
2250 area12 = area1 + area2;
2251 error12 = error1 + error2;
2261 errsum = errsum + error12 - e_i;
2262 area = area + area12 - r_i;
2264 tolerance =
GSL_MAX_DBL (epsabs, epsrel * std::abs(area));
2266 if (resasc1 != error1 && resasc2 != error2)
2268 double delta = r_i - area12;
2270 if (std::abs(delta) <= 1.0e-5 * std::abs(area12) && error12 >= 0.99 * e_i)
2281 if (iteration > 10 && error12 > e_i)
2289 if (roundoff_type1 + roundoff_type2 >= 10 || roundoff_type3 >= 20)
2294 if (roundoff_type2 >= 5)
2309 update (workspace, a1, b1, area1, error1, a2, b2, area2, error2);
2311 if (errsum <= tolerance)
2313 goto compute_result;
2321 if (iteration >= limit - 1)
2329 error_over_large_intervals = errsum;
2335 if (disallow_extrapolation)
2340 error_over_large_intervals += -last_e_i;
2342 if (current_level < workspace->maximum_level)
2344 error_over_large_intervals += error12;
2356 workspace->
nrmax = 1;
2359 if (!error_type2 && error_over_large_intervals > ertest)
2369 qelg (&table, &reseps, &abseps);
2373 if (ktmin > 5 && err_ext < 0.001 * errsum)
2378 if (abseps < err_ext)
2383 correc = error_over_large_intervals;
2384 ertest =
GSL_MAX_DBL (epsabs, epsrel * std::abs(reseps));
2385 if (err_ext <= ertest)
2393 disallow_extrapolation = 1;
2396 if (error_type == 5)
2405 error_over_large_intervals = errsum;
2408 while (iteration < limit);
2414 goto compute_result;
2416 if (error_type || error_type2)
2426 if (res_ext != 0.0 && area != 0.0)
2428 if (err_ext / std::abs(res_ext) > errsum / std::abs(area))
2429 goto compute_result;
2431 else if (err_ext > errsum)
2433 goto compute_result;
2435 else if (area == 0.0)
2444 double max_area =
GSL_MAX_DBL (std::abs(res_ext), std::abs(area));
2446 if (!positive_integrand && max_area < 0.01 * resabs0)
2451 double ratio = res_ext / area;
2453 if (ratio < 0.01 || ratio > 100.0 || errsum > std::abs(area))
2471 if (error_type == 0)
2475 else if (error_type == 1)
2479 else if (error_type == 2)
2481 GSL_ERROR (
"cannot reach tolerance because of roundoff error",
2484 else if (error_type == 3)
2486 GSL_ERROR (
"bad integrand behavior found in the integration interval",
2489 else if (error_type == 4)
2491 GSL_ERROR (
"roundoff error detected in the extrapolation table",
2494 else if (error_type == 5)
2496 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)
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
RooAbsIntegrator is the abstract interface for integrators of real-valued functions that implement th...
bool isValid() const
Is integrator in valid state.
const RooAbsFunc * _function
Pointer to function binding of integrand.
const RooAbsFunc * integrand() const
Return integrand function binding.
bool _valid
Is integrator in valid state?
RooAdaptiveGaussKronrodIntegrator1D implements the Gauss-Kronrod integration algorithm.
double integral(const double *yvec=nullptr) override
Calculate and return integral at at given parameter values.
RooAbsIntegrator * clone(const RooAbsFunc &function, const RooNumIntConfig &config) const override
Virtual constructor.
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.
RooAdaptiveGaussKronrodIntegrator1D()
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.
RooCategory is an 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.
RooNumIntConfig holds the configuration parameters of the various numeric integrators used by RooReal...
const RooArgSet & getConfigSection(const char *name) const
Retrieve configuration information specific to integrator with given name.
RooNumIntFactory is a factory to instantiate numeric integrators from a given function binding and a ...
bool storeProtoIntegrator(RooAbsIntegrator *proto, const RooArgSet &defConfig, const char *depName="")
Method accepting registration of a prototype numeric integrator along with a RooArgSet of its default...
static RooNumIntFactory & instance()
Static method returning reference to singleton instance of factory.
static constexpr int isInfinite(double x)
Return true if x is infinite by RooNumber internal specification.
RooRealVar represents a variable that can be changed from the outside.
virtual const char * ClassName() const
Returns name of class to which the object belongs.
static Roo_reg_AGKInteg1D instance
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
static void registerIntegrator()
double(* function)(double x, void *params)