Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Implements the Maxwell solver in 2D with periodic boundary conditions with PSTD method.
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... | |
|
private |
Solve ampere maxwell equation (ex,ey,hz)
[in,out] | self | maxwell equation |
ex | electric field x | |
ey | electric field y | |
hz | magnetic field z | |
dt | time step | |
jx | current x | |
jy | current y |
Definition at line 353 of file sll_m_maxwell_2d_pstd.F90.
|
private |
Solve ampere maxwell equation (hx,hy,ez)
[in,out] | self | maxwell object |
hx | magnetic field x | |
hy | magnetic field y | |
ez | electric field z | |
dt | time step | |
jz | current z |
Definition at line 314 of file sll_m_maxwell_2d_pstd.F90.
|
private |
Solve faraday equation (ex,ey,hz)
[in,out] | self | maxwell object |
[in,out] | ex | electric field x |
[in,out] | ey | electric field y |
[in,out] | hz | magnetic field z |
[in] | dt | time step |
Definition at line 280 of file sll_m_maxwell_2d_pstd.F90.
|
private |
Solve faraday equation (hx,hy,ez)
[in,out] | self | Maxwell object |
[in,out] | hx | Magnetic field x |
[in,out] | hy | Magnetic field y |
[in,out] | ez | Electric field z |
[in] | dt | time step |
Definition at line 246 of file sll_m_maxwell_2d_pstd.F90.
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.
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
self | maxwell object | |
[in] | xmin | x min |
[in] | xmax | x max |
[in] | ymin | y min |
[in] | ymax | y max |
[in] | nc_x | x cells number |
[in] | nc_y | y cells number |
[in] | polarization | TE or TM |
Definition at line 128 of file sll_m_maxwell_2d_pstd.F90.
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 | ||
) |
[in,out] | self | Maxwell object |
[in,out] | fx | field x |
[in,out] | fy | field y |
[in,out] | fz | field z |
[in] | dt | time step |
sx | source x | |
sy | source y |
Definition at line 226 of file sll_m_maxwell_2d_pstd.F90.
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.
[in,out] | self | Maxwell object |
[in,out] | fx | field x |
[in,out] | fy | field y |
[in,out] | fz | field z |
[in] | dt | time step |
Definition at line 209 of file sll_m_maxwell_2d_pstd.F90.
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.
[in,out] | self | maxwell object |
[in,out] | fx | Ex or Bx |
[in,out] | fy | Ey or By |
[in,out] | fz | Bz or Ez |
[in] | dt | time step |
Definition at line 186 of file sll_m_maxwell_2d_pstd.F90.