ccgsl 2.7.2
C++wrappersforGnuScientificLibrary
gsl::poly Namespace Reference

Various functions for handling polynomials, represented as sequences of coefficients \(a_0,\ldots,a_n\) for \(a_0+\cdots+a_nx^n\). More...

Namespaces

namespace  complex_poly
 Polynomial with complex coefficients.
 

Classes

class  complex_workspace
 Workspace for solving polynomials. More...
 

Functions

template<typename T >
double eval (T const &c, double const x)
 C++ version of gsl_poly_eval(). More...
 
template<typename T >
complex complex_eval (T const &c, gsl::complex const z)
 C++ version of gsl_poly_complex_eval(). More...
 
template<typename T , typename U >
int eval_derivs (T const &c, double const x, U &res)
 C++ version of gsl_poly_eval_derivs(). More...
 
template<typename T , typename U , typename V >
int dd_init (T &dd, U const &x, V const &y)
 C++ version of gsl_poly_dd_init(). More...
 
template<typename DD , typename Z , typename XA , typename YA , typename DYA >
int dd_hermite_init (DD &dd, Z &z, XA const &xa, YA const &ya, DYA const &dya)
 This function computes a divided-difference representation of the interpolating Hermite polynomial for the points (x, y) stored in the arrays xa and ya of length size. More...
 
template<typename T , typename U >
double dd_eval (T const &dd, U const &xa, double const x)
 C++ version of gsl_poly_dd_eval(). More...
 
template<typename C , typename DD , typename X , typename W >
int dd_taylor (C &c, double xp, DD const &dd, X const &x, W &w)
 C++ version of gsl_poly_dd_taylor(). More...
 
int solve_quadratic (double a, double b, double c, double &x0, double &x1)
 C++ version of gsl_poly_solve_quadratic(). More...
 
int complex_solve_quadratic (double a, double b, double c, complex &z0, complex &z1)
 C++ version of gsl_poly_complex_solve_quadratic(). More...
 
int solve_cubic (double a, double b, double c, double &x0, double &x1, double &x2)
 C++ version of gsl_poly_solve_cubic(). More...
 
int complex_solve_cubic (double a, double b, double c, complex &z0, complex &z1, complex &z2)
 C++ version of gsl_poly_complex_solve_cubic(). More...
 
template<typename A , typename Z >
int complex_solve (A const &a, complex_workspace &w, Z &z)
 C++ version of gsl_poly_complex_solve(). More...
 

Detailed Description

Various functions for handling polynomials, represented as sequences of coefficients \(a_0,\ldots,a_n\) for \(a_0+\cdots+a_nx^n\).

The sequences can be stored in arrays (as GSL does) or in any container that has a size() functiion and data() function returning an array of doubles of length container.size(). In particular, there are versions of the GSL functions that contain no size paramater work with a mixture of std::array, std::vector and gsl::vector objects, which have these properties. It is only possible to mix C-style arrays with containers by using the functions that require an array size and passing the remaining arrays in as container.data().

Function Documentation

◆ complex_eval()

template<typename T >
complex gsl::poly::complex_eval ( T const &  c,
gsl::complex const  z 
)
inline

C++ version of gsl_poly_complex_eval().

Works with std::array<double>, std::vector<double> and gsl::vector.

Parameters
cArray of polynomial coefficients.
zValue at which to evaluate polynomial.
Returns
Value of polynomial at z.

Definition at line 286 of file poly.hpp.

References gsl::complex::get().

◆ complex_solve()

template<typename A , typename Z >
int gsl::poly::complex_solve ( A const &  a,
complex_workspace w,
Z &  z 
)
inline

C++ version of gsl_poly_complex_solve().

Here the actual sizes of a and z should be \(n\) and \(2(n-1)\) because the complex values are stored as pairs of doubles in arrays or vectors of doubles. Works with std::array<double>, std::vector<double> and gsl::vector.

Parameters
aArray of coefficients.
wWorkspace (with same size \(n\) as a).
zRoots ( \(n-1\), return value).
Returns
Error code on failure.

Definition at line 581 of file poly.hpp.

References gsl::sf::mathieu::a(), gsl::poly::complex_workspace::get(), and gsl::rstat::n().

◆ complex_solve_cubic()

int gsl::poly::complex_solve_cubic ( double  a,
double  b,
double  c,
complex z0,
complex z1,
complex z2 
)
inline

C++ version of gsl_poly_complex_solve_cubic().

Parameters
aCoefficient of \(x^2\).
bCoefficient of \(x\).
cConstant in cubic.
z0Root (return value).
z1Root (return value).
z2Root (return value).
Returns
Error code on failure.

Definition at line 552 of file poly.hpp.

References gsl::sf::mathieu::a(), gsl::sf::mathieu::b(), and gsl::complex::get().

◆ complex_solve_quadratic()

int gsl::poly::complex_solve_quadratic ( double  a,
double  b,
double  c,
complex z0,
complex z1 
)
inline

C++ version of gsl_poly_complex_solve_quadratic().

Parameters
aCoefficient of \(x^2\).
bCoefficient of \(x\).
cConstant in quadratic.
z0Root (return value).
z1Root (return value).
Returns
Error code on failure.

Definition at line 508 of file poly.hpp.

