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

Open parametric curve fitter in \( \mathbb{R}^d \). More...

Inheritance diagram for fitpack_parametric_curves::fitpack_parametric_curve:
Collaboration diagram for fitpack_parametric_curves::fitpack_parametric_curve:

Public Member Functions

procedure destroy (this)
 Clean memory.
 
procedure new_points (this, x, u, w)
 Set new points.
 
procedure set_default_parameters (this)
 Compute default parameter values from cumulative chord lengths.
 
procedure new_fit (this, x, u, 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 parametric spline ( \( s = 0 \)) through the data points.
 
procedure least_squares (this, smoothing, reset_knots)
 Fit a least-squares parametric spline with fixed knots.
 
procedure eval_one (this, u, ierr)
 Evaluate curve at given coordinates.
 
procedure eval_many (this, u, ierr)
 Evaluate the parametric curve at multiple parameter values.
 
generic eval (this, u, ierr)
 Evaluate curve at given coordinates.
 
generic eval (this, u, ierr)
 Evaluate the parametric curve at multiple parameter values.
 
generic dfdx (this, u, order, ierr)
 Evaluate derivative at given coordinates.
 
generic dfdx (this, u, order, ierr)
 Evaluate the k-th derivative at multiple parameter values.
 
generic dfdx_all (this, u, ierr)
 Evaluate all derivatives \( \mathbf{s}^{(j)}(u) \) for \( j = 0, \ldots, k \).
 
procedure comm_size (this)
 Parallel communication interface (size/pack/expand)
 
procedure comm_pack (this, buffer)
 Pack parametric curve data into a communication buffer.
 
procedure comm_expand (this, buffer)
 Expand parametric curve data from a communication buffer.
 
type(fitpack_parametric_curve) function new_from_points (x, u, w, ierr)
 Construct a parametric 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 m = 0
 Number of points.
 
integer idim = 0
 Number of dimensions.
 
real(fp_real), dimension(:,:), allocatable x
 The data points.
 
logical has_params = .false.
 Parameter values: one for each point. They may be optional, in which case, they will be internally calculated by fitpack.
 
real(fp_real), dimension(:), allocatable u
 
integer order = 3
 Spline degree.
 
real(fp_realubegin = zero
 Interval boundaries.
 
real(fp_realuend = zero
 
real(fp_real), dimension(:), allocatable sp
 
real(fp_real), dimension(:), allocatable w
 
integer nest = 0
 
real(fp_real), dimension(:,:), allocatable dd
 
integer 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 curve_derivative (this, u, order, ierr)
 Evaluate derivative at given coordinates.
 
procedure, private curve_derivatives (this, u, order, ierr)
 Evaluate the k-th derivative at multiple parameter values.
 
procedure, private curve_all_derivatives (this, u, ierr)
 Evaluate all derivatives \( \mathbf{s}^{(j)}(u) \) for \( j = 0, \ldots, k \).
 

Detailed Description

Open parametric curve fitter in \( \mathbb{R}^d \).

Fits a parametric spline curve \( \mathbf{s}(u) = (s_1(u), \ldots, s_d(u)) \) through \( m \) data points \( \mathbf{x}_i \in \mathbb{R}^d \) associated with strictly increasing parameter values \( u_i \). If no parameter values are provided, they are computed automatically from cumulative chord lengths.

See also
Dierckx, Ch. 9, §9.1 (pp. 199–212); parcur

Member Function/Subroutine Documentation

◆ comm_expand()

procedure fitpack_parametric_curves::fitpack_parametric_curve::comm_expand ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_comm), dimension(:), intent(in) buffer )

Expand parametric curve data from a communication buffer.

◆ comm_pack()

procedure fitpack_parametric_curves::fitpack_parametric_curve::comm_pack ( class(fitpack_parametric_curve), intent(in) this,
real(fp_comm), dimension(:), intent(out) buffer )

Pack parametric curve data into a communication buffer.

◆ comm_size()

procedure fitpack_parametric_curves::fitpack_parametric_curve::comm_size ( class(fitpack_parametric_curve), intent(in) this)

Parallel communication interface (size/pack/expand)

◆ curve_all_derivatives()

procedure, private fitpack_parametric_curves::fitpack_parametric_curve::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()

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

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

procedure, private fitpack_parametric_curves::fitpack_parametric_curve::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

◆ destroy()

procedure fitpack_parametric_curves::fitpack_parametric_curve::destroy ( class(fitpack_parametric_curve), intent(inout) this)

Clean memory.

◆ dfdx() [1/2]

generic fitpack_parametric_curves::fitpack_parametric_curve::dfdx ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_real), intent(in) u,
integer, intent(in) order,
integer, intent(out), optional ierr )

Evaluate derivative at given coordinates.

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

◆ dfdx() [2/2]

generic fitpack_parametric_curves::fitpack_parametric_curve::dfdx ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_real), dimension(:), intent(in) u,
integer, intent(in) order,
integer, intent(out), optional ierr )

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

◆ dfdx_all()

generic fitpack_parametric_curves::fitpack_parametric_curve::dfdx_all ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_real), intent(in) u,
integer, intent(out), optional ierr )

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

◆ eval() [1/2]

generic fitpack_parametric_curves::fitpack_parametric_curve::eval ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_real), intent(in) u,
integer, intent(out), optional ierr )

Evaluate curve at given coordinates.

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

◆ eval() [2/2]

generic fitpack_parametric_curves::fitpack_parametric_curve::eval ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_real), dimension(:), intent(in) u,
integer, intent(out), optional ierr )

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

◆ eval_many()

procedure fitpack_parametric_curves::fitpack_parametric_curve::eval_many ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_real), dimension(:), intent(in) u,
integer, intent(out), optional ierr )

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

◆ eval_one()

procedure fitpack_parametric_curves::fitpack_parametric_curve::eval_one ( class(fitpack_parametric_curve), intent(inout) this,
real(fp_real), intent(in) u,
integer, intent(out), optional ierr )

Evaluate curve at given coordinates.

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

◆ fit()

procedure fitpack_parametric_curves::fitpack_parametric_curve::fit ( 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 )

Generate/update fitting curve, with optional smoothing.

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

See also
parcur, clocur, concur

◆ interpolate()

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

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

See also
parcur, clocur, concur

◆ least_squares()

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

Fit a least-squares parametric spline with fixed knots.

See also
parcur, clocur, concur

◆ new_fit()

procedure fitpack_parametric_curves::fitpack_parametric_curve::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 )

Generate new fit.

See also
parcur, clocur, concur

◆ new_from_points()

type(fitpack_parametric_curve) function fitpack_parametric_curves::fitpack_parametric_curve::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 )

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

◆ new_points()

procedure fitpack_parametric_curves::fitpack_parametric_curve::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 )

Set new points.

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

◆ set_default_parameters()

procedure fitpack_parametric_curves::fitpack_parametric_curve::set_default_parameters ( class(fitpack_parametric_curve), intent(inout) this)

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 \).

Member Data Documentation

◆ dd

real(fp_real), dimension(:,:), allocatable fitpack_parametric_curves::fitpack_parametric_curve::dd

◆ has_params

logical fitpack_parametric_curves::fitpack_parametric_curve::has_params = .false.

Parameter values: one for each point. They may be optional, in which case, they will be internally calculated by fitpack.

◆ idim

integer fitpack_parametric_curves::fitpack_parametric_curve::idim = 0

Number of dimensions.

◆ knots

integer fitpack_parametric_curves::fitpack_parametric_curve::knots = 0

◆ m

integer fitpack_parametric_curves::fitpack_parametric_curve::m = 0

Number of points.

◆ nest

integer fitpack_parametric_curves::fitpack_parametric_curve::nest = 0

◆ order

integer fitpack_parametric_curves::fitpack_parametric_curve::order = 3

Spline degree.

◆ sp

real(fp_real), dimension(:), allocatable fitpack_parametric_curves::fitpack_parametric_curve::sp

◆ t

real(fp_real), dimension(:), allocatable fitpack_parametric_curves::fitpack_parametric_curve::t

◆ u

real(fp_real), dimension(:), allocatable fitpack_parametric_curves::fitpack_parametric_curve::u

◆ ubegin

real(fp_real) fitpack_parametric_curves::fitpack_parametric_curve::ubegin = zero

Interval boundaries.

◆ uend

real(fp_real) fitpack_parametric_curves::fitpack_parametric_curve::uend = zero

◆ w

real(fp_real), dimension(:), allocatable fitpack_parametric_curves::fitpack_parametric_curve::w

◆ x

real(fp_real), dimension(:,:), allocatable fitpack_parametric_curves::fitpack_parametric_curve::x

The data points.


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