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

OOP wrappers for parametric curve fitting in \( \mathbb{R}^d \). More...

Data Types

type  fitpack_closed_curve
 Closed (periodic) parametric curve fitter. More...
 
type  fitpack_constrained_curve
 Parametric curve fitter with derivative constraints at the endpoints. More...
 
interface  fitpack_parametric_curve
 Open parametric curve fitter in \( \mathbb{R}^d \). More...
 

Functions/Subroutines

type(fitpack_parametric_curve) function new_from_points (x, u, w, ierr)
 Construct a parametric curve from data points and perform a default fit.
 
integer function new_fit (this, x, u, w, smoothing, order)
 Load new data points and perform a fresh fit.
 
elemental subroutine destroy (this)
 Destroy a parametric curve object and release all allocated memory.
 
elemental subroutine con_destroy (this)
 Destroy a constrained parametric curve, including endpoint constraints.
 
subroutine new_points (this, x, u, w)
 Load new data points into the parametric curve and allocate workspace.
 
elemental subroutine clean_constraints (this)
 Remove all endpoint derivative constraints.
 
subroutine set_constraints (this, ddx_begin, ddx_end, ierr)
 Set or reset endpoint derivative constraints for a constrained parametric curve.
 
real(fp_real) function, dimension(this%idim) curve_eval_one (this, u, ierr)
 Evaluate the parametric curve at a single parameter value.
 
real(fp_real) function, dimension(this%idim, size(u)) curve_eval_many (this, u, ierr)
 Evaluate the parametric curve at multiple parameter values.
 
integer function interpolating_curve (this, order, reset_knots)
 Fit an interpolating parametric spline ( \( s = 0 \)) through the data points.
 
integer function parcur_fit_least_squares (this, smoothing, reset_knots)
 Fit a least-squares parametric spline with fixed knots.
 
integer function curve_fit_automatic_knots (this, smoothing, order, keep_knots)
 Fit a smoothing parametric spline with automatic knot placement.
 
real(fp_real) function, dimension(this%idim) curve_derivative (this, u, order, ierr)
 Evaluate the k-th derivative of the parametric curve at a single parameter value.
 
real(fp_real) function, dimension(this%idim, 0:this%order) curve_all_derivatives (this, u, ierr)
 Evaluate all derivatives \( \mathbf{s}^{(j)}(u) \) for \( j = 0, \ldots, k \).
 
subroutine set_default_parameters (this)
 Compute default parameter values from cumulative chord lengths.
 
real(fp_real) function, dimension(this%idim, size(u)) curve_derivatives (this, u, order, ierr)
 Evaluate the k-th derivative at multiple parameter values.
 
elemental integer(fp_size) function parcur_comm_size (this)
 Return the communication buffer size for the parametric curve.
 
pure subroutine parcur_comm_pack (this, buffer)
 Pack parametric curve data into a communication buffer.
 
pure subroutine parcur_comm_expand (this, buffer)
 Expand parametric curve data from a communication buffer.
 
elemental integer(fp_size) function concur_comm_size (this)
 Return the communication buffer size for the constrained curve.
 
pure subroutine concur_comm_pack (this, buffer)
 Pack constrained curve data into a communication buffer.
 
pure subroutine concur_comm_expand (this, buffer)
 Expand constrained curve data from a communication buffer.
 

Variables

integer, parameter max_k = 5
 

Detailed Description

OOP wrappers for parametric curve fitting in \( \mathbb{R}^d \).

Provides three types for fitting parametric curves through data points in \( d \)-dimensional space ( \( d \geq 1 \)):

Each component \( s_j(u) \) of the curve is represented as a B-spline of degree \( k \) over a common knot vector, with the parameter \( u \) either user-supplied or automatically computed from cumulative chord lengths.

See also
Dierckx, Ch. 9 (pp. 199–228); parcur, clocur, concur

Function/Subroutine Documentation

◆ clean_constraints()

elemental subroutine fitpack_parametric_curves::clean_constraints ( class(fitpack_constrained_curve), intent(inout) this)
private

Remove all endpoint derivative constraints.

◆ con_destroy()

elemental subroutine fitpack_parametric_curves::con_destroy ( class(fitpack_constrained_curve), intent(inout) this)
private

Destroy a constrained parametric curve, including endpoint constraints.

Here is the call graph for this function:

