|
fitpack
Modern Fortran library for curve and surface fitting with splines
|
Bivariate surface fitter \( z = s(x, y) \) for gridded data. More...


Public Member Functions | |
| procedure | destroy (this) |
| Clean memory. | |
| procedure | new_points (this, x, y, z) |
| Set new points. | |
| procedure | new_fit (this, x, y, z, smoothing, order) |
| Generate new fit. | |
| procedure | fit (this, smoothing, order, keep_knots) |
| Generate/update fitting curve, with optional /Users/federico/code/fitpack/test/fitpack_curve_tests.f90smoothing. | |
| procedure | least_squares (this, smoothing, reset_knots) |
| Fit a least-squares gridded surface with fixed knots. | |
| procedure | interpolate (this, reset_knots) |
| Fit an interpolating gridded surface ( \( s = 0 \)). | |
| generic | eval (this, x, y, ierr) |
| Evaluate at scattered (x,y) points (bispeu). For gridded output, use eval_ongrid. | |
| generic | eval (this, x, y, ierr) |
| Evaluate the grid surface at scattered \( (x_i, y_i) \) points. | |
| generic | eval_ongrid (this, x, y, ierr) |
| Evaluate on a rectangular grid x(:) x y(:) (bispev). Returns f(ny,nx). | |
| generic | eval_ongrid (this, x, y, ierr) |
| Evaluate the grid surface on a rectangular evaluation grid. | |
| generic | dfdx (this, x, y, dx, dy, ierr) |
| Evaluate a partial derivative at a single \( (x, y) \) point. | |
| generic | dfdx (this, x, y, dx, dy, ierr) |
| Evaluate partial derivatives at scattered \( (x_i, y_i) \) points. | |
| generic | dfdx_ongrid (this, x, y, dx, dy, ierr) |
| Evaluate derivatives at given coordinates. | |
| procedure | integral (this, lower, upper) |
| Double integration over a rectangular domain. | |
| procedure | cross_section (this, u, along_y, ierr) |
| Extract a 1D cross-section curve from the surface. | |
| procedure | derivative_spline (this, nux, nuy, ierr) |
| Compute the derivative spline coefficients. | |
| procedure | comm_size (this) |
| Parallel communication. | |
| procedure | comm_pack (this, buffer) |
| Pack grid surface data into a communication buffer. | |
| procedure | comm_expand (this, buffer) |
| Expand grid surface data from a communication buffer. | |
| type(fitpack_grid_surface) function | surf_new_from_points (x, y, z, ierr) |
| Construct a grid surface from gridded data 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 | |
| real(fp_real), dimension(:), allocatable | x |
| The data points. | |
| real(fp_real), dimension(:), allocatable | y |
| real(fp_real), dimension(:,:), allocatable | z |
| integer(fp_size), dimension(2) | order = 3 |
| Spline degree. | |
| real(fp_real), dimension(2) | left |
| Interval boundaries. | |
| real(fp_real), dimension(2) | right |
| integer(fp_size), dimension(2) | nest = 0 |
| Node weights are not allowed. | |
| integer(fp_size) | nmax = 0 |
| integer(fp_size), dimension(2) | 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 | gridsurf_eval_one (this, x, y, ierr) |
| Evaluate at scattered (x,y) points (bispeu). For gridded output, use eval_ongrid. | |
| procedure, private | gridsurf_eval_many (this, x, y, ierr) |
| Evaluate the grid surface at scattered \( (x_i, y_i) \) points. | |
| procedure, private | gridded_eval_one (this, x, y, ierr) |
| Evaluate on a rectangular grid x(:) x y(:) (bispev). Returns f(ny,nx). | |
| procedure, private | gridded_eval_many (this, x, y, ierr) |
| Evaluate the grid surface on a rectangular evaluation grid. | |
| procedure, private | gridded_derivatives_gridded (this, x, y, dx, dy, ierr) |
| Evaluate derivatives at given coordinates. | |
| procedure, private | gridded_derivatives_many (this, x, y, dx, dy, ierr) |
| Evaluate partial derivatives at scattered \( (x_i, y_i) \) points. | |
| procedure, private | gridded_derivatives_one (this, x, y, dx, dy, ierr) |
| Evaluate a partial derivative at a single \( (x, y) \) point. | |
Bivariate surface fitter \( z = s(x, y) \) for gridded data.
Stores grid vectors \( x_i \) ( \( i = 1, \ldots, m_x \)) and \( y_j \) ( \( j = 1, \ldots, m_y \)) together with gridded function values \( z(j, i) \), the fitted tensor-product B-spline representation, and the fitting state. Uses regrid, which is significantly more efficient than surfit for gridded input.
| procedure fitpack_grid_surfaces::fitpack_grid_surface::comm_expand | ( | class(fitpack_grid_surface), intent(inout) | this, |
| real(fp_comm), dimension(:), intent(in) | buffer ) |
Expand grid surface data from a communication buffer.

