Main Page | Modules | Alphabetical List | Class List | Directories | Class Members

Polynomials


Defines

#define CPL_POLYNOMIAL_CMP
 Compare the coefficients of two polynomials.

Functions

cpl_polynomial * cpl_polynomial_new (int dim)
 Create a new cpl_polynomial.
void cpl_polynomial_delete (cpl_polynomial *p)
 Delete a cpl_polynomial.
cpl_error_code cpl_polynomial_dump (const cpl_polynomial *p, FILE *stream)
 Dump a cpl_polynomial as ASCII to a stream.
cpl_polynomial * cpl_polynomial_duplicate (const cpl_polynomial *p)
 This function duplicates an existing polynomial.
cpl_error_code cpl_polynomial_copy (cpl_polynomial *out, const cpl_polynomial *in)
 This function copies contents of a polynomial into another one.
int cpl_polynomial_get_degree (const cpl_polynomial *p)
 The degree of the polynomial.
int cpl_polynomial_get_dimension (const cpl_polynomial *p)
 The dimension of the polynomial.
double cpl_polynomial_get_coeff (const cpl_polynomial *in, const int *pow)
 Get a coefficient of the polynomial.
cpl_error_code cpl_polynomial_set_coeff (cpl_polynomial *in, const int *pow, double c)
 Set a coefficient of the polynomial.
double cpl_polynomial_eval (const cpl_polynomial *p, const cpl_vector *x)
 Evaluate the polynomial at the given point.
double cpl_polynomial_eval_1d (const cpl_polynomial *p, double x, double *pd)
 Evaluate a univariate polynomial using Horners rule.
cpl_error_code cpl_vector_fill_polynomial (cpl_vector *v, const cpl_polynomial *p, double x0, double d)
 Evaluate a 1D-polynomial on equidistant points using Horners rule.
cpl_error_code cpl_polynomial_solve_1d (const cpl_polynomial *p, double x0, double *px, int mul)
 A real solution to p(x) = 0 using Newton-Raphsons method.
cpl_error_code cpl_polynomial_shift_1d (cpl_polynomial *p, double u)
 Given p and u, modify the polynomial to p(x) := p(x+u).
cpl_polynomial * cpl_polynomial_fit_1d_create (const cpl_vector *x_pos, const cpl_vector *values, int degree, double *mse)
 Fit a 1D-polynomial to a 1D-signal in a least squares sense.
cpl_polynomial * cpl_polynomial_fit_2d_create (cpl_bivector *xy_pos, cpl_vector *values, int degree, double *mse)
 Fit a 2D-polynomial to a 2D-surface in a least squares sense.

Detailed Description

This module provides functions to handle uni- and multivariate polynomials.

Univariate polynomials use the Horner rule for evaluation, while multivariate polynomials are evaluated simply as the sum of each term.

This means that of the two polynomials

 * P1(x) = p0 + p1.x + p4.x^2
 * 
and
 * P2(x,y) = p0 + p1.x + p2.y + p3.x.y + p4.x^2 + p5.y^2
 * 
P1(x) may evaluate to more accurate results than P2(x,0), especially around the roots.

Note that a polynomial like P3(z) = p0 + p1.z + p2.z^2 + p3.z^3, z=x^4 is preferable to p4(x) = p0 + p1.x^4 + p2.x^8 + p3.x^12.


Define Documentation

#define CPL_POLYNOMIAL_CMP
 

Value:

/* Verify that it differs within tolerance */           \
                if (fabs(p1->c[i] - p2->c[j]) <= tol) {                 \
                /* Verify that the powers match */                      \
                    for (dim=0; dim < p1->dim; dim++)                   \
                        if (p1->pow[dim][i] != p2->pow[dim][j]) break;  \
                    if (dim == p1->dim) break;   /* - found it */       \
                }
Compare the coefficients of two polynomials.

