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

Bicubic spline fitter for data on a latitude-longitude grid. More...

Inheritance diagram for fitpack_gridded_sphere::fitpack_grid_sphere:
Collaboration diagram for fitpack_gridded_sphere::fitpack_grid_sphere:

Public Member Functions

procedure destroy (this)
 Clean memory.
 
procedure new_points (this, u, v, z)
 Set new points.
 
procedure new_fit (this, u, v, z, smoothing)
 Generate new fit.
 
procedure bc_north_pole (this, z0, exact, differentiable, zero_grad)
 Set pole BCs.
 
procedure bc_south_pole (this, z0, exact, differentiable, zero_grad)
 Set boundary conditions at the south pole ( \( u = \pi \)).
 
procedure fit (this, smoothing, keep_knots)
 Generate/update fitting curve, with optional smoothing.
 
procedure least_squares (this, smoothing, reset_knots)
 Fit a least-squares gridded spherical surface with fixed knots.
 
procedure interpolate (this, reset_knots)
 Fit an interpolating gridded spherical surface ( \( s = 0 \)).
 
generic eval (this, u, v, ierr)
 Evaluate gridded domain at given x,y coordinates.
 
generic eval (this, u, v, ierr)
 Evaluate the gridded spherical spline on a grid of colatitude and longitude values.
 
procedure comm_size (this)
 Parallel communication interface.
 
procedure comm_pack (this, buffer)
 Pack the gridded sphere fitter state into a communication buffer.
 
procedure comm_expand (this, buffer)
 Restore the gridded sphere fitter state from a communication buffer.
 
procedure write (this, filename)
 Write to disk.
 
type(fitpack_grid_sphere) function spgrid_new_from_points (u, v, z, ierr)
 Construct a fitpack_grid_sphere from gridded data and perform an initial 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 u
 Coordinates of the data points in grid coordinates (u,v) and domain size.
 
real(fp_real), dimension(:), allocatable v
 
real(fp_real), dimension(:,:), allocatable z
 Function values at the data points z(iv,iu) note (v,u) indices.
 
real(fp_real), dimension(2) pole_z0 = zero
 North pole BC (u=0, index=1), South pole BC (u=pi, index=2)
 
logical, dimension(2) pole_present = .false.
 
logical, dimension(2) pole_exct = .false.
 
integer, dimension(2) pole_continuity = 1
 
logical, dimension(2) pole_zero_grad = .true.
 
integer, dimension(2) nest = 0
 Node weights are not allowed.
 
integer nmax = 0
 
integer, dimension(2) 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 gridded_eval_one (this, u, v, ierr)
 Evaluate gridded domain at given x,y coordinates.
 
procedure, private gridded_eval_many (this, u, v, ierr)
 Evaluate the gridded spherical spline on a grid of colatitude and longitude values.
 

Detailed Description

Bicubic spline fitter for data on a latitude-longitude grid.

Stores the colatitude grid \( u_i \), the \( 2\pi \)-periodic longitude grid \( v_j \), gridded function values \( z(j, i) \), and the fitted B-spline representation. North and south pole boundary conditions (function value, continuity order, gradient vanishing) are independently configurable.

See also
Dierckx, Ch. 11, §11.2 (pp. 263–269)

Member Function/Subroutine Documentation

◆ bc_north_pole()

procedure fitpack_gridded_sphere::fitpack_grid_sphere::bc_north_pole ( class(fitpack_grid_sphere), intent(inout) this,
real(fp_real), intent(in), optional z0,
logical, intent(in), optional exact,
logical, intent(in), optional differentiable,
logical, intent(in), optional zero_grad )

Set pole BCs.

◆ bc_south_pole()

procedure fitpack_gridded_sphere::fitpack_grid_sphere::bc_south_pole ( class(fitpack_grid_sphere), intent(inout) this,
real(fp_real), intent(in), optional z0,
logical, intent(in), optional exact,
logical, intent(in), optional differentiable,
logical, intent(in), optional zero_grad )

Set boundary conditions at the south pole ( \( u = \pi \)).

◆ comm_expand()

procedure fitpack_gridded_sphere::fitpack_grid_sphere::comm_expand ( class(fitpack_grid_sphere), intent(inout) this,
real(fp_comm), dimension(:), intent(in) buffer )

Restore the gridded sphere fitter state from a communication buffer.

◆ comm_pack()

