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

Univariate spline curve \( y = s(x) \). More...

Inheritance diagram for fitpack_curves::fitpack_curve:
Collaboration diagram for fitpack_curves::fitpack_curve:

Public Member Functions

procedure destroy (this)
 Clean memory.
 
procedure new_points (this, x, y, w)
 Set new points.
 
procedure new_fit (this, x, y, w, smoothing, order)
 Generate new fit.
 
procedure fit (this, smoothing, order, keep_knots)
 Generate/update fitting curve, with optional smoothing.
 
procedure interpolate (this, order, reset_knots)
 Fit an interpolating spline ( \( s = 0 \)) through the data points.
 
procedure least_squares (this, smoothing, reset_knots)
 Fit a least-squares spline with fixed knots.
 
generic eval (this, x, ierr)
 Evaluate curve at given coordinates.
 
generic eval (this, x)
 Evaluate the spline at a single point (pure, returns NaN on error).
 
generic eval (this, x, ierr)
 Evaluate the spline at multiple points (with error flag).
 
generic eval (this, x)
 Evaluate the spline at multiple points (pure, returns NaN on error).
 
procedure integral (this, from, to)
 Integrate.
 
procedure fourier_coefficients (this, alpha, a, b, ierr)
 Fourier coefficients.
 
procedure zeros (this, ierr)
 Find the zeros of a spline s(x), only if it is cubic.
 
generic dfdx (this, x, order, ierr)
 Evaluate derivative at given coordinates.
 
generic dfdx (this, x, order, ierr)
 Evaluate the k-th derivative at multiple points (with error flag).
 
generic dfdx_all (this, x, ierr)
 Evaluate all derivatives \( s^{(j)}(x) \) for \( j = 0, \ldots, k \) (with error flag).
 
generic dfdx_all (this, x)
 Evaluate all derivatives \( s^{(j)}(x) \) for \( j = 0, \ldots, k \) (pure, NaN on error).
 
generic insert_knot (this, x, ierr)
 Insert knot(s) into the spline representation.
 
generic insert_knot (this, x, ierr)
 Insert multiple knots into the B-spline representation.
 
procedure comm_size (this)
 Parallel communication interface (size/pack/expand)
 
procedure comm_pack (this, buffer)
 Pack curve data into a communication buffer for parallel transfer.
 
procedure comm_expand (this, buffer)
 Expand curve data from a communication buffer.
 
