7 #ifndef WARPX_COMM_K_H_ 8 #define WARPX_COMM_K_H_ 10 #include <AMReX_FArrayBox.H> 13 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
15 amrex::Array4<amrex::Real >
const& Bxa,
16 amrex::Array4<amrex::Real const>
const& Bxf,
17 amrex::Array4<amrex::Real const>
const& Bxc)
19 using namespace amrex;
21 int const lg = amrex::coarsen(l,2);
22 int const kg = amrex::coarsen(k,2);
23 int const jg = amrex::coarsen(j,2);
25 Real
const wx = (j == jg*2) ? 0.0_rt : 0.5_rt;
26 Real
const owx = 1.0_rt - wx;
27 Bxa(j,k,l) = owx * Bxc(jg,kg,lg) + wx * Bxc(jg+1,kg,lg) + Bxf(j,k,l);
30 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
32 amrex::Array4<amrex::Real >
const& Bya,
33 amrex::Array4<amrex::Real const>
const& Byf,
34 amrex::Array4<amrex::Real const>
const& Byc)
36 using namespace amrex;
38 int const lg = amrex::coarsen(l,2);
39 int const kg = amrex::coarsen(k,2);
40 int const jg = amrex::coarsen(j,2);
44 #if (AMREX_SPACEDIM == 3) 45 Real
const wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
46 Real
const owy = 1.0_rt - wy;
47 Bya(j,k,l) = owy * Byc(jg,kg,lg) + wy * Byc(jg,kg+1,lg) + Byf(j,k,l);
49 Bya(j,k,l) = Byc(jg,kg,lg) + Byf(j,k,l);
53 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
55 amrex::Array4<amrex::Real >
const& Bza,
56 amrex::Array4<amrex::Real const>
const& Bzf,
57 amrex::Array4<amrex::Real const>
const& Bzc)
59 using namespace amrex;
61 int const lg = amrex::coarsen(l,2);
62 int const kg = amrex::coarsen(k,2);
63 int const jg = amrex::coarsen(j,2);
67 #if (AMREX_SPACEDIM == 3) 68 Real
const wz = (l == lg*2) ? 0.0_rt : 0.5_rt;
69 Real
const owz = 1.0_rt - wz;
70 Bza(j,k,l) = owz * Bzc(jg,kg,lg) + owz * Bzc(jg,kg,lg+1) + Bzf(j,k,l);
72 Real
const wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
73 Real
const owy = 1.0_rt - wy;
74 Bza(j,k,l) = owy * Bzc(jg,kg,lg) + owy * Bzc(jg,kg+1,lg) + Bzf(j,k,l);
78 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
80 amrex::Array4<amrex::Real >
const& Exa,
81 amrex::Array4<amrex::Real const>
const& Exf,
82 amrex::Array4<amrex::Real const>
const& Exc)
84 using namespace amrex;
86 int const lg = amrex::coarsen(l,2);
87 int const kg = amrex::coarsen(k,2);
88 int const jg = amrex::coarsen(j,2);
90 Real
const wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
91 Real
const owy = 1.0_rt - wy;
93 #if (AMREX_SPACEDIM == 3) 94 Real
const wz = (l == lg*2) ? 0.0_rt : 0.5_rt;
95 Real
const owz = 1.0_rt - wz;
96 Exa(j,k,l) = owy * owz * Exc(jg ,kg ,lg )
97 + wy * owz * Exc(jg ,kg+1,lg )
98 + owy * wz * Exc(jg ,kg ,lg+1)
99 + wy * wz * Exc(jg ,kg+1,lg+1)
102 Exa(j,k,l) = owy * Exc(jg,kg,lg) + wy * Exc(jg,kg+1,lg) + Exf(j,k,l);
106 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
108 amrex::Array4<amrex::Real >
const& Eya,
109 amrex::Array4<amrex::Real const>
const& Eyf,
110 amrex::Array4<amrex::Real const>
const& Eyc)
112 using namespace amrex;
114 int const lg = amrex::coarsen(l,2);
115 int const kg = amrex::coarsen(k,2);
116 int const jg = amrex::coarsen(j,2);
118 Real
const wx = (j == jg*2) ? 0.0_rt : 0.5_rt;
119 Real
const owx = 1.0_rt - wx;
121 #if (AMREX_SPACEDIM == 3) 122 Real
const wz = (l == lg*2) ? 0.0_rt : 0.5_rt;
123 Real
const owz = 1.0_rt - wz;
124 Eya(j,k,l) = owx * owz * Eyc(jg ,kg ,lg )
125 + wx * owz * Eyc(jg+1,kg ,lg )
126 + owx * wz * Eyc(jg ,kg ,lg+1)
127 + wx * wz * Eyc(jg+1,kg ,lg+1)
130 Real
const wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
131 Real
const owy = 1.0_rt - wy;
132 Eya(j,k,l) = owx * owy * Eyc(jg ,kg ,lg)
133 + wx * owy * Eyc(jg+1,kg ,lg)
134 + owx * wy * Eyc(jg ,kg+1,lg)
135 + wx * wy * Eyc(jg+1,kg+1,lg)
140 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
142 amrex::Array4<amrex::Real >
const& Eza,
143 amrex::Array4<amrex::Real const>
const& Ezf,
144 amrex::Array4<amrex::Real const>
const& Ezc)
146 using namespace amrex;
148 int const lg = amrex::coarsen(l,2);
149 int const kg = amrex::coarsen(k,2);
150 int const jg = amrex::coarsen(j,2);
152 Real
const wx = (j == jg*2) ? 0.0_rt : 0.5_rt;
153 Real
const owx = 1.0_rt - wx;
155 #if (AMREX_SPACEDIM == 3) 156 Real wy = (k == kg*2) ? 0.0_rt : 0.5_rt;
157 Real owy = 1.0_rt - wy;
158 Eza(j,k,l) = owx * owy * Ezc(jg ,kg ,lg )
159 + wx * owy * Ezc(jg+1,kg ,lg )
160 + owx * wy * Ezc(jg ,kg+1,lg )
161 + wx * wy * Ezc(jg+1,kg+1,lg )
164 Eza(j,k,l) = owx * Ezc(jg,kg,lg) + wx * Ezc(jg+1,kg,lg) + Ezf(j,k,l);
168 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
170 amrex::Array4<amrex::Real>
const& Bxa,
171 amrex::Array4<amrex::Real const>
const& Bxf,
172 amrex::Array4<amrex::Real const>
const& Bxc,
173 amrex::Array4<amrex::Real const>
const& Bxg)
175 using namespace amrex;
177 int jg = amrex::coarsen(j,2);
178 Real wx = (j == jg*2) ? 0.0 : 0.5;
181 int kg = amrex::coarsen(k,2);
182 Real wy = (k == kg*2) ? 0.0 : 0.5;
185 #if (AMREX_SPACEDIM == 2) 188 Real bg = owx * owy * Bxg(jg ,kg ,0)
189 + owx * wy * Bxg(jg ,kg+1,0)
190 + wx * owy * Bxg(jg+1,kg ,0)
191 + wx * wy * Bxg(jg+1,kg+1,0);
194 wy = 0.5-wy; owy = 1.0-wy;
195 Real bc = owx * owy * Bxc(jg ,kg ,0)
196 + owx * wy * Bxc(jg ,kg-1,0)
197 + wx * owy * Bxc(jg+1,kg ,0)
198 + wx * wy * Bxc(jg+1,kg-1,0);
201 Real bf = 0.5*(Bxf(j,k-1,0) + Bxf(j,k,0));
205 int lg = amrex::coarsen(l,2);
206 Real wz = (l == lg*2) ? 0.0 : 0.5;
210 Real bg = owx * owy * owz * Bxg(jg ,kg ,lg )
211 + wx * owy * owz * Bxg(jg+1,kg ,lg )
212 + owx * wy * owz * Bxg(jg ,kg+1,lg )
213 + wx * wy * owz * Bxg(jg+1,kg+1,lg )
214 + owx * owy * wz * Bxg(jg ,kg ,lg+1)
215 + wx * owy * wz * Bxg(jg+1,kg ,lg+1)
216 + owx * wy * wz * Bxg(jg ,kg+1,lg+1)
217 + wx * wy * wz * Bxg(jg+1,kg+1,lg+1);
220 wy = 0.5-wy; owy = 1.0-wy;
221 wz = 0.5-wz; owz = 1.0-wz;
222 Real bc = owx * owy * owz * Bxc(jg ,kg ,lg )
223 + wx * owy * owz * Bxc(jg+1,kg ,lg )
224 + owx * wy * owz * Bxc(jg ,kg-1,lg )
225 + wx * wy * owz * Bxc(jg+1,kg-1,lg )
226 + owx * owy * wz * Bxc(jg ,kg ,lg-1)
227 + wx * owy * wz * Bxc(jg+1,kg ,lg-1)
228 + owx * wy * wz * Bxc(jg ,kg-1,lg-1)
229 + wx * wy * wz * Bxc(jg+1,kg-1,lg-1);
232 Real bf = 0.25*(Bxf(j,k-1,l-1) + Bxf(j,k,l-1) + Bxf(j,k-1,l) + Bxf(j,k,l));
235 Bxa(j,k,l) = bg + (bf-bc);
238 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
240 amrex::Array4<amrex::Real>
const& Bya,
241 amrex::Array4<amrex::Real const>
const& Byf,
242 amrex::Array4<amrex::Real const>
const& Byc,
243 amrex::Array4<amrex::Real const>
const& Byg)
245 using namespace amrex;
247 int jg = amrex::coarsen(j,2);
248 Real wx = (j == jg*2) ? 0.0 : 0.5;
251 int kg = amrex::coarsen(k,2);
252 Real wy = (k == kg*2) ? 0.0 : 0.5;
255 #if (AMREX_SPACEDIM == 2) 258 Real bg = owx * owy * Byg(jg ,kg ,0)
259 + owx * wy * Byg(jg ,kg+1,0)
260 + wx * owy * Byg(jg+1,kg ,0)
261 + wx * wy * Byg(jg+1,kg+1,0);
264 wx = 0.5-wx; owx = 1.0-wx;
265 wy = 0.5-wy; owy = 1.0-wy;
266 Real bc = owx * owy * Byc(jg ,kg ,0)
267 + owx * wy * Byc(jg ,kg-1,0)
268 + wx * owy * Byc(jg-1,kg ,0)
269 + wx * wy * Byc(jg-1,kg-1,0);
272 Real bf = 0.25*(Byf(j,k,0) + Byf(j-1,k,0) + Byf(j,k-1,0) + Byf(j-1,k-1,0));
276 int lg = amrex::coarsen(l,2);
277 Real wz = (l == lg*2) ? 0.0 : 0.5;
281 Real bg = owx * owy * owz * Byg(jg ,kg ,lg )
282 + wx * owy * owz * Byg(jg+1,kg ,lg )
283 + owx * wy * owz * Byg(jg ,kg+1,lg )
284 + wx * wy * owz * Byg(jg+1,kg+1,lg )
285 + owx * owy * wz * Byg(jg ,kg ,lg+1)
286 + wx * owy * wz * Byg(jg+1,kg ,lg+1)
287 + owx * wy * wz * Byg(jg ,kg+1,lg+1)
288 + wx * wy * wz * Byg(jg+1,kg+1,lg+1);
291 wx = 0.5-wx; owx = 1.0-wx;
292 wz = 0.5-wz; owz = 1.0-wz;
293 Real bc = owx * owy * owz * Byc(jg ,kg ,lg )
294 + wx * owy * owz * Byc(jg-1,kg ,lg )
295 + owx * wy * owz * Byc(jg ,kg+1,lg )
296 + wx * wy * owz * Byc(jg-1,kg+1,lg )
297 + owx * owy * wz * Byc(jg ,kg ,lg-1)
298 + wx * owy * wz * Byc(jg-1,kg ,lg-1)
299 + owx * wy * wz * Byc(jg ,kg+1,lg-1)
300 + wx * wy * wz * Byc(jg-1,kg+1,lg-1);
303 Real bf = 0.25*(Byf(j-1,k,l-1) + Byf(j,k,l-1) + Byf(j-1,k,l) + Byf(j,k,l));
307 Bya(j,k,l) = bg + (bf-bc);
310 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
312 amrex::Array4<amrex::Real>
const& Bza,
313 amrex::Array4<amrex::Real const>
const& Bzf,
314 amrex::Array4<amrex::Real const>
const& Bzc,
315 amrex::Array4<amrex::Real const>
const& Bzg)
317 using namespace amrex;
319 int jg = amrex::coarsen(j,2);
320 Real wx = (j == jg*2) ? 0.0 : 0.5;
323 int kg = amrex::coarsen(k,2);
324 Real wy = (k == kg*2) ? 0.0 : 0.5;
327 #if (AMREX_SPACEDIM == 2) 330 Real bg = owx * owy * Bzg(jg ,kg ,0)
331 + owx * wy * Bzg(jg ,kg+1,0)
332 + wx * owy * Bzg(jg+1,kg ,0)
333 + wx * wy * Bzg(jg+1,kg+1,0);
336 wx = 0.5-wx; owx = 1.0-wx;
337 Real bc = owx * owy * Bzc(jg ,kg ,0)
338 + owx * wy * Bzc(jg ,kg+1,0)
339 + wx * owy * Bzc(jg-1,kg ,0)
340 + wx * wy * Bzc(jg-1,kg+1,0);
343 Real bf = 0.5*(Bzf(j-1,k,0) + Bzf(j,k,0));
347 int lg = amrex::coarsen(l,2);
348 Real wz = (l == lg*2) ? 0.0 : 0.5;
352 Real bg = owx * owy * owz * Bzg(jg ,kg ,lg )
353 + wx * owy * owz * Bzg(jg+1,kg ,lg )
354 + owx * wy * owz * Bzg(jg ,kg+1,lg )
355 + wx * wy * owz * Bzg(jg+1,kg+1,lg )
356 + owx * owy * wz * Bzg(jg ,kg ,lg+1)
357 + wx * owy * wz * Bzg(jg+1,kg ,lg+1)
358 + owx * wy * wz * Bzg(jg ,kg+1,lg+1)
359 + wx * wy * wz * Bzg(jg+1,kg+1,lg+1);
362 wx = 0.5-wx; owx = 1.0-wx;
363 wy = 0.5-wy; owy = 1.0-wy;
364 Real bc = owx * owy * owz * Bzc(jg ,kg ,lg )
365 + wx * owy * owz * Bzc(jg-1,kg ,lg )
366 + owx * wy * owz * Bzc(jg ,kg-1,lg )
367 + wx * wy * owz * Bzc(jg-1,kg-1,lg )
368 + owx * owy * wz * Bzc(jg ,kg ,lg+1)
369 + wx * owy * wz * Bzc(jg-1,kg ,lg+1)
370 + owx * wy * wz * Bzc(jg ,kg-1,lg+1)
371 + wx * wy * wz * Bzc(jg-1,kg-1,lg+1);
374 Real bf = 0.25*(Bzf(j-1,k-1,l) + Bzf(j,k-1,l) + Bzf(j-1,k,l) + Bzf(j,k,l));
378 Bza(j,k,l) = bg + (bf-bc);
381 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
383 amrex::Array4<amrex::Real>
const& Bxa,
384 amrex::Array4<amrex::Real const>
const& Bxf)
386 #if (AMREX_SPACEDIM == 2) 387 Bxa(j,k,0) = 0.5*(Bxf(j,k-1,0) + Bxf(j,k,0));
388 amrex::ignore_unused(l);
390 Bxa(j,k,l) = 0.25*(Bxf(j,k-1,l-1) + Bxf(j,k,l-1) + Bxf(j,k-1,l) + Bxf(j,k,l));
394 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
396 amrex::Array4<amrex::Real>
const& Bya,
397 amrex::Array4<amrex::Real const>
const& Byf)
399 #if (AMREX_SPACEDIM == 2) 400 Bya(j,k,0) = 0.25*(Byf(j,k,0) + Byf(j-1,k,0) + Byf(j,k-1,0) + Byf(j-1,k-1,0));
401 amrex::ignore_unused(l);
403 Bya(j,k,l) = 0.25*(Byf(j-1,k,l-1) + Byf(j,k,l-1) + Byf(j-1,k,l) + Byf(j,k,l));
407 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
409 amrex::Array4<amrex::Real>
const& Bza,
410 amrex::Array4<amrex::Real const>
const& Bzf)
412 #if (AMREX_SPACEDIM == 2) 413 Bza(j,k,0) = 0.5*(Bzf(j-1,k,0) + Bzf(j,k,0));
414 amrex::ignore_unused(l);
416 Bza(j,k,l) = 0.25*(Bzf(j-1,k-1,l) + Bzf(j,k-1,l) + Bzf(j-1,k,l) + Bzf(j,k,l));
420 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
422 amrex::Array4<amrex::Real>
const& Exa,
423 amrex::Array4<amrex::Real const>
const& Exf,
424 amrex::Array4<amrex::Real const>
const& Exc,
425 amrex::Array4<amrex::Real const>
const& Exg)
427 using namespace amrex;
429 int jg = amrex::coarsen(j,2);
430 Real wx = (j == jg*2) ? 0.0 : 0.5;
433 int kg = amrex::coarsen(k,2);
434 Real wy = (k == kg*2) ? 0.0 : 0.5;
437 #if (AMREX_SPACEDIM == 2) 440 Real eg = owx * owy * Exg(jg ,kg ,0)
441 + owx * wy * Exg(jg ,kg+1,0)
442 + wx * owy * Exg(jg+1,kg ,0)
443 + wx * wy * Exg(jg+1,kg+1,0);
446 wx = 0.5-wx; owx = 1.0-wx;
447 Real ec = owx * owy * Exc(jg ,kg ,0)
448 + owx * wy * Exc(jg ,kg+1,0)
449 + wx * owy * Exc(jg-1,kg ,0)
450 + wx * wy * Exc(jg-1,kg+1,0);
453 Real ef = 0.5*(Exf(j-1,k,0) + Exf(j,k,0));
457 int lg = amrex::coarsen(l,2);
458 Real wz = (l == lg*2) ? 0.0 : 0.5;
462 Real eg = owx * owy * owz * Exg(jg ,kg ,lg )
463 + wx * owy * owz * Exg(jg+1,kg ,lg )
464 + owx * wy * owz * Exg(jg ,kg+1,lg )
465 + wx * wy * owz * Exg(jg+1,kg+1,lg )
466 + owx * owy * wz * Exg(jg ,kg ,lg+1)
467 + wx * owy * wz * Exg(jg+1,kg ,lg+1)
468 + owx * wy * wz * Exg(jg ,kg+1,lg+1)
469 + wx * wy * wz * Exg(jg+1,kg+1,lg+1);
472 wx = 0.5-wx; owx = 1.0-wx;
473 Real ec = owx * owy * owz * Exc(jg ,kg ,lg )
474 + wx * owy * owz * Exc(jg-1,kg ,lg )
475 + owx * wy * owz * Exc(jg ,kg+1,lg )
476 + wx * wy * owz * Exc(jg-1,kg+1,lg )
477 + owx * owy * wz * Exc(jg ,kg ,lg+1)
478 + wx * owy * wz * Exc(jg-1,kg ,lg+1)
479 + owx * wy * wz * Exc(jg ,kg+1,lg+1)
480 + wx * wy * wz * Exc(jg-1,kg+1,lg+1);
483 Real ef = 0.5*(Exf(j-1,k,l) + Exf(j,k,l));
487 Exa(j,k,l) = eg + (ef-ec);
490 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
492 amrex::Array4<amrex::Real>
const& Eya,
493 amrex::Array4<amrex::Real const>
const& Eyf,
494 amrex::Array4<amrex::Real const>
const& Eyc,
495 amrex::Array4<amrex::Real const>
const& Eyg)
497 using namespace amrex;
499 int jg = amrex::coarsen(j,2);
500 Real wx = (j == jg*2) ? 0.0 : 0.5;
503 int kg = amrex::coarsen(k,2);
504 Real wy = (k == kg*2) ? 0.0 : 0.5;
507 #if (AMREX_SPACEDIM == 2) 510 Real eg = owx * owy * (Eyg(jg ,kg ,0) + Eyc(jg ,kg ,0))
511 + owx * wy * (Eyg(jg ,kg+1,0) + Eyc(jg ,kg+1,0))
512 + wx * owy * (Eyg(jg+1,kg ,0) + Eyc(jg+1,kg ,0))
513 + wx * wy * (Eyg(jg+1,kg+1,0) + Eyc(jg+1,kg+1,0));
517 Real ef = Eyf(j,k,0);
521 int lg = amrex::coarsen(l,2);
522 Real wz = (l == lg*2) ? 0.0 : 0.5;
526 Real eg = owx * owy * owz * Eyg(jg ,kg ,lg )
527 + wx * owy * owz * Eyg(jg+1,kg ,lg )
528 + owx * wy * owz * Eyg(jg ,kg+1,lg )
529 + wx * wy * owz * Eyg(jg+1,kg+1,lg )
530 + owx * owy * wz * Eyg(jg ,kg ,lg+1)
531 + wx * owy * wz * Eyg(jg+1,kg ,lg+1)
532 + owx * wy * wz * Eyg(jg ,kg+1,lg+1)
533 + wx * wy * wz * Eyg(jg+1,kg+1,lg+1);
536 wy = 0.5-wy; owy = 1.0-wy;
537 Real ec = owx * owy * owz * Eyc(jg ,kg ,lg )
538 + wx * owy * owz * Eyc(jg+1,kg ,lg )
539 + owx * wy * owz * Eyc(jg ,kg-1,lg )
540 + wx * wy * owz * Eyc(jg+1,kg-1,lg )
541 + owx * owy * wz * Eyc(jg ,kg ,lg+1)
542 + wx * owy * wz * Eyc(jg+1,kg ,lg+1)
543 + owx * wy * wz * Eyc(jg ,kg-1,lg+1)
544 + wx * wy * wz * Eyc(jg+1,kg-1,lg+1);
547 Real ef = 0.5*(Eyf(j,k-1,l) + Eyf(j,k,l));
551 Eya(j,k,l) = eg + (ef-ec);
554 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
556 amrex::Array4<amrex::Real>
const& Eza,
557 amrex::Array4<amrex::Real const>
const& Ezf,
558 amrex::Array4<amrex::Real const>
const& Ezc,
559 amrex::Array4<amrex::Real const>
const& Ezg)
561 using namespace amrex;
563 int jg = amrex::coarsen(j,2);
564 Real wx = (j == jg*2) ? 0.0 : 0.5;
567 int kg = amrex::coarsen(k,2);
568 Real wy = (k == kg*2) ? 0.0 : 0.5;
571 #if (AMREX_SPACEDIM == 2) 574 Real eg = owx * owy * Ezg(jg ,kg ,0)
575 + owx * wy * Ezg(jg ,kg+1,0)
576 + wx * owy * Ezg(jg+1,kg ,0)
577 + wx * wy * Ezg(jg+1,kg+1,0);
580 wy = 0.5-wy; owy = 1.0-wy;
581 Real ec = owx * owy * Ezc(jg ,kg ,0)
582 + owx * wy * Ezc(jg ,kg-1,0)
583 + wx * owy * Ezc(jg+1,kg ,0)
584 + wx * wy * Ezc(jg+1,kg-1,0);
587 Real ef = 0.5*(Ezf(j,k-1,0) + Ezf(j,k,0));
591 int lg = amrex::coarsen(l,2);
592 Real wz = (l == lg*2) ? 0.0 : 0.5;
596 Real eg = owx * owy * owz * Ezg(jg ,kg ,lg )
597 + wx * owy * owz * Ezg(jg+1,kg ,lg )
598 + owx * wy * owz * Ezg(jg ,kg+1,lg )
599 + wx * wy * owz * Ezg(jg+1,kg+1,lg )
600 + owx * owy * wz * Ezg(jg ,kg ,lg+1)
601 + wx * owy * wz * Ezg(jg+1,kg ,lg+1)
602 + owx * wy * wz * Ezg(jg ,kg+1,lg+1)
603 + wx * wy * wz * Ezg(jg+1,kg+1,lg+1);
606 wz = 0.5-wz; owz = 1.0-wz;
607 Real ec = owx * owy * owz * Ezc(jg ,kg ,lg )
608 + wx * owy * owz * Ezc(jg+1,kg ,lg )
609 + owx * wy * owz * Ezc(jg ,kg+1,lg )
610 + wx * wy * owz * Ezc(jg+1,kg+1,lg )
611 + owx * owy * wz * Ezc(jg ,kg ,lg-1)
612 + wx * owy * wz * Ezc(jg+1,kg ,lg-1)
613 + owx * wy * wz * Ezc(jg ,kg+1,lg-1)
614 + wx * wy * wz * Ezc(jg+1,kg+1,lg-1);
617 Real ef = 0.5*(Ezf(j,k,l-1) + Ezf(j,k,l));
621 Eza(j,k,l) = eg + (ef-ec);
624 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
626 amrex::Array4<amrex::Real>
const& Exa,
627 amrex::Array4<amrex::Real const>
const& Exf)
629 Exa(j,k,l) = 0.5*(Exf(j-1,k,l) + Exf(j,k,l));
632 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
634 amrex::Array4<amrex::Real>
const& Eya,
635 amrex::Array4<amrex::Real const>
const& Eyf)
637 #if (AMREX_SPACEDIM == 2) 638 Eya(j,k,0) = Eyf(j,k,0);
639 amrex::ignore_unused(l);
641 Eya(j,k,l) = 0.5*(Eyf(j,k-1,l) + Eyf(j,k,l));
645 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
647 amrex::Array4<amrex::Real>
const& Eza,
648 amrex::Array4<amrex::Real const>
const& Ezf)
650 #if (AMREX_SPACEDIM == 2) 651 Eza(j,k,0) = 0.5*(Ezf(j,k-1,0) + Ezf(j,k,0));
652 amrex::ignore_unused(l);
654 Eza(j,k,l) = 0.5*(Ezf(j,k,l-1) + Ezf(j,k,l));
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void warpx_interp_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)
Definition: WarpXComm_K.H:107
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:491
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void warpx_interp_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)
Definition: WarpXComm_K.H:54
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void warpx_interp_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)
Definition: WarpXComm_K.H:31
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void warpx_interp_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)
Definition: WarpXComm_K.H:79
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void warpx_interp_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)
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:239
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:169
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:421
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void warpx_interp_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)
Definition: WarpXComm_K.H:141
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:555
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:311