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

OOP wrapper for bivariate surface fitting to scattered data. More...

Data Types

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

Functions/Subroutines

integer(fp_flag) function surface_fit_automatic_knots (this, smoothing, order, keep_knots)
 Fit a smoothing surface with automatic knot placement.
 
integer(fp_flag) function surface_fit_interpolating (this, reset_knots)
 Fit an interpolating surface ( \( s = 0 \)) through the data points.
 
integer(fp_flag) function surface_fit_least_squares (this, smoothing, reset_knots)
 Fit a least-squares surface with fixed knots.
 
real(fp_real) function, dimension(size(x)) surface_eval_many (this, x, y, ierr)
 Evaluate the surface at a list of scattered \( (x_i, y_i) \) points.
 
real(fp_real) function surface_eval_one (this, x, y, ierr)
 Evaluate the surface at a single \( (x, y) \) point.
 
real(fp_real) function, dimension(size(y), size(x)) surface_eval_gridded (this, x, y, ierr)
 Evaluate the surface on a rectangular grid \( x_i \times y_j \).
 
elemental subroutine surf_destroy (this)
 Destroy a surface object and release all allocated memory.
 
subroutine surf_new_points (this, x, y, z, w)
 Load new scattered data points and allocate workspace for surface fitting.
 
type(fitpack_surface) function surf_new_from_points (x, y, z, w, ierr)
 Construct a surface from scattered data points and perform a default fit.
 
integer(fp_flag) function surf_new_fit (this, x, y, z, w, smoothing, order)
 Load new data points and perform a fresh surface fit.
 
real(fp_real) function, dimension(size(y), size(x)) surface_derivatives_gridded (this, x, y, dx, dy, ierr)
 Evaluate partial derivatives on a rectangular grid.
 
real(fp_real) function, dimension(size(x)) surface_derivatives_many (this, x, y, dx, dy, ierr)
 Evaluate partial derivatives at scattered \( (x_i, y_i) \) points.
 
real(fp_real) function surface_derivatives_one (this, x, y, dx, dy, ierr)
 Evaluate a partial derivative at a single \( (x, y) \) point.
 
real(fp_real) function surface_integral (this, lower, upper)
 Compute the double integral of the surface over a rectangular domain.
 
type(fitpack_curve) function surface_cross_section (this, u, along_y, ierr)
 Extract a 1D cross-section curve from the surface.
 
type(fitpack_surface) function surface_derivative_spline (this, nux, nuy, ierr)
 Compute the B-spline representation of a partial derivative surface.
 
pure real(fp_real) function, dimension(n) merge_knots (along_y, tx, ty, n)
 Select knots from the appropriate direction for cross-section extraction.
 
elemental integer(fp_size) function surf_comm_size (this)
 Return the communication buffer size for the surface.
 
pure subroutine surf_comm_pack (this, buffer)
 Pack surface data into a communication buffer.
 
pure subroutine surf_comm_expand (this, buffer)
 Expand surface data from a communication buffer.
 

Detailed Description

OOP wrapper for bivariate surface fitting to scattered data.

Provides fitpack_surface, a derived type for fitting tensor-product B-spline surfaces \( z = s(x, y) \) to scattered \( (x_i, y_i, z_i) \) data. The underlying core routine is surfit, which uses automatic knot placement with a smoothing parameter. Supports evaluation, partial derivatives, integration, and cross-section extraction.

See also
Dierckx, Ch. 5, §5.3 (pp. 85–98); surfit, bispev, parder, dblint, profil

Function/Subroutine Documentation

◆ merge_knots()

pure real(fp_real) function, dimension(n) fitpack_surfaces::merge_knots ( logical, intent(in) along_y,
real(fp_real), dimension(:), intent(in) tx,
real(fp_real), dimension(:), intent(in) ty,
integer(fp_size), intent(in) n )
private

Select knots from the appropriate direction for cross-section extraction.

◆ surf_comm_expand()

pure subroutine fitpack_surfaces::surf_comm_expand ( class(fitpack_surface), intent(inout) this,
real(fp_comm), dimension(:), intent(in) buffer )
private

Expand surface data from a communication buffer.

◆ surf_comm_pack()

pure subroutine fitpack_surfaces::surf_comm_pack ( class(fitpack_surface), intent(in) this,
real(fp_comm), dimension(:), intent(out) buffer )
private

Pack surface data into a communication buffer.

◆ surf_comm_size()

elemental integer(fp_size) function fitpack_surfaces::surf_comm_size ( class(fitpack_surface), intent(in) this)
private

Return the communication buffer size for the surface.

Here is the call graph for this function:

◆ surf_destroy()

elemental subroutine fitpack_surfaces::surf_destroy ( class(fitpack_surface), intent(inout) this)
private

Destroy a surface object and release all allocated memory.

◆ surf_new_fit()

