Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
Data Types | Modules | Macros | Functions/Subroutines
sll_m_maxwell_3d_pstd.F90 File Reference
#include "sll_fftw.h"
#include "sll_memory.h"
#include "sll_working_precision.h"
#include "sll_maxwell_solvers_macros.h"
Include dependency graph for sll_m_maxwell_3d_pstd.F90:

Go to the source code of this file.

Data Types

type  sll_t_maxwell_3d_pstd
 Maxwell solver object. More...
 

Modules

module  sll_m_maxwell_3d_pstd
 Implements the Maxwell solver in 3D with periodic boundary conditions with PSTD method.
 

Macros

#define D_DX(field)
 
#define D_DY(field)
 
#define D_DZ(field)
 

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...
 

Data Type Documentation

◆ sll_m_maxwell_3d_pstd::sll_t_maxwell_3d_pstd

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.

Collaboration diagram for sll_t_maxwell_3d_pstd:
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

Macro Definition Documentation

◆ D_DX

#define D_DX (   field)
Value:
1 self%d_dx = field; \
2 call sll_s_fft_exec_r2c_1d(self%fwx, self%d_dx, self%tmp_x); \
3 self%tmp_x(2:nc_x/2 + 1) = -cmplx(0.0_f64, self%kx(2:nc_x/2 + 1), kind=f64)*self%tmp_x(2:nc_x/2 + 1); \
4 call sll_s_fft_exec_c2r_1d(self%bwx, self%tmp_x, self%d_dx); \
5 self%d_dx = self%d_dx/nc_x

◆ D_DY

#define D_DY (   field)
Value:
1 self%d_dy = field; \
2 call sll_s_fft_exec_r2c_1d(self%fwy, self%d_dy, self%tmp_y); \
3 self%tmp_y(2:nc_y/2 + 1) = -cmplx(0.0_f64, self%ky(2:nc_y/2 + 1), kind=f64)*self%tmp_y(2:nc_y/2 + 1); \
4 call sll_s_fft_exec_c2r_1d(self%bwy, self%tmp_y, self%d_dy); \
5 self%d_dy = self%d_dy/nc_y

◆ D_DZ

#define D_DZ (   field)
Value:
1 self%d_dz = field; \
2 call sll_s_fft_exec_r2c_1d(self%fwz, self%d_dz, self%tmp_z); \
3 self%tmp_z(2:nc_z/2 + 1) = -cmplx(0.0_f64, self%kz(2:nc_z/2 + 1), kind=f64)*self%tmp_z(2:nc_z/2 + 1); \
4 call sll_s_fft_exec_c2r_1d(self%bwz, self%tmp_z, self%d_dz); \
5 self%d_dz = self%d_dz/nc_z
    Report Typos and Errors