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

Description

Provides capabilities for values and derivatives interpolation with box splines on hexagonal meshes.

Author
Laura S. Mendoza

This module contains the computation of box splines of arbitrary degree. There is a special optimized algorithm for degree 2. The box splines here are defined only on hexagonal meshes (regular equilateral triangles elements). The only boundary condition supported right now is dirichlet. Reference: http://ieeexplore.ieee.org/document/1642713/

Derived types and interfaces

type  sll_t_box_spline_2d
 Basic type for 2 dimensional box splines. More...
 
interface  sll_o_delete
 Generic sub-routine defined for 2D box spline types. Deallocates the memory associated with the given box spline object. More...
 

Functions/Subroutines

type(sll_t_box_spline_2d) function, pointer, public sll_f_new_box_spline_2d (mesh, bc_type)
 Creates a box spline element. More...
 
subroutine compute_coeff_box_spline_2d (spline, data, deg)
 Computes box splines coefficients. More...
 
subroutine compute_coeff_box_spline_2d_diri (spline, data, deg)
 Computes box splines coefficients for dirichlet BC. More...
 
subroutine compute_coeff_box_spline_2d_prdc (spline, data, deg)
 Computes box splines coefficients with periodic BC. More...
 
subroutine compute_coeff_box_spline_2d_neum (spline, data, deg)
 Computes box splines coefficients with neumann BC. More...
 
function choose (n, k)
 Computes the binomial coefficient (n, k) More...
 
function chi_gen_val (x1_in, x2_in, deg)
 Computes the box spline of degree deg on (x1, x2) More...
 
real(kind=f64) function, public sll_f_compute_box_spline (spline, x1, x2, deg)
 Computes the value of a box spline. More...
 
real(kind=f64) function change_basis_x1 (spline, x1, x2)
 1st coo. of (x1, x2) in reference hex-mesh coo. More...
 
real(kind=f64) function change_basis_x2 (spline, x1, x2)
 2nd coo. of (x1, x2) in reference hex-mesh coo. More...
 
real(kind=f64) function, public sll_f_hex_interpolate_value (mesh_geom, x1, x2, spline, deg)
 Interpolates on point of cartesian coordinates (x1, x2) More...
 
integer(kind=i32) function, dimension(3 *deg *deg) non_zeros_splines (mesh, cell_index, deg)
 Computes indices of non null splines on a given cell. More...
 
real(kind=f64) function, public sll_f_boxspline_x1_derivative (x1, x2, deg)
 Computes x-derivative on (x,y) More...
 
real(kind=f64) function, public sll_f_boxspline_x2_derivative (x1, x2, deg)
 Computes y-derivative on (x,y) More...
 
real(kind=f64) function, public sll_f_boxspline_val_der (x1, x2, deg, nderiv1, nderiv2)
 Computes the values or derivatives of box splines. More...
 
subroutine, public sll_s_write_connectivity (mesh, deg)
 Writes connectivity for CAID / DJANGO. More...
 
subroutine delete_box_spline_2d (spline)
 Generic sub-routine defined for 2D box spline types. Deallocates the memory associated with the given box spline object. More...
 

Function/Subroutine Documentation

◆ change_basis_x1()

real(kind=f64) function sll_m_box_splines::change_basis_x1 ( type(sll_t_box_spline_2d), pointer  spline,
real(kind=f64), intent(in)  x1,
real(kind=f64), intent(in)  x2 
)
private

1st coo. of (x1, x2) in reference hex-mesh coo.

This function allows to change a point of coordinates (x1, x2) on the spline basis to the mesh basis. Gives 1st coordinate.

Parameters
[in]splinebox spline who contains the reference hexagonal mesh
[in]x1real containing first coordinate of point
[in]x2real containing second coordinate of point
Returns
the first coordinate after change of coordinate sys.

Definition at line 452 of file sll_m_box_splines.F90.

Here is the caller graph for this function:

◆ change_basis_x2()

real(kind=f64) function sll_m_box_splines::change_basis_x2 ( type(sll_t_box_spline_2d), pointer  spline,
real(kind=f64), intent(in)  x1,
real(kind=f64), intent(in)  x2 
)
private

2nd coo. of (x1, x2) in reference hex-mesh coo.

This function allows to change a point of coordinates (x1, x2) on the spline basis to the mesh basis. Gives 2nd coordinate.

Parameters
[in]splinebox spline who contains the reference hexagonal mesh
[in]x1real containing first coordinate of point
[in]x2real containing second coordinate of point
Returns
the second coordinate after change of coordinate sys.

Definition at line 494 of file sll_m_box_splines.F90.

Here is the caller graph for this function:

◆ chi_gen_val()

function sll_m_box_splines::chi_gen_val ( intent(in)  x1_in,
intent(in)  x2_in,
intent(in)  deg 
)
private

Computes the box spline of degree deg on (x1, x2)

Computes the value of the box spline (chi) of degree deg on the point of cartesian coordiantes (x1, x2). The algorithm is specific to the hexagonal mesh and has been optimized for degree 2 splines. Reference : @Condat and Van De Ville (2006) "Three directional Box Splines: Characterization and Efficient Evaluation."

