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

Description

Implements the Maxwell solver in 2D with periodic boundary conditions with PSTD method.

Author
Pierre Navaro

Field derivative is made using Fast Fourier Transform.

\( \displaystyle \frac{\partial \psi}{\partial x} \Big|_i = F_x^{-1} [ -jk_xF_x(\psi)]_i, \)

For Maxwell system the scheme is

\( \displaystyle H_u\Big|^{n+1/2}_{i,j,k} = H_u\Big|^{n-1/2}_{i,j,k} - \frac{\Delta t}{\mu} \Big\{F_v^{-1}[-jk_vF_v(E_w)]|_{i,j,k}^n -F_w^{-1}[-jk_wF_w(E_v)]|_{i,j,k}^{n}\Big\}, \)

\( \displaystyle E_u\Big|^{n+1}_{i,j,k} = E_u\Big|^{n}_{i,j,k} + \frac{\Delta t}{\varepsilon} \Big\{F_v^{-1}[-jk_vF_v(H_w)]|_{i,j,k}^{n+1/2} -F_w^{-1}[-jk_wF_w(H_v)]|_{i,j,k}^{n+1/2}\Big\}, \)

where \((u,v,w) = (x,y,z),(y,z,x),(z,x,y)\)

Derived types and interfaces

type  sll_t_maxwell_2d_pstd
 Maxwell solver object. More...
 

Functions/Subroutines

subroutine, public sll_s_maxwell_2d_pstd_init (self, xmin, xmax, nc_x, ymin, ymax, nc_y, polarization)
 Solve maxwell solver 2d cartesian periodic with PSTD scheme. More...
 
subroutine, public sll_s_solve_maxwell_2d_pstd (self, fx, fy, fz, dt)
 self routine exists only for testing purpose. Use ampere and faraday in your appication. More...
 
subroutine, public sll_s_solve_faraday_2d_pstd (self, fx, fy, fz, dt)
 Solve Faraday equation. More...
 
subroutine, public sll_s_solve_ampere_2d_pstd (self, fx, fy, fz, dt, sx, sy)
 
subroutine faraday_tm_2d_pstd (self, hx, hy, ez, dt)
 Solve faraday equation (hx,hy,ez) More...
 
subroutine faraday_te_2d_pstd (self, ex, ey, hz, dt)
 Solve faraday equation (ex,ey,hz) More...
 
subroutine ampere_tm_2d_pstd (self, hx, hy, ez, dt, jz)
 Solve ampere maxwell equation (hx,hy,ez) More...
 
subroutine ampere_te_2d_pstd (self, ex, ey, hz, dt, jx, jy)
 Solve ampere maxwell equation (ex,ey,hz) More...
 
subroutine, public sll_s_free_maxwell_2d_pstd (self)
 delete maxwell solver object More...
 

Function/Subroutine Documentation

◆ ampere_te_2d_pstd()

subroutine sll_m_maxwell_2d_pstd::ampere_te_2d_pstd ( type(sll_t_maxwell_2d_pstd), intent(inout)  self,
real(kind=f64), dimension(:, :)  ex,
real(kind=f64), dimension(:, :)  ey,
real(kind=f64), dimension(:, :)  hz,
real(kind=f64)  dt,
real(kind=f64), dimension(:, :), optional  jx,
real(kind=f64), dimension(:, :), optional  jy 
)
private

Solve ampere maxwell equation (ex,ey,hz)

Parameters
[in,out]selfmaxwell equation
exelectric field x
eyelectric field y
hzmagnetic field z
dttime step
jxcurrent x
jycurrent y

Definition at line 353 of file sll_m_maxwell_2d_pstd.F90.

Here is the caller graph for this function:

◆ ampere_tm_2d_pstd()

subroutine sll_m_maxwell_2d_pstd::ampere_tm_2d_pstd ( type(sll_t_maxwell_2d_pstd), intent(inout)  self,
real(kind=f64), dimension(:, :)  hx,
real(kind=f64), dimension(:, :)  hy,
real(kind=f64), dimension(:, :)  ez,
real(kind=f64)  dt,
real(kind=f64), dimension(:, :), optional  jz 
)
private

Solve ampere maxwell equation (hx,hy,ez)

Parameters
[in,out]selfmaxwell object
hxmagnetic field x
hymagnetic field y
ezelectric field z
dttime step
jzcurrent z

Definition at line 314 of file sll_m_maxwell_2d_pstd.F90.

Here is the caller graph for this function:

◆ faraday_te_2d_pstd()

subroutine sll_m_maxwell_2d_pstd::faraday_te_2d_pstd ( type(sll_t_maxwell_2d_pstd), intent(inout)  self,
real(kind=f64), dimension(:, :), intent(inout)  ex,
real(kind=f64), dimension(:, :), intent(inout)  ey,
real(kind=f64), dimension(:, :), intent(inout)  hz,
real(kind=f64), intent(in)  dt 
)
private

Solve faraday equation (ex,ey,hz)

Parameters
[in,out]selfmaxwell object
[in,out]exelectric field x
[in,out]eyelectric field y
[in,out]hzmagnetic field z
[in]dttime step

Definition at line 280 of file sll_m_maxwell_2d_pstd.F90.

Here is the caller graph for this function:

◆ faraday_tm_2d_pstd()

subroutine sll_m_maxwell_2d_pstd::faraday_tm_2d_pstd ( type(sll_t_maxwell_2d_pstd), intent(inout)  self,
real(kind=f64), dimension(:, :), intent(inout)  hx,
real(kind=f64), dimension(:, :), intent(inout)  hy,
real(kind=f64), dimension(:, :), intent(inout)  ez,
real(kind=f64), intent(in)  dt 
)
private

Solve faraday equation (hx,hy,ez)

Parameters
[in,out]selfMaxwell object
[in,out]hxMagnetic field x
[in,out]hyMagnetic field y
[in,out]ezElectric field z
[in]dttime step

Definition at line 246 of file sll_m_maxwell_2d_pstd.F90.

Here is the caller graph for this function:

◆ sll_s_free_maxwell_2d_pstd()

subroutine, public sll_m_maxwell_2d_pstd::sll_s_free_maxwell_2d_pstd ( type(sll_t_maxwell_2d_pstd self)

delete maxwell solver object

Definition at line 394 of file sll_m_maxwell_2d_pstd.F90.

Here is the call graph for this function:

◆ sll_s_maxwell_2d_pstd_init()

subroutine, public sll_m_maxwell_2d_pstd::sll_s_maxwell_2d_pstd_init ( type(sll_t_maxwell_2d_pstd self,
real(kind=f64), intent(in)  xmin,
real(kind=f64), intent(in)  xmax,
integer(kind=i32), intent(in)  nc_x,
real(kind=f64), intent(in)  ymin,
real(kind=f64), intent(in)  ymax,
integer(kind=i32), intent(in)  nc_y,
integer(kind=i32), intent(in)  polarization 
)

Solve maxwell solver 2d cartesian periodic with PSTD scheme.

Solve ampere equation using maxwell solver 2d cartesian periodic with PSTD scheme Solve faraday equation using solver 2d cartesian periodic with PSTD scheme Delete maxwell solver 2d cartesian periodic with PSTD scheme Initialize 2d maxwell solver on cartesian mesh with PSTD scheme

Parameters
selfmaxwell object
[in]xminx min
[in]xmaxx max
[in]yminy min
[in]ymaxy max
[in]nc_xx cells number
[in]nc_yy cells number
[in]polarizationTE or TM

Definition at line 128 of file sll_m_maxwell_2d_pstd.F90.

Here is the call graph for this function:

◆ sll_s_solve_ampere_2d_pstd()

subroutine, public sll_m_maxwell_2d_pstd::sll_s_solve_ampere_2d_pstd ( type(sll_t_maxwell_2d_pstd), intent(inout)  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,
real(kind=f64), dimension(:, :), optional  sx,
real(kind=f64), dimension(:, :), optional  sy 
)
Parameters
[in,out]selfMaxwell object
[in,out]fxfield x
[in,out]fyfield y
[in,out]fzfield z
[in]dttime step
sxsource x
sysource y

Definition at line 226 of file sll_m_maxwell_2d_pstd.F90.

Here is the call graph for this function:

◆ sll_s_solve_faraday_2d_pstd()

subroutine, public sll_m_maxwell_2d_pstd::sll_s_solve_faraday_2d_pstd ( type(sll_t_maxwell_2d_pstd), intent(inout)  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 
)

Solve Faraday equation.

Parameters
[in,out]selfMaxwell object
[in,out]fxfield x
[in,out]fyfield y
[in,out]fzfield z
[in]dttime step

Definition at line 209 of file sll_m_maxwell_2d_pstd.F90.

Here is the call graph for this function:

◆ sll_s_solve_maxwell_2d_pstd()

subroutine, public sll_m_maxwell_2d_pstd::sll_s_solve_maxwell_2d_pstd ( type(sll_t_maxwell_2d_pstd), intent(inout)  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 
)

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

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

Definition at line 186 of file sll_m_maxwell_2d_pstd.F90.

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