This package contains routines for fitting (least squares) a multidimensional cubic spline to arbitrarily located data. It also contains routines for evaluating this spline (or its partial derivatives) at any point.
Coefficient calculation is performed in subroutines splcc or splcw and evaluation is performed by functions splfe or splde.
Copyright (C) 2000
University Corporation for Atmospheric Research
All Rights Reserved
The use of this Software is governed by a License Agreement.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | private, | parameter | :: | wp | = | real64 |
Real working precision if not specified [8 bytes] |
integer, | public, | parameter | :: | splpak_wp | = | wp |
Working precision |
The class contains four user entries: splcc, splcw, splfe, and splde.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | private | :: | mdim | = | 0 | ||
real(kind=wp), | private, | dimension(:), allocatable | :: | dx | |||
real(kind=wp), | private, | dimension(:), allocatable | :: | dxin | |||
integer, | private, | dimension(:), allocatable | :: | ib | |||
integer, | private, | dimension(:), allocatable | :: | ibmn | |||
integer, | private, | dimension(:), allocatable | :: | ibmx | |||
integer, | private | :: | ilast | = | 0 | ||
integer, | private | :: | isav | = | 0 | ||
integer, | private | :: | iold | = | 0 | ||
integer, | private | :: | np1 | = | 0 | ||
integer, | private | :: | l | = | 0 | ||
integer, | private | :: | il1 | = | 0 | ||
integer, | private | :: | k | = | 0 | ||
integer, | private | :: | k1 | = | 0 | ||
real(kind=wp), | private | :: | errsum | = | 0.0_wp |
generic, public :: initialize => splcc, splcw | ../../ compute the spline coefficients |
generic, public :: evaluate => splfe, splde | ../../ evaluate the spline |
procedure, public :: destroy => destroy_splpak | ../../ destory the internal class variables |
procedure, private :: splcc | |
procedure, private :: splcw | |
procedure, private :: splfe | |
procedure, private :: splde | |
procedure, private :: bascmp | |
procedure, private :: suprls |
N-dimensional cubic spline derivative evaluation.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(splpak_type), | intent(inout) | :: | me | |||
integer, | intent(in) | :: | ndim |
the dimensionality of the problem. the
spline is a function of |
||
real(kind=wp), | intent(in) | :: | x(ndim) |
an |
||
integer, | intent(in) | :: | nderiv(ndim) |
an |
||
real(kind=wp), | intent(out) | :: | coef(*) |
the array of coefficients which determine the
spline. each coefficient corresponds to a
particular basis function which in turn
corresponds to a node in the node grid. this
correspondence between the node grid and the
array |
||
real(kind=wp), | intent(in) | :: | xmin(ndim) |
a vector describing the lower extreme corner
of the node grid. a set of evenly spaced
nodes is formed along each coordinate axis
and |
||
real(kind=wp), | intent(in) | :: | xmax(ndim) |
a vector describing the upper extreme corner
of the node grid. a set of evenly spaced
nodes is formed along each coordinate axis
and |
||
integer, | intent(in) | :: | nodes(ndim) |
a vector of integers describing the number of
nodes along each axis. |
||
integer, | intent(out) | :: | ierror |
an error flag with the following meanings: |
the function value returned is the partial
derivative (indicated by nderiv
) of the
spline evaluated at x
.
N-dimensional cubic spline function evaluation.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(splpak_type), | intent(inout) | :: | me | |||
integer, | intent(in) | :: | ndim | |||
real(kind=wp), | intent(in) | :: | x(ndim) | |||
real(kind=wp), | intent(out) | :: | coef(*) | |||
real(kind=wp), | intent(in) | :: | xmin(ndim) | |||
real(kind=wp), | intent(in) | :: | xmax(ndim) | |||
integer, | intent(in) | :: | nodes(ndim) | |||
integer, | intent(out) | :: | ierror |
Destroy the internal class variables.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(splpak_type), | intent(inout) | :: | me | |||
integer, | intent(in), | optional | :: | ndim |
This routine does basis function computations for natural
splines. This routine is called by routines SPLCW and SPLDE
to compute ICOL
and BASM
, which are defined as follows:
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(splpak_type), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | x(:) | |||
integer, | intent(in) | :: | nderiv(:) | |||
real(kind=wp), | intent(in) | :: | xmin(:) | |||
integer, | intent(in) | :: | nodes(:) | |||
integer, | intent(out) | :: | icol | |||
real(kind=wp), | intent(out) | :: | basm |
To print an error number and an error message or just an error message.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | ierr |
The error number (printed only if non-zero). |
||
character(len=*), | intent(in) | :: | mess |
Message to be printed. |
N-dimensional cubic spline coefficient calculation by least squares.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(splpak_type), | intent(inout) | :: | me | |||
integer, | intent(in) | :: | ndim | |||
real(kind=wp), | intent(in) | :: | xdata(l1xdat,ndata) | |||
integer, | intent(in) | :: | l1xdat | |||
real(kind=wp), | intent(in) | :: | ydata(ndata) | |||
integer, | intent(in) | :: | ndata | |||
real(kind=wp), | intent(in) | :: | xmin(ndim) | |||
real(kind=wp), | intent(in) | :: | xmax(ndim) | |||
integer, | intent(in) | :: | nodes(ndim) | |||
real(kind=wp), | intent(in) | :: | xtrap | |||
real(kind=wp), | intent(out) | :: | coef(ncf) | |||
integer, | intent(in) | :: | ncf | |||
real(kind=wp) | :: | work(nwrk) | ||||
integer, | intent(in) | :: | nwrk | |||
integer, | intent(out) | :: | ierror |
N-dimensional cubic spline coefficient calculation by weighted least squares on arbitrarily located data.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(splpak_type), | intent(inout) | :: | me | |||
integer, | intent(in) | :: | ndim |
The dimensionality of the problem. The
spline is a function of |
||
real(kind=wp), | intent(in) | :: | xdata(l1xdat,ndata) |
A collection of locations for the data
values, i.e., points from the independent
variable space. This collection is a
2-dimensional array whose 1st dimension
indexes the |
||
integer, | intent(in) | :: | l1xdat |
The length of the 1st dimension of |
||
real(kind=wp), | intent(in) | :: | ydata(ndata) |
A collection of data values corresponding to
the points in |
||
real(kind=wp), | intent(in) | :: | wdata(:) |
A collection of weights. |
||
integer, | intent(in) | :: | ndata |
The number of data points mentioned in the above arguments. |
||
real(kind=wp), | intent(in) | :: | xmin(ndim) |
A vector describing the lower extreme corner
of the node grid. A set of evenly spaced
nodes is formed along each coordinate axis
and |
||
real(kind=wp), | intent(in) | :: | xmax(ndim) |
A vector describing the upper extreme corner
of the node grid. A set of evenly spaced
nodes is formed along each coordinate axis
and |
||
integer, | intent(in) | :: | nodes(ndim) |
A vector of integers describing the number of
nodes along each axis. |
||
real(kind=wp), | intent(in) | :: | xtrap |
A parameter to control extrapolation to data
sparse areas. The region described by |
||
real(kind=wp), | intent(out) | :: | coef(ncf) |
The array of coefficients computed by this
routine. Each coefficient corresponds to a
particular basis function which in turn
corresponds to a node in the node grid. This
correspondence between the node grid and the
array |
||
integer, | intent(in) | :: | ncf |
The length of the array |
||
real(kind=wp) | :: | work(nwrk) |
A workspace array for solving the least
squares matrix generated by this routine.
Its required size is a function of the total
number of nodes in the node grid. This
total, |
|||
integer, | intent(in) | :: | nwrk |
The length of the array |
||
integer, | intent(out) | :: | ierror |
An error flag with the following meanings: |
To determine the least squares solution of a large overdetermined linear system.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(splpak_type), | intent(inout) | :: | me | |||
integer, | intent(in) | :: | i |
the index of the row being entered. ( |
||
real(kind=wp), | intent(in) | :: | rowi(n) |
a vector which on the |
||
integer, | intent(in) | :: | n |
the length of the rows of the matrix (i.e.,
the number of columns). |
||
real(kind=wp), | intent(in) | :: | bi |
on the |
||
real(kind=wp), | intent(inout) | :: | a(nn) |
a working array which must not be changed between the successive calls to suprls. |
||
integer, | intent(in) | :: | nn |
length of scratch array |
||
real(kind=wp), | intent(out) | :: | soln(n) |
the |
||
real(kind=wp), | intent(out) | :: | err |
the euclidean norm of the residual is
returned in |
||
integer, | intent(out) | :: | ier |
error parameter. fatal errors: |