| procedure fitpack_grid_surfaces::fitpack_grid_surface::comm_pack | ( | class(fitpack_grid_surface), intent(in) | this, |
| real(fp_comm), dimension(:), intent(out) | buffer ) |
Pack grid surface data into a communication buffer.
| procedure fitpack_grid_surfaces::fitpack_grid_surface::comm_size | ( | class(fitpack_grid_surface), intent(in) | this | ) |
Parallel communication.
| procedure fitpack_grid_surfaces::fitpack_grid_surface::cross_section | ( | class(fitpack_grid_surface), intent(in) | this, |
| real(fp_real), intent(in) | u, | ||
| logical, intent(in) | along_y, | ||
| integer(fp_flag), intent(out), optional | ierr ) |
Extract a 1D cross-section curve from the surface.
If along_y=.true., returns \( f(y) = s(u, y) \); otherwise \( g(x) = s(x, u) \).
| procedure fitpack_grid_surfaces::fitpack_grid_surface::derivative_spline | ( | class(fitpack_grid_surface), intent(in) | this, |
| integer(fp_size), intent(in) | nux, | ||
| integer(fp_size), intent(in) | nuy, | ||
| integer(fp_flag), intent(out), optional | ierr ) |
Compute the derivative spline coefficients.
Returns a new grid surface with reduced degrees and trimmed knots.
| procedure fitpack_grid_surfaces::fitpack_grid_surface::destroy | ( | class(fitpack_grid_surface), intent(inout) | this | ) |
Clean memory.
| generic fitpack_grid_surfaces::fitpack_grid_surface::dfdx | ( | class(fitpack_grid_surface), intent(inout) | this, |
| real(fp_real), intent(in) | x, | ||
| real(fp_real), intent(in) | y, | ||
| integer(fp_size), intent(in) | dx, | ||
| integer(fp_size), intent(in) | dy, | ||
| integer(fp_flag), intent(out), optional | ierr ) |
Evaluate a partial derivative at a single \( (x, y) \) point.

| generic fitpack_grid_surfaces::fitpack_grid_surface::dfdx | ( | class(fitpack_grid_surface), intent(inout) | this, |
| real(fp_real), dimension(:), intent(in) | x, | ||
| real(fp_real), dimension(:), intent(in) | y, | ||
| integer(fp_size), intent(in) | dx, | ||
| integer(fp_size), intent(in) | dy, | ||
| integer(fp_flag), intent(out), optional | ierr ) |
Evaluate partial derivatives at scattered \( (x_i, y_i) \) points.
| generic fitpack_grid_surfaces::fitpack_grid_surface::dfdx_ongrid | ( | class(fitpack_grid_surface), intent(inout) | this, |
| real(fp_real), dimension(:), intent(in) | x, | ||
| real(fp_real), dimension(:), intent(in) | y, | ||
| integer(fp_size), intent(in) | dx, | ||
| integer(fp_size), intent(in) | dy, | ||
| integer(fp_flag), intent(out), optional | ierr ) |
| generic fitpack_grid_surfaces::fitpack_grid_surface::eval | ( | class(fitpack_grid_surface), intent(inout) | this, |
| real(fp_real), intent(in) | x, | ||
| real(fp_real), intent(in) | y, | ||
| integer(fp_flag), intent(out), optional | ierr ) |
Evaluate at scattered (x,y) points (bispeu). For gridded output, use eval_ongrid.

