Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Solve Maxwell's equations with boundary conditions in 1D based on spline FEM, version based on sparse matrices.
Derived types and interfaces | |
type | sll_t_maxwell_clamped_1d_fem_sm |
Functions/Subroutines | |
subroutine | sll_s_compute_e_from_b_1d_fem_sm (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_sm (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_1d_fem_sm (self, delta_t, efield, bfield, betar) |
Solve curl part of Maxwell's equations. More... | |
subroutine | sll_s_compute_e_from_rho_1d_fem_sm (self, field_in, field_out) |
compute e from rho using weak Poisson's equation ( rho = G^T M_1 G \phi, e = G \phi ) More... | |
subroutine | compute_rho_from_e_1d_fem_sm (self, field_in, field_out) |
Compute rho from Gauss law for given efield. More... | |
subroutine | compute_e_from_j_1d_fem_sm (self, current, component, E) |
Compute E_i from j_i integrated over the time interval using weak Ampere formulation. More... | |
subroutine | compute_phi_from_rho_1d_fem_sm (self, in, phi, efield) |
For model with adiabatic electrons. More... | |
subroutine | compute_phi_from_j_1d_fem_sm (self, in, phi, efield) |
For model with adiabatic electrons. More... | |
subroutine | sll_s_compute_rhs_fem_sm (self, func, degree, coefs_dofs) |
Compute the FEM right-hand-side for a given function f and clamped 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_1d_fem_sm (self, func, degree, coefs_dofs) |
Compute the L2 projection of a given function f on periodic splines of given degree. More... | |
real(kind=f64) function | l2norm_squared_1d_fem_sm (self, coefs_dofs, degree) |
Compute square of the L2norm. More... | |
real(kind=f64) function | inner_product_1d_fem_sm (self, coefs1_dofs, coefs2_dofs, degree, degree2) |
Compute inner product. More... | |
subroutine | init_1d_fem_sm (self, domain, n_cells, s_deg_0, boundary, mass_tolerance, poisson_tolerance, solver_tolerance) |
Initialization. More... | |
subroutine | init_from_file_1d_fem_sm (self, domain, n_cells, s_deg_0, boundary, nml_file) |
Initialization from nml file. More... | |
subroutine | free_1d_fem_sm (self) |
Finalization. More... | |
subroutine | multiply_g (self, in, out) |
Multiply by dicrete gradient matrix. More... | |
subroutine | multiply_gt (self, in, out) |
Multiply by transpose of dicrete gradient matrix. More... | |
subroutine | multiply_mass_1d_fem_sm (self, in, out, degree) |
Multiply by the mass matrix. More... | |
subroutine | invert_mass_1d_fem_sm (self, in, out, degree) |
Multiply by the inverse mass matrix. More... | |
subroutine | compute_field_energy (self, efield_dofs1, efield_dofs2, bfield_dofs, energy) |
Compute the field energy. More... | |
|
private |
Compute E_i from j_i integrated over the time interval using weak Ampere formulation.
self | Maxwell_Clamped 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 238 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
Compute the field energy.
self | Maxwell_Clamped solver class | |
[in] | efield_dofs1 | Ex |
[in] | efield_dofs2 | Ey |
[in] | bfield_dofs | Bz |
[out] | energy | field energy |
Definition at line 766 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
For model with adiabatic electrons.
self | Maxwell_Clamped solver class | |
[in] | in | Current integrated over time interval |
[out] | efield | E |
Definition at line 277 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
For model with adiabatic electrons.
self | Maxwell_Clamped solver class | |
[in] | in | rho |
[out] | efield | E |
Definition at line 260 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
Compute rho from Gauss law for given efield.
self | Maxwell_Clamped solver class |
Definition at line 225 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
Finalization.
self | Maxwell_Clamped solver class |
Definition at line 679 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
Initialization.
[out] | self | Maxwell_Clamped solver class |
[in] | domain | xmin, xmax |
[in] | n_cells | number of degrees of freedom (here number of cells and grid points) |
[in] | s_deg_0 | highest spline degree |
[in] | boundary | field boundary conditions |
[in] | mass_tolerance | tolerance for mass solver |
[in] | poisson_tolerance | tolerance for Poisson solver |
[in] | solver_tolerance | tolerance for Schur complement solver |
Definition at line 515 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
Initialization from nml file.
[out] | self | Maxwell_Clamped solver class |
[in] | domain | xmin, xmax |
[in] | n_cells | number of degrees of freedom (here number of cells and grid points) |
[in] | s_deg_0 | highest spline degree |
[in] | boundary | field boundary conditions |
nml_file | nml-file |
Definition at line 617 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
Compute inner product.
self | Maxwell_Clamped 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 468 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
Multiply by the inverse mass matrix.
[in,out] | self | Maxwell_Clamped 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 743 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
Compute square of the L2norm.
self | Maxwell_Clamped solver class |
coefs_dofs | Coefficient for each DoF |
degree | Specify the degree of the basis functions |
Definition at line 447 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
Compute the L2 projection of a given function f on periodic splines of given degree.
self | Maxwell_Clamped solver class | |
func | Function | |
[in] | degree | Specify the degree of the basis functions |
Definition at line 423 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
Multiply by dicrete gradient matrix.
[in] | self | Maxwell_Clamped solver class |
[in] | in | field_in |
[out] | out | G*field_in |
Definition at line 699 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
Multiply by transpose of dicrete gradient matrix.
[in] | self | Maxwell_Clamped solver class |
[in] | in | field_in |
[out] | out | G^T*field_in |
Definition at line 710 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
Multiply by the mass matrix.
[in,out] | self | Maxwell_Clamped 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 721 of file sll_m_maxwell_clamped_1d_fem_sm.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_Clamped solver class | |
[in] | delta_t | Time step |
Definition at line 155 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
Solve curl part of Maxwell's equations.
self | Maxwell_Clamped solver class |
Definition at line 169 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
compute Ey from Bz using weak Ampere formulation
self | Maxwell_Clamped solver class | |
[in] | delta_t | Time step |
[in] | field_in | Bz |
[in,out] | field_out | Ey |
Definition at line 136 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
compute e from rho using weak Poisson's equation ( rho = G^T M_1 G \phi, e = G \phi )
self | Maxwell_Clamped solver class | |
[in] | field_in | rho |
[out] | field_out | E |
Definition at line 212 of file sll_m_maxwell_clamped_1d_fem_sm.F90.
|
private |
Compute the FEM right-hand-side for a given function f and clamped 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_Clamped solver class | |
func | Function | |
[in] | degree | Specify the degree of the basis functions |
[out] | coefs_dofs | Finite Element right-hand-side |
Definition at line 300 of file sll_m_maxwell_clamped_1d_fem_sm.F90.