Parameters:
p1 the 1st polynomial
p2 the 2nd polynomial
tol The absolute (non-negative) tolerance
Returns:
0 when equal, positive when different, negative on error.
The two polynomials are considered equal iff they have identical dimensions and the absolute difference between their coefficients does not exceed the tolerance.

This means that the following pair of polynomials per definition are considered different: P1(x1,x2) = 3*x1 different from P2(x1) = 3*x1.

If all parameters are valid and p1 and p2 point to the same polynomial the functions returns 0.

Possible _cpl_error_code_ set in this function:

  • CPL_ERROR_NULL_INPUT if an input pointer is NULL
  • CPL_ERROR_ILLEGAL_INPUT if tol is negative


Function Documentation

cpl_error_code cpl_polynomial_copy cpl_polynomial *  out,
const cpl_polynomial *  in
 

This function copies contents of a polynomial into another one.

Parameters:
out Pre-allocated output cpl_polynomial
in Input cpl_polynomial
Returns:
CPL_ERROR_NONE or the relevant _cpl_error_code_
in and out must point to different polynomials.

If out already contains coefficients, they are overwritten.

This is the only function that can modify the dimension of a polynomial.

Possible _cpl_error_code_ set in this function:

  • CPL_ERROR_NULL_INPUT if an input pointer is NULL
  • CPL_ERROR_INCOMPATIBLE_INPUT if in and out point to the same polynomial

void cpl_polynomial_delete cpl_polynomial *  p  ) 
 

Delete a cpl_polynomial.

Parameters:
p cpl_polynomial to delete
Returns:
void

cpl_error_code cpl_polynomial_dump const cpl_polynomial *  p,
FILE *  stream
 

Dump a cpl_polynomial as ASCII to a stream.

Parameters:
p Input cpl_polynomial to dump
stream Output stream, accepts stdout or stderr
Returns:
CPL_ERROR_NONE or the relevant _cpl_error_code_
Each coefficient is preceded by its integer power(s) and written on a single line. If the polynomial has non-zero coefficients, only those are printed, otherwise the (zero-valued) constant term is printed.

Comment lines start with the hash character.

Possible _cpl_error_code_ set in this function:

  • CPL_ERROR_NULL_INPUT if an input pointer is NULL
  • CPL_ERROR_FILE_IO if the write operation fails

cpl_polynomial* cpl_polynomial_duplicate const cpl_polynomial *  p  ) 
 

This function duplicates an existing polynomial.

Parameters:
p The input cpl_polynomial
Returns:
A newly allocated cpl_polynomial or NULL on error
Notice that the returned object is a newly allocated cpl_polynomial that must be deallocated by the caller using cpl_polynomial_delete().

Possible _cpl_error_code_ set in this function:

  • CPL_ERROR_NULL_INPUT if an input pointer is NULL

double cpl_polynomial_eval const cpl_polynomial *  p,
const cpl_vector *  x
 

Evaluate the polynomial at the given point.

Parameters:
p The polynomial
x Point of evaluation
Returns:
The computed value or undefined on error.
The length of x must be the polynomial dimension.

A polynomial with no non-zero coefficents evaluates as 0.

For 1-dimensional polynomials the call requires 2n FLOPs where n+1 is the number of coefficients in p, see also cpl_polynomial_eval_1d(). For multivariate polynomials the call requires n*(1+dim) + d_1 + d_2 + ... + d_dim FLOPs, where dim is the dimenstion, n is the number of coefficients in p and d_i is the highest power used in dimension i, i = 1, 2, ..., dim.

Possible _cpl_error_code_ set in this function:

  • CPL_ERROR_NULL_INPUT if an input pointer is NULL
  • CPL_ERROR_INCOMPATIBLE_INPUT if the length of x differs from the dimension of the polynomial

double cpl_polynomial_eval_1d const cpl_polynomial *  p,
double  x,
double *  pd
 

Evaluate a univariate polynomial using Horners rule.

