Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Module interface to solve Maxwell's equations in 3D The linear systems are solved based on FFT diagnoalization.
Derived types and interfaces | |
type | sll_t_maxwell_3d_fem_fft |
Functions/Subroutines | |
subroutine | sll_s_compute_e_from_b_3d_fem_fft (self, delta_t, field_in, field_out) |
compute Ey from Bz using weak Ampere formulation More... | |
subroutine | sll_s_compute_b_from_e_3d_fem_fft (self, delta_t, field_in, field_out) |
Compute Bz from Ey using strong 1D Faraday equation for spline coefficients $B_z^{new}(x_j) = B_z^{old}(x_j) - \frac{\Delta t}{\Delta x} (E_y(x_j) - E_y(x_{j-1}) $. More... | |
subroutine | sll_s_compute_curl_part_3d_fem_fft (self, delta_t, efield, bfield, betar) |
Solve curl part of Maxwell's equations. More... | |
subroutine | sll_s_compute_e_from_rho_3d_fem_fft (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_fem_fft (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_fem_fft (self, current, E, component) |
Compute E_i from j_i integrated over the time interval using weak Ampere formulation. More... | |
subroutine | sll_s_compute_phi_from_rho_3d_fem_fft (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_fem_fft (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_fem_fft (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 N_i dx$ where $N_i$ is the B-spline starting at $x_i$. More... | |
subroutine | l2projection_3d_fem_fft (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_fem_fft (self, coefs, form, component) |
Compute square of the L2norm. More... | |
real(kind=f64) function | inner_product_3d_fem_fft (self, coefs1, coefs2, form, component) |
Compute inner product. More... | |
subroutine | init_3d_fem_fft (self, domain, n_dofs, s_deg_0, adiabatic_electrons, profile) |
Initialization. More... | |
subroutine | free_3d_fem_fft (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_3d_fem_fft (self, deg, coefs_in, coefs_out) |
Multiply by the mass matrix. More... | |
subroutine | multiply_mass_1form (self, coefs_in, coefs_out) |
Helper function for multiply_mass. More... | |
subroutine | multiply_mass_2form (self, coefs_in, coefs_out) |
Helper function for multiply_mass. More... | |
subroutine | multiply_mass_3dkron (self, mass_line_1, mass_line_2, mass_line_3, coefs_in, coefs_out) |
Multiply by the mass matrix. More... | |
subroutine | multiply_mass_all (self, deg, coefs_in, coefs_out) |
Multiply by the mass matrix. More... | |
subroutine | multiply_mass_inverse_all (self, form, coefs_in, coefs_out) |
Multiply by the inverse mass matrix. More... | |
|
private |
Finalization.
self | Maxwell solver class |
Definition at line 721 of file sll_m_maxwell_3d_fem_fft.F90.
|
private |
Initialization.
[out] | self | Maxwell solver class |
[in] | domain | xmin, xmax |
[in] | n_dofs | number of degrees of freedom (here number of cells and grid points) |
[in] | s_deg_0 | highest spline degree |
[in] | adiabatic_electrons | flag if adiabatic electrons are used |
[in] | profile | temperature and density profiles |
Definition at line 513 of file sll_m_maxwell_3d_fem_fft.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 |
Definition at line 451 of file sll_m_maxwell_3d_fem_fft.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 438 of file sll_m_maxwell_3d_fem_fft.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 411 of file sll_m_maxwell_3d_fem_fft.F90.
|
private |
Multiply by discrete curl matrix.
self | Maxwell solver class | |
[out] | field_out | C*field_in |
Definition at line 777 of file sll_m_maxwell_3d_fem_fft.F90.
|
private |
Multiply by transpose of discrete curl matrix.
self | Maxwell solver class | |
[out] | field_out | C^T*field_in |
Definition at line 788 of file sll_m_maxwell_3d_fem_fft.F90.
|
private |
Multiply by dicrete gradient matrix.
self | Maxwell solver class | |
[out] | field_out | G*field_in |
Definition at line 755 of file sll_m_maxwell_3d_fem_fft.F90.
|
private |
Multiply by transpose of dicrete gradient matrix.
self | Maxwell solver class | |
[out] | field_out | G^T*field_in |
Definition at line 766 of file sll_m_maxwell_3d_fem_fft.F90.
|
private |
Helper function for multiply_mass.
[in,out] | self | Maxwell solver class |
[in] | coefs_in | Coefficient for each DoF |
[out] | coefs_out | Coefficient for each DoF |
Definition at line 826 of file sll_m_maxwell_3d_fem_fft.F90.
|
private |
Helper function for multiply_mass.
[in,out] | self | Maxwell solver class |
[in] | coefs_in | Coefficient for each DoF |
[out] | coefs_out | Coefficient for each DoF |
Definition at line 853 of file sll_m_maxwell_3d_fem_fft.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 |
[in] | coefs_in | Coefficient for each DoF |
[out] | coefs_out | Coefficient for each DoF |
Definition at line 799 of file sll_m_maxwell_3d_fem_fft.F90.
|
private |
Multiply by the mass matrix.
[in,out] | self | Maxwell solver class |
[in] | mass_line_1 | massline |
[in] | mass_line_2 | massline |
[in] | mass_line_3 | massline |
[in] | coefs_in | Coefficient for each DoF |
[out] | coefs_out | Coefficient for each DoF |
Definition at line 880 of file sll_m_maxwell_3d_fem_fft.F90.
|
private |
Multiply by the mass matrix.
self | Maxwell solver class | |
[in] | deg | deg(i) specifies the degree of the 1d mass matrix in dimension i (Note: 1 for 0-form, 2 for 1-form, 3 for 0-1-form mix) |
[in] | coefs_in | Coefficient for each DoF |
[out] | coefs_out | Coefficient for each DoF |
Definition at line 939 of file sll_m_maxwell_3d_fem_fft.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) |
[in] | coefs_in | Coefficient for each DoF |
[out] | coefs_out | Coefficient for each DoF |
Definition at line 1020 of file sll_m_maxwell_3d_fem_fft.F90.
|
private |
Compute Bz from Ey using strong 1D Faraday equation for spline coefficients $B_z^{new}(x_j) = B_z^{old}(x_j) - \frac{\Delta t}{\Delta x} (E_y(x_j) - E_y(x_{j-1}) $.
self | Maxwell solver class | |
[in] | delta_t | time step |
[in] | field_in | E |
[in,out] | field_out | B |
Definition at line 179 of file sll_m_maxwell_3d_fem_fft.F90.
|
private |
Solve curl part of Maxwell's equations.
self | Maxwell solver class |
Definition at line 193 of file sll_m_maxwell_3d_fem_fft.F90.
|
private |
compute Ey from Bz 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 153 of file sll_m_maxwell_3d_fem_fft.F90.
|
private |
Compute E_i from j_i integrated over the time interval using weak Ampere formulation.
self | Maxwell solver class | |
[in] | current | Current integrated over time interval |
[in,out] | e | Updated electric field |
[in] | component | component of the Efield to be computed |
Definition at line 256 of file sll_m_maxwell_3d_fem_fft.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 231 of file sll_m_maxwell_3d_fem_fft.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 290 of file sll_m_maxwell_3d_fem_fft.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 276 of file sll_m_maxwell_3d_fem_fft.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 243 of file sll_m_maxwell_3d_fem_fft.F90.
|
private |
Compute the FEM right-hand-side for a given function f and periodic splines of given degree Its components are $\int f N_i dx$ where $N_i$ is the B-spline starting at $x_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 308 of file sll_m_maxwell_3d_fem_fft.F90.