Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Provides capabilities for values and derivatives interpolation with box splines on hexagonal meshes.
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... | |
|
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.
[in] | spline | box spline who contains the reference hexagonal mesh |
[in] | x1 | real containing first coordinate of point |
[in] | x2 | real containing second coordinate of point |
Definition at line 452 of file sll_m_box_splines.F90.
|
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.
[in] | spline | box spline who contains the reference hexagonal mesh |
[in] | x1 | real containing first coordinate of point |
[in] | x2 | real containing second coordinate of point |
Definition at line 494 of file sll_m_box_splines.F90.
|
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."
[in] | x1_in | real containing first coordinate of point where spline is to be evaluated |
[in] | x2_in | real containing second coordinate of point where spline is to be evaluated |
[in] | deg | integer containing the degree of the spline |
Definition at line 298 of file sll_m_box_splines.F90.
|
private |
Computes the binomial coefficient (n, k)
Computes the binomial coefficient (n, k)
[in] | n | integer |
[in] | k | integer |
\[ \dfrac{n!}{(k!(n-k)!} \]
Definition at line 270 of file sll_m_box_splines.F90.
|
private |
Computes box splines coefficients.
Computes box splines coefficients for a box spline of degree deg and fitted to the data vector
[in] | data | vector containing the data to be fit |
[in] | deg | integer representing the box spline degree |
[in] | spline | box spline type element, containting the mesh, bc, ... |
Definition at line 129 of file sll_m_box_splines.F90.
|
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
[in] | data | vector containing the data to be fit |
[in] | deg | integer representing the box spline degree |
[in] | spline | box spline type element, containting the mesh, bc, ... |
Definition at line 161 of file sll_m_box_splines.F90.
|
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
[in] | data | vector containing the data to be fit |
[in] | deg | integer representing the box spline degree |
[in] | spline | box spline type element, containting the mesh, bc, ... |
Definition at line 247 of file sll_m_box_splines.F90.
|
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
[in] | data | vector containing the data to be fit |
[in] | deg | integer representing the box spline degree |
[in] | spline | box spline type element, containting the mesh, bc, ... |
Definition at line 223 of file sll_m_box_splines.F90.
|
private |
Generic sub-routine defined for 2D box spline types. Deallocates the memory associated with the given box spline object.
[in,out] | spline_object. |
Definition at line 866 of file sll_m_box_splines.F90.
|
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.
[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.
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)
[in] | x1 | real containing second coordinate of point |
[in] | x2 | real containing second coordinate of point |
[in] | deg | integer with degree of splines |
[in] | nderiv1 | integer number of times to derive on the x direction |
[in] | nderiv2 | integer number of times to derive on the y direction |
no derivative to compute
derivative with respect to the second coo
Definition at line 768 of file sll_m_box_splines.F90.
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)
[in] | x1 | real containing second coordinate of point |
[in] | x2 | real containing second coordinate of point |
[in] | deg | integer with degree of splines |
Definition at line 703 of file sll_m_box_splines.F90.
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)
[in] | x1 | real containing second coordinate of point |
[in] | x2 | real containing second coordinate of point |
[in] | deg | integer with degree of splines |
Definition at line 734 of file sll_m_box_splines.F90.
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)
[in] | spline | box spline which contains the reference hexagonal mesh |
[in] | x1 | real containing first coordinate of point |
[in] | x2 | real containing second coordinate of point |
[in] | deg | real containing the degree of the spline to be computed |
Definition at line 428 of file sll_m_box_splines.F90.
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
[in] | mesh_geom | geometric mesh |
[in] | x1 | real containing 1st coordinate of point where to interpolate |
[in] | x2 | real containing 2nd coordinate of point where to interpolate |
[in] | spline | box spline with coefficients already pre computed |
[in] | deg | integer with degree of splines |
Definition at line 542 of file sll_m_box_splines.F90.
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
[in] | mesh | Hexagonal mesh on which the box spline will be defined |
[in] | bc_type | Integer defining the boundary condition, usual values are sll_p_periodic, sll_p_dirichlet, ... |
Definition at line 103 of file sll_m_box_splines.F90.
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
[in] | mesh | pointer to the hexagonal mesh |
[in] | deg | integer designing box splines degree |
Definition at line 819 of file sll_m_box_splines.F90.