7 #ifndef WARPX_COMM_K_H_
8 #define WARPX_COMM_K_H_
35 using namespace amrex;
38 const auto arr_coarse_zeropad = [arr_coarse] (
const int jj,
const int kk,
const int ll) noexcept
40 return arr_coarse.
contains(
jj,kk,ll) ? arr_coarse(
jj,kk,ll) : 0.0_rt;
47 const int rk = (AMREX_SPACEDIM == 1) ? 1 : rr[1];
48 const int rl = (AMREX_SPACEDIM <= 2) ? 1 : rr[2];
51 const int sj = arr_stag[0];
52 const int sk = (AMREX_SPACEDIM == 1) ? 0 : arr_stag[1];
53 const int sl = (AMREX_SPACEDIM <= 2) ? 0 : arr_stag[2];
66 const amrex::Real hj = (sj == 0) ? 0.5_rt : 0._rt;
67 const amrex::Real hk = (sk == 0) ? 0.5_rt : 0._rt;
68 const amrex::Real hl = (sl == 0) ? 0.5_rt : 0._rt;
70 amrex::Real res = 0.0_rt;
72 for (
int jj = 0;
jj < nj;
jj++) {
73 for (
int kk = 0; kk < nk; kk++) {
74 for (
int ll = 0; ll < nl; ll++) {
75 const amrex::Real wj = (rj - amrex::Math::abs(j + hj - (jc +
jj + hj) * rj)) /
static_cast<amrex::Real
>(rj);
76 const amrex::Real wk = (rk - amrex::Math::abs(k + hk - (kc + kk + hk) * rk)) /
static_cast<amrex::Real
>(rk);
77 const amrex::Real wl = (rl - amrex::Math::abs(l + hl - (lc + ll + hl) * rl)) /
static_cast<amrex::Real
>(rl);
78 res += wj * wk * wl * arr_coarse_zeropad(jc+
jj,kc+kk,lc+ll);
82 arr_aux(j,k,l) = arr_fine(j,k,l) + res;
111 using namespace amrex;
115 const auto arr_fine_zeropad = [arr_fine] (
const int jj,
const int kk,
const int ll) noexcept
117 return arr_fine.
contains(
jj,kk,ll) ? arr_fine(
jj,kk,ll) : 0.0_rt;
119 const auto arr_coarse_zeropad = [arr_coarse] (
const int jj,
const int kk,
const int ll) noexcept
121 return arr_coarse.
contains(
jj,kk,ll) ? arr_coarse(
jj,kk,ll) : 0.0_rt;
123 const auto arr_tmp_zeropad = [arr_tmp] (
const int jj,
const int kk,
const int ll) noexcept
125 return arr_tmp.
contains(
jj,kk,ll) ? arr_tmp(
jj,kk,ll) : 0.0_rt;
135 const int rj = rr[0];
136 #if defined(WARPX_DIM_1D_Z)
137 constexpr
int rk = 1;
138 constexpr
int rl = 1;
139 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
140 const int rk = rr[1];
141 constexpr
int rl = 1;
143 const int rk = rr[1];
144 const int rl = rr[2];
148 const int sj_fp = arr_fine_stag[0];
149 #if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
150 const int sk_fp = arr_fine_stag[1];
151 #elif defined(WARPX_DIM_3D)
152 const int sk_fp = arr_fine_stag[1];
153 const int sl_fp = arr_fine_stag[2];
157 const int sj_cp = arr_coarse_stag[0];
158 #if defined(WARPX_DIM_1D_Z)
159 constexpr
int sk_cp = 0;
160 constexpr
int sl_cp = 0;
161 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
162 const int sk_cp = arr_coarse_stag[1];
163 constexpr
int sl_cp = 0;
165 const int sk_cp = arr_coarse_stag[1];
166 const int sl_cp = arr_coarse_stag[2];
178 amrex::Real tmp = 0.0_rt;
179 amrex::Real
fine = 0.0_rt;
180 amrex::Real
coarse = 0.0_rt;
189 #if defined(WARPX_DIM_1D_Z)
192 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
203 for (
int jj = 0;
jj < nj;
jj++) {
204 for (
int kk = 0; kk < nk; kk++) {
205 for (
int ll = 0; ll < nl; ll++) {
206 wj = (rj - amrex::Math::abs(j - (jc +
jj) * rj)) /
static_cast<amrex::Real
>(rj);
207 #if (AMREX_SPACEDIM >= 2)
208 wk = (rk - amrex::Math::abs(k - (kc + kk) * rk)) /
static_cast<amrex::Real
>(rk);
210 #if (AMREX_SPACEDIM == 3)
211 wl = (rl - amrex::Math::abs(l - (lc + ll) * rl)) /
static_cast<amrex::Real
>(rl);
213 tmp += wj * wk * wl * arr_tmp_zeropad(jc+
jj,kc+kk,lc+ll);
221 #if defined(WARPX_DIM_1D_Z)
224 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
232 const int jn = (sj_cp == 1) ? j : j - rj / 2;
233 const int kn = (sk_cp == 1) ? k : k - rk / 2;
234 const int ln = (sl_cp == 1) ? l : l - rl / 2;
243 for (
int jj = 0;
jj < nj;
jj++) {
244 for (
int kk = 0; kk < nk; kk++) {
245 for (
int ll = 0; ll < nl; ll++) {
246 wj = (rj - amrex::Math::abs(jn - (jc +
jj) * rj)) /
static_cast<amrex::Real
>(rj);
247 #if (AMREX_SPACEDIM >= 2)
248 wk = (rk - amrex::Math::abs(kn - (kc + kk) * rk)) /
static_cast<amrex::Real
>(rk);
250 #if (AMREX_SPACEDIM == 3)
251 wl = (rl - amrex::Math::abs(ln - (lc + ll) * rl)) /
static_cast<amrex::Real
>(rl);
253 coarse += wj * wk * wl * arr_coarse_zeropad(jc+
jj,kc+kk,lc+ll);
260 nj = (sj_fp == 0) ? 2 : 1;
261 #if defined(WARPX_DIM_1D_Z)
264 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
265 nk = (sk_fp == 0) ? 2 : 1;
268 nk = (sk_fp == 0) ? 2 : 1;
269 nl = (sl_fp == 0) ? 2 : 1;
272 const int jm = (sj_fp == 0) ? j-1 : j;
273 #if defined(WARPX_DIM_1D_Z)
276 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
277 const int km = (sk_fp == 0) ? k-1 : k;
280 const int km = (sk_fp == 0) ? k-1 : k;
281 const int lm = (sl_fp == 0) ? l-1 : l;
284 for (
int jj = 0;
jj < nj;
jj++) {
285 for (
int kk = 0; kk < nk; kk++) {
286 for (
int ll = 0; ll < nl; ll++) {
287 wj = 1.0_rt /
static_cast<amrex::Real
>(nj);
288 wk = 1.0_rt /
static_cast<amrex::Real
>(nk);
289 wl = 1.0_rt /
static_cast<amrex::Real
>(nl);
290 fine += wj * wk * wl * arr_fine_zeropad(jm+
jj,km+kk,lm+ll);
329 amrex::Real
const* stencil_coeffs_x =
nullptr,
330 amrex::Real
const* stencil_coeffs_y =
nullptr,
331 amrex::Real
const* stencil_coeffs_z =
nullptr)
333 using namespace amrex;
337 const auto src_arr_zeropad = [src_arr] (
const int jj,
const int kk,
const int ll) noexcept
339 return src_arr.
contains(
jj,kk,ll) ? src_arr(
jj,kk,ll) : 0.0_rt;
343 #if defined(WARPX_DIM_1D_Z)
345 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
354 const int shift = (dst_nodal) ? 0 : 1;
357 const int sj = (dst_nodal) ? src_stag[0] : dst_stag[0];
358 #if (AMREX_SPACEDIM >= 2)
359 const int sk = (dst_nodal) ? src_stag[1] : dst_stag[1];
361 #if defined(WARPX_DIM_3D)
362 const int sl = (dst_nodal) ? src_stag[2] : dst_stag[2];
366 const bool interp_j = (sj == 0);
367 #if (AMREX_SPACEDIM >= 2)
368 const bool interp_k = (sk == 0);
370 #if defined(WARPX_DIM_3D)
371 const bool interp_l = (sl == 0);
374 #if defined(WARPX_DIM_1D_Z)
376 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
379 #elif defined(WARPX_DIM_3D)
386 const amrex::Real wj = (interp_j) ? 0.5_rt : 1.0_rt;
387 #if defined(WARPX_DIM_1D_Z)
388 constexpr amrex::Real wk = 1.0_rt;
389 constexpr amrex::Real wl = 1.0_rt;
390 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
391 const amrex::Real wk = (interp_k) ? 0.5_rt : 1.0_rt;
392 constexpr amrex::Real wl = 1.0_rt;
393 #elif defined(WARPX_DIM_3D)
394 const amrex::Real wk = (interp_k) ? 0.5_rt : 1.0_rt;
395 const amrex::Real wl = (interp_l) ? 0.5_rt : 1.0_rt;
399 const int jmin = (interp_j) ? j - noj/2 +
shift : j;
400 const int jmax = (interp_j) ? j + noj/2 +
shift - 1 : j;
403 #if defined(WARPX_DIM_1D_Z)
408 const int kmin = (interp_k) ? k - nok/2 +
shift : k;
409 const int kmax = (interp_k) ? k + nok/2 +
shift - 1 : k;
413 #if (AMREX_SPACEDIM <= 2)
417 #elif defined(WARPX_DIM_3D)
418 const int lmin = (interp_l) ? l - nol/2 +
shift : l;
419 const int lmax = (interp_l) ? l + nol/2 +
shift - 1 : l;
423 const int nj = jmax - jmin;
424 const int nk = kmax - kmin;
425 const int nl = lmax - lmin;
466 amrex::Real res = 0.0_rt;
468 #if defined(WARPX_DIM_1D_Z)
469 amrex::Real
const* scj = stencil_coeffs_z;
470 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
471 amrex::Real
const* scj = stencil_coeffs_x;
472 amrex::Real
const* sck = stencil_coeffs_z;
473 #elif defined(WARPX_DIM_3D)
474 amrex::Real
const* scj = stencil_coeffs_x;
475 amrex::Real
const* sck = stencil_coeffs_y;
476 amrex::Real
const* scl = stencil_coeffs_z;
479 for (
int ll = 0; ll <= nl; ll++)
481 #if defined(WARPX_DIM_3D)
482 const amrex::Real cl = (interp_l)? scl[ll] : 1.0_rt;
484 const amrex::Real cl = 1.0_rt;
486 for (
int kk = 0; kk <= nk; kk++)
488 #if (AMREX_SPACEDIM >= 2)
489 const amrex::Real ck = (interp_k)? sck[kk] : 1.0_rt;
491 const amrex::Real ck = 1.0_rt;
493 for (
int jj = 0;
jj <= nj;
jj++)
495 const amrex::Real cj = (interp_j)? scj[
jj] : 1.0_rt;
497 res += cj * ck * cl * src_arr_zeropad(jmin+
jj,kmin+kk,lmin+ll);
502 dst_arr(j,k,l) = wj * wk * wl * res;
#define AMREX_FORCE_INLINE
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void warpx_interp(int j, int k, int l, amrex::Array4< amrex::Real > const &arr_aux, amrex::Array4< amrex::Real const > const &arr_fine, amrex::Array4< amrex::Real const > const &arr_coarse, const amrex::IntVect &arr_stag, const amrex::IntVect &rr)
Interpolation function called within WarpX::UpdateAuxilaryDataSameType to interpolate data from the c...
Definition: WarpXComm_K.H:28
AMREX_GPU_HOST_DEVICE static constexpr AMREX_FORCE_INLINE IntVect TheNodeVector() noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Box coarsen(const Box &b, int ref_ratio) noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Box shift(const Box &b, int dir, int nzones) noexcept
jj
Definition: check_interp_points_and_weights.py:160
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool contains(int i, int j, int k) const noexcept