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

Univariate spline curve fitting module. More...

Data Types

interface  fitpack_curve
 Univariate spline curve \( y = s(x) \). More...
 
type  fitpack_periodic_curve
 Derived type describing a periodic curve. No changes are made to the storage, but the appropriate package functions will be called depending on the type. More...
 

Functions/Subroutines

type(fitpack_curve) function new_from_points (x, y, w, ierr)
 Construct a fitpack_curve from data points and perform a default fit.
 
integer(fp_flag) function new_fit (this, x, y, w, smoothing, order)
 Load new data points and perform a fresh fit.
 
elemental subroutine destroy (this)
 Destroy a curve object and release all allocated memory.
 
subroutine new_points (this, x, y, w)
 Load new data points into the curve, resetting all internal state.
 
real(fp_real) function curve_eval_one (this, x, ierr)
 Evaluate the spline at a single point (with error flag).
 
pure real(fp_real) function curve_eval_one_noerr (this, x)
 Evaluate the spline at a single point (pure, returns NaN on error).
 
real(fp_real) function, dimension(size(x)) curve_eval_many (this, x, ierr)
 Evaluate the spline at multiple points (with error flag).
 
pure real(fp_real) function, dimension(size(x)) curve_eval_many_pure (this, x)
 Evaluate the spline at multiple points (pure, returns NaN on error).
 
integer(fp_flag) function interpolating_curve (this, order, reset_knots)
 Fit an interpolating spline ( \( s = 0 \)) through the data points.
 
integer(fp_flag) function curve_fit_least_squares (this, smoothing, reset_knots)
 Fit a least-squares spline with fixed knots.
 
integer(fp_flag) function curve_fit_automatic_knots (this, smoothing, order, keep_knots)
 Fit a smoothing spline with automatic knot placement (curfit / percur).
 
real(fp_real) function, dimension(size(x)) curve_derivatives (this, x, order, ierr)
 Evaluate the k-th derivative at multiple points (with error flag).
 
real(fp_real) function, dimension(this%order+1) curve_all_derivatives (this, x, ierr)
 Evaluate all derivatives \( s^{(j)}(x) \) for \( j = 0, \ldots, k \) (with error flag).
 
real(fp_real) function, dimension(this%order+1) curve_all_derivatives_pure (this, x)
 Evaluate all derivatives \( s^{(j)}(x) \) for \( j = 0, \ldots, k \) (pure, NaN on error).
 
real(fp_real) function curve_derivative (this, x, order, ierr)
 Evaluate the k-th derivative at a single point (with optional error flag).
 
real(fp_real) function integral (this, from, to)
 Compute the definite integral \( \int_a^b s(x)\,dx \) via splint.
 
subroutine fourier_coefficients (this, alpha, a, b, ierr)
 Compute Fourier sine and cosine coefficients of the spline.
 
real(fp_real) function, dimension(:), allocatable zeros (this, ierr)
 Find all real zeros of a cubic spline \( s(x) = 0 \).
 
subroutine grow_knot_storage (this)
 Grow the knot and coefficient arrays when storage is insufficient.
 
subroutine curve_insert_knot_one (this, x, ierr)
 Insert a single knot into the B-spline representation.
 
subroutine curve_insert_knot_many (this, x, ierr)
 Insert multiple knots into the B-spline representation.
 
elemental integer(fp_size) function curve_comm_size (this)
 Return the communication buffer size (number of FP_COMM elements).
 
pure subroutine curve_comm_pack (this, buffer)
 Pack curve data into a communication buffer for parallel transfer.
 
pure subroutine curve_comm_expand (this, buffer)
 Expand curve data from a communication buffer.
 

Variables

integer, parameter max_k = 5
 

Detailed Description

Univariate spline curve fitting module.

Provides fitpack_curve for fitting and evaluating a spline \( y = s(x) \) to data, and fitpack_periodic_curve for periodic splines. Wraps curfit / percur for fitting, splev for evaluation, splder / spalde for derivatives, splint for integration, sproot for zero-finding, and fourco for Fourier coefficients.

See also
Dierckx, Ch. 4–5 (smoothing curves), Ch. 6 §6.1–6.2 (periodic curves)

Function/Subroutine Documentation

◆ curve_all_derivatives()

real(fp_real) function, dimension(this%order+1) fitpack_curves::curve_all_derivatives ( class(fitpack_curve), intent(inout) this,
real(fp_real), intent(in) x,
integer(fp_flag), intent(out) ierr )
private

Evaluate all derivatives \( s^{(j)}(x) \) for \( j = 0, \ldots, k \) (with error flag).

