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_spline_interpolator_1d Module Reference

Description

Interpolator for 1D splines of arbitrary degree, on uniform and non-uniform grids.

Author
Yaman Güçlü - IPP Garching
Edoardo Zoni - IPP Garching

Derived types and interfaces

type  sll_t_spline_interpolator_1d
 1D spline interpolator More...
 

Functions/Subroutines

subroutine, public sll_s_spline_1d_compute_num_cells (degree, bc_xmin, bc_xmax, nipts, ncells)
 Calculate number of cells from number of interpolation points. More...
 
subroutine s_spline_interpolator_1d__init (self, bspl, bc_xmin, bc_xmax)
 Initialize a 1D spline interpolator object. More...
 
subroutine s_spline_interpolator_1d__free (self)
 Destroy local objects and free allocated memory. More...
 
subroutine s_spline_interpolator_1d__get_interp_points (self, tau)
 Get coordinates of interpolation points (1D grid) More...
 
subroutine s_spline_interpolator_1d__compute_interpolant (self, spline, gtau, derivs_xmin, derivs_xmax)
 Compute interpolating 1D spline. More...
 
subroutine s_build_system (self, matrix)
 Private subroutine for assembling and factorizing linear system for any kind of boundary conditions at xmin and xmax. More...
 
subroutine s_compute_interpolation_points_uniform (self, tau)
 
subroutine s_compute_num_diags_uniform (self, kl, ku)
 
subroutine s_compute_interpolation_points_non_uniform (self, tau)
 
subroutine s_compute_num_diags_non_uniform (self, kl, ku)
 

Variables

integer, parameter wp = f64
 Working precision. More...
 
integer, dimension(1:3), parameter allowed_bcs = [sll_p_periodic, sll_p_hermite, sll_p_greville]
 Allowed boundary conditions. More...
 

Function/Subroutine Documentation

◆ s_build_system()

subroutine sll_m_spline_interpolator_1d::s_build_system ( class(sll_t_spline_interpolator_1d), intent(in)  self,
class(sll_c_spline_matrix), intent(inout)  matrix 
)
private

Private subroutine for assembling and factorizing linear system for any kind of boundary conditions at xmin and xmax.

Parameters
[in]self1D spline interpolator object
[in,out]matrixgeneric 'spline' matrix (dense/banded/periodic-banded)

Definition at line 330 of file sll_m_spline_interpolator_1d.F90.

Here is the caller graph for this function:

◆ s_compute_interpolation_points_non_uniform()

subroutine sll_m_spline_interpolator_1d::s_compute_interpolation_points_non_uniform ( class(sll_t_spline_interpolator_1d), intent(in)  self,
real(wp), dimension(:), intent(out), allocatable  tau 
)
private

Definition at line 546 of file sll_m_spline_interpolator_1d.F90.

Here is the caller graph for this function:

◆ s_compute_interpolation_points_uniform()

subroutine sll_m_spline_interpolator_1d::s_compute_interpolation_points_uniform ( class(sll_t_spline_interpolator_1d), intent(in)  self,
real(wp), dimension(:), intent(out), allocatable  tau 
)
private

Definition at line 427 of file sll_m_spline_interpolator_1d.F90.

Here is the caller graph for this function:

◆ s_compute_num_diags_non_uniform()

subroutine sll_m_spline_interpolator_1d::s_compute_num_diags_non_uniform ( class(sll_t_spline_interpolator_1d), intent(in)  self,
integer, intent(out)  kl,
integer, intent(out)  ku 
)
private

Definition at line 627 of file sll_m_spline_interpolator_1d.F90.

Here is the caller graph for this function:

◆ s_compute_num_diags_uniform()

subroutine sll_m_spline_interpolator_1d::s_compute_num_diags_uniform ( class(sll_t_spline_interpolator_1d), intent(in)  self,
integer, intent(out)  kl,
integer, intent(out)  ku 
)
private

Definition at line 513 of file sll_m_spline_interpolator_1d.F90.

