Common Pipeline Library Reference 7.3.2
|
Functions | |
cpl_error_code | cpl_polynomial_add (cpl_polynomial *self, const cpl_polynomial *first, const cpl_polynomial *second) |
Add two polynomials of the same dimension. | |
int | cpl_polynomial_compare (const cpl_polynomial *self, const cpl_polynomial *other, double tol) |
Compare the coefficients of two polynomials. | |
cpl_error_code | cpl_polynomial_copy (cpl_polynomial *self, const cpl_polynomial *other) |
Copy the contents of one polynomial into another one. | |
void | cpl_polynomial_delete (cpl_polynomial *self) |
Delete a cpl_polynomial. | |
cpl_error_code | cpl_polynomial_derivative (cpl_polynomial *self, cpl_size dim) |
Compute a first order partial derivative. | |
cpl_error_code | cpl_polynomial_dump (const cpl_polynomial *self, FILE *stream) |
Dump a polynomial as ASCII to a stream. | |
cpl_polynomial * | cpl_polynomial_duplicate (const cpl_polynomial *self) |
Duplicate a polynomial. | |
double | cpl_polynomial_eval (const cpl_polynomial *self, const cpl_vector *x) |
Evaluate the polynomial at the given point using Horner's rule. | |
double | cpl_polynomial_eval_1d (const cpl_polynomial *self, double x, double *pd) |
Evaluate a univariate (1D) polynomial using Horner's rule. | |
double | cpl_polynomial_eval_1d_diff (const cpl_polynomial *self, double a, double b, double *ppa) |
Evaluate p(a) - p(b) using Horner's rule. | |
double | cpl_polynomial_eval_2d (const cpl_polynomial *self, double x, double y, double *gradient) |
Evaluate a bivariate (2D) polynomial using Horner's rule and compute the derivatives if needed. | |
double | cpl_polynomial_eval_3d (const cpl_polynomial *self, double x, double y, double z, double *gradient) |
Evaluate a 3D polynomial using Horner's rule and compute the the derivatives if needed. | |
cpl_polynomial * | cpl_polynomial_extract (const cpl_polynomial *self, cpl_size dim, const cpl_polynomial *other) |
Collapse one dimension of a multi-variate polynomial by composition. | |
cpl_error_code | cpl_polynomial_fit (cpl_polynomial *self, const cpl_matrix *samppos, const cpl_boolean *sampsym, const cpl_vector *fitvals, const cpl_vector *fitsigm, cpl_boolean dimdeg, const cpl_size *mindeg, const cpl_size *maxdeg) |
Fit a polynomial to a set of samples in a least squares sense. | |
cpl_polynomial * | cpl_polynomial_fit_1d_create (const cpl_vector *x_pos, const cpl_vector *values, cpl_size degree, double *pmse) |
Fit a 1D-polynomial to a 1D-signal in a least squares sense. | |
cpl_polynomial * | cpl_polynomial_fit_2d_create (const cpl_bivector *xy_pos, const cpl_vector *values, cpl_size degree, double *pmse) |
Fit a 2D-polynomial to a 2D-surface in a least squares sense. | |
double | cpl_polynomial_get_coeff (const cpl_polynomial *self, const cpl_size *pows) |
Get a coefficient of the polynomial. | |
cpl_size | cpl_polynomial_get_degree (const cpl_polynomial *self) |
The degree of the polynomial. | |
cpl_size | cpl_polynomial_get_dimension (const cpl_polynomial *self) |
The dimension of the polynomial. | |
cpl_error_code | cpl_polynomial_multiply (cpl_polynomial *self, const cpl_polynomial *first, const cpl_polynomial *second) |
Multiply two polynomials of the same dimension. | |
cpl_error_code | cpl_polynomial_multiply_scalar (cpl_polynomial *self, const cpl_polynomial *other, double factor) |
Multiply a polynomial with a scalar. | |
cpl_error_code | cpl_polynomial_set_coeff (cpl_polynomial *self, const cpl_size *pows, double value) |
Set a coefficient of the polynomial. | |
cpl_error_code | cpl_polynomial_shift_1d (cpl_polynomial *p, cpl_size i, double u) |
Modify p, p(x0, x1, ..., xi, ...) := (x0, x1, ..., xi+u, ...) | |
cpl_error_code | cpl_polynomial_solve_1d (const cpl_polynomial *p, double x0, double *px, cpl_size mul) |
A real solution to p(x) = 0 using Newton-Raphsons method. | |
cpl_error_code | cpl_polynomial_subtract (cpl_polynomial *self, const cpl_polynomial *first, const cpl_polynomial *second) |
Subtract two polynomials of the same dimension. | |
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 Horner's rule. | |
cpl_error_code | cpl_vector_fill_polynomial_fit_residual (cpl_vector *self, const cpl_vector *fitvals, const cpl_vector *fitsigm, const cpl_polynomial *fit, const cpl_matrix *samppos, double *rechisq) |
Compute the residual of a polynomial fit. | |
This module provides functions to handle uni- and multivariate polynomials.
Comparing the two polynomials
and
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.
Polynomials are evaluated using Horner's method. For multivariate polynomials the evaluation is performed one dimension at a time, starting with the lowest dimension and proceeding upwards through the higher dimensions.
Access to a coefficient of an N-dimensional polynomial has complexity O(N).
cpl_error_code cpl_polynomial_add | ( | cpl_polynomial * | self, |
const cpl_polynomial * | first, | ||
const cpl_polynomial * | second | ||
) |
Add two polynomials of the same dimension.
self | The polynomial to hold the result |
first | The 1st polynomial to add |
second | The 2nd polynomial to add |
Possible CPL error code set in this function:
References cpl_ensure_code, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, CPL_MAX, and cpl_polynomial_get_dimension().
int cpl_polynomial_compare | ( | const cpl_polynomial * | self, |
const cpl_polynomial * | other, | ||
double | tol | ||
) |
Compare the coefficients of two polynomials.
self | The 1st polynomial |
other | The 2nd polynomial |
tol | The absolute (non-negative) tolerance |
The two polynomials are considered equal iff they have identical dimensions and the absolute difference between their coefficients does not exceed the given tolerance.
This means that the following two polynomials per definition are considered different: P1(x1) = 3*x1 different from P2(x1,x2) = 3*x1.
If all input parameters are valid and self and other point to the same polynomial the function returns 0.
If two polynomials have different dimensions, the return value is this (positive) difference.
If two 1D-polynomials differ, the return value is 1 plus the degree of the lowest order differing coefficient.
If for a higher dimension they differ, it is 1 plus the degree of a differing coefficient.
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_ILLEGAL_INPUT, and CPL_ERROR_NULL_INPUT.
cpl_error_code cpl_polynomial_copy | ( | cpl_polynomial * | self, |
const cpl_polynomial * | other | ||
) |
Copy the contents of one polynomial into another one.
self | Pre-allocated output polynomial |
other | Input polynomial |
self and other must point to different polynomials.
If self already contains coefficients, then they are overwritten.
This is the only function that can modify the dimension of a polynomial.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_wlcalib_find_best_1d().
void cpl_polynomial_delete | ( | cpl_polynomial * | self | ) |
Delete a cpl_polynomial.
self | Polynomial to deallocate |
NULL
, nothing is done and no error is set.The function deallocates the memory used by the polynomial self.
References cpl_free().
Referenced by cpl_image_fill_jacobian_polynomial(), cpl_polynomial_fit_1d_create(), cpl_polynomial_fit_2d_create(), cpl_polynomial_multiply(), cpl_ppm_match_points(), and cpl_wlcalib_find_best_1d().
cpl_error_code cpl_polynomial_derivative | ( | cpl_polynomial * | self, |
cpl_size | dim | ||
) |
Compute a first order partial derivative.
self | The polynomial to be modified in place |
dim | The dimension to differentiate (zero for first dimension) |
The dimension of the polynomial is preserved, even if the operation may cause the polynomial to become independent of the dimension dim variable.
The call requires n FLOPs, where n is the number of (non-zero) polynomial coefficients whose power in dimension dim is at least 1.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_ACCESS_OUT_OF_RANGE, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_image_fill_jacobian_polynomial().
cpl_error_code cpl_polynomial_dump | ( | const cpl_polynomial * | self, |
FILE * | stream | ||
) |
Dump a polynomial as ASCII to a stream.
self | The polynomial to process |
stream | Output stream eq. stdout or stderr |
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.
For an N-dimensional polynomial each line thus consists of N power(s) and their coefficient.
Comment lines start with the hash character.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_FILE_IO, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, cpl_free(), cpl_malloc(), and CPL_SIZE_FORMAT.
cpl_polynomial * cpl_polynomial_duplicate | ( | const cpl_polynomial * | self | ) |
Duplicate a polynomial.
self | The non-NULL polynomial to process |
Possible _cpl_error_code_ set in this function:
References cpl_ensure, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_image_fill_jacobian_polynomial().
double cpl_polynomial_eval | ( | const cpl_polynomial * | self, |
const cpl_vector * | x | ||
) |
Evaluate the polynomial at the given point using Horner's rule.
self | The polynomial to access |
x | Point of evaluation |
A polynomial with no non-zero coefficients evaluates as 0.
With n coefficients the complexity is about 2n FLOPs.
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_INCOMPATIBLE_INPUT, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_image_fill_jacobian_polynomial(), cpl_polynomial_fit_2d_create(), and cpl_vector_fill_polynomial_fit_residual().
double cpl_polynomial_eval_1d | ( | const cpl_polynomial * | self, |
double | x, | ||
double * | pd | ||
) |
Evaluate a univariate (1D) polynomial using Horner's rule.
self | The 1D-polynomial |
x | The point of evaluation |
pd | Iff pd is non-NULL, the derivative evaluated at x |
A polynomial with no non-zero coefficents evaluates to 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.
If the derivative is requested it is computed using a nested Horner rule. This requires about 4n FLOPs.
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_INVALID_TYPE, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_vector_fill_polynomial(), cpl_vector_fill_polynomial_fit_residual(), and cpl_wlcalib_find_best_1d().
double cpl_polynomial_eval_1d_diff | ( | const cpl_polynomial * | self, |
double | a, | ||
double | b, | ||
double * | ppa | ||
) |
Evaluate p(a) - p(b) using Horner's rule.
self | The 1D-polynomial |
a | The evaluation point of the minuend |
b | The evaluation point of the subtrahend |
ppa | Iff ppa is not NULL, p(a) |
The call requires about 4n FLOPs where n is the number of coefficients in self, which is the same as that required for two separate polynomial evaluations. cpl_polynomial_eval_1d_diff() is however more accurate.
ppa may be NULL. If it is not, *ppa is set to self(a), which is calculated at no extra cost.
The underlying algorithm is the same as that used in cpl_polynomial_eval_1d() when the derivative is also requested.
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_INVALID_TYPE, and CPL_ERROR_NULL_INPUT.
double cpl_polynomial_eval_2d | ( | const cpl_polynomial * | self, |
double | x, | ||
double | y, | ||
double * | gradient | ||
) |
Evaluate a bivariate (2D) polynomial using Horner's rule and compute the derivatives if needed.
self | The 2D-polynomial |
x | The x component of the evaluation point |
y | The y component of the evaluation point |
gradient | Iff gradient is non-NULL, the 2D gradient array containing the X and Y components of the gradient vector at the evaluated point |
The result is computed as p_0 + x * ( p_1 + x * ( p_2 + ... x * p_n ))
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_INVALID_TYPE, and CPL_ERROR_NULL_INPUT.
double cpl_polynomial_eval_3d | ( | const cpl_polynomial * | self, |
double | x, | ||
double | y, | ||
double | z, | ||
double * | gradient | ||
) |
Evaluate a 3D polynomial using Horner's rule and compute the the derivatives if needed.
self | The 3D-polynomial |
x | The x component of the evaluation point |
y | The y component of the evaluation point |
z | The z component of the evaluation point |
gradient | Iff gradient is non-NULL, the 3D gradient array containing the X, Y and Z components of the gradient vector at the evaluated point |
The result is computed as p_0 + x * ( p_1 + x * ( p_2 + ... x * p_n ))
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_INVALID_TYPE, and CPL_ERROR_NULL_INPUT.
cpl_polynomial * cpl_polynomial_extract | ( | const cpl_polynomial * | self, |
cpl_size | dim, | ||
const cpl_polynomial * | other | ||
) |
Collapse one dimension of a multi-variate polynomial by composition.
self | The multi-variate polynomial |
dim | The dimension to collapse (zero for first dimension) |
other | The polynomial to replace dimension dim of self |
The dimension of the polynomial self must be one greater than that of the other polynomial. Given these two polynomials the dimension dim of self is collapsed by creating a new polynomial from self(x0, x1, ..., x{dim-1}, other(x0, x1, ..., x{dim-1}, x{dim+1}, x{dim+2}, ..., x{n-1}), x{dim+1}, x{dim+2}, ..., x{n-1}).
The created polynomial thus has a dimension which is one less than the polynomial self and which is equal to that of the other polynomial. Collapsing one dimension of a 1D-polynomial is equivalent to evaluating it, which can be done with cpl_polynomial_eval_1d().
FIXME: The other polynomial must currently have a degree of zero, i.e. it must be a constant.
The collapse uses Horner's rule and requires for n coefficients requires about 2n FLOPs.
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:
References cpl_calloc(), cpl_ensure, CPL_ERROR_ACCESS_OUT_OF_RANGE, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_INVALID_TYPE, CPL_ERROR_NULL_INPUT, CPL_ERROR_UNSUPPORTED_MODE, cpl_free(), cpl_polynomial_get_coeff(), and cpl_polynomial_get_degree().
cpl_error_code cpl_polynomial_fit | ( | cpl_polynomial * | self, |
const cpl_matrix * | samppos, | ||
const cpl_boolean * | sampsym, | ||
const cpl_vector * | fitvals, | ||
const cpl_vector * | fitsigm, | ||
cpl_boolean | dimdeg, | ||
const cpl_size * | mindeg, | ||
const cpl_size * | maxdeg | ||
) |
Fit a polynomial to a set of samples in a least squares sense.
self | Polynomial of dimension d to hold the coefficients |
samppos | Matrix of p sample positions, with d rows and p columns |
sampsym | NULL, or d booleans, true iff the sampling is symmetric |
fitvals | Vector of the p values to fit |
fitsigm | Uncertainties of the sampled values, or NULL for all ones |
dimdeg | True iff there is a fitting degree per dimension |
mindeg | Pointer to 1 or d minimum fitting degree(s), or NULL |
maxdeg | Pointer to 1 or d maximum fitting degree(s), at least mindeg |
Any pre-set non-zero coefficients in self are overwritten or reset by the fit.
For 1D-polynomials N = 1 + maxdeg - mindeg coefficients are fitted. A non-zero mindeg ensures that the fitted polynomial has a fix-point at zero.
For multi-variate polynomials the fit depends on dimdeg:
If dimdeg is false, an n-degree coefficient is fitted iff mindeg <= n <= maxdeg. For a 2D-polynomial this means that N * (N + 1) / 2 coefficients are fitted. For a 3D-polynomial, it is N * (N + 1) * (N + 2) / 6.
If dimdeg is true, nci = 1 + maxdeg[i] - mindeg[i] coefficients are fitted for dimension i, i.e. for a 2D-polynomial N = nc1 * nc2 coefficients are fitted.
The number of distinct samples should exceed the number of coefficients to fit. The number of distinct samples may also equal the number of coefficients to fit, but in this case the fit has another meaning (any non-zero residual is due to rounding errors, not a fitting error). It is an error to try to fit more coefficients than there are distinct samples.
If the relative uncertainties of the sampled values are known, they may be passed via fitsigm. NULL means that all uncertainties equals one. This uncertainties can be applied to the least squares method by converting them to weights (w = 1.0 / sigma^2).
sampsym is ignored if mindeg is nonzero, otherwise the caller may use sampsym to indicate an a priori knowledge that the sampling positions are symmetric. NULL indicates no knowledge of such symmetry. sampsym[i] may be set to true iff the sampling is symmetric around u_i, where u_i is the mean of the sampling positions in dimension i.
In 1D this implies that the sampling points as pairs average u_0 (with an odd number of samples one sample must equal u_0). E.g. both x = (1, 2, 4, 6, 7) and x = (1, 6, 4, 2, 7) have sampling symmetry, while x = (1, 2, 4, 6) does not.
In 2D this implies that the sampling points are symmetric in the 2D-plane. For the first dimension sampling symmetry means that the sampling is line- symmetric around y = u_1, while for the second dimension, sampling symmetry implies line-symmetry around x = u_2. Point symmetry around (x,y) = (u_1, u_2) means that both sampsym[0] and sampsym[1] may be set to true.
Knowledge of symmetric sampling allows the fit to be both faster and eliminates certain round-off errors.
For higher order fitting the fitting problem known as "Runge's phenomenon" is minimized using the socalled "Chebyshev nodes" as sampling points. For Chebyshev nodes sampsym can be set to CPL_TRUE.
Warning: An increase in the polynomial degree will normally reduce the fitting 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 fitting error to increase. In some cases this happens with an increase of the polynomial degree from 8 to 9.
The fit is done in the following steps: 1) If fitsigm is non-NULL. The factors are applied to the values. 2) If mindeg is zero, the sampling positions are first transformed into Xhat_i = X_i - mean(X_i), i=1, .., dimension. 3) The Vandermonde matrix is formed from Xhat. If fitsigm is non-NULL, the weights are also taken into account. 4) The normal equations of the Vandermonde matrix is solved. 5) If mindeg is zero, the resulting polynomial in Xhat is transformed back to X.
For a univariate (1D) fit the call requires 6MN + N^3/3 + 7/2N^2 + O(M) FLOPs where M is the number of data points and where N is the number of polynomial coefficients to fit, N = 1 + maxdeg - mindeg.
For a bivariate fit the call requires MN^2 + N^3/3 + O(MN) FLOPs where M is the number of data points and where N is the number of polynomial coefficients to fit.
Examples of usage:
Possible _cpl_error_code_ set in this function:
References cpl_bivector_unwrap_vectors(), cpl_bivector_wrap_vectors(), cpl_ensure_code, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, CPL_ERROR_UNSUPPORTED_MODE, cpl_matrix_get_ncol_(), cpl_matrix_get_nrow(), cpl_polynomial_get_dimension(), CPL_SIZE_FORMAT, cpl_vector_get_size(), cpl_vector_unwrap(), and cpl_vector_wrap().
cpl_polynomial * cpl_polynomial_fit_1d_create | ( | const cpl_vector * | x_pos, |
const cpl_vector * | values, | ||
cpl_size | degree, | ||
double * | pmse | ||
) |
Fit a 1D-polynomial to a 1D-signal in a least squares sense.
x_pos | Vector of positions of the signal to fit. |
values | Vector of values of the signal to fit. |
degree | Non-negative polynomial degree. |
pmse | Iff pmse is not null, the mean squared error on success |
References CPL_ERROR_NONE, and cpl_polynomial_delete().
cpl_polynomial * cpl_polynomial_fit_2d_create | ( | const cpl_bivector * | xy_pos, |
const cpl_vector * | values, | ||
cpl_size | degree, | ||
double * | pmse | ||
) |
Fit a 2D-polynomial to a 2D-surface in a least squares sense.
xy_pos | Bivector positions of the surface to fit. |
values | Vector of values of the surface to fit. |
degree | Non-negative polynomial degree. |
pmse | Iff pmse is not null, the mean squared error on success |
References cpl_bivector_get_x_const(), cpl_bivector_get_y_const(), CPL_ERROR_NONE, cpl_polynomial_delete(), cpl_polynomial_eval(), cpl_vector_get_data_const(), cpl_vector_get_size(), cpl_vector_unwrap(), and cpl_vector_wrap().
double cpl_polynomial_get_coeff | ( | const cpl_polynomial * | self, |
const cpl_size * | pows | ||
) |
Get a coefficient of the polynomial.
self | The polynomial to process |
pows | The non-negative power(s) of the variable(s) |
Requesting the value of a coefficient that has not been set is allowed, in this case zero is returned.
Example of usage:
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NULL_INPUT, and CPL_SIZE_FORMAT.
Referenced by cpl_polynomial_extract().
cpl_size cpl_polynomial_get_degree | ( | const cpl_polynomial * | self | ) |
The degree of the polynomial.
self | The polynomial to access |
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:
References cpl_ensure, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_polynomial_extract(), and cpl_wlcalib_find_best_1d().
cpl_size cpl_polynomial_get_dimension | ( | const cpl_polynomial * | self | ) |
The dimension of the polynomial.
self | The polynomial to access |
Possible _cpl_error_code_ set in this function:
References cpl_ensure, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_image_fill_jacobian_polynomial(), cpl_image_fill_polynomial(), cpl_image_warp_polynomial(), cpl_polynomial_add(), cpl_polynomial_fit(), cpl_polynomial_multiply(), cpl_polynomial_multiply_scalar(), cpl_polynomial_shift_1d(), cpl_polynomial_subtract(), cpl_vector_fill_polynomial_fit_residual(), and cpl_wlcalib_find_best_1d().
cpl_error_code cpl_polynomial_multiply | ( | cpl_polynomial * | self, |
const cpl_polynomial * | first, | ||
const cpl_polynomial * | second | ||
) |
Multiply two polynomials of the same dimension.
self | The polynomial to hold the result |
first | The polynomial to multiply |
second | The polynomial to multiply |
Possible CPL error code set in this function:
References cpl_ensure_code, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, cpl_free(), cpl_malloc(), cpl_polynomial_delete(), and cpl_polynomial_get_dimension().
cpl_error_code cpl_polynomial_multiply_scalar | ( | cpl_polynomial * | self, |
const cpl_polynomial * | other, | ||
double | factor | ||
) |
Multiply a polynomial with a scalar.
self | The polynomial to hold the result |
other | The polynomial to scale of same dimension, may equal self |
factor | The factor to multiply with |
Possible CPL error code set in this function:
References cpl_ensure_code, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, and cpl_polynomial_get_dimension().
cpl_error_code cpl_polynomial_set_coeff | ( | cpl_polynomial * | self, |
const cpl_size * | pows, | ||
double | value | ||
) |
Set a coefficient of the polynomial.
self | The polynomial to modify |
pows | The non-negative power(s) of the variable(s) |
value | The coefficient |
The array pows is assumed to have the size of the polynomial dimension.
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, or if the new coefficient is zero, to decrease.
Setting the coefficient of x1^4 * x3^2 in the 4-dimensional polynomial poly4d to 12.3 would be performed by:
Setting the coefficient of x^3 in the 1-dimensional polynomial poly1d to 12.3 would be performed by:
When setting multiple coefficients there is a small efficiency gain when setting the highest powers first.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, and CPL_SIZE_FORMAT.
cpl_error_code cpl_polynomial_shift_1d | ( | cpl_polynomial * | p, |
cpl_size | i, | ||
double | u | ||
) |
Modify p, p(x0, x1, ..., xi, ...) := (x0, x1, ..., xi+u, ...)
p | The polynomial to be modified in place |
i | The dimension to shift (0 for first) |
u | The shift |
Shifting the polynomial p(x) = x^n with u = 1 will generate the binomial coefficients for n.
Shifting the coordinate system to (x,y) for the 2D-polynomium poly2d:
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_ACCESS_OUT_OF_RANGE, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, and cpl_polynomial_get_dimension().
Referenced by cpl_wlcalib_find_best_1d().
cpl_error_code cpl_polynomial_solve_1d | ( | const cpl_polynomial * | p, |
double | x0, | ||
double * | px, | ||
cpl_size | mul | ||
) |
A real solution to p(x) = 0 using Newton-Raphsons method.
p | The 1D-polynomial |
x0 | First guess of the solution |
px | The solution, on error see below |
mul | The root multiplicity (or 1 if unknown) |
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 (except in the case where the derivative at the first guess happens to be zero, see below - for an n-degree polynomial there are up to n-1 such guesses). If the constant term is zero, the solution 0 will be returned regardless of the first guess.
No solution is found 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). Either can happen due to an an actual lack of a real solution or due to an insufficiently good first guess. In these two cases *px is set to the value where the error occurred. In case of other errors *px is unmodified.
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 robustness 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.
Possible _cpl_error_code_ set in this function:
References CPL_ERROR_NONE.
cpl_error_code cpl_polynomial_subtract | ( | cpl_polynomial * | self, |
const cpl_polynomial * | first, | ||
const cpl_polynomial * | second | ||
) |
Subtract two polynomials of the same dimension.
self | The polynomial to hold the result |
first | The polynomial to subtract from, or NULL |
second | The polynomial to subtract, or NULL |
Possible CPL error code set in this function:
References cpl_ensure_code, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, CPL_MAX, and cpl_polynomial_get_dimension().
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 Horner's rule.
v | Preallocated vector to contain the result |
p | The 1D-polynomial |
x0 | The first point of evaluation |
d | The increment between points of evaluation |
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:
References cpl_ensure_code, CPL_ERROR_INVALID_TYPE, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, and cpl_polynomial_eval_1d().
cpl_error_code cpl_vector_fill_polynomial_fit_residual | ( | cpl_vector * | self, |
const cpl_vector * | fitvals, | ||
const cpl_vector * | fitsigm, | ||
const cpl_polynomial * | fit, | ||
const cpl_matrix * | samppos, | ||
double * | rechisq | ||
) |
Compute the residual of a polynomial fit.
self | Vector to hold the fitting residuals, fitvals may be used |
fitvals | Vector of the p fitted values |
fitsigm | Uncertainties of the sampled values or NULL for a uniform uncertainty |
fit | The fitted polynomial |
samppos | Matrix of p sample positions, with d rows and p columns |
rechisq | If non-NULL, the reduced chi square of the fit |
It is allowed to pass the same vector as both fitvals and as self, in which case fitvals is overwritten with the residuals.
If the relative uncertainties of the sampled values are known, they may be passed via fitsigm. NULL means that all uncertainties equal one. The uncertainties are taken into account when computing the reduced chi square value.
If rechisq is non-NULL, the reduced chi square of the fit is computed as well.
The mean square error, which was computed directly by the former CPL functions cpl_polynomial_fit_1d_create() and cpl_polynomial_fit_2d_create() can be computed from the fitting residual like this:
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_DATA_NOT_FOUND, CPL_ERROR_DIVISION_BY_ZERO, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, cpl_malloc(), cpl_matrix_get_ncol_(), cpl_polynomial_eval(), cpl_polynomial_eval_1d(), cpl_polynomial_get_dimension(), cpl_vector_delete(), cpl_vector_get_size(), cpl_vector_product(), cpl_vector_set_size(), and cpl_vector_wrap().