procedure fitpack_gridded_sphere::fitpack_grid_sphere::comm_pack ( class(fitpack_grid_sphere), intent(in) this,
real(fp_comm), dimension(:), intent(out) buffer )

Pack the gridded sphere fitter state into a communication buffer.

◆ comm_size()

procedure fitpack_gridded_sphere::fitpack_grid_sphere::comm_size ( class(fitpack_grid_sphere), intent(in) this)

Parallel communication interface.

◆ destroy()

procedure fitpack_gridded_sphere::fitpack_grid_sphere::destroy ( class(fitpack_grid_sphere), intent(inout) this)

Clean memory.

◆ eval() [1/2]

generic fitpack_gridded_sphere::fitpack_grid_sphere::eval ( class(fitpack_grid_sphere), intent(inout) this,
real(fp_real), intent(in) u,
real(fp_real), intent(in) v,
integer, intent(out), optional ierr )

Evaluate gridded domain at given x,y coordinates.

Parameters
[in]uColatitude \( u \in [0, \pi] \).
[in]vLongitude \( v \).
[out]ierrOptional error flag.
Returns
Spline value \( s(u, v) \).
See also
bispev
Here is the call graph for this function:

◆ eval() [2/2]

generic fitpack_gridded_sphere::fitpack_grid_sphere::eval ( class(fitpack_grid_sphere), intent(inout) this,
real(fp_real), dimension(:), intent(in) u,
real(fp_real), dimension(:), intent(in) v,
integer, intent(out), optional ierr )

Evaluate the gridded spherical spline on a grid of colatitude and longitude values.

Returns \( f(j,i) = s(u_i, v_j) \) for all combinations of the input vectors (tensor-product evaluation).

Parameters
[in]uColatitude evaluation points.
[in]vLongitude evaluation points.
[out]ierrOptional error flag.
Returns
Grid of spline values, dimensioned (size(v), size(u)).
See also
bispev

◆ fit()

procedure fitpack_gridded_sphere::fitpack_grid_sphere::fit ( class(fitpack_grid_sphere), intent(inout) this,
real(fp_real), intent(in), optional smoothing,
logical, intent(in), optional keep_knots )

Generate/update fitting curve, with optional smoothing.

Iterates over the smoothing schedule, calling the spgrid core routine to fit a bicubic spline on the sphere grid. Pole boundary conditions (continuity order, gradient vanishing, exact data) are applied as configured.

Parameters
[in]smoothingSmoothing factor ( \( s \ge 0 \)); default uses stored value.
[in]keep_knotsIf .true., reuse the current knot set.
Returns
Error flag.
See also
spgrid

◆ gridded_eval_many()

procedure, private fitpack_gridded_sphere::fitpack_grid_sphere::gridded_eval_many ( class(fitpack_grid_sphere), intent(inout) this,
real(fp_real), dimension(:), intent(in) u,
real(fp_real), dimension(:), intent(in) v,
integer, intent(out), optional ierr )
private

Evaluate the gridded spherical spline on a grid of colatitude and longitude values.

Returns \( f(j,i) = s(u_i, v_j) \) for all combinations of the input vectors (tensor-product evaluation).

Parameters
[in]uColatitude evaluation points.
[in]vLongitude evaluation points.
[out]ierrOptional error flag.
Returns
Grid of spline values, dimensioned (size(v), size(u)).
See also
bispev

◆ gridded_eval_one()

procedure, private fitpack_gridded_sphere::fitpack_grid_sphere::gridded_eval_one ( class(fitpack_grid_sphere), intent(inout) this,
real(fp_real), intent(in) u,
real(fp_real), intent(in) v,
integer, intent(out), optional ierr )
private

Evaluate gridded domain at given x,y coordinates.

Parameters
[in]uColatitude \( u \in [0, \pi] \).
[in]vLongitude \( v \).
[out]ierrOptional error flag.
Returns
Spline value \( s(u, v) \).
See also
bispev

◆ interpolate()

procedure fitpack_gridded_sphere::fitpack_grid_sphere::interpolate ( class(fitpack_grid_sphere), intent(inout) this,
logical, intent(in), optional reset_knots )

Fit an interpolating gridded spherical surface ( \( s = 0 \)).

Parameters
[in]reset_knotsIf .true. (default), start with fresh knot placement.
Returns
Error flag.
See also
spgrid

◆ least_squares()

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

Fit a least-squares gridded spherical surface with fixed knots.

