68struct gsl_function_struct
70 double (* function) (
double x,
void * params);
74#define GSL_FN_EVAL(F,x) (*((F)->function))(x,(F)->params)
91gsl_integration_workspace;
93gsl_integration_workspace *
101 double epsabs,
double epsrel,
size_t limit,
103 gsl_integration_workspace * workspace,
104 double *result,
double *abserr);
109 double epsabs,
double epsrel,
size_t limit,
110 gsl_integration_workspace * workspace,
111 double * result,
double * abserr) ;
115 double epsabs,
double epsrel,
size_t limit,
116 gsl_integration_workspace * workspace,
117 double *result,
double *abserr) ;
122 double epsabs,
double epsrel,
size_t limit,
123 gsl_integration_workspace * workspace,
124 double *result,
double *abserr) ;
129 double epsabs,
double epsrel,
size_t limit,
130 gsl_integration_workspace * workspace,
131 double *result,
double *abserr) ;
142 static void registerIntegrator()
149struct Roo_reg_AGKInteg1D {
150 Roo_reg_AGKInteg1D() { Roo_internal_AGKInteg1D::registerIntegrator(); }
160 RooRealVar maxSeg(
"maxSeg",
"maximum number of segments", 100);
161 RooCategory method(
"method",
"Integration method for each segment");
171 oocoutI((
TObject*)
nullptr,Integration) <<
"RooAdaptiveGaussKronrodIntegrator1D has been registered " << std::endl;
191 _epsAbs(config.epsRel()),
192 _epsRel(config.epsAbs()),
213 _epsAbs(config.epsRel()),
214 _epsRel(config.epsAbs()),
277 coutE(Integration) <<
"RooAdaptiveGaussKronrodIntegrator1D::setLimits: cannot override integrand's limits" << endl;
304 if (!infLo && !infHi) {
306 }
else if (infLo && infHi) {
308 }
else if (infLo && !infHi) {
326 return instance->integrand(instance->xvec(
x)) ;
351 double result, error;
403#define GSL_EBADTOL 13
405#define GSL_ERROR(a,b) oocoutE((TObject*)0,Integration) << "RooAdaptiveGaussKronrodIntegrator1D::integral() ERROR: " << a << endl ; return b ;
406#define GSL_DBL_MIN 2.2250738585072014e-308
407#define GSL_DBL_MAX 1.7976931348623157e+308
408#define GSL_DBL_EPSILON 2.2204460492503131e-16
411#define GSL_EMAXITER 3
414#define GSL_EDIVERGE 6
417#define GSL_ERROR_VAL(reason, gsl_errno, value) return value ;
419#define GSL_MAX(a,b) ((a) > (b) ? (a) : (b))
435#define GSL_COERCE_DBL(x) (gsl_coerce_double(x))
448 double *result,
double *abserr,
449 double *defabs,
double *resabs);
452 double *result,
double *abserr,
453 double *resabs,
double *resasc);
456 double *result,
double *abserr,
457 double *resabs,
double *resasc);
460 double *result,
double *abserr,
461 double *resabs,
double *resasc);
464 double *result,
double *abserr,
465 double *resabs,
double *resasc);
468 double *result,
double *abserr,
469 double *resabs,
double *resasc);
472 double *result,
double *abserr,
473 double *resabs,
double *resasc);
476 double *cheb12,
double *cheb24);
494 const double wg[],
const double wgk[],
495 double fv1[],
double fv2[],
497 double * result,
double * abserr,
498 double * resabs,
double * resasc);
503 double epsabs,
double epsrel,
size_t limit,
505 gsl_integration_workspace * workspace,
506 double *result,
double *abserr);
511void initialise (gsl_integration_workspace * workspace,
double a,
double b);
514void initialise (gsl_integration_workspace * workspace,
double a,
double b)
517 workspace->nrmax = 0;
519 workspace->alist[0] =
a;
520 workspace->blist[0] =
b;
521 workspace->rlist[0] = 0.0;
522 workspace->elist[0] = 0.0;
523 workspace->order[0] = 0;
524 workspace->level[0] = 0;
526 workspace->maximum_level = 0;
532 double result,
double error);
536 double result,
double error)
539 workspace->rlist[0] = result;
540 workspace->elist[0] = error;
545qpsrt (gsl_integration_workspace * workspace);
548void qpsrt (gsl_integration_workspace * workspace)
550 const size_t last = workspace->size - 1;
551 const size_t limit = workspace->limit;
553 double * elist = workspace->elist;
554 size_t * order = workspace->order;
560 size_t i_nrmax = workspace->nrmax;
561 size_t i_maxerr = order[i_nrmax] ;
569 workspace->i = i_maxerr ;
573 errmax = elist[i_maxerr] ;
580 while (i_nrmax > 0 && errmax > elist[order[i_nrmax - 1]])
582 order[i_nrmax] = order[i_nrmax - 1] ;
590 if(last < (limit/2 + 2))
596 top = limit - last + 1;
607 while (i < top && errmax < elist[order[i]])
609 order[i-1] = order[i] ;
613 order[i-1] = i_maxerr ;
617 errmin = elist[last] ;
621 while (k > i - 2 && errmin >= elist[order[k]])
623 order[k+1] = order[k] ;
631 i_maxerr = order[i_nrmax] ;
633 workspace->i = i_maxerr ;
634 workspace->nrmax = i_nrmax ;
641void update (gsl_integration_workspace * workspace,
642 double a1,
double b1,
double area1,
double error1,
643 double a2,
double b2,
double area2,
double error2);
646retrieve (
const gsl_integration_workspace * workspace,
647 double *
a,
double *
b,
double *
r,
double *
e);
652void update (gsl_integration_workspace * workspace,
653 double a1,
double b1,
double area1,
double error1,
654 double a2,
double b2,
double area2,
double error2)
656 double * alist = workspace->alist ;
657 double * blist = workspace->blist ;
658 double * rlist = workspace->rlist ;
659 double * elist = workspace->elist ;
660 size_t * level = workspace->level ;
662 const size_t i_max = workspace->i ;
663 const size_t i_new = workspace->size ;
665 const size_t new_level = workspace->level[i_max] + 1;
672 rlist[i_max] = area2;
673 elist[i_max] = error2;
674 level[i_max] = new_level;
678 rlist[i_new] = area1;
679 elist[i_new] = error1;
680 level[i_new] = new_level;
685 rlist[i_max] = area1;
686 elist[i_max] = error1;
687 level[i_max] = new_level;
691 rlist[i_new] = area2;
692 elist[i_new] = error2;
693 level[i_new] = new_level;
698 if (new_level > workspace->maximum_level)
700 workspace->maximum_level = new_level;
707retrieve (
const gsl_integration_workspace * workspace,
708 double *
a,
double *
b,
double *
r,
double *
e)
710 const size_t i = workspace->i;
711 double * alist = workspace->alist;
712 double * blist = workspace->blist;
713 double * rlist = workspace->rlist;
714 double * elist = workspace->elist;
723sum_results (
const gsl_integration_workspace * workspace);
728 const double *
const rlist = workspace->rlist ;
729 const size_t n = workspace->size;
732 double result_sum = 0;
734 for (k = 0; k <
n; k++)
736 result_sum += rlist[k];
751 double tmp = (1 + 100 *
e) * (fabs (a2) + 1000 * u);
753 int status = fabs (a1) <= tmp && fabs (b2) <= tmp;
761 const double a,
const double b,
762 const double epsabs,
const double epsrel,
764 gsl_integration_workspace * workspace,
765 double * result,
double * abserr,
771 double epsabs,
double epsrel,
size_t limit,
773 gsl_integration_workspace * workspace,
774 double * result,
double * abserr)
810 status =
qag (
f,
a,
b, epsabs, epsrel, limit,
820 const double a,
const double b,
821 const double epsabs,
const double epsrel,
823 gsl_integration_workspace * workspace,
824 double *result,
double *abserr,
828 double result0, abserr0, resabs0, resasc0;
830 size_t iteration = 0;
831 int roundoff_type1 = 0, roundoff_type2 = 0, error_type = 0;
842 if (limit > workspace->limit)
847 if (epsabs <= 0 && (epsrel < 50 *
GSL_DBL_EPSILON || epsrel < 0.5e-28))
849 GSL_ERROR (
"tolerance cannot be acheived with given epsabs and epsrel",
855 q (
f,
a,
b, &result0, &abserr0, &resabs0, &resasc0);
861 tolerance =
GSL_MAX_DBL (epsabs, epsrel * fabs (result0));
867 if (abserr0 <= round_off && abserr0 > tolerance)
872 GSL_ERROR (
"cannot reach tolerance because of roundoff error "
875 else if ((abserr0 <= tolerance && abserr0 != resasc0) || abserr0 == 0.0)
897 double a1, b1, a2, b2;
898 double a_i, b_i, r_i, e_i;
899 double area1 = 0, area2 = 0, area12 = 0;
900 double error1 = 0, error2 = 0, error12 = 0;
901 double resasc1, resasc2;
902 double resabs1, resabs2;
906 retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
909 b1 = 0.5 * (a_i + b_i);
913 q (
f, a1, b1, &area1, &error1, &resabs1, &resasc1);
914 q (
f, a2, b2, &area2, &error2, &resabs2, &resasc2);
916 area12 = area1 + area2;
917 error12 = error1 + error2;
919 errsum += (error12 - e_i);
920 area += area12 - r_i;
922 if (resasc1 != error1 && resasc2 != error2)
924 double delta = r_i - area12;
926 if (fabs (delta) <= 1.0e-5 * fabs (area12) && error12 >= 0.99 * e_i)
930 if (iteration >= 10 && error12 > e_i)
936 tolerance =
GSL_MAX_DBL (epsabs, epsrel * fabs (area));
938 if (errsum > tolerance)
940 if (roundoff_type1 >= 6 || roundoff_type2 >= 20)
954 update (workspace, a1, b1, area1, error1, a2, b2, area2, error2);
956 retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
961 while (iteration < limit && !error_type && errsum > tolerance);
966 if (errsum <= tolerance)
970 else if (error_type == 2)
972 GSL_ERROR (
"roundoff error prevents tolerance from being achieved",
975 else if (error_type == 3)
977 GSL_ERROR (
"bad integrand behavior found in the integration interval",
980 else if (iteration == limit)
990static double rescale_error (
double err,
const double result_abs,
const double result_asc) ;
997 if (result_asc != 0 && err != 0)
999 double scale =
TMath::Power((200 * err / result_asc), 1.5) ;
1003 err = result_asc * scale ;
1026 const double xgk[],
const double wg[],
const double wgk[],
1027 double fv1[],
double fv2[],
1029 double *result,
double *abserr,
1030 double *resabs,
double *resasc)
1033 const double center = 0.5 * (
a +
b);
1034 const double half_length = 0.5 * (
b -
a);
1035 const double abs_half_length = fabs (half_length);
1038 double result_gauss = 0;
1039 double result_kronrod = f_center * wgk[
n - 1];
1041 double result_abs = fabs (result_kronrod);
1042 double result_asc = 0;
1043 double mean = 0, err = 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] * (fabs (fval1) + fabs (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] * (fabs (fval1) + fabs (fval2));
1078 mean = result_kronrod * 0.5;
1080 result_asc = wgk[
n - 1] * fabs (f_center - mean);
1082 for (j = 0; j <
n - 1; j++)
1084 result_asc += wgk[j] * (fabs (fv1[j] - mean) + fabs (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)
1146 double fv1[8], fv2[8];
1148 gsl_integration_qk (8,
xgkA,
wgA,
wgkA, fv1, fv2,
f,
a,
b, result, abserr, resabs, resasc);
1157 0.995657163025808080735527280689003,
1158 0.973906528517171720077964012084452,
1159 0.930157491355708226001207180059508,
1160 0.865063366688984510732096688423493,
1161 0.780817726586416897063717578345042,
1162 0.679409568299024406234327365114874,
1163 0.562757134668604683339000099272694,
1164 0.433395394129247190799265943165784,
1165 0.294392862701460198131126603103866,
1166 0.148874338981631210884826001129720,
1167 0.000000000000000000000000000000000
1175 0.066671344308688137593568809893332,
1176 0.149451349150580593145776339657697,
1177 0.219086362515982043995534934228163,
1178 0.269266719309996355091226921569469,
1179 0.295524224714752870173892994651338
1184 0.011694638867371874278064396062192,
1185 0.032558162307964727478818972459390,
1186 0.054755896574351996031381300244580,
1187 0.075039674810919952767043140916190,
1188 0.093125454583697605535065465083366,
1189 0.109387158802297641899210590325805,
1190 0.123491976262065851077958109831074,
1191 0.134709217311473325928054001771707,
1192 0.142775938577060080797094273138717,
1193 0.147739104901338491374841515972068,
1194 0.149445554002916905664936468389821
1200 double *result,
double *abserr,
1201 double *resabs,
double *resasc)
1203 double fv1[11], fv2[11];
1205 gsl_integration_qk (11,
xgkB,
wgB,
wgkB, fv1, fv2,
f,
a,
b, result, abserr, resabs, resasc);
1214 0.998002298693397060285172840152271,
1215 0.987992518020485428489565718586613,
1216 0.967739075679139134257347978784337,
1217 0.937273392400705904307758947710209,
1218 0.897264532344081900882509656454496,
1219 0.848206583410427216200648320774217,
1220 0.790418501442465932967649294817947,
1221 0.724417731360170047416186054613938,
1222 0.650996741297416970533735895313275,
1223 0.570972172608538847537226737253911,
1224 0.485081863640239680693655740232351,
1225 0.394151347077563369897207370981045,
1226 0.299180007153168812166780024266389,
1227 0.201194093997434522300628303394596,
1228 0.101142066918717499027074231447392,
1229 0.000000000000000000000000000000000
1237 0.030753241996117268354628393577204,
1238 0.070366047488108124709267416450667,
1239 0.107159220467171935011869546685869,
1240 0.139570677926154314447804794511028,
1241 0.166269205816993933553200860481209,
1242 0.186161000015562211026800561866423,
1243 0.198431485327111576456118326443839,
1244 0.202578241925561272880620199967519
1249 0.005377479872923348987792051430128,
1250 0.015007947329316122538374763075807,
1251 0.025460847326715320186874001019653,
1252 0.035346360791375846222037948478360,
1253 0.044589751324764876608227299373280,
1254 0.053481524690928087265343147239430,
1255 0.062009567800670640285139230960803,
1256 0.069854121318728258709520077099147,
1257 0.076849680757720378894432777482659,
1258 0.083080502823133021038289247286104,
1259 0.088564443056211770647275443693774,
1260 0.093126598170825321225486872747346,
1261 0.096642726983623678505179907627589,
1262 0.099173598721791959332393173484603,
1263 0.100769845523875595044946662617570,
1264 0.101330007014791549017374792767493
1269 double *result,
double *abserr,
1270 double *resabs,
double *resasc)
1272 double fv1[16], fv2[16];
1274 gsl_integration_qk (16,
xgkC,
wgC,
wgkC, fv1, fv2,
f,
a,
b, result, abserr, resabs, resasc);
1283 0.998859031588277663838315576545863,
1284 0.993128599185094924786122388471320,
1285 0.981507877450250259193342994720217,
1286 0.963971927277913791267666131197277,
1287 0.940822633831754753519982722212443,
1288 0.912234428251325905867752441203298,
1289 0.878276811252281976077442995113078,
1290 0.839116971822218823394529061701521,
1291 0.795041428837551198350638833272788,
1292 0.746331906460150792614305070355642,
1293 0.693237656334751384805490711845932,
1294 0.636053680726515025452836696226286,
1295 0.575140446819710315342946036586425,
1296 0.510867001950827098004364050955251,
1297 0.443593175238725103199992213492640,
1298 0.373706088715419560672548177024927,
1299 0.301627868114913004320555356858592,
1300 0.227785851141645078080496195368575,
1301 0.152605465240922675505220241022678,
1302 0.076526521133497333754640409398838,
1303 0.000000000000000000000000000000000
1311 0.017614007139152118311861962351853,
1312 0.040601429800386941331039952274932,
1313 0.062672048334109063569506535187042,
1314 0.083276741576704748724758143222046,
1315 0.101930119817240435036750135480350,
1316 0.118194531961518417312377377711382,
1317 0.131688638449176626898494499748163,
1318 0.142096109318382051329298325067165,
1319 0.149172986472603746787828737001969,
1320 0.152753387130725850698084331955098
1325 0.003073583718520531501218293246031,
1326 0.008600269855642942198661787950102,
1327 0.014626169256971252983787960308868,
1328 0.020388373461266523598010231432755,
1329 0.025882133604951158834505067096153,
1330 0.031287306777032798958543119323801,
1331 0.036600169758200798030557240707211,
1332 0.041668873327973686263788305936895,
1333 0.046434821867497674720231880926108,
1334 0.050944573923728691932707670050345,
1335 0.055195105348285994744832372419777,
1336 0.059111400880639572374967220648594,
1337 0.062653237554781168025870122174255,
1338 0.065834597133618422111563556969398,
1339 0.068648672928521619345623411885368,
1340 0.071054423553444068305790361723210,
1341 0.073030690332786667495189417658913,
1342 0.074582875400499188986581418362488,
1343 0.075704497684556674659542775376617,
1344 0.076377867672080736705502835038061,
1345 0.076600711917999656445049901530102
1350 double *result,
double *abserr,
1351 double *resabs,
double *resasc)
1353 double fv1[21], fv2[21];
1355 gsl_integration_qk (21,
xgkD,
wgD,
wgkD, fv1, fv2,
f,
a,
b, result, abserr, resabs, resasc);
1364 0.999262104992609834193457486540341,
1365 0.995556969790498097908784946893902,
1366 0.988035794534077247637331014577406,
1367 0.976663921459517511498315386479594,
1368 0.961614986425842512418130033660167,
1369 0.942974571228974339414011169658471,
1370 0.920747115281701561746346084546331,
1371 0.894991997878275368851042006782805,
1372 0.865847065293275595448996969588340,
1373 0.833442628760834001421021108693570,
1374 0.797873797998500059410410904994307,
1375 0.759259263037357630577282865204361,
1376 0.717766406813084388186654079773298,
1377 0.673566368473468364485120633247622,
1378 0.626810099010317412788122681624518,
1379 0.577662930241222967723689841612654,
1380 0.526325284334719182599623778158010,
1381 0.473002731445714960522182115009192,
1382 0.417885382193037748851814394594572,
1383 0.361172305809387837735821730127641,
1384 0.303089538931107830167478909980339,
1385 0.243866883720988432045190362797452,
1386 0.183718939421048892015969888759528,
1387 0.122864692610710396387359818808037,
1388 0.061544483005685078886546392366797,
1389 0.000000000000000000000000000000000
1397 0.011393798501026287947902964113235,
1398 0.026354986615032137261901815295299,
1399 0.040939156701306312655623487711646,
1400 0.054904695975835191925936891540473,
1401 0.068038333812356917207187185656708,
1402 0.080140700335001018013234959669111,
1403 0.091028261982963649811497220702892,
1404 0.100535949067050644202206890392686,
1405 0.108519624474263653116093957050117,
1406 0.114858259145711648339325545869556,
1407 0.119455763535784772228178126512901,
1408 0.122242442990310041688959518945852,
1409 0.123176053726715451203902873079050
1414 0.001987383892330315926507851882843,
1415 0.005561932135356713758040236901066,
1416 0.009473973386174151607207710523655,
1417 0.013236229195571674813656405846976,
1418 0.016847817709128298231516667536336,
1419 0.020435371145882835456568292235939,
1420 0.024009945606953216220092489164881,
1421 0.027475317587851737802948455517811,
1422 0.030792300167387488891109020215229,
1423 0.034002130274329337836748795229551,
1424 0.037116271483415543560330625367620,
1425 0.040083825504032382074839284467076,
1426 0.042872845020170049476895792439495,
1427 0.045502913049921788909870584752660,
1428 0.047982537138836713906392255756915,
1429 0.050277679080715671963325259433440,
1430 0.052362885806407475864366712137873,
1431 0.054251129888545490144543370459876,
1432 0.055950811220412317308240686382747,
1433 0.057437116361567832853582693939506,
1434 0.058689680022394207961974175856788,
1435 0.059720340324174059979099291932562,
1436 0.060539455376045862945360267517565,
1437 0.061128509717053048305859030416293,
1438 0.061471189871425316661544131965264,
1439 0.061580818067832935078759824240066
1446 double *result,
double *abserr,
1447 double *resabs,
double *resasc)
1449 double fv1[26], fv2[26];
1451 gsl_integration_qk (26,
xgkE,
wgE,
wgkE, fv1, fv2,
f,
a,
b, result, abserr, resabs, resasc);
1460 0.999484410050490637571325895705811,
1461 0.996893484074649540271630050918695,
1462 0.991630996870404594858628366109486,
1463 0.983668123279747209970032581605663,
1464 0.973116322501126268374693868423707,
1465 0.960021864968307512216871025581798,
1466 0.944374444748559979415831324037439,
1467 0.926200047429274325879324277080474,
1468 0.905573307699907798546522558925958,
1469 0.882560535792052681543116462530226,
1470 0.857205233546061098958658510658944,
1471 0.829565762382768397442898119732502,
1472 0.799727835821839083013668942322683,
1473 0.767777432104826194917977340974503,
1474 0.733790062453226804726171131369528,
1475 0.697850494793315796932292388026640,
1476 0.660061064126626961370053668149271,
1477 0.620526182989242861140477556431189,
1478 0.579345235826361691756024932172540,
1479 0.536624148142019899264169793311073,
1480 0.492480467861778574993693061207709,
1481 0.447033769538089176780609900322854,
1482 0.400401254830394392535476211542661,
1483 0.352704725530878113471037207089374,
1484 0.304073202273625077372677107199257,
1485 0.254636926167889846439805129817805,
1486 0.204525116682309891438957671002025,
1487 0.153869913608583546963794672743256,
1488 0.102806937966737030147096751318001,
1489 0.051471842555317695833025213166723,
1490 0.000000000000000000000000000000000
1498 0.007968192496166605615465883474674,
1499 0.018466468311090959142302131912047,
1500 0.028784707883323369349719179611292,
1501 0.038799192569627049596801936446348,
1502 0.048402672830594052902938140422808,
1503 0.057493156217619066481721689402056,
1504 0.065974229882180495128128515115962,
1505 0.073755974737705206268243850022191,
1506 0.080755895229420215354694938460530,
1507 0.086899787201082979802387530715126,
1508 0.092122522237786128717632707087619,
1509 0.096368737174644259639468626351810,
1510 0.099593420586795267062780282103569,
1511 0.101762389748405504596428952168554,
1512 0.102852652893558840341285636705415
1517 0.001389013698677007624551591226760,
1518 0.003890461127099884051267201844516,
1519 0.006630703915931292173319826369750,
1520 0.009273279659517763428441146892024,
1521 0.011823015253496341742232898853251,
1522 0.014369729507045804812451432443580,
1523 0.016920889189053272627572289420322,
1524 0.019414141193942381173408951050128,
1525 0.021828035821609192297167485738339,
1526 0.024191162078080601365686370725232,
1527 0.026509954882333101610601709335075,
1528 0.028754048765041292843978785354334,
1529 0.030907257562387762472884252943092,
1530 0.032981447057483726031814191016854,
1531 0.034979338028060024137499670731468,
1532 0.036882364651821229223911065617136,
1533 0.038678945624727592950348651532281,
1534 0.040374538951535959111995279752468,
1535 0.041969810215164246147147541285970,
1536 0.043452539701356069316831728117073,
1537 0.044814800133162663192355551616723,
1538 0.046059238271006988116271735559374,
1539 0.047185546569299153945261478181099,
1540 0.048185861757087129140779492298305,
1541 0.049055434555029778887528165367238,
1542 0.049795683427074206357811569379942,
1543 0.050405921402782346840893085653585,
1544 0.050881795898749606492297473049805,
1545 0.051221547849258772170656282604944,
1546 0.051426128537459025933862879215781,
1547 0.051494729429451567558340433647099
1552 double *result,
double *abserr,
1553 double *resabs,
double *resasc)
1555 double fv1[31], fv2[31];
1557 gsl_integration_qk (31,
xgkF,
wgF,
wgkF, fv1, fv2,
f,
a,
b, result, abserr, resabs, resasc);
1560gsl_integration_workspace*
1563 gsl_integration_workspace* w ;
1567 GSL_ERROR_VAL (
"workspace length n must be positive integer",
1571 w = (gsl_integration_workspace *)
1572 malloc (
sizeof (gsl_integration_workspace));
1576 GSL_ERROR_VAL (
"failed to allocate space for workspace struct",
1580 w->alist = (
double *)
malloc (
n *
sizeof (
double));
1590 w->blist = (
double *)
malloc (
n *
sizeof (
double));
1601 w->rlist = (
double *)
malloc (
n *
sizeof (
double));
1614 w->elist = (
double *)
malloc (
n *
sizeof (
double));
1627 w->order = (
size_t *)
malloc (
n *
sizeof (
size_t));
1641 w->level = (
size_t *)
malloc (
n *
sizeof (
size_t));
1658 w->maximum_level = 0 ;
1679reset_nrmax (gsl_integration_workspace * workspace);
1684 workspace->nrmax = 0;
1685 workspace->i = workspace->order[0] ;
1701 int id = workspace->nrmax;
1704 const size_t * level = workspace->level;
1705 const size_t * order = workspace->order;
1707 size_t limit = workspace->limit ;
1708 size_t last = workspace->size - 1 ;
1710 if (last > (1 + limit / 2))
1712 jupbnd = limit + 1 - last;
1719 for (k =
id; k <= jupbnd; k++)
1721 size_t i_max = order[workspace->nrmax];
1723 workspace->i = i_max ;
1725 if (level[i_max] < workspace->maximum_level)
1739 size_t i = workspace->i ;
1740 const size_t * level = workspace->level;
1742 if (level[i] < workspace->maximum_level)
1754struct extrapolation_table
1779 table->rlist2[0] =
y;
1788 table->rlist2[
n] =
y;
1798 qelg (
struct extrapolation_table *table,
double *result,
double *abserr);
1801qelg (
struct extrapolation_table *table,
double *result,
double *abserr)
1803 double *epstab = table->rlist2;
1804 double *res3la = table->res3la;
1805 const size_t n = table->n - 1;
1807 const double current = epstab[
n];
1812 const size_t newelm =
n / 2;
1813 const size_t n_orig =
n;
1817 const size_t nres_orig = table->nres;
1829 epstab[
n + 2] = epstab[
n];
1832 for (i = 0; i < newelm; i++)
1834 double res = epstab[
n - 2 * i + 2];
1835 double e0 = epstab[
n - 2 * i - 2];
1836 double e1 = epstab[
n - 2 * i - 1];
1839 double e1abs = fabs (e1);
1840 double delta2 = e2 - e1;
1841 double err2 = fabs (delta2);
1843 double delta3 = e1 - e0;
1844 double err3 = fabs (delta3);
1847 double e3, delta1, err1, tol1, ss;
1849 if (err2 <= tol2 && err3 <= tol3)
1855 absolute = err2 + err3;
1861 e3 = epstab[
n - 2 * i];
1862 epstab[
n - 2 * i] = e1;
1864 err1 = fabs (delta1);
1870 if (err1 <= tol1 || err2 <= tol2 || err3 <= tol3)
1876 ss = (1 / delta1 + 1 / delta2) - 1 / delta3;
1882 if (fabs (ss * e1) <= 0.0001)
1892 epstab[
n - 2 * i] = res;
1895 const double error = err2 + fabs (res - e2) + err3;
1897 if (error <= *abserr)
1908 const size_t limexp = 50 - 1;
1910 if (n_final == limexp)
1912 n_final = 2 * (limexp / 2);
1916 if (n_orig % 2 == 1)
1918 for (i = 0; i <= newelm; i++)
1920 epstab[1 + i * 2] = epstab[i * 2 + 3];
1925 for (i = 0; i <= newelm; i++)
1927 epstab[i * 2] = epstab[i * 2 + 2];
1931 if (n_orig != n_final)
1933 for (i = 0; i <= n_final; i++)
1935 epstab[i] = epstab[n_orig - n_final + i];
1939 table->n = n_final + 1;
1943 res3la[nres_orig] = *result;
1948 *abserr = (fabs (*result - res3la[2]) + fabs (*result - res3la[1])
1949 + fabs (*result - res3la[0]));
1951 res3la[0] = res3la[1];
1952 res3la[1] = res3la[2];
1953 res3la[2] = *result;
1962 table->nres = nres_orig + 1;
1986 b,
const double epsabs,
const double epsrel,
const size_t limit,
1987 gsl_integration_workspace * workspace,
double *result,
double *abserr,
1993 double epsabs,
double epsrel,
size_t limit,
1994 gsl_integration_workspace * workspace,
1995 double * result,
double * abserr)
1997 int status =
qags (
f,
a,
b, epsabs, epsrel, limit,
2011static double i_transform (
double t,
void *params);
2015 double epsabs,
double epsrel,
size_t limit,
2016 gsl_integration_workspace * workspace,
2017 double *result,
double *abserr)
2024 f_transform.params =
f;
2026 status =
qags (&f_transform, 0.0, 1.0,
2027 epsabs, epsrel, limit,
2039 double x = (1 - t) / t;
2059 double epsabs,
double epsrel,
size_t limit,
2060 gsl_integration_workspace * workspace,
2061 double *result,
double *abserr)
2066 struct il_params transform_params ;
2068 transform_params.b =
b ;
2069 transform_params.f =
f ;
2072 f_transform.params = &transform_params;
2074 status =
qags (&f_transform, 0.0, 1.0,
2075 epsabs, epsrel, limit,
2086 struct il_params *p = (
struct il_params *) params;
2089 double x =
b - (1 - t) / t;
2108 double epsabs,
double epsrel,
size_t limit,
2109 gsl_integration_workspace * workspace,
2110 double *result,
double *abserr)
2115 struct iu_params transform_params ;
2117 transform_params.a =
a ;
2118 transform_params.f =
f ;
2121 f_transform.params = &transform_params;
2123 status =
qags (&f_transform, 0.0, 1.0,
2124 epsabs, epsrel, limit,
2135 struct iu_params *p = (
struct iu_params *) params;
2138 double x =
a + (1 - t) / t;
2147 const double a,
const double b,
2148 const double epsabs,
const double epsrel,
2150 gsl_integration_workspace * workspace,
2151 double *result,
double *abserr,
2154 double area, errsum;
2155 double res_ext, err_ext;
2156 double result0, abserr0, resabs0, resasc0;
2160 double error_over_large_intervals = 0;
2161 double reseps = 0, abseps = 0, correc = 0;
2163 int roundoff_type1 = 0, roundoff_type2 = 0, roundoff_type3 = 0;
2164 int error_type = 0, error_type2 = 0;
2166 size_t iteration = 0;
2168 int positive_integrand = 0;
2169 int extrapolate = 0;
2170 int disallow_extrapolation = 0;
2172 struct extrapolation_table table;
2181 if (limit > workspace->limit)
2188 if (epsabs <= 0 && (epsrel < 50 *
GSL_DBL_EPSILON || epsrel < 0.5e-28))
2190 GSL_ERROR (
"tolerance cannot be acheived with given epsabs and epsrel",
2196 q (
f,
a,
b, &result0, &abserr0, &resabs0, &resasc0);
2200 tolerance =
GSL_MAX_DBL (epsabs, epsrel * fabs (result0));
2202 if (abserr0 <= 100 * GSL_DBL_EPSILON * resabs0 && abserr0 > tolerance)
2207 GSL_ERROR (
"cannot reach tolerance because of roundoff error"
2210 else if ((abserr0 <= tolerance && abserr0 != resasc0) || abserr0 == 0.0)
2217 else if (limit == 1)
2242 size_t current_level;
2243 double a1, b1, a2, b2;
2244 double a_i, b_i, r_i, e_i;
2245 double area1 = 0, area2 = 0, area12 = 0;
2246 double error1 = 0, error2 = 0, error12 = 0;
2247 double resasc1, resasc2;
2248 double resabs1, resabs2;
2253 retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
2255 current_level = workspace->level[workspace->i] + 1;
2258 b1 = 0.5 * (a_i + b_i);
2264 q (
f, a1, b1, &area1, &error1, &resabs1, &resasc1);
2265 q (
f, a2, b2, &area2, &error2, &resabs2, &resasc2);
2267 area12 = area1 + area2;
2268 error12 = error1 + error2;
2278 errsum = errsum + error12 - e_i;
2279 area = area + area12 - r_i;
2281 tolerance =
GSL_MAX_DBL (epsabs, epsrel * fabs (area));
2283 if (resasc1 != error1 && resasc2 != error2)
2285 double delta = r_i - area12;
2287 if (fabs (delta) <= 1.0e-5 * fabs (area12) && error12 >= 0.99 * e_i)
2298 if (iteration > 10 && error12 > e_i)
2306 if (roundoff_type1 + roundoff_type2 >= 10 || roundoff_type3 >= 20)
2311 if (roundoff_type2 >= 5)
2326 update (workspace, a1, b1, area1, error1, a2, b2, area2, error2);
2328 if (errsum <= tolerance)
2330 goto compute_result;
2338 if (iteration >= limit - 1)
2346 error_over_large_intervals = errsum;
2352 if (disallow_extrapolation)
2357 error_over_large_intervals += -last_e_i;
2359 if (current_level < workspace->maximum_level)
2361 error_over_large_intervals += error12;
2373 workspace->nrmax = 1;
2376 if (!error_type2 && error_over_large_intervals > ertest)
2386 qelg (&table, &reseps, &abseps);
2390 if (ktmin > 5 && err_ext < 0.001 * errsum)
2395 if (abseps < err_ext)
2400 correc = error_over_large_intervals;
2401 ertest =
GSL_MAX_DBL (epsabs, epsrel * fabs (reseps));
2402 if (err_ext <= ertest)
2410 disallow_extrapolation = 1;
2413 if (error_type == 5)
2422 error_over_large_intervals = errsum;
2425 while (iteration < limit);
2431 goto compute_result;
2433 if (error_type || error_type2)
2443 if (res_ext != 0.0 && area != 0.0)
2445 if (err_ext / fabs (res_ext) > errsum / fabs (area))
2446 goto compute_result;
2448 else if (err_ext > errsum)
2450 goto compute_result;
2452 else if (area == 0.0)
2461 double max_area =
GSL_MAX_DBL (fabs (res_ext), fabs (area));
2463 if (!positive_integrand && max_area < 0.01 * resabs0)
2468 double ratio = res_ext / area;
2470 if (ratio < 0.01 || ratio > 100.0 || errsum > fabs (area))
2488 if (error_type == 0)
2492 else if (error_type == 1)
2496 else if (error_type == 2)
2498 GSL_ERROR (
"cannot reach tolerance because of roundoff error",
2501 else if (error_type == 3)
2503 GSL_ERROR (
"bad integrand behavior found in the integration interval",
2506 else if (error_type == 4)
2508 GSL_ERROR (
"roundoff error detected in the extrapolation table",
2511 else if (error_type == 5)
2513 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)
struct gsl_function_struct gsl_function
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)
Int_t getCatIndex(const char *name, Int_t defVal=0, Bool_t verbose=kFALSE) const
Get index value of a RooAbsCategory stored in set with given name.
Double_t getRealValue(const char *name, Double_t defVal=0, Bool_t verbose=kFALSE) const
Get value of a RooAbsReal stored in set with given name.
Abstract interface for evaluating a real-valued function of one real variable and performing numerica...
virtual Double_t getMinLimit(UInt_t dimension) const =0
virtual Double_t getMaxLimit(UInt_t dimension) const =0
UInt_t getDimension() const
RooAbsIntegrator is the abstract interface for integrators of real-valued functions that implement th...
const RooAbsFunc * _function
const RooAbsFunc * integrand() const
RooAdaptiveGaussKronrodIntegrator1D implements the Gauss-Kronrod integration algorithm.
virtual RooAbsIntegrator * clone(const RooAbsFunc &function, const RooNumIntConfig &config) const
Virtual constructor.
Double_t _xmax
Lower integration bound.
Double_t _epsAbs
Current coordinate.
virtual Bool_t checkLimits() const
Check that our integration range is finite and otherwise return kFALSE.
Bool_t setLimits(Double_t *xmin, Double_t *xmax)
Change our integration limits.
RooAdaptiveGaussKronrodIntegrator1D()
coverity[UNINIT_CTOR] Default constructor
friend double RooAdaptiveGaussKronrodIntegrator1D_GSL_GlueFunction(double x, void *data)
Glue function interacing to GSL code.
Bool_t initialize()
Initialize integrator allocate buffers and setup GSL workspace.
virtual ~RooAdaptiveGaussKronrodIntegrator1D()
Destructor.
virtual Double_t integral(const Double_t *yvec=0)
Calculate and return integral at at given parameter values.
static void registerIntegrator(RooNumIntFactory &fact)
Register this class with RooNumIntConfig as a possible choice of numeric integrator for one-dimension...
Bool_t _useIntegrandLimits
RooArgSet is a container object that can hold multiple RooAbsArg objects.
RooCategory is an object to represent discrete states.
bool defineType(const std::string &label)
Define a state with given name.
virtual Bool_t setIndex(Int_t index, bool printError=true) override
Set value by specifying the index code of the desired state.
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_t 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 Int_t isInfinite(Double_t x)
Return true if x is infinite by RooNumBer internal specification.
RooRealVar represents a variable that can be changed from the outside.
Mother of all ROOT objects.
virtual const char * GetName() const
Returns name of object.
static Roo_reg_AGKInteg1D instance
LongDouble_t Power(LongDouble_t x, LongDouble_t y)