Merge pull request #258 from p15-git-acc/backlund_s_bound

add backlund_s_bound
This commit is contained in:
Fredrik Johansson 2019-02-08 02:17:31 +01:00 committed by GitHub
commit 2a5650a006
No known key found for this signature in database
GPG key ID: 4AEE18F83AFDEB23
5 changed files with 143 additions and 0 deletions

View file

@ -166,6 +166,7 @@ void _acb_dirichlet_hardy_z_series(acb_ptr res, acb_srcptr s, slong slen, const
void acb_dirichlet_hardy_z_series(acb_poly_t res, const acb_poly_t s, const dirichlet_group_t G, const dirichlet_char_t chi, slong len, slong prec);
void acb_dirichlet_gram_point(arb_t res, const fmpz_t n, const dirichlet_group_t G, const dirichlet_char_t chi, slong prec);
void acb_dirichlet_backlund_s_bound(mag_t res, const arb_t t);
/* Discrete Fourier Transform */

View file

@ -0,0 +1,64 @@
/*
Copyright (C) 2019 D.H.J Polymath
This file is part of Arb.
Arb is free software: you can redistribute it and/or modify it under
the terms of the GNU Lesser General Public License (LGPL) as published
by the Free Software Foundation; either version 2.1 of the License, or
(at your option) any later version. See <http://www.gnu.org/licenses/>.
*/
#include "acb_dirichlet.h"
static void
_mag_div_ui_ui(mag_t res, ulong a, ulong b)
{
mag_set_ui(res, a);
mag_div_ui(res, res, b);
}
void
acb_dirichlet_backlund_s_bound(mag_t res, const arb_t t)
{
if (!arb_is_nonnegative(t))
{
mag_inf(res);
}
else
{
mag_t m;
mag_init(m);
arb_get_mag(m, t);
if (mag_cmp_2exp_si(m, 8) < 0) /* 2^8 < 280 */
{
mag_one(res);
}
else if (mag_cmp_2exp_si(m, 22) < 0) /* 2^22 < 6.8*10^6 */
{
mag_set_ui(res, 2);
}
else if (mag_cmp_2exp_si(m, 29) < 0) /* 2^29 < 5.45*10^8 */
{
_mag_div_ui_ui(res, 231366, 100000);
}
else
{
/* |S(t)| <= 0.112*log(t) + 0.278*log(log(t)) + 2.51 */
mag_t c, logm;
mag_init(c);
mag_init(logm);
mag_log(logm, m);
_mag_div_ui_ui(c, 278, 1000);
mag_log(res, logm);
mag_mul(res, res, c);
_mag_div_ui_ui(c, 112, 1000);
mag_addmul(res, c, logm);
_mag_div_ui_ui(c, 251, 100);
mag_add(res, res, c);
mag_clear(c);
mag_clear(logm);
}
mag_clear(m);
}
}

View file

@ -0,0 +1,71 @@
/*
Copyright (C) 2019 D.H.J Polymath
This file is part of Arb.
Arb is free software: you can redistribute it and/or modify it under
the terms of the GNU Lesser General Public License (LGPL) as published
by the Free Software Foundation; either version 2.1 of the License, or
(at your option) any later version. See <http://www.gnu.org/licenses/>.
*/
#include "acb_dirichlet.h"
int main()
{
slong iter;
flint_rand_t state;
flint_printf("backlund_s_bound....");
fflush(stdout);
flint_randinit(state);
for (iter = 0; iter < 500 * arb_test_multiplier(); iter++)
{
arb_t a, b;
mag_t u, v;
slong aprec, bprec;
slong abits, bbits;
aprec = 2 + n_randint(state, 1000);
bprec = 2 + n_randint(state, 1000);
abits = 2 + n_randint(state, 100);
bbits = 2 + n_randint(state, 100);
arb_init(a);
arb_init(b);
mag_init(u);
mag_init(v);
arb_randtest(a, state, aprec, abits);
arb_randtest(b, state, bprec, bbits);
if (arb_is_nonnegative(a) && arb_is_nonnegative(b))
{
acb_dirichlet_backlund_s_bound(u, a);
acb_dirichlet_backlund_s_bound(v, b);
if ((arb_lt(a, b) && mag_cmp(u, v) > 0) ||
(arb_gt(a, b) && mag_cmp(u, v) < 0))
{
flint_printf("FAIL: increasing on t >= 0\n\n");
flint_printf("a = "); arb_print(a); flint_printf("\n\n");
flint_printf("b = "); arb_print(b); flint_printf("\n\n");
flint_printf("u = "); mag_print(u); flint_printf("\n\n");
flint_printf("v = "); mag_print(v); flint_printf("\n\n");
flint_abort();
}
}
arb_clear(a);
arb_clear(b);
mag_clear(u);
mag_clear(v);
}
flint_randclear(state);
flint_cleanup();
flint_printf("PASS\n");
return EXIT_SUCCESS;
}

View file

@ -620,3 +620,8 @@ Currently, these methods require *chi* to be a primitive character.
Riemann zeta function are supported and *G* and *chi* must both be set to
*NULL*. Requires `n \ge 0`.
.. function:: void acb_dirichlet_backlund_s_bound(mag_t res, const arb_t t)
Computes an upper bound for `|S(t)|` quickly. Theorem 1
and the bounds in (1.2) in [Tru2014]_ are used.

View file

@ -250,3 +250,5 @@ Bibliography
.. [Tak2000] \D. Takahashi, "A fast algorithm for computing large Fibonacci numbers", Information Processing Letters 75 (2000) 243-246, http://www.ii.uni.wroc.pl/~lorys/IPL/article75-6-1.pdf
.. [Tre2008] \L. N. Trefethen, "Is Gauss Quadrature Better than Clenshaw-Curtis?", SIAM Review, 50:1 (2008), 67-87, https://doi.org/10.1137/060659831
.. [Tru2014] \T. S. Trudgian, "An improved upper bound for the argument of the Riemann zeta-function on the critical line II", Journal of Number Theory 134 (2014), 280-292, https://doi.org/10.1016/j.jnt.2013.07.017