Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
Data Types | Modules | Functions/Subroutines
sll_m_maxwell_2d_fem_fft.F90 File Reference

Go to the source code of this file.

Data Types

type  sll_t_maxwell_2d_fem_fft
 
interface  sll_i_function_2d_real64
 2d real function More...
 

Modules

module  sll_m_maxwell_2d_fem_fft
 Module interface to solve Maxwell's equations in 2D The linear systems are solved based on FFT diagnoalization.
 

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...
 
    Report Typos and Errors