/* bound (1 + 1/m)^n, m > 0, n >= 0 */ void mag_binpow_uiui(mag_t b, ulong m, ulong n) { mag_t t; if (m == 0) { mag_inf(b); return; } mag_init(t); /* bound by exp(n/m) <= 1 + (n/m) + (n/m)^2 */ if (m > n) { mag_set_ui(t, n); /* x = n/m */ mag_div_ui(t, t, m); mag_mul(b, t, t); /* x^2 */ mag_add(b, b, t); /* x */ mag_one(t); mag_add(b, b, t); /* 1 */ } else { mag_one(b); mag_div_ui(b, b, m); mag_one(t); mag_add(t, t, b); mag_pow_ui(b, t, n); } mag_clear(t); }
void mag_geom_series(mag_t res, const mag_t x, ulong n) { if (mag_is_zero(x)) { if (n == 0) mag_one(res); else mag_zero(res); } else if (mag_is_inf(x)) { mag_inf(res); } else { mag_t t; mag_init(t); mag_one(t); mag_sub_lower(t, t, x); if (mag_is_zero(t)) { mag_inf(res); } else { mag_pow_ui(res, x, n); mag_div(res, res, t); } mag_clear(t); } }
void acb_hypgeom_mag_chi(mag_t chi, ulong n) { mag_t p, q; ulong k; mag_init(p); mag_init(q); if (n % 2 == 0) { mag_one(p); mag_one(q); } else { /* upper bound for pi/2 */ mag_set_ui_2exp_si(p, 843314857, -28); mag_one(q); } for (k = n; k >= 2; k -= 2) { mag_mul_ui(p, p, k); mag_mul_ui_lower(q, q, k - 1); } mag_div(chi, p, q); mag_clear(p); mag_clear(q); }
void arb_sin_cos_pi(arb_t s, arb_t c, const arb_t x, long prec) { arb_t t; arb_t u; fmpz_t v; if (arf_cmpabs_2exp_si(arb_midref(x), FLINT_MAX(65536, (4*prec))) > 0) { arf_zero(arb_midref(s)); mag_one(arb_radref(s)); arf_zero(arb_midref(c)); mag_one(arb_radref(c)); return; } arb_init(t); arb_init(u); fmpz_init(v); arb_mul_2exp_si(t, x, 1); arf_get_fmpz(v, arb_midref(t), ARF_RND_NEAR); arb_sub_fmpz(t, t, v, prec); arb_const_pi(u, prec); arb_mul(t, t, u, prec); arb_mul_2exp_si(t, t, -1); switch (fmpz_fdiv_ui(v, 4)) { case 0: arb_sin_cos(s, c, t, prec); break; case 1: arb_sin_cos(c, s, t, prec); arb_neg(c, c); break; case 2: arb_sin_cos(s, c, t, prec); arb_neg(s, s); arb_neg(c, c); break; default: arb_sin_cos(c, s, t, prec); arb_neg(s, s); break; } fmpz_clear(v); arb_clear(t); arb_clear(u); }
void mag_expinv(mag_t res, const mag_t x) { if (mag_is_zero(x)) { mag_one(res); } else if (mag_is_inf(x)) { mag_zero(res); } else if (fmpz_sgn(MAG_EXPREF(x)) <= 0) { mag_one(res); } else if (fmpz_cmp_ui(MAG_EXPREF(x), 2 * MAG_BITS) > 0) { fmpz_t t; fmpz_init(t); /* If x > 2^60, exp(-x) < 2^(-2^60 / log(2)) */ /* -1/log(2) < -369/256 */ fmpz_set_si(t, -369); fmpz_mul_2exp(t, t, 2 * MAG_BITS - 8); mag_one(res); mag_mul_2exp_fmpz(res, res, t); fmpz_clear(t); } else { fmpz_t t; slong e = MAG_EXP(x); fmpz_init(t); fmpz_set_ui(t, MAG_MAN(x)); if (e >= MAG_BITS) fmpz_mul_2exp(t, t, e - MAG_BITS); else fmpz_tdiv_q_2exp(t, t, MAG_BITS - e); /* upper bound for 1/e */ mag_set_ui_2exp_si(res, 395007543, -30); mag_pow_fmpz(res, res, t); fmpz_clear(t); } }
static void arb_sqrt1pm1_tiny(arb_t r, const arb_t z, slong prec) { mag_t b, c; arb_t t; mag_init(b); mag_init(c); arb_init(t); /* if |z| < 1, then |(sqrt(1+z)-1) - (z/2-z^2/8)| <= |z|^3/(1-|z|)/16 */ arb_get_mag(b, z); mag_one(c); mag_sub_lower(c, c, b); mag_pow_ui(b, b, 3); mag_div(b, b, c); mag_mul_2exp_si(b, b, -4); arb_mul(t, z, z, prec); arb_mul_2exp_si(t, t, -2); arb_sub(r, z, t, prec); arb_mul_2exp_si(r, r, -1); if (mag_is_finite(b)) arb_add_error_mag(r, b); else arb_indeterminate(r); mag_clear(b); mag_clear(c); arb_clear(t); }
static void acb_log1p_tiny(acb_t r, const acb_t z, slong prec) { mag_t b, c; acb_t t; int real; mag_init(b); mag_init(c); acb_init(t); real = acb_is_real(z); /* if |z| < 1, then |log(1+z) - [z - z^2/2]| <= |z|^3/(1-|z|) */ acb_get_mag(b, z); mag_one(c); mag_sub_lower(c, c, b); mag_pow_ui(b, b, 3); mag_div(b, b, c); acb_mul(t, z, z, prec); acb_mul_2exp_si(t, t, -1); acb_sub(r, z, t, prec); if (real && mag_is_finite(b)) arb_add_error_mag(acb_realref(r), b); else acb_add_error_mag(r, b); mag_clear(b); mag_clear(c); acb_clear(t); }
void _hypgeom_precompute(hypgeom_t hyp, const fmpz_poly_t P, const fmpz_poly_t Q) { slong k; fmpz_t t; fmpz_init(t); hyp->r = fmpz_poly_degree(Q) - fmpz_poly_degree(P); hyp->boundC = hypgeom_root_norm(P); hyp->boundD = hypgeom_root_norm(Q); hyp->boundK = 1 + FLINT_MAX(hyp->boundC, 2 * hyp->boundD); mag_one(hyp->MK); for (k = 1; k <= hyp->boundK; k++) { fmpz_poly_evaluate_si(t, P, k); mag_mul_fmpz(hyp->MK, hyp->MK, t); fmpz_poly_evaluate_si(t, Q, k); mag_div_fmpz(hyp->MK, hyp->MK, t); } fmpz_clear(t); }
/* Differential equation for F(a,b,c,y+z): (y+z)(y-1+z) F''(z) + ((y+z)(a+b+1) - c) F'(z) + a b F(z) = 0 Coefficients in the Taylor series are bounded by A * binomial(N+k, k) * nu^k using the Cauchy-Kovalevskaya majorant method. See J. van der Hoeven, "Fast evaluation of holonomic functions near and in regular singularities" */ static void bound(mag_t A, mag_t nu, mag_t N, const acb_t a, const acb_t b, const acb_t c, const acb_t y, const acb_t f0, const acb_t f1) { mag_t M0, M1, t, u; acb_t d; acb_init(d); mag_init(M0); mag_init(M1); mag_init(t); mag_init(u); /* nu = max(1/|y-1|, 1/|y|) = 1/min(|y-1|, |y|) */ acb_get_mag_lower(t, y); acb_sub_ui(d, y, 1, MAG_BITS); acb_get_mag_lower(u, d); mag_min(t, t, u); mag_one(u); mag_div(nu, u, t); /* M0 = 2 nu |ab| */ acb_get_mag(t, a); acb_get_mag(u, b); mag_mul(M0, t, u); mag_mul(M0, M0, nu); mag_mul_2exp_si(M0, M0, 1); /* M1 = 2 nu |(a+b+1)y-c| + 2|a+b+1| */ acb_add(d, a, b, MAG_BITS); acb_add_ui(d, d, 1, MAG_BITS); acb_get_mag(t, d); acb_mul(d, d, y, MAG_BITS); acb_sub(d, d, c, MAG_BITS); acb_get_mag(u, d); mag_mul(u, u, nu); mag_add(M1, t, u); mag_mul_2exp_si(M1, M1, 1); /* N = max(sqrt(2 M0), 2 M1) / nu */ mag_mul_2exp_si(M0, M0, 1); mag_sqrt(M0, M0); mag_mul_2exp_si(M1, M1, 1); mag_max(N, M0, M1); mag_div(N, N, nu); /* A = max(|f0|, |f1| / (nu (N+1)) */ acb_get_mag(t, f0); acb_get_mag(u, f1); mag_div(u, u, nu); mag_div(u, u, N); /* upper bound for dividing by N+1 */ mag_max(A, t, u); acb_clear(d); mag_clear(M0); mag_clear(M1); mag_clear(t); mag_clear(u); }
static void acb_hypgeom_mag_Cn(mag_t Cn, int R, const mag_t nu, const mag_t sigma, ulong n) { if (R == 1) { mag_one(Cn); } else { acb_hypgeom_mag_chi(Cn, n); if (R == 3) { mag_t tmp; mag_init(tmp); mag_mul(tmp, nu, nu); mag_mul(tmp, tmp, sigma); if (n != 1) mag_mul_ui(tmp, tmp, n); mag_add(Cn, Cn, tmp); mag_pow_ui(tmp, nu, n); mag_mul(Cn, Cn, tmp); mag_clear(tmp); } } }
void acb_hypgeom_erf_asymp(acb_t res, const acb_t z, slong prec, slong prec2) { acb_t a, t, u; acb_init(a); acb_init(t); acb_init(u); acb_one(a); acb_mul_2exp_si(a, a, -1); acb_mul(t, z, z, prec2); acb_hypgeom_u_asymp(u, a, a, t, -1, prec2); acb_neg(t, t); acb_exp(t, t, prec2); acb_mul(u, u, t, prec2); acb_const_pi(t, prec2); acb_sqrt(t, t, prec2); acb_mul(t, t, z, prec2); acb_div(u, u, t, prec2); /* branch cut term: -1 or 1 */ if (arb_contains_zero(acb_realref(z))) { arb_zero(acb_imagref(t)); arf_zero(arb_midref(acb_realref(t))); mag_one(arb_radref(acb_realref(t))); } else { acb_set_si(t, arf_sgn(arb_midref(acb_realref(z)))); } acb_sub(t, t, u, prec); if (arb_is_zero(acb_imagref(z))) arb_zero(acb_imagref(t)); else if (arb_is_zero(acb_realref(z))) arb_zero(acb_realref(t)); acb_set(res, t); acb_clear(a); acb_clear(t); acb_clear(u); }
void arb_sinc(arb_t z, const arb_t x, slong prec) { mag_t c, r; mag_init(c); mag_init(r); mag_set_ui_2exp_si(c, 5, -1); arb_get_mag_lower(r, x); if (mag_cmp(c, r) < 0) { /* x is not near the origin */ _arb_sinc_direct(z, x, prec); } else if (mag_cmp_2exp_si(arb_radref(x), 1) < 0) { /* determine error magnitude using the derivative bound */ if (arb_is_exact(x)) { mag_zero(c); } else { _arb_sinc_derivative_bound(r, x); mag_mul(c, arb_radref(x), r); } /* evaluate sinc at the midpoint of x */ if (arf_is_zero(arb_midref(x))) { arb_one(z); } else { arb_get_mid_arb(z, x); _arb_sinc_direct(z, z, prec); } /* add the error */ mag_add(arb_radref(z), arb_radref(z), c); } else { /* x has a large radius and includes points near the origin */ arf_zero(arb_midref(z)); mag_one(arb_radref(z)); } mag_clear(c); mag_clear(r); }
void _arb_sinc_derivative_bound(mag_t d, const arb_t x) { /* |f'(x)| < min(arb_get_mag(x), 1) / 2 */ mag_t r, one; mag_init(r); mag_init(one); arb_get_mag(r, x); mag_one(one); mag_min(d, r, one); mag_mul_2exp_si(d, d, -1); mag_clear(r); mag_clear(one); }
void arb_root_ui_algebraic(arb_t res, const arb_t x, ulong k, slong prec) { mag_t r, msubr, m1k, t; if (arb_is_exact(x)) { arb_root_arf(res, arb_midref(x), k, prec); return; } if (!arb_is_nonnegative(x)) { arb_indeterminate(res); return; } mag_init(r); mag_init(msubr); mag_init(m1k); mag_init(t); /* x = [m-r, m+r] */ mag_set(r, arb_radref(x)); /* m - r */ arb_get_mag_lower(msubr, x); /* m^(1/k) */ arb_root_arf(res, arb_midref(x), k, prec); /* bound for m^(1/k) */ arb_get_mag(m1k, res); /* C = min(1, log(1+r/(m-r))/k) */ mag_div(t, r, msubr); mag_log1p(t, t); mag_div_ui(t, t, k); if (mag_cmp_2exp_si(t, 0) > 0) mag_one(t); /* C m^(1/k) */ mag_mul(t, m1k, t); mag_add(arb_radref(res), arb_radref(res), t); mag_clear(r); mag_clear(msubr); mag_clear(m1k); mag_clear(t); }
void arb_bernoulli_fmpz(arb_t res, const fmpz_t n, slong prec) { if (fmpz_cmp_ui(n, UWORD_MAX) <= 0) { if (fmpz_sgn(n) >= 0) arb_bernoulli_ui(res, fmpz_get_ui(n), prec); else arb_zero(res); } else if (fmpz_is_odd(n)) { arb_zero(res); } else { arb_t t; slong wp; arb_init(t); wp = prec + 2 * fmpz_bits(n); /* zeta(n) ~= 1 */ arf_one(arb_midref(res)); mag_one(arb_radref(res)); mag_mul_2exp_si(arb_radref(res), arb_radref(res), WORD_MIN); /* |B_n| = 2 * n! / (2*pi)^n * zeta(n) */ arb_gamma_fmpz(t, n, wp); arb_mul_fmpz(t, t, n, wp); arb_mul(res, res, t, wp); arb_const_pi(t, wp); arb_mul_2exp_si(t, t, 1); arb_pow_fmpz(t, t, n, wp); arb_div(res, res, t, prec); arb_mul_2exp_si(res, res, 1); if (fmpz_fdiv_ui(n, 4) == 0) arb_neg(res, res); arb_clear(t); } }
slong renf_set_embeddings_fmpz_poly(renf * nf, fmpz_poly_t pol, slong lim, slong prec) { slong i, n, n_exact, n_interval; fmpq_poly_t p2; arb_t a; fmpz * c; slong * k; n = fmpz_poly_num_real_roots_upper_bound(pol); c = _fmpz_vec_init(n); k = (slong *) flint_malloc(n * sizeof(slong)); fmpz_poly_isolate_real_roots(NULL, &n_exact, c, k, &n_interval, pol); if (n_exact) { fprintf(stderr, "ERROR (fmpz_poly_real_embeddings): rational roots\n"); abort(); } arb_init(a); fmpq_poly_init(p2); fmpz_one(fmpq_poly_denref(p2)); fmpq_poly_fit_length(p2, pol->length); _fmpz_vec_set(p2->coeffs, pol->coeffs, pol->length); p2->length = pol->length; for (i = 0; i < FLINT_MIN(lim, n_interval); i++) { arb_set_fmpz(a, c + i); arb_mul_2exp_si(a, a, 1); arb_add_si(a, a, 1, prec); mag_one(arb_radref(a)); arb_mul_2exp_si(a, a, k[i] - 1); renf_init(nf + i, p2, a, prec); } arb_clear(a); fmpq_poly_clear(p2); _fmpz_vec_clear(c, n); flint_free(k); return n_interval; }
static void phase(acb_t res, const arb_t re, const arb_t im) { if (arb_is_nonnegative(re) || arb_is_negative(im)) { acb_one(res); } else if (arb_is_negative(re) && arb_is_nonnegative(im)) { acb_set_si(res, -3); } else { arb_zero(acb_imagref(res)); /* -1 +/- 2 */ arf_set_si(arb_midref(acb_realref(res)), -1); mag_one(arb_radref(acb_realref(res))); mag_mul_2exp_si(arb_radref(acb_realref(res)), arb_radref(acb_realref(res)), 1); } }
/* f(z) = sin(1/z), assume on real interval */ int f_essing(acb_ptr res, const acb_t z, void * param, slong order, slong prec) { if (order > 1) flint_abort(); /* Would be needed for Taylor method. */ if ((order == 0) && acb_is_real(z) && arb_contains_zero(acb_realref(z))) { /* todo: arb_zero_pm_one, arb_unit_interval? */ acb_zero(res); mag_one(arb_radref(acb_realref(res))); } else { acb_inv(res, z, prec); acb_sin(res, res, prec); } return 0; }
void arb_atan(arb_t z, const arb_t x, slong prec) { if (arb_is_exact(x)) { arb_atan_arf(z, arb_midref(x), prec); } else { mag_t t, u; mag_init(t); mag_init(u); arb_get_mag_lower(t, x); if (mag_is_zero(t)) { mag_set(t, arb_radref(x)); } else { mag_mul_lower(t, t, t); mag_one(u); mag_add_lower(t, t, u); mag_div(t, arb_radref(x), t); } if (mag_cmp_2exp_si(t, 0) > 0) { mag_const_pi(u); mag_min(t, t, u); } arb_atan_arf(z, arb_midref(x), prec); mag_add(arb_radref(z), arb_radref(z), t); mag_clear(t); mag_clear(u); } }
void mag_pow_ui_lower(mag_t z, const mag_t x, ulong e) { if (e <= 2) { if (e == 0) mag_one(z); else if (e == 1) mag_set(z, x); else mag_mul_lower(z, x, x); } else if (mag_is_inf(x)) { mag_inf(z); } else { mag_t y; int i, bits; mag_init_set(y, x); bits = FLINT_BIT_COUNT(e); for (i = bits - 2; i >= 0; i--) { mag_mul_lower(y, y, y); if (e & (1UL << i)) mag_mul_lower(y, y, x); } mag_swap(z, y); mag_clear(y); } }
int main() { long iter; flint_rand_t state; printf("const_glaisher...."); fflush(stdout); flint_randinit(state); for (iter = 0; iter < 250; iter++) { arb_t r, s, t; fmpz_t v; long accuracy, prec; prec = 2 + n_randint(state, 2000); arb_init(r); arb_init(s); arb_init(t); fmpz_init(v); arb_const_glaisher(r, prec); arb_const_glaisher(s, prec + 100); if (!arb_overlaps(r, s)) { printf("FAIL: containment\n\n"); printf("prec = %ld\n", prec); printf("r = "); arb_printd(r, prec / 3.33); printf("\n\n"); abort(); } accuracy = arb_rel_accuracy_bits(r); if (accuracy < prec - 4) { printf("FAIL: poor accuracy\n\n"); printf("prec = %ld\n", prec); printf("r = "); arb_printd(r, prec / 3.33); printf("\n\n"); abort(); } if (n_randint(state, 30) == 0) { flint_cleanup(); } fmpz_set_str(v, "128242712910062263687534256886979172776768892732500", 10); arb_set_fmpz(t, v); mag_one(arb_radref(t)); fmpz_ui_pow_ui(v, 10, 50); arb_div_fmpz(t, t, v, 170); if (!arb_overlaps(r, t)) { printf("FAIL: reference value\n\n"); printf("prec = %ld\n", prec); printf("r = "); arb_printd(r, prec / 3.33); printf("\n\n"); abort(); } arb_clear(r); arb_clear(s); arb_clear(t); fmpz_clear(v); } flint_randclear(state); flint_cleanup(); printf("PASS\n"); return EXIT_SUCCESS; }
void _arb_bell_sum_taylor(arb_t res, const fmpz_t n, const fmpz_t a, const fmpz_t b, const fmpz_t mmag, long tol) { fmpz_t m, r, R, tmp; mag_t B, C, D, bound; arb_t t, u; long wp, k, N; if (_fmpz_sub_small(b, a) < 5) { arb_bell_sum_bsplit(res, n, a, b, mmag, tol); return; } fmpz_init(m); fmpz_init(r); fmpz_init(R); fmpz_init(tmp); /* r = max(m - a, b - m) */ /* m = a + (b - a) / 2 */ fmpz_sub(r, b, a); fmpz_cdiv_q_2exp(r, r, 1); fmpz_add(m, a, r); fmpz_mul_2exp(R, r, RADIUS_BITS); mag_init(B); mag_init(C); mag_init(D); mag_init(bound); arb_init(t); arb_init(u); if (fmpz_cmp(R, m) >= 0) { mag_inf(C); mag_inf(D); } else { /* C = exp(R * |F'(m)| + (1/2) R^2 * (n/(m-R)^2 + 1/(m-R))) */ /* C = exp(R * (|F'(m)| + (1/2) R * (n/(m-R) + 1)/(m-R))) */ /* D = (1/2) R * (n/(m-R) + 1)/(m-R) */ fmpz_sub(tmp, m, R); mag_set_fmpz(D, n); mag_div_fmpz(D, D, tmp); mag_one(C); mag_add(D, D, C); mag_div_fmpz(D, D, tmp); mag_mul_fmpz(D, D, R); mag_mul_2exp_si(D, D, -1); /* C = |F'(m)| */ wp = 20 + 1.05 * fmpz_bits(n); arb_set_fmpz(t, n); arb_div_fmpz(t, t, m, wp); fmpz_add_ui(tmp, m, 1); arb_set_fmpz(u, tmp); arb_digamma(u, u, wp); arb_sub(t, t, u, wp); arb_get_mag(C, t); /* C = exp(R * (C + D)) */ mag_add(C, C, D); mag_mul_fmpz(C, C, R); mag_exp(C, C); } if (mag_cmp_2exp_si(C, tol / 4 + 2) > 0) { _arb_bell_sum_taylor(res, n, a, m, mmag, tol); _arb_bell_sum_taylor(t, n, m, b, mmag, tol); arb_add(res, res, t, 2 * tol); } else { arb_ptr mx, ser1, ser2, ser3; /* D = T(m) */ wp = 20 + 1.05 * fmpz_bits(n); arb_set_fmpz(t, m); arb_pow_fmpz(t, t, n, wp); fmpz_add_ui(tmp, m, 1); arb_gamma_fmpz(u, tmp, wp); arb_div(t, t, u, wp); arb_get_mag(D, t); /* error bound: (b-a) * C * D * B^N / (1 - B), B = r/R */ /* ((b-a) * C * D * 2) * 2^(-N*RADIUS_BITS) */ /* ((b-a) * C * D * 2) */ mag_mul(bound, C, D); mag_mul_2exp_si(bound, bound, 1); fmpz_sub(tmp, b, a); mag_mul_fmpz(bound, bound, tmp); /* N = (tol + log2((b-a)*C*D*2) - mmag) / RADIUS_BITS */ if (mmag == NULL) { /* estimate D ~= 2^mmag */ fmpz_add_ui(tmp, MAG_EXPREF(C), tol); fmpz_cdiv_q_ui(tmp, tmp, RADIUS_BITS); } else { fmpz_sub(tmp, MAG_EXPREF(bound), mmag); fmpz_add_ui(tmp, tmp, tol); fmpz_cdiv_q_ui(tmp, tmp, RADIUS_BITS); } if (fmpz_cmp_ui(tmp, 5 * tol / 4) > 0) N = 5 * tol / 4; else if (fmpz_cmp_ui(tmp, 2) < 0) N = 2; else N = fmpz_get_ui(tmp); /* multiply by 2^(-N*RADIUS_BITS) */ mag_mul_2exp_si(bound, bound, -N * RADIUS_BITS); mx = _arb_vec_init(2); ser1 = _arb_vec_init(N); ser2 = _arb_vec_init(N); ser3 = _arb_vec_init(N); /* estimate (this should work for moderate n and tol) */ wp = 1.1 * tol + 1.05 * fmpz_bits(n) + 5; /* increase precision until convergence */ while (1) { /* (m+x)^n / gamma(m+1+x) */ arb_set_fmpz(mx, m); arb_one(mx + 1); _arb_poly_log_series(ser1, mx, 2, N, wp); for (k = 0; k < N; k++) arb_mul_fmpz(ser1 + k, ser1 + k, n, wp); arb_add_ui(mx, mx, 1, wp); _arb_poly_lgamma_series(ser2, mx, 2, N, wp); _arb_vec_sub(ser1, ser1, ser2, N, wp); _arb_poly_exp_series(ser3, ser1, N, N, wp); /* t = a - m, u = b - m */ arb_set_fmpz(t, a); arb_sub_fmpz(t, t, m, wp); arb_set_fmpz(u, b); arb_sub_fmpz(u, u, m, wp); arb_power_sum_vec(ser1, t, u, N, wp); arb_zero(res); for (k = 0; k < N; k++) arb_addmul(res, ser3 + k, ser1 + k, wp); if (mmag != NULL) { if (_fmpz_sub_small(MAG_EXPREF(arb_radref(res)), mmag) <= -tol) break; } else { if (arb_rel_accuracy_bits(res) >= tol) break; } wp = 2 * wp; } /* add the series truncation bound */ arb_add_error_mag(res, bound); _arb_vec_clear(mx, 2); _arb_vec_clear(ser1, N); _arb_vec_clear(ser2, N); _arb_vec_clear(ser3, N); } mag_clear(B); mag_clear(C); mag_clear(D); mag_clear(bound); arb_clear(t); arb_clear(u); fmpz_clear(m); fmpz_clear(r); fmpz_clear(R); fmpz_clear(tmp); }
void mag_polylog_tail(mag_t u, const mag_t z, long sigma, ulong d, ulong N) { mag_t TN, UN, t; if (N < 2) { mag_inf(u); return; } mag_init(TN); mag_init(UN); mag_init(t); if (mag_cmp_2exp_si(z, 0) >= 0) { mag_inf(u); } else { /* Bound T(N) */ mag_pow_ui(TN, z, N); /* multiply by log(N)^d */ if (d > 0) { mag_log_ui(t, N); mag_pow_ui(t, t, d); mag_mul(TN, TN, t); } /* multiply by 1/k^s */ if (sigma > 0) { mag_set_ui_lower(t, N); mag_pow_ui_lower(t, t, sigma); mag_div(TN, TN, t); } else if (sigma < 0) { mag_set_ui(t, N); mag_pow_ui(t, t, -sigma); mag_mul(TN, TN, t); } /* Bound U(N) */ mag_set(UN, z); /* multiply by (1 + 1/N)**S */ if (sigma < 0) { mag_binpow_uiui(t, N, -sigma); mag_mul(UN, UN, t); } /* multiply by (1 + 1/(N log(N)))^d */ if (d > 0) { ulong nl; /* rounds down */ nl = mag_d_log_lower_bound(N) * N * (1 - 1e-13); mag_binpow_uiui(t, nl, d); mag_mul(UN, UN, t); } /* T(N) / (1 - U(N)) */ if (mag_cmp_2exp_si(UN, 0) >= 0) { mag_inf(u); } else { mag_one(t); mag_sub_lower(t, t, UN); mag_div(u, TN, t); } } mag_clear(TN); mag_clear(UN); mag_clear(t); }
void acb_rising_ui_get_mag(mag_t bound, const acb_t s, ulong n) { if (n == 0) { mag_one(bound); return; } if (n == 1) { acb_get_mag(bound, s); return; } if (!acb_is_finite(s)) { mag_inf(bound); return; } if (arf_sgn(arb_midref(acb_realref(s))) >= 0) { acb_rising_get_mag2_right(bound, acb_realref(s), acb_imagref(s), n); } else { arb_t a; long k; mag_t bound2, t, u; arb_init(a); mag_init(bound2); mag_init(t); mag_init(u); arb_get_mag(u, acb_imagref(s)); mag_mul(u, u, u); mag_one(bound); for (k = 0; k < n; k++) { arb_add_ui(a, acb_realref(s), k, MAG_BITS); if (arf_sgn(arb_midref(a)) >= 0) { acb_rising_get_mag2_right(bound2, a, acb_imagref(s), n - k); mag_mul(bound, bound, bound2); break; } else { arb_get_mag(t, a); mag_mul(t, t, t); mag_add(t, t, u); mag_mul(bound, bound, t); } } arb_clear(a); mag_clear(bound2); mag_clear(t); mag_clear(u); } mag_sqrt(bound, bound); }
/* computes the factors that are independent of n (all are upper bounds) */ void acb_hypgeom_u_asymp_bound_factors(int * R, mag_t alpha, mag_t nu, mag_t sigma, mag_t rho, mag_t zinv, const acb_t a, const acb_t b, const acb_t z) { mag_t r, u, zre, zim, zlo, sigma_prime; acb_t t; mag_init(r); mag_init(u); mag_init(zre); mag_init(zim); mag_init(zlo); mag_init(sigma_prime); acb_init(t); /* lower bounds for |re(z)|, |im(z)|, |z| */ arb_get_mag_lower(zre, acb_realref(z)); arb_get_mag_lower(zim, acb_imagref(z)); acb_get_mag_lower(zlo, z); /* todo: hypot */ /* upper bound for 1/|z| */ mag_one(u); mag_div(zinv, u, zlo); /* upper bound for r = |b - 2a| */ acb_mul_2exp_si(t, a, 1); acb_sub(t, b, t, MAG_BITS); acb_get_mag(r, t); /* determine region */ *R = 0; if (mag_cmp(zlo, r) >= 0) { int znonneg = arb_is_nonnegative(acb_realref(z)); if (znonneg && mag_cmp(zre, r) >= 0) { *R = 1; } else if (mag_cmp(zim, r) >= 0 || znonneg) { *R = 2; } else { mag_mul_2exp_si(u, r, 1); if (mag_cmp(zlo, u) >= 0) *R = 3; } } if (R == 0) { mag_inf(alpha); mag_inf(nu); mag_inf(sigma); mag_inf(rho); } else { /* sigma = |(b-2a)/z| */ mag_mul(sigma, r, zinv); /* nu = (1/2 + 1/2 sqrt(1-4 sigma^2))^(-1/2) <= 1 + 2 sigma^2 */ if (mag_cmp_2exp_si(sigma, -1) <= 0) { mag_mul(nu, sigma, sigma); mag_mul_2exp_si(nu, nu, 1); mag_one(u); mag_add(nu, nu, u); } else { mag_inf(nu); } /* modified sigma for alpha, beta, rho when in R3 */ if (*R == 3) mag_mul(sigma_prime, sigma, nu); else mag_set(sigma_prime, sigma); /* alpha = 1/(1-sigma') */ mag_one(alpha); mag_sub_lower(alpha, alpha, sigma_prime); mag_one(u); mag_div(alpha, u, alpha); /* rho = |2a^2-2ab+b|/2 + sigma'*(1+sigma'/4)/(1-sigma')^2 */ mag_mul_2exp_si(rho, sigma_prime, -2); mag_one(u); mag_add(rho, rho, u); mag_mul(rho, rho, sigma_prime); mag_mul(rho, rho, alpha); mag_mul(rho, rho, alpha); acb_sub(t, a, b, MAG_BITS); acb_mul(t, t, a, MAG_BITS); acb_mul_2exp_si(t, t, 1); acb_add(t, t, b, MAG_BITS); acb_get_mag(u, t); mag_mul_2exp_si(u, u, -1); mag_add(rho, rho, u); } mag_clear(r); mag_clear(u); mag_clear(zre); mag_clear(zim); mag_clear(zlo); mag_clear(sigma_prime); acb_clear(t); }
slong hypgeom_bound(mag_t error, int r, slong A, slong B, slong K, const mag_t TK, const mag_t z, slong tol_2exp) { mag_t Tn, t, u, one, tol, num, den; slong n, m; mag_init(Tn); mag_init(t); mag_init(u); mag_init(one); mag_init(tol); mag_init(num); mag_init(den); mag_one(one); mag_set_ui_2exp_si(tol, UWORD(1), -tol_2exp); /* approximate number of needed terms */ n = hypgeom_estimate_terms(z, r, tol_2exp); /* required for 1 + O(1/k) part to be decreasing */ n = FLINT_MAX(n, K + 1); /* required for z^k / (k!)^r to be decreasing */ m = hypgeom_root_bound(z, r); n = FLINT_MAX(n, m); /* We now have |R(k)| <= G(k) where G(k) is monotonically decreasing, and can bound the tail using a geometric series as soon as soon as G(k) < 1. */ /* bound T(n-1) */ hypgeom_term_bound(Tn, TK, K, A, B, r, z, n-1); while (1) { /* bound R(n) */ mag_mul_ui(num, z, n); mag_mul_ui(num, num, n - B); mag_set_ui_lower(den, n - A); mag_mul_ui_lower(den, den, n - 2*B); if (r != 0) { mag_set_ui_lower(u, n); mag_pow_ui_lower(u, u, r); mag_mul_lower(den, den, u); } mag_div(t, num, den); /* multiply bound for T(n-1) by bound for R(n) to bound T(n) */ mag_mul(Tn, Tn, t); /* geometric series termination check */ /* u = max(1-t, 0), rounding down [lower bound] */ mag_sub_lower(u, one, t); if (!mag_is_zero(u)) { mag_div(u, Tn, u); if (mag_cmp(u, tol) < 0) { mag_set(error, u); break; } } /* move on to next term */ n++; } mag_clear(Tn); mag_clear(t); mag_clear(u); mag_clear(one); mag_clear(tol); mag_clear(num); mag_clear(den); return n; }
void acb_hypgeom_pfq_sum_rs(acb_t res, acb_t term, acb_srcptr a, slong p, acb_srcptr b, slong q, const acb_t z, slong n, slong prec) { acb_ptr zpow; acb_t s, t, u; slong i, j, k, m; mag_t B, C; if (n == 0) { acb_zero(res); acb_one(term); return; } if (n < 0) abort(); m = n_sqrt(n); m = FLINT_MIN(m, 150); mag_init(B); mag_init(C); acb_init(s); acb_init(t); acb_init(u); zpow = _acb_vec_init(m + 1); _acb_vec_set_powers(zpow, z, m + 1, prec); mag_one(B); for (k = n; k >= 0; k--) { j = k % m; if (k < n) acb_add(s, s, zpow + j, prec); if (k > 0) { if (p > 0) { acb_add_ui(u, a, k - 1, prec); for (i = 1; i < p; i++) { acb_add_ui(t, a + i, k - 1, prec); acb_mul(u, u, t, prec); } if (k < n) acb_mul(s, s, u, prec); acb_get_mag(C, u); mag_mul(B, B, C); } if (q > 0) { acb_add_ui(u, b, k - 1, prec); for (i = 1; i < q; i++) { acb_add_ui(t, b + i, k - 1, prec); acb_mul(u, u, t, prec); } if (k < n) acb_div(s, s, u, prec); acb_get_mag_lower(C, u); mag_div(B, B, C); } if (j == 0 && k < n) { acb_mul(s, s, zpow + m, prec); } } } acb_get_mag(C, z); mag_pow_ui(C, C, n); mag_mul(B, B, C); acb_zero(term); if (_acb_vec_is_real(a, p) && _acb_vec_is_real(b, q) && acb_is_real(z)) arb_add_error_mag(acb_realref(term), B); else acb_add_error_mag(term, B); acb_set(res, s); mag_clear(B); mag_clear(C); acb_clear(s); acb_clear(t); acb_clear(u); _acb_vec_clear(zpow, m + 1); }