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_1d_base Module Reference

Description

Module interface to solve Maxwell's equations in 1D.

Contains the abstract class to create a Maxwell solver in 1D.

Author
Katharina Kormann

Derived types and interfaces

type  sll_c_maxwell_1d_base
 
interface  empty
 
interface  norm_squared
 
interface  signature_mass
 
interface  signature_inner_product
 
interface  sll_i_function_1d_real64
 1d real function More...
 
interface  update_dofs_function
 
interface  compute_field1_from_field2
 
interface  signature_compute_field_from_field
 
interface  signature_compute_E_from_j_1d
 
interface  signature_compute_phi_from_field
 
interface  multiply_derivative
 

Functions/Subroutines

subroutine, public sll_s_plot_two_fields_1d (fname, n1, f1, f2, iplot, time)
 write files to visualize 1d fields with gnuplot More...
 
subroutine transform_dofs (self, in, out, degree)
 Transformation of the dofs (identity by default) More...
 
subroutine compute_curl_part (self, delta_t, efield, bfield, betar)
 Solve source-free Maxwell's equations: Default implementation by solving first B_from_E then E_from_B. More...
 

Function/Subroutine Documentation

◆ compute_curl_part()

subroutine sll_m_maxwell_1d_base::compute_curl_part ( class(sll_c_maxwell_1d_base self,
real(kind=f64), intent(in)  delta_t,
real(kind=f64), dimension(:), intent(inout)  efield,
real(kind=f64), dimension(:), intent(inout)  bfield,
real(kind=f64), optional  betar 
)
private

Solve source-free Maxwell's equations: Default implementation by solving first B_from_E then E_from_B.

Parameters
[in]delta_tTime step
[in,out]efieldEy
[in,out]bfieldBz
betar1/beta

Definition at line 268 of file sll_m_maxwell_1d_base.F90.

◆ sll_s_plot_two_fields_1d()

subroutine, public sll_m_maxwell_1d_base::sll_s_plot_two_fields_1d ( character(len=*), intent(in)  fname,
intent(in)  n1,
intent(in)  f1,
intent(in)  f2,
intent(in)  iplot,
intent(in)  time 
)

write files to visualize 1d fields with gnuplot

Parameters
[in]fnameoutput file name

Definition at line 216 of file sll_m_maxwell_1d_base.F90.

Here is the call graph for this function:

◆ transform_dofs()

subroutine sll_m_maxwell_1d_base::transform_dofs ( class(sll_c_maxwell_1d_base), intent(inout)  self,
dimension(:), intent(in)  in,
dimension(:), intent(out)  out,
intent(in)  degree 
)
private

Transformation of the dofs (identity by default)

Definition at line 256 of file sll_m_maxwell_1d_base.F90.

    Report Typos and Errors