References gsl::sf::mathieu::a(), gsl::sf::mathieu::b(), and gsl::complex::get().

◆ dd_eval()

template<typename T , typename U >
double gsl::poly::dd_eval ( T const &  dd,
U const &  xa,
double const  x 
)
inline

C++ version of gsl_poly_dd_eval().

Works with std::array<double>, std::vector<double> and gsl::vector.

Parameters
ddAn array of length size representing the divided differences (input).
xaAn array of x values of length size representing the divided differences (input).
xValue at which to evaluate polynomial.
Returns
Value of inrterpolating polynomial at x.

Definition at line 436 of file poly.hpp.

References gsl::cheb::size().

◆ dd_hermite_init()

template<typename DD , typename Z , typename XA , typename YA , typename DYA >
int gsl::poly::dd_hermite_init ( DD &  dd,
Z &  z,
XA const &  xa,
YA const &  ya,
DYA const &  dya 
)
inline

This function computes a divided-difference representation of the interpolating Hermite polynomial for the points (x, y) stored in the arrays xa and ya of length size.

Hermite interpolation constructs polynomials which also match first derivatives dy/dx which are provided in the array dya also of length size. The first derivatives can be incorported into the usual divided-difference algorithm by forming a new dataset \(z = \{x_0,x_0,x_1,x_1,...\}\), which is stored in the array za of length 2*size on output. On output the divided-differences of the Hermite representation are stored in the array dd, also of length 2*size. Using the notation above, \(dd[k] = [z_0,z_1,...,z_k]\). The resulting Hermite polynomial can be evaluated by calling gsl_poly_dd_eval and using za for the input argument xa.

Parameters
ddAn array of length size representing the divided differences (return value).
zThe first derivatives
xaThe x-coördinates of the interpolating Hermite polynomial
yaThe y-coördinates of the interpolating Hermite polynomial
dyaFirst derivatives
Returns
Error code on failure.

Definition at line 393 of file poly.hpp.

References gsl::cheb::size().

◆ dd_init()

template<typename T , typename U , typename V >
int gsl::poly::dd_init ( T &  dd,
U const &  x,
V const &  y 
)
inline

C++ version of gsl_poly_dd_init().

Works with std::array<double>, std::vector<double> and gsl::vector.

Parameters
ddAn array of length size representing the divided differences (return value).
xAn array of x values of length size (input).
yAn array of y values (input).
Returns
Error code on failure.

Definition at line 342 of file poly.hpp.

References gsl::cheb::size().

◆ dd_taylor()

template<typename C , typename DD , typename X , typename W >
int gsl::poly::dd_taylor ( C &  c,
double  xp,
DD const &  dd,
X const &  x,
W &  w 
)
inline

C++ version of gsl_poly_dd_taylor().

Works with std::array<double>, std::vector<double> and gsl::vector.

Parameters
cArray of Taylor coefficients about xp.
xpPoint at which Taylor expansion is computed.
ddAn array of length size of c representing the divided differences.
xAn array of x values of length size of c.
wA workspace in the form of an array of length size of c.
Returns
Error code on failure.

Definition at line 467 of file poly.hpp.

References gsl::cheb::size().

◆ eval()

template<typename T >
double gsl::poly::eval ( T const &  c,
double const  x 
)
inline

C++ version of gsl_poly_eval().

Works with std::array<double>, std::vector<double> and gsl::vector.

Parameters
cArray of polynomial coefficients.
xValue at which to evaluate polynomial.
Returns
Value of polynomial at x.

Definition at line 263 of file poly.hpp.

◆ eval_derivs()

template<typename T , typename U >
int gsl::poly::eval_derivs ( T const &  c,
double const  x,
U &  res 
)
inline

C++ version of gsl_poly_eval_derivs().

Works with std::array<double>, std::vector<double> and gsl::vector.

Parameters
cArray of polynomial coefficients.
xValue at which to evaluate polynomial.
resAn array containing 0th, 1st, 2nd, ... derivatives of polynomial evaluated at x: res.size() is the largest value of \(k\) at which to evaluate \(k\)th derivative of polynomial.
Returns
Error code on failure.

Definition at line 317 of file poly.hpp.

◆ solve_cubic()

int gsl::poly::solve_cubic ( double  a,
double  b,
double  c,
double &  x0,
double &  x1,
double &  x2 
)
inline

C++ version of gsl_poly_solve_cubic().

Cubic with coefficient of \(x^3=1\).

Parameters
aCoefficient of \(x^2\).
bCoefficient of \(x\).
cConstant in cubic.
x0Real root (return value).
x1Real root (return value).
x2Real root (return value).
Returns
The number of real roots.

Definition at line 538 of file poly.hpp.

References gsl::sf::mathieu::a(), and gsl::sf::mathieu::b().

◆ solve_quadratic()

int gsl::poly::solve_quadratic ( double  a,
double  b,
double  c,
double &  x0,
double &  x1 
)
inline

C++ version of gsl_poly_solve_quadratic().

Parameters
aCoefficient of \(x^2\).
bCoefficient of \(x\).
cConstant in quadratic.
x0Real root (return value).
x1Real root (return value).
Returns
The number of real roots.

Definition at line 496 of file poly.hpp.