..function:: void _fmpcb_poly_mullow_classical(fmpcb_struct * res, const fmpcb_struct * poly1, long len1, const fmpcb_struct * poly2, long len2, long n, long prec)
..function:: void fmpcb_poly_mullow_classical(fmpcb_poly_t res, const fmpcb_poly_t poly1, const fmpcb_poly_t poly2, long n, long prec)
..function:: void _fmpcb_poly_mullow_transpose(fmpcb_struct * res, const fmpcb_struct * poly1, long len1, const fmpcb_struct * poly2, long len2, long n, long prec)
..function:: void fmpcb_poly_mullow_transpose(fmpcb_poly_t res, const fmpcb_poly_t poly1, const fmpcb_poly_t poly2, long n, long prec)
..function:: void _fmpcb_poly_mullow(fmpcb_struct * res, const fmpcb_struct * poly1, long len1, const fmpcb_struct * poly2, long len2, long n, long prec)
..function:: void fmpcb_poly_mullow(fmpcb_poly_t res, const fmpcb_poly_t poly1, const fmpcb_poly_t poly2, long n, long prec)
..function:: void _fmpcb_poly_mul(fmpcb_struct * C, const fmpcb_struct * A, long lenA, const fmpcb_struct * B, long lenB, long prec)
Given any complex number `m`, and a nonconstant polynomial `f` and its
derivative `f'`, sets *r* to an interval centered on `m` that is
guaranteed to contain at least one root of `f`.
Such an interval is obtained by taking a ball of radius `|f(m)/f'(m)| n`
where `n` is the degree of `f`.
..function:: long _fmpcb_poly_validate_roots(fmpcb_struct * roots, const fmpcb_struct * poly, long len, long prec)
Given a list of approximate roots of the input polynomial, this
function sets a rigorous bounding interval for each root, and determines
which roots are isolated from all the other roots.
It then rearranges the list of roots so that the isolated roots
are at the front of the list, and returns the count of isolated roots.
If the return value equals the degree of the polynomial, then all
roots have been found. If the return value is smaller, all the
remaining output intervals are guaranteed to contain roots, but
it is possible that not all of the polynomial's roots are contained
among them.
..function:: void _fmpcb_poly_refine_roots_durand_kerner(fmpcb_struct * roots, const fmpcb_struct * poly, long len, long prec)
Refines the given roots simultaneously using a single iteration
of the Durand-Kerner method. The radius of each root is set to an
approximation of the correction, giving a rough estimate of its error (not
a rigorous bound).
..function:: long _fmpcb_poly_find_roots(fmpcb_struct * roots, const fmpcb_struct * poly, const fmpcb_struct * initial, long len, long maxiter, long prec)
..function:: long fmpcb_poly_find_roots(fmpcb_struct * roots, const fmpcb_poly_t poly, const fmpcb_struct * initial, long maxiter, long prec)
Attempts to compute all the roots of the given nonzero polynomial *poly*
using a working precision of *prec* bits. If *n* denotes the degree of *poly*,
the function writes *n* approximate roots with rigorous error bounds to
the preallocated array *roots*, and returns the number of
roots that are isolated.
If the return value equals the degree of the polynomial, then all
roots have been found. If the return value is smaller, all the output
intervals are guaranteed to contain roots, but it is possible that
not all of the polynomial's roots are contained among them.
The roots are computed numerically by performing several steps with
the Durand-Kerner method and terminating if the estimated accuracy of
the roots approaches the working precision or if the number
of steps exceeds *maxiter*, which can be set to zero in order to use
a default value. Finally, the approximate roots are validated rigorously.
Initial values for the iteration can be provided as the array *initial*.
If *initial* is set to *NULL*, default values `(0.4+0.9i)^k` are used.
The polynomial is assumed to be squarefree. If there are repeated
roots, the iteration is likely to find them (with low numerical accuracy),
but the error bounds will not converge as the precision increases.