FastPolyEval  1.0
Fast Evaluation of Real and Complex Polynomials
eval.h
Go to the documentation of this file.
1 //
2 // eval.h
3 //
4 // Authors: Nicolae Mihalache & François Vigneron
5 //
6 // This software is released under BSD licence, with an attribution clause (see Licence file).
7 // Please cite the reference below if you use or distribute this software.
8 //
9 // • [1] R. Anton, N. Mihalache & F. Vigneron. Fast evaluation of real and complex polynomials. 2022.
10 // https://hal.archives-ouvertes.fr/hal-03820369
11 //
12 // Copyright 2022 Univ. Paris-Est Créteil & Univ. de Reims Champagne-Ardenne.
13 //
14 
20 #ifndef eval_h
21 #define eval_h
22 
23 #include <mpfr.h>
24 
25 #include "ntypes.h"
26 #include "comp.h"
27 #include "concave.h"
28 #include "poly.h"
29 #include "polyr.h"
30 #include "pows.h"
31 #include "powsr.h"
32 
33 // MARK: constants &data types
34 
36 #define MIN_EVAL_PREC 8
37 #define MIN_EVAL_PREC_STR "8"
38 
41 typedef struct {
42  poly P;
43  polyr Q;
44  bool real;
47  pows zn;
54  mpfr_t br;
55 } eval_struct;
56 
60 typedef eval_struct eval_t[1];
61 
63 typedef eval_struct *eval;
64 
65 // MARK: macros & functions prototypes
66 
68 #define eval_lastValError(ev) (ev->valErr)
69 
71 #define eval_lastDerError(ev) (ev->derErr)
72 
80 
88 
95 bool eval_free(eval ev);
96 
107 bool eval_val(comp v, eval ev, comp z);
108 
119 bool eval_der(comp d, eval ev, comp z);
120 
133 bool eval_val_der(comp v, comp d, eval ev, comp z);
134 
146 bool eval_newton(comp nt, eval ev, comp z);
147 
148 // MARK: Bare-metal variants
149 
160 
171 
183 bool eval_val_cc(comp v, eval ev, comp z);
184 
196 bool eval_der_cc(comp d, eval ev, comp z);
197 
211 bool eval_val_der_cc(comp v, comp d, eval ev, comp z);
212 
225 bool eval_newton_cc(comp nt, eval ev, comp z);
226 
238 bool eval_val_cr(comp v, eval ev, mpfr_t x);
239 
251 bool eval_der_cr(comp d, eval ev, mpfr_t x);
252 
266 bool eval_val_der_cr(comp v, comp d, eval ev, mpfr_t x);
267 
280 bool eval_newton_cr(comp nt, eval ev, mpfr_t x);
281 
293 bool eval_val_rc(comp v, eval ev, comp z);
294 
306 bool eval_der_rc(comp d, eval ev, comp z);
307 
321 bool eval_val_der_rc(comp v, comp d, eval ev, comp z);
322 
335 bool eval_newton_rc(comp nt, eval ev, comp z);
336 
348 bool eval_val_rr(mpfr_t v, eval ev, mpfr_t x);
349 
361 bool eval_der_rr(mpfr_t d, eval ev, mpfr_t x);
362 
376 bool eval_val_der_rr(mpfr_t v, mpfr_t d, eval ev, mpfr_t x);
377 
390 bool eval_newton_rr(mpfr_t nt, eval ev, mpfr_t x);
391 
392 #endif /* eval_h */
Definition of MPFR complex numbers.
comp_struct comp[1]
Practical wrapper for comp_struct.
Definition: comp.h:39
Definition of a concave function computed from the coefficients of some polynomial.
bool eval_val_der_rr(mpfr_t v, mpfr_t d, eval ev, mpfr_t x)
Evaluates ev->Q(x) and ev->Q'(x) using the method described in [1].
bool eval_val_cc(comp v, eval ev, comp z)
Evaluates ev->P(x) using the method described in [1].
bool eval_analyse(eval ev)
Analyses the complex polynomial ev->P, after some of its coefficients have been changed.
bool eval_val_der_cr(comp v, comp d, eval ev, mpfr_t x)
Evaluates ev->P(x) and ev->P'(x) using the method described in [1].
bool eval_analyse_r(eval ev)
Analyses the real polynomial ev->Q, after some of its coefficients have been changed.
bool eval_val_rr(mpfr_t v, eval ev, mpfr_t x)
Evaluates the real polynomial ev->Q(x) using the method described in [1].
eval eval_new_r(polyr Q, prec_t prec)
Returns a new evaluator of the real polynomial Q.
bool eval_der(comp d, eval ev, comp z)
Evaluates ev->P'(z) (or ev->Q'(z)) using the method described in [1].
bool eval_val_der(comp v, comp d, eval ev, comp z)
Evaluates ev->P(z) and ev->P'(z) (or ev->Q(z) and ev->Q'(z)) using the method described in [1].
eval_struct eval_t[1]
Practical wrapper for eval_struct.
Definition: eval.h:60
bool eval_val_rc(comp v, eval ev, comp z)
Evaluates ev->Q(x) using the method described in [1].
bool eval_der_cr(comp d, eval ev, mpfr_t x)
Evaluates ev->P'(x) using the method described in [1].
bool eval_newton_cc(comp nt, eval ev, comp z)
Computes the Newthon method step of ev->P using the method described in [1].
bool eval_val_der_rc(comp v, comp d, eval ev, comp z)
Evaluates ev->Q(x) and ev->Q'(x) using the method described in [1].
eval eval_new(poly P, prec_t prec)
Returns a new evaluator of the complex polynomial P.
bool eval_val_der_cc(comp v, comp d, eval ev, comp z)
Evaluates ev->P(x) and ev->P'(x) using the method described in [1].
bool eval_newton_rr(mpfr_t nt, eval ev, mpfr_t x)
Computes the Newthon method step of the real polynomial ev->Q using the method described in [1].
eval_struct * eval
Convenience pointer to eval_struct.
Definition: eval.h:63
bool eval_newton_rc(comp nt, eval ev, comp z)
Computes the Newthon method step of ev->Q using the method described in [1].
bool eval_newton_cr(comp nt, eval ev, mpfr_t x)
Computes the Newthon method step of ev->P using the method described in [1].
bool eval_der_cc(comp d, eval ev, comp z)
Evaluates ev->P'(x) using the method described in [1].
bool eval_der_rc(comp d, eval ev, comp z)
Evaluates ev->Q'(x) using the method described in [1].
bool eval_val_cr(comp v, eval ev, mpfr_t x)
Evaluates ev->P(x) using the method described in [1].
bool eval_der_rr(mpfr_t d, eval ev, mpfr_t x)
Evaluates the derivative of real polynomial ev->Q'(x) using the method described in [1].
bool eval_free(eval ev)
Frees all the memory used by the evaluator ev, assuming the struct has been allocated with malloc(),...
bool eval_newton(comp nt, eval ev, comp z)
Computes the Newthon method step of ev->P (or ev->Q) using the method described in [1].
bool eval_val(comp v, eval ev, comp z)
Evaluates ev->P(z) (or ev->Q(z)) using the method described in [1].
Definition of basic types.
double real_t
The machine number type to use for polynomial analysis and preconditionning.
Definition: ntypes.h:163
ulong deg_t
The integer number type to use for polynomial degrees and indexes.
Definition: ntypes.h:128
ulong prec_t
The integer number type to use for polynomial degrees and indexes.
Definition: ntypes.h:188
Definition of complex polynomials with arbitary precision coefficients.
Definition of real polynomials with arbitary precision coefficients.
Definition of a buffer for pre-computed powers of a complex number with arbitary precision.
Definition of a buffer for pre-computed powers of a real number with arbitary precision.
Description of a concave function computed from the coefficients of some polynomial.
Definition: concave.h:35
Evaluator of polynomials with multi-precision floating point coefficients.
Definition: eval.h:41
concave f
concave cover of the scales of coefficients
Definition: eval.h:46
real_t derErr
the [approximative] upper bound for the absolute error of the last derivative evaluation,...
Definition: eval.h:50
bool real
the type of polynomial to evaluate
Definition: eval.h:44
prec_t prec
the precision of evaluations, in bits
Definition: eval.h:45
pows zn
powers of a complex argument
Definition: eval.h:47
poly P
the polynomial, with complex coefficients
Definition: eval.h:42
real_t ntErr
the [approximative] upper bound for the absolute error of the last Newton term evaluation,...
Definition: eval.h:51
deg_t terms
the number of polynomial terms computed by the last operation
Definition: eval.h:52
powsr xn
powers of a real argument
Definition: eval.h:48
polyr Q
the polynomial, with real coefficients
Definition: eval.h:43
real_t valErr
the [approximative] upper bound for the absolute error of the last evaluation, in bits
Definition: eval.h:49
mpfr_t br
a real buffer for internal computations
Definition: eval.h:54
comp buf
a buffer for internal computations
Definition: eval.h:53
Polynomial with multi-precision floating point complex coefficients.
Definition: poly.h:31
Polynomial with multi-precision floating point complex coefficients.
Definition: polyr.h:31
The powers of the complex number z using multi-precision floating point numbers.
Definition: pows.h:30
The powers of the real number x using multi-precision floating point numbers.
Definition: powsr.h:37