8 #ifndef WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_YEE_H_
9 #define WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_YEE_H_
35 using namespace amrex;
37 stencil_coefs_x.resize(1);
39 stencil_coefs_y.resize(1);
41 stencil_coefs_z.resize(1);
49 using namespace amrex::literals;
50 amrex::Real
const delta_t = 1._rt / ( std::sqrt(
AMREX_D_TERM(
51 1._rt / (
dx[0]*
dx[0]),
52 + 1._rt / (
dx[1]*
dx[1]),
53 + 1._rt / (
dx[2]*
dx[2])
71 amrex::Real
const *
const coefs_x,
int const ,
72 int const i,
int const j,
int const k,
int const ncomp=0 ) {
74 using namespace amrex;
75 #if (defined WARPX_DIM_1D_Z)
79 amrex::Real
const inv_dx = coefs_x[0];
80 return inv_dx*( F(
i+1,j,k,ncomp) - F(
i,j,k,ncomp) );
86 template<
typename T_Field>
90 amrex::Real
const *
const coefs_x,
int const ,
91 int const i,
int const j,
int const k,
int const ncomp=0 ) {
93 using namespace amrex;
94 #if (defined WARPX_DIM_1D_Z)
98 amrex::Real
const inv_dx = coefs_x[0];
99 return inv_dx*( F(
i,j,k,ncomp) - F(
i-1,j,k,ncomp) );
108 amrex::Real
const *
const coefs_y,
int const n_coefs_y,
109 int const i,
int const j,
int const k,
int const ncomp=0 ) {
111 using namespace amrex;
112 #if defined WARPX_DIM_3D
113 Real
const inv_dy = coefs_y[0];
115 return inv_dy*( F(
i,j+1,k,ncomp) - F(
i,j,k,ncomp) );
116 #elif (defined WARPX_DIM_XZ || WARPX_DIM_1D_Z)
128 template<
typename T_Field>
132 amrex::Real
const *
const coefs_y,
int const n_coefs_y,
133 int const i,
int const j,
int const k,
int const ncomp=0 ) {
135 using namespace amrex;
136 #if defined WARPX_DIM_3D
137 Real
const inv_dy = coefs_y[0];
139 return inv_dy*( F(
i,j,k,ncomp) - F(
i,j-1,k,ncomp) );
140 #elif (defined WARPX_DIM_XZ || WARPX_DIM_1D_Z)
155 amrex::Real
const *
const coefs_z,
int const ,
156 int const i,
int const j,
int const k,
int const ncomp=0 ) {
158 using namespace amrex;
159 Real
const inv_dz = coefs_z[0];
160 #if defined WARPX_DIM_3D
161 return inv_dz*( F(
i,j,k+1,ncomp) - F(
i,j,k,ncomp) );
162 #elif (defined WARPX_DIM_XZ)
163 return inv_dz*( F(
i,j+1,k,ncomp) - F(
i,j,k,ncomp) );
164 #elif (defined WARPX_DIM_1D_Z)
165 return inv_dz*( F(
i+1,j,k,ncomp) - F(
i,j,k,ncomp) );
171 template<
typename T_Field>
175 amrex::Real
const *
const coefs_z,
int const ,
176 int const i,
int const j,
int const k,
int const ncomp=0 ) {
178 using namespace amrex;
179 Real
const inv_dz = coefs_z[0];
180 #if defined WARPX_DIM_3D
181 return inv_dz*( F(
i,j,k,ncomp) - F(
i,j,k-1,ncomp) );
182 #elif (defined WARPX_DIM_XZ)
183 return inv_dz*( F(
i,j,k,ncomp) - F(
i,j-1,k,ncomp) );
184 #elif (defined WARPX_DIM_1D_Z)
185 return inv_dz*( F(
i,j,k,ncomp) - F(
i-1,j,k,ncomp) );
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
#define AMREX_D_TERM(a, b, c)
#define AMREX_D_DECL(a, b, c)
static constexpr auto c
vacuum speed of light [m/s]
Definition: constant.H:44
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
i
Definition: check_interp_points_and_weights.py:174
cell_size
Definition: compute_domain.py:37
tuple dx
lab frame
Definition: stencil.py:429
Definition: CartesianYeeAlgorithm.H:27
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDx(T_Field const &F, amrex::Real const *const coefs_x, int const, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianYeeAlgorithm.H:88
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real UpwardDz(amrex::Array4< amrex::Real > const &F, amrex::Real const *const coefs_z, int const, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianYeeAlgorithm.H:153
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real UpwardDy(amrex::Array4< amrex::Real > const &F, amrex::Real const *const coefs_y, int const n_coefs_y, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianYeeAlgorithm.H:106
static amrex::IntVect GetMaxGuardCell()
Returns maximum number of guard cells required by the field-solve.
Definition: CartesianYeeAlgorithm.H:61
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real UpwardDx(amrex::Array4< amrex::Real > const &F, amrex::Real const *const coefs_x, int const, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianYeeAlgorithm.H:69
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDy(T_Field const &F, amrex::Real const *const coefs_y, int const n_coefs_y, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianYeeAlgorithm.H:130
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDz(T_Field const &F, amrex::Real const *const coefs_z, int const, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianYeeAlgorithm.H:173
static void InitializeStencilCoefficients(std::array< amrex::Real, 3 > &cell_size, amrex::Vector< amrex::Real > &stencil_coefs_x, amrex::Vector< amrex::Real > &stencil_coefs_y, amrex::Vector< amrex::Real > &stencil_coefs_z)
Definition: CartesianYeeAlgorithm.H:29
static amrex::Real ComputeMaxDt(amrex::Real const *const dx)
Definition: CartesianYeeAlgorithm.H:48