Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Solve Maxwell's equations in 1D.
This version corresponds to conforming spline finite elements
Derived types and interfaces | |
type | sll_t_maxwell_1d_fem |
Maxwell solver class. More... | |
Functions/Subroutines | |
subroutine | sll_s_compute_e_from_b_1d_fem (self, delta_t, field_in, field_out) |
compute Ey from Bz using weak Ampere formulation More... | |
subroutine | sll_s_compute_b_from_e_1d_fem (self, delta_t, field_in, field_out) |
subroutine | sll_s_compute_curl_part_1d_fem (self, delta_t, efield, bfield, betar) |
subroutine | weak_curl (self, delta_t, field_in, field_out) |
compute Ey from Bz using weak Ampere formulation More... | |
subroutine | strong_curl (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_e_from_rho_1d_fem (self, field_in, field_out) |
subroutine | compute_rho_from_e_1d_fem (self, field_in, field_out) |
Compute rho from Gauss law for given efield. More... | |
subroutine | choose_interpolation (self, current, component, E) |
Choose between compute_E_from_j_1d_fem and compute_E_from_j_1d_fem_shape. More... | |
subroutine | compute_e_from_j_1d_fem (self, current, component, E) |
Compute E_i from j_i integrated over the time interval using weak Ampere formulation More... | |
subroutine | compute_e_from_j_1d_fem_shape (self, current, component, E) |
subroutine | multiply_interpolation_inverse_transpose (self, in, out) |
subroutine | solve_e_b_1d_fem (self, delta_t, bfield, efield) |
Solves for the efield part in a implicit curl part solve. More... | |
subroutine | solve_circulant (self, eigvals, rhs, res) |
subroutine | compute_phi_from_rho_1d_fem (self, in, phi, efield) |
For model with adiabatic electrons. More... | |
subroutine | compute_phi_from_j_1d_fem (self, in, phi, efield) |
For model with adiabatic electrons. More... | |
subroutine | sll_s_compute_fem_rhs (self, func, degree, 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... | |
real(kind=f64) function | l2norm_squared_1d_fem (self, coefs_dofs, degree) |
Compute square of the L2norm. More... | |
real(kind=f64) function | inner_product_1d_fem (self, coefs1_dofs, coefs2_dofs, degree, degree2) |
subroutine | l2projection_1d_fem (self, func, degree, coefs_dofs) |
Compute the L2 projection of a given function f on periodic splines of given degree. More... | |
subroutine | init_1d_fem (self, domain, n_dofs, s_deg_0, delta_t, strong_ampere, solver_tolerance, force_sign, adiabatic_electrons) |
subroutine | free_1d_fem (self) |
subroutine | multiply_g (self, in, out) |
subroutine | multiply_gt (self, in, out) |
subroutine | multiply_mass_1d_fem (self, in, out, degree) |
subroutine | invert_mass_1d_fem (self, in, out, degree) |
Invert the mass matrix. More... | |
subroutine | transform_dofs_1d_fem (self, in, out, degree) |
Invert the mass matrix. More... | |
|
private |
Choose between compute_E_from_j_1d_fem and compute_E_from_j_1d_fem_shape.
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 326 of file sll_m_maxwell_1d_fem.F90.
|
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 342 of file sll_m_maxwell_1d_fem.F90.
|
private |
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 363 of file sll_m_maxwell_1d_fem.F90.
|
private |
For model with adiabatic electrons.
self | Maxwell solver class | |
[in] | in | Current integrated over time interval |
[out] | efield | E |
Definition at line 455 of file sll_m_maxwell_1d_fem.F90.
|
private |
For model with adiabatic electrons.
self | Maxwell solver class | |
[in] | in | rho |
[out] | efield | E |
Definition at line 438 of file sll_m_maxwell_1d_fem.F90.
|
private |
Compute rho from Gauss law for given efield.
self | Maxwell solver class | |
[out] | field_out | rho |
Definition at line 294 of file sll_m_maxwell_1d_fem.F90.
|
private |
self | Maxwell solver class |
Definition at line 863 of file sll_m_maxwell_1d_fem.F90.
|
private |
[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] | delta_t | Time step |
[in] | strong_ampere | flag to switch between strong and weak Ampere formulation |
[in] | solver_tolerance | tolerance for Schur complement solver |
[in] | force_sign | sign of particle force |
[in] | adiabatic_electrons | flag if adiabatic electrons are used |
Definition at line 615 of file sll_m_maxwell_1d_fem.F90.
|
private |
self | Maxwell solver class |
coefs1_dofs | Coefficient for each DoF |
coefs2_dofs | Coefficient for each DoF |
degree | Specify the degree of the basis functions |
degree2 | Specify the degree of the basis functions |
Definition at line 541 of file sll_m_maxwell_1d_fem.F90.
|
private |
Invert the mass matrix.
[in,out] | self | Maxwell solver class |
[in] | degree | Specify the degree of the basis functions |
[in] | in | spline coefficients of projection |
[out] | out | spline coefficients of projection |
Definition at line 928 of file sll_m_maxwell_1d_fem.F90.
|
private |
Compute square of the L2norm.
self | Maxwell solver class |
coefs_dofs | Coefficient for each DoF |
degree | Specify the degree of the basis functions |
Definition at line 516 of file sll_m_maxwell_1d_fem.F90.
|
private |
Compute the L2 projection of a given function f on periodic splines of given degree.
self | Maxwell solver class | |
func | function | |
[in] | degree | Specify the degree of the basis functions |
Definition at line 590 of file sll_m_maxwell_1d_fem.F90.
|
private |
[in] | self | Maxwell_Clamped solver class |
[in] | in | field_in |
[out] | out | G*field_in |
Definition at line 881 of file sll_m_maxwell_1d_fem.F90.
|
private |
[in] | self | Maxwell_Clamped solver class |
[in] | in | field_in |
[out] | out | G^T*field_in |
Definition at line 891 of file sll_m_maxwell_1d_fem.F90.
|
private |
self | Maxwell solver class | |
[in] | in | input |
[in,out] | out | output |
Definition at line 378 of file sll_m_maxwell_1d_fem.F90.
|
private |
[in,out] | self | Maxwell solver class |
[in] | in | Coefficient for each DoF |
[out] | out | Coefficient for each DoF |
[in] | degree | Specify the degree of the basis functions |
Definition at line 901 of file sll_m_maxwell_1d_fem.F90.
|
private |
self | Maxwell solver class | |
[in] | delta_t | Time step |
[in] | field_in | Ey |
[in,out] | field_out | Bz |
Definition at line 169 of file sll_m_maxwell_1d_fem.F90.
|
private |
self | Maxwell solver class | |
[in] | delta_t | Time step |
[in,out] | efield | Ey |
[in,out] | bfield | Bz |
betar | 1/beta |
Definition at line 184 of file sll_m_maxwell_1d_fem.F90.
|
private |
compute Ey from Bz using weak Ampere formulation
self | Maxwell solver class | |
[in] | delta_t | Time step |
[in] | field_in | Bz |
[in,out] | field_out | Ey |
Definition at line 153 of file sll_m_maxwell_1d_fem.F90.
|
private |
self | Maxwell solver class | |
[in] | field_in | rho |
[out] | field_out | E |
Definition at line 265 of file sll_m_maxwell_1d_fem.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 | |
func | function | |
[in] | degree | Specify the degree of the basis functions |
Definition at line 477 of file sll_m_maxwell_1d_fem.F90.
|
private |
self | Maxwell solver class | |
[in] | eigvals | eigenvalues of circulant matrix |
[out] | res | result |
Definition at line 405 of file sll_m_maxwell_1d_fem.F90.
|
private |
Solves for the efield part in a implicit curl part solve.
self | Maxwell solver class | |
[in] | delta_t | Time step |
[in] | bfield | Component component of the current integrated over time interval |
[in,out] | efield | Updated electric field |
Definition at line 390 of file sll_m_maxwell_1d_fem.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 | Ey |
[in,out] | field_out | Bz |
Definition at line 246 of file sll_m_maxwell_1d_fem.F90.
|
private |
Invert the mass matrix.
[in,out] | self | Maxwell solver class |
[in] | degree | Specify the degree of the basis functions ! this is 0 for 0 form and 1 for 1 form |
[in] | in | spline coefficients of projection |
[out] | out | spline coefficients of projection |
Definition at line 950 of file sll_m_maxwell_1d_fem.F90.
|
private |
compute Ey from Bz using weak Ampere formulation
self | Maxwell solver class | |
[in] | delta_t | Time step |
[in] | field_in | input field |
[in,out] | field_out | output_field |
Definition at line 226 of file sll_m_maxwell_1d_fem.F90.