WarpX
WarpXComm_K.H
Go to the documentation of this file.
1 /* Copyright 2019 Axel Huebl, Weiqun Zhang
2  *
3  * This file is part of WarpX.
4  *
5  * License: BSD-3-Clause-LBNL
6  */
7 #ifndef WARPX_COMM_K_H_
8 #define WARPX_COMM_K_H_
9 
10 #include <AMReX.H>
11 #include <AMReX_FArrayBox.H>
12 
28 void warpx_interp (int j, int k, int l,
29  amrex::Array4<amrex::Real > const& arr_aux,
30  amrex::Array4<amrex::Real const> const& arr_fine,
31  amrex::Array4<amrex::Real const> const& arr_coarse,
32  const amrex::IntVect& arr_stag,
33  const amrex::IntVect& rr)
34 {
35  using namespace amrex;
36 
37  // Pad arr_coarse with zeros beyond ghost cells for out-of-bound accesses
38  const auto arr_coarse_zeropad = [arr_coarse] (const int jj, const int kk, const int ll) noexcept
39  {
40  return arr_coarse.contains(jj,kk,ll) ? arr_coarse(jj,kk,ll) : 0.0_rt;
41  };
42 
43  // NOTE Indices (j,k,l) in the following refer to (z,-,-) in 1D, (x,z,-) in 2D, and (x,y,z) in 3D
44 
45  // Refinement ratio
46  const int rj = rr[0];
47  const int rk = (AMREX_SPACEDIM == 1) ? 1 : rr[1];
48  const int rl = (AMREX_SPACEDIM <= 2) ? 1 : rr[2];
49 
50  // Staggering (0: cell-centered; 1: nodal)
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];
54 
55  // Number of points used for interpolation from coarse grid to fine grid
56  const int nj = 2;
57  const int nk = 2;
58  const int nl = 2;
59 
60  const int jc = (sj == 0) ? amrex::coarsen(j - rj/2, rj) : amrex::coarsen(j, rj);
61  const int kc = (sk == 0) ? amrex::coarsen(k - rk/2, rk) : amrex::coarsen(k, rk);
62  const int lc = (sl == 0) ? amrex::coarsen(l - rl/2, rl) : amrex::coarsen(l, rl);
63 
64  // Interpolate from coarse grid to fine grid using 2 points
65  // with weights depending on the distance, for both nodal and cell-centered grids
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;
69 
70  amrex::Real res = 0.0_rt;
71 
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);
79  }
80  }
81  }
82  arr_aux(j,k,l) = arr_fine(j,k,l) + res;
83 }
84 
102 void warpx_interp (int j, int k, int l,
103  amrex::Array4<amrex::Real > const& arr_aux,
104  amrex::Array4<amrex::Real const> const& arr_fine,
105  amrex::Array4<amrex::Real const> const& arr_coarse,
106  amrex::Array4<amrex::Real const> const& arr_tmp,
107  const amrex::IntVect& arr_fine_stag,
108  const amrex::IntVect& arr_coarse_stag,
109  const amrex::IntVect& rr)
110 {
111  using namespace amrex;
112 
113  // Pad input arrays with zeros beyond ghost cells
114  // for out-of-bound accesses due to large-stencil operations
115  const auto arr_fine_zeropad = [arr_fine] (const int jj, const int kk, const int ll) noexcept
116  {
117  return arr_fine.contains(jj,kk,ll) ? arr_fine(jj,kk,ll) : 0.0_rt;
118  };
119  const auto arr_coarse_zeropad = [arr_coarse] (const int jj, const int kk, const int ll) noexcept
120  {
121  return arr_coarse.contains(jj,kk,ll) ? arr_coarse(jj,kk,ll) : 0.0_rt;
122  };
123  const auto arr_tmp_zeropad = [arr_tmp] (const int jj, const int kk, const int ll) noexcept
124  {
125  return arr_tmp.contains(jj,kk,ll) ? arr_tmp(jj,kk,ll) : 0.0_rt;
126  };
127 
128  // NOTE Indices (j,k,l) in the following refer to:
129  // - (z,-,-) in 1D
130  // - (x,z,-) in 2D
131  // - (r,z,-) in RZ
132  // - (x,y,z) in 3D
133 
134  // Refinement ratio
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;
142 #else
143  const int rk = rr[1];
144  const int rl = rr[2];
145 #endif
146 
147  // Staggering of fine array (0: cell-centered; 1: nodal)
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];
154 #endif
155 
156  // Staggering of coarse array (0: cell-centered; 1: nodal)
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;
164 #else
165  const int sk_cp = arr_coarse_stag[1];
166  const int sl_cp = arr_coarse_stag[2];
167 #endif
168 
169  // Number of points used for interpolation from coarse grid to fine grid
170  int nj;
171  int nk;
172  int nl;
173 
174  int jc = amrex::coarsen(j, rj);
175  int kc = amrex::coarsen(k, rk);
176  int lc = amrex::coarsen(l, rl);
177 
178  amrex::Real tmp = 0.0_rt;
179  amrex::Real fine = 0.0_rt;
180  amrex::Real coarse = 0.0_rt;
181 
182  amrex::Real wj;
183  amrex::Real wk;
184  amrex::Real wl;
185 
186  // 1) Interpolation from coarse nodal to fine nodal
187 
188  nj = 2;
189 #if defined(WARPX_DIM_1D_Z)
190  nk = 1;
191  nl = 1;
192 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
193  nk = 2;
194  nl = 1;
195 #else
196  nk = 2;
197  nl = 2;
198 #endif
199 
200  wj = 1.0_rt;
201  wk = 1.0_rt;
202  wl = 1.0_rt;
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);
209 #endif
210 #if (AMREX_SPACEDIM == 3)
211  wl = (rl - amrex::Math::abs(l - (lc + ll) * rl)) / static_cast<amrex::Real>(rl);
212 #endif
213  tmp += wj * wk * wl * arr_tmp_zeropad(jc+jj,kc+kk,lc+ll);
214  }
215  }
216  }
217 
218  // 2) Interpolation from coarse staggered to fine nodal
219 
220  nj = 2;
221 #if defined(WARPX_DIM_1D_Z)
222  nk = 1;
223  nl = 1;
224 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
225  nk = 2;
226  nl = 1;
227 #else
228  nk = 2;
229  nl = 2;
230 #endif
231 
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;
235 
236  jc = amrex::coarsen(jn, rj);
237  kc = amrex::coarsen(kn, rk);
238  lc = amrex::coarsen(ln, rl);
239 
240  wj = 1.0_rt;
241  wk = 1.0_rt;
242  wl = 1.0_rt;
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);
249 #endif
250 #if (AMREX_SPACEDIM == 3)
251  wl = (rl - amrex::Math::abs(ln - (lc + ll) * rl)) / static_cast<amrex::Real>(rl);
252 #endif
253  coarse += wj * wk * wl * arr_coarse_zeropad(jc+jj,kc+kk,lc+ll);
254  }
255  }
256  }
257 
258  // 3) Interpolation from fine staggered to fine nodal
259 
260  nj = (sj_fp == 0) ? 2 : 1;
261 #if defined(WARPX_DIM_1D_Z)
262  nk = 1;
263  nl = 1;
264 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
265  nk = (sk_fp == 0) ? 2 : 1;
266  nl = 1;
267 #else
268  nk = (sk_fp == 0) ? 2 : 1;
269  nl = (sl_fp == 0) ? 2 : 1;
270 #endif
271 
272  const int jm = (sj_fp == 0) ? j-1 : j;
273 #if defined(WARPX_DIM_1D_Z)
274  const int km = k;
275  const int lm = l;
276 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
277  const int km = (sk_fp == 0) ? k-1 : k;
278  const int lm = l;
279 #else
280  const int km = (sk_fp == 0) ? k-1 : k;
281  const int lm = (sl_fp == 0) ? l-1 : l;
282 #endif
283 
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);
291  }
292  }
293  }
294 
295  // Final result
296  arr_aux(j,k,l) = tmp + (fine - coarse);
297 }
298 
319 void warpx_interp (const int j,
320  const int k,
321  const int l,
322  amrex::Array4<amrex::Real > const& dst_arr,
323  amrex::Array4<amrex::Real const> const& src_arr,
324  const amrex::IntVect& dst_stag,
325  const amrex::IntVect& src_stag,
326  const int nox = 2,
327  const int noy = 2,
328  const int noz = 2,
329  amrex::Real const* stencil_coeffs_x = nullptr,
330  amrex::Real const* stencil_coeffs_y = nullptr,
331  amrex::Real const* stencil_coeffs_z = nullptr)
332 {
333  using namespace amrex;
334 
335  // Pad input array with zeros beyond ghost cells
336  // for out-of-bound accesses due to large-stencil operations
337  const auto src_arr_zeropad = [src_arr] (const int jj, const int kk, const int ll) noexcept
338  {
339  return src_arr.contains(jj,kk,ll) ? src_arr(jj,kk,ll) : 0.0_rt;
340  };
341 
342  // Avoid compiler warnings
343 #if defined(WARPX_DIM_1D_Z)
344  amrex::ignore_unused(nox, noy, stencil_coeffs_x, stencil_coeffs_y);
345 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
346  amrex::ignore_unused(noy, stencil_coeffs_y);
347 #endif
348 
349  // If dst_nodal = true , we are centering from a staggered grid to a nodal grid
350  // If dst_nodal = false, we are centering from a nodal grid to a staggered grid
351  const bool dst_nodal = (dst_stag == amrex::IntVect::TheNodeVector());
352 
353  // See 1D examples below to understand the meaning of this integer shift
354  const int shift = (dst_nodal) ? 0 : 1;
355 
356  // Staggering (s = 0 if cell-centered, s = 1 if nodal)
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];
360 #endif
361 #if defined(WARPX_DIM_3D)
362  const int sl = (dst_nodal) ? src_stag[2] : dst_stag[2];
363 #endif
364 
365  // Interpolate along j,k,l only if source MultiFab is staggered along j,k,l
366  const bool interp_j = (sj == 0);
367 #if (AMREX_SPACEDIM >= 2)
368  const bool interp_k = (sk == 0);
369 #endif
370 #if defined(WARPX_DIM_3D)
371  const bool interp_l = (sl == 0);
372 #endif
373 
374 #if defined(WARPX_DIM_1D_Z)
375  const int noj = noz;
376 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
377  const int noj = nox;
378  const int nok = noz;
379 #elif defined(WARPX_DIM_3D)
380  const int noj = nox;
381  const int nok = noy;
382  const int nol = noz;
383 #endif
384 
385  // Additional normalization factor
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;
396 #endif
397 
398  // Min and max for interpolation loop along j
399  const int jmin = (interp_j) ? j - noj/2 + shift : j;
400  const int jmax = (interp_j) ? j + noj/2 + shift - 1 : j;
401 
402  // Min and max for interpolation loop along k
403 #if defined(WARPX_DIM_1D_Z)
404  // k = 0 always
405  const int kmin = k;
406  const int kmax = k;
407 #else
408  const int kmin = (interp_k) ? k - nok/2 + shift : k;
409  const int kmax = (interp_k) ? k + nok/2 + shift - 1 : k;
410 #endif
411 
412  // Min and max for interpolation loop along l
413 #if (AMREX_SPACEDIM <= 2)
414  // l = 0 always
415  const int lmin = l;
416  const int lmax = l;
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;
420 #endif
421 
422  // Number of interpolation points
423  const int nj = jmax - jmin;
424  const int nk = kmax - kmin;
425  const int nl = lmax - lmin;
426 
427  // Example of 1D centering from nodal grid to nodal grid (simple copy):
428  //
429  // j
430  // --o-----o-----o-- result(j) = f(j)
431  // --o-----o-----o--
432  // j-1 j j+1
433  //
434  // Example of 1D linear centering from staggered grid to nodal grid:
435  //
436  // j
437  // --o-----o-----o-- result(j) = (f(j-1) + f(j)) / 2
438  // -----x-----x-----
439  // j-1 j
440  //
441  // Example of 1D linear centering from nodal grid to staggered grid:
442  // (note the shift of +1 in the indices with respect to the case above, see variable "shift")
443  //
444  // j
445  // --x-----x-----x-- result(j) = (f(j) + f(j+1)) / 2
446  // -----o-----o-----
447  // j j+1
448  //
449  // Example of 1D finite-order centering from staggered grid to nodal grid:
450  //
451  // j
452  // --o-----o-----o-----o-----o-----o-----o-- result(j) = c_0 * (f(j-1) + f(j) ) / 2
453  // -----x-----x-----x-----x-----x-----x----- + c_1 * (f(j-2) + f(j+1)) / 2
454  // j-3 j-2 j-1 j j+1 j+2 + c_2 * (f(j-3) + f(j+2)) / 2
455  // c_2 c_1 c_0 c_0 c_1 c_2 + ...
456  //
457  // Example of 1D finite-order centering from nodal grid to staggered grid:
458  // (note the shift of +1 in the indices with respect to the case above, see variable "shift")
459  //
460  // j
461  // --x-----x-----x-----x-----x-----x-----x-- result(j) = c_0 * (f(j) + f(j+1)) / 2
462  // -----o-----o-----o-----o-----o-----o----- + c_1 * (f(j-1) + f(j+2)) / 2
463  // j-2 j-1 j j+1 j+2 j+3 + c_2 * (f(j-2) + f(j+3)) / 2
464  // c_2 c_1 c_0 c_0 c_1 c_2 + ...
465 
466  amrex::Real res = 0.0_rt;
467 
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;
477 #endif
478 
479  for (int ll = 0; ll <= nl; ll++)
480  {
481 #if defined(WARPX_DIM_3D)
482  const amrex::Real cl = (interp_l)? scl[ll] : 1.0_rt;
483 #else
484  const amrex::Real cl = 1.0_rt;
485 #endif
486  for (int kk = 0; kk <= nk; kk++)
487  {
488 #if (AMREX_SPACEDIM >= 2)
489  const amrex::Real ck = (interp_k)? sck[kk] : 1.0_rt;
490 #else
491  const amrex::Real ck = 1.0_rt;
492 #endif
493  for (int jj = 0; jj <= nj; jj++)
494  {
495  const amrex::Real cj = (interp_j)? scj[jj] : 1.0_rt;
496 
497  res += cj * ck * cl * src_arr_zeropad(jmin+jj,kmin+kk,lmin+ll);
498  }
499  }
500  }
501 
502  dst_arr(j,k,l) = wj * wk * wl * res;
503 }
504 
505 #endif
#define AMREX_FORCE_INLINE
#define AMREX_GPU_DEVICE
Array4< Real > fine
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