integer(fp_flag) function fitpack_surfaces::surf_new_fit ( class(fitpack_surface), intent(inout) this,
real(fp_real), dimension(:), intent(in) x,
real(fp_real), dimension(size(x)), intent(in) y,
real(fp_real), dimension(size(x)), intent(in) z,
real(fp_real), dimension(size(x)), intent(in), optional w,
real(fp_real), intent(in), optional smoothing,
integer, intent(in), optional order )
private

Load new data points and perform a fresh surface fit.

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

◆ surf_new_from_points()

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

Construct a surface from scattered data points and perform a default fit.

Parameters
[in]xX coordinates.
[in]yY coordinates, same size as x.
[in]zFunction values, same size as x.
[in]wOptional weights.
[out]ierrOptional error flag.
Returns
Fitted surface object.
See also
surfit

◆ surf_new_points()

subroutine fitpack_surfaces::surf_new_points ( class(fitpack_surface), intent(inout) this,
real(fp_real), dimension(:), intent(in) x,
real(fp_real), dimension(size(x)), intent(in) y,
real(fp_real), dimension(size(x)), intent(in) z,
real(fp_real), dimension(size(x)), intent(in), optional w )
private

Load new scattered data points and allocate workspace for surface fitting.

Parameters
[in,out]thisThe surface object (destroyed and reinitialized).
[in]xX coordinates of scattered data points.
[in]yY coordinates, same size as x.
[in]zFunction values, same size as x.
[in]wOptional positive weights, same size as x.
See also
surfit

◆ surface_cross_section()

type(fitpack_curve) function fitpack_surfaces::surface_cross_section ( class(fitpack_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 surface.

If along_y=.true., returns \( f(y) = s(u, y) \); otherwise \( g(x) = s(x, u) \). The result is a fitpack_curve with the appropriate knots and coefficients.

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

◆ surface_derivative_spline()

type(fitpack_surface) function fitpack_surfaces::surface_derivative_spline ( class(fitpack_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 surface representing \( \partial^{n_x+n_y} s / \partial x^{n_x} \partial y^{n_y} \), with reduced degrees \( (k_x - n_x, k_y - n_y) \) and trimmed knots.

See also
pardtc

◆ surface_derivatives_gridded()

real(fp_real) function, dimension(size(y),size(x)) fitpack_surfaces::surface_derivatives_gridded ( class(fitpack_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.

Returns f(j,i) = \( \partial^{dx+dy} s / \partial x^{dx} \partial y^{dy} \) at \( (x_i, y_j) \).

See also
parder

◆ surface_derivatives_many()

real(fp_real) function, dimension(size(x)) fitpack_surfaces::surface_derivatives_many ( class(fitpack_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

◆ surface_derivatives_one()

real(fp_real) function fitpack_surfaces::surface_derivatives_one ( class(fitpack_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:

◆ surface_eval_gridded()

real(fp_real) function, dimension(size(y),size(x)) fitpack_surfaces::surface_eval_gridded ( class(fitpack_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 surface on a rectangular grid \( x_i \times y_j \).

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

See also
bispev

◆ surface_eval_many()

real(fp_real) function, dimension(size(x)) fitpack_surfaces::surface_eval_many ( class(fitpack_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 surface at a list of scattered \( (x_i, y_i) \) points.

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

◆ surface_eval_one()

real(fp_real) function fitpack_surfaces::surface_eval_one ( class(fitpack_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 surface at a single \( (x, y) \) point.

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

◆ surface_fit_automatic_knots()

integer(fp_flag) function fitpack_surfaces::surface_fit_automatic_knots ( class(fitpack_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 surface with automatic knot placement.

Uses the surfit core routine. The smoothing parameter controls the trade-off between closeness of fit and smoothness.

Parameters
[in,out]thisThe surface object.
[in]smoothingOptional smoothing factor.
[in]orderOptional spline degree (applied to both directions).
[in]keep_knotsIf .true., continue from current knot set.
Returns
Error flag.
See also
surfit

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

Here is the call graph for this function:

◆ surface_fit_interpolating()

integer(fp_flag) function fitpack_surfaces::surface_fit_interpolating ( class(fitpack_surface), intent(inout) this,
logical, intent(in), optional reset_knots )
private

Fit an interpolating surface ( \( s = 0 \)) through the data points.

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

◆ surface_fit_least_squares()

integer(fp_flag) function fitpack_surfaces::surface_fit_least_squares ( class(fitpack_surface), intent(inout) this,
real(fp_real), intent(in), optional smoothing,
logical, intent(in), optional reset_knots )
private

Fit a least-squares surface with fixed knots.

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

◆ surface_integral()

real(fp_real) function fitpack_surfaces::surface_integral ( class(fitpack_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 surface over a rectangular domain.

Computes \( \int_{x_a}^{x_b} \int_{y_a}^{y_b} s(x,y)\,dy\,dx \).

Parameters
[in]thisThe fitted surface.
[in]lowerLower bounds \( (x_a, y_a) \).
[in]upperUpper bounds \( (x_b, y_b) \).
Returns
The integral value.
See also
dblint
Here is the call graph for this function: