7 #ifndef ABLASTR_NODALFIELDGATHER_H_
8 #define ABLASTR_NODALFIELDGATHER_H_
36 template <amrex::IndexType::CellIndex IdxType>
39 const amrex::ParticleReal yp,
40 const amrex::ParticleReal zp,
43 int&
i,
int& j,
int& k, amrex::Real W[AMREX_SPACEDIM][2]) noexcept
45 using namespace amrex::literals;
47 constexpr
auto half_if_cell_centered = [](){
48 if constexpr (IdxType == amrex::IndexType::CellIndex::NODE){
56 #if (defined WARPX_DIM_3D)
57 const amrex::Real
x = (xp - plo[0]) * dxi[0] - half_if_cell_centered;
58 const amrex::Real y = (yp - plo[1]) * dxi[1] - half_if_cell_centered;
59 const amrex::Real z = (zp - plo[2]) * dxi[2] - half_if_cell_centered;
61 i =
static_cast<int>(amrex::Math::floor(
x));
62 j =
static_cast<int>(amrex::Math::floor(y));
63 k =
static_cast<int>(amrex::Math::floor(z));
69 W[0][0] = 1.0_rt - W[0][1];
70 W[1][0] = 1.0_rt - W[1][1];
71 W[2][0] = 1.0_rt - W[2][1];
73 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
75 # if (defined WARPX_DIM_XZ)
76 const amrex::Real
x = (xp - plo[0]) * dxi[0] - half_if_cell_centered;
78 i =
static_cast<int>(amrex::Math::floor(
x));
80 # elif (defined WARPX_DIM_RZ)
81 const amrex::Real
r = (std::sqrt(xp*xp+yp*yp) - plo[0]) * dxi[0] - half_if_cell_centered;
82 i =
static_cast<int>(amrex::Math::floor(
r));
86 const amrex::Real z = (zp - plo[1]) * dxi[1] - half_if_cell_centered;
87 j =
static_cast<int>(amrex::Math::floor(z));
90 W[0][0] = 1.0_rt - W[0][1];
91 W[1][0] = 1.0_rt - W[1][1];
95 const amrex::Real z = (zp - plo[0]) * dxi[0] - half_if_cell_centered;
97 i =
static_cast<int>(amrex::Math::floor(z));
100 W[0][0] = 1.0_rt - W[0][1];
116 const amrex::Real W[AMREX_SPACEDIM][2],
119 amrex::Real
value = 0;
120 #if (defined WARPX_DIM_3D)
121 value += scalar_field(
i, j , k ) * W[0][0] * W[1][0] * W[2][0];
122 value += scalar_field(
i+1, j , k ) * W[0][1] * W[1][0] * W[2][0];
123 value += scalar_field(
i, j+1, k ) * W[0][0] * W[1][1] * W[2][0];
124 value += scalar_field(
i+1, j+1, k ) * W[0][1] * W[1][1] * W[2][0];
125 value += scalar_field(
i, j , k+1) * W[0][0] * W[1][0] * W[2][1];
126 value += scalar_field(
i+1, j , k+1) * W[0][1] * W[1][0] * W[2][1];
127 value += scalar_field(
i , j+1, k+1) * W[0][0] * W[1][1] * W[2][1];
128 value += scalar_field(
i+1, j+1, k+1) * W[0][1] * W[1][1] * W[2][1];
129 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
130 value += scalar_field(
i, j , k) * W[0][0] * W[1][0];
131 value += scalar_field(
i+1, j , k) * W[0][1] * W[1][0];
132 value += scalar_field(
i, j+1, k) * W[0][0] * W[1][1];
133 value += scalar_field(
i+1, j+1, k) * W[0][1] * W[1][1];
135 value += scalar_field(
i, j , k) * W[0][0];
136 value += scalar_field(
i+1, j , k) * W[0][1];
152 const amrex::ParticleReal yp,
153 const amrex::ParticleReal zp,
160 amrex::Real W[AMREX_SPACEDIM][2];
161 compute_weights<amrex::IndexType::NODE>(xp, yp, zp, lo, dxi,
ii,
jj, kk, W);
178 const amrex::ParticleReal yp,
179 const amrex::ParticleReal zp,
188 amrex::Real W[AMREX_SPACEDIM][2];
189 compute_weights<amrex::IndexType::NODE>(xp, yp, zp, lo, dxi,
ii,
jj, kk, W);
#define AMREX_GPU_HOST_DEVICE
Definition: DepositCharge.H:21
AMREX_GPU_HOST_DEVICE AMREX_INLINE amrex::Real interp_field_nodal(int i, int j, int k, const amrex::Real W[AMREX_SPACEDIM][2], amrex::Array4< const amrex::Real > const &scalar_field) noexcept
Interpolate nodal field value based on surrounding indices and weights.
Definition: NodalFieldGather.H:115
AMREX_GPU_HOST_DEVICE AMREX_INLINE void compute_weights(const amrex::ParticleReal xp, const amrex::ParticleReal yp, const amrex::ParticleReal zp, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const &plo, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const &dxi, int &i, int &j, int &k, amrex::Real W[AMREX_SPACEDIM][2]) noexcept
Compute weight of each surrounding node (or cell-centered nodes) in interpolating a nodal ((or a cell...
Definition: NodalFieldGather.H:38
AMREX_GPU_HOST_DEVICE AMREX_INLINE amrex::Real doGatherScalarFieldNodal(const amrex::ParticleReal xp, const amrex::ParticleReal yp, const amrex::ParticleReal zp, amrex::Array4< const amrex::Real > const &scalar_field, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const &dxi, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const &lo) noexcept
Scalar field gather for a single particle. The field has to be defined at the cell nodes (see https:/...
Definition: NodalFieldGather.H:151
AMREX_GPU_HOST_DEVICE AMREX_INLINE amrex::GpuArray< amrex::Real, 3 > doGatherVectorFieldNodal(const amrex::ParticleReal xp, const amrex::ParticleReal yp, const amrex::ParticleReal zp, amrex::Array4< const amrex::Real > const &vector_field_x, amrex::Array4< const amrex::Real > const &vector_field_y, amrex::Array4< const amrex::Real > const &vector_field_z, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const &dxi, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const &lo) noexcept
Vector field gather for a single particle. The field has to be defined at the cell nodes (see https:/...
Definition: NodalFieldGather.H:177
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
i
Definition: check_interp_points_and_weights.py:174
ii
Definition: check_interp_points_and_weights.py:148
jj
Definition: check_interp_points_and_weights.py:160
value
Definition: updateAMReX.py:141