Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
Derived types and interfaces | Functions/Subroutines | Variables
sll_m_low_level_bsplines Module Reference

Description

Low level arbitrary degree splines.

This module defines low level algorithms for arbitrary degree splines. For non periodic boundary conditions, knots defined from a grid are either duplicated (open knot sequence) or mirror outside boundary. It is a selalib implementation of the classical algorithms found in the de Boor book "A Practical Guide to Splines" or the NURBS book by Piegl and Tiller.

Derived types and interfaces

type  sll_t_bsplines
 Information for evaluation of B-splines on non-uniform grid. More...
 

Functions/Subroutines

subroutine, public sll_s_bsplines_init_from_grid (basis, degree, grid, bc_xmin, bc_xmax)
 Build new sll_t_bsplines object based on a on a grid of strictly increasing points including the last point also for periodic domains. More...
 
subroutine, public sll_s_bsplines_init_from_knots (basis, degree, n, knots)
 Build new sll_t_bsplines object. More...
 
pure integer function, public sll_f_find_cell (basis, x)
 return index i of grid cell such that: basisknots(i) <= x <= basisknots(i+1). More...
 
subroutine, public sll_s_bsplines_free (spline)
 Evaluates B-spline values at a point x in a given cell. More...
 

Variables

integer, parameter wp = f64
 Working precision. More...
 
integer, dimension(1:3), parameter allowed_bcs = [sll_p_periodic, sll_p_open, sll_p_mirror]
 List of allowed boundary conditions. More...
 

Function/Subroutine Documentation

◆ sll_f_find_cell()

pure integer function, public sll_m_low_level_bsplines::sll_f_find_cell ( type(sll_t_bsplines), intent(in)  basis,
real(wp), intent(in)  x 
)

return index i of grid cell such that: basisknots(i) <= x <= basisknots(i+1).

@detail If x is not between basisknots(1) and basisknots(basisnum_pts), then the value -1 is returned.

Definition at line 339 of file sll_m_low_level_bsplines.F90.

◆ sll_s_bsplines_free()

subroutine, public sll_m_low_level_bsplines::sll_s_bsplines_free ( type(sll_t_bsplines), intent(inout)  spline)

Evaluates B-spline values at a point x in a given cell.

sll_s_bsplines_eval_basis( basis, cell, x, splines_at_x ) computes the values of all the splines which have support in 'cell' and evaluates them at point 'x', which is supposed to be in cell. The spline object should have already been initialized and will contain information on the spline degree to use and the type of boundary condition desired. The algorithm implemented is numerically stable and known as The Cox - de Boor algorithm, which is a generalisation to splines of the de Casteljau algorithm for Bezier curves.

Returns
b_spline_at_x B-spline values

returns first derivative values at x of all b-splines with support in cell

sll_s_bsplines_eval_deriv returns an array with the derivative values of the B-splines of a requested order that are supported in 'cell' and evaluated at 'x'. Algorithm derived from algorithm A3.2 of NURBS book The return value has the format:

\[ B'[deg,i-deg](x), B'[deg,i-deg+1](x), ..., B'[deg,i](x) \]

where 'deg' is the degree of the spline.

Parameters
[out]bsdxB-spline derivatives

returns splines and first derivatives

See sll_s_bsplines_eval_deriv and sll_f_splines_at_x

Returns
b_spline_and_derivs_at_x B-spline values and derivatives

returns splines and first derivatives

See sll_s_bsplines_eval_deriv and sll_f_splines_at_x

Parameters
[in]basisbspline object
[in]icellcell where bsplines are to be computed
[in]xvalue of point where bsplines are to be computed
[in]nnumber of derivatives to be computed
[out]bsdxB-spline values and the first n derivatives

Alternative direct implentation of recursion formula.

@detail This provides an evaluation of B-splines directly based on the recurrence formula. It is 10% faster than the classical Cox - de Boor formula that is implented in sll_f_splines_at_x, but can have numerical stability issues. For this reason the Cox - de Boor formula should be the default implementation

Alternative direct implentation of recursion formula.

@detail This provides an evaluation of B-splines directly based on the recurrence formula. It is about 80% faster than the classical Cox - de Boor formula that is implented in sll_f_splines_at_x, but can have numerical stability issues. For this reason the Cox - de Boor formula should be the default implementation

Definition at line 842 of file sll_m_low_level_bsplines.F90.

◆ sll_s_bsplines_init_from_grid()

subroutine, public sll_m_low_level_bsplines::sll_s_bsplines_init_from_grid ( type(sll_t_bsplines), intent(out)  basis,
integer, intent(in)  degree,
real(wp), dimension(:), intent(in)  grid,
integer, intent(in)  bc_xmin,
integer, intent(in)  bc_xmax 
)

Build new sll_t_bsplines object based on a on a grid of strictly increasing points including the last point also for periodic domains.

Parameters
[in,out]basisarbitrary_degree_spline_1d object
[in]degreespline_degree
[in]gridx coordinates of grid points
[in]num_ptsnumber of grid points
[in]bc_xminboundary condition at xmin [periodic/open/mirror]
[in]bc_xmaxboundary condition at xmax [periodic/open/mirror]

The logic behind the periodic boundary condition is the following. The given grid array has minimum (grid(1)) and maximum (grid(n)) values at either end. This defines a length 'L'. If interpreted as a periodic space, this is also the period. Thus, as we extend the number of knots at both ends of the given array, we use the periodicity condition to fill out the new values:

               .
               .
               .
      knots(-1) = knots(n-2) - L
      knots( 0) = knots(n-1) - L
               .
               .
               .
      knots(n+1) = knots(1) + L
      knots(n+2) = knots(2) + L
               .
               .
               .

The 'open' boundary condition simply extends the new values of the local array at both ends with repeated endpoint values. That is

... = knots(-2) = knots(-1) = knots(0) = knots(1)

and

knots(n+1) = knots(n+2) = knots(n+3) = ... = knots(n)

The mirror boundary condition mirrors the knot values on each side of the grid of the local array at both ends with repeated endpoint values. That is

... =  = knots(-1) = knots(0) = knots(1)

and

knots(n+1) = knots(n+2) = knots(n+3) = ... = knots(n)

Definition at line 145 of file sll_m_low_level_bsplines.F90.

◆ sll_s_bsplines_init_from_knots()

subroutine, public sll_m_low_level_bsplines::sll_s_bsplines_init_from_knots ( type(sll_t_bsplines), intent(out)  basis,
integer, intent(in)  degree,
integer, intent(in)  n,
real(wp), dimension(:), intent(in)  knots 
)

Build new sll_t_bsplines object.

based on a given array of knots which is a non decreasing array

Parameters
[in,out]basisspline object
[in]degreespline degree
[in]ndimension of spline space
[in]knotsarray of give knots

Definition at line 289 of file sll_m_low_level_bsplines.F90.

Variable Documentation

◆ allowed_bcs

integer, dimension(1:3), parameter allowed_bcs = [sll_p_periodic, sll_p_open, sll_p_mirror]
private

List of allowed boundary conditions.

Definition at line 72 of file sll_m_low_level_bsplines.F90.

◆ wp

integer, parameter wp = f64
private

Working precision.

Definition at line 69 of file sll_m_low_level_bsplines.F90.

    Report Typos and Errors