About News Download Benchmark Documentation FAQ Links

Fortran 77 interface

The Fortran 77 wrappers provide a very similar interface to the einspline library as in C, with a few exceptions:
Boundary condition type Code Meaning
PERIODIC 0 Use periodic boundary conditions. The value, first derivative and second derivative at the left boundary match those at the right boundary.
DERIV1 1 The value of the first derivative is specified in lVal or rVal.
DERIV2 2 The value of the second derivative is specified in lVal or rVal.
FLAT 3 The value of the first derivative is set to zero at the boundary.
NATURAL 4 The value of the second derivative is set to zero at the boundary.
Please see the documentation for the C routines for more information. The subroutine names and parameters are given below:

Nonuniform grid creation routines

The general grid is a catch-all grid type which can be constructed from any monotonic sequence of grid points. To create it, we pass an array with the grid points, the number of points to the grid FCREATE_GENERAL_GRID. A handle to the grid object is returned in the parameter grid.

SUBROUTINE FCREATE_GENERAL_GRID (points, num, grid) 
Argument Intent Type Description
points IN REAL*8 Monotonic set of points for grid with num elements
num IN REAL*8 Number of elements in points array
grid OUT INTEGER*8 Grid handle

The center grid is a grid in which the points are clustered more closely around the center than at the outsides. It is useful for splines representing functions which have rapid oscillations near the origin and are more smooth further out. Atomic wave functions have this property. A create_center_grid function takes 5 arguments, as described below. The ratio argument is the ratio of the largest spacing to the smallest spacing. Thus with a ratio of 10, the points at the center would be spaced ten times closer together than the points at the left and right extremes of the grid.
SUBROUTINE FCREATE_CENTER_GRID (start, end, ratio, num, grid) 
Argument Intent Type Description
start IN REAL*8 First grid point
end IN REAL*8 Last grid point
ratio IN REAL*8 Largest spacing divided by smallest spacing
num IN REAL*8 Number of points in the grid
grid OUT INTEGER*8 Grid handle

Once the grids are created, their handles are passed to the nonuniform spline create routines documented below.


Nonuniform grid destruction routine

The grid objects must be destroyed once they are no longer needed. If they have been passed to a spline creation routine, they should not be destroyed until after the spline has been destroyed. The following routine can be used to destroy either of the nonuniform grids.
SUBROUTINE FDESTROY_GRID (grid)
Argument Intent Type Description
grid IN INTEGER*8 Handle of grid object

Nonuniform spline creation routines

One-dimensional:

Single-precision real

SUBROUTINE FCREATE_NUBSPLINE_1D_S (x_grid, x0_code, x0_val, x1_code, x1_val, data, spline)
Argument Intent Type Description
x_grid IN INTEGER*8 Handle of x grid object
x0_code IN INTEGER Left BC type
x0_val IN REAL*4 Left BC value
x1_code IN INTEGER Right BC type
x1_val IN REAL*4 Right BC value
data IN REAL*4 Data to interpolate
spline OUT INTEGER*8 Handle for spline object

Double-precision real

SUBROUTINE FCREATE_NUBSPLINE_1D_D (x_grid, x0_code, x0_val, x1_code, x1_val, data, spline)
Argument Intent Type Description
x_grid IN INTEGER*8 Handle of x grid object
x0_code IN INTEGER Left BC type
x0_val IN REAL*8 Left BC value
x1_code IN INTEGER Right BC type
x1_val IN REAL*8 Right BC value
data IN REAL*8 Data to interpolate
spline OUT INTEGER*8 Handle for spline object

Single-precision complex

SUBROUTINE FCREATE_NUBSPLINE_1D_C (x_grid, x0_code, x0_val, x1_code, x1_val, data, spline)
Argument Intent Type Description
x_grid IN INTEGER*8 Handle of x grid object
x0_code IN INTEGER Left BC type
x0_val IN COMPLEX*8 Left BC value
x1_code IN INTEGER Right BC type
x1_val IN COMPLEX*8 Right BC value
data IN COMPLEX*8 Data to interpolate
spline OUT INTEGER*8 Handle for spline object

Double-precision complex

