8 #ifndef WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_NODAL_H_ 9 #define WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_NODAL_H_ 14 #include <AMReX_Array4.H> 15 #include <AMReX_Gpu.H> 16 #include <AMReX_REAL.H> 30 amrex::Vector<amrex::Real>& stencil_coefs_x,
31 amrex::Vector<amrex::Real>& stencil_coefs_y,
32 amrex::Vector<amrex::Real>& stencil_coefs_z ) {
34 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(
52 + 1._rt/(dx[1]*dx[1]),
63 return (amrex::IntVect(AMREX_D_DECL(1,1,1)));
70 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
72 amrex::Array4<amrex::Real>
const& F,
73 amrex::Real
const *
const coefs_x,
int const ,
74 int const i,
int const j,
int const k,
int const ncomp=0 ) {
76 using namespace amrex;
77 Real
const inv_dx = coefs_x[0];
78 return 0.5_rt*inv_dx*( F(i+1,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_x,
int const n_coefs_x,
89 int const i,
int const j,
int const k,
int const ncomp=0 ) {
91 return UpwardDx( F, coefs_x, n_coefs_x, i, j, k ,ncomp);
99 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
101 amrex::Array4<amrex::Real>
const& F,
102 amrex::Real
const *
const coefs_y,
int const ,
103 int const i,
int const j,
int const k,
int const ncomp=0 ) {
105 using namespace amrex;
106 #if defined WARPX_DIM_3D 107 Real
const inv_dy = coefs_y[0];
108 return 0.5_rt*inv_dy*( F(i,j+1,k,ncomp) - F(i,j-1,k,ncomp) );
109 #elif (defined WARPX_DIM_XZ) 110 ignore_unused(i, j, k, coefs_y, ncomp, F);
119 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
121 amrex::Array4<amrex::Real>
const& F,
122 amrex::Real
const *
const coefs_y,
int const n_coefs_y,
123 int const i,
int const j,
int const k,
int const ncomp=0 ) {
125 return UpwardDy( F, coefs_y, n_coefs_y, i, j, k ,ncomp);
133 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
135 amrex::Array4<amrex::Real>
const& F,
136 amrex::Real
const *
const coefs_z,
int const ,
137 int const i,
int const j,
int const k,
int const ncomp=0 ) {
139 using namespace amrex;
140 Real
const inv_dz = coefs_z[0];
141 #if defined WARPX_DIM_3D 142 return 0.5_rt*inv_dz*( F(i,j,k+1,ncomp) - F(i,j,k-1,ncomp) );
143 #elif (defined WARPX_DIM_XZ) 144 return 0.5_rt*inv_dz*( F(i,j+1,k,ncomp) - F(i,j-1,k,ncomp) );
152 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
154 amrex::Array4<amrex::Real>
const& F,
155 amrex::Real
const *
const coefs_z,
int const n_coefs_z,
156 int const i,
int const j,
int const k,
int const ncomp=0 ) {
158 return UpwardDz( F, coefs_z, n_coefs_z, i, j, k ,ncomp);
164 #endif // WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_NODAL_H_ AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDx(amrex::Array4< amrex::Real > const &F, amrex::Real const *const coefs_x, int const n_coefs_x, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianNodalAlgorithm.H:86
Definition: CartesianNodalAlgorithm.H:26
static amrex::IntVect GetMaxGuardCell()
Returns maximum number of guard cells required by the field-solve.
Definition: CartesianNodalAlgorithm.H:61
int dx
Definition: compute_domain.py:35
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, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianNodalAlgorithm.H:100
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: CartesianNodalAlgorithm.H:134
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDy(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: CartesianNodalAlgorithm.H:120
cell_size
Definition: compute_domain.py:37
i
Definition: check_interp_points_and_weights.py:171
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDz(amrex::Array4< amrex::Real > const &F, amrex::Real const *const coefs_z, int const n_coefs_z, int const i, int const j, int const k, int const ncomp=0)
Definition: CartesianNodalAlgorithm.H:153
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: CartesianNodalAlgorithm.H:28
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: CartesianNodalAlgorithm.H:71
static amrex::Real ComputeMaxDt(amrex::Real const *const dx)
Definition: CartesianNodalAlgorithm.H:48
Definition: BreitWheelerEngineWrapper.H:35