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

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

Inheritance diagram for fitpack_surfaces::fitpack_surface:
Collaboration diagram for fitpack_surfaces::fitpack_surface:

Public Member Functions

procedure destroy (this)
 Clean memory.
 
procedure new_points (this, x, y, z, w)
 Set new points.
 
procedure new_fit (this, x, y, z, w, smoothing, order)
 Generate new fit.
 
procedure fit (this, smoothing, order, keep_knots)
 Generate/update fitting curve, with optional smoothing.
 
procedure interpolate (this, reset_knots)
 Fit an interpolating surface ( \( s = 0 \)) through the data points.
 
procedure least_squares (this, smoothing, reset_knots)
 Fit a least-squares surface with fixed knots.
 
generic eval (this, x, y, ierr)
 Evaluate surface at given coordinates.
 
generic eval (this, x, y, ierr)
 Evaluate the surface at a list of scattered \( (x_i, y_i) \) points.
 
generic eval_ongrid (this, x, y, ierr)
 Evaluate surface on a grid.
 
generic dfdx (this, x, y, dx, dy, ierr)
 Evaluate a partial derivative at a single \( (x, y) \) point.
 
generic dfdx (this, x, y, dx, dy, ierr)
 Evaluate partial derivatives at scattered \( (x_i, y_i) \) points.
 
generic dfdx_ongrid (this, x, y, dx, dy, ierr)
 Evaluate derivatives at given coordinates.
 
procedure integral (this, lower, upper)
 Double integration over a rectangular domain.
 
procedure cross_section (this, u, along_y, ierr)
 Extract a 1D cross-section curve from the surface.
 
procedure derivative_spline (this, nux, nuy, ierr)
 Compute the derivative spline coefficients.
 
procedure comm_size (this)
 Parallel communication interface.
 
procedure comm_pack (this, buffer)
 Pack surface data into a communication buffer.
 
procedure comm_expand (this, buffer)
 Expand surface data from a communication buffer.
 
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.
 
- Public Member Functions inherited from fitpack_fitters::fitpack_fitter
procedure, non_overridable mse (this)
 MSE accessor (shared by all types)
 
procedure, non_overridable core_comm_size (this)
 Base field helpers for comm (non-overridable, called by subtypes)
 
procedure, non_overridable core_comm_pack (this, buffer)
 Pack base fields into communication buffer.
 
procedure, non_overridable core_comm_expand (this, buffer)
 Expand base fields from communication buffer.
 
procedure, non_overridable destroy_base (this)
 Base field reset (called by subtype destroy methods)
 
procedure(comm_size_if), deferred comm_size (this)
 Deferred communication interface.
 
procedure(comm_pack_if), deferred comm_pack (this, buffer)
 
procedure(comm_expand_if), deferred comm_expand (this, buffer)
 

Public Attributes

integer m = 0
 The input data points.
 
real(fp_real), dimension(:), allocatable x
 
real(fp_real), dimension(:), allocatable y
 
real(fp_real), dimension(:), allocatable z
 
integer, dimension(2) order = 3
 Spline degree.
 
real(fp_real), dimension(2) left
 Interval boundaries.
 
real(fp_real), dimension(2) right
 
real(fp_real), dimension(:), allocatable w
 
integer, dimension(2) nest = 0
 
integer nmax = 0
 
integer lwrk2 = 0
 
real(fp_real), dimension(:), allocatable wrk2
 
integer bc = OUTSIDE_NEAREST_BND
 
integer, dimension(2) knots = 0
 
real(fp_real), dimension(:,:), allocatable t
 
