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

Bicubic spline fitter for data on a gridded polar disc. More...

Inheritance diagram for fitpack_gridded_polar::fitpack_grid_polar:
Collaboration diagram for fitpack_gridded_polar::fitpack_grid_polar:

Public Member Functions

procedure destroy (this)
 Clean memory.
 
procedure new_points (this, u, v, r, z, z0)
 Set new points.
 
procedure new_fit (this, u, v, r, z, z0, smoothing)
 Generate new fit.
 
procedure set_origin_bc (this, z0, exact, differentiable)
 Set origin BC.
 
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 polar surface with fixed knots.
 
procedure interpolate (this, reset_knots)
 Fit an interpolating gridded polar 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 polar surface on a rectangular \( u \times v \) grid.
 
procedure comm_size (this)
 Parallel communication interface.
 
procedure comm_pack (this, buffer)
 Pack gridded polar data into a communication buffer.
 
procedure comm_expand (this, buffer)
 Expand gridded polar data from a communication buffer.
 
procedure write (this, filename)
 Write to disk.
 
type(fitpack_grid_polar) function surf_new_from_points (u, v, r, z, z0, ierr)
 Construct a gridded polar surface from 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 u
 Coordinates of the data points in grid coordinates (u,v) and domain size.
 
real(fp_real), dimension(:), allocatable v
 
