ARTEMIS
NodalFieldGather.H
Go to the documentation of this file.
1 /* Copyright 2019-2022 Modern Electron, Axel Huebl, Remi Lehe
2  *
3  * This file is part of WarpX.
4  *
5  * License: BSD-3-Clause-LBNL
6  */
7 #ifndef ABLASTR_NODALFIELDGATHER_H_
8 #define ABLASTR_NODALFIELDGATHER_H_
9 
10 #include <AMReX_Array.H>
11 #include <AMReX_Extension.H>
12 #include <AMReX_GpuQualifiers.H>
13 #include <AMReX_Math.H>
14 #include <AMReX_REAL.H>
15 
16 
17 namespace ablastr::particles
18 {
32 void compute_weights_nodal (const amrex::ParticleReal xp,
33  const amrex::ParticleReal yp,
34  const amrex::ParticleReal zp,
37  int& i, int& j, int& k, amrex::Real W[AMREX_SPACEDIM][2]) noexcept
38 {
39  using namespace amrex::literals;
40 #if (defined WARPX_DIM_3D)
41  amrex::Real x = (xp - plo[0]) * dxi[0];
42  amrex::Real y = (yp - plo[1]) * dxi[1];
43  amrex::Real z = (zp - plo[2]) * dxi[2];
44 
45  i = static_cast<int>(amrex::Math::floor(x));
46  j = static_cast<int>(amrex::Math::floor(y));
47  k = static_cast<int>(amrex::Math::floor(z));
48 
49  W[0][1] = x - i;
50  W[1][1] = y - j;
51  W[2][1] = z - k;
52 
53  W[0][0] = 1.0_rt - W[0][1];
54  W[1][0] = 1.0_rt - W[1][1];
55  W[2][0] = 1.0_rt - W[2][1];
56 
57 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
58 
59 # if (defined WARPX_DIM_XZ)
60  amrex::Real x = (xp - plo[0]) * dxi[0];
62  i = static_cast<int>(amrex::Math::floor(x));
63  W[0][1] = x - i;
64 # elif (defined WARPX_DIM_RZ)
65  amrex::Real r = (std::sqrt(xp*xp+yp*yp) - plo[0]) * dxi[0];
66  i = static_cast<int>(amrex::Math::floor(r));
67  W[0][1] = r - i;
68 # endif
69 
70  amrex::Real z = (zp - plo[1]) * dxi[1];
71  j = static_cast<int>(amrex::Math::floor(z));
72  W[1][1] = z - j;
73 
74  W[0][0] = 1.0_rt - W[0][1];
75  W[1][0] = 1.0_rt - W[1][1];
76 
77  k = 0;
78 #else
79  amrex::ignore_unused(xp, yp, zp, plo, dxi, i, j, k, W);
80  amrex::Abort("Error: compute_weights not yet implemented in 1D");
81 #endif
82 }
83 
92 amrex::Real interp_field_nodal (int i, int j, int k,
93  const amrex::Real W[AMREX_SPACEDIM][2],
94  amrex::Array4<const amrex::Real> const& scalar_field) noexcept
95 {
96  amrex::Real value = 0;
97 #if (defined WARPX_DIM_3D)
98  value += scalar_field(i, j , k ) * W[0][0] * W[1][0] * W[2][0];
99  value += scalar_field(i+1, j , k ) * W[0][1] * W[1][0] * W[2][0];
100  value += scalar_field(i, j+1, k ) * W[0][0] * W[1][1] * W[2][0];
101  value += scalar_field(i+1, j+1, k ) * W[0][1] * W[1][1] * W[2][0];
102  value += scalar_field(i, j , k+1) * W[0][0] * W[1][0] * W[2][1];
103  value += scalar_field(i+1, j , k+1) * W[0][1] * W[1][0] * W[2][1];
104  value += scalar_field(i , j+1, k+1) * W[0][0] * W[1][1] * W[2][1];
105  value += scalar_field(i+1, j+1, k+1) * W[0][1] * W[1][1] * W[2][1];
106 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
107  value += scalar_field(i, j , k) * W[0][0] * W[1][0];
108  value += scalar_field(i+1, j , k) * W[0][1] * W[1][0];
109  value += scalar_field(i, j+1, k) * W[0][0] * W[1][1];
110  value += scalar_field(i+1, j+1, k) * W[0][1] * W[1][1];
111 #else
112  amrex::ignore_unused(i, j, k, W, scalar_field);
113  amrex::Abort("Error: interp_field not yet implemented in 1D");
114 #endif
115  return value;
116 }
117 
128 amrex::Real doGatherScalarFieldNodal (const amrex::ParticleReal xp,
129  const amrex::ParticleReal yp,
130  const amrex::ParticleReal zp,
131  amrex::Array4<const amrex::Real> const& scalar_field,
134 {
135  // first find the weight of surrounding nodes to use during interpolation
136  int ii, jj, kk;
137  amrex::Real W[AMREX_SPACEDIM][2];
138  compute_weights_nodal(xp, yp, zp, lo, dxi, ii, jj, kk, W);
139 
140  return interp_field_nodal(ii, jj, kk, W, scalar_field);
141 }
142 
154 doGatherVectorFieldNodal (const amrex::ParticleReal xp,
155  const amrex::ParticleReal yp,
156  const amrex::ParticleReal zp,
157  amrex::Array4<const amrex::Real> const& vector_field_x,
158  amrex::Array4<const amrex::Real> const& vector_field_y,
159  amrex::Array4<const amrex::Real> const& vector_field_z,
162 {
163  // first find the weight of surrounding nodes to use during interpolation
164  int ii, jj, kk;
165  amrex::Real W[AMREX_SPACEDIM][2];
166  compute_weights_nodal(xp, yp, zp, lo, dxi, ii, jj, kk, W);
167 
168  amrex::GpuArray<amrex::Real, 3> const field_interp = {
169  interp_field_nodal(ii, jj, kk, W, vector_field_x),
170  interp_field_nodal(ii, jj, kk, W, vector_field_y),
171  interp_field_nodal(ii, jj, kk, W, vector_field_z)
172  };
173 
174  return field_interp;
175 }
176 
177 } // namespace ablastr::particles
178 
179 #endif // ABLASTR_NODALFIELDGATHER_H_
#define AMREX_INLINE
#define AMREX_GPU_HOST_DEVICE
def y
Definition: Excitation_Flag_Generator.py:76
def z
Definition: Excitation_Flag_Generator.py:77
Definition: DepositCharge.H:22
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:92
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:32
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:128
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:154
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
void Abort(const std::string &msg)
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