- Public Attributes inherited from fitpack_fitters::fitpack_fitter
integer(fp_flagiopt = IOPT_NEW_SMOOTHING
 Fitting state flag.
 
real(fp_realsmoothing = 1000.0_FP_REAL
 Smoothing parameter.
 
real(fp_realfp = zero
 Weighted sum of squared residuals.
 
real(fp_real), dimension(:), allocatable c
 B-spline coefficients.
 
integer(fp_sizelwrk = 0
 Real workspace and its size.
 
real(fp_real), dimension(:), allocatable wrk
 
integer(fp_sizeliwrk = 0
 Integer workspace and its size.
 
integer(fp_size), dimension(:), allocatable iwrk
 

Private Member Functions

procedure, private surface_eval_one (this, x, y, ierr)
 Evaluate surface at given coordinates.
 
procedure, private surface_eval_many (this, x, y, ierr)
 Evaluate the surface at a list of scattered \( (x_i, y_i) \) points.
 
procedure, private surface_eval_gridded (this, x, y, ierr)
 Evaluate surface on a grid.
 
procedure, private surface_derivatives_gridded (this, x, y, dx, dy, ierr)
 Evaluate derivatives at given coordinates.
 
procedure, private surface_derivatives_many (this, x, y, dx, dy, ierr)
 Evaluate partial derivatives at scattered \( (x_i, y_i) \) points.
 
procedure, private surface_derivatives_one (this, x, y, dx, dy, ierr)
 Evaluate a partial derivative at a single \( (x, y) \) point.
 

Detailed Description

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

Stores the input data points \( (x_i, y_i, z_i) \), optional weights \( w_i \), the fitted tensor-product B-spline representation (knots and coefficients), and the fitting state. The smoothing parameter \( s \) controls the trade-off between closeness of fit and smoothness.

See also
Dierckx, Ch. 5, §5.3 (pp. 85–98)

Member Function/Subroutine Documentation

◆ comm_expand()

procedure fitpack_surfaces::fitpack_surface::comm_expand ( class(fitpack_surface), intent(inout) this,
real(fp_comm), dimension(:), intent(in) buffer )

Expand surface data from a communication buffer.

Here is the call graph for this function:

◆ comm_pack()

procedure fitpack_surfaces::fitpack_surface::comm_pack ( class(fitpack_surface), intent(in) this,
real(fp_comm), dimension(:), intent(out) buffer )

Pack surface data into a communication buffer.

◆ comm_size()

procedure fitpack_surfaces::fitpack_surface::comm_size ( class(fitpack_surface), intent(in) this)

Parallel communication interface.

◆ cross_section()

procedure fitpack_surfaces::fitpack_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 )

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

◆ derivative_spline()

procedure fitpack_surfaces::fitpack_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 )

Compute the derivative spline coefficients.

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

◆ destroy()

procedure fitpack_surfaces::fitpack_surface::destroy ( class(fitpack_surface), intent(inout) this)

Clean memory.

◆ dfdx() [1/2]

generic fitpack_surfaces::fitpack_surface::dfdx ( 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 )

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

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

◆ dfdx() [2/2]

generic fitpack_surfaces::fitpack_surface::dfdx ( 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 )

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

See also
pardeu

◆ dfdx_ongrid()

generic fitpack_surfaces::fitpack_surface::dfdx_ongrid ( 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 )

Evaluate derivatives at given coordinates.

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

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

◆ eval() [1/2]

generic fitpack_surfaces::fitpack_surface::eval ( 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 )

Evaluate surface at given coordinates.

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

◆ eval() [2/2]

generic fitpack_surfaces::fitpack_surface::eval ( 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 )

Evaluate the surface at a list of scattered \( (x_i, y_i) \) points.

See also
bispeu

◆ eval_ongrid()

generic fitpack_surfaces::fitpack_surface::eval_ongrid ( 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 )

Evaluate surface on a grid.

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

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

◆ fit()

procedure fitpack_surfaces::fitpack_surface::fit ( class(fitpack_surface), intent(inout) this,
real(fp_real), intent(in), optional smoothing,
integer, intent(in), optional order,
logical, intent(in), optional keep_knots )

Generate/update fitting curve, with optional smoothing.

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

◆ integral()

procedure fitpack_surfaces::fitpack_surface::integral ( class(fitpack_surface), intent(in) this,
real(fp_real), dimension(2), intent(in) lower,
real(fp_real), dimension(2), intent(in) upper )

Double integration 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

◆ interpolate()

procedure fitpack_surfaces::fitpack_surface::interpolate ( class(fitpack_surface), intent(inout) this,
logical, intent(in), optional reset_knots )

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

See also
surfit

◆ least_squares()

procedure fitpack_surfaces::fitpack_surface::least_squares ( class(fitpack_surface), intent(inout) this,
real(fp_real), intent(in), optional smoothing,
logical, intent(in), optional reset_knots )

Fit a least-squares surface with fixed knots.

See also
surfit

◆ new_fit()

procedure fitpack_surfaces::fitpack_surface::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 )

Generate new fit.

See also
surfit

◆ new_points()

procedure fitpack_surfaces::fitpack_surface::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 )

Set new points.

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

◆ surf_new_from_points()

type(fitpack_surface) function fitpack_surfaces::fitpack_surface::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 )

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

◆ surface_derivatives_gridded()

procedure, private fitpack_surfaces::fitpack_surface::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 derivatives at given coordinates.

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()

procedure, private fitpack_surfaces::fitpack_surface::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()

procedure, private fitpack_surfaces::fitpack_surface::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

◆ surface_eval_gridded()

procedure, private fitpack_surfaces::fitpack_surface::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 surface on a grid.

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

See also
bispev

◆ surface_eval_many()

procedure, private fitpack_surfaces::fitpack_surface::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

◆ surface_eval_one()

procedure, private fitpack_surfaces::fitpack_surface::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 surface at given coordinates.

See also
bispeu

Member Data Documentation

◆ bc

integer fitpack_surfaces::fitpack_surface::bc = OUTSIDE_NEAREST_BND

◆ knots

integer, dimension(2) fitpack_surfaces::fitpack_surface::knots = 0

◆ left

real(fp_real), dimension(2) fitpack_surfaces::fitpack_surface::left

Interval boundaries.

◆ lwrk2

integer fitpack_surfaces::fitpack_surface::lwrk2 = 0

◆ m

integer fitpack_surfaces::fitpack_surface::m = 0

The input data points.

◆ nest

integer, dimension(2) fitpack_surfaces::fitpack_surface::nest = 0

◆ nmax

integer fitpack_surfaces::fitpack_surface::nmax = 0

◆ order

integer, dimension(2) fitpack_surfaces::fitpack_surface::order = 3

Spline degree.

◆ right

real(fp_real), dimension(2) fitpack_surfaces::fitpack_surface::right

◆ t

real(fp_real), dimension(:,:), allocatable fitpack_surfaces::fitpack_surface::t

◆ w

real(fp_real), dimension(:), allocatable fitpack_surfaces::fitpack_surface::w

◆ wrk2

real(fp_real), dimension(:), allocatable fitpack_surfaces::fitpack_surface::wrk2

◆ x

real(fp_real), dimension(:), allocatable fitpack_surfaces::fitpack_surface::x

◆ y

real(fp_real), dimension(:), allocatable fitpack_surfaces::fitpack_surface::y

◆ z

real(fp_real), dimension(:), allocatable fitpack_surfaces::fitpack_surface::z

The documentation for this interface was generated from the following file: