fitpack
Modern Fortran library for curve and surface fitting with splines
Loading...
Searching...
No Matches
fitpack_grid_surfaces Module Reference

OOP wrapper for bivariate surface fitting to data on a rectangular grid. More...

Data Types

interface  fitpack_grid_surface
 Bivariate surface fitter \( z = s(x, y) \) for gridded data. More...
 

Functions/Subroutines

integer function surface_fit_least_squares (this, smoothing, reset_knots)
 Fit a least-squares gridded surface with fixed knots.
 
integer function surface_fit_interpolating (this, reset_knots)
 Fit an interpolating gridded surface ( \( s = 0 \)).
 
integer(fp_flag) function surface_fit_automatic_knots (this, smoothing, order, keep_knots)
 Fit a smoothing gridded surface with automatic knot placement.
 
elemental subroutine surf_destroy (this)
 Destroy a grid surface object and release all allocated memory.
 
subroutine surf_new_points (this, x, y, z)
 Load new gridded data and allocate workspace.
 
type(fitpack_grid_surface) function surf_new_from_points (x, y, z, ierr)
 Construct a grid surface from gridded data and perform a default fit.
 
integer function surf_new_fit (this, x, y, z, smoothing, order)
 Load new gridded data and perform a fresh fit.
 
real(fp_real) function, dimension(size(x)) gridsurf_eval_many (this, x, y, ierr)
 Evaluate the grid surface at scattered \( (x_i, y_i) \) points.
 
real(fp_real) function gridsurf_eval_one (this, x, y, ierr)
 Evaluate the grid surface at a single \( (x, y) \) point.
 
real(fp_real) function, dimension(size(y), size(x)) gridded_eval_many (this, x, y, ierr)
 Evaluate the grid surface on a rectangular evaluation grid.
 
real(fp_real) function gridded_eval_one (this, x, y, ierr)
 Evaluate the grid surface at a single grid point via bispev.
 
real(fp_real) function, dimension(size(y), size(x)) gridded_derivatives_gridded (this, x, y, dx, dy, ierr)
 Evaluate partial derivatives on a rectangular grid.
 
real(fp_real) function, dimension(size(x)) gridded_derivatives_many (this, x, y, dx, dy, ierr)
 Evaluate partial derivatives at scattered \( (x_i, y_i) \) points.
 
real(fp_real) function gridded_derivatives_one (this, x, y, dx, dy, ierr)
 Evaluate a partial derivative at a single \( (x, y) \) point.
 
real(fp_real) function gridsurf_integral (this, lower, upper)
 Compute the double integral of the grid surface over a rectangular domain.
 
type(fitpack_curve) function gridsurf_cross_section (this, u, along_y, ierr)
 Extract a 1D cross-section curve from the grid surface.
 
type(fitpack_grid_surface) function gridsurf_derivative_spline (this, nux, nuy, ierr)
 Compute the B-spline representation of a partial derivative surface.
 
elemental integer(fp_size) function gridsurf_comm_size (this)
 Return the communication buffer size for the grid surface.
 
pure subroutine gridsurf_comm_pack (this, buffer)
 Pack grid surface data into a communication buffer.
 
pure subroutine gridsurf_comm_expand (this, buffer)
 Expand grid surface data from a communication buffer.
 

Detailed Description

OOP wrapper for bivariate surface fitting to data on a rectangular grid.

Provides fitpack_grid_surface, a derived type for fitting tensor-product B-spline surfaces \( z = s(x, y) \) to data values given on a rectangular grid \( (x_i, y_j) \). The underlying core routine is regrid, which exploits the grid structure for faster fitting than surfit. Supports evaluation, partial derivatives, integration, cross-section extraction, and derivative-spline computation.

See also
Dierckx, Ch. 5, §5.4 (pp. 98–103); regrid, bispev, parder, pardeu, dblint, profil

Function/Subroutine Documentation

◆ gridded_derivatives_gridded()

real(fp_real) function, dimension(size(y),size(x)) fitpack_grid_surfaces::gridded_derivatives_gridded ( class(fitpack_grid_surface), intent(inout) this,
real(fp_real), dimension(:), intent(in) x,
real(fp_real), dimension(:), intent(in) y,
integer(fp_size), intent(in) dx,
integer(fp_size), intent(in) dy,
integer(fp_flag), intent(out), optional ierr )
private