◆ concur_comm_expand()

pure subroutine fitpack_parametric_curves::concur_comm_expand ( class(fitpack_constrained_curve), intent(inout) this,
real(fp_comm), dimension(:), intent(in) buffer )
private

Expand constrained curve data from a communication buffer.

Here is the call graph for this function:

◆ concur_comm_pack()

pure subroutine fitpack_parametric_curves::concur_comm_pack ( class(fitpack_constrained_curve), intent(in) this,
real(fp_comm), dimension(:), intent(out) buffer )
private

Pack constrained curve data into a communication buffer.

Here is the call graph for this function:

◆ concur_comm_size()

elemental integer(fp_size) function fitpack_parametric_curves::concur_comm_size ( class(fitpack_constrained_curve), intent(in) this)
private

Return the communication buffer size for the constrained curve.

Here is the call graph for this function:

◆ curve_all_derivatives()

real(fp_real) function, dimension(this%idim,0:this%order) fitpack_parametric_curves::curve_all_derivatives ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_real), intent(in) u,
integer, intent(out), optional ierr )
private

Evaluate all derivatives \( \mathbf{s}^{(j)}(u) \) for \( j = 0, \ldots, k \).

Parameters
[in,out]thisThe fitted parametric curve.
[in]uParameter value.
[out]ierrOptional error flag.
Returns
Array \( (d \times (k+1)) \) of derivatives at u.
See also
cualde

◆ curve_derivative()

real(fp_real) function, dimension(this%idim) fitpack_parametric_curves::curve_derivative ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_real), intent(in) u,
integer, intent(in) order,
integer, intent(out), optional ierr )
private

Evaluate the k-th derivative of the parametric curve at a single parameter value.

Parameters
[in,out]thisThe fitted parametric curve.
[in]uParameter value.
[in]orderDerivative order ( \( 0 \leq \text{order} \leq k \)).
[out]ierrOptional error flag.
Returns
Derivative vector \( \mathbf{s}^{(\text{order})}(u) \in \mathbb{R}^d \).
See also
cualde

◆ curve_derivatives()

real(fp_real) function, dimension(this%idim,size(u)) fitpack_parametric_curves::curve_derivatives ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_real), dimension(:), intent(in) u,
integer, intent(in) order,
integer, intent(out), optional ierr )
private

Evaluate the k-th derivative at multiple parameter values.

Parameters
[in,out]thisThe fitted parametric curve.
[in]uArray of parameter values.
[in]orderDerivative order.
[out]ierrOptional error flag.
Returns
Array \( (d \times m) \) of derivative vectors.
See also
cualde
Here is the call graph for this function:

◆ curve_eval_many()

real(fp_real) function, dimension(this%idim,size(u)) fitpack_parametric_curves::curve_eval_many ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_real), dimension(:), intent(in) u,
integer, intent(out), optional ierr )
private

Evaluate the parametric curve at multiple parameter values.

Parameters
[in,out]thisThe fitted parametric curve.
[in]uArray of parameter values.
[out]ierrOptional error flag.
Returns
Points \( \mathbf{s}(u_i) \) as columns of a \( d \times m \) array.
See also
curev

◆ curve_eval_one()

real(fp_real) function, dimension(this%idim) fitpack_parametric_curves::curve_eval_one ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_real), intent(in) u,
integer, intent(out), optional ierr )
private

Evaluate the parametric curve at a single parameter value.

Parameters
[in,out]thisThe fitted parametric curve.
[in]uParameter value.
[out]ierrOptional error flag.
Returns
Point \( \mathbf{s}(u) \in \mathbb{R}^d \).
See also
curev
Here is the call graph for this function:

◆ curve_fit_automatic_knots()

