8 #ifndef WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_YEE_H_
9 #define WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_YEE_H_
38 using namespace amrex;
40 stencil_coefs_x.resize(1);
42 stencil_coefs_y.resize(1);
44 stencil_coefs_z.resize(1);
52 using namespace amrex::literals;
53 amrex::Real
const delta_t = 1._rt / ( std::sqrt(
AMREX_D_TERM(
54 1._rt / (
dx[0]*
dx[0]),
55 + 1._rt / (
dx[1]*
dx[1]),
56 + 1._rt / (
dx[2]*
dx[2])
74 amrex::Real
const *
const coefs_x,
int const ,
75 int const i,
int const j,
int const k,
int const ncomp=0 ) {
77 using namespace amrex;
78 #if (defined WARPX_DIM_1D_Z)
82 amrex::Real
const inv_dx = coefs_x[0];
83 return inv_dx*(
F(
i+1,j,k,ncomp) -
F(
i,j,k,ncomp) );
89 template<
typename T_Field>
93 amrex::Real
const *
const coefs_x,
int const ,
94 int const i,
int const j,
int const k,
int const ncomp=0 ) {
96 using namespace amrex;
97 #if (defined WARPX_DIM_1D_Z)
101 amrex::Real
const inv_dx = coefs_x[0];
102 return inv_dx*(
F(
i,j,k,ncomp) -
F(
i-1,j,k,ncomp) );
111 amrex::Real
const *
const coefs_y,
int const n_coefs_y,
112 int const i,
int const j,
int const k,
int const ncomp=0 ) {
114 using namespace amrex;
115 #if defined WARPX_DIM_3D
116 Real
const inv_dy = coefs_y[0];
118 return inv_dy*(
F(
i,j+1,k,ncomp) -
F(
i,j,k,ncomp) );
119 #elif (defined WARPX_DIM_XZ || WARPX_DIM_1D_Z)
132 template<
typename T_Field>
136 amrex::Real
const *
const coefs_y,
int const n_coefs_y,
137 int const i,
int const j,
int const k,
int const ncomp=0 ) {
139 using namespace amrex;
140 #if defined WARPX_DIM_3D
141 Real
const inv_dy = coefs_y[0];
143 return inv_dy*(
F(
i,j,k,ncomp) -
F(
i,j-1,k,ncomp) );
144 #elif (defined WARPX_DIM_XZ || WARPX_DIM_1D_Z)
160 amrex::Real
const *
const coefs_z,
int const ,
161 int const i,
int const j,
int const k,
int const ncomp=0 ) {
163 using namespace amrex;
164 Real
const inv_dz = coefs_z[0];
165 #if defined WARPX_DIM_3D
166 return inv_dz*(
F(
i,j,k+1,ncomp) -
F(
i,j,k,ncomp) );
167 #elif (defined WARPX_DIM_XZ)
168 return inv_dz*(
F(
i,j+1,k,ncomp) -
F(
i,j,k,ncomp) );
169 #elif (defined WARPX_DIM_1D_Z)
170 return inv_dz*(
F(
i+1,j,k,ncomp) -
F(
i,j,k,ncomp) );
176 template<
typename T_Field>
180 amrex::Real
const *
const coefs_z,
int const ,
181 int const i,
int const j,
int const k,
int const ncomp=0 ) {
183 using namespace amrex;
184 Real
const inv_dz = coefs_z[0];
185 #if defined WARPX_DIM_3D
186 return inv_dz*(
F(
i,j,k,ncomp) -
F(
i,j,k-1,ncomp) );
187 #elif (defined WARPX_DIM_XZ)
188 return inv_dz*(
F(
i,j,k,ncomp) -
F(
i,j-1,k,ncomp) );
189 #elif (defined WARPX_DIM_1D_Z)
190 return inv_dz*(
F(
i,j,k,ncomp) -
F(
i-1,j,k,ncomp) );
199 static amrex::Real LaplacianDx (
201 amrex::Real
const *
const coefs_x,
202 int const i,
int const j,
int const k,
int const ncomp=0 ) {
204 amrex::Real
const inv_dx = coefs_x[0];
205 return inv_dx*inv_dx*(
F(
i+1,j,k,ncomp) - 2.0*
F(
i,j,k,ncomp) +
F(
i-1,j,k,ncomp));
211 static amrex::Real LaplacianDy (
213 amrex::Real
const *
const coefs_y,
214 int const i,
int const j,
int const k,
int const ncomp=0 ) {
216 amrex::Real
const inv_dy = coefs_y[0];
217 return inv_dy*inv_dy*(
F(
i,j+1,k,ncomp) - 2.0*
F(
i,j,k,ncomp) +
F(
i,j-1,k,ncomp));
223 static amrex::Real LaplacianDz (
225 amrex::Real
const *
const coefs_z,
226 int const i,
int const j,
int const k,
int const ncomp=0 ) {
228 amrex::Real
const inv_dz = coefs_z[0];
229 return inv_dz*inv_dz*(
F(
i,j,k+1,ncomp) - 2.0*
F(
i,j,k,ncomp) +
F(
i,j,k-1,ncomp));
235 static amrex::Real Laplacian (
237 amrex::Real
const *
const coefs_x, amrex::Real
const *
const coefs_y, amrex::Real
const *
const coefs_z,
238 int const i,
int const j,
int const k,
int const ncomp=0 ) {
240 return LaplacianDx(F, coefs_x,
i, j, k, ncomp) + LaplacianDy(F, coefs_y,
i, j, k, ncomp) + LaplacianDz(F, coefs_z,
i, j, k, ncomp);
250 static amrex::Real LaplacianDx_Mag (
252 amrex::Real
const *
const coefs_x,
int const n_coefs_x, amrex::Real
const Ms_lo_x, amrex::Real
const Ms_hi_x,
253 int const i,
int const j,
int const k,
int const ncomp=0,
int const nodality=0) {
255 amrex::Real
const inv_dx = coefs_x[0];
258 return 0.5 * inv_dx * inv_dx * (8. *
F(
i-1, j, k, ncomp) -
F(
i-2, j, k, ncomp) - 7. *
F(
i, j, k, ncomp));
259 }
else if (Ms_lo_x == 0){
260 return 0.5 * inv_dx * inv_dx * (8. *
F(
i+1, j, k, ncomp) -
F(
i+2, j, k, ncomp) - 7. *
F(
i, j, k, ncomp));
262 return inv_dx*(
UpwardDx(F, coefs_x, n_coefs_x,
i, j, k, ncomp) -
DownwardDx(F, coefs_x, n_coefs_x,
i, j, k, ncomp));
266 return inv_dx*(0. -
DownwardDx(F, coefs_x, n_coefs_x,
i, j, k, ncomp));
267 }
else if (Ms_lo_x == 0.){
268 return inv_dx*(
UpwardDx(F, coefs_x, n_coefs_x,
i, j, k, ncomp) - 0.);
270 return inv_dx*(
UpwardDx(F, coefs_x, n_coefs_x,
i, j, k, ncomp) -
DownwardDx(F, coefs_x, n_coefs_x,
i, j, k, ncomp));
278 static amrex::Real LaplacianDy_Mag (
280 amrex::Real
const *
const coefs_y,
int const n_coefs_y, amrex::Real
const Ms_lo_y, amrex::Real
const Ms_hi_y,
281 int const i,
int const j,
int const k,
int const ncomp=0,
int const nodality=0) {
283 amrex::Real
const inv_dy = coefs_y[0];
286 return 0.5 * inv_dy * inv_dy * (8. *
F(
i, j-1, k, ncomp) -
F(
i, j-2, k, ncomp) - 7. *
F(
i, j, k, ncomp));
287 }
else if (Ms_lo_y == 0.) {
288 return 0.5 * inv_dy * inv_dy * (8. *
F(
i, j+1, k, ncomp) -
F(
i, j+2, k, ncomp) - 7. *
F(
i, j, k, ncomp));
290 return inv_dy*(
UpwardDy(F, coefs_y, n_coefs_y,
i, j, k, ncomp) -
DownwardDy(F, coefs_y, n_coefs_y,
i, j, k, ncomp));
295 return inv_dy*(0. -
DownwardDy(F, coefs_y, n_coefs_y,
i, j, k, ncomp));
296 }
else if (Ms_lo_y == 0.) {
297 return inv_dy*(
UpwardDy(F, coefs_y, n_coefs_y,
i, j, k, ncomp) - 0.);
299 return inv_dy*(
UpwardDy(F, coefs_y, n_coefs_y,
i, j, k, ncomp) -
DownwardDy(F, coefs_y, n_coefs_y,
i, j, k, ncomp));
307 static amrex::Real LaplacianDz_Mag (
309 amrex::Real
const *
const coefs_z,
int const n_coefs_z, amrex::Real
const Ms_lo_z, amrex::Real
const Ms_hi_z,
310 int const i,
int const j,
int const k,
int const ncomp=0,
int const nodality=0) {
312 amrex::Real
const inv_dz = coefs_z[0];
314 if ( Ms_hi_z == 0.) {
315 return 0.5 * inv_dz * inv_dz * (8. *
F(
i, j, k-1, ncomp) -
F(
i, j, k-2, ncomp) - 7. *
F(
i, j, k, ncomp));
316 }
else if ( Ms_lo_z == 0.){
317 return 0.5 * inv_dz * inv_dz * (8. *
F(
i, j, k+1, ncomp) -
F(
i, j, k+2, ncomp) - 7. *
F(
i, j, k, ncomp));
319 return inv_dz*(
UpwardDz(F, coefs_z, n_coefs_z,
i, j, k, ncomp) -
DownwardDz(F, coefs_z, n_coefs_z,
i, j, k, ncomp));
323 if ( Ms_hi_z == 0.) {
324 return inv_dz*(0. -
DownwardDz(F, coefs_z, n_coefs_z,
i, j, k, ncomp));
325 }
else if ( Ms_lo_z == 0.){
326 return inv_dz*(
UpwardDz(F, coefs_z, n_coefs_z,
i, j, k, ncomp) - 0.);
328 return inv_dz*(
UpwardDz(F, coefs_z, n_coefs_z,
i, j, k, ncomp) -
DownwardDz(F, coefs_z, n_coefs_z,
i, j, k, ncomp));
336 static amrex::Real Laplacian_Mag (
338 amrex::Real
const *
const coefs_x, amrex::Real
const *
const coefs_y, amrex::Real
const *
const coefs_z,
339 int const n_coefs_x,
int const n_coefs_y,
int const n_coefs_z,
340 amrex::Real
const Ms_lo_x, amrex::Real
const Ms_hi_x, amrex::Real
const Ms_lo_y, amrex::Real
const Ms_hi_y, amrex::Real
const Ms_lo_z, amrex::Real
const Ms_hi_z,
341 int const i,
int const j,
int const k,
int const ncomp=0,
int const nodality=0) {
343 return LaplacianDx_Mag(F, coefs_x, n_coefs_x, Ms_lo_x, Ms_hi_x,
i, j, k, ncomp, nodality) + LaplacianDy_Mag(F, coefs_y, n_coefs_y, Ms_lo_y, Ms_hi_y,
i, j, k, ncomp, nodality) + LaplacianDz_Mag(F, coefs_z, n_coefs_z, Ms_lo_z, Ms_hi_z,
i, j, k, ncomp, nodality);
#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
int dx
Definition: stencil.py:436
Definition: CartesianYeeAlgorithm.H:30
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:91
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:158
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:109
static amrex::IntVect GetMaxGuardCell()
Returns maximum number of guard cells required by the field-solve.
Definition: CartesianYeeAlgorithm.H:64
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:72
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:134
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:178
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:32
static amrex::Real ComputeMaxDt(amrex::Real const *const dx)
Definition: CartesianYeeAlgorithm.H:51