| generic fitpack_grid_surfaces::fitpack_grid_surface::eval | ( | class(fitpack_grid_surface), intent(inout) | this, |
| real(fp_real), dimension(:), intent(in) | x, | ||
| real(fp_real), dimension(size(x)), intent(in) | y, | ||
| integer(fp_flag), intent(out), optional | ierr ) |
Evaluate the grid surface at scattered \( (x_i, y_i) \) points.
| generic fitpack_grid_surfaces::fitpack_grid_surface::eval_ongrid | ( | class(fitpack_grid_surface), intent(inout) | this, |
| real(fp_real), intent(in) | x, | ||
| real(fp_real), intent(in) | y, | ||
| integer(fp_flag), intent(out), optional | ierr ) |
| generic fitpack_grid_surfaces::fitpack_grid_surface::eval_ongrid | ( | class(fitpack_grid_surface), intent(inout) | this, |
| real(fp_real), dimension(:), intent(in) | x, | ||
| real(fp_real), dimension(:), intent(in) | y, | ||
| integer(fp_flag), intent(out), optional | ierr ) |
Evaluate the grid surface on a rectangular evaluation grid.
Returns f(j,i) = \( s(x_i, y_j) \).
| procedure fitpack_grid_surfaces::fitpack_grid_surface::fit | ( | class(fitpack_grid_surface), intent(inout) | this, |
| real(fp_real), intent(in), optional | smoothing, | ||
| integer, intent(in), optional | order, | ||
| logical, intent(in), optional | keep_knots ) |
Generate/update fitting curve, with optional /Users/federico/code/fitpack/test/fitpack_curve_tests.f90smoothing.
Uses the regrid core routine, which exploits the rectangular grid structure for efficiency.
|
private |
Evaluate derivatives at given coordinates.
|
private |
Evaluate partial derivatives at scattered \( (x_i, y_i) \) points.
|
private |
Evaluate a partial derivative at a single \( (x, y) \) point.
|
private |
Evaluate the grid surface on a rectangular evaluation grid.
Returns f(j,i) = \( s(x_i, y_j) \).
|
private |
|
private |
Evaluate the grid surface at scattered \( (x_i, y_i) \) points.
|
private |
Evaluate at scattered (x,y) points (bispeu). For gridded output, use eval_ongrid.
| procedure fitpack_grid_surfaces::fitpack_grid_surface::integral | ( | class(fitpack_grid_surface), intent(in) | this, |
| real(fp_real), dimension(2), intent(in) | lower, | ||
| real(fp_real), dimension(2), intent(in) | upper ) |
Double integration over a rectangular domain.
| procedure fitpack_grid_surfaces::fitpack_grid_surface::interpolate | ( | class(fitpack_grid_surface), intent(inout) | this, |
| logical, intent(in), optional | reset_knots ) |
Fit an interpolating gridded surface ( \( s = 0 \)).
| procedure fitpack_grid_surfaces::fitpack_grid_surface::least_squares | ( | class(fitpack_grid_surface), intent(inout) | this, |
| real(fp_real), intent(in), optional | smoothing, | ||
| logical, intent(in), optional | reset_knots ) |
Fit a least-squares gridded surface with fixed knots.
| procedure fitpack_grid_surfaces::fitpack_grid_surface::new_fit | ( | class(fitpack_grid_surface), intent(inout) | this, |
| real(fp_real), dimension(:), intent(in) | x, | ||
| real(fp_real), dimension(:), intent(in) | y, | ||
| real(fp_real), dimension(size(y),size(x)), intent(in) | z, | ||
| real(fp_real), intent(in), optional | smoothing, | ||
| integer, intent(in), optional | order ) |
Generate new fit.
| type(fitpack_grid_surface) function fitpack_grid_surfaces::fitpack_grid_surface::surf_new_from_points | ( | real(fp_real), dimension(:), intent(in) | x, |
| real(fp_real), dimension(:), intent(in) | y, | ||
| real(fp_real), dimension(size(y),size(x)), intent(in) | z, | ||
| integer(fp_flag), intent(out), optional | ierr ) |
Construct a grid surface from gridded data and perform a default fit.

| integer(fp_size), dimension(2) fitpack_grid_surfaces::fitpack_grid_surface::knots = 0 |
| real(fp_real), dimension(2) fitpack_grid_surfaces::fitpack_grid_surface::left |
Interval boundaries.
| integer(fp_size), dimension(2) fitpack_grid_surfaces::fitpack_grid_surface::nest = 0 |
Node weights are not allowed.
| integer(fp_size) fitpack_grid_surfaces::fitpack_grid_surface::nmax = 0 |
| integer(fp_size), dimension(2) fitpack_grid_surfaces::fitpack_grid_surface::order = 3 |
Spline degree.
| real(fp_real), dimension(2) fitpack_grid_surfaces::fitpack_grid_surface::right |
| real(fp_real), dimension(:,:), allocatable fitpack_grid_surfaces::fitpack_grid_surface::t |
| real(fp_real), dimension(:), allocatable fitpack_grid_surfaces::fitpack_grid_surface::x |
The data points.
| real(fp_real), dimension(:), allocatable fitpack_grid_surfaces::fitpack_grid_surface::y |
| real(fp_real), dimension(:,:), allocatable fitpack_grid_surfaces::fitpack_grid_surface::z |