hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_interp_nd_scat_shep (e01zm)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_interp_nd_scat_shep (e01zm) generates a multidimensional interpolant to a set of scattered data points, using a modified Shepard method. When the number of dimensions is no more than five, there are corresponding functions in Chapter E01 which are specific to the given dimensionality. nag_interp_2d_scat_shep (e01sg) generates the two-dimensional interpolant, while nag_interp_3d_scat_shep (e01tg), nag_interp_4d_scat_shep (e01tk) and nag_interp_5d_scat_shep (e01tm) generate the three-, four- and five-dimensional interpolants respectively.

Syntax

[iq, rq, ifail] = e01zm(x, f, 'd', d, 'm', m, 'nw', nw, 'nq', nq)
[iq, rq, ifail] = nag_interp_nd_scat_shep(x, f, 'd', d, 'm', m, 'nw', nw, 'nq', nq)
Note: the interface to this routine has changed since earlier releases of the toolbox:
At Mark 24: nw and nq were made optional

Description

nag_interp_nd_scat_shep (e01zm) constructs a smooth function Q x , xd which interpolates a set of m scattered data points xr,fr , for r=1,2,,m, using a modification of Shepard's method. The surface is continuous and has continuous first partial derivatives.
The basic Shepard method, which is a generalization of the two-dimensional method described in Shepard (1968), interpolates the input data with the weighted mean
Q x = r=1 m wr x qr r=1 m wr x ,  
where qr = fr , wr x = 1 x - xr 2 2 .
The basic method is global in that the interpolated value at any point depends on all the data, but nag_interp_nd_scat_shep (e01zm) uses a modification (see Franke and Nielson (1980) and Renka (1988a)), whereby the method becomes local by adjusting each wr x  to be zero outside a hypersphere with centre xr  and some radius Rw. Also, to improve the performance of the basic method, each qr above is replaced by a function qr x , which is a quadratic fitted by weighted least squares to data local to xr  and forced to interpolate xr,fr . In this context, a point x is defined to be local to another point if it lies within some distance Rq of it.
The efficiency of nag_interp_nd_scat_shep (e01zm) is enhanced by using a cell method for nearest neighbour searching due to Bentley and Friedman (1979) with a cell density of 3.
The radii Rw and Rq are chosen to be just large enough to include Nw and Nq data points, respectively, for user-supplied constants Nw and Nq. Default values of these parameters are provided, and advice on alternatives is given in Choice of and .
nag_interp_nd_scat_shep (e01zm) is derived from the new implementation of QSHEP3 described by Renka (1988b). It uses the modification for high-dimensional interpolation described by Berry and Minser (1999).
Values of the interpolant Q x  generated by nag_interp_nd_scat_shep (e01zm), and its first partial derivatives, can subsequently be evaluated for points in the domain of the data by a call to nag_interp_nd_scat_shep_eval (e01zn).

References

Bentley J L and Friedman J H (1979) Data structures for range searching ACM Comput. Surv. 11 397–409
Berry M W, Minser K S (1999) Algorithm 798: high-dimensional interpolation using the modified Shepard method ACM Trans. Math. Software 25 353–366
Franke R and Nielson G (1980) Smooth interpolation of large sets of scattered data Internat. J. Num. Methods Engrg. 15 1691–1704
Renka R J (1988a) Multivariate interpolation of large sets of scattered data ACM Trans. Math. Software 14 139–148
Renka R J (1988b) Algorithm 661: QSHEP3D: Quadratic Shepard method for trivariate interpolation of scattered data ACM Trans. Math. Software 14 151–152
Shepard D (1968) A two-dimensional interpolation function for irregularly spaced data Proc. 23rd Nat. Conf. ACM 517–523 Brandon/Systems Press Inc., Princeton

Parameters

Compulsory Input Parameters

1:     xdm – double array
x1:dr must be set to the Cartesian coordinates of the data point xr, for r=1,2,,m.
Constraint: these coordinates must be distinct, and must not all lie on the same d-1-dimensional hypersurface.
2:     fm – double array
fr must be set to the data value fr, for r=1,2,,m.

