7 #ifndef ABLASTR_NODALFIELDGATHER_H_
8 #define ABLASTR_NODALFIELDGATHER_H_
34 const amrex::ParticleReal yp,
35 const amrex::ParticleReal zp,
38 int&
i,
int& j,
int& k, amrex::Real W[AMREX_SPACEDIM][2]) noexcept
40 using namespace amrex::literals;
41 #if (defined WARPX_DIM_3D)
42 const amrex::Real x = (xp - plo[0]) * dxi[0];
43 const amrex::Real y = (yp - plo[1]) * dxi[1];
44 const amrex::Real z = (zp - plo[2]) * dxi[2];
46 i =
static_cast<int>(amrex::Math::floor(x));
47 j =
static_cast<int>(amrex::Math::floor(y));
48 k =
static_cast<int>(amrex::Math::floor(z));
54 W[0][0] = 1.0_rt - W[0][1];
55 W[1][0] = 1.0_rt - W[1][1];
56 W[2][0] = 1.0_rt - W[2][1];
58 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
60 # if (defined WARPX_DIM_XZ)
61 const amrex::Real x = (xp - plo[0]) * dxi[0];
63 i =
static_cast<int>(amrex::Math::floor(x));
65 # elif (defined WARPX_DIM_RZ)
66 const amrex::Real r = (std::sqrt(xp*xp+yp*yp) - plo[0]) * dxi[0];
67 i =
static_cast<int>(amrex::Math::floor(r));
71 const amrex::Real z = (zp - plo[1]) * dxi[1];
72 j =
static_cast<int>(amrex::Math::floor(z));
75 W[0][0] = 1.0_rt - W[0][1];
76 W[1][0] = 1.0_rt - W[1][1];
94 const amrex::Real W[AMREX_SPACEDIM][2],
97 #if (defined WARPX_DIM_3D)
98 amrex::Real
value = 0;
99 value += scalar_field(
i, j , k ) * W[0][0] * W[1][0] * W[2][0];
100 value += scalar_field(
i+1, j , k ) * W[0][1] * W[1][0] * W[2][0];
101 value += scalar_field(
i, j+1, k ) * W[0][0] * W[1][1] * W[2][0];
102 value += scalar_field(
i+1, j+1, k ) * W[0][1] * W[1][1] * W[2][0];
103 value += scalar_field(
i, j , k+1) * W[0][0] * W[1][0] * W[2][1];
104 value += scalar_field(
i+1, j , k+1) * W[0][1] * W[1][0] * W[2][1];
105 value += scalar_field(
i , j+1, k+1) * W[0][0] * W[1][1] * W[2][1];
106 value += scalar_field(
i+1, j+1, k+1) * W[0][1] * W[1][1] * W[2][1];
107 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
108 amrex::Real
value = 0;
109 value += scalar_field(
i, j , k) * W[0][0] * W[1][0];
110 value += scalar_field(
i+1, j , k) * W[0][1] * W[1][0];
111 value += scalar_field(
i, j+1, k) * W[0][0] * W[1][1];
112 value += scalar_field(
i+1, j+1, k) * W[0][1] * W[1][1];
114 const amrex::Real
value = 0;
132 const amrex::ParticleReal yp,
133 const amrex::ParticleReal zp,
140 amrex::Real W[AMREX_SPACEDIM][2];
158 const amrex::ParticleReal yp,
159 const amrex::ParticleReal zp,
168 amrex::Real W[AMREX_SPACEDIM][2];
#define AMREX_GPU_HOST_DEVICE
#define ABLASTR_ABORT_WITH_MESSAGE(MSG)
Definition: TextMsg.H:77
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:93
AMREX_GPU_HOST_DEVICE AMREX_INLINE void compute_weights_nodal(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 in interpolating a nodal field to the given coordinates.
Definition: NodalFieldGather.H:33
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:131
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:157
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