static int Integrate(This *t, real *integral, real *error, real *prob) { TYPEDEFREGION; typedef struct pool { struct pool *next; Region region[POOLSIZE]; } Pool; count dim, comp, ncur, ipool, npool; int fail; Totals totals[NCOMP]; Pool *cur = NULL, *pool; Region *region; if( VERBOSE > 1 ) { char s[256]; sprintf(s, "Cuhre input parameters:\n" " ndim " COUNT "\n ncomp " COUNT "\n" " epsrel " REAL "\n epsabs " REAL "\n" " flags %d\n mineval " NUMBER "\n maxeval " NUMBER "\n" " key " COUNT, t->ndim, t->ncomp, t->epsrel, t->epsabs, t->flags, t->mineval, t->maxeval, t->key); Print(s); } if( BadComponent(t) ) return -2; if( BadDimension(t) ) return -1; t->epsabs = Max(t->epsabs, NOTZERO); RuleAlloc(t); t->mineval = IMax(t->mineval, t->rule.n + 1); FrameAlloc(t, ShmRm(t)); ForkCores(t); if( (fail = setjmp(t->abort)) ) goto abort; Alloc(cur, 1); cur->next = NULL; ncur = 1; region = cur->region; region->div = 0; for( dim = 0; dim < t->ndim; ++dim ) { Bounds *b = ®ion->bounds[dim]; b->lower = 0; b->upper = 1; } Sample(t, region); for( comp = 0; comp < t->ncomp; ++comp ) { Totals *tot = &totals[comp]; Result *r = ®ion->result[comp]; tot->avg = tot->lastavg = tot->guess = r->avg; tot->err = tot->lasterr = r->err; tot->weightsum = 1/Max(Sq(r->err), NOTZERO); tot->avgsum = tot->weightsum*r->avg; tot->chisq = tot->chisqsum = tot->chisum = 0; } for( t->nregions = 1; ; ++t->nregions ) { count maxcomp, bisectdim; real maxratio, maxerr; Result result[NCOMP]; Region *regionL, *regionR; Bounds *bL, *bR; if( VERBOSE ) { char s[128 + 128*NCOMP], *p = s; p += sprintf(p, "\n" "Iteration " COUNT ": " NUMBER " integrand evaluations so far", t->nregions, t->neval); for( comp = 0; comp < t->ncomp; ++comp ) { cTotals *tot = &totals[comp]; p += sprintf(p, "\n[" COUNT "] " REAL " +- " REAL " \tchisq " REAL " (" COUNT " df)", comp + 1, tot->avg, tot->err, tot->chisq, t->nregions - 1); } Print(s); } maxratio = -INFTY; maxcomp = 0; for( comp = 0; comp < t->ncomp; ++comp ) { creal ratio = totals[comp].err/MaxErr(totals[comp].avg); if( ratio > maxratio ) { maxratio = ratio; maxcomp = comp; } } if( maxratio <= 1 && t->neval >= t->mineval ) break; if( t->neval >= t->maxeval ) { fail = 1; break; } maxerr = -INFTY; regionL = cur->region; npool = ncur; for( pool = cur; pool; npool = POOLSIZE, pool = pool->next ) for( ipool = 0; ipool < npool; ++ipool ) { Region *region = &pool->region[ipool]; creal err = region->result[maxcomp].err; if( err > maxerr ) { maxerr = err; regionL = region; } } if( ncur == POOLSIZE ) { Pool *prev = cur; Alloc(cur, 1); cur->next = prev; ncur = 0; } regionR = &cur->region[ncur++]; regionR->div = ++regionL->div; FCopy(result, regionL->result); XCopy(regionR->bounds, regionL->bounds); bisectdim = result[maxcomp].bisectdim; bL = ®ionL->bounds[bisectdim]; bR = ®ionR->bounds[bisectdim]; bL->upper = bR->lower = .5*(bL->upper + bL->lower); Sample(t, regionL); Sample(t, regionR); for( comp = 0; comp < t->ncomp; ++comp ) { cResult *r = &result[comp]; Result *rL = ®ionL->result[comp]; Result *rR = ®ionR->result[comp]; Totals *tot = &totals[comp]; real diff, err, w, avg, sigsq; tot->lastavg += diff = rL->avg + rR->avg - r->avg; diff = fabs(.25*diff); err = rL->err + rR->err; if( err > 0 ) { creal c = 1 + 2*diff/err; rL->err *= c; rR->err *= c; } rL->err += diff; rR->err += diff; tot->lasterr += rL->err + rR->err - r->err; tot->weightsum += w = 1/Max(Sq(tot->lasterr), NOTZERO); sigsq = 1/tot->weightsum; tot->avgsum += w*tot->lastavg; avg = sigsq*tot->avgsum; tot->chisum += w *= tot->lastavg - tot->guess; tot->chisqsum += w*tot->lastavg; tot->chisq = tot->chisqsum - avg*tot->chisum; if( LAST ) { tot->avg = tot->lastavg; tot->err = tot->lasterr; } else { tot->avg = avg; tot->err = sqrt(sigsq); } } } for( comp = 0; comp < t->ncomp; ++comp ) { cTotals *tot = &totals[comp]; integral[comp] = tot->avg; error[comp] = tot->err; prob[comp] = ChiSquare(tot->chisq, t->nregions - 1); } #ifdef MLVERSION if( REGIONS ) { MLPutFunction(stdlink, "List", 2); MLPutFunction(stdlink, "List", t->nregions); npool = ncur; for( pool = cur; pool; npool = POOLSIZE, pool = pool->next ) for( ipool = 0; ipool < npool; ++ipool ) { Region const *region = &pool->region[ipool]; real lower[NDIM], upper[NDIM]; for( dim = 0; dim < t->ndim; ++dim ) { cBounds *b = ®ion->bounds[dim]; lower[dim] = b->lower; upper[dim] = b->upper; } MLPutFunction(stdlink, "Cuba`Cuhre`region", 3); MLPutRealList(stdlink, lower, t->ndim); MLPutRealList(stdlink, upper, t->ndim); MLPutFunction(stdlink, "List", t->ncomp); for( comp = 0; comp < t->ncomp; ++comp ) { cResult *r = ®ion->result[comp]; real res[] = {r->avg, r->err}; MLPutRealList(stdlink, res, Elements(res)); } } } #endif abort: while( (pool = cur) ) { cur = cur->next; free(pool); } WaitCores(t); FrameFree(t); RuleFree(t); return fail; }
static real LocalSearch(This *t, ccount nfree, ccount *ifree, cBounds *b, creal *x, creal fx, real *z) { Vector(real, y, NDIM); Vector(real, p, NDIM); real delta, smax, sopp, spmax, snmax; real fy, fz, ftest; int sign; count i; /* Choose a direction p along which to move away from the present x. We choose the direction which leads farthest away from all borders. */ smax = INFTY; for( i = 0; i < nfree; ++i ) { ccount dim = ifree[i]; creal sp = b[dim].upper - x[dim]; creal sn = x[dim] - b[dim].lower; if( sp < sn ) { smax = Min(smax, sn); p[i] = -1; } else { smax = Min(smax, sp); p[i] = 1; } } smax *= .9; /* Move along p until the integrand changes appreciably or we come close to a border. */ XCopy(y, x); ftest = SUFTOL*(1 + fabsx(fx)); delta = RTDELTA/5; do { delta = Min(5*delta, smax); for( i = 0; i < nfree; ++i ) { ccount dim = ifree[i]; y[dim] = x[dim] + delta*p[i]; } fy = Sample(t, y); if( fabsx(fy - fx) > ftest ) break; } while( delta != smax ); /* Construct a second direction p' orthogonal to p, i.e. p.p' = 0. We let pairs of coordinates cancel in the dot product, i.e. we choose p'[0] = p[0], p'[1] = -p[1], etc. (It should really be 1/p and -1/p, but p consists of 1's and -1's.) For odd nfree, we let the last three components cancel by choosing p'[nfree - 3] = p[nfree - 3], p'[nfree - 2] = -1/2 p[nfree - 2], and p'[nfree - 1] = -1/2 p[nfree - 1]. */ sign = (nfree <= 1 && fy > fx) ? 1 : -1; spmax = snmax = INFTY; for( i = 0; i < nfree; ++i ) { ccount dim = ifree[i]; real sp, sn; p[i] *= (nfree & 1 && nfree - i <= 2) ? -.5*sign : (sign = -sign); sp = (b[dim].upper - y[dim])/p[i]; sn = (y[dim] - b[dim].lower)/p[i]; if( p[i] > 0 ) { spmax = Min(spmax, sp); snmax = Min(snmax, sn); } else { spmax = Min(spmax, -sn); snmax = Min(snmax, -sp); } } smax = .9*spmax; sopp = .9*snmax; if( nfree > 1 && smax < snmax ) { real tmp = smax; smax = sopp; sopp = tmp; for( i = 0; i < nfree; ++i ) p[i] = -p[i]; } /* Move along p' until the integrand changes appreciably or we come close to a border. */ XCopy(z, y); ftest = SUFTOL*(1 + fabsx(fy)); delta = RTDELTA/5; do { delta = Min(5*delta, smax); for( i = 0; i < nfree; ++i ) { ccount dim = ifree[i]; z[dim] = y[dim] + delta*p[i]; } fz = Sample(t, z); if( fabsx(fz - fy) > ftest ) break; } while( delta != smax ); if( fy != fz ) { real pleneps, grad, range, step; Point low; if( fy > fz ) { grad = (fz - fy)/delta; range = smax/.9; step = Min(delta + delta, smax); } else { grad = (fy - fz)/delta; range = sopp/.9 + delta; step = Min(delta + delta, sopp); XCopy(y, z); fy = fz; for( i = 0; i < nfree; ++i ) p[i] = -p[i]; } pleneps = Length(nfree, p) + RTEPS; low = LineSearch(t, nfree, ifree, p, y, fy, z, step, range, grad, RTEPS/pleneps, 0., RTEPS); fz = low.f; } if( fz != fx ) { real pleneps, grad, range, step; Point low; spmax = snmax = INFTY; for( i = 0; i < nfree; ++i ) { ccount dim = ifree[i]; p[i] = z[dim] - x[dim]; if( p[i] != 0 ) { creal sp = (b[dim].upper - x[dim])/p[i]; creal sn = (x[dim] - b[dim].lower)/p[i]; if( p[i] > 0 ) { spmax = Min(spmax, sp); snmax = Min(snmax, sn); } else { spmax = Min(spmax, -sn); snmax = Min(snmax, -sp); } } } grad = fz - fx; range = spmax; step = Min(.9*spmax, 2.); pleneps = Length(nfree, p) + RTEPS; if( fz > fx ) { delta = Min(.9*snmax, RTDELTA/pleneps); for( i = 0; i < nfree; ++i ) { ccount dim = ifree[i]; z[dim] = x[dim] - delta*p[i]; } fz = Sample(t, z); if( fz < fx ) { grad = (fz - fx)/delta; range = snmax; step = Min(.9*snmax, delta + delta); for( i = 0; i < nfree; ++i ) p[i] = -p[i]; } else if( delta < 1 ) grad = (fx - fz)/delta; } low = LineSearch(t, nfree, ifree, p, x, fx, z, step, range, grad, RTEPS/pleneps, 0., RTEPS); fz = low.f; } return fz; }
static real FindMinimum(This *t, cBounds *b, real *xmin, real fmin) { Vector(real, hessian, NDIM*NDIM); Vector(real, gfree, NDIM); Vector(real, p, NDIM); Vector(real, tmp, NDIM); Vector(count, ifree, NDIM); Vector(count, ifix, NDIM); real ftmp, fini = fmin; ccount maxeval = t->neval + 50*t->ndim; count nfree, nfix; count dim, local; Clear(hessian, t->ndim*t->ndim); for( dim = 0; dim < t->ndim; ++dim ) Hessian(dim, dim) = 1; /* Step 1: - classify the variables as "fixed" (sufficiently close to a border) and "free", - if the integrand is flat in the direction of the gradient w.r.t. the free dimensions, perform a local search. */ for( local = 0; local < 2; ++local ) { bool resample = false; nfree = nfix = 0; for( dim = 0; dim < t->ndim; ++dim ) { if( xmin[dim] < b[dim].lower + (1 + fabsx(b[dim].lower))*QEPS ) { xmin[dim] = b[dim].lower; ifix[nfix++] = dim; resample = true; } else if( xmin[dim] > b[dim].upper - (1 + fabsx(b[dim].upper))*QEPS ) { xmin[dim] = b[dim].upper; ifix[nfix++] = Tag(dim); resample = true; } else ifree[nfree++] = dim; } if( resample ) fini = fmin = Sample(t, xmin); if( nfree == 0 ) goto releasebounds; Gradient(t, nfree, ifree, b, xmin, fmin, gfree); if( local || Length(nfree, gfree) > GTOL ) break; ftmp = LocalSearch(t, nfree, ifree, b, xmin, fmin, tmp); if( ftmp > fmin - (1 + fabsx(fmin))*RTEPS ) goto releasebounds; fmin = ftmp; XCopy(xmin, tmp); } while( t->neval <= maxeval ) { /* Step 2a: perform a quasi-Newton iteration on the free variables only. */ if( nfree > 0 ) { real plen, pleneps; real minstep; count i, mini = 0, minfix = 0; Point low; LinearSolve(t, nfree, hessian, gfree, p); plen = Length(nfree, p); pleneps = plen + RTEPS; minstep = INFTY; for( i = 0; i < nfree; ++i ) { count dim = Untag(ifree[i]); if( fabsx(p[i]) > EPS ) { real step; count fix; if( p[i] < 0 ) { step = (b[dim].lower - xmin[dim])/p[i]; fix = dim; } else { step = (b[dim].upper - xmin[dim])/p[i]; fix = Tag(dim); } if( step < minstep ) { minstep = step; mini = i; minfix = fix; } } } if( minstep*pleneps <= DELTA ) { fixbound: ifix[nfix++] = minfix; if( mini < --nfree ) { creal diag = Hessian(mini, mini); Clear(tmp, mini); for( i = mini; i < nfree; ++i ) tmp[i] = Hessian(i + 1, mini); for( i = mini; i < nfree; ++i ) { Move(&Hessian(i, 0), &Hessian(i + 1, 0), i); Hessian(i, i) = Hessian(i + 1, i + 1); } RenormalizeCholesky(t, nfree, hessian, tmp, diag); Move(&ifree[mini], &ifree[mini + 1], nfree - mini); Move(&gfree[mini], &gfree[mini + 1], nfree - mini); } continue; } low = LineSearch(t, nfree, ifree, p, xmin, fmin, tmp, Min(minstep, 1.), Min(minstep, 100.), Dot(nfree, gfree, p), RTEPS/pleneps, DELTA/pleneps, .2); if( low.dx > 0 ) { real fdiff; fmin = low.f; XCopy(xmin, tmp); Gradient(t, nfree, ifree, b, xmin, fmin, tmp); BFGS(t, nfree, hessian, tmp, gfree, p, low.dx); XCopy(gfree, tmp); if( fabsx(low.dx - minstep) < QEPS*minstep ) goto fixbound; fdiff = fini - fmin; fini = fmin; if( fdiff > (1 + fabsx(fmin))*FTOL || low.dx*plen > (1 + Length(t->ndim, xmin))*FTOL ) continue; } } /* Step 2b: check whether freeing any fixed variable will lead to a reduction in f. */ releasebounds: if( nfix > 0 ) { real mingrad = INFTY; count i, mini = 0; bool repeat = false; Gradient(t, nfix, ifix, b, xmin, fmin, tmp); for( i = 0; i < nfix; ++i ) { creal grad = Sign(ifix[i])*tmp[i]; if( grad < -RTEPS ) { repeat = true; if( grad < mingrad ) { mingrad = grad; mini = i; } } } if( repeat ) { gfree[nfree] = tmp[mini]; ifree[nfree] = Untag(ifix[mini]); Clear(&Hessian(nfree, 0), nfree); Hessian(nfree, nfree) = 1; ++nfree; --nfix; Move(&ifix[mini], &ifix[mini + 1], nfix - mini); continue; } } break; } return fmin; }
static Point LineSearch(This *t, ccount nfree, ccount *ifree, creal *p, creal *xini, real fini, real *x, real step, creal range, creal grad, creal ftol, creal xtol, creal gtol) { real tol = ftol, tol2 = tol + tol; Point cur = {0, fini}; XCopy(x, xini); /* don't even try if a) we'd walk backwards, b) the range to explore is too small, c) the gradient is positive, i.e. we'd move uphill */ if( step > 0 && range > tol2 && grad <= 0 ) { creal eps = RTEPS*fabsx(range) + ftol; creal mingrad = -1e-4*grad, maxgrad = -gtol*grad; real end = range + eps; real maxstep = range - eps/(1 + RTEPS); Point min = cur, v = cur, w = cur; Point a = cur, b = {end, 0}; real a1, b1 = end; /* distmin: distance along p from xini to the minimum, u: second-lowest point, v: third-lowest point, a, b: interval in which the minimum is sought. */ real distmin = 0, dist, mid, q, r, s; count i; int shift; bool first; for( first = true; ; first = false ) { if( step >= maxstep ) { step = maxstep; maxstep = maxstep*(1 + .75*RTEPS) + .75*tol; } cur.dx = (fabsx(step) >= tol) ? step : (step > 0) ? tol : -tol; dist = distmin + cur.dx; for( i = 0; i < nfree; ++i ) { ccount dim = ifree[i]; x[dim] = xini[dim] + dist*p[i]; } cur.f = Sample(t, x); if( cur.f <= min.f ) { v = w; w = min; min.f = cur.f; distmin = dist; /* shift everything to the new minimum position */ maxstep -= cur.dx; v.dx -= cur.dx; w.dx -= cur.dx; a.dx -= cur.dx; b.dx -= cur.dx; if( cur.dx < 0 ) b = w; else a = w; tol = RTEPS*fabsx(distmin) + ftol; tol2 = tol + tol; } else { if( cur.dx < 0 ) a = cur; else b = cur; if( cur.f <= w.f || w.dx == 0 ) v = w, w = cur; else if( cur.f <= v.f || v.dx == 0 || v.dx == w.dx ) v = cur; } if( distmin + b.dx <= xtol ) break; if( min.f < fini && a.f - min.f <= fabsx(a.dx)*maxgrad && (fabsx(distmin - range) > tol || maxstep < b.dx) ) break; mid = .5*(a.dx + b.dx); if( fabsx(mid) <= tol2 - .5*(b.dx - a.dx) ) break; r = q = s = 0; if( fabsx(end) > tol ) { if( first ) { creal s1 = w.dx*grad; creal s2 = w.f - min.f; s = (s1 - ((distmin == 0) ? 0 : 2*s2))*w.dx; q = 2*(s2 - s1); } else { creal s1 = w.dx*(v.f - min.f); creal s2 = v.dx*(w.f - min.f); s = s1*w.dx - s2*v.dx; q = 2*(s2 - s1); } if( q > 0 ) s = -s; q = fabsx(q); r = end; if( step != b1 || b.dx <= maxstep ) end = step; } if( distmin == a.dx ) step = mid; else if( b.dx > maxstep ) step = (step < b.dx) ? -4*a.dx : maxstep; else { real num = a.dx, den = b.dx; if( fabsx(b.dx) <= tol || (w.dx > 0 && fabsx(a.dx) > tol) ) num = b.dx, den = a.dx; num /= -den; step = (num < 1) ? .5*den*sqrtx(num) : 5/11.*den*(.1 + 1/num); } if( step > 0 ) a1 = a.dx, b1 = step; else a1 = step, b1 = b.dx; if( fabsx(s) < fabsx(.5*q*r) && s > q*a1 && s < q*b1 ) { step = s/q; if( step - a.dx < tol2 || b.dx - step < tol2 ) step = (mid > 0) ? tol : -tol; } else end = (mid > 0) ? b.dx : a.dx; } first = true; if( fabsx(distmin - range) < tol ) { distmin = range; if( maxstep > b.dx ) first = false; } for( cur.dx = distmin, cur.f = min.f, shift = -1; ; cur.dx = Max(ldexp(distmin, shift), ftol), shift <<= 1 ) { for( i = 0; i < nfree; ++i ) { ccount dim = ifree[i]; x[dim] = xini[dim] + cur.dx*p[i]; } if( !first ) cur.f = Sample(t, x); if( cur.dx + b.dx <= xtol ) { cur.dx = 0; break; } if( fini - cur.f > cur.dx*mingrad ) break; if( cur.dx <= ftol ) { cur.dx = 0; break; } first = false; } } return cur; }