Parameters
[in,out]thisThe fitted curve object.
[in]xScalar evaluation point.
[out]ierrError flag.
Returns
Array of size \( k+1 \) containing \( s(x), s'(x), \ldots, s^{(k)}(x) \).
See also
spalde

◆ curve_all_derivatives_pure()

real(fp_real) function, dimension(this%order+1) fitpack_curves::curve_all_derivatives_pure ( class(fitpack_curve), intent(in) this,
real(fp_real), intent(in) x )
private

Evaluate all derivatives \( s^{(j)}(x) \) for \( j = 0, \ldots, k \) (pure, NaN on error).

Parameters
[in]thisThe fitted curve object.
[in]xScalar evaluation point.
Returns
Array of size \( k+1 \) containing \( s(x), s'(x), \ldots, s^{(k)}(x) \).
See also
spalde

◆ curve_comm_expand()

pure subroutine fitpack_curves::curve_comm_expand ( class(fitpack_curve), intent(inout) this,
real(fp_comm), dimension(:), intent(in) buffer )
private

Expand curve data from a communication buffer.

◆ curve_comm_pack()

pure subroutine fitpack_curves::curve_comm_pack ( class(fitpack_curve), intent(in) this,
real(fp_comm), dimension(:), intent(out) buffer )
private

Pack curve data into a communication buffer for parallel transfer.

◆ curve_comm_size()

elemental integer(fp_size) function fitpack_curves::curve_comm_size ( class(fitpack_curve), intent(in) this)
private

Return the communication buffer size (number of FP_COMM elements).

Here is the call graph for this function:

◆ curve_derivative()

real(fp_real) function fitpack_curves::curve_derivative ( class(fitpack_curve), intent(inout) this,
real(fp_real), intent(in) x,
integer, intent(in) order,
integer(fp_flag), intent(out), optional ierr )
private

Evaluate the k-th derivative at a single point (with optional error flag).

Parameters
[in,out]thisThe fitted curve object.
[in]xScalar evaluation point.
[in]orderDerivative order ( \( 0 \leq \text{order} \leq k \)).
[out]ierrOptional error flag.
Returns
Scalar derivative value \( s^{(\text{order})}(x) \).
See also
splder
Here is the call graph for this function:

◆ curve_derivatives()

real(fp_real) function, dimension(size(x)) fitpack_curves::curve_derivatives ( class(fitpack_curve), intent(inout) this,
real(fp_real), dimension(:), intent(in) x,
integer, intent(in) order,
integer(fp_flag), intent(out), optional ierr )
private

Evaluate the k-th derivative at multiple points (with error flag).

Parameters
[in,out]thisThe fitted curve object.
[in]xArray of evaluation points.
[in]orderDerivative order ( \( 0 \leq \text{order} \leq k \)).
[out]ierrOptional error flag.
Returns
Array of derivative values \( s^{(\text{order})}(x_i) \).
See also
splder

◆ curve_eval_many()

real(fp_real) function, dimension(size(x)) fitpack_curves::curve_eval_many ( class(fitpack_curve), intent(inout) this,
real(fp_real), dimension(:), intent(in) x,
integer(fp_flag), intent(out) ierr )
private

Evaluate the spline at multiple points (with error flag).

◆ curve_eval_many_pure()

pure real(fp_real) function, dimension(size(x)) fitpack_curves::curve_eval_many_pure ( class(fitpack_curve), intent(in) this,
real(fp_real), dimension(:), intent(in) x )
private

Evaluate the spline at multiple points (pure, returns NaN on error).

◆ curve_eval_one()

real(fp_real) function fitpack_curves::curve_eval_one ( class(fitpack_curve), intent(inout) this,
real(fp_real), intent(in) x,
integer(fp_flag), intent(out) ierr )
private

Evaluate the spline at a single point (with error flag).

Parameters
[in,out]thisThe fitted curve object.
[in]xEvaluation point.
[out]ierrError flag.
Returns
Spline value \( s(x) \).
See also
splev
Here is the call graph for this function:

◆ curve_eval_one_noerr()

pure real(fp_real) function fitpack_curves::curve_eval_one_noerr ( class(fitpack_curve), intent(in) this,
real(fp_real), intent(in) x )
private

Evaluate the spline at a single point (pure, returns NaN on error).

Parameters
[in]thisThe fitted curve object.
[in]xEvaluation point.
Returns
Spline value \( s(x) \), or NaN on error.
See also
splev
Here is the call graph for this function:

◆ curve_fit_automatic_knots()

integer(fp_flag) function fitpack_curves::curve_fit_automatic_knots ( class(fitpack_curve), intent(inout) this,
real(fp_real), intent(in), optional smoothing,
integer(fp_size), intent(in), optional order,
logical, intent(in), optional keep_knots )
private

Fit a smoothing spline with automatic knot placement (curfit / percur).

Get smoothing trajectory

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

Set/update order

◆ curve_fit_least_squares()

integer(fp_flag) function fitpack_curves::curve_fit_least_squares ( class(fitpack_curve), intent(inout) this,
real(fp_real), intent(in), optional smoothing,
logical, intent(in), optional reset_knots )
private

Fit a least-squares spline with fixed knots.

◆ curve_insert_knot_many()

subroutine fitpack_curves::curve_insert_knot_many ( class(fitpack_curve), intent(inout) this,
real(fp_real), dimension(:), intent(in) x,
integer(fp_flag), intent(out), optional ierr )
private

Insert multiple knots into the B-spline representation.

Parameters
[in,out]thisThe fitted curve.
[in]xArray of knot positions to insert.
[out]ierrOptional error flag.
See also
insert_inplace

◆ curve_insert_knot_one()

subroutine fitpack_curves::curve_insert_knot_one ( class(fitpack_curve), intent(inout) this,
real(fp_real), intent(in) x,
integer(fp_flag), intent(out), optional ierr )
private

Insert a single knot into the B-spline representation.

The spline value is unchanged; only the knot vector and coefficients are updated. Workspace state is reset after insertion.

Parameters
[in,out]thisThe fitted curve.
[in]xKnot position to insert.
[out]ierrOptional error flag.
See also
insert_inplace
Here is the call graph for this function:

◆ destroy()

elemental subroutine fitpack_curves::destroy ( class(fitpack_curve), intent(inout) this)
private

Destroy a curve object and release all allocated memory.

◆ fourier_coefficients()

subroutine fitpack_curves::fourier_coefficients ( class(fitpack_curve), intent(inout) this,
real(fp_real), dimension(:), intent(in) alpha,
real(fp_real), dimension(size(alpha)), intent(out) a,
real(fp_real), dimension(size(alpha)), intent(out) b,
integer(fp_flag), intent(out), optional ierr )
private

Compute Fourier sine and cosine coefficients of the spline.

Computes integrals over the interior knot range \( [t_4, t_{n-3}] \):

\[ A_i = \int_{t_4}^{t_{n-3}} s(x) \sin(\alpha_i x) \, dx, \quad B_i = \int_{t_4}^{t_{n-3}} s(x) \cos(\alpha_i x) \, dx \]

Requires at least 10 knots.

Parameters
[in,out]thisThe fitted curve (must have \( \geq 10 \) knots).
[in]alphaArray of frequency values \( \alpha_i \).
[out]ASine coefficients, same size as alpha.
[out]BCosine coefficients, same size as alpha.
[out]ierrOptional error flag.
See also
fourco

◆ grow_knot_storage()

subroutine fitpack_curves::grow_knot_storage ( class(fitpack_curve), intent(inout) this)
private

Grow the knot and coefficient arrays when storage is insufficient.

◆ integral()

real(fp_real) function fitpack_curves::integral ( class(fitpack_curve), intent(inout) this,
real(fp_real), intent(in) from,
real(fp_real), intent(in) to )
private

Compute the definite integral \( \int_a^b s(x)\,dx \) via splint.

Here is the call graph for this function:

◆ interpolating_curve()

integer(fp_flag) function fitpack_curves::interpolating_curve ( class(fitpack_curve), intent(inout) this,
integer(fp_size), intent(in), optional order,
logical, intent(in), optional reset_knots )
private

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

Here is the call graph for this function:

◆ new_fit()

integer(fp_flag) function fitpack_curves::new_fit ( class(fitpack_curve), 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), optional w,
real(fp_real), intent(in), optional smoothing,
integer(fp_size), intent(in), optional order )
private

Load new data points and perform a fresh fit.

Here is the call graph for this function:

◆ new_from_points()

type(fitpack_curve) function fitpack_curves::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), optional w,
integer(fp_flag), intent(out), optional ierr )
private

