|
fitpack
Modern Fortran library for curve and surface fitting with splines
|
Univariate spline curve \( y = s(x) \). More...


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_size) | m = 0 |
| The data points. | |
| real(fp_real), dimension(:), allocatable | x |
| real(fp_real), dimension(:), allocatable | y |
| integer(fp_size) | order = 3 |
| Spline degree. | |
| real(fp_real) | xleft |
| Interval boundaries. | |
| real(fp_real) | xright |
| real(fp_real), dimension(:), allocatable | sp |
| real(fp_real), dimension(:), allocatable | w |
| integer(fp_size) | nest = 0 |
| real(fp_real), dimension(:,:), allocatable | wrk_fou |
| integer(fp_flag) | bc = OUTSIDE_NEAREST_BND |
| integer(fp_size) | knots = 0 |
| real(fp_real), dimension(:), allocatable | t |
Public Attributes inherited from fitpack_fitters::fitpack_fitter | |
| integer(fp_flag) | iopt = IOPT_NEW_SMOOTHING |
| Fitting state flag. | |
| real(fp_real) | smoothing = 1000.0_FP_REAL |
| Smoothing parameter. | |
| real(fp_real) | fp = zero |
| Weighted sum of squared residuals. | |
| real(fp_real), dimension(:), allocatable | c |
| B-spline coefficients. | |
| integer(fp_size) | lwrk = 0 |
| Real workspace and its size. | |
| real(fp_real), dimension(:), allocatable | wrk |
| integer(fp_size) | liwrk = 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. | |
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:
| 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.
| 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.
| procedure fitpack_curves::fitpack_curve::comm_size | ( | class(fitpack_curve), intent(in) | this | ) |
Parallel communication interface (size/pack/expand)
|
private |
Evaluate all derivatives \( s^{(j)}(x) \) for \( j = 0, \ldots, k \) (with error flag).
| [in,out] | this | The fitted curve object. |
| [in] | x | Scalar evaluation point. |
| [out] | ierr | Error flag. |
|
private |
Evaluate all derivatives \( s^{(j)}(x) \) for \( j = 0, \ldots, k \) (pure, NaN on error).
| [in] | this | The fitted curve object. |
| [in] | x | Scalar evaluation point. |
|
private |
Evaluate derivative at given coordinates.
| [in,out] | this | The fitted curve object. |
| [in] | x | Scalar evaluation point. |
| [in] | order | Derivative order ( \( 0 \leq \text{order} \leq k \)). |
| [out] | ierr | Optional error flag. |
|
private |
Evaluate the k-th derivative at multiple points (with error flag).
| [in,out] | this | The fitted curve object. |
| [in] | x | Array of evaluation points. |
| [in] | order | Derivative order ( \( 0 \leq \text{order} \leq k \)). |
| [out] | ierr | Optional error flag. |
|
private |
Evaluate the spline at multiple points (with error flag).
|
private |
Evaluate the spline at multiple points (pure, returns NaN on error).
|
private |
Evaluate curve at given coordinates.
| [in,out] | this | The fitted curve object. |
| [in] | x | Evaluation point. |
| [out] | ierr | Error flag. |
|
private |
Evaluate the spline at a single point (pure, returns NaN on error).
| [in] | this | The fitted curve object. |
| [in] | x | Evaluation point. |
|
private |
Insert multiple knots into the B-spline representation.
| [in,out] | this | The fitted curve. |
| [in] | x | Array of knot positions to insert. |
| [out] | ierr | Optional error flag. |
|
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.
| [in,out] | this | The fitted curve. |
| [in] | x | Knot position to insert. |
| [out] | ierr | Optional error flag. |
| procedure fitpack_curves::fitpack_curve::destroy | ( | class(fitpack_curve), intent(inout) | this | ) |
Clean memory.
| 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.
| [in,out] | this | The fitted curve object. |
| [in] | x | Scalar evaluation point. |
| [in] | order | Derivative order ( \( 0 \leq \text{order} \leq k \)). |
| [out] | ierr | Optional error flag. |

| 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).
| [in,out] | this | The fitted curve object. |
| [in] | x | Array of evaluation points. |
| [in] | order | Derivative order ( \( 0 \leq \text{order} \leq k \)). |
| [out] | ierr | Optional error flag. |
| 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).
| [in] | this | The fitted curve object. |
| [in] | x | Scalar evaluation point. |
| 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).
| [in,out] | this | The fitted curve object. |
| [in] | x | Scalar evaluation point. |
| [out] | ierr | Error flag. |

| 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).
| [in] | this | The fitted curve object. |
| [in] | x | Evaluation point. |
| 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).
| 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.
| [in,out] | this | The fitted curve object. |
| [in] | x | Evaluation point. |
| [out] | ierr | Error flag. |

| 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).
| 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.
| 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.
| [in,out] | this | The fitted curve (must have \( \geq 10 \) knots). |
| [in] | alpha | Array of frequency values \( \alpha_i \). |
| [out] | A | Sine coefficients, same size as alpha. |
| [out] | B | Cosine coefficients, same size as alpha. |
| [out] | ierr | Optional error flag. |
| 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.
| [in,out] | this | The fitted curve. |
| [in] | x | Knot position to insert. |
| [out] | ierr | Optional error flag. |

| 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.
| [in,out] | this | The fitted curve. |
| [in] | x | Array of knot positions to insert. |
| [out] | ierr | Optional error flag. |
| 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.
| 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.
| 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.
| 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.
| 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.
| [in] | x | Independent variable values (must be sortable). |
| [in] | y | Dependent variable values, same size as x. |
| [in] | w | Optional weights (positive), same size as x. |
| [out] | ierr | Optional error flag. |

| 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.
| [in,out] | this | The curve object (destroyed and reinitialized). |
| [in] | x | Independent variable values (unsorted OK; will be sorted internally). |
| [in] | y | Dependent variable values, same size as x. |
| [in] | w | Optional positive weights, same size as x. Default: uniform. |
| 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.
| [in] | this | The fitted cubic curve. |
| [out] | ierr | Optional error flag. |
| integer(fp_flag) fitpack_curves::fitpack_curve::bc = OUTSIDE_NEAREST_BND |
| integer(fp_size) fitpack_curves::fitpack_curve::knots = 0 |
| integer(fp_size) fitpack_curves::fitpack_curve::m = 0 |
The data points.
| integer(fp_size) fitpack_curves::fitpack_curve::nest = 0 |
| integer(fp_size) fitpack_curves::fitpack_curve::order = 3 |
Spline degree.
| real(fp_real), dimension(:), allocatable fitpack_curves::fitpack_curve::sp |
| real(fp_real), dimension(:), allocatable fitpack_curves::fitpack_curve::t |
| real(fp_real), dimension(:), allocatable fitpack_curves::fitpack_curve::w |
| real(fp_real), dimension(:,:), allocatable fitpack_curves::fitpack_curve::wrk_fou |
| real(fp_real), dimension(:), allocatable fitpack_curves::fitpack_curve::x |
| real(fp_real) fitpack_curves::fitpack_curve::xleft |
Interval boundaries.
| real(fp_real) fitpack_curves::fitpack_curve::xright |
| real(fp_real), dimension(:), allocatable fitpack_curves::fitpack_curve::y |