8#ifndef WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_NODAL_H_
9#define WARPX_FINITE_DIFFERENCE_ALGORITHM_CARTESIAN_NODAL_H_
29 std::array<amrex::Real,3>& cell_size,
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];
49 using namespace amrex::literals;
50 amrex::Real
const delta_t = 1._rt / ( std::sqrt(
AMREX_D_TERM(
52 + 1._rt/(dx[1]*dx[1]),
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#if (defined WARPX_DIM_1D_Z)
81 Real
const inv_dx = coefs_x[0];
82 return 0.5_rt*inv_dx*( F(i+1,j,k,ncomp) - F(i-1,j,k,ncomp) );
93 amrex::Real
const *
const coefs_x,
int const n_coefs_x,
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 return UpwardDx( F, coefs_x, n_coefs_x, i, j, k ,ncomp);
108 template<
typename T_Field>
112 amrex::Real
const *
const coefs_x,
int const ,
113 int const i,
int const j,
int const k,
int const ncomp=0 ) {
115 using namespace amrex;
116#if (defined WARPX_DIM_1D_Z)
120 amrex::Real
const inv_dx2 = coefs_x[0]*coefs_x[0];
121 return inv_dx2*( F(i-1,j,k,ncomp) - 2._rt*F(i,j,k,ncomp) + F(i+1,j,k,ncomp) );
132 amrex::Real
const *
const coefs_y,
int const ,
133 int const i,
int const j,
int const k,
int const ncomp=0 ) {
135 using namespace amrex;
136#if defined WARPX_DIM_3D
137 Real
const inv_dy = coefs_y[0];
138 return 0.5_rt*inv_dy*( F(i,j+1,k,ncomp) - F(i,j-1,k,ncomp) );
139#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_1D_Z)
152 amrex::Real
const *
const coefs_y,
int const n_coefs_y,
153 int const i,
int const j,
int const k,
int const ncomp=0 ) {
155 return UpwardDy( F, coefs_y, n_coefs_y, i, j, k ,ncomp);
161 template<
typename T_Field>
165 amrex::Real
const *
const coefs_y,
int const ,
166 int const i,
int const j,
int const k,
int const ncomp=0 ) {
168 using namespace amrex;
169#if defined WARPX_DIM_3D
170 Real
const inv_dy2 = coefs_y[0]*coefs_y[0];
171 return inv_dy2*( F(i,j-1,k,ncomp) - 2._rt*F(i,j,k,ncomp) + F(i,j+1,k,ncomp) );
172#elif (defined WARPX_DIM_XZ || WARPX_DIM_1D_Z)
185 amrex::Real
const *
const coefs_z,
int const ,
186 int const i,
int const j,
int const k,
int const ncomp=0 ) {
188 using namespace amrex;
189 Real
const inv_dz = coefs_z[0];
190#if defined WARPX_DIM_3D
191 return 0.5_rt*inv_dz*( F(i,j,k+1,ncomp) - F(i,j,k-1,ncomp) );
192#elif (defined WARPX_DIM_XZ)
193 return 0.5_rt*inv_dz*( F(i,j+1,k,ncomp) - F(i,j-1,k,ncomp) );
194#elif (defined WARPX_DIM_1D_Z)
195 return 0.5_rt*inv_dz*( F(i+1,j,k,ncomp) - F(i-1,j,k,ncomp) );
206 amrex::Real
const *
const coefs_z,
int const n_coefs_z,
207 int const i,
int const j,
int const k,
int const ncomp=0 ) {
209 return UpwardDz( F, coefs_z, n_coefs_z, i, j, k ,ncomp);
216 template<
typename T_Field>
220 amrex::Real
const *
const coefs_z,
int const ,
221 int const i,
int const j,
int const k,
int const ncomp=0 ) {
223 using namespace amrex;
224 Real
const inv_dz2 = coefs_z[0]*coefs_z[0];
225#if defined WARPX_DIM_3D
226 return inv_dz2*( F(i,j,k-1,ncomp) - 2._rt*F(i,j,k,ncomp) + F(i,j,k+1,ncomp) );
227#elif (defined WARPX_DIM_XZ)
228 return inv_dz2*( F(i,j-1,k,ncomp) - 2._rt*F(i,j,k,ncomp) + F(i,j+1,k,ncomp) );
229#elif (defined WARPX_DIM_1D_Z)
230 return inv_dz2*( F(i-1,j,k,ncomp) - 2._rt*F(i,j,k,ncomp) + F(i+1,j,k,ncomp) );
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
#define AMREX_D_TERM(a, b, c)
#define AMREX_D_DECL(a, b, c)
constexpr auto c
vacuum speed of light [m/s]
Definition constant.H:149
__host__ __device__ void ignore_unused(const Ts &...)
Definition CartesianNodalAlgorithm.H:26
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real UpwardDy(amrex::Array4< amrex::Real const > 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:130
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDz(amrex::Array4< amrex::Real const > 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:204
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real UpwardDz(amrex::Array4< amrex::Real const > 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:183
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real Dzz(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 CartesianNodalAlgorithm.H:218
static amrex::Real ComputeMaxDt(amrex::Real const *const dx)
Definition CartesianNodalAlgorithm.H:48
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real Dxx(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 CartesianNodalAlgorithm.H:110
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 > 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
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDx(amrex::Array4< amrex::Real const > 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:91
static amrex::IntVect GetMaxGuardCell()
Returns maximum number of guard cells required by the field-solve.
Definition CartesianNodalAlgorithm.H:61
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDy(amrex::Array4< amrex::Real const > 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:150
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real Dyy(T_Field 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:163