Construct a fitpack_curve from data points and perform a default fit.

Parameters
[in]xIndependent variable values (must be sortable).
[in]yDependent variable values, same size as x.
[in]wOptional weights (positive), same size as x.
[out]ierrOptional error flag.
Returns
Fitted curve object.

◆ new_points()

subroutine fitpack_curves::new_points ( class(fitpack_curve), 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), optional w )
private

Load new data points into the curve, resetting all internal state.

Sorts data by x, allocates knot/coefficient storage, and prepares workspaces. A subsequent call to fit, interpolate, or least_squares is required to compute the spline representation.

Parameters
[in,out]thisThe curve object (destroyed and reinitialized).
[in]xIndependent variable values (unsorted OK; will be sorted internally).
[in]yDependent variable values, same size as x.
[in]wOptional positive weights, same size as x. Default: uniform.
See also
curfit, percur

◆ zeros()

real(fp_real) function, dimension(:), allocatable fitpack_curves::zeros ( class(fitpack_curve), intent(in) this,
integer(fp_flag), intent(out), optional ierr )
private

Find all real zeros of a cubic spline \( s(x) = 0 \).

Only valid for cubic splines ( \( k = 3 \)). Iteratively grows storage until all zeros are found.

Parameters
[in]thisThe fitted cubic curve.
[out]ierrOptional error flag.
Returns
Allocatable array of zero locations (empty if none found).
See also
sproot

Always return an allocated array

Here is the call graph for this function:

Variable Documentation

◆ max_k

integer, parameter fitpack_curves::max_k = 5
private