Parameters:
p The 1D-polynomial
x The point of evaluation
pd Iff pd is not null, the derivative p`(x) evaluated at x
Returns:
The result or undefined on error.
A polynomial with no non-zero coefficents evaluates as 0 with a derivative that does likewise.

The result is computed as p_0 + x * ( p_1 + x * ( p_2 + ... x * p_n )) and requires 2n FLOPs where n+1 is the number of coefficients in p.

If the derivative is requested it is computed along with p(x), using a nested Horner rule. This requires 4n FLOPs in total.

Possible _cpl_error_code_ set in this function:

  • CPL_ERROR_NULL_INPUT if an input pointer is NULL
  • CPL_ERROR_INVALID_TYPE if the polynomial has the wrong dimension

cpl_polynomial* cpl_polynomial_fit_1d_create const cpl_vector *  x_pos,
const cpl_vector *  values,
int  degree,
double *  mse
 

Fit a 1D-polynomial to a 1D-signal in a least squares sense.

Parameters:
x_pos Vector of positions of the signal to fit.
values Vector of values of the signal to fit.
degree Non-negative polynomial degree.
mse Iff mse is not null, the mean squared error on success
Returns:
The fitted polynomial or NULL on error
This function fits a 1D-polynomial to a 1D-signal in a least squares sense.

The input signal is given in x_pos and values which must be of the same size.

The size of x_pos must be greater than degree.

mse may be NULL. If it is not, *mse is set to the mean squared error on success, while it is unchanged on error.

The fit is done in the following steps: 1) x_pos is first transformed into xhat = x_pos - mean(x_pos). 2) The Vandermonde matrix is formed from xhat. 3) The normal equations of the Vandermonde matrix is solved. 4) The resulting polynomial in xhat is transformed to back x_pos.

Warning: An increase in the polynomial degree will normally reduce the mean squared error. However, due to rounding errors and the limited accuracy of the solver of the normal equations, an increase in the polynomial degree may at some point cause the mse to _increase_. In some cases this happens with an increase of the polynomial degree from 9 to 10.

Possible _cpl_error_code_ set in this function:

  • CPL_ERROR_NULL_INPUT if an input pointer is NULL
  • CPL_ERROR_ILLEGAL_INPUT if the size of x_pos is lower than degree, or if degree is negative.
  • CPL_ERROR_INCOMPATIBLE_INPUT if x_pos and values have different sizes
  • CPL_ERROR_DIVISION_BY_ZERO if all values in x_pos are identical.
  • CPL_ERROR_SINGULAR_MATRIX if the fit does not succeed

cpl_polynomial* cpl_polynomial_fit_2d_create cpl_bivector *  xy_pos,
cpl_vector *  values,
int  degree,
double *  mse
 

Fit a 2D-polynomial to a 2D-surface in a least squares sense.

Parameters:
xy_pos Bivector positions of the surface to fit.
values Vector of values of the surface to fit.
degree Non-negative polynomial degree.
mse Iff mse is not null, the mean squared error on success
Returns:
The fitted polynomial or NULL on error
See also:
cpl_polynomial_fit_1d_create

cpl_polynomial_get_degree

This function fits a 2D-polynomial to a 2D-surface in a least squares sense.

The input signal is given in xy_pos and values which must be of the same size.

The size of xy_pos must be at least (degree+1)*(degree+2)/2, which is the number of polynomial coefficients to be determined.

mse may be NULL. If it is not, *mse is set to the mean squared error on success, while it is unchanged on error.

Example: For degree=3, the following terms will be computed:

1 x x^2 x^3 y x.y x^2.y y^2 x.y^2 y^3

The fit is done in the following steps: 1) The x-positions are first transformed into xhat = x - mean(x), and the y-positions are first transformed into yhat = y - mean(y). 2) The Vandermonde matrix is formed from xhat and yhat. 3) The normal equations of the Vandermonde matrix is solved. 4) The resulting polynomial in (xhat, yhat) is transformed to back (x,y).

Warning: An increase in the polynomial degree will normally reduce the mean squared error. However, due to rounding errors and the limited accuracy of the solver of the normal equations, an increase in the polynomial degree may at some point cause the mse to _increase_. In some cases this happens with an increase of the polynomial degree from 8 to 9.

Possible _cpl_error_code_ set in this function:

  • CPL_ERROR_NULL_INPUT if an input pointer is NULL
  • CPL_ERROR_ILLEGAL_INPUT if the size of xy_pos is lower than the number of coefficients to be determined, or if degree is negative.
  • CPL_ERROR_INCOMPATIBLE_INPUT if x_pos and values have different sizes
  • CPL_ERROR_DIVISION_BY_ZERO if all x-values or all y-values in xy_pos are identical.
  • CPL_ERROR_SINGULAR_MATRIX if the fit does not succeed
  • CPL_ERROR_NULL_INPUT if an input pointer is NULL

double cpl_polynomial_get_coeff const cpl_polynomial *  in,
const int *  pow
 

Get a coefficient of the polynomial.

Parameters:
in the input polynomial
pow the powers of the different variables
Returns:
The coefficient or undefined on error
The array pow must have the size of the polynomial dimension and have non-negative elements.

It is allowed to specify a (set of) power(s) for which no coefficient has previously been set. In this case functions returns zero.

Possible _cpl_error_code_ set in this function:

  • CPL_ERROR_NULL_INPUT if an input pointer is NULL
  • CPL_ERROR_ILLEGAL_INPUT if pow contains negative values

int cpl_polynomial_get_degree const cpl_polynomial *  p  ) 
 

The degree of the polynomial.

Parameters:
p the polynomial
Returns:
The degree or negative on error
The degree is the highest sum of exponents (with a non-zero coefficient).

If there are no non-zero coefficients the degree is zero.

Possible _cpl_error_code_ set in this function:

  • CPL_ERROR_NULL_INPUT if an input pointer is NULL

int cpl_polynomial_get_dimension const cpl_polynomial *  p  ) 
 

The dimension of the polynomial.

Parameters:
p the polynomial
Returns:
The dimension or negative on error
Possible _cpl_error_code_ set in this function:
  • CPL_ERROR_NULL_INPUT if an input pointer is NULL

cpl_polynomial* cpl_polynomial_new int  dim  ) 
 

Create a new cpl_polynomial.

Parameters:
dim The positive polynomial dimension (number of variables)
Returns:
1 newly allocated cpl_polynomial, or NULL on error
The returned object must be deallocated using cpl_polynomial_delete().

A newly created polynomial has degree 0 and evaluates as 0.

Possible _cpl_error_code_ set in this function:

  • CPL_ERROR_ILLEGAL_INPUT if dim is negative or zero

cpl_error_code cpl_polynomial_set_coeff cpl_polynomial *  in,
const int *  pow,
double  c
 

Set a coefficient of the polynomial.

Parameters:
in the input polynomial
pow the powers of the different variables
c the coefficient
Returns:
CPL_ERROR_NONE or the relevant _cpl_error_code_
The array pow must have the size of the polynomial dimension and have non-negative elements.

If the coefficient is already there, it is overwritten, if not, a new coefficient is added to the polynomial. This may cause the degree of the polynomial to be increased.

Setting the coefficient of x1^4 * x3^2 in the 4 dimensional polynomial p to 12.3 would be performed by:

cpl_polynomial_set_coeff(p, pow, 12.3); where pow is the integer array [4, 0, 2, 0].

For efficiency reasons the coefficients of a 1d-polynomial are best inserted with the leading coefficient first.

Possible _cpl_error_code_ set in this function:

  • CPL_ERROR_NULL_INPUT if an input pointer is NULL
  • CPL_ERROR_ILLEGAL_INPUT if pow contains negative values

cpl_error_code cpl_polynomial_shift_1d cpl_polynomial *  p,
double  u
 

Given p and u, modify the polynomial to p(x) := p(x+u).

Parameters:
p The 1D-polynomial to be modified in place
u The shift
Returns:
CPL_ERROR_NONE or the relevant _cpl_error_code_
The call requires n*(n+1) FLOPs, where n is the degree of p.

The transformation p(x) := (x+1)^n will generate the binomial coefficients, p_i = p_{n-i} = ( n ), i =0, 1, ..., n. ( i )

Transformation with u = -p_{n-1}/n/p_n will (in absence of rounding errors) yield a polynomial with p_{n-1} = 0 and roots that have the sum zero.

Possible _cpl_error_code_ set in this function:

  • CPL_ERROR_NULL_INPUT if an input pointer is NULL
  • CPL_ERROR_INVALID_TYPE if the polynomial has the wrong dimension

cpl_error_code cpl_polynomial_solve_1d const cpl_polynomial *  p,
double  x0,
double *  px,
int  mul
 

A real solution to p(x) = 0 using Newton-Raphsons method.

Parameters:
p The 1D-polynomial
x0 First guess of the solution
px The solution or undefined on error
mul The root multiplicity (or 1 if unknown)
Returns:
CPL_ERROR_NONE or the relevant _cpl_error_code_
Even if a real solution exists, it may not be found if the first guess is too far from the solution. But a solution is guaranteed to be found if all roots of p are real. If the constant coefficient is zero, the solution 0 will returned regardless of the first guess.

No solution is found and *px is undefined when the iterative process stops because: 1) It can not proceed because p`(x) = 0 (CPL_ERROR_DIVISION_BY_ZERO). 2) Only a finite number of iterations are allowed. (CPL_ERROR_CONTINUE). Both cases may be due to lack of a real solution or a bad first guess.

