Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Implements the Maxwell solver in 3D 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_3d_pstd |
Maxwell solver object. More... | |
Functions/Subroutines | |
subroutine | sll_s_maxwell_3d_pstd_init (self, xmin, xmax, nc_x, ymin, ymax, nc_y, zmin, zmax, nc_z) |
Initialize 3d maxwell solver on cartesian mesh with PSTD scheme. More... | |
subroutine | sll_s_maxwell_3d_pstd_faraday (self, ex, ey, ez, hx, hy, hz, dt) |
Solve faraday equation (hx,hy,ez) More... | |
subroutine | sll_s_maxwell_3d_pstd_ampere (self, hx, hy, hz, ex, ey, ez, dt, jx, jy, jz) |
Solve ampere maxwell equation (hx,hy,ez) More... | |
subroutine | sll_s_maxwell_3d_pstd_free (self) |
delete maxwell solver object More... | |
type sll_m_maxwell_3d_pstd::sll_t_maxwell_3d_pstd |
Maxwell solver object.
Definition at line 78 of file sll_m_maxwell_3d_pstd.F90.
Class Members | ||
---|---|---|
type(sll_t_fft) | bwx | backward fft plan along x |
type(sll_t_fft) | bwy | backward fft plan along y |
type(sll_t_fft) | bwz | backward fft plan along y |
real(kind=f64), dimension(:), pointer | d_dx | field x derivative |
real(kind=f64), dimension(:), pointer | d_dy | field y derivative |
real(kind=f64), dimension(:), pointer | d_dz | field y derivative |
real(kind=f64) | e_0 | electric conductivity |
type(sll_t_fft) | fwx | forward fft plan along x |
type(sll_t_fft) | fwy | forward fft plan along y |
type(sll_t_fft) | fwz | forward fft plan along y |
real(kind=f64), dimension(:), pointer | kx | x wave number |
real(kind=f64), dimension(:), pointer | ky | y wave number |
real(kind=f64), dimension(:), pointer | kz | z wave number |
real(kind=f64) | mu_0 | magnetic permeability |
integer(kind=i32) | nc_x | x cells number |
integer(kind=i32) | nc_y | y cells number |
integer(kind=i32) | nc_z | z cells number |
integer(kind=i32) | polarization | TE or TM. |
complex(kind=f64), dimension(:), pointer | tmp_x | x fft transform |
complex(kind=f64), dimension(:), pointer | tmp_y | y fft transform |
complex(kind=f64), dimension(:), pointer | tmp_z | y fft transform |
subroutine sll_m_maxwell_3d_pstd::sll_s_maxwell_3d_pstd_ampere | ( | type(sll_t_maxwell_3d_pstd), intent(inout) | self, |
real(kind=f64), dimension(:, :, :) | hx, | ||
real(kind=f64), dimension(:, :, :) | hy, | ||
real(kind=f64), dimension(:, :, :) | hz, | ||
real(kind=f64), dimension(:, :, :) | ex, | ||
real(kind=f64), dimension(:, :, :) | ey, | ||
real(kind=f64), dimension(:, :, :) | ez, | ||
real(kind=f64) | dt, | ||
real(kind=f64), dimension(:, :, :), optional | jx, | ||
real(kind=f64), dimension(:, :, :), optional | jy, | ||
real(kind=f64), dimension(:, :, :), optional | jz | ||
) |
Solve ampere maxwell equation (hx,hy,ez)
[in,out] | self | maxwell object |
hx | magnetic field x | |
hy | magnetic field y | |
hz | magnetic field z | |
ex | electric field x | |
ey | electric field y | |
ez | electric field z | |
dt | time step | |
jx | current x | |
jy | current y | |
jz | current z |
Definition at line 243 of file sll_m_maxwell_3d_pstd.F90.
subroutine sll_m_maxwell_3d_pstd::sll_s_maxwell_3d_pstd_faraday | ( | type(sll_t_maxwell_3d_pstd), intent(inout) | self, |
real(kind=f64), dimension(:, :, :), intent(inout) | ex, | ||
real(kind=f64), dimension(:, :, :), intent(inout) | ey, | ||
real(kind=f64), dimension(:, :, :), intent(inout) | ez, | ||
real(kind=f64), dimension(:, :, :), intent(inout) | hx, | ||
real(kind=f64), dimension(:, :, :), intent(inout) | hy, | ||
real(kind=f64), dimension(:, :, :), intent(inout) | hz, | ||
real(kind=f64), intent(in) | dt | ||
) |
Solve faraday equation (hx,hy,ez)
[in,out] | self | Maxwell object |
[in,out] | ex | Electric field x |
[in,out] | ey | Electric field y |
[in,out] | ez | Electric field z |
[in,out] | hx | Magnetic field x |
[in,out] | hy | Magnetic field y |
[in,out] | hz | Magnetic field z |
[in] | dt | time step |
Definition at line 181 of file sll_m_maxwell_3d_pstd.F90.
subroutine sll_m_maxwell_3d_pstd::sll_s_maxwell_3d_pstd_free | ( | type(sll_t_maxwell_3d_pstd) | self | ) |
delete maxwell solver object
Definition at line 315 of file sll_m_maxwell_3d_pstd.F90.
subroutine sll_m_maxwell_3d_pstd::sll_s_maxwell_3d_pstd_init | ( | type(sll_t_maxwell_3d_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, | ||
real(kind=f64), intent(in) | zmin, | ||
real(kind=f64), intent(in) | zmax, | ||
integer(kind=i32), intent(in) | nc_z | ||
) |
Initialize 3d 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] | zmin | z min |
[in] | zmax | z max |
[in] | nc_x | x cells number |
[in] | nc_y | y cells number |
[in] | nc_z | z cells number |
Definition at line 106 of file sll_m_maxwell_3d_pstd.F90.