8 #ifndef WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_NODAL_H_ 9 #define WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_NODAL_H_ 13 #include <AMReX_REAL.H> 14 #include <AMReX_Array4.H> 15 #include <AMReX_Gpu.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]),
62 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
64 amrex::Array4<amrex::Real>
const& F,
65 amrex::Real
const *
const coefs_x,
int const ,
66 int const i,
int const j,
int const k,
int const ncomp=0 ) {
68 using namespace amrex;
69 Real
const inv_dx = coefs_x[0];
70 return 0.5_rt*inv_dx*( F(i+1,j,k,ncomp) - F(i-1,j,k,ncomp) );
77 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
79 amrex::Array4<amrex::Real>
const& F,
80 amrex::Real
const *
const coefs_x,
int const n_coefs_x,
81 int const i,
int const j,
int const k,
int const ncomp=0 ) {
83 return UpwardDx( F, coefs_x, n_coefs_x, i, j, k ,ncomp);
91 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
93 amrex::Array4<amrex::Real>
const& F,
94 amrex::Real
const *
const coefs_y,
int const ,
95 int const i,
int const j,
int const k,
int const ncomp=0 ) {
97 using namespace amrex;
98 #if defined WARPX_DIM_3D 99 Real
const inv_dy = coefs_y[0];
100 return 0.5_rt*inv_dy*( F(i,j+1,k,ncomp) - F(i,j-1,k,ncomp) );
101 #elif (defined WARPX_DIM_XZ) 102 ignore_unused(i, j, k, coefs_y, ncomp, F);
111 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
113 amrex::Array4<amrex::Real>
const& F,
114 amrex::Real
const *
const coefs_y,
int const n_coefs_y,
115 int const i,
int const j,
int const k,
int const ncomp=0 ) {
117 return UpwardDy( F, coefs_y, n_coefs_y, i, j, k ,ncomp);
125 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
127 amrex::Array4<amrex::Real>
const& F,
128 amrex::Real
const *
const coefs_z,
int const ,
129 int const i,
int const j,
int const k,
int const ncomp=0 ) {
131 using namespace amrex;
132 Real
const inv_dz = coefs_z[0];
133 #if defined WARPX_DIM_3D 134 return 0.5_rt*inv_dz*( F(i,j,k+1,ncomp) - F(i,j,k-1,ncomp) );
135 #elif (defined WARPX_DIM_XZ) 136 return 0.5_rt*inv_dz*( F(i,j+1,k,ncomp) - F(i,j-1,k,ncomp) );
144 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
146 amrex::Array4<amrex::Real>
const& F,
147 amrex::Real
const *
const coefs_z,
int const n_coefs_z,
148 int const i,
int const j,
int const k,
int const ncomp=0 ) {
150 return UpwardDz( F, coefs_z, n_coefs_z, i, j, k ,ncomp);
156 #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:78
Definition: CartesianNodalAlgorithm.H:26
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:92
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:126
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:112
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:145
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:63
static amrex::Real ComputeMaxDt(amrex::Real const *const dx)
Definition: CartesianNodalAlgorithm.H:48