Exemple #1
0
/* assumes no aliasing */
slong
acb_lambertw_initial(acb_t res, const acb_t z, const acb_t ez1, const fmpz_t k, slong prec)
{
    /* Handle z very close to 0 on the principal branch. */
    if (fmpz_is_zero(k) && 
            (arf_cmpabs_2exp_si(arb_midref(acb_realref(z)), -20) <= 0 &&
             arf_cmpabs_2exp_si(arb_midref(acb_imagref(z)), -20) <= 0))
    {
        acb_set(res, z);
        acb_submul(res, res, res, prec);
        return 40;  /* could be tightened... */
    }

    /* For moderate input not close to the branch point, compute a double
       approximation as the initial value. */
    if (fmpz_is_zero(k) &&
        arf_cmpabs_2exp_si(arb_midref(acb_realref(z)), 400) < 0 &&
        arf_cmpabs_2exp_si(arb_midref(acb_imagref(z)), 400) < 0 &&
          (arf_cmp_d(arb_midref(acb_realref(z)), -0.37) < 0 ||
           arf_cmp_d(arb_midref(acb_realref(z)), -0.36) > 0 ||
           arf_cmpabs_d(arb_midref(acb_imagref(z)), 0.01) > 0))
    {
        acb_lambertw_principal_d(res, z);
        return 48;
    }

    /* Check if we are close to the branch point at -1/e. */
    if ((fmpz_is_zero(k) || (fmpz_is_one(k) && arb_is_negative(acb_imagref(z)))
                         || (fmpz_equal_si(k, -1) && arb_is_nonnegative(acb_imagref(z))))
        && ((arf_cmpabs_2exp_si(arb_midref(acb_realref(ez1)), -2) <= 0 &&
             arf_cmpabs_2exp_si(arb_midref(acb_imagref(ez1)), -2) <= 0)))
    {
        acb_t t;
        acb_init(t);
        acb_mul_2exp_si(t, ez1, 1);
        mag_zero(arb_radref(acb_realref(t)));
        mag_zero(arb_radref(acb_imagref(t)));
        acb_mul_ui(t, t, 3, prec);
        acb_sqrt(t, t, prec);
        if (!fmpz_is_zero(k))
            acb_neg(t, t);
        acb_lambertw_branchpoint_series(res, t, 0, prec);
        acb_clear(t);
        return 1;  /* todo: estimate */
    }

    acb_lambertw_initial_asymp(res, z, k, prec);
    return 1;  /* todo: estimate */
}
Exemple #2
0
/* f(z) = sech(10(x-0.2))^2 + sech(100(x-0.4))^4 + sech(1000(x-0.6))^6 */
int
f_spike(acb_ptr res, const acb_t z, void * param, slong order, slong prec)
{
    acb_t a, b, c;

    if (order > 1)
        flint_abort();  /* Would be needed for Taylor method. */

    acb_init(a);
    acb_init(b);
    acb_init(c);

    acb_mul_ui(a, z, 10, prec);
    acb_sub_ui(a, a, 2, prec);
    acb_sech(a, a, prec);
    acb_pow_ui(a, a, 2, prec);

    acb_mul_ui(b, z, 100, prec);
    acb_sub_ui(b, b, 40, prec);
    acb_sech(b, b, prec);
    acb_pow_ui(b, b, 4, prec);

    acb_mul_ui(c, z, 1000, prec);
    acb_sub_ui(c, c, 600, prec);
    acb_sech(c, c, prec);
    acb_pow_ui(c, c, 6, prec);

    acb_add(res, a, b, prec);
    acb_add(res, res, c, prec);

    acb_clear(a);
    acb_clear(b);
    acb_clear(c);

    return 0;
}
Exemple #3
0
void
acb_modular_eisenstein(acb_ptr r, const acb_t tau, slong len, slong prec)
{
    psl2z_t g;
    arf_t one_minus_eps;
    acb_t tau_prime, t1, t2, t3, t4, q;
    slong m, n;

    if (len < 1)
        return;

    psl2z_init(g);
    arf_init(one_minus_eps);
    acb_init(tau_prime);
    acb_init(t1);
    acb_init(t2);
    acb_init(t3);
    acb_init(t4);
    acb_init(q);

    arf_set_ui_2exp_si(one_minus_eps, 63, -6);
    acb_modular_fundamental_domain_approx(tau_prime, g, tau,
        one_minus_eps, prec);

    acb_exp_pi_i(q, tau_prime, prec);
    acb_modular_theta_const_sum(t2, t3, t4, q, prec);

    /* fourth powers of the theta functions (a, b, c) */
    acb_mul(t2, t2, t2, prec);
    acb_mul(t2, t2, t2, prec);
    acb_mul(t2, t2, q, prec);

    acb_mul(t3, t3, t3, prec);
    acb_mul(t3, t3, t3, prec);

    acb_mul(t4, t4, t4, prec);
    acb_mul(t4, t4, t4, prec);

    /* c2 = pi^4 * (a^8 + b^8 + c^8) / 30 */
    /* c3 = pi^6 * (b^12 + c^12 - 3a^8 * (b^4+c^4)) / 180 */

    /* r = a^8 */
    acb_mul(r, t2, t2, prec);

    if (len > 1)
    {
        /* r[1] = -3 a^8 * (b^4 + c^4) */
        acb_add(r + 1, t3, t4, prec);
        acb_mul(r + 1, r + 1, r, prec);
        acb_mul_si(r + 1, r + 1, -3, prec);
    }

    /* b^8 */
    acb_mul(t1, t3, t3, prec);
    acb_add(r, r, t1, prec);

    /* b^12 */
    if (len > 1)
        acb_addmul(r + 1, t1, t3, prec);

    /* c^8 */
    acb_mul(t1, t4, t4, prec);
    acb_add(r, r, t1, prec);

    /* c^12 */
    if (len > 1)
        acb_addmul(r + 1, t1, t4, prec);

    acb_const_pi(t1, prec);
    acb_mul(t1, t1, t1, prec);
    acb_mul(t2, t1, t1, prec);
    acb_mul(r, r, t2, prec);
    acb_div_ui(r, r, 30, prec);

    if (len > 1)
    {
        acb_mul(t2, t2, t1, prec);
        acb_mul(r + 1, r + 1, t2, prec);
        acb_div_ui(r + 1, r + 1, 189, prec);
    }

    /* apply modular transformation */
    if (!fmpz_is_zero(&g->c))
    {

        acb_mul_fmpz(t1, tau, &g->c, prec);
        acb_add_fmpz(t1, t1, &g->d, prec);
        acb_inv(t1, t1, prec);
        acb_mul(t1, t1, t1, prec);
        acb_mul(t2, t1, t1, prec);
        acb_mul(r, r, t2, prec);

        if (len > 1)
        {
            acb_mul(t2, t1, t2, prec);
            acb_mul(r + 1, r + 1, t2, prec);
        }
    }

    /* compute more coefficients using recurrence */
    for (n = 4; n < len + 2; n++)
    {
        acb_zero(r + n - 2);

        m = 2;
        for (m = 2; m * 2 < n; m++)
            acb_addmul(r + n - 2, r + m - 2, r + n - m - 2, prec);

        acb_mul_2exp_si(r + n - 2, r + n - 2, 1);

        if (n % 2 == 0)
            acb_addmul(r + n - 2, r + n / 2 - 2, r + n / 2 - 2, prec);

        acb_mul_ui(r + n - 2, r + n - 2, 3, prec);
        acb_div_ui(r + n - 2, r + n - 2, (2 * n + 1) * (n - 3), prec);
    }

    /* convert c's to G's */
    for (n = 0; n < len; n++)
        acb_div_ui(r + n, r + n, 2 * n + 3, prec);

    psl2z_clear(g);
    arf_clear(one_minus_eps);
    acb_clear(tau_prime);
    acb_clear(t1);
    acb_clear(t2);
    acb_clear(t3);
    acb_clear(t4);
    acb_clear(q);
}
Exemple #4
0
int main(int argc, char *argv[])
{
    acb_t s, t, a, b;
    mag_t tol;
    slong prec, goal;
    slong N;
    ulong k;
    int integral, ifrom, ito;
    int i, twice, havegoal, havetol;
    acb_calc_integrate_opt_t options;

    ifrom = ito = -1;

    for (i = 1; i < argc; i++)
    {
        if (!strcmp(argv[i], "-i"))
        {
            if (!strcmp(argv[i+1], "all"))
            {
                ifrom = 0;
                ito = NUM_INTEGRALS - 1;
            }
            else
            {
                ifrom = ito = atol(argv[i+1]);
                if (ito < 0 || ito >= NUM_INTEGRALS)
                    flint_abort();
            }
        }
    }

    if (ifrom == -1)
    {
        flint_printf("Compute integrals using acb_calc_integrate.\n");
        flint_printf("Usage: integrals -i n [-prec p] [-tol eps] [-twice] [...]\n\n");
        flint_printf("-i n       - compute integral n (0 <= n <= %d), or \"-i all\"\n", NUM_INTEGRALS - 1);
        flint_printf("-prec p    - precision in bits (default p = 64)\n");
        flint_printf("-goal p    - approximate relative accuracy goal (default p)\n");
        flint_printf("-tol eps   - approximate absolute error goal (default 2^-p)\n");
        flint_printf("-twice     - run twice (to see overhead of computing nodes)\n");
        flint_printf("-heap      - use heap for subinterval queue\n");
        flint_printf("-verbose   - show information\n");
        flint_printf("-verbose2  - show more information\n");
        flint_printf("-deg n     - use quadrature degree up to n\n");
        flint_printf("-eval n    - limit number of function evaluations to n\n");
        flint_printf("-depth n   - limit subinterval queue size to n\n\n");
        flint_printf("Implemented integrals:\n");
        for (integral = 0; integral < NUM_INTEGRALS; integral++)
            flint_printf("I%d = %s\n", integral, descr[integral]);
        flint_printf("\n");
        return 1;
    }

    acb_calc_integrate_opt_init(options);

    prec = 64;
    twice = 0;
    goal = 0;
    havetol = havegoal = 0;

    acb_init(a);
    acb_init(b);
    acb_init(s);
    acb_init(t);
    mag_init(tol);

    for (i = 1; i < argc; i++)
    {
        if (!strcmp(argv[i], "-prec"))
        {
            prec = atol(argv[i+1]);
        }
        else if (!strcmp(argv[i], "-twice"))
        {
            twice = 1;
        }
        else if (!strcmp(argv[i], "-goal"))
        {
            goal = atol(argv[i+1]);
            if (goal < 0)
            {
                flint_printf("expected goal >= 0\n");
                return 1;
            }
            havegoal = 1;
        }
        else if (!strcmp(argv[i], "-tol"))
        {
            arb_t x;
            arb_init(x);
            arb_set_str(x, argv[i+1], 10);
            arb_get_mag(tol, x);
            arb_clear(x);
            havetol = 1;
        }
        else if (!strcmp(argv[i], "-deg"))
        {
            options->deg_limit = atol(argv[i+1]);
        }
        else if (!strcmp(argv[i], "-eval"))
        {
            options->eval_limit = atol(argv[i+1]);
        }
        else if (!strcmp(argv[i], "-depth"))
        {
            options->depth_limit = atol(argv[i+1]);
        }
        else if (!strcmp(argv[i], "-verbose"))
        {
            options->verbose = 1;
        }
        else if (!strcmp(argv[i], "-verbose2"))
        {
            options->verbose = 2;
        }
        else if (!strcmp(argv[i], "-heap"))
        {
            options->use_heap = 1;
        }
    }

    if (!havegoal)
        goal = prec;

    if (!havetol)
        mag_set_ui_2exp_si(tol, 1, -prec);

    for (integral = ifrom; integral <= ito; integral++)
    {
        flint_printf("I%d = %s ...\n", integral, descr[integral]);

        for (i = 0; i < 1 + twice; i++)
        {
            TIMEIT_ONCE_START
            switch (integral)
            {
            case 0:
                acb_set_d(a, 0);
                acb_set_d(b, 100);
                acb_calc_integrate(s, f_sin, NULL, a, b, goal, tol, options, prec);
                break;

            case 1:
                acb_set_d(a, 0);
                acb_set_d(b, 1);
                acb_calc_integrate(s, f_atanderiv, NULL, a, b, goal, tol, options, prec);
                acb_mul_2exp_si(s, s, 2);
                break;

            case 2:
                acb_set_d(a, 0);
                acb_one(b);
                acb_mul_2exp_si(b, b, goal);
                acb_calc_integrate(s, f_atanderiv, NULL, a, b, goal, tol, options, prec);
                arb_add_error_2exp_si(acb_realref(s), -goal);
                acb_mul_2exp_si(s, s, 1);
                break;

            case 3:
                acb_set_d(a, 0);
                acb_set_d(b, 1);
                acb_calc_integrate(s, f_circle, NULL, a, b, goal, tol, options, prec);
                acb_mul_2exp_si(s, s, 2);
                break;

            case 4:
                acb_set_d(a, 0);
                acb_set_d(b, 8);
                acb_calc_integrate(s, f_rump, NULL, a, b, goal, tol, options, prec);
                break;

            case 5:
                acb_set_d(a, 1);
                acb_set_d(b, 101);
                acb_calc_integrate(s, f_floor, NULL, a, b, goal, tol, options, prec);
                break;

            case 6:
                acb_set_d(a, 0);
                acb_set_d(b, 1);
                acb_calc_integrate(s, f_helfgott, NULL, a, b, goal, tol, options, prec);
                break;

            case 7:
                acb_zero(s);

                acb_set_d_d(a, -1.0, -1.0);
                acb_set_d_d(b, 2.0, -1.0);
                acb_calc_integrate(t, f_zeta, NULL, a, b, goal, tol, options, prec);
                acb_add(s, s, t, prec);

                acb_set_d_d(a, 2.0, -1.0);
                acb_set_d_d(b, 2.0, 1.0);
                acb_calc_integrate(t, f_zeta, NULL, a, b, goal, tol, options, prec);
                acb_add(s, s, t, prec);

                acb_set_d_d(a, 2.0, 1.0);
                acb_set_d_d(b, -1.0, 1.0);
                acb_calc_integrate(t, f_zeta, NULL, a, b, goal, tol, options, prec);
                acb_add(s, s, t, prec);

                acb_set_d_d(a, -1.0, 1.0);
                acb_set_d_d(b, -1.0, -1.0);
                acb_calc_integrate(t, f_zeta, NULL, a, b, goal, tol, options, prec);
                acb_add(s, s, t, prec);

                acb_const_pi(t, prec);
                acb_div(s, s, t, prec);
                acb_mul_2exp_si(s, s, -1);
                acb_div_onei(s, s);
                break;

            case 8:
                acb_set_d(a, 0);
                acb_set_d(b, 1);
                acb_calc_integrate(s, f_essing, NULL, a, b, goal, tol, options, prec);
                break;

            case 9:
                acb_set_d(a, 0);
                acb_set_d(b, 1);
                acb_calc_integrate(s, f_essing2, NULL, a, b, goal, tol, options, prec);
                break;

            case 10:
                acb_set_d(a, 0);
                acb_set_d(b, 10000);
                acb_calc_integrate(s, f_factorial1000, NULL, a, b, goal, tol, options, prec);
                break;

            case 11:
                acb_set_d_d(a, 1.0, 0.0);
                acb_set_d_d(b, 1.0, 1000.0);
                acb_calc_integrate(s, f_gamma, NULL, a, b, goal, tol, options, prec);
                break;

            case 12:
                acb_set_d(a, -10.0);
                acb_set_d(b, 10.0);
                acb_calc_integrate(s, f_sin_plus_small, NULL, a, b, goal, tol, options, prec);
                break;

            case 13:
                acb_set_d(a, -1020.0);
                acb_set_d(b, -1010.0);
                acb_calc_integrate(s, f_exp, NULL, a, b, goal, tol, options, prec);
                break;

            case 14:
                acb_set_d(a, 0);
                acb_set_d(b, ceil(sqrt(goal * 0.693147181) + 1.0));
                acb_calc_integrate(s, f_gaussian, NULL, a, b, goal, tol, options, prec);
                acb_mul(b, b, b, prec);
                acb_neg(b, b);
                acb_exp(b, b, prec);
                arb_add_error(acb_realref(s), acb_realref(b));
                break;

            case 15:
                acb_set_d(a, 0.0);
                acb_set_d(b, 1.0);
                acb_calc_integrate(s, f_spike, NULL, a, b, goal, tol, options, prec);
                break;

            case 16:
                acb_set_d(a, 0.0);
                acb_set_d(b, 8.0);
                acb_calc_integrate(s, f_monster, NULL, a, b, goal, tol, options, prec);
                break;

            case 17:
                acb_set_d(a, 0);
                acb_set_d(b, ceil(goal * 0.693147181 + 1.0));
                acb_calc_integrate(s, f_sech, NULL, a, b, goal, tol, options, prec);
                acb_neg(b, b);
                acb_exp(b, b, prec);
                acb_mul_2exp_si(b, b, 1);
                arb_add_error(acb_realref(s), acb_realref(b));
                break;

            case 18:
                acb_set_d(a, 0);
                acb_set_d(b, ceil(goal * 0.693147181 / 3.0 + 2.0));
                acb_calc_integrate(s, f_sech3, NULL, a, b, goal, tol, options, prec);
                acb_neg(b, b);
                acb_mul_ui(b, b, 3, prec);
                acb_exp(b, b, prec);
                acb_mul_2exp_si(b, b, 3);
                acb_div_ui(b, b, 3, prec);
                arb_add_error(acb_realref(s), acb_realref(b));
                break;

            case 19:
                if (goal < 0)
                    abort();
                /* error bound 2^-N (1+N) when truncated at 2^-N */
                N = goal + FLINT_BIT_COUNT(goal);
                acb_one(a);
                acb_mul_2exp_si(a, a, -N);
                acb_one(b);
                acb_calc_integrate(s, f_log_div1p, NULL, a, b, goal, tol, options, prec);
                acb_set_ui(b, N + 1);
                acb_mul_2exp_si(b, b, -N);
                arb_add_error(acb_realref(s), acb_realref(b));
                break;

           case 20:
                if (goal < 0)
                    abort();
                /* error bound (N+1) exp(-N) when truncated at N */
                N = goal + FLINT_BIT_COUNT(goal);
                acb_zero(a);
                acb_set_ui(b, N);
                acb_calc_integrate(s, f_log_div1p_transformed, NULL, a, b, goal, tol, options, prec);
                acb_neg(b, b);
                acb_exp(b, b, prec);
                acb_mul_ui(b, b, N + 1, prec);
                arb_add_error(acb_realref(s), acb_realref(b));
                break;

            case 21:

                acb_zero(s);

                N = 10;

                acb_set_d_d(a, 0.5, -0.5);
                acb_set_d_d(b, 0.5, 0.5);
                acb_calc_integrate(t, f_elliptic_p_laurent_n, &N, a, b, goal, tol, options, prec);
                acb_add(s, s, t, prec);

                acb_set_d_d(a, 0.5, 0.5);
                acb_set_d_d(b, -0.5, 0.5);
                acb_calc_integrate(t, f_elliptic_p_laurent_n, &N, a, b, goal, tol, options, prec);
                acb_add(s, s, t, prec);

                acb_set_d_d(a, -0.5, 0.5);
                acb_set_d_d(b, -0.5, -0.5);
                acb_calc_integrate(t, f_elliptic_p_laurent_n, &N, a, b, goal, tol, options, prec);
                acb_add(s, s, t, prec);

                acb_set_d_d(a, -0.5, -0.5);
                acb_set_d_d(b, 0.5, -0.5);
                acb_calc_integrate(t, f_elliptic_p_laurent_n, &N, a, b, goal, tol, options, prec);
                acb_add(s, s, t, prec);

                acb_const_pi(t, prec);
                acb_div(s, s, t, prec);
                acb_mul_2exp_si(s, s, -1);
                acb_div_onei(s, s);
                break;

            case 22:

                acb_zero(s);

                N = 1000;

                acb_set_d_d(a, 100.0, 0.0);
                acb_set_d_d(b, 100.0, N);
                acb_calc_integrate(t, f_zeta_frac, NULL, a, b, goal, tol, options, prec);
                acb_add(s, s, t, prec);

                acb_set_d_d(a, 100, N);
                acb_set_d_d(b, 0.5, N);
                acb_calc_integrate(t, f_zeta_frac, NULL, a, b, goal, tol, options, prec);
                acb_add(s, s, t, prec);

                acb_div_onei(s, s);
                arb_zero(acb_imagref(s));

                acb_set_ui(t, N);
                acb_dirichlet_hardy_theta(t, t, NULL, NULL, 1, prec);
                acb_add(s, s, t, prec);

                acb_const_pi(t, prec);
                acb_div(s, s, t, prec);
                acb_add_ui(s, s, 1, prec);
                break;

            case 23:
                acb_set_d(a, 0.0);
                acb_set_d(b, 1000.0);
                acb_calc_integrate(s, f_lambertw, NULL, a, b, goal, tol, options, prec);
                break;

            case 24:
                acb_set_d(a, 0.0);
                acb_const_pi(b, prec);
                acb_calc_integrate(s, f_max_sin_cos, NULL, a, b, goal, tol, options, prec);
                break;

            case 25:
                acb_set_si(a, -1);
                acb_set_si(b, 1);
                acb_calc_integrate(s, f_erf_bent, NULL, a, b, goal, tol, options, prec);
                break;

            case 26:
                acb_set_si(a, -10);
                acb_set_si(b, 10);
                acb_calc_integrate(s, f_airy_ai, NULL, a, b, goal, tol, options, prec);
                break;

            case 27:
                acb_set_si(a, 0);
                acb_set_si(b, 10);
                acb_calc_integrate(s, f_horror, NULL, a, b, goal, tol, options, prec);
                break;

            case 28:
                acb_set_d_d(a, -1, -1);
                acb_set_d_d(b, -1, 1);
                acb_calc_integrate(s, f_sqrt, NULL, a, b, goal, tol, options, prec);
                break;

            case 29:
                acb_set_d(a, 0);
                acb_set_d(b, ceil(sqrt(goal * 0.693147181) + 1.0));
                acb_calc_integrate(s, f_gaussian_twist, NULL, a, b, goal, tol, options, prec);
                acb_mul(b, b, b, prec);
                acb_neg(b, b);
                acb_exp(b, b, prec);
                arb_add_error(acb_realref(s), acb_realref(b));
                arb_add_error(acb_imagref(s), acb_realref(b));
                break;

            case 30:
                acb_set_d(a, 0);
                acb_set_d(b, ceil(goal * 0.693147181 + 1.0));
                acb_calc_integrate(s, f_exp_airy, NULL, a, b, goal, tol, options, prec);
                acb_neg(b, b);
                acb_exp(b, b, prec);
                acb_mul_2exp_si(b, b, 1);
                arb_add_error(acb_realref(s), acb_realref(b));
                break;

            case 31:
                acb_zero(a);
                acb_const_pi(b, prec);
                acb_calc_integrate(s, f_sin_cos_frac, NULL, a, b, goal, tol, options, prec);
                break;

            case 32:
                acb_zero(a);
                acb_set_ui(b, 3);
                acb_calc_integrate(s, f_sin_near_essing, NULL, a, b, goal, tol, options, prec);
                break;

            case 33:
                acb_zero(a);
                acb_zero(b);
                k = 3;
                scaled_bessel_select_N(acb_realref(b), k, prec);
                acb_calc_integrate(s, f_scaled_bessel, &k, a, b, goal, tol, options, prec);
                scaled_bessel_tail_bound(acb_realref(a), k, acb_realref(b), prec);
                arb_add_error(acb_realref(s), acb_realref(a));
                break;

            case 34:
                acb_zero(a);
                acb_zero(b);
                k = 15;
                scaled_bessel_select_N(acb_realref(b), k, prec);
                acb_calc_integrate(s, f_scaled_bessel, &k, a, b, goal, tol, options, prec);
                scaled_bessel_tail_bound(acb_realref(a), k, acb_realref(b), prec);
                arb_add_error(acb_realref(s), acb_realref(a));
                break;

            case 35:
                acb_set_d_d(a, -1, -1);
                acb_set_d_d(b, -1, 1);
                acb_calc_integrate(s, f_rsqrt, NULL, a, b, goal, tol, options, prec);
                break;

            default:
                abort();
            }

            TIMEIT_ONCE_STOP
        }
        flint_printf("I%d = ", integral);
        acb_printn(s, 3.333 * prec, 0);
        flint_printf("\n\n");
    }

    acb_clear(a);
    acb_clear(b);
    acb_clear(s);
    acb_clear(t);
    mag_clear(tol);

    flint_cleanup();
    return 0;
}
Exemple #5
0
static void
evaluate(acb_poly_t A, acb_srcptr a, slong p, const acb_t z, slong n, slong prec)
{
    acb_poly_fit_length(A, p + 1);

    if (p == 1)
    {
        acb_add_ui(A->coeffs, a, n, prec);
        if (z != NULL)
            acb_mul(A->coeffs, A->coeffs, z, prec);
    }
    else if (p == 2)
    {
        acb_add(A->coeffs, a + 0, a + 1, prec);
        acb_add_ui(A->coeffs + 1, A->coeffs, 2 * n, prec);
        acb_add_ui(A->coeffs, A->coeffs, n, prec);
        acb_mul_ui(A->coeffs, A->coeffs, n, prec);
        acb_addmul(A->coeffs, a + 0, a + 1, prec);
        if (z != NULL)
        {
            acb_mul(A->coeffs, A->coeffs, z, prec);
            acb_mul(A->coeffs + 1, A->coeffs + 1, z, prec);
        }
    }
    else if (p == 3)
    {
        acb_t t, u;
        acb_init(t);
        acb_init(u);

        acb_add(t, a + 0, a + 1, prec);
        acb_add(t, t, a + 2, prec);

        acb_mul(u, a + 0, a + 1, prec);
        acb_mul(A->coeffs, u, a + 2, prec);

        acb_addmul(u, a + 0, a + 2, prec);
        acb_addmul(u, a + 1, a + 2, prec);

        /*
        (a0 + n)(a1 + n)(a2 + n) = a0 a1 a2 + (a0 a1 + a0 a2 + a1 a2) n + (a0 + a1 + a2) n^2 + n^3
        (a0 a1 + a0 a2 + a1 a2) + 2 (a0 + a1 + a2) n + 3 n^2
        (a0 + a1 + a2) + 3n
        1
        */

        acb_addmul_ui(A->coeffs, u, n, prec);
        acb_addmul_ui(A->coeffs, t, n * n, prec);
        acb_add_ui(A->coeffs, A->coeffs, n * n * n, prec);

        acb_set(A->coeffs + 1, u);
        acb_addmul_ui(A->coeffs + 1, t, 2 * n, prec);
        acb_add_ui(A->coeffs + 1, A->coeffs + 1, 3 * n * n, prec);

        acb_add_ui(A->coeffs + 2, t, 3 * n, prec);

        if (z != NULL)
        {
            acb_mul(A->coeffs + 0, A->coeffs + 0, z, prec);
            acb_mul(A->coeffs + 1, A->coeffs + 1, z, prec);
            acb_mul(A->coeffs + 2, A->coeffs + 2, z, prec);
        }

        acb_clear(t);
        acb_clear(u);
    }
    else if (p != 0)
    {
        flint_abort();
    }

    if (z != NULL)
        acb_set(A->coeffs + p, z);
    else
        acb_one(A->coeffs + p);

    _acb_poly_set_length(A, p + 1);
    _acb_poly_normalise(A);
}
int main()
{
    long iter;
    flint_rand_t state;

    printf("elliptic_p_zpx....");
    fflush(stdout);

    flint_randinit(state);

    /* Test differential equation */
    for (iter = 0; iter < 5000; iter++)
    {
        acb_t tau, z;
        acb_ptr g, wp, wp3, wpd, wpd2;
        long prec, len, i;

        len = 1 + n_randint(state, 15);
        prec = 2 + n_randint(state, 1000);

        acb_init(tau);
        acb_init(z);
        g = _acb_vec_init(2);
        wp = _acb_vec_init(len + 1);
        wp3 = _acb_vec_init(len);
        wpd = _acb_vec_init(len);
        wpd2 = _acb_vec_init(len);

        acb_randtest(tau, state, prec, 10);
        acb_randtest(z, state, prec, 10);

        acb_modular_elliptic_p_zpx(wp, z, tau, len + 1, prec);
        acb_modular_eisenstein(g, tau, 2, prec);
        acb_mul_ui(g, g, 60, prec);
        acb_mul_ui(g + 1, g + 1, 140, prec);

        _acb_poly_derivative(wpd, wp, len + 1, prec);
        _acb_poly_mullow(wpd2, wpd, len, wpd, len, len, prec);
        _acb_poly_pow_ui_trunc_binexp(wp3, wp, len, 3, len, prec);
        _acb_vec_scalar_mul_ui(wp3, wp3, len, 4, prec);
        _acb_vec_scalar_submul(wp3, wp, len, g, prec);
        acb_sub(wp3, wp3, g + 1, prec);

        for (i = 0; i < len; i++)
        {
            if (!acb_overlaps(wpd2 + i, wp3 + i))
            {
                printf("FAIL (overlap)\n");
                printf("i = %ld  len = %ld  prec = %ld\n\n", i, len, prec);
                printf("z = "); acb_printd(z, 15); printf("\n\n");
                printf("tau = "); acb_printd(tau, 15); printf("\n\n");
                printf("wp = "); acb_printd(wp + i, 15); printf("\n\n");
                printf("wpd = "); acb_printd(wpd + i, 15); printf("\n\n");
                printf("wp3 = "); acb_printd(wp3 + i, 15); printf("\n\n");
                abort();
            }
        }

        acb_clear(tau);
        acb_clear(z);
        _acb_vec_clear(g, 2);
        _acb_vec_clear(wp, len + 1);
        _acb_vec_clear(wp3, len);
        _acb_vec_clear(wpd, len);
        _acb_vec_clear(wpd2, len);
    }

    /* Consistency test */
    for (iter = 0; iter < 5000; iter++)
    {
        acb_t tau, z;
        acb_ptr wp1, wp2;
        long prec1, prec2, len1, len2, i;

        len1 = n_randint(state, 15);
        len2 = n_randint(state, 15);
        prec1 = 2 + n_randint(state, 1000);
        prec2 = 2 + n_randint(state, 1000);

        acb_init(tau);
        acb_init(z);
        wp1 = _acb_vec_init(len1);
        wp2 = _acb_vec_init(len2);

        acb_randtest(tau, state, prec1, 10);
        acb_randtest(z, state, prec1, 10);

        acb_modular_elliptic_p_zpx(wp1, z, tau, len1, prec1);
        acb_modular_elliptic_p_zpx(wp2, z, tau, len2, prec2);

        for (i = 0; i < FLINT_MIN(len1, len2); i++)
        {
            if (!acb_overlaps(wp1 + i, wp2 + i))
            {
                printf("FAIL (overlap)\n");
                printf("i = %ld len1 = %ld len2 = %ld\n\n", i, len1, len2);
                printf("tau = "); acb_printd(tau, 15); printf("\n\n");
                printf("z = "); acb_printd(z, 15); printf("\n\n");
                printf("wp1 = "); acb_printd(wp1 + i, 15); printf("\n\n");
                printf("wp2 = "); acb_printd(wp2 + i, 15); printf("\n\n");
                abort();
            }
        }

        acb_clear(tau);
        acb_clear(z);
        _acb_vec_clear(wp1, len1);
        _acb_vec_clear(wp2, len2);
    }

    flint_randclear(state);
    flint_cleanup();
    printf("PASS\n");
    return EXIT_SUCCESS;
}
Exemple #7
0
int main()
{
    slong iter;
    flint_rand_t state;

    flint_printf("bernoulli_poly_ui....");
    fflush(stdout);

    flint_randinit(state);

    /* test multiplication theorem */
    for (iter = 0; iter < 1000 * arb_test_multiplier(); iter++)
    {
        acb_t x, t, res1, res2;
        ulong n, m, k;
        slong prec;

        n = n_randint(state, 50);
        m = 1 + n_randint(state, 5);
        prec = 2 + n_randint(state, 200);

        acb_init(x);
        acb_init(t);
        acb_init(res1);
        acb_init(res2);

        acb_randtest(x, state, 2 + n_randint(state, 200), 20);
        acb_randtest(res1, state, 2 + n_randint(state, 200), 20);

        acb_mul_ui(t, x, m, prec);
        acb_bernoulli_poly_ui(res1, n, t, prec);

        acb_zero(res2);
        for (k = 0; k < m; k++)
        {
            acb_set_ui(t, k);
            acb_div_ui(t, t, m, prec);
            acb_add(t, t, x, prec);
            acb_bernoulli_poly_ui(t, n, t, prec);
            acb_add(res2, res2, t, prec);
        }

        if (n > 0)
        {
            arb_ui_pow_ui(acb_realref(t), m, n - 1, prec);
            acb_mul_arb(res2, res2, acb_realref(t), prec);
        }
        else
        {
            acb_div_ui(res2, res2, m, prec);
        }

        if (!acb_overlaps(res1, res2))
        {
            flint_printf("FAIL: overlap\n\n");
            flint_printf("n = %wu, m = %wu\n\n", n, m);
            flint_printf("x = "); acb_printd(x, 15); flint_printf("\n\n");
            flint_printf("res1 = "); acb_printd(res1, 15); flint_printf("\n\n");
            flint_printf("res2 = "); acb_printd(res2, 15); flint_printf("\n\n");
            abort();
        }

        acb_clear(x);
        acb_clear(t);
        acb_clear(res1);
        acb_clear(res2);
    }

    flint_randclear(state);
    flint_cleanup();
    flint_printf("PASS\n");
    return EXIT_SUCCESS;
}
Exemple #8
0
void
acb_modular_elliptic_k_cpx(acb_ptr w, const acb_t m, slong len, slong prec)
{
    acb_t t, u, msub1m, m2sub1;
    slong k, n;

    if (len < 1)
        return;

    if (len == 1)
    {
        acb_modular_elliptic_k(w, m, prec);
        return;
    }

    if (acb_is_zero(m))
    {
        acb_const_pi(w, prec);
        acb_mul_2exp_si(w, w, -1);

        for (k = 1; k < len; k++)
        {
            acb_mul_ui(w + k, w + k - 1, (2 * k - 1) * (2 * k - 1), prec);
            acb_div_ui(w + k, w + k, 4 * k * k, prec);
        }

        return;
    }

    acb_init(t);
    acb_init(u);
    acb_init(msub1m);
    acb_init(m2sub1);

    acb_sub_ui(msub1m, m, 1, prec);
    acb_neg(t, msub1m);
    acb_sqrt(t, t, prec);
    acb_mul(msub1m, msub1m, m, prec);

    acb_mul_2exp_si(m2sub1, m, 1);
    acb_sub_ui(m2sub1, m2sub1, 1, prec);

    acb_agm1_cpx(w, t, 2, prec);

    /* pi M'(t) / (4 t M(t)^2) */
    acb_mul(u, w, w, prec);
    acb_mul(t, t, u, prec);
    acb_div(w + 1, w + 1, t, prec);

    acb_const_pi(u, prec);
    acb_mul(w + 1, w + 1, u, prec);
    acb_mul_2exp_si(w + 1, w + 1, -2);

    /* pi / (2 M(t)) */
    acb_const_pi(u, prec);
    acb_div(w, u, w, prec);
    acb_mul_2exp_si(w, w, -1);

    acb_inv(t, msub1m, prec);

    for (k = 2; k < len; k++)
    {
        n = k - 2;

        acb_mul_ui(w + k, w + n, (2 * n + 1) * (2 * n + 1), prec);

        acb_mul(u, w + n + 1, m2sub1, prec);
        acb_addmul_ui(w + k, u, (n + 1) * (n + 1) * 4, prec);

        acb_mul(w + k, w + k, t, prec);
        acb_div_ui(w + k, w + k, 4 * (n + 1) * (n + 2), prec);
        acb_neg(w + k, w + k);
    }

    acb_clear(t);
    acb_clear(u);
    acb_clear(msub1m);
    acb_clear(m2sub1);
}
Exemple #9
0
/* 
   F(x)  = c0   +     c1 x    +     c2 x^2   +     c3 x^3    +    [...]
   F'(x) = c1   +   2 c2 x    +   3 c3 x^2   +   4 c4 x^3    +    [...]
*/
static void
evaluate_sum(acb_t res, acb_t res1,
    const acb_t a, const acb_t b, const acb_t c, const acb_t y,
    const acb_t x, const acb_t f0, const acb_t f1, long num, long prec)
{
    acb_t s, s2, w, d, e, xpow, ck, cknext;
    long k;

    acb_init(s);
    acb_init(s2);
    acb_init(w);
    acb_init(d);
    acb_init(e);
    acb_init(xpow);
    acb_init(ck);
    acb_init(cknext);

    /* d = (y-1)*y */
    acb_sub_ui(d, y, 1, prec);
    acb_mul(d, d, y, prec);
    acb_one(xpow);

    for (k = 0; k < num; k++)
    {
        if (k == 0)
        {
            acb_set(ck, f0);
            acb_set(cknext, f1);
        }
        else
        {
            acb_add_ui(w, b, k-1, prec);
            acb_mul(w, w, ck, prec);
            acb_add_ui(e, a, k-1, prec);
            acb_mul(w, w, e, prec);

            acb_add(e, a, b, prec);
            acb_add_ui(e, e, 2*(k+1)-3, prec);
            acb_mul(e, e, y, prec);
            acb_sub(e, e, c, prec);
            acb_sub_ui(e, e, k-1, prec);
            acb_mul_ui(e, e, k, prec);
            acb_addmul(w, e, cknext, prec);

            acb_mul_ui(e, d, k+1, prec);
            acb_mul_ui(e, e, k, prec);
            acb_div(w, w, e, prec);
            acb_neg(w, w);

            acb_set(ck, cknext);
            acb_set(cknext, w);
        }

        acb_addmul(s, ck, xpow, prec);
        acb_mul_ui(w, cknext, k+1, prec);
        acb_addmul(s2, w, xpow, prec);

        acb_mul(xpow, xpow, x, prec);
    }

    acb_set(res, s);
    acb_set(res1, s2);

    acb_clear(s);
    acb_clear(s2);
    acb_clear(w);
    acb_clear(d);
    acb_clear(e);
    acb_clear(xpow);
    acb_clear(ck);
    acb_clear(cknext);
}
int
acb_calc_integrate_taylor(acb_t res,
    acb_calc_func_t func, void * param,
    const acb_t a, const acb_t b,
    const arf_t inner_radius,
    const arf_t outer_radius,
    long accuracy_goal, long prec)
{
    long num_steps, step, N, bp;
    int result;

    acb_t delta, m, x, y1, y2, sum;
    acb_ptr taylor_poly;
    arf_t err;

    acb_init(delta);
    acb_init(m);
    acb_init(x);
    acb_init(y1);
    acb_init(y2);
    acb_init(sum);
    arf_init(err);

    acb_sub(delta, b, a, prec);

    /* precision used for bounds calculations */
    bp = MAG_BITS;

    /* compute the number of steps */
    {
        arf_t t;
        arf_init(t);
        acb_get_abs_ubound_arf(t, delta, bp);
        arf_div(t, t, inner_radius, bp, ARF_RND_UP);
        arf_mul_2exp_si(t, t, -1);
        num_steps = (long) (arf_get_d(t, ARF_RND_UP) + 1.0);
        /* make sure it's not something absurd */
        num_steps = FLINT_MIN(num_steps, 10 * prec);
        num_steps = FLINT_MAX(num_steps, 1);
        arf_clear(t);
    }

    result = ARB_CALC_SUCCESS;

    acb_zero(sum);

    for (step = 0; step < num_steps; step++)
    {
        /* midpoint of subinterval */
        acb_mul_ui(m, delta, 2 * step + 1, prec);
        acb_div_ui(m, m, 2 * num_steps, prec);
        acb_add(m, m, a, prec);

        if (arb_calc_verbose)
        {
            printf("integration point %ld/%ld: ", 2 * step + 1, 2 * num_steps);
            acb_printd(m, 15); printf("\n");
        }

        /* evaluate at +/- x */
        /* TODO: exactify m, and include error in x? */
        acb_div_ui(x, delta, 2 * num_steps, prec);

        /* compute bounds and number of terms to use */
        {
            arb_t cbound, xbound, rbound;
            arf_t C, D, R, X, T;
            double DD, TT, NN;

            arb_init(cbound);
            arb_init(xbound);
            arb_init(rbound);
            arf_init(C);
            arf_init(D);
            arf_init(R);
            arf_init(X);
            arf_init(T);

            /* R is the outer radius */
            arf_set(R, outer_radius);

            /* X = upper bound for |x| */
            acb_get_abs_ubound_arf(X, x, bp);
            arb_set_arf(xbound, X);

            /* Compute C(m,R). Important subtlety: due to rounding when
               computing m, we will in general be farther than R away from
               the integration path. But since acb_calc_cauchy_bound
               actually integrates over the area traced by a complex
               interval, it will catch any extra singularities (giving
               an infinite bound). */
            arb_set_arf(rbound, outer_radius);
            acb_calc_cauchy_bound(cbound, func, param, m, rbound, 8, bp);
            arf_set_mag(C, arb_radref(cbound));
            arf_add(C, arb_midref(cbound), C, bp, ARF_RND_UP);

            /* Sanity check: we need C < inf and R > X */
            if (arf_is_finite(C) && arf_cmp(R, X) > 0)
            {
                /* Compute upper bound for D = C * R * X / (R - X) */
                arf_mul(D, C, R, bp, ARF_RND_UP);
                arf_mul(D, D, X, bp, ARF_RND_UP);
                arf_sub(T, R, X, bp, ARF_RND_DOWN);
                arf_div(D, D, T, bp, ARF_RND_UP);

                /* Compute upper bound for T = (X / R) */
                arf_div(T, X, R, bp, ARF_RND_UP);

                /* Choose N */
                /* TODO: use arf arithmetic to avoid overflow */
                /* TODO: use relative accuracy (look at |f(m)|?) */
                DD = arf_get_d(D, ARF_RND_UP);
                TT = arf_get_d(T, ARF_RND_UP);
                NN = -(accuracy_goal * 0.69314718055994530942 + log(DD)) / log(TT);
                N = NN + 0.5;
                N = FLINT_MIN(N, 100 * prec);
                N = FLINT_MAX(N, 1);

                /* Tail bound: D / (N + 1) * T^N */
                {
                    mag_t TT;
                    mag_init(TT);
                    arf_get_mag(TT, T);
                    mag_pow_ui(TT, TT, N);
                    arf_set_mag(T, TT);
                    mag_clear(TT);
                }
                arf_mul(D, D, T, bp, ARF_RND_UP);
                arf_div_ui(err, D, N + 1, bp, ARF_RND_UP);
            }
            else
            {
                N = 1;
                arf_pos_inf(err);
                result = ARB_CALC_NO_CONVERGENCE;
            }

            if (arb_calc_verbose)
            {
                printf("N = %ld; bound: ", N); arf_printd(err, 15); printf("\n");
                printf("R: "); arf_printd(R, 15); printf("\n");
                printf("C: "); arf_printd(C, 15); printf("\n");
                printf("X: "); arf_printd(X, 15); printf("\n");
            }

            arb_clear(cbound);
            arb_clear(xbound);
            arb_clear(rbound);
            arf_clear(C);
            arf_clear(D);
            arf_clear(R);
            arf_clear(X);
            arf_clear(T);
        }

        /* evaluate Taylor polynomial */
        taylor_poly = _acb_vec_init(N + 1);
        func(taylor_poly, m, param, N, prec);
        _acb_poly_integral(taylor_poly, taylor_poly, N + 1, prec);
        _acb_poly_evaluate(y2, taylor_poly, N + 1, x, prec);
        acb_neg(x, x);
        _acb_poly_evaluate(y1, taylor_poly, N + 1, x, prec);
        acb_neg(x, x);

        /* add truncation error */
        arb_add_error_arf(acb_realref(y1), err);
        arb_add_error_arf(acb_imagref(y1), err);
        arb_add_error_arf(acb_realref(y2), err);
        arb_add_error_arf(acb_imagref(y2), err);

        acb_add(sum, sum, y2, prec);
        acb_sub(sum, sum, y1, prec);

        if (arb_calc_verbose)
        {
            printf("values:  ");
            acb_printd(y1, 15); printf("  ");
            acb_printd(y2, 15); printf("\n");
        }

        _acb_vec_clear(taylor_poly, N + 1);

        if (result == ARB_CALC_NO_CONVERGENCE)
            break;
    }

    acb_set(res, sum);

    acb_clear(delta);
    acb_clear(m);
    acb_clear(x);
    acb_clear(y1);
    acb_clear(y2);
    acb_clear(sum);
    arf_clear(err);

    return result;
}