Here is the caller graph for this function:

◆ s_spline_interpolator_1d__compute_interpolant()

subroutine sll_m_spline_interpolator_1d::s_spline_interpolator_1d__compute_interpolant ( class(sll_t_spline_interpolator_1d), intent(in)  self,
type(sll_t_spline_1d), intent(inout)  spline,
real(wp), dimension(:), intent(in)  gtau,
real(wp), dimension(:), intent(in), optional  derivs_xmin,
real(wp), dimension(:), intent(in), optional  derivs_xmax 
)
private

Compute interpolating 1D spline.

Compute coefficients of 1D spline that interpolates function values on grid. If Hermite BCs are used, function derivatives at appropriate boundary are also needed.

Parameters
[in]self1D spline interpolator
[in,out]spline1D spline
[in]gtaufunction values of interpolation points
[in]derivs_xmin(optional) array with boundary conditions at xmin
[in]derivs_xmax(optional) array with boundary conditions at xmax

Definition at line 243 of file sll_m_spline_interpolator_1d.F90.

◆ s_spline_interpolator_1d__free()

subroutine sll_m_spline_interpolator_1d::s_spline_interpolator_1d__free ( class(sll_t_spline_interpolator_1d), intent(inout)  self)
private

Destroy local objects and free allocated memory.

Parameters
[in,out]self1D spline interpolator

Definition at line 202 of file sll_m_spline_interpolator_1d.F90.

◆ s_spline_interpolator_1d__get_interp_points()

subroutine sll_m_spline_interpolator_1d::s_spline_interpolator_1d__get_interp_points ( class(sll_t_spline_interpolator_1d), intent(in)  self,
real(wp), dimension(:), intent(out), allocatable  tau 
)
private

Get coordinates of interpolation points (1D grid)

Parameters
[in]self1D spline interpolator
[out]taux coordinates of interpolation points

Definition at line 221 of file sll_m_spline_interpolator_1d.F90.

◆ s_spline_interpolator_1d__init()

subroutine sll_m_spline_interpolator_1d::s_spline_interpolator_1d__init ( class(sll_t_spline_interpolator_1d), intent(out)  self,
class(sll_c_bsplines), intent(in), target  bspl,
integer, intent(in)  bc_xmin,
integer, intent(in)  bc_xmax 
)
private

Initialize a 1D spline interpolator object.

Parameters
[out]self1D spline interpolator
[in]bsplB-splines (basis)
[in]bc_xminboundary condition at xmin
[in]bc_xmaxboundary condition at xmax

Definition at line 124 of file sll_m_spline_interpolator_1d.F90.

Here is the call graph for this function:

◆ sll_s_spline_1d_compute_num_cells()

subroutine, public sll_m_spline_interpolator_1d::sll_s_spline_1d_compute_num_cells ( integer, intent(in)  degree,
integer, intent(in)  bc_xmin,
integer, intent(in)  bc_xmax,
integer, intent(in)  nipts,
integer, intent(out)  ncells 
)

Calculate number of cells from number of interpolation points.

Important for parallelization: for given spline degree and BCs, calculate the number of grid cells that yields the desired number of interpolation points

Parameters
[in]degreespline degree
[in]bc_xminboundary condition type at left boundary (x=xmin)
[in]bc_xmaxboundary condition type at right boundary (x=xmax)
[in]niptsdesired number of interpolation points
[out]ncellscalculated number of cells in domain

Definition at line 84 of file sll_m_spline_interpolator_1d.F90.

Here is the caller graph for this function:

Variable Documentation

◆ allowed_bcs

integer, dimension(1:3), parameter allowed_bcs = [sll_p_periodic, sll_p_hermite, sll_p_greville]
private

Allowed boundary conditions.

Definition at line 39 of file sll_m_spline_interpolator_1d.F90.

◆ wp

integer, parameter wp = f64
private

Working precision.

Definition at line 36 of file sll_m_spline_interpolator_1d.F90.

    Report Typos and Errors