E02DHF (PDF version)
E02 Chapter Contents
E02 Chapter Introduction
NAG Library Manual

NAG Library Routine Document

E02DHF

Note:  before using this routine, please read the Users' Note for your implementation to check the interpretation of bold italicised terms and other implementation-dependent details.

+ Contents

    1  Purpose
    7  Accuracy

1  Purpose

E02DHF computes the partial derivative (of order νx, νy), of a bicubic spline approximation to a set of data values, from its B-spline representation, at points on a rectangular grid in the x-y plane. This routine may be used to calculate derivatives of a bicubic spline given in the form produced by E01DAF, E02DAF, E02DCF and E02DDF.

2  Specification

SUBROUTINE E02DHF ( MX, MY, PX, PY, X, Y, LAMDA, MU, C, NUX, NUY, Z, IFAIL)
INTEGER  MX, MY, PX, PY, NUX, NUY, IFAIL
REAL (KIND=nag_wp)  X(MX), Y(MY), LAMDA(PX), MU(PY), C((PX-4)*(PY-4)), Z(MX*MY)

3  Description

E02DHF determines the partial derivative νx + νy x νx y νy  of a smooth bicubic spline approximation sx,y at the set of data points xq,yr.
The spline is given in the B-spline representation
sx,y = i=1 nx-4 j=1 ny-4 cij Mix Njy , (1)
where Mix and Njy denote normalized cubic B-splines, the former defined on the knots λi to λi+4 and the latter on the knots μj to μj+4, with nx and ny the total numbers of knots of the computed spline with respect to the x and y variables respectively. For further details, see Hayes and Halliday (1974) for bicubic splines and de Boor (1972) for normalized B-splines. This routine is suitable for B-spline representations returned by E01DAF, E02DAF, E02DCF and E02DDF.
The partial derivatives can be up to order 2 in each direction; thus the highest mixed derivative available is 4 x2 y2 .
The points in the grid are defined by coordinates xq, for q=1,2,,mx, along the x axis, and coordinates yr, for r=1,2,,my, along the y axis.

4  References

de Boor C (1972) On calculating with B-splines J. Approx. Theory 6 50–62
Dierckx P (1981) An improved algorithm for curve fitting with spline functions Report TW54 Department of Computer Science, Katholieke Univerciteit Leuven
Dierckx P (1982) A fast algorithm for smoothing data on a rectangular grid while using spline functions SIAM J. Numer. Anal. 19 1286–1304
Hayes J G and Halliday J (1974) The least-squares fitting of cubic spline surfaces to general data sets J. Inst. Math. Appl. 14 89–103
Reinsch C H (1967) Smoothing by spline functions Numer. Math. 10 177–183

5  Parameters

1:     MX – INTEGERInput
On entry: mx, the number of grid points along the x axis.
Constraint: MX1.
2:     MY – INTEGERInput
On entry: my, the number of grid points along the y axis.
Constraint: MY1.
3:     PX – INTEGERInput
On entry: the total number of knots in the x-direction of the bicubic spline approximation, e.g., the value NX as returned by E02DCF.
4:     PY – INTEGERInput
On entry: the total number of knots in the y-direction of the bicubic spline approximation, e.g., the value NY as returned by E02DCF.
5:     X(MX) – REAL (KIND=nag_wp) arrayInput
On entry: Xq must be set to xq, the x coordinate of the qth grid point along the x axis, for q=1,2,,mx, on which values of the partial derivative are sought.
Constraint: x1<x2<<xmx.
6:     Y(MY) – REAL (KIND=nag_wp) arrayInput
On entry: Yr must be set to yr, the y coordinate of the rth grid point along the y axis, for r=1,2,,my on which values of the partial derivative are sought.
Constraint: y1<y2<<ymy.
7:     LAMDA(PX) – REAL (KIND=nag_wp) arrayInput
On entry: contains the position of the knots in the x-direction of the bicubic spline approximation to be differentiated, e.g., LAMDA as returned by E02DCF.
8:     MU(PY) – REAL (KIND=nag_wp) arrayInput
On entry: contains the position of the knots in the y-direction of the bicubic spline approximation to be differentiated, e.g., MU as returned by E02DCF.
9:     C(PX-4×PY-4) – REAL (KIND=nag_wp) arrayInput
On entry: the coefficients of the bicubic spline approximation to be differentiated, e.g., C as returned by E02DCF.
10:   NUX – INTEGERInput
On entry: specifies the order, νx of the partial derivative in the x-direction.
Constraint: 0NUX2.
11:   NUY – INTEGERInput
On entry: specifies the order, νy of the partial derivative in the y-direction.
Constraint: 0NUY2.
12:   Z(MX×MY) – REAL (KIND=nag_wp) arrayOutput
On exit: Zmy×q-1+r contains the derivative νx+νy xνx yνy s xq,yr , for q=1,2,,mx and r=1,2,,my.
13:   IFAIL – INTEGERInput/Output
On entry: IFAIL must be set to 0, -1​ or ​1. If you are unfamiliar with this parameter you should refer to Section 3.3 in the Essential Introduction for details.
For environments where it might be inappropriate to halt program execution when an error is detected, the value -1​ or ​1 is recommended. If the output of error messages is undesirable, then the value 1 is recommended. Otherwise, if you are not familiar with this parameter, the recommended value is 0. When the value -1​ or ​1 is used it is essential to test the value of IFAIL on exit.
On exit: IFAIL=0 unless the routine detects an error or a warning has been flagged (see Section 6).

6  Error Indicators and Warnings

If on entry IFAIL=0 or -1, explanatory error messages are output on the current error message unit (as defined by X04AAF).
Errors or warnings detected by the routine:
IFAIL=1
On entry,NUX<0,
orNUX>2.
IFAIL=2
On entry,NUY<0,
orNUY>2.
IFAIL=3
On entry,MX<1.
IFAIL=4
On entry,MY<1.
IFAIL=5
The x coordinates defining the mesh are not well-ordered, that is, xi-1>xi, for at least one i.
IFAIL=6
The y coordinates defining the mesh are not well-ordered, that is, yi-1>yi, for at least one i.

7  Accuracy

On successful exit, the partial derivatives on the given mesh are accurate to machine precision with respect to the supplied bicubic spline. Please refer to Section 7 in E01DAF, E02DAF, E02DCF and E02DDF of the routine document for the respective routine which calculated the spline approximant for details on the accuracy of that approximation.

8  Further Comments

None.

9  Example

This example reads in values of mx, my, xq, for q=1,2,,mx, and yr, for r=1,2,,my, followed by values of the ordinates fq,r defined at the grid points xq,yr. It then calls E02DCF to compute a bicubic spline approximation for one specified value of S. Finally it evaluates the spline and its first x derivative at a small sample of points on a rectangular grid by calling E02DHF.

9.1  Program Text

Program Text (e02dhfe.f90)

9.2  Program Data

Program Data (e02dhfe.d)

9.3  Program Results

Program Results (e02dhfe.r)


E02DHF (PDF version)
E02 Chapter Contents
E02 Chapter Introduction
NAG Library Manual

© The Numerical Algorithms Group Ltd, Oxford, UK. 2012