SUBROUTINE FCREATE_NUBSPLINE_1D_Z (x_grid, x0_code, x0_val, x1_code, x1_val, data, spline)
Argument Intent Type Description
x_grid IN INTEGER*8 Handle of x grid object
x0_code IN INTEGER Left BC type
x0_val IN COMPLEX*16 Left BC value
x1_code IN INTEGER Right BC type
x1_val IN COMPLEX*16 Right BC value
data IN COMPLEX*16 Data to interpolate
spline OUT INTEGER*8 Handle for spline object

Two-dimensional:

Single-precision real

SUBROUTINE FCREATE_NUBSPLINE_2D_S (x_grid, y_grid, 
                                  x0_code, x0_val, x1_code, x1_val, 
				  y0_code, y0_val, y1_code, y1_val, 
                                  data, spline)
Argument Intent Type Description
x_grid IN INTEGER*8 Handle of x grid object
y_grid IN INTEGER*8 Handle of y grid object
x0_code IN INTEGER Left x BC type
x0_val IN REAL*4 Left x BC value
x1_code IN INTEGER Right x BC type
x1_val IN REAL*4 Right x BC value
y0_code IN INTEGER Left y BC type
y0_val IN REAL*4 Left y BC value
y1_code IN INTEGER Right y BC type
y1_val IN REAL*4 Right y BC value
data IN REAL*4 Data to interpolate
spline OUT INTEGER*8 Handle for spline object

Double-precision real

SUBROUTINE FCREATE_NUBSPLINE_2D_D (x_grid, y_grid, 
                                  x0_code, x0_val, x1_code, x1_val, 
				  y0_code, y0_val, y1_code, y1_val, 
                                  data, spline)
Argument Intent Type Description
x_grid IN INTEGER*8 Handle of x grid object
y_grid IN INTEGER*8 Handle of y grid object
x0_code IN INTEGER Left x BC type
x0_val IN REAL*8 Left x BC value
x1_code IN INTEGER Right x BC type
x1_val IN REAL*8 Right x BC value
y0_code IN INTEGER Left y BC type
y0_val IN REAL*8 Left y BC value
y1_code IN INTEGER Right y BC type
y1_val IN REAL*8 Right y BC value
data IN REAL*8 Data to interpolate
spline OUT INTEGER*8 Handle for spline object

Single-precision complex

SUBROUTINE FCREATE_NUBSPLINE_2D_C (x_grid, y_grid, 
                                  x0_code, x0_val, x1_code, x1_val, 
				  y0_code, y0_val, y1_code, y1_val, 
                                  data, spline)
Argument Intent Type Description
x_grid IN INTEGER*8 Handle of x grid object
y_grid IN INTEGER*8 Handle of y grid object
x0_code IN INTEGER Left x BC type
x0_val IN COMPLEX*8 Left x BC value
x1_code IN INTEGER Right x BC type
x1_val IN COMPLEX*8 Right x BC value
y0_code IN INTEGER Left y BC type
y0_val IN COMPLEX*8 Left y BC value
y1_code IN INTEGER Right y BC type
y1_val IN COMPLEX*8 Right y BC value
data IN COMPLEX*8 Data to interpolate
spline OUT INTEGER*8 Handle for spline object

Double-precision complex

SUBROUTINE FCREATE_NUBSPLINE_2D_Z (x_grid, y_grid, 
                                  x0_code, x0_val, x1_code, x1_val, 
				  y0_code, y0_val, y1_code, y1_val, 
                                  data, spline)
Argument Intent Type Description
x_grid IN INTEGER*8 Handle of x grid object
y_grid IN INTEGER*8 Handle of y grid object
x0_code IN INTEGER Left x BC type
x0_val IN COMPLEX*16 Left x BC value
x1_code IN INTEGER Right x BC type
x1_val IN COMPLEX*16 Right x BC value
y0_code IN INTEGER Left y BC type
y0_val IN COMPLEX*16 Left y BC value
y1_code IN INTEGER Right y BC type
y1_val IN COMPLEX*16 Right y BC value
data IN COMPLEX*16 Data to interpolate
spline OUT INTEGER*8 Handle for spline object

Three-dimensional:

Single-precision real

SUBROUTINE FCREATE_NUBSPLINE_3D_S (x_grid, y_grid, z_grid,
                                  x0_code, x0_val, x1_code, x1_val, 
				  y0_code, y0_val, y1_code, y1_val, 
				  z0_code, z0_val, z1_code, z1_val, 
                                  data, spline)
Argument Intent Type Description
x_grid IN INTEGER*8 Handle of x grid object
y_grid IN INTEGER*8 Handle of y grid object
z_grid IN INTEGER*8 Handle of z grid object
num_z IN INTEGER # of z grid points
x0_code IN INTEGER Left x BC type
x0_val IN REAL*4 Left x BC value
x1_code IN INTEGER Right x BC type
x1_val IN REAL*4 Right x BC value
y0_code IN INTEGER Left y BC type
y0_val IN REAL*4 Left y BC value
y1_code IN INTEGER Right y BC type
y1_val IN REAL*4 Right y BC value
z0_code IN INTEGER Left z BC type
z0_val IN REAL*4 Left z BC value
z1_code IN INTEGER Right z BC type
z1_val IN REAL*4 Right z BC value
data IN REAL*4 Data to interpolate
spline OUT INTEGER*8 Handle for spline object

Double-precision real

SUBROUTINE FCREATE_NUBSPLINE_3D_D (x_grid, y_grid, z_grid,
                                  x0_code, x0_val, x1_code, x1_val, 
				  y0_code, y0_val, y1_code, y1_val, 
				  z0_code, z0_val, z1_code, z1_val, 
                                  data, spline)
Argument Intent Type Description
x_grid IN INTEGER*8 Handle of x grid object
y_grid IN INTEGER*8 Handle of y grid object
z_grid IN INTEGER*8 Handle of z grid object
x0_code IN INTEGER Left x BC type
x0_val IN REAL*8 Left x BC value
x1_code IN INTEGER Right x BC type
x1_val IN REAL*8 Right x BC value
y0_code IN INTEGER Left y BC type
y0_val IN REAL*8 Left y BC value
y1_code IN INTEGER Right y BC type
y1_val IN REAL*8 Right y BC value
z0_code IN INTEGER Left z BC type
z0_val IN REAL*8 Left z BC value
z1_code IN INTEGER Right z BC type
z1_val IN REAL*8 Right z BC value
data IN REAL*8 Data to interpolate
spline OUT INTEGER*8 Handle for spline object

Single-precision complex

SUBROUTINE FCREATE_NUBSPLINE_3D_C (x_grid, y_grid, z_grid,
                                  x0_code, x0_val, x1_code, x1_val, 
				  y0_code, y0_val, y1_code, y1_val, 
				  z0_code, z0_val, z1_code, z1_val, 
                                  data, spline)
Argument Intent Type Description
x_grid IN INTEGER*8 Handle of x grid object
y_grid IN INTEGER*8 Handle of y grid object
z_grid IN INTEGER*8 Handle of z grid object
x0_code IN INTEGER Left x BC type
x0_val IN COMPLEX*8 Left x BC value
x1_code IN INTEGER Right x BC type
x1_val IN COMPLEX*8 Right x BC value
y0_code IN INTEGER Left y BC type
y0_val IN COMPLEX*8 Left y BC value
y1_code IN INTEGER Right y BC type
y1_val IN COMPLEX*8 Right y BC value
z0_code IN INTEGER Left z BC type
z0_val IN COMPLEX*8 Left z BC value
z1_code IN INTEGER Right z BC type
z1_val IN COMPLEX*8 Right z BC value
data IN COMPLEX*8 Data to interpolate
spline OUT INTEGER*8 Handle for spline object

Double-precision complex

SUBROUTINE FCREATE_NUBSPLINE_3D_Z (x_grid, y_grid, z_grid,
                                  x0_code, x0_val, x1_code, x1_val, 
				  y0_code, y0_val, y1_code, y1_val, 
				  z0_code, z0_val, z1_code, z1_val, 
                                  data, spline)
Argument Intent Type Description
x_grid IN INTEGER*8 Handle of x grid object
y_grid IN INTEGER*8 Handle of y grid object
z_grid IN INTEGER*8 Handle of z grid object
x0_code IN INTEGER Left x BC type
x0_val IN COMPLEX*16 Left x BC value
x1_code IN INTEGER Right x BC type
x1_val IN COMPLEX*16 Right x BC value
y0_code IN INTEGER Left y BC type
y0_val IN COMPLEX*16 Left y BC value
y1_code IN INTEGER Right y BC type
y1_val IN COMPLEX*16 Right y BC value
z0_code IN INTEGER Left z BC type
z0_val IN COMPLEX*16 Left z BC value
z1_code IN INTEGER Right z BC type
z1_val IN COMPLEX*16 Right z BC value
data IN COMPLEX*16 Data to interpolate
spline OUT INTEGER*8 Handle for spline object

