8 #ifndef WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_YEE_H_ 9 #define WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_YEE_H_ 15 #include <AMReX_REAL.H> 16 #include <AMReX_Array4.H> 17 #include <AMReX_Gpu.H> 31 amrex::Vector<amrex::Real>& stencil_coefs_x,
32 amrex::Vector<amrex::Real>& stencil_coefs_y,
33 amrex::Vector<amrex::Real>& stencil_coefs_z ) {
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];
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])
60 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
62 amrex::Array4<amrex::Real>
const& F,
63 amrex::Real
const *
const coefs_x,
int const ,
64 int const i,
int const j,
int const k,
int const ncomp=0 ) {
66 amrex::Real
const inv_dx = coefs_x[0];
67 return inv_dx*( F(i+1,j,k,ncomp) - F(i,j,k,ncomp) );
72 template<
typename T_Field>
73 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
76 amrex::Real
const *
const coefs_x,
int const ,
77 int const i,
int const j,
int const k,
int const ncomp=0 ) {
79 amrex::Real
const inv_dx = coefs_x[0];
80 return inv_dx*( F(i,j,k,ncomp) - F(i-1,j,k,ncomp) );
85 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
87 amrex::Array4<amrex::Real>
const& F,
88 amrex::Real
const *
const coefs_y,
int const n_coefs_y,
89 int const i,
int const j,
int const k,
int const ncomp=0 ) {
91 using namespace amrex;
92 #if defined WARPX_DIM_3D 93 Real
const inv_dy = coefs_y[0];
94 return inv_dy*( F(i,j+1,k,ncomp) - F(i,j,k,ncomp) );
95 amrex::ignore_unused(n_coefs_y);
96 #elif (defined WARPX_DIM_XZ) 97 amrex::ignore_unused(F, coefs_y, n_coefs_y,
101 amrex::ignore_unused(F, coefs_y, n_coefs_y,
108 template<
typename T_Field>
109 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
112 amrex::Real
const *
const coefs_y,
int const n_coefs_y,
113 int const i,
int const j,
int const k,
int const ncomp=0 ) {
115 using namespace amrex;
116 #if defined WARPX_DIM_3D 117 Real
const inv_dy = coefs_y[0];
118 return inv_dy*( F(i,j,k,ncomp) - F(i,j-1,k,ncomp) );
119 amrex::ignore_unused(n_coefs_y);
120 #elif (defined WARPX_DIM_XZ) 121 amrex::ignore_unused(F, coefs_y, n_coefs_y,
125 amrex::ignore_unused(F, coefs_y, n_coefs_y,
132 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
134 amrex::Array4<amrex::Real>
const& F,
135 amrex::Real
const *
const coefs_z,
int const ,
136 int const i,
int const j,
int const k,
int const ncomp=0 ) {
138 using namespace amrex;
139 Real
const inv_dz = coefs_z[0];
140 #if defined WARPX_DIM_3D 141 return inv_dz*( F(i,j,k+1,ncomp) - F(i,j,k,ncomp) );
142 #elif (defined WARPX_DIM_XZ) 143 return inv_dz*( F(i,j+1,k,ncomp) - F(i,j,k,ncomp) );
149 template<
typename T_Field>
150 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
153 amrex::Real
const *
const coefs_z,
int const ,
154 int const i,
int const j,
int const k,
int const ncomp=0 ) {
156 using namespace amrex;
157 Real
const inv_dz = coefs_z[0];
158 #if defined WARPX_DIM_3D 159 return inv_dz*( F(i,j,k,ncomp) - F(i,j,k-1,ncomp) );
160 #elif (defined WARPX_DIM_XZ) 161 return inv_dz*( F(i,j,k,ncomp) - F(i,j-1,k,ncomp) );
167 #endif // WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_YEE_H_ 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
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:86
int dx
Definition: compute_domain.py:35
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:61
cell_size
Definition: compute_domain.py:37
Definition: CartesianYeeAlgorithm.H:27
i
Definition: check_interp_points_and_weights.py:171
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:151
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:110
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:74
static amrex::Real ComputeMaxDt(amrex::Real const *const dx)
Definition: CartesianYeeAlgorithm.H:48
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:133