8#ifndef WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_YEE_H_
9#define WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_YEE_H_
30 std::array<amrex::Real,3>& cell_size,
35 using namespace amrex;
37 stencil_coefs_x.resize(1);
38 stencil_coefs_x[0] = 1._rt/cell_size[0];
39 stencil_coefs_y.resize(1);
40 stencil_coefs_y[0] = 1._rt/cell_size[1];
41 stencil_coefs_z.resize(1);
42 stencil_coefs_z[0] = 1._rt/cell_size[2];
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) );
105 template<
typename T_Field>
109 amrex::Real
const *
const coefs_x,
int const ,
110 int const i,
int const j,
int const k,
int const ncomp=0 ) {
112 using namespace amrex;
113#if (defined WARPX_DIM_1D_Z)
117 amrex::Real
const inv_dx2 = coefs_x[0]*coefs_x[0];
118 return inv_dx2*( F(i-1,j,k,ncomp) - 2._rt*F(i,j,k,ncomp) + F(i+1,j,k,ncomp) );
127 amrex::Real
const *
const coefs_y,
int const n_coefs_y,
128 int const i,
int const j,
int const k,
int const ncomp=0 ) {
130 using namespace amrex;
131#if defined WARPX_DIM_3D
132 Real
const inv_dy = coefs_y[0];
134 return inv_dy*( F(i,j+1,k,ncomp) - F(i,j,k,ncomp) );
135#elif (defined WARPX_DIM_XZ || WARPX_DIM_1D_Z)
147 template<
typename T_Field>
151 amrex::Real
const *
const coefs_y,
int const n_coefs_y,
152 int const i,
int const j,
int const k,
int const ncomp=0 ) {
154 using namespace amrex;
155#if defined WARPX_DIM_3D
156 Real
const inv_dy = coefs_y[0];
158 return inv_dy*( F(i,j,k,ncomp) - F(i,j-1,k,ncomp) );
159#elif (defined WARPX_DIM_XZ || WARPX_DIM_1D_Z)
171 template<
typename T_Field>
175 amrex::Real
const *
const coefs_y,
int const ,
176 int const i,
int const j,
int const k,
int const ncomp=0 ) {
178 using namespace amrex;
179#if defined WARPX_DIM_3D
180 Real
const inv_dy2 = coefs_y[0]*coefs_y[0];
181 return inv_dy2*( F(i,j-1,k,ncomp) - 2._rt*F(i,j,k,ncomp) + F(i,j+1,k,ncomp) );
182#elif (defined WARPX_DIM_XZ || WARPX_DIM_1D_Z)
193 amrex::Real
const *
const coefs_z,
int const ,
194 int const i,
int const j,
int const k,
int const ncomp=0 ) {
196 using namespace amrex;
197 Real
const inv_dz = coefs_z[0];
198#if defined WARPX_DIM_3D
199 return inv_dz*( F(i,j,k+1,ncomp) - F(i,j,k,ncomp) );
200#elif (defined WARPX_DIM_XZ)
201 return inv_dz*( F(i,j+1,k,ncomp) - F(i,j,k,ncomp) );
202#elif (defined WARPX_DIM_1D_Z)
203 return inv_dz*( F(i+1,j,k,ncomp) - F(i,j,k,ncomp) );
209 template<
typename T_Field>
213 amrex::Real
const *
const coefs_z,
int const ,
214 int const i,
int const j,
int const k,
int const ncomp=0 ) {
216 using namespace amrex;
217 Real
const inv_dz = coefs_z[0];
218#if defined WARPX_DIM_3D
219 return inv_dz*( F(i,j,k,ncomp) - F(i,j,k-1,ncomp) );
220#elif (defined WARPX_DIM_XZ)
221 return inv_dz*( F(i,j,k,ncomp) - F(i,j-1,k,ncomp) );
222#elif (defined WARPX_DIM_1D_Z)
223 return inv_dz*( F(i,j,k,ncomp) - F(i-1,j,k,ncomp) );
229 template<
typename T_Field>
233 amrex::Real
const *
const coefs_z,
int const ,
234 int const i,
int const j,
int const k,
int const ncomp=0 ) {
236 using namespace amrex;
237 Real
const inv_dz2 = coefs_z[0]*coefs_z[0];
238#if defined WARPX_DIM_3D
239 return inv_dz2*( F(i,j,k-1,ncomp) - 2._rt*F(i,j,k,ncomp) + F(i,j,k+1,ncomp) );
240#elif (defined WARPX_DIM_XZ)
241 return inv_dz2*( F(i,j-1,k,ncomp) - 2._rt*F(i,j,k,ncomp) + F(i,j+1,k,ncomp) );
242#elif (defined WARPX_DIM_1D_Z)
243 return inv_dz2*( F(i-1,j,k,ncomp) - 2._rt*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)
constexpr auto c
vacuum speed of light [m/s]
Definition constant.H:149
__host__ __device__ void ignore_unused(const Ts &...)
Definition CartesianYeeAlgorithm.H:27
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real UpwardDz(amrex::Array4< amrex::Real const > 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:191
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real Dxx(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:107
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 UpwardDx(amrex::Array4< amrex::Real const > 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
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 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:149
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real UpwardDy(amrex::Array4< amrex::Real const > 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:125
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real Dzz(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:231
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:211
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real Dyy(T_Field const &F, amrex::Real const *const coefs_y, 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