Spline destruction routine

The following subroutine can be used to deallocate the memory for any Bspline object. Note that in the nonuniform case, the grid objects must be destroyed after the splines that refer to them.

SUBROUTINE FDESTROY_BSPLINE (spline) 
Argument Intent Type Description
spline IN INTEGER*8 Spline object handle

Uniform spline evaulation routines

One-dimensional

SUBROUTINE FEVAL_NUBSPLINE_1D_S     (spline, x, val)
SUBROUTINE FEVAL_NUBSPLINE_1D_D     (spline, x, val)
SUBROUTINE FEVAL_NUBSPLINE_1D_C     (spline, x, val)
SUBROUTINE FEVAL_NUBSPLINE_1D_Z     (spline, x, val)

SUBROUTINE FEVAL_NUBSPLINE_1D_S_VG  (spline, x, val, grad)
SUBROUTINE FEVAL_NUBSPLINE_1D_D_VG  (spline, x, val, grad)
SUBROUTINE FEVAL_NUBSPLINE_1D_C_VG  (spline, x, val, grad)
SUBROUTINE FEVAL_NUBSPLINE_1D_Z_VG  (spline, x, val, grad)

SUBROUTINE FEVAL_NUBSPLINE_1D_S_VGL (spline, x, val, grad, lapl)
SUBROUTINE FEVAL_NUBSPLINE_1D_D_VGL (spline, x, val, grad, lapl)
SUBROUTINE FEVAL_NUBSPLINE_1D_C_VGL (spline, x, val, grad, lapl)
SUBROUTINE FEVAL_NUBSPLINE_1D_Z_VGL (spline, x, val, grad, lapl)

SUBROUTINE FEVAL_NUBSPLINE_1D_S_VGH (spline, x, val, grad, hess)
SUBROUTINE FEVAL_NUBSPLINE_1D_D_VGH (spline, x, val, grad, hess)
SUBROUTINE FEVAL_NUBSPLINE_1D_C_VGH (spline, x, val, grad, hess)
SUBROUTINE FEVAL_NUBSPLINE_1D_Z_VGH (spline, x, val, grad, hess)
Argument Intent _S Type _D Type _C Type _Z Type Description
spline IN INTEGER*8 INTEGER*8 INTEGER*8 INTEGER*8 Spline handle
x IN REAL*8 REAL*8 REAL*8 REAL*8 Interpolation position
val OUT REAL*4 REAL*8 COMPLEX*8 COMPLEX*16 Interpolated value
grad OUT REAL*4 REAL*8 COMPLEX*8 COMPLEX*16 Interpolated first derivative
lapl OUT REAL*4 REAL*8 COMPLEX*8 COMPLEX*16 Interpolated second derivative
hess OUT REAL*4 REAL*8 COMPLEX*8 COMPLEX*16 Interpolated second derivative

Two-dimensional

SUBROUTINE FEVAL_NUBSPLINE_2D_S     (spline, x, y, val)
SUBROUTINE FEVAL_NUBSPLINE_2D_D     (spline, x, y, val)
SUBROUTINE FEVAL_NUBSPLINE_2D_C     (spline, x, y, val)
SUBROUTINE FEVAL_NUBSPLINE_2D_Z     (spline, x, y, val)

SUBROUTINE FEVAL_NUBSPLINE_2D_S_VG  (spline, x, y, val, grad)
SUBROUTINE FEVAL_NUBSPLINE_2D_D_VG  (spline, x, y, val, grad)
SUBROUTINE FEVAL_NUBSPLINE_2D_C_VG  (spline, x, y, val, grad)
SUBROUTINE FEVAL_NUBSPLINE_2D_Z_VG  (spline, x, y, val, grad)

SUBROUTINE FEVAL_NUBSPLINE_2D_S_VGL (spline, x, y, val, grad, lapl)
SUBROUTINE FEVAL_NUBSPLINE_2D_D_VGL (spline, x, y, val, grad, lapl)
SUBROUTINE FEVAL_NUBSPLINE_2D_C_VGL (spline, x, y, val, grad, lapl)
SUBROUTINE FEVAL_NUBSPLINE_2D_Z_VGL (spline, x, y, val, grad, lapl)