Parameters
[in]smoothingOptional smoothing factor for the reset pass.
[in]reset_knotsIf .true., recompute knots via a smoothing fit first.
Returns
Error flag.
See also
spgrid

◆ new_fit()

procedure fitpack_gridded_sphere::fitpack_grid_sphere::new_fit ( class(fitpack_grid_sphere), intent(inout) this,
real(fp_real), dimension(:), intent(in) u,
real(fp_real), dimension(:), intent(in) v,
real(fp_real), dimension(size(v),size(u)), intent(in) z,
real(fp_real), intent(in), optional smoothing )

Generate new fit.

Parameters
[in]uColatitude grid.
[in]vLongitude grid.
[in]zGridded function values.
[in]smoothingOptional smoothing factor.
Returns
Error flag.
See also
spgrid

◆ new_points()

procedure fitpack_gridded_sphere::fitpack_grid_sphere::new_points ( class(fitpack_grid_sphere), intent(inout) this,
real(fp_real), dimension(:), intent(in) u,
real(fp_real), dimension(:), intent(in) v,
real(fp_real), dimension(size(v),size(u)), intent(in) z )

Set new points.

Parameters
[in]uColatitude grid \( u_i \in [0, \pi] \).
[in]vLongitude grid \( v_j \) ( \( 2\pi \)-periodic).
[in]zFunction values \( z(j, i) \), dimensioned (size(v), size(u)).

◆ spgrid_new_from_points()

type(fitpack_grid_sphere) function fitpack_gridded_sphere::fitpack_grid_sphere::spgrid_new_from_points ( real(fp_real), dimension(:), intent(in) u,
real(fp_real), dimension(:), intent(in) v,
real(fp_real), dimension(size(v),size(u)), intent(in) z,
integer, intent(out), optional ierr )

Construct a fitpack_grid_sphere from gridded data and perform an initial fit.

Parameters
[in]uColatitude grid \( u_i \in [0, \pi] \).
[in]vLongitude grid \( v_j \) ( \( 2\pi \)-periodic).
[in]zFunction values \( z(j, i) \), dimensioned (size(v), size(u)).
[out]ierrOptional error flag; if absent, halts on error.
See also
spgrid
Here is the call graph for this function:

◆ write()

procedure fitpack_gridded_sphere::fitpack_grid_sphere::write ( class(fitpack_grid_sphere), intent(inout) this,
character(*), intent(in) filename )

Write to disk.

Outputs a table with colatitude grid values in rows and longitude grid values in columns, suitable for external plotting or inspection.

Parameters
[in]fileNameOutput file path.
Here is the call graph for this function:

Member Data Documentation

◆ knots

integer, dimension(2) fitpack_gridded_sphere::fitpack_grid_sphere::knots = 0

◆ nest

integer, dimension(2) fitpack_gridded_sphere::fitpack_grid_sphere::nest = 0

Node weights are not allowed.

◆ nmax

integer fitpack_gridded_sphere::fitpack_grid_sphere::nmax = 0

◆ pole_continuity

integer, dimension(2) fitpack_gridded_sphere::fitpack_grid_sphere::pole_continuity = 1

◆ pole_exct

logical, dimension(2) fitpack_gridded_sphere::fitpack_grid_sphere::pole_exct = .false.

◆ pole_present

logical, dimension(2) fitpack_gridded_sphere::fitpack_grid_sphere::pole_present = .false.

◆ pole_z0

real(fp_real), dimension(2) fitpack_gridded_sphere::fitpack_grid_sphere::pole_z0 = zero

North pole BC (u=0, index=1), South pole BC (u=pi, index=2)

◆ pole_zero_grad

logical, dimension(2) fitpack_gridded_sphere::fitpack_grid_sphere::pole_zero_grad = .true.

◆ t

real(fp_real), dimension(:,:), allocatable fitpack_gridded_sphere::fitpack_grid_sphere::t

◆ u

real(fp_real), dimension(:), allocatable fitpack_gridded_sphere::fitpack_grid_sphere::u

Coordinates of the data points in grid coordinates (u,v) and domain size.

◆ v

real(fp_real), dimension(:), allocatable fitpack_gridded_sphere::fitpack_grid_sphere::v

◆ z

real(fp_real), dimension(:,:), allocatable fitpack_gridded_sphere::fitpack_grid_sphere::z

Function values at the data points z(iv,iu) note (v,u) indices.


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