![]() |
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Module interface to solve Maxwell's equations in 2D The linear systems are solved based on FFT diagnoalization.
Derived types and interfaces | |
| type | sll_t_maxwell_2d_fem_fft |
| interface | sll_i_function_2d_real64 |
| 2d real function More... | |
Functions/Subroutines | |
| subroutine | compute_rho_from_e_2d_fem (self, efield, rho) |
| compute rho from e using weak Gauss law ( rho = G^T M_1 e ) More... | |
| subroutine | sll_s_compute_e_from_b_2d_fem (self, delta_t, field_in, field_out) |
| compute Ey from Bz using weak Ampere formulation More... | |
| subroutine | sll_s_compute_b_from_e_2d_fem (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 | compute_e_from_j_2d_fem (self, current, component, E) |
| Compute E_i from j_i integrated over the time interval using weak Ampere formulation. More... | |
| subroutine | sll_s_compute_e_from_rho_2d_fem (self, rho, E) |
| subroutine | sll_s_compute_fem_rhs (self, func, component, form, coefs_dofs) |
| 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_2d_fem (self, func, component, form, coefs_dofs) |
| Compute the L2 projection of a given function f on periodic splines of given degree. More... | |
| real(kind=f64) function | l2norm_squared_2d_fem (self, coefs_dofs, component, form) |
| Compute square of the L2norm. More... | |
| real(kind=f64) function | inner_product_2d_fem (self, coefs1_dofs, coefs2_dofs, component, form) |
| subroutine | init_2d_fem (self, domain, n_dofs, s_deg_0) |
| subroutine | free_2d_fem (self) |
| 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_gt (self, field_in, field_out) |
| Multiply by transpose of dicrete gradient matrix. More... | |
| subroutine | multiply_mass_1form (self, coefs_in, coefs_out) |
| subroutine | multiply_mass_2form (self, coefs_in, coefs_out) |
| subroutine | multiply_mass_2dkron (self, mass_line_1, mass_line_2, 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... | |
|
private |
Compute E_i from j_i integrated over the time interval using weak Ampere formulation.
| self | Maxwell solver class | |
| [in] | current | Component component of the current integrated over time interval |
| [in] | component | Component of the Efield to be computed |
| [in,out] | e | Updated electric field |
Definition at line 178 of file sll_m_maxwell_2d_fem_fft.F90.
| subroutine sll_m_maxwell_2d_fem_fft::compute_rho_from_e_2d_fem | ( | class(sll_t_maxwell_2d_fem_fft) | self, |
| real(kind=f64), dimension(:), intent(in) | efield, | ||
| real(kind=f64), dimension(:), intent(inout) | rho | ||
| ) |
compute rho from e using weak Gauss law ( rho = G^T M_1 e )
Definition at line 119 of file sll_m_maxwell_2d_fem_fft.F90.
|
private |
Definition at line 479 of file sll_m_maxwell_2d_fem_fft.F90.
|
private |
| [out] | self | solver object |
Definition at line 398 of file sll_m_maxwell_2d_fem_fft.F90.
|
private |
| self | Maxwell solver object |
| coefs1_dofs | Coefficient for each DoF |
| coefs2_dofs | Coefficient for each DoF |
| component | Specify the component |
| form | Specify 0,1,2 or 3 form |
Definition at line 337 of file sll_m_maxwell_2d_fem_fft.F90.
|
private |
Compute square of the L2norm.
| self | Maxwell solver object |
| coefs_dofs | Coefficient for each DoF |
| component | Specify the component |
| form | Specify 0,1,2 or 3 form |
Definition at line 323 of file sll_m_maxwell_2d_fem_fft.F90.
|
private |
Compute the L2 projection of a given function f on periodic splines of given degree.
| component | Specify the component |
| form | Specify 0,1,2 or 3 form |
Definition at line 294 of file sll_m_maxwell_2d_fem_fft.F90.
|
private |
Multiply by discrete curl matrix.
Definition at line 490 of file sll_m_maxwell_2d_fem_fft.F90.
|
private |
Multiply by transpose of discrete curl matrix.
| self | Maxwell solver object | |
| [in] | field_in | Matrix to be multiplied |
| [in,out] | field_out | C*field_in |
Definition at line 601 of file sll_m_maxwell_2d_fem_fft.F90.
|
private |
Multiply by transpose of dicrete gradient matrix.
| self | Maxwell solver object | |
| [in] | field_in | Matrix to be multiplied |
| [in,out] | field_out | C*field_in |
Definition at line 712 of file sll_m_maxwell_2d_fem_fft.F90.
|
private |
Definition at line 768 of file sll_m_maxwell_2d_fem_fft.F90.
|
private |
Multiply by the mass matrix.
Definition at line 824 of file sll_m_maxwell_2d_fem_fft.F90.
|
private |
Definition at line 795 of file sll_m_maxwell_2d_fem_fft.F90.
|
private |
Multiply by the mass matrix.
| [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) |
Definition at line 863 of file sll_m_maxwell_2d_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}) $.
Definition at line 163 of file sll_m_maxwell_2d_fem_fft.F90.
|
private |
compute Ey from Bz using weak Ampere formulation
| [in] | delta_t | Time step |
| [in] | field_in | B |
| [in,out] | field_out | E |
Definition at line 137 of file sll_m_maxwell_2d_fem_fft.F90.
|
private |
Definition at line 194 of file sll_m_maxwell_2d_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$.
| component | Specify the component |
| form | Specify 0,1,2 or 3 form |
Definition at line 206 of file sll_m_maxwell_2d_fem_fft.F90.
1.9.1