8 #ifndef WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_YEE_H_ 9 #define WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_YEE_H_ 15 #include <AMReX_Array4.H> 16 #include <AMReX_Gpu.H> 17 #include <AMReX_REAL.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])
63 return (amrex::IntVect(AMREX_D_DECL(1,1,1)));
68 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
70 amrex::Array4<amrex::Real>
const& F,
71 amrex::Real
const *
const coefs_x,
int const ,
72 int const i,
int const j,
int const k,
int const ncomp=0 ) {
74 amrex::Real
const inv_dx = coefs_x[0];
75 return inv_dx*( F(i+1,j,k,ncomp) - F(i,j,k,ncomp) );
80 template<
typename T_Field>
81 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
84 amrex::Real
const *
const coefs_x,
int const ,
85 int const i,
int const j,
int const k,
int const ncomp=0 ) {
87 amrex::Real
const inv_dx = coefs_x[0];
88 return inv_dx*( F(i,j,k,ncomp) - F(i-1,j,k,ncomp) );
93 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
95 amrex::Array4<amrex::Real>
const& F,
96 amrex::Real
const *
const coefs_y,
int const n_coefs_y,
97 int const i,
int const j,
int const k,
int const ncomp=0 ) {
99 using namespace amrex;
100 #if defined WARPX_DIM_3D 101 Real
const inv_dy = coefs_y[0];
102 amrex::ignore_unused(n_coefs_y);
103 return inv_dy*( F(i,j+1,k,ncomp) - F(i,j,k,ncomp) );
104 #elif (defined WARPX_DIM_XZ) 105 amrex::ignore_unused(F, coefs_y, n_coefs_y,
109 amrex::ignore_unused(F, coefs_y, n_coefs_y,
116 template<
typename T_Field>
117 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
120 amrex::Real
const *
const coefs_y,
int const n_coefs_y,
121 int const i,
int const j,
int const k,
int const ncomp=0 ) {
123 using namespace amrex;
124 #if defined WARPX_DIM_3D 125 Real
const inv_dy = coefs_y[0];
126 amrex::ignore_unused(n_coefs_y);
127 return inv_dy*( F(i,j,k,ncomp) - F(i,j-1,k,ncomp) );
128 #elif (defined WARPX_DIM_XZ) 129 amrex::ignore_unused(F, coefs_y, n_coefs_y,
133 amrex::ignore_unused(F, coefs_y, n_coefs_y,
140 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
142 amrex::Array4<amrex::Real>
const& F,
143 amrex::Real
const *
const coefs_z,
int const ,
144 int const i,
int const j,
int const k,
int const ncomp=0 ) {
146 using namespace amrex;
147 Real
const inv_dz = coefs_z[0];
148 #if defined WARPX_DIM_3D 149 return inv_dz*( F(i,j,k+1,ncomp) - F(i,j,k,ncomp) );
150 #elif (defined WARPX_DIM_XZ) 151 return inv_dz*( F(i,j+1,k,ncomp) - F(i,j,k,ncomp) );
157 template<
typename T_Field>
158 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
161 amrex::Real
const *
const coefs_z,
int const ,
162 int const i,
int const j,
int const k,
int const ncomp=0 ) {
164 using namespace amrex;
165 Real
const inv_dz = coefs_z[0];
166 #if defined WARPX_DIM_3D 167 return inv_dz*( F(i,j,k,ncomp) - F(i,j,k-1,ncomp) );
168 #elif (defined WARPX_DIM_XZ) 169 return inv_dz*( F(i,j,k,ncomp) - F(i,j-1,k,ncomp) );
175 #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:94
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:69
cell_size
Definition: compute_domain.py:37
Definition: CartesianYeeAlgorithm.H:27
i
Definition: check_interp_points_and_weights.py:171
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 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:159
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:118
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:82
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:141
Definition: BreitWheelerEngineWrapper.H:35