Parameters
[in]x1_inreal containing first coordinate of point where spline is to be evaluated
[in]x2_inreal containing second coordinate of point where spline is to be evaluated
[in]deginteger containing the degree of the spline
Returns
chi_gen_val value of the box spline

Definition at line 298 of file sll_m_box_splines.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ choose()

function sll_m_box_splines::choose ( intent(in)  n,
intent(in)  k 
)
private

Computes the binomial coefficient (n, k)

Computes the binomial coefficient (n, k)

Parameters
[in]ninteger
[in]kinteger
Returns
res binomial coefficient

\[ \dfrac{n!}{(k!(n-k)!} \]

Definition at line 270 of file sll_m_box_splines.F90.

Here is the caller graph for this function:

◆ compute_coeff_box_spline_2d()

subroutine sll_m_box_splines::compute_coeff_box_spline_2d ( class(sll_t_box_spline_2d), intent(inout)  spline,
real(kind=f64), dimension(:), intent(in), target  data,
integer(kind=i32), intent(in)  deg 
)
private

Computes box splines coefficients.

Computes box splines coefficients for a box spline of degree deg and fitted to the data vector

Parameters
[in]datavector containing the data to be fit
[in]deginteger representing the box spline degree
[in]splinebox spline type element, containting the mesh, bc, ...

Definition at line 129 of file sll_m_box_splines.F90.

◆ compute_coeff_box_spline_2d_diri()

subroutine sll_m_box_splines::compute_coeff_box_spline_2d_diri ( class(sll_t_box_spline_2d), intent(inout)  spline,
real(kind=f64), dimension(:), intent(in), target  data,
integer(kind=i32), intent(in)  deg 
)
private

Computes box splines coefficients for dirichlet BC.

Computes box splines coefficients for a box spline of degree deg and fitted to the data vector on a hexmesh with dirichlet BC

Parameters
[in]datavector containing the data to be fit
[in]deginteger representing the box spline degree
[in]splinebox spline type element, containting the mesh, bc, ...

Definition at line 161 of file sll_m_box_splines.F90.

Here is the call graph for this function:

◆ compute_coeff_box_spline_2d_neum()

subroutine sll_m_box_splines::compute_coeff_box_spline_2d_neum ( class(sll_t_box_spline_2d), intent(inout)  spline,
target  data,
intent(in)  deg 
)
private

Computes box splines coefficients with neumann BC.

NOT YET IMPLEMENTED. Computes neuman box splines coefficients for a box spline of degree deg and fitted to the data vector

Parameters
[in]datavector containing the data to be fit
[in]deginteger representing the box spline degree
[in]splinebox spline type element, containting the mesh, bc, ...

Definition at line 247 of file sll_m_box_splines.F90.

◆ compute_coeff_box_spline_2d_prdc()

subroutine sll_m_box_splines::compute_coeff_box_spline_2d_prdc ( class(sll_t_box_spline_2d), intent(inout)  spline,
target  data,
intent(in)  deg 
)
private

Computes box splines coefficients with periodic BC.

NOT YET IMPLEMENTED. Computes periodic box splines coefficients for a box spline of degree deg and fitted to the data vector

Parameters
[in]datavector containing the data to be fit
[in]deginteger representing the box spline degree
[in]splinebox spline type element, containting the mesh, bc, ...

Definition at line 223 of file sll_m_box_splines.F90.

◆ delete_box_spline_2d()

subroutine sll_m_box_splines::delete_box_spline_2d ( type(sll_t_box_spline_2d), intent(inout), pointer  spline)
private

Generic sub-routine defined for 2D box spline types. Deallocates the memory associated with the given box spline object.

Parameters
[in,out]spline_object.

Definition at line 866 of file sll_m_box_splines.F90.

◆ non_zeros_splines()

integer(kind=i32) function, dimension(3*deg*deg) sll_m_box_splines::non_zeros_splines ( type(sll_t_hex_mesh_2d), intent(in), pointer  mesh,
integer(kind=i32), intent(in)  cell_index,
integer(kind=i32), intent(in)  deg 
)
private

Computes indices of non null splines on a given cell.

The function returns for a given cell and a certain degree the indices of the splines of that degree that are different to 0 on that cell.

Parameters
[IN]mesh hexagonal mesh, the domain
[IN]cell_index index of the cell where we wish to know the indices of the non vanishing splines
[IN]deg integer of the degree of the splines
[OUT]index_nZ vector of size 3*deg*deg containing the global indices of all non-zero splines

Definition at line 632 of file sll_m_box_splines.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sll_f_boxspline_val_der()

real(kind=f64) function, public sll_m_box_splines::sll_f_boxspline_val_der ( real(kind=f64), intent(in)  x1,
real(kind=f64), intent(in)  x2,
integer(kind=i32), intent(in)  deg,
integer(kind=i32), intent(in)  nderiv1,
integer(kind=i32), intent(in)  nderiv2 
)

Computes the values or derivatives of box splines.

Depending on nderiv1 and nderiv2 will compute box spline or derivative with respect to x (nderiv1 > 0) or/and to y (nderiv2 > 0)

Parameters
[in]x1real containing second coordinate of point
[in]x2real containing second coordinate of point
[in]deginteger with degree of splines
[in]nderiv1integer number of times to derive on the x direction
[in]nderiv2integer number of times to derive on the y direction
Returns
real nderiv-derivatives of boxspline

no derivative to compute

derivative with respect to the second coo

Definition at line 768 of file sll_m_box_splines.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sll_f_boxspline_x1_derivative()

real(kind=f64) function, public sll_m_box_splines::sll_f_boxspline_x1_derivative ( real(kind=f64), intent(in)  x1,
real(kind=f64), intent(in)  x2,
integer(kind=i32), intent(in)  deg 
)

Computes x-derivative on (x,y)

Using the 5 point stencil computes the x-derivative on (x,y)

Parameters
[in]x1real containing second coordinate of point
[in]x2real containing second coordinate of point
[in]deginteger with degree of splines
Returns
real derivative on x1 of box spline

Definition at line 703 of file sll_m_box_splines.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sll_f_boxspline_x2_derivative()

real(kind=f64) function, public sll_m_box_splines::sll_f_boxspline_x2_derivative ( real(kind=f64), intent(in)  x1,
real(kind=f64), intent(in)  x2,
integer(kind=i32), intent(in)  deg 
)

Computes y-derivative on (x,y)

Using the 5 point stencil computes the y-derivative on (x,y)

Parameters
[in]x1real containing second coordinate of point
[in]x2real containing second coordinate of point
[in]deginteger with degree of splines
Returns
real derivative on x2 of box spline

Definition at line 734 of file sll_m_box_splines.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sll_f_compute_box_spline()

real(kind=f64) function, public sll_m_box_splines::sll_f_compute_box_spline ( type(sll_t_box_spline_2d), intent(in), pointer  spline,
real(kind=f64), intent(in)  x1,
real(kind=f64), intent(in)  x2,
integer(kind=i32), intent(in)  deg 
)

Computes the value of a box spline.

This function computes the value of a box spline of degree deg at the point (x1,x2)

Parameters
[in]splinebox spline which contains the reference hexagonal mesh
[in]x1real containing first coordinate of point
[in]x2real containing second coordinate of point
[in]degreal containing the degree of the spline to be computed
Returns
the value of the box spline at (x1,x2)

Definition at line 428 of file sll_m_box_splines.F90.

Here is the call graph for this function:

◆ sll_f_hex_interpolate_value()

real(kind=f64) function, public sll_m_box_splines::sll_f_hex_interpolate_value ( type(sll_t_hex_mesh_2d), intent(in), pointer  mesh_geom,
real(kind=f64), intent(in)  x1,
real(kind=f64), intent(in)  x2,
type(sll_t_box_spline_2d), pointer  spline,
integer(kind=i32), intent(in)  deg 
)

Interpolates on point of cartesian coordinates (x1, x2)

Interpolation with box splines of degree deg on the point of cartesian coordinates (x1, x2) on the mesh mesh_geom

Parameters
[in]mesh_geomgeometric mesh
[in]x1real containing 1st coordinate of point where to interpolate
[in]x2real containing 2nd coordinate of point where to interpolate
[in]splinebox spline with coefficients already pre computed
[in]deginteger with degree of splines
Returns
real Interpolated value

Definition at line 542 of file sll_m_box_splines.F90.

Here is the call graph for this function:

◆ sll_f_new_box_spline_2d()

type(sll_t_box_spline_2d) function, pointer, public sll_m_box_splines::sll_f_new_box_spline_2d ( type(sll_t_hex_mesh_2d), pointer  mesh,
integer(kind=i32), intent(in)  bc_type 
)

Creates a box spline element.

Creates a box spline element to the associated mesh and with the given boundary conditions. Also allocates memory for the coefficients

Parameters
[in]meshHexagonal mesh on which the box spline will be defined
[in]bc_typeInteger defining the boundary condition, usual values are sll_p_periodic, sll_p_dirichlet, ...
Returns
box spline element

Definition at line 103 of file sll_m_box_splines.F90.

Here is the caller graph for this function:

◆ sll_s_write_connectivity()

subroutine, public sll_m_box_splines::sll_s_write_connectivity ( type(sll_t_hex_mesh_2d), pointer  mesh,
integer(kind=i32), intent(in)  deg 
)

Writes connectivity for CAID / DJANGO.

write connectivity info for CAID/DJANGO. This function was intented to couple Pigasus poisson solver to the hex-mesh. Output file : boxsplines_connectivity.txt

Parameters
[in]meshpointer to the hexagonal mesh
[in]deginteger designing box splines degree

Definition at line 819 of file sll_m_box_splines.F90.

Here is the call graph for this function:
Here is the caller graph for this function:
    Report Typos and Errors