Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
Data Types | Modules | Macros | Functions/Subroutines
sll_m_box_splines.F90 File Reference
#include "sll_memory.h"
#include "sll_working_precision.h"
Include dependency graph for sll_m_box_splines.F90:

Go to the source code of this file.

Data Types

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...
 

Modules

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

Macros

#define MAKE_GET_SLOT_FUNCTION(fname, datatype, slot, ret_type)
 

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...
 

Macro Definition Documentation

◆ MAKE_GET_SLOT_FUNCTION

#define MAKE_GET_SLOT_FUNCTION (   fname,
  datatype,
  slot,
  ret_type 
)
Value:
1  function fname(spline_obj) result(val); \
2  type(datatype), pointer :: spline_obj; \
3  ret_type :: val; \
4  val = spline_obj%slot; \
5  end function fname
    Report Typos and Errors