Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Module interface to solve Maxwell's equations with coordinate transformation in 3D.
Derived types and interfaces | |
type | sll_t_maxwell_3d_trafo |
Functions/Subroutines | |
subroutine | sll_s_compute_e_from_b_3d_trafo (self, delta_t, field_in, field_out) |
compute E from B using weak Ampere formulation More... | |
subroutine | sll_s_compute_b_from_e_3d_trafo (self, delta_t, field_in, field_out) |
Compute B from E using strong 3D Faraday equation for spline coefficients. More... | |
subroutine | sll_s_compute_curl_part_3d_trafo (self, delta_t, efield, bfield, betar) |
Solve curl part of Maxwell's equations. More... | |
subroutine | sll_s_compute_e_from_rho_3d_trafo (self, field_in, field_out) |
Compute E_i from rho_i integrated over the time interval using weak Poisson's equation. More... | |
subroutine | sll_s_compute_rho_from_e_3d_trafo (self, field_in, field_out) |
compute rho from e using weak Gauss law ( rho = G^T M_1 e ) More... | |
subroutine | sll_s_compute_e_from_j_3d_trafo (self, current, E, component) |
Compute E_i from j_i integrated over the time interval using weak Ampere formulation, delta_t is already included. More... | |
subroutine | sll_s_compute_phi_from_rho_3d_trafo (self, field_in, field_out, efield_dofs) |
Compute phi from rho_i integrated over the time interval. More... | |
subroutine | sll_s_compute_phi_from_j_3d_trafo (self, field_in, field_out, efield_dofs) |
Compute phi from j_i integrated over the time interval, delta_t is already included. More... | |
subroutine | sll_s_compute_rhs_trafo (self, form, component, coefs_dofs, func1, func2, func3) |
Compute the FEM right-hand-side for a given function f and periodic splines of given degree Its components are $\int f(F(\xi)) N_i(F(\xi)) J_F d\xi$ where $N_i$ is the B-spline starting at $\xi_i$. More... | |
subroutine | rhs_zeroform (self, deg, func, coefs_dofs) |
Compute $\int f(F(\xi)) N_i(\xi) J_F(\xi) d\xi$ for scalar function f, where $N_i$ is the B-spline starting at $\xi_i$, replace modulo(a,b) by a-floor(a/b)*b. More... | |
subroutine | rhs_oneform (self, deg, func1, func2, func3, component, coefs_dofs) |
Compute $\int f(F(\xi))\cdot DF^{-T}(\xi) N_i(\xi) J_F(\xi) d\xi$ for vector function f, where $N_i$ is the B-spline starting at $\xi_i$. More... | |
subroutine | rhs_twoform (self, deg, func1, func2, func3, component, coefs_dofs) |
Compute $\int f(F(\xi))\cdot DF(\xi) N_i(\xi) d\xi$ for vector function f, where $N_i$ is the B-spline starting at $\xi_i$. More... | |
subroutine | rhs_threeform (self, deg, func, coefs_dofs) |
Compute $\int f(F(\xi)) N_i(\xi) d\xi$ for scalar function f, where $N_i$ is the B-spline starting at $\xi_i$. More... | |
subroutine | l2projection_3d_trafo (self, form, component, coefs_dofs, func1, func2, func3) |
Compute the L2 projection of a given function f on periodic splines of given degree. More... | |
real(kind=f64) function | l2norm_squared_3d_trafo (self, coefs, form, component) |
Compute square of the L2norm. More... | |
real(kind=f64) function | inner_product_3d_trafo (self, coefs1, coefs2, form, component) |
Compute inner product. More... | |
subroutine | init_3d_trafo (self, domain, n_dofs, s_deg_0, map, mass_tolerance, poisson_tolerance, solver_tolerance, adiabatic_electrons, profile) |
Initialization. More... | |
subroutine | init_from_file_3d_trafo (self, domain, n_dofs, s_deg_0, map, nml_file, adiabatic_electrons, profile) |
Initialization from nml file. More... | |
subroutine | free_3d_trafo (self) |
Finalization. More... | |
subroutine | multiply_g (self, field_in, field_out) |
Multiply by dicrete gradient matrix. More... | |
subroutine | multiply_gt (self, field_in, field_out) |
Multiply by transpose of dicrete gradient matrix. More... | |
subroutine | multiply_c (self, field_in, field_out) |
Multiply by discrete curl matrix. More... | |
subroutine | multiply_ct (self, field_in, field_out) |
Multiply by transpose of discrete curl matrix. More... | |
subroutine | multiply_mass_trafo (self, deg, coefs_in, coefs_out) |
Multiply by the mass matrix. More... | |
subroutine | multiply_mass_inverse_trafo (self, form, coefs_in, coefs_out) |
Multiply by the inverse mass matrix. More... | |
|
private |
Finalization.
self | Maxwell solver class |
Definition at line 1052 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Initialization.
[in,out] | self | Maxwell solver class |
[in] | n_dofs | number of degrees of freedom (here number of cells and grid points) |
[in] | s_deg_0 | highest spline degree |
[in,out] | map | coordinate transformation |
[in] | mass_tolerance | tolerance for mass solver |
[in] | poisson_tolerance | tolerance for Poisson solver |
[in] | solver_tolerance | tolerance for Schur complement solver |
[in] | adiabatic_electrons | flag if adiabatic electrons are used |
[in] | profile | temperature and density profiles |
Definition at line 749 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Initialization from nml file.
[in,out] | self | Maxwell solver class |
[in] | n_dofs | number of degrees of freedom (here number of cells and grid points) |
[in] | s_deg_0 | highest spline degree |
[in,out] | map | coordinate transformation |
[in] | nml_file | nml-file |
[in] | adiabatic_electrons | flag if adiabatic electrons are used |
[in] | profile | temperature and density profiles |
Definition at line 971 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Compute inner product.
self | Maxwell solver class |
coefs1 | Coefficient for each DoF |
coefs2 | Coefficient for each DoF |
form | Specify 0,1,2 or 3-form |
component | Specify the component of the form |
Definition at line 719 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Compute square of the L2norm.
self | Maxwell solver class |
coefs | Coefficient for each DoF |
form | Specify 0,1,2 or 3-form |
component | Specify the component of the form |
Definition at line 706 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Compute the L2 projection of a given function f on periodic splines of given degree.
self | Maxwell solver class | |
[in] | form | Specify if the function is a0,1,2 or 3-form |
[in] | component | Specify the component of the function |
[out] | coefs_dofs | Finite Element right-hand-side |
func1 | Function first component | |
func2 | Function second component | |
func3 | Function third component |
Definition at line 668 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Multiply by discrete curl matrix.
self | Maxwell solver class | |
[out] | field_out | C*field_in |
Definition at line 1096 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Multiply by transpose of discrete curl matrix.
self | Maxwell solver class | |
[out] | field_out | C^T*field_in |
Definition at line 1107 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Multiply by dicrete gradient matrix.
self | Maxwell solver class | |
[out] | field_out | G*field_in |
Definition at line 1074 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Multiply by transpose of dicrete gradient matrix.
self | Maxwell solver class | |
[out] | field_out | G^T*field_in |
Definition at line 1085 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Multiply by the inverse mass matrix.
self | Maxwell solver class | |
[in] | form | form specifies the form (Note: 0 for 0-form, 1 for 1-form, 2 for 2-form, 3 for 0-1-form mix) |
[out] | coefs_out | M^{-1}*coefs_in |
Definition at line 1145 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Multiply by the mass matrix.
self | Maxwell solver class | |
[in] | deg | deg/form specifies if we multiply the mass to a 1- or 2-form or a mix of both |
[out] | coefs_out | M*coefs_in |
Definition at line 1118 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Compute $\int f(F(\xi))\cdot DF^{-T}(\xi) N_i(\xi) J_F(\xi) d\xi$ for vector function f, where $N_i$ is the B-spline starting at $\xi_i$.
self | Maxwell solver class | |
[in] | deg | spline degree |
func1 | Function first component | |
func2 | Function second component | |
func3 | Function third component | |
[in] | component | Specify the component of the function |
[out] | coefs_dofs | Finite Element right-hand-side |
Definition at line 409 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Compute $\int f(F(\xi)) N_i(\xi) d\xi$ for scalar function f, where $N_i$ is the B-spline starting at $\xi_i$.
self | Maxwell solver class | |
[in] | deg | spline degree |
func | Function | |
[out] | coefs_dofs | Finite Element right-hand-side |
Definition at line 589 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Compute $\int f(F(\xi))\cdot DF(\xi) N_i(\xi) d\xi$ for vector function f, where $N_i$ is the B-spline starting at $\xi_i$.
self | Maxwell solver class | |
[in] | deg | spline degree |
func1 | Function first component | |
func2 | Function second component | |
func3 | Function third component | |
[in] | component | Specify the component of the function |
[out] | coefs_dofs | Finite Element right-hand-side |
Definition at line 499 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Compute $\int f(F(\xi)) N_i(\xi) J_F(\xi) d\xi$ for scalar function f, where $N_i$ is the B-spline starting at $\xi_i$, replace modulo(a,b) by a-floor(a/b)*b.
self | Maxwell solver class | |
[in] | deg | spline degree |
func | function | |
[out] | coefs_dofs | Finite Element right-hand-side |
Definition at line 326 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Compute B from E using strong 3D Faraday equation for spline coefficients.
self | Maxwell solver class | |
[in] | delta_t | Time step |
[in] | field_in | E |
[in,out] | field_out | B |
Definition at line 168 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Solve curl part of Maxwell's equations.
self | Maxwell solver class |
Definition at line 182 of file sll_m_maxwell_3d_trafo.F90.
|
private |
compute E from B using weak Ampere formulation
self | Maxwell solver class | |
[in] | delta_t | Time step |
[in] | field_in | B |
[in,out] | field_out | E |
Definition at line 150 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Compute E_i from j_i integrated over the time interval using weak Ampere formulation, delta_t is already included.
self | Maxwell solver class | |
[in] | current | Current integrated over time interval |
[in,out] | e | Updated electric field |
Definition at line 249 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Compute E_i from rho_i integrated over the time interval using weak Poisson's equation.
self | Maxwell solver class | |
[in] | field_in | rho |
[out] | field_out | E |
Definition at line 221 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Compute phi from j_i integrated over the time interval, delta_t is already included.
self | Maxwell solver class | |
[in] | field_in | Current integrated over time interval |
[in,out] | field_out | phi |
[out] | efield_dofs | E |
Definition at line 276 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Compute phi from rho_i integrated over the time interval.
self | Maxwell solver class | |
[in] | field_in | rho |
[in,out] | field_out | phi |
[out] | efield_dofs | E |
Definition at line 262 of file sll_m_maxwell_3d_trafo.F90.
|
private |
compute rho from e using weak Gauss law ( rho = G^T M_1 e )
self | Maxwell solver class | |
[in] | field_in | E |
[out] | field_out | rho |
Definition at line 234 of file sll_m_maxwell_3d_trafo.F90.
|
private |
Compute the FEM right-hand-side for a given function f and periodic splines of given degree Its components are $\int f(F(\xi)) N_i(F(\xi)) J_F d\xi$ where $N_i$ is the B-spline starting at $\xi_i$.
self | Maxwell solver class | |
[in] | form | Specify if the function is a0,1,2 or 3-form |
[in] | component | Specify the component of the function |
[out] | coefs_dofs | Finite Element right-hand-side |
func1 | Function first component | |
func2 | Function second component | |
func3 | Function third component |
Definition at line 294 of file sll_m_maxwell_3d_trafo.F90.