WarpX
Loading...
Searching...
No Matches
FieldPoyntingFlux.H
Go to the documentation of this file.
1/* Copyright 2019-2020
2 *
3 * This file is part of WarpX.
4 *
5 * License: BSD-3-Clause-LBNL
6 */
7
8#ifndef WARPX_DIAGNOSTICS_REDUCEDDIAGS_FIELDPOYTINGFLUX_H_
9#define WARPX_DIAGNOSTICS_REDUCEDDIAGS_FIELDPOYTINGFLUX_H_
10
11#include "ReducedDiags.H"
12
13#include <AMReX_Array4.H>
14#include <AMReX_GpuControl.H>
15#include <AMReX_Box.H>
16#include <AMReX_REAL.H>
17#include <AMReX_Reduce.H>
18
19#include <string>
20
22 // This requires the arrays to have Yee centering
23 // The procedure (in 3D) is to interpolate the B (which is centered on faces normal
24 // to the plane of the calculaton) to the edge in the plane of the calculation to match
25 // the location of the E. Then the product E*B is interpolated to the face center in
26 // the plane of the calculation.
27 // The procedure is mapped to lower dimensions, and in some cases the second interpolation
28 // is not needed.
29
31 static amrex::Real EyBx(int i, int j, int k, int comp,
33 amrex::Array4<const amrex::Real> const & Bx_arr) {
34 using namespace amrex::literals;
35#if defined(WARPX_DIM_3D)
36 const amrex::Real EyBx_dn = Ey_arr(i,j,k,comp)*0.5_rt*(Bx_arr(i,j,k-1,comp) + Bx_arr(i,j,k,comp));
37 const amrex::Real EyBx_up = Ey_arr(i+1,j,k,comp)*0.5_rt*(Bx_arr(i+1,j,k-1,comp) + Bx_arr(i+1,j,k,comp));
38 return 0.5_rt*(EyBx_dn + EyBx_up);
39#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
40 const amrex::Real EyBx_dn = Ey_arr(i,j,k,comp)*0.5_rt*(Bx_arr(i,j-1,k,comp) + Bx_arr(i,j,k,comp));
41 const amrex::Real EyBx_up = Ey_arr(i+1,j,k,comp)*0.5_rt*(Bx_arr(i+1,j-1,k,comp) + Bx_arr(i+1,j,k,comp));
42 return 0.5_rt*(EyBx_dn + EyBx_up);
43#elif defined(WARPX_DIM_1D_Z)
44 return Ey_arr(i,j,k,comp)*0.5_rt*(Bx_arr(i-1,j,k,comp) + Bx_arr(i,j,k,comp));
45#elif defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
46 amrex::ignore_unused(i, j, k, comp, Ey_arr, Bx_arr);
47 return 0._rt;
48#endif
49 }
50
52 static amrex::Real ExBy(int i, int j, int k, int comp,
54 amrex::Array4<const amrex::Real> const & By_arr) {
55 using namespace amrex::literals;
56#if defined(WARPX_DIM_3D)
57 const amrex::Real ExBy_dn = Ex_arr(i,j,k,comp)*0.5_rt*(By_arr(i,j,k-1,comp) + By_arr(i,j,k,comp));
58 const amrex::Real ExBy_up = Ex_arr(i,j+1,k,comp)*0.5_rt*(By_arr(i,j+1,k-1,comp) + By_arr(i,j+1,k,comp));
59 return 0.5_rt*(ExBy_dn + ExBy_up);
60#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) || defined(WARPX_DIM_1D_Z)
61 return Ex_arr(i,j,k,comp)*0.5_rt*(By_arr(i-1,j,k,comp) + By_arr(i,j,k,comp));
62#elif defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
63 amrex::ignore_unused(i, j, k, comp, Ex_arr, By_arr);
64 return 0._rt;
65#endif
66 }
67
69 static amrex::Real EzBx(int i, int j, int k, int comp,
71 amrex::Array4<const amrex::Real> const & Bx_arr) {
72 using namespace amrex::literals;
73#if defined(WARPX_DIM_3D)
74 const amrex::Real EzBx_dn = Ez_arr(i,j,k,comp)*0.5_rt*(Bx_arr(i,j-1,k,comp) + Bx_arr(i,j,k,comp));
75 const amrex::Real EzBx_up = Ez_arr(i+1,j,k,comp)*0.5_rt*(Bx_arr(i+1,j-1,k,comp) + Bx_arr(i+1,j,k,comp));
76 return 0.5_rt*(EzBx_dn + EzBx_up);
77#else
78 amrex::ignore_unused(i, j, k, comp, Ez_arr, Bx_arr);
79 return 0._rt;
80#endif
81 }
82
84 static amrex::Real ExBz(int i, int j, int k, int comp,
86 amrex::Array4<const amrex::Real> const & Bz_arr) {
87 using namespace amrex::literals;
88#if defined(WARPX_DIM_3D)
89 const amrex::Real ExBz_dn = Ex_arr(i,j,k,comp)*0.5_rt*(Bz_arr(i,j-1,k,comp) + Bz_arr(i,j,k,comp));
90 const amrex::Real ExBz_up = Ex_arr(i,j,k+1,comp)*0.5_rt*(Bz_arr(i,j-1,k+1,comp) + Bz_arr(i,j,k+1,comp));
91 return 0.5_rt*(ExBz_dn + ExBz_up);
92#else
93 amrex::ignore_unused(i, j, k, comp, Ex_arr, Bz_arr);
94 return 0._rt;
95#endif
96 }
97
99 static amrex::Real EyBz(int i, int j, int k, int comp,
101 amrex::Array4<const amrex::Real> const & Bz_arr) {
102 using namespace amrex::literals;
103#if defined(WARPX_DIM_3D)
104 const amrex::Real EyBz_dn = Ey_arr(i,j,k,comp)*0.5_rt*(Bz_arr(i-1,j,k,comp) + Bz_arr(i,j,k,comp));
105 const amrex::Real EyBz_up = Ey_arr(i,j,k+1,comp)*0.5_rt*(Bz_arr(i-1,j,k+1,comp) + Bz_arr(i,j,k+1,comp));
106 return 0.5_rt*(EyBz_dn + EyBz_up);
107#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
108 const amrex::Real EyBz_dn = Ey_arr(i,j,k,comp)*0.5_rt*(Bz_arr(i-1,j,k,comp) + Bz_arr(i,j,k,comp));
109 const amrex::Real EyBz_up = Ey_arr(i,j+1,k,comp)*0.5_rt*(Bz_arr(i-1,j+1,k,comp) + Bz_arr(i,j+1,k,comp));
110 return 0.5_rt*(EyBz_dn + EyBz_up);
111#elif defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
112 return Ey_arr(i,j,k,comp)*0.5_rt*(Bz_arr(i-1,j,k,comp) + Bz_arr(i,j,k,comp));
113#elif defined(WARPX_DIM_1D_Z)
114 amrex::ignore_unused(i, j, k, comp, Ey_arr, Bz_arr);
115 return 0._rt;
116#endif
117 }
118
120 static amrex::Real EzBy(int i, int j, int k, int comp,
122 amrex::Array4<const amrex::Real> const & By_arr) {
123 using namespace amrex::literals;
124#if defined(WARPX_DIM_3D)
125 const amrex::Real EzBy_dn = Ez_arr(i,j,k,comp)*0.5_rt*(By_arr(i-1,j,k,comp) + By_arr(i,j,k,comp));
126 const amrex::Real EzBy_up = Ez_arr(i,j+1,k,comp)*0.5_rt*(By_arr(i-1,j+1,k,comp) + By_arr(i,j+1,k,comp));
127 return 0.5_rt*(EzBy_dn + EzBy_up);
128#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
129 return Ez_arr(i,j,k,comp)*0.5_rt*(By_arr(i-1,j,k,comp) + By_arr(i,j,k,comp));
130#elif defined(WARPX_DIM_1D_Z)
131 amrex::ignore_unused(i, j, k, comp, Ez_arr, By_arr);
132 return 0._rt;
133#endif
134 }
135
136};
137
139
141 static amrex::Real EyBx(int i, int j, int k, int comp,
143 amrex::Array4<const amrex::Real> const & Bx_arr) {
144 using namespace amrex::literals;
145#if defined(WARPX_DIM_3D) || defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) || defined(WARPX_DIM_1D_Z)
146 return 0.5_rt*(Ey_arr(i,j,k-1,comp)*Bx_arr(i,j,k-1,comp)
147 + Ey_arr(i,j,k,comp)*Bx_arr(i,j,k,comp));
148#elif defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
149 amrex::ignore_unused(i, j, k, comp, Ey_arr, Bx_arr);
150 return 0._rt;
151#endif
152 }
153
155 static amrex::Real ExBy(int i, int j, int k, int comp,
157 amrex::Array4<const amrex::Real> const & By_arr) {
158 using namespace amrex::literals;
159#if defined(WARPX_DIM_3D) || defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) || defined(WARPX_DIM_1D_Z)
160 return 0.5_rt*(Ex_arr(i,j,k-1,comp)*By_arr(i,j,k-1,comp)
161 + Ex_arr(i,j,k,comp)*By_arr(i,j,k,comp));
162#elif defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
163 amrex::ignore_unused(i, j, k, comp, Ex_arr, By_arr);
164 return 0._rt;
165#endif
166 }
167
169 static amrex::Real EzBx(int i, int j, int k, int comp,
171 amrex::Array4<const amrex::Real> const & Bx_arr) {
172 using namespace amrex::literals;
173#if defined(WARPX_DIM_3D)
174 return 0.5_rt*(Ez_arr(i,j-1,k,comp)*Bx_arr(i,j-1,k,comp)
175 + Ez_arr(i,j,k,comp)*Bx_arr(i,j,k,comp));
176#else
177 amrex::ignore_unused(i, j, k, comp, Ez_arr, Bx_arr);
178 return 0._rt;
179#endif
180 }
181
183 static amrex::Real ExBz(int i, int j, int k, int comp,
185 amrex::Array4<const amrex::Real> const & Bz_arr) {
186 using namespace amrex::literals;
187#if defined(WARPX_DIM_3D)
188 return 0.5_rt*(Ex_arr(i,j-1,k,comp)*Bz_arr(i,j-1,k,comp)
189 + Ex_arr(i,j,k,comp)*Bz_arr(i,j,k,comp));
190#else
191 amrex::ignore_unused(i, j, k, comp, Ex_arr, Bz_arr);
192 return 0._rt;
193#endif
194 }
195
197 static amrex::Real EyBz(int i, int j, int k, int comp,
199 amrex::Array4<const amrex::Real> const & Bz_arr) {
200 using namespace amrex::literals;
201#if defined(WARPX_DIM_3D) || defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
202 return 0.5_rt*(Ey_arr(i-1,j,k,comp)*Bz_arr(i-1,j,k,comp)
203 + Ey_arr(i,j,k,comp)*Bz_arr(i,j,k,comp));
204#elif defined(WARPX_DIM_1D_Z)
205 amrex::ignore_unused(i, j, k, comp, Ey_arr, Bz_arr);
206 return 0._rt;
207#endif
208 }
209
211 static amrex::Real EzBy(int i, int j, int k, int comp,
213 amrex::Array4<const amrex::Real> const & By_arr) {
214 using namespace amrex::literals;
215#if defined(WARPX_DIM_3D) || defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
216 return 0.5_rt*(Ez_arr(i-1,j,k,comp)*By_arr(i-1,j,k,comp)
217 + Ez_arr(i,j,k,comp)*By_arr(i,j,k,comp));
218#elif defined(WARPX_DIM_1D_Z)
219 amrex::ignore_unused(i, j, k, comp, Ez_arr, By_arr);
220 return 0._rt;
221#endif
222 }
223
224};
225
227
229 static amrex::Real EyBx(int i, int j, int k, int comp,
231 amrex::Array4<const amrex::Real> const & Bx_arr) {
232 using namespace amrex::literals;
233#if defined(WARPX_DIM_3D)
234 return 0.25_rt*(Ey_arr(i,j,k,comp)*Bx_arr(i,j,k,comp)
235 + Ey_arr(i+1,j,k,comp)*Bx_arr(i+1,j,k,comp)
236 + Ey_arr(i,j+1,k,comp)*Bx_arr(i,j+1,k,comp)
237 + Ey_arr(i+1,j+1,k,comp)*Bx_arr(i+1,j+1,k,comp));
238#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
239 return 0.5_rt*(Ey_arr(i,j,k,comp)*Bx_arr(i,j,k,comp) +
240 Ey_arr(i+1,j,k,comp)*Bx_arr(i+1,j,k,comp));
241#elif defined(WARPX_DIM_1D_Z)
242 return Ey_arr(i,j,k,comp)*Bx_arr(i,j,k,comp);
243#elif defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
244 amrex::ignore_unused(i, j, k, comp, Ey_arr, Bx_arr);
245 return 0._rt;
246#endif
247 }
248
250 static amrex::Real ExBy(int i, int j, int k, int comp,
252 amrex::Array4<const amrex::Real> const & By_arr) {
253 using namespace amrex::literals;
254#if defined(WARPX_DIM_3D)
255 return 0.25_rt*(Ex_arr(i,j,k,comp)*By_arr(i,j,k,comp)
256 + Ex_arr(i+1,j,k,comp)*By_arr(i+1,j,k,comp)
257 + Ex_arr(i,j+1,k,comp)*By_arr(i,j+1,k,comp)
258 + Ex_arr(i+1,j+1,k,comp)*By_arr(i+1,j+1,k,comp));
259#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
260 return 0.5_rt*(Ex_arr(i,j,k,comp)*By_arr(i,j,k,comp)
261 + Ex_arr(i+1,j,k,comp)*By_arr(i+1,j,k,comp));
262#elif defined(WARPX_DIM_1D_Z)
263 return Ex_arr(i,j,k,comp)*By_arr(i,j,k,comp);
264#elif defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
265 amrex::ignore_unused(i, j, k, comp, Ex_arr, By_arr);
266 return 0._rt;
267#endif
268 }
269
271 static amrex::Real EzBx(int i, int j, int k, int comp,
273 amrex::Array4<const amrex::Real> const & Bx_arr) {
274 using namespace amrex::literals;
275#if defined(WARPX_DIM_3D)
276 return 0.25_rt*(Ez_arr(i,j,k,comp)*Bx_arr(i,j,k,comp)
277 + Ez_arr(i+1,j,k,comp)*Bx_arr(i+1,j,k,comp)
278 + Ez_arr(i,j,k+1,comp)*Bx_arr(i,j,k+1,comp)
279 + Ez_arr(i+1,j,k+1,comp)*Bx_arr(i+1,j,k+1,comp));
280#else
281 amrex::ignore_unused(i, j, k, comp, Ez_arr, Bx_arr);
282 return 0._rt;
283#endif
284 }
285
287 static amrex::Real ExBz(int i, int j, int k, int comp,
289 amrex::Array4<const amrex::Real> const & Bz_arr) {
290 using namespace amrex::literals;
291#if defined(WARPX_DIM_3D)
292 return 0.25_rt*(Ex_arr(i,j,k,comp)*Bz_arr(i,j,k,comp)
293 + Ex_arr(i+1,j,k,comp)*Bz_arr(i+1,j,k,comp)
294 + Ex_arr(i,j,k+1,comp)*Bz_arr(i,j,k+1,comp)
295 + Ex_arr(i+1,j,k+1,comp)*Bz_arr(i+1,j,k+1,comp));
296#else
297 amrex::ignore_unused(i, j, k, comp, Ex_arr, Bz_arr);
298 return 0._rt;
299#endif
300 }
301
303 static amrex::Real EyBz(int i, int j, int k, int comp,
305 amrex::Array4<const amrex::Real> const & Bz_arr) {
306 using namespace amrex::literals;
307#if defined(WARPX_DIM_3D)
308 return 0.25_rt*(Ey_arr(i,j,k,comp)*Bz_arr(i,j,k,comp)
309 + Ey_arr(i,j+1,k,comp)*Bz_arr(i,j+1,k,comp)
310 + Ey_arr(i,j,k+1,comp)*Bz_arr(i,j,k+1,comp)
311 + Ey_arr(i,j+1,k+1,comp)*Bz_arr(i,j+1,k+1,comp));
312#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
313 return 0.5_rt*(Ey_arr(i,j,k,comp)*Bz_arr(i,j,k,comp)
314 + Ey_arr(i,j+1,k,comp)*Bz_arr(i,j+1,k,comp));
315#elif defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
316 return Ey_arr(i,j,k,comp)*Bz_arr(i,j,k,comp);
317#elif defined(WARPX_DIM_1D_Z)
318 amrex::ignore_unused(i, j, k, comp, Ey_arr, Bz_arr);
319 return 0._rt;
320#endif
321 }
322
324 static amrex::Real EzBy(int i, int j, int k, int comp,
326 amrex::Array4<const amrex::Real> const & By_arr) {
327 using namespace amrex::literals;
328#if defined(WARPX_DIM_3D)
329 return 0.25_rt*(Ez_arr(i,j,k,comp)*By_arr(i,j,k,comp)
330 + Ez_arr(i,j+1,k,comp)*By_arr(i,j+1,k,comp)
331 + Ez_arr(i,j,k+1,comp)*By_arr(i,j,k+1,comp)
332 + Ez_arr(i,j+1,k+1,comp)*By_arr(i,j+1,k+1,comp));
333#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
334 return 0.5_rt*(Ez_arr(i,j,k,comp)*By_arr(i,j,k,comp)
335 + Ez_arr(i,j+1,k,comp)*By_arr(i,j+1,k,comp));
336#elif defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
337 return Ez_arr(i,j,k,comp)*By_arr(i,j,k,comp);
338#elif defined(WARPX_DIM_1D_Z)
339 amrex::ignore_unused(i, j, k, comp, Ez_arr, By_arr);
340 return 0._rt;
341#endif
342 }
343
344};
345
346namespace Poynting {
347
348 template<int normal_dir, typename T_Algo, typename AreaFunc>
356 AreaFunc area_factor)
357 {
358
360 amrex::ReduceData<amrex::Real> reduce_data(reduce_ops);
361
362 reduce_ops.eval(box, reduce_data,
363 [=] AMREX_GPU_DEVICE (int i, int j, int k) -> amrex::GpuTuple<amrex::Real>
364 {
365
366 constexpr int comp = 0;
367
368 // This is needed by the GPU compilers where captured variables cannot first appear
369 // in a constexpr-if context
370 amrex::ignore_unused(Ez_arr, By_arr, Ex_arr, Bz_arr, Ey_arr, Bx_arr);
371
372 amrex::Real const af = area_factor(i,j,k);
373 if constexpr (normal_dir == 0) {
374 return af*(T_Algo::EyBz(i,j,k,comp,Ey_arr,Bz_arr) - T_Algo::EzBy(i,j,k,comp,Ez_arr,By_arr));
375 }
376
377 else if constexpr (normal_dir == 1) {
378 return af*(T_Algo::EzBx(i,j,k,comp,Ez_arr,Bx_arr) - T_Algo::ExBz(i,j,k,comp,Ex_arr,Bz_arr));
379 }
380
381 else if constexpr (normal_dir == 2) {
382 return af*(T_Algo::ExBy(i,j,k,comp,Ex_arr,By_arr) - T_Algo::EyBx(i,j,k,comp,Ey_arr,Bx_arr));
383 }
384
385 });
386
387 auto r = reduce_data.value();
388 return amrex::get<0>(r);
389 }
390}
391
397{
398public:
399
405 explicit FieldPoyntingFlux (const std::string& rd_name);
406
412 void ComputeDiags (int step) final;
413
419 void ComputeDiagsMidStep (int step) final;
420
426 void ComputePoyntingFlux ();
427
428 void WriteCheckpointData (std::string const & dir) final;
429
430 void ReadCheckpointData (std::string const & dir) final;
431
432private:
433
434 bool use_mid_step_value = false;
435
436};
437
438#endif
#define AMREX_INLINE
#define AMREX_GPU_DEVICE
#define AMREX_GPU_HOST_DEVICE
void WriteCheckpointData(std::string const &dir) final
Write out checkpoint related data.
Definition FieldPoyntingFlux.cpp:312
FieldPoyntingFlux(const std::string &rd_name)
Constructor.
Definition FieldPoyntingFlux.cpp:43
void ComputePoyntingFlux()
This function computes the electromagnetic Poynting flux, obtained by integrating the electromagnetic...
Definition FieldPoyntingFlux.cpp:117
void ReadCheckpointData(std::string const &dir) final
Read in checkpoint related data.
Definition FieldPoyntingFlux.cpp:328
bool use_mid_step_value
Definition FieldPoyntingFlux.H:434
void ComputeDiagsMidStep(int step) final
Call the routine to compute the Poynting flux at the mid step time level.
Definition FieldPoyntingFlux.cpp:110
void ComputeDiags(int step) final
Call the routine to compute the Poynting flux if needed.
Definition FieldPoyntingFlux.cpp:101
ReducedDiags(const std::string &rd_name)
Definition ReducedDiags.cpp:26
std::enable_if_t< IsFabArray< MF >::value > eval(MF const &mf, IntVect const &nghost, D &reduce_data, F &&f)
amrex_real Real
ArrayND< T, 4, true > Array4
Definition FieldPoyntingFlux.H:346
amrex::Real Kernel(amrex::Box const &box, amrex::Array4< const amrex::Real > const &Ex_arr, amrex::Array4< const amrex::Real > const &Ey_arr, amrex::Array4< const amrex::Real > const &Ez_arr, amrex::Array4< const amrex::Real > const &Bx_arr, amrex::Array4< const amrex::Real > const &By_arr, amrex::Array4< const amrex::Real > const &Bz_arr, AreaFunc area_factor)
Definition FieldPoyntingFlux.H:349
__host__ __device__ void ignore_unused(const Ts &...)
BoxND< 3 > Box
__host__ __device__ constexpr int get(IntVectND< dim > const &iv) noexcept
Definition FieldPoyntingFlux.H:138
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real EzBy(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ez_arr, amrex::Array4< const amrex::Real > const &By_arr)
Definition FieldPoyntingFlux.H:211
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real EyBx(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ey_arr, amrex::Array4< const amrex::Real > const &Bx_arr)
Definition FieldPoyntingFlux.H:141
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real EzBx(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ez_arr, amrex::Array4< const amrex::Real > const &Bx_arr)
Definition FieldPoyntingFlux.H:169
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real ExBy(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ex_arr, amrex::Array4< const amrex::Real > const &By_arr)
Definition FieldPoyntingFlux.H:155
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real EyBz(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ey_arr, amrex::Array4< const amrex::Real > const &Bz_arr)
Definition FieldPoyntingFlux.H:197
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real ExBz(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ex_arr, amrex::Array4< const amrex::Real > const &Bz_arr)
Definition FieldPoyntingFlux.H:183
Definition FieldPoyntingFlux.H:226
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real ExBz(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ex_arr, amrex::Array4< const amrex::Real > const &Bz_arr)
Definition FieldPoyntingFlux.H:287
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real EyBz(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ey_arr, amrex::Array4< const amrex::Real > const &Bz_arr)
Definition FieldPoyntingFlux.H:303
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real EzBy(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ez_arr, amrex::Array4< const amrex::Real > const &By_arr)
Definition FieldPoyntingFlux.H:324
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real EzBx(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ez_arr, amrex::Array4< const amrex::Real > const &Bx_arr)
Definition FieldPoyntingFlux.H:271
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real ExBy(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ex_arr, amrex::Array4< const amrex::Real > const &By_arr)
Definition FieldPoyntingFlux.H:250
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real EyBx(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ey_arr, amrex::Array4< const amrex::Real > const &Bx_arr)
Definition FieldPoyntingFlux.H:229
Definition FieldPoyntingFlux.H:21
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real ExBy(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ex_arr, amrex::Array4< const amrex::Real > const &By_arr)
Definition FieldPoyntingFlux.H:52
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real ExBz(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ex_arr, amrex::Array4< const amrex::Real > const &Bz_arr)
Definition FieldPoyntingFlux.H:84
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real EzBx(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ez_arr, amrex::Array4< const amrex::Real > const &Bx_arr)
Definition FieldPoyntingFlux.H:69
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real EyBz(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ey_arr, amrex::Array4< const amrex::Real > const &Bz_arr)
Definition FieldPoyntingFlux.H:99
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real EzBy(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ez_arr, amrex::Array4< const amrex::Real > const &By_arr)
Definition FieldPoyntingFlux.H:120
AMREX_GPU_HOST_DEVICE static AMREX_INLINE amrex::Real EyBx(int i, int j, int k, int comp, amrex::Array4< const amrex::Real > const &Ey_arr, amrex::Array4< const amrex::Real > const &Bx_arr)
Definition FieldPoyntingFlux.H:31