real(fp_realr = zero
 
real(fp_real), dimension(:,:), allocatable z
 Function values at the points.
 
real(fp_realz0 = zero
 Function value at the origin.
 
logical z0_present = .false.
 Fit spline to z0 exactly at the origin.
 
logical z0_exact = .false.
 
logical z0_zero_gradient = .false.
 Decide if the origin should have zero-gradient BCs in both u and v.
 
integer, dimension(2) nest = 0
 Node weights are not allowed.
 
integer nmax = 0
 
integer, dimension(2) knots = 0
 
real(fp_real), dimension(:,:), allocatable t
 
integer bc_continuity_origin = 1
 
integer bc_boundary = OUTSIDE_EXTRAPOLATE
 
- 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 polar surface on a rectangular \( u \times v \) grid.
 

Detailed Description

Bicubic spline fitter for data on a gridded polar disc.

Stores the grid coordinates \( u_i \), \( v_j \), the constant boundary radius \( r \), gridded function values \( z(i, j) \), and the fitted B-spline representation. The origin value \( z_0 \) may optionally be prescribed (exactly or approximately), and zero-gradient boundary conditions at the origin can be enforced.

See also
Dierckx, Ch. 11, §11.1 (pp. 255–263)

Member Function/Subroutine Documentation

◆ comm_expand()

procedure fitpack_gridded_polar::fitpack_grid_polar::comm_expand ( class(fitpack_grid_polar), intent(inout) this,
real(fp_comm), dimension(:), intent(in) buffer )

Expand gridded polar data from a communication buffer.

◆ comm_pack()

procedure fitpack_gridded_polar::fitpack_grid_polar::comm_pack ( class(fitpack_grid_polar), intent(in) this,
real(fp_comm), dimension(:), intent(out) buffer )

Pack gridded polar data into a communication buffer.

◆ comm_size()

procedure fitpack_gridded_polar::fitpack_grid_polar::comm_size ( class(fitpack_grid_polar), intent(in) this)

Parallel communication interface.

◆ destroy()

procedure fitpack_gridded_polar::fitpack_grid_polar::destroy ( class(fitpack_grid_polar), intent(inout) this)

Clean memory.

◆ eval() [1/2]

generic fitpack_gridded_polar::fitpack_grid_polar::eval ( class(fitpack_grid_polar), 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.

See also
bispev
Here is the call graph for this function:

◆ eval() [2/2]

generic fitpack_gridded_polar::fitpack_grid_polar::eval ( class(fitpack_grid_polar), 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 polar surface on a rectangular \( u \times v \) grid.

Returns f(j,i) = \( s(u_i, v_j) \) in normalized polar coordinates.

See also
bispev

◆ fit()

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

Generate/update fitting curve, with optional smoothing.

Uses the pogrid core routine. Origin continuity constraints and optional function value at the origin are configured via set_origin_BC.

See also
pogrid

◆ gridded_eval_many()

procedure, private fitpack_gridded_polar::fitpack_grid_polar::gridded_eval_many ( class(fitpack_grid_polar), 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 polar surface on a rectangular \( u \times v \) grid.

Returns f(j,i) = \( s(u_i, v_j) \) in normalized polar coordinates.

See also
bispev

◆ gridded_eval_one()

procedure, private fitpack_gridded_polar::fitpack_grid_polar::gridded_eval_one ( class(fitpack_grid_polar), 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.

See also
bispev

◆ interpolate()

procedure fitpack_gridded_polar::fitpack_grid_polar::interpolate ( class(fitpack_grid_polar), intent(inout) this,
logical, intent(in), optional reset_knots )

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

See also
pogrid

◆ least_squares()

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

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

See also
pogrid

◆ new_fit()

procedure fitpack_gridded_polar::fitpack_grid_polar::new_fit ( class(fitpack_grid_polar), intent(inout) this,
real(fp_real), dimension(:), intent(in) u,
real(fp_real), dimension(:), intent(in) v,
real(fp_real), intent(in) r,
real(fp_real), dimension(size(v),size(u)), intent(in) z,
real(fp_real), intent(in), optional z0,
real(fp_real), intent(in), optional smoothing )

Generate new fit.

See also
pogrid

◆ new_points()

procedure fitpack_gridded_polar::fitpack_grid_polar::new_points ( class(fitpack_grid_polar), intent(inout) this,
real(fp_real), dimension(:), intent(in) u,
real(fp_real), dimension(:), intent(in) v,
real(fp_real), intent(in) r,
real(fp_real), dimension(size(v),size(u)), intent(in) z,
real(fp_real), intent(in), optional z0 )

Set new points.

Parameters
[in,out]thisThe grid polar surface (destroyed and reinitialized).
[in]uRadial grid values \( u_i \in [0, 1] \).
[in]vAngular grid values \( v_j \in [-\pi, \pi] \).
[in]rConstant boundary radius.
[in]zGridded function values z(j,i).
[in]z0Optional function value at the origin.
See also
pogrid

◆ set_origin_bc()

procedure fitpack_gridded_polar::fitpack_grid_polar::set_origin_bc ( class(fitpack_grid_polar), intent(inout) this,
real(fp_real), intent(in), optional z0,
logical, intent(in), optional exact,
logical, intent(in), optional differentiable )

Set origin BC.

Parameters
[in,out]thisThe grid polar surface.
[in]z0Optional function value at the origin.
[in]exactIf .true., fit the origin value exactly.
[in]differentiableIf .true., enforce \( C^1 \) continuity at the origin.

◆ surf_new_from_points()

type(fitpack_grid_polar) function fitpack_gridded_polar::fitpack_grid_polar::surf_new_from_points ( real(fp_real), dimension(:), intent(in) u,
real(fp_real), dimension(:), intent(in) v,
real(fp_real), intent(in) r,
real(fp_real), dimension(size(v),size(u)), intent(in) z,
real(fp_real), intent(in), optional z0,
integer, intent(out), optional ierr )

Construct a gridded polar surface from data and perform a default fit.

See also
pogrid
Here is the call graph for this function:

◆ write()

procedure fitpack_gridded_polar::fitpack_grid_polar::write ( class(fitpack_grid_polar), intent(inout) this,
character(*), intent(in) filename )

Write to disk.

Here is the call graph for this function:

Member Data Documentation

◆ bc_boundary

integer fitpack_gridded_polar::fitpack_grid_polar::bc_boundary = OUTSIDE_EXTRAPOLATE

◆ bc_continuity_origin

integer fitpack_gridded_polar::fitpack_grid_polar::bc_continuity_origin = 1

◆ knots

integer, dimension(2) fitpack_gridded_polar::fitpack_grid_polar::knots = 0

◆ nest

integer, dimension(2) fitpack_gridded_polar::fitpack_grid_polar::nest = 0

Node weights are not allowed.

◆ nmax

integer fitpack_gridded_polar::fitpack_grid_polar::nmax = 0

◆ r

real(fp_real) fitpack_gridded_polar::fitpack_grid_polar::r = zero

◆ t

real(fp_real), dimension(:,:), allocatable fitpack_gridded_polar::fitpack_grid_polar::t

◆ u

real(fp_real), dimension(:), allocatable fitpack_gridded_polar::fitpack_grid_polar::u

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

◆ v

real(fp_real), dimension(:), allocatable fitpack_gridded_polar::fitpack_grid_polar::v

◆ z

real(fp_real), dimension(:,:), allocatable fitpack_gridded_polar::fitpack_grid_polar::z

Function values at the points.

◆ z0

real(fp_real) fitpack_gridded_polar::fitpack_grid_polar::z0 = zero

Function value at the origin.

◆ z0_exact

logical fitpack_gridded_polar::fitpack_grid_polar::z0_exact = .false.

◆ z0_present

logical fitpack_gridded_polar::fitpack_grid_polar::z0_present = .false.

Fit spline to z0 exactly at the origin.

◆ z0_zero_gradient

logical fitpack_gridded_polar::fitpack_grid_polar::z0_zero_gradient = .false.

Decide if the origin should have zero-gradient BCs in both u and v.


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