SUBROUTINE FEVAL_NUBSPLINE_2D_S_VGH (spline, x, y, val, grad, hess)
SUBROUTINE FEVAL_NUBSPLINE_2D_D_VGH (spline, x, y, val, grad, hess)
SUBROUTINE FEVAL_NUBSPLINE_2D_C_VGH (spline, x, y, val, grad, hess)
SUBROUTINE FEVAL_NUBSPLINE_2D_Z_VGH (spline, x, y, val, grad, hess)
Argument Intent _S Type _D Type _C Type _Z Type Description
spline IN INTEGER*8 INTEGER*8 INTEGER*8 INTEGER*8 Spline handle
x IN REAL*8 REAL*8 REAL*8 REAL*8 x coordinate for interpolation
y IN REAL*8 REAL*8 REAL*8 REAL*8 y coordinate for interpolation
val OUT REAL*4 REAL*8 COMPLEX*8 COMPLEX*16 Interpolated value
grad OUT REAL*4 REAL*8 COMPLEX*8 COMPLEX*16 Interpolated gradient (2 elements)
lapl OUT REAL*4 REAL*8 COMPLEX*8 COMPLEX*16 Interpolated Laplacian
hess OUT REAL*4 REAL*8 COMPLEX*8 COMPLEX*16 Interpolated Hessian (4 elements)

Three-dimensional

SUBROUTINE FEVAL_NUBSPLINE_3D_S     (spline, x, y, z, val)
SUBROUTINE FEVAL_NUBSPLINE_3D_D     (spline, x, y, z, val)
SUBROUTINE FEVAL_NUBSPLINE_3D_C     (spline, x, y, z, val)
SUBROUTINE FEVAL_NUBSPLINE_3D_Z     (spline, x, y, z, val)

SUBROUTINE FEVAL_NUBSPLINE_3D_S_VG  (spline, x, y, z, val, grad)
SUBROUTINE FEVAL_NUBSPLINE_3D_D_VG  (spline, x, y, z, val, grad)
SUBROUTINE FEVAL_NUBSPLINE_3D_C_VG  (spline, x, y, z, val, grad)
SUBROUTINE FEVAL_NUBSPLINE_3D_Z_VG  (spline, x, y, z, val, grad)

SUBROUTINE FEVAL_NUBSPLINE_3D_S_VGL (spline, x, y, z, val, grad, lapl)
SUBROUTINE FEVAL_NUBSPLINE_3D_D_VGL (spline, x, y, z, val, grad, lapl)
SUBROUTINE FEVAL_NUBSPLINE_3D_C_VGL (spline, x, y, z, val, grad, lapl)
SUBROUTINE FEVAL_NUBSPLINE_3D_Z_VGL (spline, x, y, z, val, grad, lapl)

SUBROUTINE FEVAL_NUBSPLINE_3D_S_VGH (spline, x, y, z, val, grad, hess)
SUBROUTINE FEVAL_NUBSPLINE_3D_D_VGH (spline, x, y, z, val, grad, hess)
SUBROUTINE FEVAL_NUBSPLINE_3D_C_VGH (spline, x, y, z, val, grad, hess)
SUBROUTINE FEVAL_NUBSPLINE_3D_Z_VGH (spline, x, y, z, val, grad, hess)
Argument Intent _S Type _D Type _C Type _Z Type Description
spline IN INTEGER*8 INTEGER*8 INTEGER*8 INTEGER*8 Spline handle
x IN REAL*8 REAL*8 REAL*8 REAL*8 x coordinate for interpolation
y IN REAL*8 REAL*8 REAL*8 REAL*8 y coordinate for interpolation
z IN REAL*8 REAL*8 REAL*8 REAL*8 z coordinate for interpolation
val OUT REAL*4 REAL*8 COMPLEX*8 COMPLEX*16 Interpolated value
grad OUT REAL*4 REAL*8 COMPLEX*8 COMPLEX*16 Interpolated gradient (3 elements)
lapl OUT REAL*4 REAL*8 COMPLEX*8 COMPLEX*16 Interpolated Laplacian
hess OUT REAL*4 REAL*8 COMPLEX*8 COMPLEX*16 Interpolated Hessian (9 elements)


SourceForge.net Logo