7 #ifndef WARPX_COMM_K_H_
8 #define WARPX_COMM_K_H_
21 using namespace amrex;
24 const auto arr_coarse_zeropad = [arr_coarse] (
const int jj,
const int kk,
const int ll) noexcept
26 return arr_coarse.
contains(
jj,kk,ll) ? arr_coarse(
jj,kk,ll) : 0.0_rt;
33 const int rk = (AMREX_SPACEDIM == 1) ? 1 : rr[1];
34 const int rl = (AMREX_SPACEDIM <= 2) ? 1 : rr[2];
37 const int sj = arr_stag[0];
38 const int sk = (AMREX_SPACEDIM == 1) ? 0 : arr_stag[1];
39 const int sl = (AMREX_SPACEDIM <= 2) ? 0 : arr_stag[2];
56 const amrex::Real hj = (sj == 0) ? 0.5_rt : 0._rt;
57 const amrex::Real hk = (sk == 0) ? 0.5_rt : 0._rt;
58 const amrex::Real hl = (sl == 0) ? 0.5_rt : 0._rt;
60 amrex::Real res = 0.0_rt;
62 for (
int jj = 0;
jj < nj;
jj++) {
63 for (
int kk = 0; kk < nk; kk++) {
64 for (
int ll = 0; ll < nl; ll++) {
65 wj = (rj - amrex::Math::abs(j + hj - (jc +
jj + hj) * rj)) /
static_cast<amrex::Real
>(rj);
66 wk = (rk - amrex::Math::abs(k + hk - (kc + kk + hk) * rk)) /
static_cast<amrex::Real
>(rk);
67 wl = (rl - amrex::Math::abs(l + hl - (lc + ll + hl) * rl)) /
static_cast<amrex::Real
>(rl);
68 res += wj * wk * wl * arr_coarse_zeropad(jc+
jj,kc+kk,lc+ll);
72 arr_aux(j,k,l) = arr_fine(j,k,l) + res;
82 using namespace amrex;
85 const auto Bxf_zeropad = [Bxf] (
const int jj,
const int kk,
const int ll) noexcept
91 const Real wx = (j == jg*2) ? 0.0_rt : 0.5_rt;
92 const Real owx = 1.0_rt-wx;
97 wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
100 #if defined(WARPX_DIM_1D_Z)
103 const Real bg = owx * Bxg(jg ,0,0)
104 + wx * Bxg(jg+1,0,0);
107 const Real bc = owx * Bxc(jg ,0,0)
108 + wx * Bxc(jg+1,0,0);
111 const Real bf = 0.5_rt*(Bxf_zeropad(j-1,0,0) + Bxf_zeropad(j,0,0));
114 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
117 const Real bg = owx * owy * Bxg(jg ,kg ,0)
118 + owx * wy * Bxg(jg ,kg+1,0)
119 + wx * owy * Bxg(jg+1,kg ,0)
120 + wx * wy * Bxg(jg+1,kg+1,0);
123 wy = 0.5_rt-wy; owy = 1.0_rt-wy;
124 const Real bc = owx * owy * Bxc(jg ,kg ,0)
125 + owx * wy * Bxc(jg ,kg-1,0)
126 + wx * owy * Bxc(jg+1,kg ,0)
127 + wx * wy * Bxc(jg+1,kg-1,0);
130 const Real bf = 0.5_rt*(Bxf_zeropad(j,k-1,0) + Bxf_zeropad(j,k,0));
135 Real wz = (l == lg*2) ? 0.0_rt : 0.5_rt;
136 Real owz = 1.0_rt-wz;
139 const Real bg = owx * owy * owz * Bxg(jg ,kg ,lg )
140 + wx * owy * owz * Bxg(jg+1,kg ,lg )
141 + owx * wy * owz * Bxg(jg ,kg+1,lg )
142 + wx * wy * owz * Bxg(jg+1,kg+1,lg )
143 + owx * owy * wz * Bxg(jg ,kg ,lg+1)
144 + wx * owy * wz * Bxg(jg+1,kg ,lg+1)
145 + owx * wy * wz * Bxg(jg ,kg+1,lg+1)
146 + wx * wy * wz * Bxg(jg+1,kg+1,lg+1);
149 wy = 0.5_rt-wy; owy = 1.0_rt-wy;
150 wz = 0.5_rt-wz; owz = 1.0_rt-wz;
151 const Real bc = owx * owy * owz * Bxc(jg ,kg ,lg )
152 + wx * owy * owz * Bxc(jg+1,kg ,lg )
153 + owx * wy * owz * Bxc(jg ,kg-1,lg )
154 + wx * wy * owz * Bxc(jg+1,kg-1,lg )
155 + owx * owy * wz * Bxc(jg ,kg ,lg-1)
156 + wx * owy * wz * Bxc(jg+1,kg ,lg-1)
157 + owx * wy * wz * Bxc(jg ,kg-1,lg-1)
158 + wx * wy * wz * Bxc(jg+1,kg-1,lg-1);
161 const Real bf = 0.25_rt*(Bxf_zeropad(j,k-1,l-1) + Bxf_zeropad(j,k,l-1) + Bxf_zeropad(j,k-1,l) + Bxf_zeropad(j,k,l));
164 Bxa(j,k,l) = bg + (bf-bc);
174 using namespace amrex;
177 const auto Byf_zeropad = [Byf] (
const int jj,
const int kk,
const int ll) noexcept
184 wx = (j == jg*2) ? 0.0_rt : 0.5_rt;
190 wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
193 #if defined(WARPX_DIM_1D_Z)
196 const Real bg = owx * Byg(jg ,0,0)
197 + wx * Byg(jg+1,0,0);
200 const Real bc = owx * Byc(jg ,0,0)
201 + wx * Byc(jg+1,0,0);
204 const Real bf = 0.5_rt*(Byf_zeropad(j-1,0,0) + Byf_zeropad(j,0,0));
207 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
210 const Real bg = owx * owy * Byg(jg ,kg ,0)
211 + owx * wy * Byg(jg ,kg+1,0)
212 + wx * owy * Byg(jg+1,kg ,0)
213 + wx * wy * Byg(jg+1,kg+1,0);
216 wx = 0.5_rt-wx; owx = 1.0_rt-wx;
217 wy = 0.5_rt-wy; owy = 1.0_rt-wy;
218 const Real bc = owx * owy * Byc(jg ,kg ,0)
219 + owx * wy * Byc(jg ,kg-1,0)
220 + wx * owy * Byc(jg-1,kg ,0)
221 + wx * wy * Byc(jg-1,kg-1,0);
224 const Real bf = 0.25_rt*(Byf_zeropad(j,k,0) + Byf_zeropad(j-1,k,0) + Byf_zeropad(j,k-1,0) + Byf_zeropad(j-1,k-1,0));
229 Real wz = (l == lg*2) ? 0.0_rt : 0.5_rt;
230 Real owz = 1.0_rt-wz;
233 const Real bg = owx * owy * owz * Byg(jg ,kg ,lg )
234 + wx * owy * owz * Byg(jg+1,kg ,lg )
235 + owx * wy * owz * Byg(jg ,kg+1,lg )
236 + wx * wy * owz * Byg(jg+1,kg+1,lg )
237 + owx * owy * wz * Byg(jg ,kg ,lg+1)
238 + wx * owy * wz * Byg(jg+1,kg ,lg+1)
239 + owx * wy * wz * Byg(jg ,kg+1,lg+1)
240 + wx * wy * wz * Byg(jg+1,kg+1,lg+1);
243 wx = 0.5_rt-wx; owx = 1.0_rt-wx;
244 wz = 0.5_rt-wz; owz = 1.0_rt-wz;
245 const Real bc = owx * owy * owz * Byc(jg ,kg ,lg )
246 + wx * owy * owz * Byc(jg-1,kg ,lg )
247 + owx * wy * owz * Byc(jg ,kg+1,lg )
248 + wx * wy * owz * Byc(jg-1,kg+1,lg )
249 + owx * owy * wz * Byc(jg ,kg ,lg-1)
250 + wx * owy * wz * Byc(jg-1,kg ,lg-1)
251 + owx * wy * wz * Byc(jg ,kg+1,lg-1)
252 + wx * wy * wz * Byc(jg-1,kg+1,lg-1);
255 const Real bf = 0.25_rt*(Byf_zeropad(j-1,k,l-1) + Byf_zeropad(j,k,l-1) + Byf_zeropad(j-1,k,l) + Byf_zeropad(j,k,l));
259 Bya(j,k,l) = bg + (bf-bc);
269 using namespace amrex;
272 const auto Bzf_zeropad = [Bzf] (
const int jj,
const int kk,
const int ll) noexcept
279 wx = (j == jg*2) ? 0.0_rt : 0.5_rt;
285 wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
288 #if defined(WARPX_DIM_1D_Z)
291 const Real bg = owx * Bzg(jg ,0,0)
292 + wx * Bzg(jg+1,0,0);
295 const Real bc = owx * Bzc(jg ,0,0)
296 + wx * Bzc(jg+1,0,0);
299 const Real bf = Bzf_zeropad(j,0,0);
302 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
305 const Real bg = owx * owy * Bzg(jg ,kg ,0)
306 + owx * wy * Bzg(jg ,kg+1,0)
307 + wx * owy * Bzg(jg+1,kg ,0)
308 + wx * wy * Bzg(jg+1,kg+1,0);
311 wx = 0.5_rt-wx; owx = 1.0_rt-wx;
312 const Real bc = owx * owy * Bzc(jg ,kg ,0)
313 + owx * wy * Bzc(jg ,kg+1,0)
314 + wx * owy * Bzc(jg-1,kg ,0)
315 + wx * wy * Bzc(jg-1,kg+1,0);
318 const Real bf = 0.5_rt*(Bzf_zeropad(j-1,k,0) + Bzf_zeropad(j,k,0));
323 const Real wz = (l == lg*2) ? 0.0_rt : 0.5_rt;
324 const Real owz = 1.0_rt-wz;
327 const Real bg = owx * owy * owz * Bzg(jg ,kg ,lg )
328 + wx * owy * owz * Bzg(jg+1,kg ,lg )
329 + owx * wy * owz * Bzg(jg ,kg+1,lg )
330 + wx * wy * owz * Bzg(jg+1,kg+1,lg )
331 + owx * owy * wz * Bzg(jg ,kg ,lg+1)
332 + wx * owy * wz * Bzg(jg+1,kg ,lg+1)
333 + owx * wy * wz * Bzg(jg ,kg+1,lg+1)
334 + wx * wy * wz * Bzg(jg+1,kg+1,lg+1);
337 wx = 0.5_rt-wx; owx = 1.0_rt-wx;
338 wy = 0.5_rt-wy; owy = 1.0_rt-wy;
339 const Real bc = owx * owy * owz * Bzc(jg ,kg ,lg )
340 + wx * owy * owz * Bzc(jg-1,kg ,lg )
341 + owx * wy * owz * Bzc(jg ,kg-1,lg )
342 + wx * wy * owz * Bzc(jg-1,kg-1,lg )
343 + owx * owy * wz * Bzc(jg ,kg ,lg+1)
344 + wx * owy * wz * Bzc(jg-1,kg ,lg+1)
345 + owx * wy * wz * Bzc(jg ,kg-1,lg+1)
346 + wx * wy * wz * Bzc(jg-1,kg-1,lg+1);
349 const Real bf = 0.25_rt*(Bzf_zeropad(j-1,k-1,l) + Bzf_zeropad(j,k-1,l) + Bzf_zeropad(j-1,k,l) + Bzf_zeropad(j,k,l));
353 Bza(j,k,l) = bg + (bf-bc);
363 using namespace amrex;
366 const auto Exf_zeropad = [Exf] (
const int jj,
const int kk,
const int ll) noexcept
373 wx = (j == jg*2) ? 0.0_rt : 0.5_rt;
379 wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
382 #if defined(WARPX_DIM_1D_Z)
385 const Real eg = owx * Exg(jg ,0,0)
386 + wx * Exg(jg+1,0,0);
389 const Real ec = owx * Exc(jg ,0,0)
390 + wx * Exc(jg+1,0,0);
393 const Real ef = Exf_zeropad(j,0,0);
396 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
399 const Real eg = owx * owy * Exg(jg ,kg ,0)
400 + owx * wy * Exg(jg ,kg+1,0)
401 + wx * owy * Exg(jg+1,kg ,0)
402 + wx * wy * Exg(jg+1,kg+1,0);
405 wx = 0.5_rt-wx; owx = 1.0_rt-wx;
406 const Real ec = owx * owy * Exc(jg ,kg ,0)
407 + owx * wy * Exc(jg ,kg+1,0)
408 + wx * owy * Exc(jg-1,kg ,0)
409 + wx * wy * Exc(jg-1,kg+1,0);
412 const Real ef = 0.5_rt*(Exf_zeropad(j-1,k,0) + Exf_zeropad(j,k,0));
417 const Real wz = (l == lg*2) ? 0.0 : 0.5;
418 const Real owz = 1.0_rt-wz;
421 const Real eg = owx * owy * owz * Exg(jg ,kg ,lg )
422 + wx * owy * owz * Exg(jg+1,kg ,lg )
423 + owx * wy * owz * Exg(jg ,kg+1,lg )
424 + wx * wy * owz * Exg(jg+1,kg+1,lg )
425 + owx * owy * wz * Exg(jg ,kg ,lg+1)
426 + wx * owy * wz * Exg(jg+1,kg ,lg+1)
427 + owx * wy * wz * Exg(jg ,kg+1,lg+1)
428 + wx * wy * wz * Exg(jg+1,kg+1,lg+1);
431 wx = 0.5_rt-wx; owx = 1.0_rt-wx;
432 const Real ec = owx * owy * owz * Exc(jg ,kg ,lg )
433 + wx * owy * owz * Exc(jg-1,kg ,lg )
434 + owx * wy * owz * Exc(jg ,kg+1,lg )
435 + wx * wy * owz * Exc(jg-1,kg+1,lg )
436 + owx * owy * wz * Exc(jg ,kg ,lg+1)
437 + wx * owy * wz * Exc(jg-1,kg ,lg+1)
438 + owx * wy * wz * Exc(jg ,kg+1,lg+1)
439 + wx * wy * wz * Exc(jg-1,kg+1,lg+1);
442 const Real ef = 0.5_rt*(Exf_zeropad(j-1,k,l) + Exf_zeropad(j,k,l));
446 Exa(j,k,l) = eg + (ef-ec);
456 using namespace amrex;
459 const auto Eyf_zeropad = [Eyf] (
const int jj,
const int kk,
const int ll) noexcept
465 const Real wx = (j == jg*2) ? 0.0_rt : 0.5_rt;
466 const Real owx = 1.0_rt-wx;
471 wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
474 #if defined(WARPX_DIM_1D_Z)
477 const Real eg = owx * Eyg(jg ,0,0)
478 + wx * Eyg(jg+1,0,0);
481 const Real ec = owx * Eyc(jg ,0,0)
482 + wx * Eyc(jg+1,0,0);
485 const Real ef = Eyf_zeropad(j,0,0);
488 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
491 const Real eg = owx * owy * Eyg(jg ,kg ,0)
492 + owx * wy * Eyg(jg ,kg+1,0)
493 + wx * owy * Eyg(jg+1,kg ,0)
494 + wx * wy * Eyg(jg+1,kg+1,0);
497 const Real ec = owx * owy * Eyc(jg ,kg ,0)
498 + owx * wy * Eyc(jg ,kg+1,0)
499 + wx * owy * Eyc(jg+1,kg ,0)
500 + wx * wy * Eyc(jg+1,kg+1,0);
503 const Real ef = Eyf_zeropad(j,k,0);
508 const Real wz = (l == lg*2) ? 0.0 : 0.5;
509 const Real owz = 1.0_rt-wz;
512 const Real eg = owx * owy * owz * Eyg(jg ,kg ,lg )
513 + wx * owy * owz * Eyg(jg+1,kg ,lg )
514 + owx * wy * owz * Eyg(jg ,kg+1,lg )
515 + wx * wy * owz * Eyg(jg+1,kg+1,lg )
516 + owx * owy * wz * Eyg(jg ,kg ,lg+1)
517 + wx * owy * wz * Eyg(jg+1,kg ,lg+1)
518 + owx * wy * wz * Eyg(jg ,kg+1,lg+1)
519 + wx * wy * wz * Eyg(jg+1,kg+1,lg+1);
522 wy = 0.5_rt-wy; owy = 1.0_rt-wy;
523 const Real ec = owx * owy * owz * Eyc(jg ,kg ,lg )
524 + wx * owy * owz * Eyc(jg+1,kg ,lg )
525 + owx * wy * owz * Eyc(jg ,kg-1,lg )
526 + wx * wy * owz * Eyc(jg+1,kg-1,lg )
527 + owx * owy * wz * Eyc(jg ,kg ,lg+1)
528 + wx * owy * wz * Eyc(jg+1,kg ,lg+1)
529 + owx * wy * wz * Eyc(jg ,kg-1,lg+1)
530 + wx * wy * wz * Eyc(jg+1,kg-1,lg+1);
533 const Real ef = 0.5_rt*(Eyf_zeropad(j,k-1,l) + Eyf_zeropad(j,k,l));
537 Eya(j,k,l) = eg + (ef-ec);
547 using namespace amrex;
550 const auto Ezf_zeropad = [Ezf] (
const int jj,
const int kk,
const int ll) noexcept
556 const Real wx = (j == jg*2) ? 0.0_rt : 0.5_rt;
557 const Real owx = 1.0_rt-wx;
562 wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
565 #if defined(WARPX_DIM_1D_Z)
568 const Real eg = owx * Ezg(jg ,0,0)
569 + wx * Ezg(jg+1,0,0);
572 const Real ec = owx * Ezc(jg ,0,0)
573 + wx * Ezc(jg+1,0,0);
576 const Real ef = 0.5_rt*(Ezf_zeropad(j-1,0,0) + Ezf_zeropad(j,0,0));
579 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
582 const Real eg = owx * owy * Ezg(jg ,kg ,0)
583 + owx * wy * Ezg(jg ,kg+1,0)
584 + wx * owy * Ezg(jg+1,kg ,0)
585 + wx * wy * Ezg(jg+1,kg+1,0);
588 wy = 0.5_rt-wy; owy = 1.0_rt-wy;
589 const Real ec = owx * owy * Ezc(jg ,kg ,0)
590 + owx * wy * Ezc(jg ,kg-1,0)
591 + wx * owy * Ezc(jg+1,kg ,0)
592 + wx * wy * Ezc(jg+1,kg-1,0);
595 const Real ef = 0.5_rt*(Ezf_zeropad(j,k-1,0) + Ezf_zeropad(j,k,0));
600 Real wz = (l == lg*2) ? 0.0_rt : 0.5_rt;
601 Real owz = 1.0_rt-wz;
604 const Real eg = owx * owy * owz * Ezg(jg ,kg ,lg )
605 + wx * owy * owz * Ezg(jg+1,kg ,lg )
606 + owx * wy * owz * Ezg(jg ,kg+1,lg )
607 + wx * wy * owz * Ezg(jg+1,kg+1,lg )
608 + owx * owy * wz * Ezg(jg ,kg ,lg+1)
609 + wx * owy * wz * Ezg(jg+1,kg ,lg+1)
610 + owx * wy * wz * Ezg(jg ,kg+1,lg+1)
611 + wx * wy * wz * Ezg(jg+1,kg+1,lg+1);
614 wz = 0.5_rt-wz; owz = 1.0_rt-wz;
615 const Real ec = owx * owy * owz * Ezc(jg ,kg ,lg )
616 + wx * owy * owz * Ezc(jg+1,kg ,lg )
617 + owx * wy * owz * Ezc(jg ,kg+1,lg )
618 + wx * wy * owz * Ezc(jg+1,kg+1,lg )
619 + owx * owy * wz * Ezc(jg ,kg ,lg-1)
620 + wx * owy * wz * Ezc(jg+1,kg ,lg-1)
621 + owx * wy * wz * Ezc(jg ,kg+1,lg-1)
622 + wx * wy * wz * Ezc(jg+1,kg+1,lg-1);
625 const Real ef = 0.5_rt*(Ezf_zeropad(j,k,l-1) + Ezf_zeropad(j,k,l));
629 Eza(j,k,l) = eg + (ef-ec);
662 amrex::Real
const* stencil_coeffs_x =
nullptr,
663 amrex::Real
const* stencil_coeffs_y =
nullptr,
664 amrex::Real
const* stencil_coeffs_z =
nullptr)
666 using namespace amrex;
670 const auto src_arr_zeropad = [src_arr] (
const int jj,
const int kk,
const int ll) noexcept
672 return src_arr.
contains(
jj,kk,ll) ? src_arr(
jj,kk,ll) : 0.0_rt;
676 #if defined(WARPX_DIM_1D_Z)
678 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
687 const int shift = (dst_nodal) ? 0 : 1;
690 const int sj = (dst_nodal) ? src_stag[0] : dst_stag[0];
691 #if (AMREX_SPACEDIM >= 2)
692 const int sk = (dst_nodal) ? src_stag[1] : dst_stag[1];
694 #if defined(WARPX_DIM_3D)
695 const int sl = (dst_nodal) ? src_stag[2] : dst_stag[2];
699 const bool interp_j = (sj == 0);
700 #if (AMREX_SPACEDIM >= 2)
701 const bool interp_k = (sk == 0);
703 #if defined(WARPX_DIM_3D)
704 const bool interp_l = (sl == 0);
707 #if defined(WARPX_DIM_1D_Z)
709 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
712 #elif defined(WARPX_DIM_3D)
719 const amrex::Real wj = (interp_j) ? 0.5_rt : 1.0_rt;
720 #if defined(WARPX_DIM_1D_Z)
721 constexpr amrex::Real wk = 1.0_rt;
722 constexpr amrex::Real wl = 1.0_rt;
723 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
724 const amrex::Real wk = (interp_k) ? 0.5_rt : 1.0_rt;
725 constexpr amrex::Real wl = 1.0_rt;
726 #elif defined(WARPX_DIM_3D)
727 const amrex::Real wk = (interp_k) ? 0.5_rt : 1.0_rt;
728 const amrex::Real wl = (interp_l) ? 0.5_rt : 1.0_rt;
732 const int jmin = (interp_j) ? j - noj/2 +
shift : j;
733 const int jmax = (interp_j) ? j + noj/2 +
shift - 1 : j;
736 #if defined(WARPX_DIM_1D_Z)
741 const int kmin = (interp_k) ? k - nok/2 +
shift : k;
742 const int kmax = (interp_k) ? k + nok/2 +
shift - 1 : k;
746 #if (AMREX_SPACEDIM <= 2)
750 #elif defined(WARPX_DIM_3D)
751 const int lmin = (interp_l) ? l - nol/2 +
shift : l;
752 const int lmax = (interp_l) ? l + nol/2 +
shift - 1 : l;
756 const int nj = jmax - jmin;
757 const int nk = kmax - kmin;
758 const int nl = lmax - lmin;
799 amrex::Real res = 0.0_rt;
801 #if defined(WARPX_DIM_1D_Z)
802 amrex::Real
const* scj = stencil_coeffs_z;
803 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
804 amrex::Real
const* scj = stencil_coeffs_x;
805 amrex::Real
const* sck = stencil_coeffs_z;
806 #elif defined(WARPX_DIM_3D)
807 amrex::Real
const* scj = stencil_coeffs_x;
808 amrex::Real
const* sck = stencil_coeffs_y;
809 amrex::Real
const* scl = stencil_coeffs_z;
812 for (
int ll = 0; ll <= nl; ll++)
814 #if defined(WARPX_DIM_3D)
815 const amrex::Real cl = (interp_l)? scl[ll] : 1.0_rt;
817 const amrex::Real cl = 1.0_rt;
819 for (
int kk = 0; kk <= nk; kk++)
821 #if (AMREX_SPACEDIM >= 2)
822 const amrex::Real ck = (interp_k)? sck[kk] : 1.0_rt;
824 const amrex::Real ck = 1.0_rt;
826 for (
int jj = 0;
jj <= nj;
jj++)
828 const amrex::Real cj = (interp_j)? scj[
jj] : 1.0_rt;
830 res += cj * ck * cl * src_arr_zeropad(jmin+
jj,kmin+kk,lmin+ll);
835 dst_arr(j,k,l) = wj * wk * wl * res;
#define AMREX_FORCE_INLINE
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void warpx_interp_nd_efield_x(int j, int k, int l, amrex::Array4< amrex::Real > const &Exa, amrex::Array4< amrex::Real const > const &Exf, amrex::Array4< amrex::Real const > const &Exc, amrex::Array4< amrex::Real const > const &Exg)
Definition: WarpXComm_K.H:357
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void warpx_interp_nd_bfield_x(int j, int k, int l, amrex::Array4< amrex::Real > const &Bxa, amrex::Array4< amrex::Real const > const &Bxf, amrex::Array4< amrex::Real const > const &Bxc, amrex::Array4< amrex::Real const > const &Bxg)
Definition: WarpXComm_K.H:76
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void warpx_interp_nd_efield_y(int j, int k, int l, amrex::Array4< amrex::Real > const &Eya, amrex::Array4< amrex::Real const > const &Eyf, amrex::Array4< amrex::Real const > const &Eyc, amrex::Array4< amrex::Real const > const &Eyg)
Definition: WarpXComm_K.H:450
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)
Definition: WarpXComm_K.H:14
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void warpx_interp_nd_bfield_y(int j, int k, int l, amrex::Array4< amrex::Real > const &Bya, amrex::Array4< amrex::Real const > const &Byf, amrex::Array4< amrex::Real const > const &Byc, amrex::Array4< amrex::Real const > const &Byg)
Definition: WarpXComm_K.H:168
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void warpx_interp_nd_bfield_z(int j, int k, int l, amrex::Array4< amrex::Real > const &Bza, amrex::Array4< amrex::Real const > const &Bzf, amrex::Array4< amrex::Real const > const &Bzc, amrex::Array4< amrex::Real const > const &Bzg)
Definition: WarpXComm_K.H:263
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void warpx_interp_nd_efield_z(int j, int k, int l, amrex::Array4< amrex::Real > const &Eza, amrex::Array4< amrex::Real const > const &Ezf, amrex::Array4< amrex::Real const > const &Ezc, amrex::Array4< amrex::Real const > const &Ezg)
Definition: WarpXComm_K.H:541
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
int nox
Definition: stencil.py:442
int noy
Definition: stencil.py:443
int noz
Definition: stencil.py:444
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool contains(int i, int j, int k) const noexcept