55 #ifndef O2SCL_CUBATURE_H 56 #define O2SCL_CUBATURE_H 60 #include <boost/numeric/ublas/vector.hpp> 62 #ifndef O2SCL_CLENCURT_H 63 #define O2SCL_CLENCURT_H 64 #include <o2scl/clencurt.h> 66 #include <o2scl/err_hnd.h> 67 #include <o2scl/vector.h> 132 template<
class func_t>
138 typedef boost::numeric::ublas::vector_range<ubvector> ubvector_range;
168 double errMax(
const std::vector<esterr> &ee) {
170 for (
size_t k = 0; k < ee.size(); ++k) {
171 if (ee[k].err > errmax) errmax = ee[k].err;
208 for (
size_t i = 0; i < h.
dim; ++i) {
216 template<
class vec_t>
221 h.
data.resize(dim*2);
223 for (
size_t i = 0; i < dim; ++i) {
224 h.
data[i] = center[i];
225 h.
data[i + dim] = halfwidth[i];
227 h.
vol = compute_vol(h);
237 h.
data.resize(dim*2);
239 for (
size_t i = 0; i < dim; ++i) {
241 h.
data[i + dim] = dat[i+dim];
243 h.
vol = compute_vol(h);
249 template<
class vec_t>
250 void make_hypercube_range
251 (
size_t dim,
const vec_t &xmin,
const vec_t &xmax,
hypercube &h) {
253 make_hypercube(dim,xmin,xmax,h);
254 for (
size_t i = 0; i < dim; ++i) {
255 h.
data[i] = 0.5 * (xmin[i] + xmax[i]);
256 h.
data[i + dim] = 0.5 * (xmax[i] - xmin[i]);
258 h.
vol = compute_vol(h);
302 std::vector<esterr>
ee;
326 R.
h.
data[d + dim] *= 0.5;
328 make_hypercube2(dim, R.
h.
data,R2.
h);
392 r.
pts.resize((num_regions
424 if (rule15gauss_evalError(r, R[0].fdim, f, nR, R)) {
428 if (rule75genzmalik_evalError(r, R[0].fdim, f, nR, R)) {
432 for (iR = 0; iR < nR; ++iR) {
433 R[iR].errmax = errMax(R[iR].ee);
447 #if defined(__GNUC__) && \ 448 ((__GNUC__ == 3 && __GNUC_MINOR__ >= 4) || __GNUC__ > 3) 449 return __builtin_ctz(~n);
451 const size_t bits[256] = {
452 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 4,
453 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 5,
454 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 4,
455 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 6,
456 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 4,
457 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 5,
458 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 4,
459 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 7,
460 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 4,
461 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 5,
462 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 4,
463 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 6,
464 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 4,
465 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 5,
466 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 4,
467 0, 1, 0, 2, 0, 1, 0, 3, 0, 1, 0, 2, 0, 1, 0, 8,
470 while ((n & 0xff) == 0xff) {
474 return bit + bits[n & 0xff];
486 std::vector<double> &p,
size_t p_ix,
487 const std::vector<double> &c,
size_t c_ix,
488 const std::vector<double> &r,
size_t r_ix) {
496 for (i = 0; i < dim; ++i) {
497 p[p_ix+i] = c[c_ix+i] + r[r_ix+i];
504 for(
size_t k=0;k<dim;k++) pts[pts_ix+k]=p[p_ix+k];
515 p[p_ix+d] = (signs & mask) ? c[c_ix+d] - r[r_ix+d] :
516 c[c_ix+d] + r[r_ix+d];
524 std::vector<double> &p,
size_t p_ix,
525 const std::vector<double> &c,
size_t c_ix,
526 const std::vector<double> &r,
size_t r_ix) {
528 for (
size_t i = 0; i < dim - 1; ++i) {
529 p[p_ix+i] = c[c_ix+i] - r[r_ix+i];
530 for (
size_t j = i + 1; j < dim; ++j) {
531 p[p_ix+j] = c[c_ix+j] - r[r_ix+j];
532 for(
size_t k=0;k<dim;k++) pts[pts_ix+k]=p[p_ix+k];
534 p[p_ix+i] = c[c_ix+i] + r[r_ix+i];
535 for(
size_t k=0;k<dim;k++) pts[pts_ix+k]=p[p_ix+k];
537 p[p_ix+j] = c[c_ix+j] + r[r_ix+j];
538 for(
size_t k=0;k<dim;k++) pts[pts_ix+k]=p[p_ix+k];
540 p[p_ix+i] = c[c_ix+i] - r[r_ix+i];
541 for(
size_t k=0;k<dim;k++) pts[pts_ix+k]=p[p_ix+k];
545 p[p_ix+j] = c[c_ix+j];
548 p[p_ix+i] = c[c_ix+i];
556 (ubvector &pts,
size_t pts_ix,
size_t dim,
557 std::vector<double> &p,
size_t p_ix,
558 const std::vector<double> &c,
size_t c_ix,
559 const std::vector<double> &r1,
size_t r1_ix,
560 const std::vector<double> &r2,
size_t r2_ix) {
562 for(
size_t k=0;k<dim;k++) pts[pts_ix+k]=p[p_ix+k];
565 for (
size_t i = 0; i < dim; i++) {
566 p[p_ix+i] = c[c_ix+i] - r1[r1_ix+i];
567 for(
size_t k=0;k<dim;k++) pts[pts_ix+k]=p[p_ix+k];
570 p[p_ix+i] = c[c_ix+i] + r1[r1_ix+i];
571 for(
size_t k=0;k<dim;k++) pts[pts_ix+k]=p[p_ix+k];
574 p[p_ix+i] = c[c_ix+i] - r2[r2_ix+i];
575 for(
size_t k=0;k<dim;k++) pts[pts_ix+k]=p[p_ix+k];
578 p[p_ix+i] = c[c_ix+i] + r2[r2_ix+i];
579 for(
size_t k=0;k<dim;k++) pts[pts_ix+k]=p[p_ix+k];
582 p[p_ix+i] = c[c_ix+i];
644 std::vector<double>
p;
658 #ifdef O2SCL_NEVER_DEFINED 664 int rule75genzmalik_evalError
665 (
rule &runder,
size_t fdim, func_t &f,
size_t nR,
666 std::vector<region> &R) {
669 const double lambda2 = 0.3585685828003180919906451539079374954541;
670 const double lambda4 = 0.9486832980505137995996680633298155601160;
671 const double lambda5 = 0.6882472016116852977216287342936235251269;
672 const double weight2 = 980. / 6561.;
673 const double weight4 = 200. / 19683.;
674 const double weightE2 = 245. / 486.;
675 const double weightE4 = 25. / 729.;
676 const double ratio = (lambda2 * lambda2) / (lambda4 * lambda4);
679 size_t i, j, iR, dim = runder.
dim;
682 ubvector &pts2=runder.
pts;
684 alloc_rule_pts(runder, nR);
687 for (iR = 0; iR < nR; ++iR) {
688 std::vector<double> ¢er2=R[iR].h.data;
690 for (i = 0; i < dim; ++i) {
691 r->
p[i] = center2[i];
694 for (i = 0; i < dim; ++i) {
695 r->
p[i+2*dim] = center2[i+dim] * lambda2;
697 for (i = 0; i < dim; ++i) {
698 r->
p[i+dim] = center2[i+dim] * lambda4;
703 evalR0_0fs4d2(runder.
pts,npts*dim, dim, r->
p,0,center2,0,
704 r->
p,2*dim,r->
p,dim);
705 npts += num0_0(dim) + 2 * numR0_0fs(dim);
708 evalRR0_0fs2(runder.
pts,npts*dim,dim,r->
p,0,center2,0,r->
p,dim);
709 npts += numRR0_0fs(dim);
712 for (i = 0; i < dim; ++i) {
713 r->
p[i+dim] = center2[i+dim] * lambda5;
715 evalR_Rfs2(runder.
pts,npts*dim,dim,r->
p,0,center2,0,r->
p,dim);
716 npts += numR_Rfs(dim);
720 if (f(dim, npts, &(pts2[0]), fdim, vals)) {
728 for (i = 0; i < dim * nR; ++i) pts2[i] = 0;
730 for (j = 0; j < fdim; ++j) {
732 const double *v = vals + j;
734 for (iR = 0; iR < nR; ++iR) {
735 double result, res5th;
736 double val0, sum2=0, sum3=0, sum4=0, sum5=0;
746 for (k = 0; k < dim; ++k) {
747 double v0 = v[fdim*(k0 + 4*k)];
748 double v1 = v[fdim*((k0 + 4*k) + 1)];
749 double v2 = v[fdim*((k0 + 4*k) + 2)];
750 double v3 = v[fdim*((k0 + 4*k) + 3)];
755 pts2[iR * dim + k] +=
756 fabs(v0 + v1 - 2*val0 - ratio * (v2 + v3 - 2*val0));
758 #ifdef O2SCL_NEVER_DEFINED 763 for (k = 0; k < numRR0_0fs(dim); ++k) {
764 sum4 += v[fdim*(k0 + k)];
768 for (k = 0; k < numR_Rfs(dim); ++k) {
769 sum5 += v[fdim*(k0 + k)];
773 result = R[iR].h.vol * (r->
weight1 * val0 + weight2 * sum2 +
774 r->
weight3 * sum3 + weight4 * sum4 +
776 res5th = R[iR].h.vol * (r->
weightE1 * val0 + weightE2 * sum2 +
777 r->
weightE3 * sum3 + weightE4 * sum4);
779 R[iR].ee[j].val = result;
780 R[iR].ee[j].err = fabs(res5th - result);
787 for (iR = 0; iR < nR; ++iR) {
789 size_t dimDiffMax = 0;
791 for (i = 0; i < dim; ++i) {
792 if (pts2[iR*dim + i] > maxdiff) {
793 maxdiff = pts2[iR*dim + i];
797 R[iR].splitDim = dimDiffMax;
802 #ifdef O2SCL_NEVER_DEFINED 812 O2SCL_ERR(
"this rule does not support 1d integrals",
821 if (dim >=
sizeof(
size_t) * 8) {
822 O2SCL_ERR(
"this rule does not support large dims",
827 make_rule(dim, fdim,num0_0(dim) + 2 * numR0_0fs(dim)
828 + numRR0_0fs(dim) + numR_Rfs(dim),r);
830 r.
weight1=(12824.0-9120.0*dim+400.0*dim*dim)/19683.0;
831 r.
weight3=(1820.0-400.0*dim)/19683.0;
832 r.
weight5=6859.0/19683.0/((double)(1U << dim));
833 r.
weightE1=(729.0-950.0*dim+50.0*dim*dim)/729.0;
834 r.
weightE3=(265.0-100.0*dim)/1458.0;
844 int rule15gauss_evalError
845 (
rule &r,
size_t fdim, func_t &f,
size_t nR, std::vector<region> &R) {
847 static const double cub_dbl_min=std::numeric_limits<double>::min();
848 static const double cub_dbl_eps=std::numeric_limits<double>::epsilon();
854 const double xgk[8] = {
855 0.991455371120812639206854697526329,
856 0.949107912342758524526189684047851,
857 0.864864423359769072789712788640926,
858 0.741531185599394439863864773280788,
859 0.586087235467691130294144838258730,
860 0.405845151377397166906606412076961,
861 0.207784955007898467600689403773245,
862 0.000000000000000000000000000000000
866 static const double wg[4] = {
867 0.129484966168869693270611432679082,
868 0.279705391489276667901467771423780,
869 0.381830050505118944950369775488975,
870 0.417959183673469387755102040816327
872 static const double wgk[8] = {
873 0.022935322010529224963732008058970,
874 0.063092092629978553290700663189204,
875 0.104790010322250183839876322541518,
876 0.140653259715525918745189590510238,
877 0.169004726639267902826583426598550,
878 0.190350578064785409913256402421014,
879 0.204432940075298892414161999234649,
880 0.209482141084727828012999174891714
886 alloc_rule_pts(r, nR);
890 for (iR = 0; iR < nR; ++iR) {
891 const double center = R[iR].h.data[0];
892 const double halfwidth = R[iR].h.data[1];
894 pts[npts++] = center;
896 for (j = 0; j < (n - 1) / 2; ++j) {
898 double w = halfwidth * xgk[j2];
899 pts[npts++] = center - w;
900 pts[npts++] = center + w;
902 for (j = 0; j < n/2; ++j) {
904 double w = halfwidth * xgk[j2];
905 pts[npts++] = center - w;
906 pts[npts++] = center + w;
912 if (f(1, npts, pts, fdim, vals)) {
916 for (k = 0; k < fdim; ++k) {
917 const double *vk = vals + k;
918 for (iR = 0; iR < nR; ++iR) {
919 const double halfwidth = R[iR].h.data[1];
920 double result_gauss = vk[0] * wg[n/2 - 1];
921 double result_kronrod = vk[0] * wgk[n - 1];
922 double result_abs = fabs(result_kronrod);
923 double result_asc, mean, err;
927 for (j = 0; j < (n - 1) / 2; ++j) {
929 double v = vk[fdim*npts] + vk[fdim*npts+fdim];
930 result_gauss += wg[j] * v;
931 result_kronrod += wgk[j2] * v;
932 result_abs += wgk[j2] * (fabs(vk[fdim*npts])
933 + fabs(vk[fdim*npts+fdim]));
936 for (j = 0; j < n/2; ++j) {
938 result_kronrod += wgk[j2] * (vk[fdim*npts]
939 + vk[fdim*npts+fdim]);
940 result_abs += wgk[j2] * (fabs(vk[fdim*npts])
941 + fabs(vk[fdim*npts+fdim]));
946 R[iR].ee[k].val = result_kronrod * halfwidth;
952 mean = result_kronrod * 0.5;
953 result_asc = wgk[n - 1] * fabs(vk[0] - mean);
955 for (j = 0; j < (n - 1) / 2; ++j) {
957 result_asc += wgk[j2] * (fabs(vk[fdim*npts]-mean)
958 + fabs(vk[fdim*npts+fdim]-mean));
961 for (j = 0; j < n/2; ++j) {
963 result_asc += wgk[j2] * (fabs(vk[fdim*npts]-mean)
964 + fabs(vk[fdim*npts+fdim]-mean));
967 err = fabs(result_kronrod - result_gauss) * halfwidth;
968 result_abs *= halfwidth;
969 result_asc *= halfwidth;
970 if (result_asc != 0 && err != 0) {
971 double scale = pow((200 * err / result_asc), 1.5);
972 err = (scale < 1) ? result_asc * scale : result_asc;
974 if (result_abs > cub_dbl_min / (50 * cub_dbl_eps)) {
975 double min_err = 50 * cub_dbl_eps * result_abs;
976 if (min_err > err) err = min_err;
978 R[iR].ee[k].err = err;
995 return make_rule(dim,fdim,15,r);
1051 h.
items.resize(nalloc);
1065 for (
size_t i = 0; i < fdim; ++i) {
1069 heap_resize(h, nalloc);
1089 size_t fdim = h.
fdim;
1091 for (
size_t i = 0; i < fdim; ++i) {
1092 h.
ee[i].val += hi.
ee[i].val;
1093 h.
ee[i].err += hi.
ee[i].err;
1097 heap_resize(h, h.
n * 2);
1101 int parent = (insert - 1) / 2;
1108 h.
items[insert] = hi;
1115 for (
size_t i = 0; i < ni; ++i) {
1129 O2SCL_ERR(
"Attempted to pop an empty heap in cubature.",
1136 while ((child = i * 2 + 1) < n) {
1141 if (h.
items[child].errmax <= h.
items[i].errmax) {
1147 if (++child < n && h.
items[largest].errmax <
1148 h.
items[child].errmax) {
1156 h.
items[i = largest] = swap;
1160 size_t i, fdim = h.
fdim;
1161 for (i = 0; i < fdim; ++i) {
1162 h.
ee[i].val -= ret.
ee[i].val;
1163 h.
ee[i].err -= ret.
ee[i].err;
1173 double reqAbsError,
double reqRelError,
1182 for (j = 0; j < fdim; ++j) {
1183 if (ee[j].err > reqAbsError && ee[j].err >
1184 fabs(ee[j].val)*reqRelError) {
1192 for (j = 0; j+1 < fdim; j += 2) {
1193 double maxerr, serr, err, maxval, sval, val;
1195 maxerr = ee[j].err > ee[j+1].err ? ee[j].err : ee[j+1].err;
1196 maxval = ee[j].val > ee[j+1].val ? ee[j].val : ee[j+1].val;
1197 serr = maxerr > 0 ? 1/maxerr : 1;
1198 sval = maxval > 0 ? 1/maxval : 1;
1199 err=std::hypot(ee[j].err*serr,ee[j+1].err*serr)*maxerr;
1200 val=std::hypot(ee[j].val*sval,ee[j+1].val*sval)*maxval;
1201 if (err > reqAbsError && err > val*reqRelError) {
1208 if (ee[j].err > reqAbsError && ee[j].err >
1209 fabs(ee[j].val)*reqRelError) {
1218 double err = 0, val = 0;
1219 for (j = 0; j < fdim; ++j) {
1221 val += fabs(ee[j].val);
1223 return err <= reqAbsError || err <= val*reqRelError;
1228 double err = 0, val = 0;
1229 for (j = 0; j < fdim; ++j) {
1230 double absval = fabs(ee[j].val);
1231 if (ee[j].err > err) err = ee[j].err;
1232 if (absval > val) val = absval;
1234 return err <= reqAbsError || err <= val*reqRelError;
1239 double maxerr = 0, maxval = 0, serr, sval, err = 0, val = 0;
1241 for (j = 0; j < fdim; ++j) {
1242 double absval = fabs(ee[j].val);
1243 if (ee[j].err > maxerr) maxerr = ee[j].err;
1244 if (absval > maxval) maxval = absval;
1246 serr = maxerr > 0 ? 1/maxerr : 1;
1247 sval = maxval > 0 ? 1/maxval : 1;
1248 for (j = 0; j < fdim; ++j) {
1249 err += (ee[j].err * serr)*(ee[j].err * serr);
1250 val += fabs(ee[j].val) * sval*fabs(ee[j].val) * sval;
1252 err = sqrt(err) * maxerr;
1253 val = sqrt(val) * maxval;
1254 return err <= reqAbsError || err <= val*reqRelError;
1263 template<
class vec_t>
1266 double reqAbsError,
double reqRelError,
1268 vec_t &err,
int parallel) {
1274 std::vector<region> R;
1275 size_t nR_alloc = 0;
1276 std::vector<esterr> ee(fdim);
1283 regions = heap_alloc(1, fdim);
1287 make_region(h, fdim, R[0]);
1288 if (eval_regions(1, R, f, r) || heap_push(regions, R[0])) {
1295 while (numEval < maxEval || !maxEval) {
1297 if (converged(fdim, regions.
ee, reqAbsError, reqRelError, norm)) {
1330 for (j = 0; j < fdim; ++j) ee[j] = regions.
ee[j];
1333 if (nR + 2 > nR_alloc) {
1334 nR_alloc = (nR + 2) * 2;
1337 R[nR] = heap_pop(regions);
1338 for (j = 0; j < fdim; ++j) ee[j].err -= R[nR].ee[j].err;
1341 if (converged(fdim, ee, reqAbsError, reqRelError, norm)) {
1345 }
while (regions.
n > 0 && (numEval < maxEval || !maxEval));
1347 if (eval_regions(nR, R, f, r)
1348 || heap_push_many(regions, nR, &(R[0]))) {
1360 R[0] = heap_pop(regions);
1361 if (cut_region(R[0], R[1]) || eval_regions(2, R, f, r)
1362 || heap_push_many(regions, 2, &(R[0]))) {
1372 for (j = 0; j < fdim; ++j) {
1376 for (i = 0; i < regions.n; ++i) {
1377 for (j = 0; j < fdim; ++j) {
1378 val[j] += regions.items[i].ee[j].val;
1379 err[j] += regions.items[i].ee[j].err;
1382 destroy_hypercube(regions.items[i].h);
1383 regions.items[i].ee.clear();
1395 template<
class vec_t>
1396 int cubature(
size_t fdim, func_t &f,
size_t dim,
const vec_t &xmin,
1397 const vec_t &xmax,
size_t maxEval,
double reqAbsError,
1398 double reqRelError,
error_norm norm, vec_t &val,
1399 vec_t &err,
int parallel) {
1409 if (f(0, 1, &(xmin[0]), fdim, &(val[0]))) {
1412 for (
size_t i = 0; i < fdim; ++i) err[i] = 0;
1418 make_rule15gauss(dim,fdim,r);
1419 make_hypercube_range(dim,xmin,xmax,h);
1420 status = rulecubature(r,fdim,f,h,maxEval,reqAbsError,
1421 reqRelError,norm,val,err,parallel);
1424 make_rule75genzmalik(dim,fdim,r);
1425 make_hypercube_range(dim,xmin,xmax,h);
1426 status = rulecubature(r,fdim,f,h,maxEval,reqAbsError,
1427 reqRelError,norm,val,err,parallel);
1429 destroy_hypercube(h);
1445 template<
class vec_t>
1446 int integ(
size_t fdim, func_t &f,
size_t dim,
1447 const vec_t &xmin,
const vec_t &xmax,
size_t maxEval,
1448 double reqAbsError,
double reqRelError,
error_norm norm,
1449 vec_t &val, vec_t &err) {
1455 return cubature(fdim,f,dim,xmin,xmax,maxEval,reqAbsError,
1456 reqRelError,norm,val,err,use_parallel);
1462 #ifdef O2SCL_NEVER_DEFINED 1481 template<
class func_t,
class vec_t,
class vec_crange_t,
class vec_range_t>
1488 static const size_t MAXDIM=20;
1523 std::vector<size_t>
m;
1537 size_t &vali,
size_t fdim, func_t &f,
size_t dim,
1538 size_t id, std::vector<double> &p,
const vec_t &xmin,
1539 const vec_t &xmax, vec_t &buf,
size_t nbuf,
1544 for(
size_t k=0;k<dim;k++) {
1545 buf[k+ibuf*dim]=p[k];
1550 (((
const vec_t &)buf),0,buf.size());
1553 if (f(dim, nbuf, &(buf[0]), fdim, &(val[0]) + vali)) {
1556 vali += ibuf * fdim;
1562 double c = (xmin[id] + xmax[id]) * 0.5;
1563 double r = (xmax[id] - xmin[id]) * 0.5;
1564 const double *x = clencurt_x
1565 + ((
id == mi) ? (m[
id] ? (1 << (m[
id] - 1)) : 0) : 0);
1567 size_t nx = (
id == mi ? (m[id] ? (1 << (m[id] - 1)) : 1)
1571 if (compute_cacheval(m, mi, val, vali, fdim, f,
1573 xmin, xmax, buf, nbuf, ibuf)) {
1577 for (i = 0; i < nx; ++i) {
1578 p[id] = c + r * x[i];
1579 if (compute_cacheval(m, mi, val, vali, fdim, f,
1581 xmin, xmax, buf, nbuf, ibuf)) {
1584 p[id] = c - r * x[i];
1585 if (compute_cacheval(m, mi, val, vali, fdim, f,
1587 xmin, xmax, buf, nbuf, ibuf)) {
1601 for (
size_t i = 0; i < dim; ++i) {
1603 if (m[i+i_shift]==0) {
1606 nval*= (1 << (m[i+i_shift]));
1609 nval *= (1 << (m[i+i_shift] + 1)) + 1;
1618 const std::vector<size_t> &m,
1619 size_t mi,
size_t fdim, func_t &f,
size_t dim,
1620 const vec_t &xmin,
const vec_t &xmax, vec_t &buf,
1623 size_t ic = vc.size();
1624 size_t nval, vali = 0, ibuf = 0;
1625 std::vector<double> p(MAXDIM);
1631 for(
size_t j=0;j<dim;j++) {
1634 nval = fdim * num_cacheval(m,mi,dim,0);
1635 vc[ic].val.resize(nval);
1637 if (compute_cacheval(m,mi,vc[ic].val,vali,fdim,f,
1638 dim,0,p,xmin,xmax,buf,nbuf,ibuf)) {
1645 (((
const vec_t &)buf),0,buf.size());
1647 return f(dim, ibuf, &(buf[0]), fdim, &((vc[ic].val)[vali]));
1661 size_t eval(
const std::vector<size_t> &cm,
size_t cmi,
1662 vec_t &cval,
const std::vector<size_t> &m,
size_t md,
1663 size_t fdim,
size_t dim,
size_t id,
1664 double weight, vec_t &val,
size_t voff2) {
1669 for (i = 0; i < fdim; ++i) {
1670 val[i] += cval[i+voff2] * weight;
1674 }
else if (m[
id] == 0 &&
id == md) {
1677 voff = eval(cm, cmi, cval, m, md, fdim, dim,
id+1, weight*2, val,voff2);
1678 voff += fdim * (1 << cm[id]) * 2
1679 * num_cacheval(cm, cmi - (
id+1), dim - (
id+1),
id+1);
1685 size_t mid = m[id] - (
id == md);
1686 const double *w = clencurt_w + mid + (1 << mid) - 1
1687 + (
id == cmi ? (cm[
id] ? 1 + (1 << (cm[
id]-1)) : 1) : 0);
1688 size_t cnx = (
id == cmi ? (cm[id] ? (1 << (cm[id]-1)) : 1)
1690 size_t nx = cm[id] <= mid ? cnx : (1 << mid);
1693 voff = eval(cm, cmi, cval, m, md, fdim, dim,
id + 1,
1694 weight * w[0], val,voff2);
1697 for (i = 0; i < nx; ++i) {
1698 voff += eval(cm, cmi, cval, m, md, fdim, dim,
id + 1,
1699 weight * w[i], val,voff+voff2);
1700 voff += eval(cm, cmi, cval, m, md, fdim, dim,
id + 1,
1701 weight * w[i], val,voff+voff2);
1704 voff += (cnx - nx) * fdim * 2
1705 * num_cacheval(cm, cmi - (
id+1), dim - (
id+1),
id+1);
1715 void evals(std::vector<cache> &vc,
const std::vector<size_t> &m,
1716 size_t md,
size_t fdim,
size_t dim,
double V, vec_t &val) {
1718 for(
size_t k=0;k<fdim;k++) {
1721 for (
size_t i = 0; i < vc.size(); ++i) {
1722 if (vc[i].mi >= dim ||
1723 vc[i].m[vc[i].mi] + (vc[i].mi == md) <= m[vc[i].mi]) {
1724 eval(vc[i].m, vc[i].mi, vc[i].val,
1725 m, md, fdim, dim, 0, V, val,0);
1739 const std::vector<size_t> &m,
1740 size_t fdim,
size_t dim,
double V,
1741 size_t &mi, vec_t &val, vec_t &err, vec_t &val1) {
1746 evals(vc, m, dim, fdim, dim, V, val);
1751 for(
size_t j=0;j<fdim;j++) {
1755 for (i = 0; i < dim; ++i) {
1757 evals(vc, m, i, fdim, dim, V, val1);
1758 for (j = 0; j < fdim; ++j) {
1759 double e = fabs(val[j] - val1[j]);
1760 if (e > emax) emax = e;
1761 if (e > err[j]) err[j] = e;
1763 if (emax > maxerr) {
1773 template<
class vec2_t>
1774 int converged(
size_t fdim,
const vec2_t &vals,
const vec2_t &errs,
1775 double reqAbsError,
double reqRelError,
error_norm norm) {
1781 for (
size_t j = 0; j < fdim; ++j) {
1782 if (errs[j] > reqAbsError && errs[j] > fabs(vals[j])*reqRelError) {
1793 for (j = 0; j+1 < fdim; j += 2) {
1794 double maxerr, serr, err, maxval, sval, val;
1796 maxerr = errs[j] > errs[j+1] ? errs[j] : errs[j+1];
1797 maxval = vals[j] > vals[j+1] ? vals[j] : vals[j+1];
1798 serr = maxerr > 0 ? 1/maxerr : 1;
1799 sval = maxval > 0 ? 1/maxval : 1;
1800 err=std::hypot(errs[j]*serr,errs[j+1]*serr)*maxerr;
1801 val=std::hypot(vals[j]*sval,vals[j+1]*sval)*maxval;
1802 if (err > reqAbsError && err > val*reqRelError) {
1808 if (errs[j] > reqAbsError && errs[j] > fabs(vals[j])*reqRelError) {
1818 double err = 0, val = 0;
1819 for (
size_t j = 0; j < fdim; ++j) {
1821 val += fabs(vals[j]);
1823 return err <= reqAbsError || err <= val*reqRelError;
1829 double err = 0, val = 0;
1830 for (
size_t j = 0; j < fdim; ++j) {
1831 double absval = fabs(vals[j]);
1832 if (errs[j] > err) err = errs[j];
1833 if (absval > val) val = absval;
1835 return err <= reqAbsError || err <= val*reqRelError;
1841 double maxerr = 0, maxval = 0, serr, sval, err = 0, val = 0;
1843 for (
size_t j = 0; j < fdim; ++j) {
1844 double absval = fabs(vals[j]);
1845 if (errs[j] > maxerr) maxerr = errs[j];
1846 if (absval > maxval) maxval = absval;
1848 serr = maxerr > 0 ? 1/maxerr : 1;
1849 sval = maxval > 0 ? 1/maxval : 1;
1850 for (
size_t j = 0; j < fdim; ++j) {
1851 err += (errs[j] * serr)*(errs[j] * serr);
1852 val += (fabs(vals[j]) * sval)*(fabs(vals[j]) * sval);
1854 err = sqrt(err) * maxerr;
1855 val = sqrt(val) * maxval;
1856 return err <= reqAbsError || err <= val*reqRelError;
1860 O2SCL_ERR(
"Improper value of 'norm' in cubature::converged().",
1880 int integ_v_buf(
size_t fdim, func_t &f,
size_t dim,
const vec_t &xmin,
1881 const vec_t &xmax,
size_t maxEval,
double reqAbsError,
1882 double reqRelError,
error_norm norm, std::vector<size_t> &m,
1883 vec_t &buf,
size_t &nbuf,
size_t max_nbuf, vec_t &val,
1888 size_t numEval = 0, new_nbuf;
1890 std::vector<cache> vc;
1908 for (i = 0; i < fdim; ++i) err[i] = 0;
1912 for (i = 0; i < fdim; ++i) {
1917 for (i = 0; i < dim; ++i) {
1919 V *= (xmax[i] - xmin[i]) * 0.5;
1922 new_nbuf = num_cacheval(m,dim,dim,0);
1924 if (max_nbuf < 1) max_nbuf = 1;
1925 if (new_nbuf > max_nbuf) new_nbuf = max_nbuf;
1926 if (nbuf < new_nbuf) {
1927 buf.resize((nbuf=new_nbuf)*dim);
1931 if (add_cacheval(vc,m, dim, fdim, f, dim, xmin, xmax,
1938 eval_integral(vc,m,fdim,dim,V,mi,val,err,val1);
1940 if (converged(fdim,val,err,reqAbsError, reqRelError, norm) ||
1941 (numEval > maxEval && maxEval)) {
1946 if (m[mi] > clencurt_M) {
1951 new_nbuf = num_cacheval(m, mi, dim,0);
1952 if (new_nbuf > nbuf && nbuf < max_nbuf) {
1954 if (nbuf > max_nbuf) nbuf = max_nbuf;
1955 buf.resize(nbuf*dim);
1958 if (add_cacheval(vc,m, mi, fdim, f,
1963 numEval += new_nbuf;
1971 static const size_t DEFAULT_MAX_NBUF=(1U << 20);
1975 int integ(
size_t fdim, func_t &f,
size_t dim,
1976 const vec_t &xmin,
const vec_t &xmax,
size_t maxEval,
1977 double reqAbsError,
double reqRelError,
error_norm norm,
1978 vec_t &val, vec_t &err) {
1982 std::vector<size_t> m(dim);
1986 ret = integ_v_buf(fdim,f,dim,xmin,xmax,
1987 maxEval,reqAbsError,reqRelError,norm,
1988 m,buf,nbuf,16,val,err);
size_t num_points
The number of evaluation points.
size_t eval(const std::vector< size_t > &cm, size_t cmi, vec_t &cval, const std::vector< size_t > &m, size_t md, size_t fdim, size_t dim, size_t id, double weight, vec_t &val, size_t voff2)
Desc.
Base class for integration routines from the Cubature library.
size_t num0_0(size_t dim)
Desc.
The main O<span style='position: relative; top: 0.3em; font-size: 0.8em'>2</span>scl O$_2$scl names...
double vol
cache volume = product of widths
int rulecubature(rule &r, size_t fdim, func_t &f, const hypercube &h, size_t maxEval, double reqAbsError, double reqRelError, error_norm norm, vec_t &val, vec_t &err, int parallel)
Desc.
int heap_push(heap &h, heap_item &hi)
Desc.
int cubature(size_t fdim, func_t &f, size_t dim, const vec_t &xmin, const vec_t &xmax, size_t maxEval, double reqAbsError, double reqRelError, error_norm norm, vec_t &val, vec_t &err, int parallel)
Desc.
sanity check failed - shouldn't happen
std::vector< double > p
Desc.
invalid argument supplied by user
void make_rule15gauss(size_t dim, size_t fdim, rule &r)
Desc.
void evals(std::vector< cache > &vc, const std::vector< size_t > &m, size_t md, size_t fdim, size_t dim, double V, vec_t &val)
Desc.
void make_rule(size_t dim, size_t fdim, size_t num_points, rule &r)
Desc.
size_t num_regions
The max number of regions evaluated at once.
size_t numR0_0fs(size_t dim)
Desc.
std::vector< heap_item > items
Desc.
int converged(size_t fdim, const std::vector< esterr > &ee, double reqAbsError, double reqRelError, error_norm norm)
Desc.
size_t vals_ix
num_regions * num_points * fdim
void eval_integral(std::vector< cache > &vc, const std::vector< size_t > &m, size_t fdim, size_t dim, double V, size_t &mi, vec_t &val, vec_t &err, vec_t &val1)
Desc.
size_t fdim
The number of functions.
int integ(size_t fdim, func_t &f, size_t dim, const vec_t &xmin, const vec_t &xmax, size_t maxEval, double reqAbsError, double reqRelError, error_norm norm, vec_t &val, vec_t &err)
Desc.
int heap_push_many(heap &h, size_t ni, heap_item *hi)
Desc.
double weight1
dimension-dependent constants
void make_region(const hypercube &h, size_t fdim, region &R)
Desc.
void heap_free(heap &h)
Note that heap_free does not deallocate anything referenced by the items.
int integ_v_buf(size_t fdim, func_t &f, size_t dim, const vec_t &xmin, const vec_t &xmax, size_t maxEval, double reqAbsError, double reqRelError, error_norm norm, std::vector< size_t > &m, vec_t &buf, size_t &nbuf, size_t max_nbuf, vec_t &val, vec_t &err)
Desc.
ubvector pts
points to eval: num_regions * num_points * dim
individual relerr criteria in each component
int compute_cacheval(const std::vector< size_t > &m, size_t mi, vec_t &val, size_t &vali, size_t fdim, func_t &f, size_t dim, size_t id, std::vector< double > &p, const vec_t &xmin, const vec_t &xmax, vec_t &buf, size_t nbuf, size_t &ibuf)
Desc.
dat_t * vector_range(dat_t *v, size_t start, size_t last)
Vector range function for pointers.
void destroy_hypercube(hypercube &h)
Desc.
error_norm
Different ways of measuring the absolute and relative error.
void make_hypercube2(size_t dim, const std::vector< double > &dat, hypercube &h)
Desc.
std::vector< size_t > m
Desc.
Cache of the values for the m[dim] grid.
void evalRR0_0fs2(ubvector &pts, size_t pts_ix, size_t dim, std::vector< double > &p, size_t p_ix, const std::vector< double > &c, size_t c_ix, const std::vector< double > &r, size_t r_ix)
Desc.
int eval_regions(size_t nR, std::vector< region > &R, func_t &f, rule &r)
Desc.
size_t numR_Rfs(size_t dim)
Desc.
void evalR_Rfs2(ubvector &pts, size_t pts_ix, size_t dim, std::vector< double > &p, size_t p_ix, const std::vector< double > &c, size_t c_ix, const std::vector< double > &r, size_t r_ix)
Evaluate the integration points for all points (+/-r,...+/-r)
#define O2SCL_ERR(d, n)
Set an error with message d and code n.
abserr is L_2 norm |e|, and relerr is |e|/|v|
abserr is L_1 norm |e|, and relerr is |e|/|v|
Adaptive multidimensional integration on hyper-rectangles using cubature rules from the Cubature libr...
void make_hypercube(size_t dim, const vec_t ¢er, const vec_t &halfwidth, hypercube &h)
Desc.
double errMax(const std::vector< esterr > &ee)
Return the maximum error from ee.
heap heap_alloc(size_t nalloc, size_t fdim)
Desc.
int integ(size_t fdim, func_t &f, size_t dim, const vec_t &xmin, const vec_t &xmax, size_t maxEval, double reqAbsError, double reqRelError, error_norm norm, vec_t &val, vec_t &err)
Desc.
Integration by p-adaptive cubature from the Cubature library.
double errmax
max ee[k].err
abserr is norm |e|, and relerr is |e|/|v|
heap_item heap_pop(heap &h)
Desc.
std::vector< double > data
length 2*dim = center followed by half-widths
std::vector< esterr > ee
array of length fdim
void alloc_rule_pts(rule &r, size_t num_regions)
Desc.
size_t fdim
dimensionality of vector integrand
size_t numRR0_0fs(size_t dim)
Desc.
size_t num_cacheval(const std::vector< size_t > &m, size_t mi, size_t dim, size_t i_shift)
Desc.
size_t ls0(size_t n)
Functions to loop over points in a hypercube.
int add_cacheval(std::vector< cache > &vc, const std::vector< size_t > &m, size_t mi, size_t fdim, func_t &f, size_t dim, const vec_t &xmin, const vec_t &xmax, vec_t &buf, size_t nbuf)
Desc.
int converged(size_t fdim, const vec2_t &vals, const vec2_t &errs, double reqAbsError, double reqRelError, error_norm norm)
Desc.
double compute_vol(const hypercube &h)
Desc.
void make_rule75genzmalik(size_t dim, size_t fdim, rule75genzmalik &r)
Desc.
const dat_t * const_vector_range(const dat_t *v, size_t start, size_t last)
Vector range function for const pointers.
void heap_resize(heap &h, size_t nalloc)
Desc.
paired L2 norms of errors in each component, mainly for integrating vectors of complex numbers ...
size_t dim
The dimensionality.
int cut_region(region &R, region &R2)
Desc.