42 template<
class Float,
class integr_fn>
static void rdqagie(integr_fn f,
void *ex,
43 Float *,
int *, Float * , Float *,
int *,
44 Float *, Float *,
int *,
45 int *, Float *, Float *, Float *, Float *,
48 template<
class Float,
class integr_fn>
static void rdqk15i(integr_fn f,
void *ex,
49 Float *,
int *, Float * , Float *,
50 Float *, Float *, Float *, Float *);
52 template<
class Float,
class integr_fn>
static void rdqagse(integr_fn f,
void *ex, Float *, Float *,
53 Float *, Float *,
int *, Float *, Float *,
54 int *,
int *, Float *, Float *, Float *,
55 Float *,
int *,
int *);
57 template<
class Float,
class integr_fn>
static void rdqk21(integr_fn f,
void *ex,
58 Float *, Float *, Float *, Float *, Float *, Float *);
60 template<
class Float>
static void rdqpsrt(
int *,
int *,
int *, Float *, Float *,
int *,
int *);
62 template<
class Float>
static void rdqelg(
int *, Float *, Float *, Float *, Float *,
int *);
64 template<
class Float,
class integr_fn>
65 void Rdqagi(integr_fn f,
void *ex, Float *bound,
int *inf,
66 Float *epsabs, Float *epsrel,
67 Float *result, Float *abserr,
int *neval,
int *ier,
68 int *limit,
int *lenw,
int *last,
69 int *iwork, Float *work)
235 if (*limit < 1 || *lenw < *limit << 2)
return;
241 rdqagie(f, ex, bound, inf, epsabs, epsrel, limit, result, abserr, neval, ier,
242 work, &work[l1], &work[l2], &work[l3], iwork, last);
247 template<
class Float,
class integr_fn>
static 248 void rdqagie(integr_fn f,
void *ex, Float *bound,
int *inf, Float *
249 epsabs, Float *epsrel,
int *limit, Float *result,
250 Float *abserr,
int *neval,
int *ier, Float *alist,
251 Float *blist, Float *rlist, Float *elist,
int *
259 Float area1, area2, area12;
261 Float small = 0.0, erro12;
263 Float a1, a2, b1, b2, defab1, defab2, oflow;
267 int iroff1, iroff2, iroff3;
268 Float res3la[3], error1, error2;
272 Float defabs, epmach, erlarg = 0.0, abseps, correc = 0.0, errbnd, resabs;
274 Float erlast, errmax;
278 Float ertest = 0.0, errsum;
489 epmach = DBL_EPSILON;
504 if (*epsabs <= 0. && (*epsrel < fmax2(epmach * 50., 5e-29))) *ier = 6;
505 if (*ier == 6)
return;
521 static Float c_b6 = 0.;
522 static Float c_b7 = 1.;
524 rdqk15i(f, ex, &boun, inf, &c_b6, &c_b7, result, abserr, &defabs, &resabs);
532 dres = fabs(*result);
533 errbnd = fmax2(*epsabs, *epsrel * dres);
534 if (*abserr <= epmach * 100. * defabs && *abserr > errbnd) *ier = 2;
535 if (*limit == 1) *ier = 1;
536 if (*ier != 0 || (*abserr <= errbnd && *abserr != resabs)
537 || *abserr == 0.)
goto L130;
561 if (dres >= (1. - epmach * 50.) * defabs) {
568 for (*last = 2; *last <= *limit; ++(*last)) {
573 b1 = (alist[maxerr] + blist[maxerr]) * .5;
577 rdqk15i(f, ex, &boun, inf, &a1, &b1, &area1, &error1, &resabs, &defab1);
578 rdqk15i(f, ex, &boun, inf, &a2, &b2, &area2, &error2, &resabs, &defab2);
583 area12 = area1 + area2;
584 erro12 = error1 + error2;
585 errsum = errsum + erro12 - errmax;
586 area = area + area12 - rlist[maxerr];
587 if (!(defab1 == error1 || defab2 == error2)) {
588 if (fabs(rlist[maxerr] - area12) <= fabs(area12) * 1e-5 &&
589 erro12 >= errmax * .99) {
595 if (*last > 10 && erro12 > errmax)
599 rlist[maxerr] = area1;
600 rlist[*last] = area2;
601 errbnd = fmax2(*epsabs, *epsrel * fabs(area));
605 if (iroff1 + iroff2 >= 10 || iroff3 >= 20)
619 if (fmax2(fabs(a1), fabs(b2)) <=
620 (epmach * 100. + 1.) * (fabs(a2) + uflow * 1e3))
627 if (error2 <= error1) {
631 elist[maxerr] = error1;
632 elist[*last] = error2;
638 rlist[maxerr] = area2;
639 rlist[*last] = area1;
640 elist[maxerr] = error2;
641 elist[*last] = error1;
648 rdqpsrt(limit, last, &maxerr, &errmax, &elist[1], &iord[1], &nrmax);
649 if (errsum <= errbnd) {
652 if (*ier != 0)
break;
657 rlist2[1] = area;
continue;
662 if (fabs(b1 - a1) > small) {
670 if (fabs(blist[maxerr] - alist[maxerr]) > small) {
677 if (ierro != 3 && erlarg > ertest) {
685 if (*last > *limit / 2 + 2) {
686 jupbnd = *limit + 3 - *last;
688 for (k =
id; k <= jupbnd; ++k) {
689 maxerr = iord[nrmax];
690 errmax = elist[maxerr];
691 if (fabs(blist[maxerr] - alist[maxerr]) > small) {
700 rlist2[numrl2 - 1] = area;
701 rdqelg(&numrl2, rlist2, &reseps, &abseps, res3la, &nres);
703 if (ktmin > 5 && *abserr < errsum * .001) {
706 if (abseps >= *abserr) {
713 ertest = fmax2(*epsabs, *epsrel * fabs(reseps));
714 if (*abserr <= ertest) {
728 errmax = elist[maxerr];
740 if (*abserr == oflow) {
743 if (*ier + ierro == 0) {
752 if (*result == 0. || area == 0.) {
753 if (*abserr > errsum)
760 if (*abserr / fabs(*result) > errsum / fabs(area)) {
767 if (ksgn == -1 && fmax2(fabs(*result), fabs(area)) <= defabs * .01) {
770 if (.01 > *result / area || *result / area > 100. || errsum > fabs(area)) {
779 for (k = 1; k <= *last; ++k)
784 *neval = *last * 30 - 15;
794 template<
class Float,
class integr_fn>
795 void Rdqags(integr_fn f,
void *ex, Float *a, Float *b,
796 Float *epsabs, Float *epsrel,
797 Float *result, Float *abserr,
int *neval,
int *ier,
798 int *limit,
int *lenw,
int *last,
int *iwork, Float *work)
961 if (*limit < 1 || *lenw < *limit *4)
return;
969 rdqagse(f, ex, a, b, epsabs, epsrel, limit, result, abserr, neval, ier,
970 work, &work[l1], &work[l2], &work[l3], iwork, last);
975 template<
class Float,
class integr_fn>
static 976 void rdqagse(integr_fn f,
void *ex, Float *a, Float *b, Float *
977 epsabs, Float *epsrel,
int *limit, Float *result,
978 Float *abserr,
int *neval,
int *ier, Float *alist,
979 Float *blist, Float *rlist, Float *elist,
int *
987 int iroff1, iroff2, iroff3;
994 Float abseps, area, area1, area2, area12, dres, epmach;
995 Float a1, a2, b1, b2, defabs, defab1, defab2, oflow, uflow, resabs, reseps;
996 Float error1, error2, erro12, errbnd, erlast, errmax, errsum;
998 Float correc = 0.0, erlarg = 0.0, ertest = 0.0, small = 0.0;
1210 epmach = DBL_EPSILON;
1223 if (*epsabs <= 0. && *epsrel < fmax2(epmach * 50., 5e-29)) {
1234 rdqk21(f, ex, a, b, result, abserr, &defabs, &resabs);
1238 dres = fabs(*result);
1239 errbnd = fmax2(*epsabs, *epsrel * dres);
1244 if (*abserr <= epmach * 100. * defabs && *abserr > errbnd)
1248 if (*ier != 0 || (*abserr <= errbnd && *abserr != resabs)
1249 || *abserr == 0.)
goto L140;
1254 rlist2[0] = *result;
1270 if (dres >= (1. - epmach * 50.) * defabs) {
1277 for (*last = 2; *last <= *limit; ++(*last)) {
1282 b1 = (alist[maxerr] + blist[maxerr]) * .5;
1286 rdqk21(f, ex, &a1, &b1, &area1, &error1, &resabs, &defab1);
1287 rdqk21(f, ex, &a2, &b2, &area2, &error2, &resabs, &defab2);
1292 area12 = area1 + area2;
1293 erro12 = error1 + error2;
1294 errsum = errsum + erro12 - errmax;
1295 area = area + area12 - rlist[maxerr];
1296 if (!(defab1 == error1 || defab2 == error2)) {
1298 if (fabs(rlist[maxerr] - area12) <= fabs(area12) * 1e-5 &&
1299 erro12 >= errmax * .99) {
1305 if (*last > 10 && erro12 > errmax)
1308 rlist[maxerr] = area1;
1309 rlist[*last] = area2;
1310 errbnd = fmax2(*epsabs, *epsrel * fabs(area));
1314 if (iroff1 + iroff2 >= 10 || iroff3 >= 20)
1320 if (*last == *limit)
1326 if (fmax2(fabs(a1), fabs(b2)) <=
1327 (epmach * 100. + 1.) * (fabs(a2) + uflow * 1e3)) {
1333 if (error2 > error1) {
1337 rlist[maxerr] = area2;
1338 rlist[*last] = area1;
1339 elist[maxerr] = error2;
1340 elist[*last] = error1;
1345 elist[maxerr] = error1;
1346 elist[*last] = error2;
1354 rdqpsrt(limit, last, &maxerr, &errmax, &elist[1], &iord[1], &nrmax);
1356 if (errsum <= errbnd)
goto L115;
1357 if (*ier != 0)
break;
1359 small = fabs(*b - *a) * .375;
1362 rlist2[1] = area;
continue;
1364 if (noext)
continue;
1367 if (fabs(b1 - a1) > small) {
1375 if (fabs(blist[maxerr] - alist[maxerr]) > small) {
1382 if (ierro != 3 && erlarg > ertest) {
1390 if (*last > *limit / 2 + 2) {
1391 jupbnd = *limit + 3 - *last;
1393 for (k =
id; k <= jupbnd; ++k) {
1394 maxerr = iord[nrmax];
1395 errmax = elist[maxerr];
1396 if (fabs(blist[maxerr] - alist[maxerr]) > small) {
1406 rlist2[numrl2 - 1] = area;
1407 rdqelg(&numrl2, rlist2, &reseps, &abseps, res3la, &nres);
1409 if (ktmin > 5 && *abserr < errsum * .001) {
1412 if (abseps < *abserr) {
1417 ertest = fmax2(*epsabs, *epsrel * fabs(reseps));
1418 if (*abserr <= ertest) {
1432 errmax = elist[maxerr];
1445 if (*abserr == oflow)
goto L115;
1446 if (*ier + ierro == 0)
goto L110;
1451 if (*result == 0. || area == 0.) {
1452 if (*abserr > errsum)
goto L115;
1453 if (area == 0.)
goto L130;
1456 if (*abserr / fabs(*result) > errsum / fabs(area))
1461 if (ksgn == -1 && fmax2(fabs(*result), fabs(area)) <= defabs * .01) {
1464 if (.01 > *result / area || *result / area > 100. || errsum > fabs(area)) {
1471 for (k = 1; k <= *last; ++k)
1472 *result += rlist[k];
1477 *neval = *last * 42 - 21;
1482 template<
class Float,
class integr_fn>
static void rdqk15i(integr_fn f,
void *ex,
1483 Float *boun,
int *inf, Float *a, Float *b,
1485 Float *abserr, Float *resabs, Float *resasc)
1489 static double wg[8] = {
1490 0., .129484966168869693270611432679082,
1491 0., .27970539148927666790146777142378,
1492 0., .381830050505118944950369775488975,
1493 0., .417959183673469387755102040816327 };
1494 static double xgk[8] = {
1495 .991455371120812639206854697526329,
1496 .949107912342758524526189684047851,
1497 .864864423359769072789712788640926,
1498 .741531185599394439863864773280788,
1499 .58608723546769113029414483825873,
1500 .405845151377397166906606412076961,
1501 .207784955007898467600689403773245, 0. };
1502 static double wgk[8] = {
1503 .02293532201052922496373200805897,
1504 .063092092629978553290700663189204,
1505 .104790010322250183839876322541518,
1506 .140653259715525918745189590510238,
1507 .16900472663926790282658342659855,
1508 .190350578064785409913256402421014,
1509 .204432940075298892414161999234649,
1510 .209482141084727828012999174891714 };
1513 Float absc, dinf, resg, resk, fsum, absc1, absc2, fval1, fval2;
1515 Float hlgth, centr, reskh, uflow;
1516 Float tabsc1, tabsc2, fc, epmach;
1517 Float fv1[7], fv2[7], vec[15], vec2[15];
1631 epmach = DBL_EPSILON;
1633 dinf = (double) imin2(1, *inf);
1635 centr = (*a + *b) * .5;
1636 hlgth = (*b - *a) * .5;
1637 tabsc1 = *boun + dinf * (1. - centr) / centr;
1642 for (j = 1; j <= 7; ++j) {
1643 absc = hlgth * xgk[j - 1];
1644 absc1 = centr - absc;
1645 absc2 = centr + absc;
1646 tabsc1 = *boun + dinf * (1. - absc1) / absc1;
1647 tabsc2 = *boun + dinf * (1. - absc2) / absc2;
1648 vec[(j << 1) - 1] = tabsc1;
1649 vec[j * 2] = tabsc2;
1651 vec2[(j << 1) - 1] = -tabsc1;
1652 vec2[j * 2] = -tabsc2;
1657 if (*inf == 2) f(vec2, 15, ex);
1659 if (*inf == 2) fval1 += vec2[0];
1660 fc = fval1 / centr / centr;
1667 *resabs = fabs(resk);
1668 for (j = 1; j <= 7; ++j) {
1669 absc = hlgth * xgk[j - 1];
1670 absc1 = centr - absc;
1671 absc2 = centr + absc;
1672 tabsc1 = *boun + dinf * (1. - absc1) / absc1;
1673 tabsc2 = *boun + dinf * (1. - absc2) / absc2;
1674 fval1 = vec[(j << 1) - 1];
1677 fval1 += vec2[(j << 1) - 1];
1680 fval2 += vec2[j * 2];
1682 fval1 = fval1 / absc1 / absc1;
1683 fval2 = fval2 / absc2 / absc2;
1686 fsum = fval1 + fval2;
1687 resg += wg[j - 1] * fsum;
1688 resk += wgk[j - 1] * fsum;
1689 *resabs += wgk[j - 1] * (fabs(fval1) + fabs(fval2));
1693 *resasc = wgk[7] * fabs(fc - reskh);
1694 for (j = 1; j <= 7; ++j) {
1695 *resasc += wgk[j - 1] * (fabs(fv1[j - 1] - reskh) +
1696 fabs(fv2[j - 1] - reskh));
1699 *result = resk * hlgth;
1702 *abserr = fabs((resk - resg) * hlgth);
1703 if (*resasc != 0. && *abserr != 0.) {
1704 *abserr = *resasc * fmin2(1., pow(*abserr * 200. / *resasc, 1.5));
1706 if (*resabs > uflow / (epmach * 50.)) {
1707 *abserr = fmax2(epmach * 50. * *resabs, *abserr);
1712 template<
class Float>
static void rdqelg(
int *n, Float *epstab, Float *
1713 result, Float *abserr, Float *res3la,
int *nres)
1716 int i__, indx, ib, ib2, ie, k1, k2, k3, num, newelm, limexp;
1717 Float delta1, delta2, delta3, e0, e1, e1abs, e2, e3, epmach, epsinf;
1718 Float oflow, ss, res;
1719 Float errA, err1, err2, err3, tol1, tol2, tol3;
1799 epmach = DBL_EPSILON;
1803 *result = epstab[*n];
1808 epstab[*n + 2] = epstab[*n];
1809 newelm = (*n - 1) / 2;
1813 for (i__ = 1; i__ <= newelm; ++i__) {
1816 res = epstab[k1 + 2];
1822 err2 = fabs(delta2);
1823 tol2 = fmax2(fabs(e2), e1abs) * epmach;
1825 err3 = fabs(delta3);
1826 tol3 = fmax2(e1abs, fabs(e0)) * epmach;
1827 if (err2 <= tol2 && err3 <= tol3) {
1831 *abserr = err2 + err3;
1839 err1 = fabs(delta1);
1840 tol1 = fmax2(e1abs, fabs(e3)) * epmach;
1845 if (err1 > tol1 && err2 > tol2 && err3 > tol3) {
1846 ss = 1. / delta1 + 1. / delta2 - 1. / delta3;
1847 epsinf = fabs(ss * e1);
1852 if (epsinf > 1e-4) {
1866 errA = err2 + fabs(res - e2) + err3;
1867 if (errA <= *abserr) {
1877 *n = (limexp / 2 << 1) - 1;
1880 if (num / 2 << 1 == num) ib = 2;
else ib = 1;
1882 for (i__ = 1; i__ <= ie; ++i__) {
1884 epstab[ib] = epstab[ib2];
1888 indx = num - *n + 1;
1889 for (i__ = 1; i__ <= *n; ++i__) {
1890 epstab[i__] = epstab[indx];
1897 *abserr = fabs(*result - res3la[3]) +
1898 fabs(*result - res3la[2]) +
1899 fabs(*result - res3la[1]);
1900 res3la[1] = res3la[2];
1901 res3la[2] = res3la[3];
1902 res3la[3] = *result;
1904 res3la[*nres] = *result;
1909 *abserr = fmax2(*abserr, epmach * 5. * fabs(*result));
1913 template<
class Float,
class integr_fn>
static void rdqk21(integr_fn f,
void *ex, Float *a, Float *b, Float *result,
1914 Float *abserr, Float *resabs, Float *resasc)
1918 static double wg[5] = { .066671344308688137593568809893332,
1919 .149451349150580593145776339657697,
1920 .219086362515982043995534934228163,
1921 .269266719309996355091226921569469,
1922 .295524224714752870173892994651338 };
1923 static double xgk[11] = { .995657163025808080735527280689003,
1924 .973906528517171720077964012084452,
1925 .930157491355708226001207180059508,
1926 .865063366688984510732096688423493,
1927 .780817726586416897063717578345042,
1928 .679409568299024406234327365114874,
1929 .562757134668604683339000099272694,
1930 .433395394129247190799265943165784,
1931 .294392862701460198131126603103866,
1932 .14887433898163121088482600112972,0. };
1933 static double wgk[11] = { .011694638867371874278064396062192,
1934 .03255816230796472747881897245939,
1935 .05475589657435199603138130024458,
1936 .07503967481091995276704314091619,
1937 .093125454583697605535065465083366,
1938 .109387158802297641899210590325805,
1939 .123491976262065851077958109831074,
1940 .134709217311473325928054001771707,
1941 .142775938577060080797094273138717,
1942 .147739104901338491374841515972068,
1943 .149445554002916905664936468389821 };
1947 Float fv1[10], fv2[10], vec[21];
1948 Float absc, resg, resk, fsum, fval1, fval2;
1949 Float hlgth, centr, reskh, uflow;
1950 Float fc, epmach, dhlgth;
2048 epmach = DBL_EPSILON;
2051 centr = (*a + *b) * .5;
2052 hlgth = (*b - *a) * .5;
2053 dhlgth = fabs(hlgth);
2060 for (j = 1; j <= 5; ++j) {
2062 absc = hlgth * xgk[jtw - 1];
2063 vec[(j << 1) - 1] = centr - absc;
2065 vec[j * 2] = centr + absc;
2067 for (j = 1; j <= 5; ++j) {
2068 jtwm1 = (j << 1) - 1;
2069 absc = hlgth * xgk[jtwm1 - 1];
2070 vec[(j << 1) + 9] = centr - absc;
2071 vec[(j << 1) + 10] = centr + absc;
2075 resk = wgk[10] * fc;
2076 *resabs = fabs(resk);
2077 for (j = 1; j <= 5; ++j) {
2079 absc = hlgth * xgk[jtw - 1];
2080 fval1 = vec[(j << 1) - 1];
2082 fv1[jtw - 1] = fval1;
2083 fv2[jtw - 1] = fval2;
2084 fsum = fval1 + fval2;
2085 resg += wg[j - 1] * fsum;
2086 resk += wgk[jtw - 1] * fsum;
2087 *resabs += wgk[jtw - 1] * (fabs(fval1) + fabs(fval2));
2090 for (j = 1; j <= 5; ++j) {
2091 jtwm1 = (j << 1) - 1;
2092 absc = hlgth * xgk[jtwm1 - 1];
2093 fval1 = vec[(j << 1) + 9];
2094 fval2 = vec[(j << 1) + 10];
2095 fv1[jtwm1 - 1] = fval1;
2096 fv2[jtwm1 - 1] = fval2;
2097 fsum = fval1 + fval2;
2098 resk += wgk[jtwm1 - 1] * fsum;
2099 *resabs += wgk[jtwm1 - 1] * (fabs(fval1) + fabs(fval2));
2103 *resasc = wgk[10] * fabs(fc - reskh);
2104 for (j = 1; j <= 10; ++j) {
2105 *resasc += wgk[j - 1] * (fabs(fv1[j - 1] - reskh) +
2106 fabs(fv2[j - 1] - reskh));
2109 *result = resk * hlgth;
2112 *abserr = fabs((resk - resg) * hlgth);
2113 if (*resasc != 0. && *abserr != 0.) {
2114 *abserr = *resasc * fmin2(1., pow(*abserr * 200. / *resasc, 1.5));
2116 if (*resabs > uflow / (epmach * 50.)) {
2117 *abserr = fmax2(epmach * 50. * *resabs, *abserr);
2122 template<
class Float>
static void rdqpsrt(
int *limit,
int *last,
int *maxerr,
2123 Float *ermax, Float *elist,
int *iord,
int *nrmax)
2126 int i, j, k, ido, jbnd, isucc, jupbn;
2127 Float errmin, errmax;
2202 errmax = elist[*maxerr];
2205 for (i = 1; i <= ido; ++i) {
2206 isucc = iord[*nrmax - 1];
2207 if (errmax <= elist[isucc])
2209 iord[*nrmax] = isucc;
2218 if (*last > *limit / 2 + 2)
2219 jupbn = *limit + 3 - *last;
2223 errmin = elist[*last];
2229 for (i = *nrmax + 1; i <= jbnd; ++i) {
2231 if (errmax >= elist[isucc]) {
2233 iord[i - 1] = *maxerr;
2234 for (j = i, k = jbnd; j <= jbnd; j++, k--) {
2236 if (errmin < elist[isucc]) {
2238 iord[k + 1] = *last;
2241 iord[k + 1] = isucc;
2246 iord[i - 1] = isucc;
2249 iord[jbnd] = *maxerr;
2250 iord[jupbn] = *last;
2254 *maxerr = iord[*nrmax];
2255 *ermax = elist[*maxerr];