acb_poly.h – polynomials over the complex numbers¶
An acb_poly_t represents a polynomial over the complex numbers,
implemented as an array of coefficients of type acb_struct.
Most functions are provided in two versions: an underscore method which operates directly on pre-allocated arrays of coefficients and generally has some restrictions (such as requiring the lengths to be nonzero and not supporting aliasing of the input and output arrays), and a non-underscore method which performs automatic memory management and handles degenerate cases.
Types, macros and constants¶
-
type
acb_poly_struct¶
-
type
acb_poly_t¶ Contains a pointer to an array of coefficients (coeffs), the used length (length), and the allocated size of the array (alloc).
An acb_poly_t is defined as an array of length one of type acb_poly_struct, permitting an acb_poly_t to be passed by reference.
Memory management¶
-
void
acb_poly_init(acb_poly_t poly)¶ Initializes the polynomial for use, setting it to the zero polynomial.
-
void
acb_poly_clear(acb_poly_t poly)¶ Clears the polynomial, deallocating all coefficients and the coefficient array.
-
void
acb_poly_fit_length(acb_poly_t poly, slong len)¶ Makes sure that the coefficient array of the polynomial contains at least len initialized coefficients.
-
void
_acb_poly_set_length(acb_poly_t poly, slong len)¶ Directly changes the length of the polynomial, without allocating or deallocating coefficients. The value should not exceed the allocation length.
-
void
_acb_poly_normalise(acb_poly_t poly)¶ Strips any trailing coefficients which are identical to zero.
-
void
acb_poly_swap(acb_poly_t poly1, acb_poly_t poly2)¶ Swaps poly1 and poly2 efficiently.
-
slong
acb_poly_allocated_bytes(const acb_poly_t x)¶ Returns the total number of bytes heap-allocated internally by this object. The count excludes the size of the structure itself. Add
sizeof(acb_poly_struct)to get the size of the object as a whole.
Basic properties and manipulation¶
-
slong
acb_poly_length(const acb_poly_t poly)¶ Returns the length of poly, i.e. zero if poly is identically zero, and otherwise one more than the index of the highest term that is not identically zero.
-
slong
acb_poly_degree(const acb_poly_t poly)¶ Returns the degree of poly, defined as one less than its length. Note that if one or several leading coefficients are balls containing zero, this value can be larger than the true degree of the exact polynomial represented by poly, so the return value of this function is effectively an upper bound.
-
int
acb_poly_is_zero(const acb_poly_t poly)¶
-
int
acb_poly_is_one(const acb_poly_t poly)¶
-
int
acb_poly_is_x(const acb_poly_t poly)¶ Returns 1 if poly is exactly the polynomial 0, 1 or x respectively. Returns 0 otherwise.
-
void
acb_poly_zero(acb_poly_t poly)¶ Sets poly to the zero polynomial.
-
void
acb_poly_one(acb_poly_t poly)¶ Sets poly to the constant polynomial 1.
-
void
acb_poly_set(acb_poly_t dest, const acb_poly_t src)¶ Sets dest to a copy of src.
-
void
acb_poly_set_round(acb_poly_t dest, const acb_poly_t src, slong prec)¶ Sets dest to a copy of src, rounded to prec bits.
-
void
acb_poly_set_trunc(acb_poly_t dest, const acb_poly_t src, slong n)¶
-
void
acb_poly_set_trunc_round(acb_poly_t dest, const acb_poly_t src, slong n, slong prec)¶ Sets dest to a copy of src, truncated to length n and rounded to prec bits.
-
void
acb_poly_set_coeff_si(acb_poly_t poly, slong n, slong c)¶
-
void
acb_poly_set_coeff_acb(acb_poly_t poly, slong n, const acb_t c)¶ Sets the coefficient with index n in poly to the value c. We require that n is nonnegative.
-
void
acb_poly_get_coeff_acb(acb_t v, const acb_poly_t poly, slong n)¶ Sets v to the value of the coefficient with index n in poly. We require that n is nonnegative.
-
acb_poly_get_coeff_ptr(poly, n)¶ Given \(n \ge 0\), returns a pointer to coefficient n of poly, or NULL if n exceeds the length of poly.
-
void
_acb_poly_shift_right(acb_ptr res, acb_srcptr poly, slong len, slong n)¶
-
void
acb_poly_shift_right(acb_poly_t res, const acb_poly_t poly, slong n)¶ Sets res to poly divided by \(x^n\), throwing away the lower coefficients. We require that n is nonnegative.
-
void
_acb_poly_shift_left(acb_ptr res, acb_srcptr poly, slong len, slong n)¶
-
void
acb_poly_shift_left(acb_poly_t res, const acb_poly_t poly, slong n)¶ Sets res to poly multiplied by \(x^n\). We require that n is nonnegative.
-
void
acb_poly_truncate(acb_poly_t poly, slong n)¶ Truncates poly to have length at most n, i.e. degree strictly smaller than n. We require that n is nonnegative.
-
slong
acb_poly_valuation(const acb_poly_t poly)¶ Returns the degree of the lowest term that is not exactly zero in poly. Returns -1 if poly is the zero polynomial.
Input and output¶
-
void
acb_poly_printd(const acb_poly_t poly, slong digits)¶ Prints the polynomial as an array of coefficients, printing each coefficient using acb_printd.
-
void
acb_poly_fprintd(FILE *file, const acb_poly_t poly, slong digits)¶ Prints the polynomial as an array of coefficients to the stream file, printing each coefficient using acb_fprintd.
Random generation¶
-
void
acb_poly_randtest(acb_poly_t poly, flint_rand_t state, slong len, slong prec, slong mag_bits)¶ Creates a random polynomial with length at most len.
Comparisons¶
-
int
acb_poly_equal(const acb_poly_t A, const acb_poly_t B)¶ Returns nonzero iff A and B are identical as interval polynomials.
-
int
acb_poly_contains(const acb_poly_t poly1, const acb_poly_t poly2)¶
-
int
acb_poly_contains_fmpz_poly(const acb_poly_t poly1, const fmpz_poly_t poly2)¶
-
int
acb_poly_contains_fmpq_poly(const acb_poly_t poly1, const fmpq_poly_t poly2)¶ Returns nonzero iff poly2 is contained in poly1.
-
int
_acb_poly_overlaps(acb_srcptr poly1, slong len1, acb_srcptr poly2, slong len2)¶
-
int
acb_poly_overlaps(const acb_poly_t poly1, const acb_poly_t poly2)¶ Returns nonzero iff poly1 overlaps with poly2. The underscore function requires that len1 ist at least as large as len2.
-
int
acb_poly_get_unique_fmpz_poly(fmpz_poly_t z, const acb_poly_t x)¶ If x contains a unique integer polynomial, sets z to that value and returns nonzero. Otherwise (if x represents no integers or more than one integer), returns zero, possibly partially modifying z.
-
int
acb_poly_is_real(const acb_poly_t poly)¶ Returns nonzero iff all coefficients in poly have zero imaginary part.
Conversions¶
-
void
acb_poly_set_fmpz_poly(acb_poly_t poly, const fmpz_poly_t re, slong prec)¶
-
void
acb_poly_set2_fmpz_poly(acb_poly_t poly, const fmpz_poly_t re, const fmpz_poly_t im, slong prec)¶
-
void
acb_poly_set_arb_poly(acb_poly_t poly, const arb_poly_t re)¶
-
void
acb_poly_set2_arb_poly(acb_poly_t poly, const arb_poly_t re, const arb_poly_t im)¶
-
void
acb_poly_set_fmpq_poly(acb_poly_t poly, const fmpq_poly_t re, slong prec)¶
-
void
acb_poly_set2_fmpq_poly(acb_poly_t poly, const fmpq_poly_t re, const fmpq_poly_t im, slong prec)¶ Sets poly to the given real part re plus the imaginary part im, both rounded to prec bits.
-
void
acb_poly_set_acb(acb_poly_t poly, const acb_t src)¶
-
void
acb_poly_set_si(acb_poly_t poly, slong src)¶ Sets poly to src.
Bounds¶
-
void
_acb_poly_majorant(arb_ptr res, acb_srcptr poly, slong len, slong prec)¶
-
void
acb_poly_majorant(arb_poly_t res, const acb_poly_t poly, slong prec)¶ Sets res to an exact real polynomial whose coefficients are upper bounds for the absolute values of the coefficients in poly, rounded to prec bits.
Arithmetic¶
-
void
_acb_poly_add(acb_ptr C, acb_srcptr A, slong lenA, acb_srcptr B, slong lenB, slong prec)¶ Sets {C, max(lenA, lenB)} to the sum of {A, lenA} and {B, lenB}. Allows aliasing of the input and output operands.
-
void
acb_poly_add(acb_poly_t C, const acb_poly_t A, const acb_poly_t B, slong prec)¶
-
void
acb_poly_add_si(acb_poly_t C, const acb_poly_t A, slong B, slong prec)¶ Sets C to the sum of A and B.
-
void
_acb_poly_sub(acb_ptr C, acb_srcptr A, slong lenA, acb_srcptr B, slong lenB, slong prec)¶ Sets {C, max(lenA, lenB)} to the difference of {A, lenA} and {B, lenB}. Allows aliasing of the input and output operands.
-
void
acb_poly_sub(acb_poly_t C, const acb_poly_t A, const acb_poly_t B, slong prec)¶ Sets C to the difference of A and B.
-
void
acb_poly_add_series(acb_poly_t C, const acb_poly_t A, const acb_poly_t B, slong len, slong prec)¶ Sets C to the sum of A and B, truncated to length len.
-
void
acb_poly_sub_series(acb_poly_t C, const acb_poly_t A, const acb_poly_t B, slong len, slong prec)¶ Sets C to the difference of A and B, truncated to length len.
-
void
acb_poly_neg(acb_poly_t C, const acb_poly_t A)¶ Sets C to the negation of A.
-
void
acb_poly_scalar_mul_2exp_si(acb_poly_t C, const acb_poly_t A, slong c)¶ Sets C to A multiplied by \(2^c\).
-
void
acb_poly_scalar_mul(acb_poly_t C, const acb_poly_t A, const acb_t c, slong prec)¶ Sets C to A multiplied by c.
-
void
acb_poly_scalar_div(acb_poly_t C, const acb_poly_t A, const acb_t c, slong prec)¶ Sets C to A divided by c.
-
void
_acb_poly_mullow_classical(acb_ptr C, acb_srcptr A, slong lenA, acb_srcptr B, slong lenB, slong n, slong prec)¶
-
void
_acb_poly_mullow_transpose(acb_ptr C, acb_srcptr A, slong lenA, acb_srcptr B, slong lenB, slong n, slong prec)¶
-
void
_acb_poly_mullow_transpose_gauss(acb_ptr C, acb_srcptr A, slong lenA, acb_srcptr B, slong lenB, slong n, slong prec)¶
-
void
_acb_poly_mullow(acb_ptr C, acb_srcptr A, slong lenA, acb_srcptr B, slong lenB, slong n, slong prec)¶ Sets {C, n} to the product of {A, lenA} and {B, lenB}, truncated to length n. The output is not allowed to be aliased with either of the inputs. We require \(\mathrm{lenA} \ge \mathrm{lenB} > 0\), \(n > 0\), \(\mathrm{lenA} + \mathrm{lenB} - 1 \ge n\).
The classical version uses a plain loop.
The transpose version evaluates the product using four real polynomial multiplications (via
_arb_poly_mullow()).The transpose_gauss version evaluates the product using three real polynomial multiplications. This is almost always faster than transpose, but has worse numerical stability when the coefficients vary in magnitude.
The default function
_acb_poly_mullow()automatically switches been classical and transpose multiplication.If the input pointers are identical (and the lengths are the same), they are assumed to represent the same polynomial, and its square is computed.
-
void
acb_poly_mullow_classical(acb_poly_t C, const acb_poly_t A, const acb_poly_t B, slong n, slong prec)¶
-
void
acb_poly_mullow_transpose(acb_poly_t C, const acb_poly_t A, const acb_poly_t B, slong n, slong prec)¶
-
void
acb_poly_mullow_transpose_gauss(acb_poly_t C, const acb_poly_t A, const acb_poly_t B, slong n, slong prec)¶
-
void
acb_poly_mullow(acb_poly_t C, const acb_poly_t A, const acb_poly_t B, slong n, slong prec)¶ Sets C to the product of A and B, truncated to length n. If the same variable is passed for A and B, sets C to the square of A truncated to length n.
-
void
_acb_poly_mul(acb_ptr C, acb_srcptr A, slong lenA, acb_srcptr B, slong lenB, slong prec)¶ Sets {C, lenA + lenB - 1} to the product of {A, lenA} and {B, lenB}. The output is not allowed to be aliased with either of the inputs. We require \(\mathrm{lenA} \ge \mathrm{lenB} > 0\). This function is implemented as a simple wrapper for
_acb_poly_mullow().If the input pointers are identical (and the lengths are the same), they are assumed to represent the same polynomial, and its square is computed.
-
void
acb_poly_mul(acb_poly_t C, const acb_poly_t A1, const acb_poly_t B2, slong prec)¶ Sets C to the product of A and B. If the same variable is passed for A and B, sets C to the square of A.
-
void
_acb_poly_inv_series(acb_ptr Qinv, acb_srcptr Q, slong Qlen, slong len, slong prec)¶ Sets {Qinv, len} to the power series inverse of {Q, Qlen}. Uses Newton iteration.
-
void
acb_poly_inv_series(acb_poly_t Qinv, const acb_poly_t Q, slong n, slong prec)¶ Sets Qinv to the power series inverse of Q.
-
void
_acb_poly_div_series(acb_ptr Q, acb_srcptr A, slong Alen, acb_srcptr B, slong Blen, slong n, slong prec)¶ Sets {Q, n} to the power series quotient of {A, Alen} by {B, Blen}. Uses Newton iteration followed by multiplication.
-
void
acb_poly_div_series(acb_poly_t Q, const acb_poly_t A, const acb_poly_t B, slong n, slong prec)¶ Sets Q to the power series quotient A divided by B, truncated to length n.
-
void
_acb_poly_div(acb_ptr Q, acb_srcptr A, slong lenA, acb_srcptr B, slong lenB, slong prec)¶
-
void
_acb_poly_rem(acb_ptr R, acb_srcptr A, slong lenA, acb_srcptr B, slong lenB, slong prec)¶
-
void
_acb_poly_divrem(acb_ptr Q, acb_ptr R, acb_srcptr A, slong lenA, acb_srcptr B, slong lenB, slong prec)¶
-
int
acb_poly_divrem(acb_poly_t Q, acb_poly_t R, const acb_poly_t A, const acb_poly_t B, slong prec)¶ Performs polynomial division with remainder, computing a quotient \(Q\) and a remainder \(R\) such that \(A = BQ + R\). The implementation reverses the inputs and performs power series division.
If the leading coefficient of \(B\) contains zero (or if \(B\) is identically zero), returns 0 indicating failure without modifying the outputs. Otherwise returns nonzero.
Composition¶
-
void
acb_poly_taylor_shift_horner(acb_poly_t g, const acb_poly_t f, const acb_t c, slong prec)¶
-
void
acb_poly_taylor_shift_divconquer(acb_poly_t g, const acb_poly_t f, const acb_t c, slong prec)¶
-
void
acb_poly_taylor_shift_convolution(acb_poly_t g, const acb_poly_t f, const acb_t c, slong prec)¶
-
void
acb_poly_taylor_shift(acb_poly_t g, const acb_poly_t f, const acb_t c, slong prec)¶ Sets g to the Taylor shift \(f(x+c)\), computed respectively using an optimized form of Horner’s rule, divide-and-conquer, a single convolution, and an automatic choice between the three algorithms.
The underscore methods act in-place on g = f which has length n.
-
void
_acb_poly_compose_horner(acb_ptr res, acb_srcptr poly1, slong len1, acb_srcptr poly2, slong len2, slong prec)¶
-
void
acb_poly_compose_horner(acb_poly_t res, const acb_poly_t poly1, const acb_poly_t poly2, slong prec)¶
-
void
_acb_poly_compose_divconquer(acb_ptr res, acb_srcptr poly1, slong len1, acb_srcptr poly2, slong len2, slong prec)¶
-
void
acb_poly_compose_divconquer(acb_poly_t res, const acb_poly_t poly1, const acb_poly_t poly2, slong prec)¶
-
void
_acb_poly_compose(acb_ptr res, acb_srcptr poly1, slong len1, acb_srcptr poly2, slong len2, slong prec)¶
-
void
acb_poly_compose(acb_poly_t res, const acb_poly_t poly1, const acb_poly_t poly2, slong prec)¶ Sets res to the composition \(h(x) = f(g(x))\) where \(f\) is given by poly1 and \(g\) is given by poly2, respectively using Horner’s rule, divide-and-conquer, and an automatic choice between the two algorithms.
The default algorithm also handles special-form input \(g = ax^n + c\) efficiently by performing a Taylor shift followed by a rescaling.
The underscore methods do not support aliasing of the output with either input polynomial.
-
void
_acb_poly_compose_series_horner(acb_ptr res, acb_srcptr poly1, slong len1, acb_srcptr poly2, slong len2, slong n, slong prec)¶
-
void
acb_poly_compose_series_horner(acb_poly_t res, const acb_poly_t poly1, const acb_poly_t poly2, slong n, slong prec)¶
-
void
_acb_poly_compose_series_brent_kung(acb_ptr res, acb_srcptr poly1, slong len1, acb_srcptr poly2, slong len2, slong n, slong prec)¶
-
void
acb_poly_compose_series_brent_kung(acb_poly_t res, const acb_poly_t poly1, const acb_poly_t poly2, slong n, slong prec)¶
-
void
_acb_poly_compose_series(acb_ptr res, acb_srcptr poly1, slong len1, acb_srcptr poly2, slong len2, slong n, slong prec)¶
-
void
acb_poly_compose_series(acb_poly_t res, const acb_poly_t poly1, const acb_poly_t poly2, slong n, slong prec)¶ Sets res to the power series composition \(h(x) = f(g(x))\) truncated to order \(O(x^n)\) where \(f\) is given by poly1 and \(g\) is given by poly2, respectively using Horner’s rule, the Brent-Kung baby step-giant step algorithm, and an automatic choice between the two algorithms.
The default algorithm also handles special-form input \(g = ax^n\) efficiently.
We require that the constant term in \(g(x)\) is exactly zero. The underscore methods do not support aliasing of the output with either input polynomial.
-
void
acb_poly_revert_series_lagrange(acb_poly_t h, const acb_poly_t f, slong n, slong prec)¶
-
void
acb_poly_revert_series_newton(acb_poly_t h, const acb_poly_t f, slong n, slong prec)¶
-
void
_acb_poly_revert_series_lagrange_fast(acb_ptr h, acb_srcptr f, slong flen, slong n, slong prec)¶
-
void
acb_poly_revert_series_lagrange_fast(acb_poly_t h, const acb_poly_t f, slong n, slong prec)¶
-
void
acb_poly_revert_series(acb_poly_t h, const acb_poly_t f, slong n, slong prec)¶ Sets \(h\) to the power series reversion of \(f\), i.e. the expansion of the compositional inverse function \(f^{-1}(x)\), truncated to order \(O(x^n)\), using respectively Lagrange inversion, Newton iteration, fast Lagrange inversion, and a default algorithm choice.
We require that the constant term in \(f\) is exactly zero and that the linear term is nonzero. The underscore methods assume that flen is at least 2, and do not support aliasing.
Evaluation¶
-
void
acb_poly_evaluate_horner(acb_t y, const acb_poly_t f, const acb_t x, slong prec)¶
-
void
acb_poly_evaluate_rectangular(acb_t y, const acb_poly_t f, const acb_t x, slong prec)¶
-
void
acb_poly_evaluate(acb_t y, const acb_poly_t f, const acb_t x, slong prec)¶ Sets \(y = f(x)\), evaluated respectively using Horner’s rule, rectangular splitting, and an automatic algorithm choice.
-
void
_acb_poly_evaluate2_horner(acb_t y, acb_t z, acb_srcptr f, slong len, const acb_t x, slong prec)¶
-
void
_acb_poly_evaluate2_rectangular(acb_t y, acb_t z, acb_srcptr f, slong len, const acb_t x, slong prec)¶
-
void
acb_poly_evaluate2_rectangular(acb_t y, acb_t z, const acb_poly_t f, const acb_t x, slong prec)¶
-
void
acb_poly_evaluate2(acb_t y, acb_t z, const acb_poly_t f, const acb_t x, slong prec)¶ Sets \(y = f(x), z = f'(x)\), evaluated respectively using Horner’s rule, rectangular splitting, and an automatic algorithm choice.
When Horner’s rule is used, the only advantage of evaluating the function and its derivative simultaneously is that one does not have to generate the derivative polynomial explicitly. With the rectangular splitting algorithm, the powers can be reused, making simultaneous evaluation slightly faster.
Product trees¶
-
void
_acb_poly_product_roots(acb_ptr poly, acb_srcptr xs, slong n, slong prec)¶
-
void
acb_poly_product_roots(acb_poly_t poly, acb_srcptr xs, slong n, slong prec)¶ Generates the polynomial \((x-x_0)(x-x_1)\cdots(x-x_{n-1})\).
-
acb_ptr *
_acb_poly_tree_alloc(slong len)¶ Returns an initialized data structured capable of representing a remainder tree (product tree) of len roots.
-
void
_acb_poly_tree_free(acb_ptr *tree, slong len)¶ Deallocates a tree structure as allocated using _acb_poly_tree_alloc.
-
void
_acb_poly_tree_build(acb_ptr *tree, acb_srcptr roots, slong len, slong prec)¶ Constructs a product tree from a given array of len roots. The tree structure must be pre-allocated to the specified length using
_acb_poly_tree_alloc().
Multipoint evaluation¶
-
void
_acb_poly_evaluate_vec_iter(acb_ptr ys, acb_srcptr poly, slong plen, acb_srcptr xs, slong n, slong prec)¶
-
void
acb_poly_evaluate_vec_iter(acb_ptr ys, const acb_poly_t poly, acb_srcptr xs, slong n, slong prec)¶ Evaluates the polynomial simultaneously at n given points, calling
_acb_poly_evaluate()repeatedly.
-
void
_acb_poly_evaluate_vec_fast_precomp(acb_ptr vs, acb_srcptr poly, slong plen, acb_ptr *tree, slong len, slong prec)¶
-
void
_acb_poly_evaluate_vec_fast(acb_ptr ys, acb_srcptr poly, slong plen, acb_srcptr xs, slong n, slong prec)¶
-
void
acb_poly_evaluate_vec_fast(acb_ptr ys, const acb_poly_t poly, acb_srcptr xs, slong n, slong prec)¶ Evaluates the polynomial simultaneously at n given points, using fast multipoint evaluation.
Interpolation¶
-
void
_acb_poly_interpolate_newton(acb_ptr poly, acb_srcptr xs, acb_srcptr ys, slong n, slong prec)¶
-
void
acb_poly_interpolate_newton(acb_poly_t poly, acb_srcptr xs, acb_srcptr ys, slong n, slong prec)¶ Recovers the unique polynomial of length at most n that interpolates the given x and y values. This implementation first interpolates in the Newton basis and then converts back to the monomial basis.
-
void
_acb_poly_interpolate_barycentric(acb_ptr poly, acb_srcptr xs, acb_srcptr ys, slong n, slong prec)¶
-
void
acb_poly_interpolate_barycentric(acb_poly_t poly, acb_srcptr xs, acb_srcptr ys, slong n, slong prec)¶ Recovers the unique polynomial of length at most n that interpolates the given x and y values. This implementation uses the barycentric form of Lagrange interpolation.
-
void
_acb_poly_interpolate_fast_precomp(acb_ptr poly, acb_srcptr ys, acb_ptr *tree, acb_srcptr weights, slong len, slong prec)¶
-
void
_acb_poly_interpolate_fast(acb_ptr poly, acb_srcptr xs, acb_srcptr ys, slong len, slong prec)¶
-
void
acb_poly_interpolate_fast(acb_poly_t poly, acb_srcptr xs, acb_srcptr ys, slong n, slong prec)¶ Recovers the unique polynomial of length at most n that interpolates the given x and y values, using fast Lagrange interpolation. The precomp function takes a precomputed product tree over the x values and a vector of interpolation weights as additional inputs.
Differentiation¶
-
void
_acb_poly_derivative(acb_ptr res, acb_srcptr poly, slong len, slong prec)¶ Sets {res, len - 1} to the derivative of {poly, len}. Allows aliasing of the input and output.
-
void
acb_poly_derivative(acb_poly_t res, const acb_poly_t poly, slong prec)¶ Sets res to the derivative of poly.
-
void
_acb_poly_integral(acb_ptr res, acb_srcptr poly, slong len, slong prec)¶ Sets {res, len} to the integral of {poly, len - 1}. Allows aliasing of the input and output.
-
void
acb_poly_integral(acb_poly_t res, const acb_poly_t poly, slong prec)¶ Sets res to the integral of poly.
Transforms¶
-
void
_acb_poly_borel_transform(acb_ptr res, acb_srcptr poly, slong len, slong prec)¶
-
void
acb_poly_borel_transform(acb_poly_t res, const acb_poly_t poly, slong prec)¶ Computes the Borel transform of the input polynomial, mapping \(\sum_k a_k x^k\) to \(\sum_k (a_k / k!) x^k\). The underscore method allows aliasing.
-
void
_acb_poly_inv_borel_transform(acb_ptr res, acb_srcptr poly, slong len, slong prec)¶
-
void
acb_poly_inv_borel_transform(acb_poly_t res, const acb_poly_t poly, slong prec)¶ Computes the inverse Borel transform of the input polynomial, mapping \(\sum_k a_k x^k\) to \(\sum_k a_k k! x^k\). The underscore method allows aliasing.
-
void
_acb_poly_binomial_transform_basecase(acb_ptr b, acb_srcptr a, slong alen, slong len, slong prec)¶
-
void
acb_poly_binomial_transform_basecase(acb_poly_t b, const acb_poly_t a, slong len, slong prec)¶
-
void
_acb_poly_binomial_transform_convolution(acb_ptr b, acb_srcptr a, slong alen, slong len, slong prec)¶
-
void
acb_poly_binomial_transform_convolution(acb_poly_t b, const acb_poly_t a, slong len, slong prec)¶
-
void
acb_poly_binomial_transform(acb_poly_t b, const acb_poly_t a, slong len, slong prec)¶ Computes the binomial transform of the input polynomial, truncating the output to length len. See
arb_poly_binomial_transform()for details.The underscore methods do not support aliasing, and assume that the lengths are nonzero.
Elementary functions¶
-
void
_acb_poly_pow_ui_trunc_binexp(acb_ptr res, acb_srcptr f, slong flen, ulong exp, slong len, slong prec)¶ Sets {res, len} to {f, flen} raised to the power exp, truncated to length len. Requires that len is no longer than the length of the power as computed without truncation (i.e. no zero-padding is performed). Does not support aliasing of the input and output, and requires that flen and len are positive. Uses binary expontiation.
-
void
acb_poly_pow_ui_trunc_binexp(acb_poly_t res, const acb_poly_t poly, ulong exp, slong len, slong prec)¶ Sets res to poly raised to the power exp, truncated to length len. Uses binary exponentiation.
-
void
_acb_poly_pow_ui(acb_ptr res, acb_srcptr f, slong flen, ulong exp, slong prec)¶ Sets res to {f, flen} raised to the power exp. Does not support aliasing of the input and output, and requires that flen is positive.
-
void
acb_poly_pow_ui(acb_poly_t res, const acb_poly_t poly, ulong exp, slong prec)¶ Sets res to poly raised to the power exp.
-
void
_acb_poly_pow_series(acb_ptr h, acb_srcptr f, slong flen, acb_srcptr g, slong glen, slong len, slong prec)¶ Sets {h, len} to the power series \(f(x)^{g(x)} = \exp(g(x) \log f(x))\) truncated to length len. This function detects special cases such as g being an exact small integer or \(\pm 1/2\), and computes such powers more efficiently. This function does not support aliasing of the output with either of the input operands. It requires that all lengths are positive, and assumes that flen and glen do not exceed len.
-
void
acb_poly_pow_series(acb_poly_t h, const acb_poly_t f, const acb_poly_t g, slong len, slong prec)¶ Sets h to the power series \(f(x)^{g(x)} = \exp(g(x) \log f(x))\) truncated to length len. This function detects special cases such as g being an exact small integer or \(\pm 1/2\), and computes such powers more efficiently.
-
void
_acb_poly_pow_acb_series(acb_ptr h, acb_srcptr f, slong flen, const acb_t g, slong len, slong prec)¶ Sets {h, len} to the power series \(f(x)^g = \exp(g \log f(x))\) truncated to length len. This function detects special cases such as g being an exact small integer or \(\pm 1/2\), and computes such powers more efficiently. This function does not support aliasing of the output with either of the input operands. It requires that all lengths are positive, and assumes that flen does not exceed len.
-
void
acb_poly_pow_acb_series(acb_poly_t h, const acb_poly_t f, const acb_t g, slong len, slong prec)¶ Sets h to the power series \(f(x)^g = \exp(g \log f(x))\) truncated to length len.
-
void
acb_poly_sqrt_series(acb_poly_t g, const acb_poly_t h, slong n, slong prec)¶ Sets g to the power series square root of h, truncated to length n. Uses division-free Newton iteration for the reciprocal square root, followed by a multiplication.
The underscore method does not support aliasing of the input and output arrays. It requires that hlen and n are greater than zero.
-
void
acb_poly_rsqrt_series(acb_poly_t g, const acb_poly_t h, slong n, slong prec)¶ Sets g to the reciprocal power series square root of h, truncated to length n. Uses division-free Newton iteration.
The underscore method does not support aliasing of the input and output arrays. It requires that hlen and n are greater than zero.
-
void
acb_poly_log_series(acb_poly_t res, const acb_poly_t f, slong n, slong prec)¶ Sets res to the power series logarithm of f, truncated to length n. Uses the formula \(\log(f(x)) = \int f'(x) / f(x) dx\), adding the logarithm of the constant term in f as the constant of integration.
The underscore method supports aliasing of the input and output arrays. It requires that flen and n are greater than zero.
-
void
acb_poly_log1p_series(acb_poly_t res, const acb_poly_t f, slong n, slong prec)¶ Computes the power series \(\log(1+f)\), with better accuracy when the constant term of f is small.
-
void
acb_poly_atan_series(acb_poly_t res, const acb_poly_t f, slong n, slong prec)¶ Sets res the power series inverse tangent of f, truncated to length n.
Uses the formula
\[\tan^{-1}(f(x)) = \int f'(x) / (1+f(x)^2) dx,\]adding the function of the constant term in f as the constant of integration.
The underscore method supports aliasing of the input and output arrays. It requires that flen and n are greater than zero.
-
void
acb_poly_exp_series_basecase(acb_poly_t f, const acb_poly_t h, slong n, slong prec)¶
-
void
acb_poly_exp_series(acb_poly_t f, const acb_poly_t h, slong n, slong prec)¶ Sets \(f\) to the power series exponential of \(h\), truncated to length \(n\).
The basecase version uses a simple recurrence for the coefficients, requiring \(O(nm)\) operations where \(m\) is the length of \(h\).
The main implementation uses Newton iteration, starting from a small number of terms given by the basecase algorithm. The complexity is \(O(M(n))\). Redundant operations in the Newton iteration are avoided by using the scheme described in [HZ2004].
The underscore methods support aliasing and allow the input to be shorter than the output, but require the lengths to be nonzero.
-
void
acb_poly_exp_pi_i_series(acb_poly_t f, const acb_poly_t h, slong n, slong prec)¶ Sets f to the power series \(\exp(\pi i h)\) truncated to length n. The underscore method supports aliasing and allows the input to be shorter than the output, but requires the lengths to be nonzero.
-
void
_acb_poly_sin_cos_series_basecase(acb_ptr s, acb_ptr c, acb_srcptr h, slong hlen, slong n, slong prec, int times_pi)¶
-
void
acb_poly_sin_cos_series_basecase(acb_poly_t s, acb_poly_t c, const acb_poly_t h, slong n, slong prec, int times_pi)¶
-
void
_acb_poly_sin_cos_series_tangent(acb_ptr s, acb_ptr c, acb_srcptr h, slong hlen, slong n, slong prec, int times_pi)¶
-
void
acb_poly_sin_cos_series_tangent(acb_poly_t s, acb_poly_t c, const acb_poly_t h, slong n, slong prec, int times_pi)¶