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

Description

Implements the Maxwell solver in 2D with FDTD method.

Equation

\( \displaystyle \frac{\partial E_x}{\partial t} = {c^2} \frac{\partial B_z}{\partial y} \),

\(\displaystyle \frac{\partial E_y}{\partial t} = -{c^2} \frac{\partial B_z}{\partial x} \)

\(\displaystyle \frac{\partial B_z}{\partial t} = \frac{\partial E_x}{\partial y} - \frac{\partial E_y}{\partial x} \).

FDTD scheme

\(\displaystyle B_{z}^{n+1/2} = B_z^{n-1/2} - \Delta t \big( \frac{\partial E_y^n}{\partial x} - \frac{\partial E_x^n}{\partial y} \big)\)

\(\displaystyle E_x^{n+1} = E_x^{n} + c^2\Delta t \frac{\partial B_z^{n+1/2}}{\partial y} \)

\(\displaystyle E_y^{n+1} = E_y^{n} - c^2\Delta t \frac{\partial B_z^{n+1/2}}{\partial x} \)

Derived types and interfaces

interface  sll_o_create
 Initialize maxwell solver 2d with FDTD scheme. More...
 
interface  sll_o_solve
 Solve maxwell solver 2d with FDTD scheme. More...
 
interface  sll_solve_ampere
 Solve Ampere-Maxwell equation. More...
 
interface  sll_solve_faraday
 Solve Faraday equation. More...
 
type  sll_t_maxwell_2d_fdtd
 Object with data to solve Maxwell equation Maxwell in TE mode: (Ex,Ey,Bz) More...
 

Functions/Subroutines

subroutine initialize_maxwell_2d_fdtd_alt (self, x1, x2, nc_x, y1, y2, nc_y, polarization)
 Initilialize the maxwell solver. More...
 
subroutine initialize_maxwell_2d_fdtd (self, i1, j1, i2, j2, dx, dy, polarization)
 Initilialize the maxwell solver. More...
 
subroutine solve_maxwell_2d_fdtd (self, fx, fy, fz, dt)
 self routine exists only for testing purpose. Use ampere and faraday in your appication. More...
 
subroutine faraday_2d_fdtd (self, fx, fy, fz, dt)
 Solve Faraday equation. More...
 
subroutine ampere_2d_fdtd (self, fx, fy, fz, dt, jx, jy)
 Solve ampere-maxwell equation with FDTD scheme. More...
 
subroutine bc_periodic_2d_fdtd (self, fx, fy, fz, dt)
 Set boundary conditions. More...
 

Function/Subroutine Documentation

◆ ampere_2d_fdtd()

subroutine sll_m_maxwell_2d_fdtd::ampere_2d_fdtd ( type(sll_t_maxwell_2d_fdtd self,
real(kind=f64), dimension(:, :), intent(inout), target  fx,
real(kind=f64), dimension(:, :), intent(inout), target  fy,
real(kind=f64), dimension(:, :), intent(inout), target  fz,
real(kind=f64), intent(in)  dt,
real(kind=f64), dimension(:, :), intent(in), optional  jx,
real(kind=f64), dimension(:, :), intent(in), optional  jy 
)
private

Solve ampere-maxwell equation with FDTD scheme.

Parameters
selfMaxwell object
[in,out]fxEx or Bx
[in,out]fyEy or By
[in,out]fzBz or Ez
[in]jxJx current
[in]jyJy current
[in]dttime step

Definition at line 241 of file sll_m_maxwell_2d_fdtd.F90.

Here is the caller graph for this function:

◆ bc_periodic_2d_fdtd()

subroutine sll_m_maxwell_2d_fdtd::bc_periodic_2d_fdtd ( type(sll_t_maxwell_2d_fdtd self,
real(kind=f64), dimension(:, :), intent(inout)  fx,
real(kind=f64), dimension(:, :), intent(inout)  fy,
real(kind=f64), dimension(:, :), intent(inout)  fz,
real(kind=f64), intent(in)  dt 
)
private

Set boundary conditions.

Parameters
selfmaxwell solver object
[in,out]fxEx or Bx
[in,out]fyEy or By
[in,out]fzBz or Ez
[in]dttime step

Definition at line 328 of file sll_m_maxwell_2d_fdtd.F90.

Here is the caller graph for this function:

◆ faraday_2d_fdtd()

subroutine sll_m_maxwell_2d_fdtd::faraday_2d_fdtd ( type(sll_t_maxwell_2d_fdtd self,
real(kind=f64), dimension(:, :), target  fx,
real(kind=f64), dimension(:, :), target  fy,
real(kind=f64), dimension(:, :), target  fz,
real(kind=f64), intent(in)  dt 
)
private

Solve Faraday equation.

Parameters
selfMaxwell object
fxEx or Bx
fyEy or By
fzBz or Ez
[in]dttime step

Definition at line 174 of file sll_m_maxwell_2d_fdtd.F90.

Here is the caller graph for this function:

◆ initialize_maxwell_2d_fdtd()

subroutine sll_m_maxwell_2d_fdtd::initialize_maxwell_2d_fdtd ( type(sll_t_maxwell_2d_fdtd self,
integer(kind=i32)  i1,
integer(kind=i32)  j1,
integer(kind=i32)  i2,
integer(kind=i32)  j2,
real(kind=f64)  dx,
real(kind=f64)  dy,
integer(kind=i32)  polarization 
)
private

Initilialize the maxwell solver.

Parameters
selfmaxwell solver object
i1first incidice along x
j1last indice along x
i2first indice along y
j2last indice along y
dxsize step along x
dysize step along y
polarizationTE or TM

Definition at line 127 of file sll_m_maxwell_2d_fdtd.F90.

◆ initialize_maxwell_2d_fdtd_alt()

subroutine sll_m_maxwell_2d_fdtd::initialize_maxwell_2d_fdtd_alt ( type(sll_t_maxwell_2d_fdtd self,
real(kind=f64)  x1,
real(kind=f64)  x2,
integer(kind=i32)  nc_x,
real(kind=f64)  y1,
real(kind=f64)  y2,
integer(kind=i32)  nc_y,
integer(kind=i32)  polarization 
)
private

Initilialize the maxwell solver.

Parameters
selfmaxwell solver object
x1first incidice along x
y1last indice along x
x2first indice along y
y2last indice along y
nc_xsize step along y
nc_ysize step along y
polarizationTE or TM

Definition at line 100 of file sll_m_maxwell_2d_fdtd.F90.

◆ solve_maxwell_2d_fdtd()

subroutine sll_m_maxwell_2d_fdtd::solve_maxwell_2d_fdtd ( type(sll_t_maxwell_2d_fdtd self,
real(kind=f64), dimension(:, :)  fx,
real(kind=f64), dimension(:, :)  fy,
real(kind=f64), dimension(:, :)  fz,
real(kind=f64), intent(in)  dt 
)
private

self routine exists only for testing purpose. Use ampere and faraday in your appication.

Parameters
selfmaxwell object
fxEx or Bx
fyEy or By
fzBz or Ez
[in]dttime step

Definition at line 154 of file sll_m_maxwell_2d_fdtd.F90.

    Report Typos and Errors