The accuracy and robustness deteriorates with increasing multiplicity of the solution. This is also the case with numerical multiplicity, i.e. when multiple solutions are located close together.

mul is assumed to be the multiplicity of the solution. Knowledge of the root multiplicity often improves the robustnes and accuracy. If there is no knowledge of the root multiplicity mul should be 1. Setting mul to a too high value should be avoided.

Reverse order of the coefficients: Given x such that p(x) = 0 (p having non-zero constant and leading coefficient) then q(1/x) = 0, where q is obtained by reversing the order of the coefficients of p.

Possible _cpl_error_code_ set in this function:

  • CPL_ERROR_NULL_INPUT if an input pointer is NULL
  • CPL_ERROR_INVALID_TYPE if the polynomial has the wrong dimension
  • CPL_ERROR_ILLEGAL_INPUT if the multiplicity is non-positive
  • CPL_ERROR_DIVISION_BY_ZERO if a division by zero occurs
  • CPL_ERROR_CONTINUE if the algorithm does not converge

cpl_error_code cpl_vector_fill_polynomial cpl_vector *  v,
const cpl_polynomial *  p,
double  x0,
double  d
 

Evaluate a 1D-polynomial on equidistant points using Horners rule.

Parameters:
v Preallocated vector to contain the result
p The 1D-polynomial
x0 The first point of evaluation
d The increment between points of evaluation
Returns:
CPL_ERROR_NONE or the relevant _cpl_error_code_
See also:
cpl_vector_fill
The evaluation points are x_i = x0 + i * d, i=0, 1, ..., n-1, where n is the length of the vector.

If d is zero it is preferable to simply use cpl_vector_fill(v, cpl_polynomial_eval_1d(p, x0, NULL)).

The call requires about 2nm FLOPs, where m+1 is the number of coefficients in p.

Possible _cpl_error_code_ set in this function:

  • CPL_ERROR_NULL_INPUT if an input pointer is NULL
  • CPL_ERROR_INVALID_TYPE if the polynomial has the wrong dimension


Generated on Mon Sep 26 14:38:17 2005 for Common Pipeline Library Reference Manual by  doxygen 1.4.1