Optional Input Parameters

1:     d int64int32nag_int scalar
Default: the first dimension of the array x.
d, the number of dimensions.
Constraint: d2.
2:     m int64int32nag_int scalar
Default: the dimension of the array f and the second dimension of the array x. (An error is raised if these dimensions are not equal.)
m, the number of data points.
Note: on the basis of experimental results reported in Berry and Minser (1999), when d5 it is recommended to use m4000.
Constraint: m d+1× d+2/ 2+2.
3:     nw int64int32nag_int scalar
Default: nw=-1
The number Nw of data points that determines each radius of influence Rw, appearing in the definition of each of the weights wr, for r=1,2,,m (see Description). Note that Rw is different for each weight. If nw0 the default value nw= min 2× d+1× d+2 ,m-1  is used instead.
Constraint: nwm-1.
4:     nq int64int32nag_int scalar
Default: nq=-1
The number Nq of data points to be used in the least squares fit for coefficients defining the quadratic functions qr x  (see Description). If nq0 the default value nq=mind+1×d+2×6/5,m-1 is used instead.
Constraint: nq0 or d+1×d+2/2-1nqm-1.

Output Parameters

1:     iq2×m+1 int64int32nag_int array
Integer data defining the interpolant Qx.
2:     rq* – double array
The dimension of the array rq will be d+1×d+2/2×m+2×d+1
Real data defining the interpolant Qx.
3:     ifail int64int32nag_int scalar
ifail=0 unless the function detects an error (see Error Indicators and Warnings).

Error Indicators and Warnings

Errors or warnings detected by the function:
   ifail=1
Constraint: d2.
Constraint: md+1×d+2/2+2.
Constraint: nq0 or nqd+1×d+2/2-1.
Constraint: nqm-1.
Constraint: nwm-1.
On entry, d+1×d+2/2×m+2×d+1 exceeds the largest machine integer.
   ifail=2
There are duplicate nodes in the dataset.
   ifail=3
On entry, all the data points lie on the same hypersurface. No unique solution exists.
   ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
   ifail=-399
Your licence key may have expired or may not have been installed correctly.
   ifail=-999
Dynamic memory allocation failed.

Accuracy

In experiments undertaken by Berry and Minser (1999), the accuracies obtained for a conditional function resulting in sharp functional transitions were of the order of 10-1 at best. In other cases in these experiments, the function generated interpolates the input data with maximum absolute error of the order of 10-2.

Further Comments

Timing

The time taken for a call to nag_interp_nd_scat_shep (e01zm) will depend in general on the distribution of the data points and on the choice of Nw and Nq parameters. If the data points are uniformly randomly distributed, then the time taken should be Om. At worst Om2 time will be required.

Choice of Nw and Nq

Default values of the parameters Nw and Nq may be selected by calling nag_interp_nd_scat_shep (e01zm) with nw0 and nq0. These default values may well be satisfactory for many applications.
If non-default values are required they must be supplied to nag_interp_nd_scat_shep (e01zm) through positive values of nw and nq. Increasing these argument values makes the method less local. This may increase the accuracy of the resulting interpolant at the expense of increased computational cost. The default values nw = min 2× d+1× d+2 ,m-1  and nq = min d+1× d+2× 6/5,m-1  have been chosen on the basis of experimental results reported in Renka (1988a) and Berry and Minser (1999). For further advice on the choice of these arguments see Renka (1988a) and Berry and Minser (1999).

Example

This program reads in a set of 30 data points and calls nag_interp_nd_scat_shep (e01zm) to construct an interpolating function Q x . It then calls nag_interp_nd_scat_shep_eval (e01zn) to evaluate the interpolant at a set of points.
Note that this example is not typical of a realistic problem: the number of data points would normally be very much larger.
See also Example in nag_interp_nd_scat_shep_eval (e01zn).
function e01zm_example


fprintf('e01zm example results\n\n');