Evaluate partial derivatives on a rectangular grid.

See also
parder

◆ gridded_derivatives_many()

real(fp_real) function, dimension(size(x)) fitpack_grid_surfaces::gridded_derivatives_many ( class(fitpack_grid_surface), intent(inout) this,
real(fp_real), dimension(:), intent(in) x,
real(fp_real), dimension(:), intent(in) y,
integer(fp_size), intent(in) dx,
integer(fp_size), intent(in) dy,
integer(fp_flag), intent(out), optional ierr )
private

Evaluate partial derivatives at scattered \( (x_i, y_i) \) points.

See also
pardeu

◆ gridded_derivatives_one()

real(fp_real) function fitpack_grid_surfaces::gridded_derivatives_one ( class(fitpack_grid_surface), intent(inout) this,
real(fp_real), intent(in) x,
real(fp_real), intent(in) y,
integer(fp_size), intent(in) dx,
integer(fp_size), intent(in) dy,
integer(fp_flag), intent(out), optional ierr )
private

Evaluate a partial derivative at a single \( (x, y) \) point.

See also
pardeu
Here is the call graph for this function:

◆ gridded_eval_many()

real(fp_real) function, dimension(size(y),size(x)) fitpack_grid_surfaces::gridded_eval_many ( class(fitpack_grid_surface), intent(inout) this,
real(fp_real), dimension(:), intent(in) x,
real(fp_real), dimension(:), intent(in) y,
integer(fp_flag), intent(out), optional ierr )
private

Evaluate the grid surface on a rectangular evaluation grid.

Returns f(j,i) = \( s(x_i, y_j) \).

See also
bispev

◆ gridded_eval_one()

real(fp_real) function fitpack_grid_surfaces::gridded_eval_one ( class(fitpack_grid_surface), intent(inout) this,
real(fp_real), intent(in) x,
real(fp_real), intent(in) y,
integer(fp_flag), intent(out), optional ierr )
private

Evaluate the grid surface at a single grid point via bispev.

Here is the call graph for this function:

◆ gridsurf_comm_expand()

pure subroutine fitpack_grid_surfaces::gridsurf_comm_expand ( class(fitpack_grid_surface), intent(inout) this,
real(fp_comm), dimension(:), intent(in) buffer )
private

Expand grid surface data from a communication buffer.

◆ gridsurf_comm_pack()

pure subroutine fitpack_grid_surfaces::gridsurf_comm_pack ( class(fitpack_grid_surface), intent(in) this,
real(fp_comm), dimension(:), intent(out) buffer )
private

Pack grid surface data into a communication buffer.

◆ gridsurf_comm_size()

elemental integer(fp_size) function fitpack_grid_surfaces::gridsurf_comm_size ( class(fitpack_grid_surface), intent(in) this)
private

Return the communication buffer size for the grid surface.

Here is the call graph for this function:

◆ gridsurf_cross_section()

type(fitpack_curve) function fitpack_grid_surfaces::gridsurf_cross_section ( class(fitpack_grid_surface), intent(in) this,
real(fp_real), intent(in) u,
logical, intent(in) along_y,
integer(fp_flag), intent(out), optional ierr )
private

Extract a 1D cross-section curve from the grid surface.

If along_y=.true., returns \( f(y) = s(u, y) \); otherwise \( g(x) = s(x, u) \).

See also
profil

◆ gridsurf_derivative_spline()

type(fitpack_grid_surface) function fitpack_grid_surfaces::gridsurf_derivative_spline ( class(fitpack_grid_surface), intent(in) this,
integer(fp_size), intent(in) nux,
integer(fp_size), intent(in) nuy,
integer(fp_flag), intent(out), optional ierr )
private

Compute the B-spline representation of a partial derivative surface.

Returns a new grid surface with reduced degrees and trimmed knots.

See also
pardtc

◆ gridsurf_eval_many()