type(fitpack_curve) function new_from_points (x, y, w, ierr)
 Construct a fitpack_curve from 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(fp_sizem = 0
 The data points.
 
real(fp_real), dimension(:), allocatable x
 
real(fp_real), dimension(:), allocatable y
 
integer(fp_sizeorder = 3
 Spline degree.
 
real(fp_realxleft
 Interval boundaries.
 
real(fp_realxright
 
real(fp_real), dimension(:), allocatable sp
 
real(fp_real), dimension(:), allocatable w
 
integer(fp_sizenest = 0
 
real(fp_real), dimension(:,:), allocatable wrk_fou
 
integer(fp_flagbc = OUTSIDE_NEAREST_BND
 
integer(fp_sizeknots = 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 curve_eval_one (this, x, ierr)
 Evaluate curve at given coordinates.
 
procedure, private curve_eval_one_noerr (this, x)
 Evaluate the spline at a single point (pure, returns NaN on error).
 
procedure, private curve_eval_many (this, x, ierr)
 Evaluate the spline at multiple points (with error flag).
 
procedure, private curve_eval_many_pure (this, x)
 Evaluate the spline at multiple points (pure, returns NaN on error).
 
procedure, private curve_derivative (this, x, order, ierr)
 Evaluate derivative at given coordinates.
 
procedure, private curve_derivatives (this, x, order, ierr)
 Evaluate the k-th derivative at multiple points (with error flag).
 
procedure, private curve_all_derivatives (this, x, ierr)
 Evaluate all derivatives \( s^{(j)}(x) \) for \( j = 0, \ldots, k \) (with error flag).
 
procedure, private curve_all_derivatives_pure (this, x)
 Evaluate all derivatives \( s^{(j)}(x) \) for \( j = 0, \ldots, k \) (pure, NaN on error).
 
procedure, private curve_insert_knot_one (this, x, ierr)
 Insert knot(s) into the spline representation.
 
procedure, private curve_insert_knot_many (this, x, ierr)
 Insert multiple knots into the B-spline representation.
 

Detailed Description

Univariate spline curve \( y = s(x) \).

Default constructor: create a curve from data points and fit it.

Fits a smoothing or interpolating spline of degree \( k \) (default 3 = cubic) to weighted data points \( (x_i, y_i, w_i) \). After fitting, the curve can be evaluated, differentiated, integrated, and its zeros found.

Typical usage:

type(fitpack_curve) :: curve
curve = fitpack_curve(x, y) ! construct and fit
y_eval = curve%eval(x_new) ! evaluate at new points
dy = curve%dfdx(x_new, order=1) ! first derivative
area = curve%integral(a, b) ! definite integral

Member Function/Subroutine Documentation

◆ comm_expand()

procedure fitpack_curves::fitpack_curve::comm_expand ( class(fitpack_curve), intent(inout) this,
real(fp_comm), dimension(:), intent(in) buffer )

Expand curve data from a communication buffer.

◆ comm_pack()

procedure fitpack_curves::fitpack_curve::comm_pack ( class(fitpack_curve), intent(in) this,
real(fp_comm), dimension(:), intent(out) buffer )

Pack curve data into a communication buffer for parallel transfer.

◆ comm_size()

procedure fitpack_curves::fitpack_curve::comm_size ( class(fitpack_curve), intent(in) this)

Parallel communication interface (size/pack/expand)

◆ curve_all_derivatives()

procedure, private fitpack_curves::fitpack_curve::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()

procedure, private fitpack_curves::fitpack_curve::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_derivative()

procedure, private fitpack_curves::fitpack_curve::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 derivative at given coordinates.

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

◆ curve_derivatives()

procedure, private fitpack_curves::fitpack_curve::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()

procedure, private fitpack_curves::fitpack_curve::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()

procedure, private fitpack_curves::fitpack_curve::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()

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

Evaluate curve at given coordinates.

Parameters
[in,out]thisThe fitted curve object.
[in]xEvaluation point.
[out]ierrError flag.
Returns
Spline value \( s(x) \).
See also
splev

◆ curve_eval_one_noerr()

procedure, private fitpack_curves::fitpack_curve::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

◆ curve_insert_knot_many()

procedure, private fitpack_curves::fitpack_curve::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()

procedure, private fitpack_curves::fitpack_curve::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 knot(s) into the 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

◆ destroy()

procedure fitpack_curves::fitpack_curve::destroy ( class(fitpack_curve), intent(inout) this)

Clean memory.

◆ dfdx() [1/2]

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

Evaluate derivative at given coordinates.

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:

◆ dfdx() [2/2]

generic fitpack_curves::fitpack_curve::dfdx ( class(fitpack_curve), intent(inout) this,
real(fp_real), dimension(:), intent(in) x,
integer, intent(in) order,
integer(fp_flag), intent(out), optional ierr )

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

◆ dfdx_all() [1/2]

generic fitpack_curves::fitpack_curve::dfdx_all ( class(fitpack_curve), intent(in) this,
real(fp_real), intent(in) x )

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

◆ dfdx_all() [2/2]

generic fitpack_curves::fitpack_curve::dfdx_all ( class(fitpack_curve), intent(inout) this,
real(fp_real), intent(in) x,
integer(fp_flag), intent(out) ierr )

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
Here is the call graph for this function:

◆ eval() [1/4]

generic fitpack_curves::fitpack_curve::eval ( class(fitpack_curve), intent(in) this,
real(fp_real), intent(in) x )

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

◆ eval() [2/4]

generic fitpack_curves::fitpack_curve::eval ( class(fitpack_curve), intent(in) this,
real(fp_real), dimension(:), intent(in) x )

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

◆ eval() [3/4]

generic fitpack_curves::fitpack_curve::eval ( class(fitpack_curve), intent(inout) this,
real(fp_real), intent(in) x,
integer(fp_flag), intent(out) ierr )

Evaluate curve at given coordinates.

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:

◆ eval() [4/4]

generic fitpack_curves::fitpack_curve::eval ( class(fitpack_curve), intent(inout) this,
real(fp_real), dimension(:), intent(in) x,
integer(fp_flag), intent(out) ierr )

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

◆ fit()

procedure fitpack_curves::fitpack_curve::fit ( 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 )

Generate/update fitting curve, with optional smoothing.

◆ fourier_coefficients()

procedure fitpack_curves::fitpack_curve::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 )

Fourier coefficients.

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

◆ insert_knot() [1/2]

generic fitpack_curves::fitpack_curve::insert_knot ( class(fitpack_curve), intent(inout) this,
real(fp_real), intent(in) x,
integer(fp_flag), intent(out), optional ierr )

Insert knot(s) into the 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:

◆ insert_knot() [2/2]

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

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

◆ integral()

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

Integrate.

◆ interpolate()

procedure fitpack_curves::fitpack_curve::interpolate ( class(fitpack_curve), intent(inout) this,
integer(fp_size), intent(in), optional order,
logical, intent(in), optional reset_knots )

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

◆ least_squares()

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

Fit a least-squares spline with fixed knots.

◆ new_fit()

procedure fitpack_curves::fitpack_curve::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 )

Generate new fit.

◆ new_from_points()

type(fitpack_curve) function fitpack_curves::fitpack_curve::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 )

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.
Here is the call graph for this function:

◆ new_points()

procedure fitpack_curves::fitpack_curve::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 )

Set new points.

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

procedure fitpack_curves::fitpack_curve::zeros ( class(fitpack_curve), intent(in) this,
integer(fp_flag), intent(out), optional ierr )

Find the zeros of a spline s(x), only if it is cubic.

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

Member Data Documentation

◆ bc

integer(fp_flag) fitpack_curves::fitpack_curve::bc = OUTSIDE_NEAREST_BND

◆ knots

integer(fp_size) fitpack_curves::fitpack_curve::knots = 0

◆ m

integer(fp_size) fitpack_curves::fitpack_curve::m = 0

The data points.

◆ nest

integer(fp_size) fitpack_curves::fitpack_curve::nest = 0

◆ order

integer(fp_size) fitpack_curves::fitpack_curve::order = 3

Spline degree.

◆ sp

real(fp_real), dimension(:), allocatable fitpack_curves::fitpack_curve::sp

◆ t

real(fp_real), dimension(:), allocatable fitpack_curves::fitpack_curve::t

◆ w

real(fp_real), dimension(:), allocatable fitpack_curves::fitpack_curve::w

◆ wrk_fou

real(fp_real), dimension(:,:), allocatable fitpack_curves::fitpack_curve::wrk_fou

◆ x

real(fp_real), dimension(:), allocatable fitpack_curves::fitpack_curve::x

◆ xleft

real(fp_real) fitpack_curves::fitpack_curve::xleft

Interval boundaries.

◆ xright

real(fp_real) fitpack_curves::fitpack_curve::xright

◆ y

real(fp_real), dimension(:), allocatable fitpack_curves::fitpack_curve::y

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