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

Description

Utilites for 3D Maxwell solvers with spline finite elements.

Author
Benedikt Perse

Derived types and interfaces

interface  sll_i_profile_function
 

Functions/Subroutines

subroutine, public sll_s_spline_fem_mass3d (n_cells, deg, component, matrix, profile, n_cells_min, n_cells_max)
 Set up 3d mass matrix for specific spline degrees and profile function. More...
 
subroutine, public sll_s_spline_fem_mass_line (q, deg, profile, row, n_cells, component, mass_line, xw_gauss_d1, xw_gauss_d2, xw_gauss_d3, bspl_d1, bspl_d2, bspl_d3)
 Computes the mass line for a mass matrix with degree splines. More...
 
subroutine, public sll_s_spline_fem_mixedmass3d (n_cells, deg1, deg2, component, matrix, profile, n_cells_min, n_cells_max)
 Set up 3d mixed mass matrix for specific spline degree and profile function. More...
 
subroutine, public sll_s_spline_fem_mixedmass_line (q, deg1, deg2, profile, row, n_cells, component, mass_line, xw_gauss_d1, xw_gauss_d2, xw_gauss_d3, bspl_d1a, bspl_d2a, bspl_d3a, bspl_d1b, bspl_d2b, bspl_d3b)
 

Function/Subroutine Documentation

◆ sll_s_spline_fem_mass3d()

subroutine, public sll_m_spline_fem_utilities_3d::sll_s_spline_fem_mass3d ( integer(kind=i32), dimension(3), intent(in)  n_cells,
integer(kind=i32), dimension(3), intent(in)  deg,
integer(kind=i32), intent(in)  component,
type(sll_t_matrix_csr), intent(out)  matrix,
procedure(sll_i_profile_function profile,
integer(kind=i32), dimension(3), intent(in), optional  n_cells_min,
integer(kind=i32), dimension(3), intent(in), optional  n_cells_max 
)

Set up 3d mass matrix for specific spline degrees and profile function.

Parameters
[in]n_cellsnumber of cells (and grid points)
[in]deghighest spline degree
[in]componentSpecify the component
[out]matrixsparse mass matrix
profileprofile function
[in]n_cells_minminimal cell number for mpi process
[in]n_cells_maxmaximal cell number for mpi process

Definition at line 55 of file sll_m_spline_fem_utilities_3d.F90.

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

◆ sll_s_spline_fem_mass_line()

subroutine, public sll_m_spline_fem_utilities_3d::sll_s_spline_fem_mass_line ( integer(kind=i32), dimension(3), intent(in)  q,
integer(kind=i32), dimension(3), intent(in)  deg,
procedure(sll_i_profile_function profile,
integer(kind=i32), dimension(3), intent(in)  row,
integer(kind=i32), dimension(3), intent(in)  n_cells,
integer(kind=i32), intent(in)  component,
real(kind=f64), dimension((2*deg(1)+1)*(2*deg(2)+1)*(2*deg(3)+1)), intent(out)  mass_line,
real(kind=f64), dimension(2, q(1)), intent(in)  xw_gauss_d1,
real(kind=f64), dimension(2, q(2)), intent(in)  xw_gauss_d2,
real(kind=f64), dimension(2, q(3)), intent(in)  xw_gauss_d3,
real(kind=f64), dimension(deg(1)+1, q(1)), intent(in)  bspl_d1,
real(kind=f64), dimension(deg(2)+1, q(2)), intent(in)  bspl_d2,
real(kind=f64), dimension(deg(3)+1, q(3)), intent(in)  bspl_d3 
)

Computes the mass line for a mass matrix with degree splines.

Parameters
[in]qnumber of quadrature points
[in]degspline degree in every direction
profileprofile function
[in]rowcurrent row in the matrix
[in]n_cellsgridpoints in every direction
[in]componentspecify component of the form for which the massline is computed
[out]mass_linemassline for one row of the sparse mass matrix
[in]xw_gauss_d3quadrature gauss points
[in]bspl_d3spline values

Definition at line 129 of file sll_m_spline_fem_utilities_3d.F90.

Here is the caller graph for this function:

◆ sll_s_spline_fem_mixedmass3d()

subroutine, public sll_m_spline_fem_utilities_3d::sll_s_spline_fem_mixedmass3d ( integer(kind=i32), dimension(3), intent(in)  n_cells,
integer(kind=i32), dimension(3), intent(in)  deg1,
integer(kind=i32), dimension(3), intent(in)  deg2,
integer(kind=i32), dimension(2), intent(in)  component,
type(sll_t_matrix_csr), intent(out)  matrix,
procedure(sll_i_profile_function profile,
integer(kind=i32), dimension(3), intent(in), optional  n_cells_min,
integer(kind=i32), dimension(3), intent(in), optional  n_cells_max 
)

Set up 3d mixed mass matrix for specific spline degree and profile function.

Parameters
[in]n_cellsnumber of cells (and grid points)
[in]deg2spline degrees
[in]componentSpecify the component
profileprofile function
[out]matrixsparse mass matrix
[in]n_cells_minminimal cell number for mpi process
[in]n_cells_maxmaximal cell number for mpi process

Definition at line 239 of file sll_m_spline_fem_utilities_3d.F90.

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

◆ sll_s_spline_fem_mixedmass_line()

subroutine, public sll_m_spline_fem_utilities_3d::sll_s_spline_fem_mixedmass_line ( integer(kind=i32), dimension(3), intent(in)  q,
integer(kind=i32), dimension(3), intent(in)  deg1,
integer(kind=i32), dimension(3), intent(in)  deg2,
procedure(sll_i_profile_function profile,
integer(kind=i32), dimension(3), intent(in)  row,
integer(kind=i32), dimension(3), intent(in)  n_cells,
integer(kind=i32), dimension(2), intent(in)  component,
real(kind=f64), dimension((deg1(1)+deg2(1)+1)*(deg1(2)+deg2(2)+1)*(deg1(3)+deg2(3)+1)), intent(out)  mass_line,
real(kind=f64), dimension(2, q(1)), intent(in)  xw_gauss_d1,
real(kind=f64), dimension(2, q(2)), intent(in)  xw_gauss_d2,
real(kind=f64), dimension(2, q(3)), intent(in)  xw_gauss_d3,
real(kind=f64), dimension(deg1(1)+1, q(1)), intent(in)  bspl_d1a,
real(kind=f64), dimension(deg1(2)+1, q(2)), intent(in)  bspl_d2a,
real(kind=f64), dimension(deg1(3)+1, q(3)), intent(in)  bspl_d3a,
real(kind=f64), dimension(deg2(1)+1, q(1)), intent(in)  bspl_d1b,
real(kind=f64), dimension(deg2(2)+1, q(2)), intent(in)  bspl_d2b,
real(kind=f64), dimension(deg2(3)+1, q(3)), intent(in)  bspl_d3b 
)
Parameters
[in]qnumber of quadrature points
[in]deg2spline degrees in every direction
profileprofile function
[in]rowcurrent row in the mass matrix
[in]n_cellsgridpoints in every direction
[in]componentspecify component of the form for which the massline is computed
[out]mass_linemassline for one row of the sparse mass matrix
[in]xw_gauss_d3quadrature gauss points
[in]bspl_d3aspline values
[in]bspl_d3bspline values

Definition at line 333 of file sll_m_spline_fem_utilities_3d.F90.

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