Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
Derived types and interfaces | Functions/Subroutines
sll_m_maxwell_2d_fem_fft Module Reference

Description

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

Author
Katharina Kormann

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...
 

Function/Subroutine Documentation

◆ compute_e_from_j_2d_fem()

subroutine sll_m_maxwell_2d_fem_fft::compute_e_from_j_2d_fem ( class(sll_t_maxwell_2d_fem_fft self,
real(kind=f64), dimension(:), intent(in)  current,
integer(kind=i32), intent(in)  component,
real(kind=f64), dimension(:), intent(inout)  E 
)
private

Compute E_i from j_i integrated over the time interval using weak Ampere formulation.

Parameters
selfMaxwell solver class
[in]currentComponent component of the current integrated over time interval
[in]componentComponent of the Efield to be computed
[in,out]eUpdated electric field

Definition at line 178 of file sll_m_maxwell_2d_fem_fft.F90.

◆ compute_rho_from_e_2d_fem()

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.

Here is the call graph for this function:

◆ free_2d_fem()

subroutine sll_m_maxwell_2d_fem_fft::free_2d_fem ( class(sll_t_maxwell_2d_fem_fft self)
private

Definition at line 479 of file sll_m_maxwell_2d_fem_fft.F90.

◆ init_2d_fem()

subroutine sll_m_maxwell_2d_fem_fft::init_2d_fem ( class(sll_t_maxwell_2d_fem_fft), intent(out)  self,
real(kind=f64), dimension(2,2), intent(in)  domain,
integer(kind=i32), dimension(2), intent(in)  n_dofs,
integer(kind=i32), intent(in)  s_deg_0 
)
private
Parameters
[out]selfsolver object

Definition at line 398 of file sll_m_maxwell_2d_fem_fft.F90.

Here is the call graph for this function:

◆ inner_product_2d_fem()

real(kind=f64) function sll_m_maxwell_2d_fem_fft::inner_product_2d_fem ( class(sll_t_maxwell_2d_fem_fft self,
real(kind=f64), dimension(:)  coefs1_dofs,
real(kind=f64), dimension(:)  coefs2_dofs,
integer(kind=i32)  component,
integer(kind=i32)  form 
)
private
Parameters
selfMaxwell solver object
coefs1_dofsCoefficient for each DoF
coefs2_dofsCoefficient for each DoF
componentSpecify the component
formSpecify 0,1,2 or 3 form
Returns
Result: squared L2 norm

Definition at line 337 of file sll_m_maxwell_2d_fem_fft.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ l2norm_squared_2d_fem()

real(kind=f64) function sll_m_maxwell_2d_fem_fft::l2norm_squared_2d_fem ( class(sll_t_maxwell_2d_fem_fft self,
real(kind=f64), dimension(:)  coefs_dofs,
integer(kind=i32)  component,
integer(kind=i32)  form 
)
private

Compute square of the L2norm.

Parameters
selfMaxwell solver object
coefs_dofsCoefficient for each DoF
componentSpecify the component
formSpecify 0,1,2 or 3 form
Returns
Result: squared L2 norm

Definition at line 323 of file sll_m_maxwell_2d_fem_fft.F90.

Here is the call graph for this function:

◆ l2projection_2d_fem()

subroutine sll_m_maxwell_2d_fem_fft::l2projection_2d_fem ( class(sll_t_maxwell_2d_fem_fft self,
procedure(sll_i_function_2d_real64 func,
integer(kind=i32)  component,
integer(kind=i32)  form,
real(kind=f64), dimension(:), intent(out)  coefs_dofs 
)
private

Compute the L2 projection of a given function f on periodic splines of given degree.

Parameters
componentSpecify the component
formSpecify 0,1,2 or 3 form

Definition at line 294 of file sll_m_maxwell_2d_fem_fft.F90.

Here is the call graph for this function:

◆ multiply_c()

subroutine sll_m_maxwell_2d_fem_fft::multiply_c ( class(sll_t_maxwell_2d_fem_fft self,
real(kind=f64), dimension(:), intent(in)  field_in,
real(kind=f64), dimension(:), intent(inout)  field_out 
)
private

Multiply by discrete curl matrix.

Definition at line 490 of file sll_m_maxwell_2d_fem_fft.F90.

Here is the caller graph for this function:

◆ multiply_ct()

subroutine sll_m_maxwell_2d_fem_fft::multiply_ct ( class(sll_t_maxwell_2d_fem_fft self,
real(kind=f64), dimension(:), intent(in)  field_in,
real(kind=f64), dimension(:), intent(inout)  field_out 
)
private

Multiply by transpose of discrete curl matrix.

Parameters
selfMaxwell solver object
[in]field_inMatrix to be multiplied
[in,out]field_outC*field_in

Definition at line 601 of file sll_m_maxwell_2d_fem_fft.F90.

Here is the caller graph for this function:

◆ multiply_gt()

subroutine sll_m_maxwell_2d_fem_fft::multiply_gt ( class(sll_t_maxwell_2d_fem_fft self,
real(kind=f64), dimension(:), intent(in)  field_in,
real(kind=f64), dimension(:), intent(inout)  field_out 
)
private

Multiply by transpose of dicrete gradient matrix.

Parameters
selfMaxwell solver object
[in]field_inMatrix to be multiplied
[in,out]field_outC*field_in

Definition at line 712 of file sll_m_maxwell_2d_fem_fft.F90.

Here is the caller graph for this function:

◆ multiply_mass_1form()

subroutine sll_m_maxwell_2d_fem_fft::multiply_mass_1form ( class(sll_t_maxwell_2d_fem_fft), intent(inout)  self,
real(kind=f64), dimension(:), intent(in)  coefs_in,
real(kind=f64), dimension(:), intent(out)  coefs_out 
)
private

Definition at line 768 of file sll_m_maxwell_2d_fem_fft.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ multiply_mass_2dkron()

subroutine sll_m_maxwell_2d_fem_fft::multiply_mass_2dkron ( class(sll_t_maxwell_2d_fem_fft), intent(inout)  self,
real(kind=f64), dimension(:), intent(in)  mass_line_1,
real(kind=f64), dimension(:), intent(in)  mass_line_2,
real(kind=f64), dimension(:), intent(in)  coefs_in,
real(kind=f64), dimension(:), intent(out)  coefs_out 
)
private

Multiply by the mass matrix.

Definition at line 824 of file sll_m_maxwell_2d_fem_fft.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ multiply_mass_2form()

subroutine sll_m_maxwell_2d_fem_fft::multiply_mass_2form ( class(sll_t_maxwell_2d_fem_fft), intent(inout)  self,
real(kind=f64), dimension(:), intent(in)  coefs_in,
real(kind=f64), dimension(:), intent(out)  coefs_out 
)
private

Definition at line 795 of file sll_m_maxwell_2d_fem_fft.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ multiply_mass_all()

subroutine sll_m_maxwell_2d_fem_fft::multiply_mass_all ( class(sll_t_maxwell_2d_fem_fft), intent(inout)  self,
integer(kind=i32), dimension(2), intent(in)  deg,
real(kind=f64), dimension(:), intent(in)  coefs_in,
real(kind=f64), dimension(:), intent(out)  coefs_out 
)
private

Multiply by the mass matrix.

Parameters
[in]degdeg(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.

Here is the call graph for this function:

◆ sll_s_compute_b_from_e_2d_fem()

subroutine sll_m_maxwell_2d_fem_fft::sll_s_compute_b_from_e_2d_fem ( class(sll_t_maxwell_2d_fem_fft self,
real(kind=f64), intent(in)  delta_t,
real(kind=f64), dimension(:), intent(in)  field_in,
real(kind=f64), dimension(:), intent(inout)  field_out 
)
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.

Here is the call graph for this function:

◆ sll_s_compute_e_from_b_2d_fem()

subroutine sll_m_maxwell_2d_fem_fft::sll_s_compute_e_from_b_2d_fem ( class(sll_t_maxwell_2d_fem_fft self,
real(kind=f64), intent(in)  delta_t,
real(kind=f64), dimension(:), intent(in)  field_in,
real(kind=f64), dimension(:), intent(inout)  field_out 
)
private

compute Ey from Bz using weak Ampere formulation

Parameters
[in]delta_tTime step
[in]field_inB
[in,out]field_outE

Definition at line 137 of file sll_m_maxwell_2d_fem_fft.F90.

Here is the call graph for this function:

◆ sll_s_compute_e_from_rho_2d_fem()

subroutine sll_m_maxwell_2d_fem_fft::sll_s_compute_e_from_rho_2d_fem ( class(sll_t_maxwell_2d_fem_fft self,
real(kind=f64), dimension(:), intent(in)  rho,
real(kind=f64), dimension(:), intent(out)  E 
)
private

Definition at line 194 of file sll_m_maxwell_2d_fem_fft.F90.

◆ sll_s_compute_fem_rhs()

subroutine sll_m_maxwell_2d_fem_fft::sll_s_compute_fem_rhs ( class(sll_t_maxwell_2d_fem_fft self,
procedure(sll_i_function_2d_real64 func,
integer(kind=i32)  component,
integer(kind=i32)  form,
real(kind=f64), dimension(:), intent(out)  coefs_dofs 
)
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$.

Parameters
componentSpecify the component
formSpecify 0,1,2 or 3 form

Definition at line 206 of file sll_m_maxwell_2d_fem_fft.F90.

Here is the call graph for this function:
Here is the caller graph for this function:
    Report Typos and Errors