integer function fitpack_parametric_curves::curve_fit_automatic_knots ( class(fitpack_parametric_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 parametric spline with automatic knot placement.

Dispatches to parcur (open), clocur (closed), or concur (constrained) depending on the dynamic type.

See also
parcur, clocur, concur

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

◆ destroy()

elemental subroutine fitpack_parametric_curves::destroy ( class(fitpack_parametric_curve), intent(inout) this)
private

Destroy a parametric curve object and release all allocated memory.

◆ interpolating_curve()

integer function fitpack_parametric_curves::interpolating_curve ( class(fitpack_parametric_curve), intent(inout) this,
integer(fp_size), intent(in), optional order,
logical, intent(in), optional reset_knots )
private

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

See also
parcur, clocur, concur
Here is the call graph for this function:

◆ new_fit()

integer function fitpack_parametric_curves::new_fit ( class(fitpack_parametric_curve), intent(inout) this,
real (fp_real), dimension(:,:), intent(in) x,
real (fp_real), dimension(size(x,2)), intent(in), optional u,
real (fp_real), dimension(size(x,2)), 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.

See also
parcur, clocur, concur
Here is the call graph for this function:

◆ new_from_points()

type(fitpack_parametric_curve) function fitpack_parametric_curves::new_from_points ( real(fp_real), dimension(:,:), intent(in) x,
real(fp_real), dimension(size(x,2)), intent(in), optional u,
real(fp_real), dimension(size(x,2)), intent(in), optional w,
integer, intent(out), optional ierr )
private

Construct a parametric curve from data points and perform a default fit.

Parameters
[in]xData points \( (d \times m) \): each column is a point in \( \mathbb{R}^d \).
[in]uOptional parameter values (must be strictly increasing). If omitted, computed from cumulative chord lengths.
[in]wOptional positive weights, one per data point.
[out]ierrOptional error flag.
Returns
Fitted parametric curve object.
See also
parcur

◆ new_points()

subroutine fitpack_parametric_curves::new_points ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_real), dimension(:,:), intent(in) x,
real(fp_real), dimension(size(x,2)), intent(in), optional u,
real(fp_real), dimension(size(x,2)), intent(in), optional w )
private

Load new data points into the parametric curve and allocate workspace.

Parameters
[in,out]thisThe parametric curve (destroyed and reinitialized).
[in]xData points \( (d \times m) \): each column is a point.
[in]uOptional parameter values. If omitted, computed from chord lengths.
[in]wOptional positive weights.
See also
parcur, clocur, concur
Here is the call graph for this function:

◆ parcur_comm_expand()

pure subroutine fitpack_parametric_curves::parcur_comm_expand ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_comm), dimension(:), intent(in) buffer )
private

Expand parametric curve data from a communication buffer.

◆ parcur_comm_pack()

pure subroutine fitpack_parametric_curves::parcur_comm_pack ( class(fitpack_parametric_curve), intent(in) this,
real(fp_comm), dimension(:), intent(out) buffer )
private

Pack parametric curve data into a communication buffer.

◆ parcur_comm_size()

elemental integer(fp_size) function fitpack_parametric_curves::parcur_comm_size ( class(fitpack_parametric_curve), intent(in) this)
private

Return the communication buffer size for the parametric curve.

Here is the call graph for this function:

◆ parcur_fit_least_squares()

integer function fitpack_parametric_curves::parcur_fit_least_squares ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_real), intent(in), optional smoothing,
logical, intent(in), optional reset_knots )
private

Fit a least-squares parametric spline with fixed knots.

See also
parcur, clocur, concur

◆ set_constraints()

subroutine fitpack_parametric_curves::set_constraints ( class(fitpack_constrained_curve), intent(inout) this,
real(fp_real), dimension(:,0:), intent(in), optional ddx_begin,
real(fp_real), dimension (:,0:), intent(in), optional ddx_end,
integer, intent(out), optional ierr )
private

Set or reset endpoint derivative constraints for a constrained parametric curve.

Resets all constraints, then applies the supplied left and/or right boundary conditions. Omitting ddx_begin or ddx_end removes constraints at that endpoint.

Parameters
[in,out]thisThe constrained curve object.
[in]ddx_beginLeft endpoint constraints: column 0 = function value, column \( \ell \) = \( \ell \)-th derivative ( \( d \times (ib) \)).
[in]ddx_endRight endpoint constraints (same layout as ddx_begin).
[out]ierrError flag (optional).
Here is the call graph for this function:

◆ set_default_parameters()

subroutine fitpack_parametric_curves::set_default_parameters ( class(fitpack_parametric_curve), intent(inout) this)
private

Compute default parameter values from cumulative chord lengths.

Sets \( u_i \) to normalized cumulative Euclidean distances between consecutive data points, so that \( u_1 = 0 \) and \( u_m = 1 \).

Variable Documentation

◆ max_k

integer, parameter fitpack_parametric_curves::max_k = 5
private