Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Utilites for Maxwell solver's with spline finite elements.
Derived types and interfaces | |
interface | sll_i_profile_function_1d |
Functions/Subroutines | |
subroutine, public | sll_s_spline_fem_mass_line (degree, mass_line) |
Computes the mass line for a mass matrix with degree splines. More... | |
subroutine, public | sll_s_spline_fem_mass_line_boundary (degree, spline, mass_line) |
Function computing the boundary part for a mass matrix with clamped splines of degree. More... | |
subroutine, public | sll_s_spline_fem_mass_line_boundary_full (deg, profile, spline, row, n_cells, mass_line) |
Function computing the boundary part for a mass matrix with clamped splines of degree (full version without symmetry part) More... | |
subroutine, public | sll_s_spline_fem_mass_line_full (deg, profile, mass_line, row, n_cells) |
Helper function to find line of mass matrix (full version without symmetry part) More... | |
subroutine, public | sll_s_spline_fem_mixedmass_line_boundary (deg1, deg2, profile, spline1, spline2, row, n_cells, mass_line) |
Function computing the boundary part for a mass matrix with clamped splines of degree (full version without symmetry part) More... | |
subroutine, public | sll_s_spline_fem_mixedmass_line_full (deg1, deg2, profile, mass_line, cell, n_cells) |
Helper function to find line of mass matrix ( N_i^p N_j^{p-1})) More... | |
subroutine, public | sll_s_spline_fem_mixedmass_line (deg, mass_line) |
Helper function to find line of mass matrix ( N_i^p N_j^{p-1})) More... | |
subroutine, public | sll_s_spline_fem_compute_mass_eig (n_cells, degree, mass_line, eig_values_mass) |
Compute eigenvalues of mass matrix with given line mass_line. More... | |
subroutine, public | sll_s_spline_fem_multiply_mass (n_cells, degree, mass, invec, outvec) |
Multiply the vector invec with the spline FEM mass matrix with mass line mass. More... | |
subroutine, public | sll_s_spline_fem_multiply_massmixed (n_cells, degree, mass, invec, outvec) |
Multiplication of the mix mass matrix given by a mass line mass. More... | |
subroutine, public | sll_s_spline_fem_interpolation_eigenvalues (degree, ndofs, eig) |
Compute the eigenvalues of the interpolation matrix for splines of degree degree (with first grid point as starting point of the first spline) The interpolation is for the Greville points, i.e. the grid points for odd degree and the midpoints for even degree. The grid points are numbered as they appear in the grid and the first mid point is the one in the first cell. More... | |
subroutine, public | sll_s_multiply_g_1d (n_dofs, delta_x, in, out) |
Multiplication of the input vector in by the derivative matrix G. More... | |
subroutine, public | sll_s_multiply_gt_1d (n_dofs, delta_x, in, out) |
Multiplication of the input vector in by the transposed derivative matrix G^T More... | |
subroutine, public | sll_s_multiply_g (n_dofs, delta_x, field_in, field_out) |
Multiply by dicrete gradient matrix (3d version) More... | |
subroutine, public | sll_s_multiply_gt (n_dofs, delta_x, field_in, field_out) |
Multiply by transpose of dicrete gradient matrix (3d version) More... | |
subroutine, public | sll_s_multiply_c (n_dofs, delta_x, field_in, field_out) |
Multiply by discrete curl matrix. More... | |
subroutine, public | sll_s_multiply_ct (n_dofs, delta_x, field_in, field_out) |
Multiply by transpose of discrete curl matrix. More... | |
subroutine, public sll_m_spline_fem_utilities::sll_s_multiply_c | ( | dimension(3), intent(in) | n_dofs, |
dimension(3), intent(in) | delta_x, | ||
dimension(:), intent(in) | field_in, | ||
dimension(:), intent(out) | field_out | ||
) |
Multiply by discrete curl matrix.
Definition at line 927 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_multiply_ct | ( | dimension(3), intent(in) | n_dofs, |
dimension(3), intent(in) | delta_x, | ||
dimension(:), intent(in) | field_in, | ||
dimension(:), intent(out) | field_out | ||
) |
Multiply by transpose of discrete curl matrix.
Definition at line 1057 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_multiply_g | ( | dimension(3), intent(in) | n_dofs, |
dimension(3), intent(in) | delta_x, | ||
dimension(:), intent(in) | field_in, | ||
dimension(:), intent(out) | field_out | ||
) |
Multiply by dicrete gradient matrix (3d version)
Definition at line 771 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_multiply_g_1d | ( | intent(in) | n_dofs, |
intent(in) | delta_x, | ||
dimension(:), intent(in) | in, | ||
dimension(:), intent(out) | out | ||
) |
Multiplication of the input vector in by the derivative matrix G.
Definition at line 733 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_multiply_gt | ( | dimension(3), intent(in) | n_dofs, |
dimension(3), intent(in) | delta_x, | ||
dimension(:), intent(in) | field_in, | ||
dimension(:), intent(out) | field_out | ||
) |
Multiply by transpose of dicrete gradient matrix (3d version)
Definition at line 845 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_multiply_gt_1d | ( | intent(in) | n_dofs, |
intent(in) | delta_x, | ||
dimension(:), intent(in) | in, | ||
dimension(:), intent(out) | out | ||
) |
Multiplication of the input vector in by the transposed derivative matrix G^T
Definition at line 752 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_spline_fem_compute_mass_eig | ( | intent(in) | n_cells, |
intent(in) | degree, | ||
dimension(0:degree), intent(in) | mass_line, | ||
dimension(n_cells), intent(out) | eig_values_mass | ||
) |
Compute eigenvalues of mass matrix with given line mass_line.
Definition at line 571 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_spline_fem_interpolation_eigenvalues | ( | intent(in) | degree, |
intent(in) | ndofs, | ||
dimension(ndofs), intent(out) | eig | ||
) |
Compute the eigenvalues of the interpolation matrix for splines of degree degree (with first grid point as starting point of the first spline) The interpolation is for the Greville points, i.e. the grid points for odd degree and the midpoints for even degree. The grid points are numbered as they appear in the grid and the first mid point is the one in the first cell.
Definition at line 689 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_spline_fem_mass_line | ( | intent(in) | degree, |
dimension(degree+1), intent(out) | mass_line | ||
) |
Computes the mass line for a mass matrix with degree splines.
Definition at line 64 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_spline_fem_mass_line_boundary | ( | intent(in) | degree, |
type(sll_t_spline_pp_1d), intent(in) | spline, | ||
dimension(:), intent(out) | mass_line | ||
) |
Function computing the boundary part for a mass matrix with clamped splines of degree.
[in] | spline | pp spline class |
Definition at line 98 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_spline_fem_mass_line_boundary_full | ( | intent(in) | deg, |
procedure(sll_i_profile_function_1d) | profile, | ||
type(sll_t_spline_pp_1d), intent(in) | spline, | ||
intent(in) | row, | ||
intent(in) | n_cells, | ||
dimension((7*deg**2-deg-2)/2), intent(out) | mass_line | ||
) |
Function computing the boundary part for a mass matrix with clamped splines of degree (full version without symmetry part)
profile | profile function | |
[in] | spline | pp spline class |
Definition at line 234 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_spline_fem_mass_line_full | ( | intent(in) | deg, |
procedure(sll_i_profile_function_1d) | profile, | ||
dimension(2*deg+1), intent(out) | mass_line, | ||
intent(in) | row, | ||
intent(in) | n_cells | ||
) |
Helper function to find line of mass matrix (full version without symmetry part)
profile | profile function |
Definition at line 309 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_spline_fem_mixedmass_line | ( | intent(in) | deg, |
dimension(deg*2), intent(out) | mass_line | ||
) |
Helper function to find line of mass matrix ( N_i^p N_j^{p-1}))
Definition at line 525 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_spline_fem_mixedmass_line_boundary | ( | intent(in) | deg1, |
intent(in) | deg2, | ||
procedure(sll_i_profile_function_1d) | profile, | ||
type(sll_t_spline_pp_1d), intent(in) | spline1, | ||
type(sll_t_spline_pp_1d), intent(in) | spline2, | ||
intent(in) | row, | ||
intent(in) | n_cells, | ||
dimension( (deg1+deg2)**2+(2*deg2+deg1-deg1**2)/2 ), intent(out) | mass_line | ||
) |
Function computing the boundary part for a mass matrix with clamped splines of degree (full version without symmetry part)
profile | profile function | |
[in] | spline1 | 1D pp-spline |
[in] | spline2 | 1D pp-spline |
Definition at line 365 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_spline_fem_mixedmass_line_full | ( | intent(in) | deg1, |
intent(in) | deg2, | ||
procedure(sll_i_profile_function_1d) | profile, | ||
dimension(deg1+deg2+1), intent(out) | mass_line, | ||
intent(in) | cell, | ||
intent(in) | n_cells | ||
) |
Helper function to find line of mass matrix ( N_i^p N_j^{p-1}))
profile | profile function |
Definition at line 465 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_spline_fem_multiply_mass | ( | intent(in) | n_cells, |
intent(in) | degree, | ||
dimension(degree+1), intent(in) | mass, | ||
dimension(:), intent(in) | invec, | ||
dimension(:), intent(out) | outvec | ||
) |
Multiply the vector invec with the spline FEM mass matrix with mass line mass.
Definition at line 593 of file sll_m_spline_fem_utilities.F90.
subroutine, public sll_m_spline_fem_utilities::sll_s_spline_fem_multiply_massmixed | ( | intent(in) | n_cells, |
intent(in) | degree, | ||
dimension(degree*2), intent(in) | mass, | ||
dimension(:), intent(in) | invec, | ||
dimension(:), intent(out) | outvec | ||
) |
Multiplication of the mix mass matrix given by a mass line mass.
Definition at line 643 of file sll_m_spline_fem_utilities.F90.