WarpX
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 "ablastr/utils/TextMsg.H"
11 
12 #include <AMReX_Array.H>
13 #include <AMReX_Extension.H>
14 #include <AMReX_GpuQualifiers.H>
15 #include <AMReX_Math.H>
16 #include <AMReX_REAL.H>
17 
18 namespace ablastr::particles
19 {
33 void compute_weights_nodal (const amrex::ParticleReal xp,
34  const amrex::ParticleReal yp,
35  const amrex::ParticleReal zp,
38  int& i, int& j, int& k, amrex::Real W[AMREX_SPACEDIM][2]) noexcept
39 {
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];
45 
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));
49 
50  W[0][1] = x - i;
51  W[1][1] = y - j;
52  W[2][1] = z - k;
53 
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];
57 
58 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
59 
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));
64  W[0][1] = x - i;
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));
68  W[0][1] = r - i;
69 # endif
70 
71  const amrex::Real z = (zp - plo[1]) * dxi[1];
72  j = static_cast<int>(amrex::Math::floor(z));
73  W[1][1] = z - j;
74 
75  W[0][0] = 1.0_rt - W[0][1];
76  W[1][0] = 1.0_rt - W[1][1];
77 
78  k = 0;
79 #else
80  amrex::ignore_unused(xp, yp, zp, plo, dxi, i, j, k, W);
81  ABLASTR_ABORT_WITH_MESSAGE("Error: compute_weights not yet implemented in 1D");
82 #endif
83 }
84 
93 amrex::Real interp_field_nodal (int i, int j, int k,
94  const amrex::Real W[AMREX_SPACEDIM][2],
95  amrex::Array4<const amrex::Real> const& scalar_field) noexcept
96 {
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];
113 #else
114  const amrex::Real value = 0;
115  amrex::ignore_unused(i, j, k, W, scalar_field);
116  ABLASTR_ABORT_WITH_MESSAGE("Error: interp_field not yet implemented in 1D");
117 #endif
118  return value;
119 }
120 
131 amrex::Real doGatherScalarFieldNodal (const amrex::ParticleReal xp,
132  const amrex::ParticleReal yp,
133  const amrex::ParticleReal zp,
134  amrex::Array4<const amrex::Real> const& scalar_field,
137 {
138  // first find the weight of surrounding nodes to use during interpolation
139  int ii, jj, kk;
140  amrex::Real W[AMREX_SPACEDIM][2];
141  compute_weights_nodal(xp, yp, zp, lo, dxi, ii, jj, kk, W);
142 
143  return interp_field_nodal(ii, jj, kk, W, scalar_field);
144 }
145 
157 doGatherVectorFieldNodal (const amrex::ParticleReal xp,
158  const amrex::ParticleReal yp,
159  const amrex::ParticleReal zp,
160  amrex::Array4<const amrex::Real> const& vector_field_x,
161  amrex::Array4<const amrex::Real> const& vector_field_y,
162  amrex::Array4<const amrex::Real> const& vector_field_z,
165 {
166  // first find the weight of surrounding nodes to use during interpolation
167  int ii, jj, kk;
168  amrex::Real W[AMREX_SPACEDIM][2];
169  compute_weights_nodal(xp, yp, zp, lo, dxi, ii, jj, kk, W);
170 
171  amrex::GpuArray<amrex::Real, 3> const field_interp = {
172  interp_field_nodal(ii, jj, kk, W, vector_field_x),
173  interp_field_nodal(ii, jj, kk, W, vector_field_y),
174  interp_field_nodal(ii, jj, kk, W, vector_field_z)
175  };
176 
177  return field_interp;
178 }
179 
180 } // namespace ablastr::particles
181 
182 #endif // ABLASTR_NODALFIELDGATHER_H_
#define AMREX_INLINE
#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