NAG CPP Interface
nagcpp::fit::dim1_spline_eval (e02bb)

1 Purpose

dim1_spline_eval evaluates a cubic spline from its B-spline representation.

2 Specification

#include "e02/nagcpp_e02bb.hpp"
template <typename LAMDA, typename C>

void function dim1_spline_eval(const LAMDA &lamda, const C &c, const double x, double &s, OptionalE02BB opt)
template <typename LAMDA, typename C>

void function dim1_spline_eval(const LAMDA &lamda, const C &c, const double x, double &s)

3 Description

dim1_spline_eval evaluates the cubic spline s(x) at a prescribed argument x from its augmented knot set λi, for i=1,2,,n+7, (see e02baf (no CPP interface)) and from the coefficients ci, for i=1,2,,qin its B-spline representation
s(x)=i=1qciNi(x).  
Here q=n¯+3, where n¯ is the number of intervals of the spline, and Ni(x) denotes the normalized B-spline of degree 3 defined upon the knots λi,λi+1,,λi+4. The prescribed argument x must satisfy λ4xλn¯+4.
It is assumed that λjλj-1, for j=2,3,,n¯+7, and λn¯+4>λ4.
If x is a point at which 4 knots coincide, s(x) is discontinuous at x; in this case, s contains the value defined as x is approached from the right.
The method employed is that of evaluation by taking convex combinations due to de Boor (1972). For further details of the algorithm and its use see Cox (1972) and Cox and Hayes (1973).
It is expected that a common use of dim1_spline_eval will be the evaluation of the cubic spline approximations produced by e02baf (no CPP interface). A generalization of dim1_spline_eval which also forms the derivative of s(x) is e02bcf (no CPP interface). e02bcf (no CPP interface) takes about 50% longer than dim1_spline_eval.

4 References

Cox M G (1972) The numerical evaluation of B-splines J. Inst. Math. Appl. 10 134–149
Cox M G (1978) The numerical evaluation of a spline from its B-spline representation J. Inst. Math. Appl. 21 135–143
Cox M G and Hayes J G (1973) Curve fitting: a guide and suite of algorithms for the non-specialist user NPL Report NAC26 National Physical Laboratory
de Boor C (1972) On calculating with B-splines J. Approx. Theory 6 50–62

5 Arguments

1: lamda(ncap7) double array Input
On entry: lamda(j-1) must be set to the value of the jth member of the complete set of knots, λj, for j=1,2,,n¯+7.
Constraint: the lamda(j-1) must be in nondecreasing order with lamda(ncap7-4)> lamda(3).
2: c(ncap7) double array Input
On entry: the coefficient ci of the B-spline Ni(x), for i=1,2,,n¯+3. The remaining elements of the array are not referenced.
3: x double Input
On entry: the argument x at which the cubic spline is to be evaluated.
Constraint: lamda(3)xlamda(ncap7-4).
4: s double Output
On exit: the value of the spline, s(x).
5: opt OptionalE02BB Input/Output
Optional parameter container, derived from Optional.

5.1Additional Quantities

1: ncap7
n-+7, where n- is the number of intervals (one greater than the number of interior knots, i.e., the knots strictly within the range λ4 to λn-+4) over which the spline is defined.

6 Exceptions and Warnings

Errors or warnings detected by the function:
All errors and warnings have an associated numeric error code field, errorid, stored either as a member of the thrown exception object (see errorid), or as a member of opt.ifail, depending on how errors and warnings are being handled (see Error Handling for more details).
Raises: ErrorException
errorid=1
On entry, x=value, ncap7 = value
and lamda[ncap7-4]=value.
Constraint: xlamda[ncap7-4].
errorid=1
On entry, x=value and
lamda[3]=value.
Constraint: xlamda[3].
errorid=2
On entry, ncap7 = value.
Constraint: ncap78.
errorid=10601
On entry, argument value must be a vector of size value array.
Supplied argument has value dimensions.
errorid=10601
On entry, argument value must be a vector of size value array.
Supplied argument was a vector of size value.
errorid=10601
On entry, argument value must be a vector of size value array.
The size for the supplied array could not be ascertained.
errorid=10602
On entry, the raw data component of value is null.
errorid=10603
On entry, unable to ascertain a value for value.
errorid=-99
An unexpected error has been triggered by this routine.
errorid=-399
Your licence key may have expired or may not have been installed correctly.
errorid=-999
Dynamic memory allocation failed.

7 Accuracy

The computed value of s(x) has negligible error in most practical situations. Specifically, this value has an absolute error bounded in modulus by 18×cmax×machine precision, where cmax is the largest in modulus of cj,cj+1,cj+2 and cj+3, and j is an integer such that λj+3xλj+4. If cj,cj+1,cj+2 and cj+3 are all of the same sign, then the computed value of s(x) has a relative error not exceeding 20×machine precision in modulus. For further details see Cox (1978).

8 Parallelism and Performance

Please see the description for the underlying computational routine in this section of the FL Interface documentation.

9 Further Comments

The time taken is approximately c×(1+0.1×log(n¯+7)) seconds, where c is a machine-dependent constant.
Note:  the function does not test all the conditions on the knots given in the description of lamda in Section 5, since to do this would result in a computation time approximately linear in n¯+7 instead of log(n¯+7). All the conditions are tested in e02baf (no CPP interface), however.

10 Example

Evaluate at nine equally-spaced points in the interval 1.0x9.0 the cubic spline with (augmented) knots 1.0, 1.0, 1.0, 1.0, 3.0, 6.0, 8.0, 9.0, 9.0, 9.0, 9.0 and normalized cubic B-spline coefficients 1.0, 2.0, 4.0, 7.0, 6.0, 4.0, 3.0.
The example program is written in a general form that will enable a cubic spline with n¯ intervals, in its normalized cubic B-spline form, to be evaluated at m equally-spaced points in the interval lamda(3)xlamda(n¯+3). The program is self-starting in that any number of datasets may be supplied.
Source FileDataResults
ex_e02bb.cppNoneex_e02bb.r