real(fp_real) function, dimension(size(x)) fitpack_grid_surfaces::gridsurf_eval_many ( class(fitpack_grid_surface), intent(inout) this,
real(fp_real), dimension(:), intent(in) x,
real(fp_real), dimension(size(x)), intent(in) y,
integer(fp_flag), intent(out), optional ierr )
private

Evaluate the grid surface at scattered \( (x_i, y_i) \) points.

See also
bispeu

◆ gridsurf_eval_one()

real(fp_real) function fitpack_grid_surfaces::gridsurf_eval_one ( class(fitpack_grid_surface), intent(inout) this,
real(fp_real), intent(in) x,
real(fp_real), intent(in) y,
integer(fp_flag), intent(out), optional ierr )
private

Evaluate the grid surface at a single \( (x, y) \) point.

See also
bispeu
Here is the call graph for this function:

◆ gridsurf_integral()

real(fp_real) function fitpack_grid_surfaces::gridsurf_integral ( class(fitpack_grid_surface), intent(in) this,
real(fp_real), dimension(2), intent(in) lower,
real(fp_real), dimension(2), intent(in) upper )
private

Compute the double integral of the grid surface over a rectangular domain.

See also
dblint
Here is the call graph for this function:

◆ surf_destroy()

elemental subroutine fitpack_grid_surfaces::surf_destroy ( class(fitpack_grid_surface), intent(inout) this)
private

Destroy a grid surface object and release all allocated memory.

◆ surf_new_fit()

integer function fitpack_grid_surfaces::surf_new_fit ( class(fitpack_grid_surface), intent(inout) this,
real(fp_real), dimension(:), intent(in) x,
real(fp_real), dimension(:), intent(in) y,
real(fp_real), dimension(size(y),size(x)), intent(in) z,
real(fp_real), intent(in), optional smoothing,
integer, intent(in), optional order )
private

Load new gridded data and perform a fresh fit.

See also
regrid
Here is the call graph for this function:

◆ surf_new_from_points()

type(fitpack_grid_surface) function fitpack_grid_surfaces::surf_new_from_points ( real(fp_real), dimension(:), intent(in) x,
real(fp_real), dimension(:), intent(in) y,
real(fp_real), dimension(size(y),size(x)), intent(in) z,
integer(fp_flag), intent(out), optional ierr )
private

Construct a grid surface from gridded data and perform a default fit.

See also
regrid

◆ surf_new_points()

subroutine fitpack_grid_surfaces::surf_new_points ( class(fitpack_grid_surface), intent(inout) this,
real(fp_real), dimension(:), intent(in) x,
real(fp_real), dimension(:), intent(in) y,
real(fp_real), dimension(size(y),size(x)), intent(in) z )
private

Load new gridded data and allocate workspace.

Parameters
[in,out]thisThe grid surface (destroyed and reinitialized).
[in]xGrid coordinates in the x direction.
[in]yGrid coordinates in the y direction.
[in]zFunction values z(j,i) = \( f(x_i, y_j) \).
See also
regrid

◆ surface_fit_automatic_knots()

integer(fp_flag) function fitpack_grid_surfaces::surface_fit_automatic_knots ( class(fitpack_grid_surface), intent(inout) this,
real(fp_real), intent(in), optional smoothing,
integer, intent(in), optional order,
logical, intent(in), optional keep_knots )
private

Fit a smoothing gridded surface with automatic knot placement.

Uses the regrid core routine, which exploits the rectangular grid structure for efficiency.

See also
regrid

Ensure we start with new knots (unless caller wants to keep them)

Here is the call graph for this function:

◆ surface_fit_interpolating()

integer function fitpack_grid_surfaces::surface_fit_interpolating ( class(fitpack_grid_surface), intent(inout) this,
logical, intent(in), optional reset_knots )
private

Fit an interpolating gridded surface ( \( s = 0 \)).

See also
regrid
Here is the call graph for this function:

◆ surface_fit_least_squares()

integer function fitpack_grid_surfaces::surface_fit_least_squares ( class(fitpack_grid_surface), intent(inout) this,
real(fp_real), intent(in), optional smoothing,
logical, intent(in), optional reset_knots )
private

Fit a least-squares gridded surface with fixed knots.

See also
regrid
Here is the call graph for this function: