![]() |
FastPolyEval
1.0
Fast Evaluation of Real and Complex Polynomials
|
Definition of real polynomials with arbitary precision coefficients. More...
Go to the source code of this file.
Data Structures | |
struct | polyr_struct |
Polynomial with multi-precision floating point complex coefficients. More... | |
Typedefs | |
typedef polyr_struct | polyr_t[1] |
Practical wrapper for polyr_struct . More... | |
typedef polyr_struct * | polyr |
Convenience pointer to polyr_struct . | |
Functions | |
polyr | polyr_new (deg_t degree, prec_t prec) |
Returns a new real polynomial of given degree , with coefficients of precision prec . More... | |
bool | polyr_free (polyr P) |
Frees all the memory used by the polynomial P , assuming the struct has been allocated with malloc() , for example with polyr_new() . More... | |
bool | polyr_set (polyr P, mpfr_t coeff, deg_t ind) |
Sets the coefficient of the polynomial P corresponding to the power ind to coeff . More... | |
bool | polyr_seti (polyr P, long coeff, deg_t ind) |
Sets the coefficient of the polynomial P corresponding to the power ind to coeff . More... | |
bool | polyr_eval_c (comp res, polyr P, comp z) |
Evaluates P(z) using Horner's method. More... | |
bool | polyr_eval (mpfr_t res, polyr P, mpfr_t x) |
Evaluates P(x) using Horner's method. More... | |
polyr | polyr_derivative (polyr P) |
Computes the derivative of P . More... | |
polyr | polyr_sum (polyr P, polyr Q) |
Computes P+Q . More... | |
polyr | polyr_diff (polyr P, polyr Q) |
Computes P-Q . More... | |
polyr | polyr_prod (polyr P, polyr Q) |
Computes P*Q . More... | |
polyr | polyr_sqr (polyr P) |
Computes the square of P . More... | |
polyr | poly_hyp (int n, prec_t prec) |
Computes the n-th hyperbolic polynomial, the n-th image of 0 under the iteration of z->z^2+c . It is a polynomial of degree 2^{n-1} in c . More... | |
polyr | poly_cheb (int n, prec_t prec) |
Computes the Chebyshev polynomial of degree n . More... | |
polyr | poly_leg (int n, prec_t prec) |
Computes the Legendre polynomial of degree n . More... | |
polyr | poly_her (int n, prec_t prec) |
Computes the Hermite polynomial of degree n . More... | |
polyr | poly_lag (int n, prec_t prec) |
Computes the Laguerre polynomial of degree n . More... | |
Definition of real polynomials with arbitary precision coefficients.
Definition in file polyr.h.
typedef polyr_struct polyr_t[1] |
Practical wrapper for polyr_struct
.
To avoid the constant use *
and &
the type polyr_t
is a pointer.
Computes the Chebyshev polynomial of degree n
.
n | the degree of the polynomial |
prec | the precision of the coefficients |
NULL
if some error occurred. Computes the Hermite polynomial of degree n
.
n | the degree of the polynomial |
prec | the precision of the coefficients |
NULL
if some error occurred. Computes the n-th
hyperbolic polynomial, the n-th
image of 0
under the iteration of z->z^2+c
. It is a polynomial of degree 2^{n-1}
in c
.
n | the order of the polynomial |
prec | the precision of the coefficients |
NULL
if some error occurred. Computes the Laguerre polynomial of degree n
.
n | the degree of the polynomial |
prec | the precision of the coefficients |
NULL
if some error occurred. Computes the Legendre polynomial of degree n
.
n | the degree of the polynomial |
prec | the precision of the coefficients |
NULL
if some error occurred. Computes the derivative of P
.
P | the polynomial |
NULL
if some error occurred. Computes P-Q
.
P | a polynomial |
Q | another polynomial |
NULL
if some error occurred. Frees all the memory used by the polynomial P
, assuming the struct has been allocated with malloc()
, for example with polyr_new()
.
P | the polynomial |
Returns a new real polynomial of given degree
, with coefficients of precision prec
.
prec
must be at least precf
and the degree at most MAX_DEG
.degree | the degree of the polynomial |
prec | the precision of the coefficients, in bits |
NULL
if some error occurred. Computes P*Q
.
P | a polynomial |
Q | another polynomial |
NULL
if some error occurred. Computes the square of P
.
P | the polynomial |
NULL
if some error occurred.