% Data: 30 6-dimensional points and function values
nd = 6;
npts = 30;
x = zeros(npts,nd);
x = [0.81, 0.15, 0.44, 0.83, 0.21, 0.64;
     0.91, 0.96, 0.00, 0.09, 0.98, 0.37;
     0.13, 0.88, 0.22, 0.21, 0.73, 1.00;
     0.91, 0.49, 0.39, 0.79, 0.47, 0.71;
     0.63, 0.41, 0.72, 0.68, 0.65, 0.83;
     0.10, 0.13, 0.77, 0.47, 0.22, 0.09;
     0.28, 0.93, 0.24, 0.90, 0.96, 0.21;
     0.55, 0.01, 0.04, 0.41, 0.26, 0.79;
     0.96, 0.19, 0.95, 0.66, 0.99, 0.68;
     0.96, 0.32, 0.53, 0.96, 0.84, 0.47;
     0.16, 0.05, 0.16, 0.30, 0.58, 0.90;
     0.97, 0.14, 0.36, 0.72, 0.78, 0.06;
     0.96, 0.73, 0.28, 0.75, 0.28, 0.68;
     0.49, 0.48, 0.58, 0.19, 0.25, 0.67;
     0.80, 0.34, 0.64, 0.57, 0.08, 0.13;
     0.14, 0.24, 0.12, 0.06, 0.63, 0.89;
     0.42, 0.45, 0.03, 0.68, 0.66, 0.17;
     0.92, 0.19, 0.48, 0.67, 0.28, 0.54;
     0.79, 0.32, 0.15, 0.13, 0.40, 0.03;
     0.96, 0.26, 0.93, 0.89, 0.61, 0.81;
     0.66, 0.83, 0.41, 0.17, 0.09, 0.60;
     0.04, 0.70, 0.40, 0.54, 0.37, 0.41;
     0.85, 0.33, 0.15, 0.03, 0.36, 5.77;
     0.93, 0.58, 0.88, 0.81, 0.40, 0.66;
     0.68, 0.29, 0.88, 0.60, 0.47, 0.96;
     0.76, 0.26, 0.09, 0.41, 0.14, 0.30;
     0.74, 0.26, 0.33, 0.64, 0.36, 0.72;
     0.39, 0.68, 0.69, 0.37, 0.12, 0.75;
     0.66, 0.52, 0.17, 1.00, 0.43, 0.19;
     0.17, 0.08, 0.35, 0.71, 0.17, 0.57];
x=x';
f = [6.39;  2.50; 9.34; 7.52; 6.91; 4.68; 45.40; 5.48;  2.75;  7.43; ...
     6.05;  0.41; 8.68; 2.38; 3.70; 1.34; 15.18; 4.35;  1.50;  3.43; ...
     3.10; 14.33; 0.35; 4.30; 3.77; 4.16;  6.75; 5.22; 16.23; 10.62];

% Interpolation points
nip = 6;
xe = zeros(nd,nip);
for i=1:npts
  xe(1:nd,i) = 0.1*i;
end

% Generate the interpolant
[iq, rq, ifail] = e01zm(x, f);

% Evaluate the interpolant
[q, qx, ifail] = e01zn(x, f, iq, rq, xe);

T = array2table(double([xe' q]));
T.Properties.VariableNames = {'x_1','x_2','x_3','x_4','x_5','x_6','q'};

fprintf('Interpolated Values\n\n');
disp(T(1:nip,1:nd+1));


e01zm example results

Interpolated Values

    x_1    x_2    x_3    x_4    x_5    x_6       q    
    ___    ___    ___    ___    ___    ___    ________

    0.1    0.1    0.1    0.1    0.1    0.1     -6.5071
    0.2    0.2    0.2    0.2    0.2    0.2     -3.5708
    0.3    0.3    0.3    0.3    0.3    0.3    -0.80057
    0.4    0.4    0.4    0.4    0.4    0.4      1.7257
    0.5    0.5    0.5    0.5    0.5    0.5      3.9278
    0.6    0.6    0.6    0.6    0.6    0.6      5.8951


PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015