36 template <
int depos_order,
int galerkin_
interpolation>
37 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
39 const amrex::ParticleReal yp,
40 const amrex::ParticleReal zp,
41 amrex::ParticleReal& Exp,
42 amrex::ParticleReal& Eyp,
43 amrex::ParticleReal& Ezp,
44 amrex::ParticleReal& Bxp,
45 amrex::ParticleReal& Byp,
46 amrex::ParticleReal& Bzp,
47 amrex::Array4<amrex::Real const>
const& ex_arr,
48 amrex::Array4<amrex::Real const>
const& ey_arr,
49 amrex::Array4<amrex::Real const>
const& ez_arr,
50 amrex::Array4<amrex::Real const>
const& bx_arr,
51 amrex::Array4<amrex::Real const>
const& by_arr,
52 amrex::Array4<amrex::Real const>
const& bz_arr,
53 const amrex::IndexType ex_type,
54 const amrex::IndexType ey_type,
55 const amrex::IndexType ez_type,
56 const amrex::IndexType bx_type,
57 const amrex::IndexType by_type,
58 const amrex::IndexType bz_type,
59 const amrex::GpuArray<amrex::Real, 3>&
dx,
60 const amrex::GpuArray<amrex::Real, 3>& xyzmin,
61 const amrex::Dim3& lo,
62 const long n_rz_azimuthal_modes)
64 using namespace amrex;
66 #if defined(WARPX_DIM_XZ) 67 amrex::ignore_unused(yp);
71 amrex::ignore_unused(n_rz_azimuthal_modes);
74 const amrex::Real dxi = 1.0/dx[0];
75 const amrex::Real dzi = 1.0/dx[2];
76 #if (AMREX_SPACEDIM == 3) 77 const amrex::Real dyi = 1.0/dx[1];
80 const amrex::Real xmin = xyzmin[0];
81 #if (AMREX_SPACEDIM == 3) 82 const amrex::Real ymin = xyzmin[1];
84 const amrex::Real zmin = xyzmin[2];
86 constexpr
int zdir = (AMREX_SPACEDIM - 1);
88 constexpr
int CELL = amrex::IndexType::CELL;
94 const amrex::Real rp = std::sqrt(xp*xp + yp*yp);
95 const amrex::Real
x = (rp - xmin)*dxi;
97 const amrex::Real x = (xp-xmin)*dxi;
104 amrex::Real sx_node[depos_order + 1];
105 amrex::Real sx_cell[depos_order + 1];
106 amrex::Real sx_node_galerkin[depos_order + 1 - galerkin_interpolation];
107 amrex::Real sx_cell_galerkin[depos_order + 1 - galerkin_interpolation];
113 Compute_shape_factor<depos_order - galerkin_interpolation >
const compute_shape_factor_galerkin;
114 if ((ey_type[0] == NODE) || (ez_type[0] ==
NODE) || (bx_type[0] == NODE)) {
115 j_node = compute_shape_factor(sx_node, x);
117 if ((ey_type[0] == CELL) || (ez_type[0] == CELL) || (bx_type[0] == CELL)) {
118 j_cell = compute_shape_factor(sx_cell, x - 0.5_rt);
120 if ((ex_type[0] == NODE) || (by_type[0] == NODE) || (bz_type[0] == NODE)) {
121 j_node_v = compute_shape_factor_galerkin(sx_node_galerkin, x);
123 if ((ex_type[0] == CELL) || (by_type[0] == CELL) || (bz_type[0] == CELL)) {
124 j_cell_v = compute_shape_factor_galerkin(sx_cell_galerkin, x - 0.5_rt);
126 const amrex::Real (&sx_ex)[depos_order + 1 - galerkin_interpolation] = ((ex_type[0] ==
NODE) ? sx_node_galerkin : sx_cell_galerkin);
127 const amrex::Real (&sx_ey)[depos_order + 1 ] = ((ey_type[0] ==
NODE) ? sx_node : sx_cell );
128 const amrex::Real (&sx_ez)[depos_order + 1 ] = ((ez_type[0] ==
NODE) ? sx_node : sx_cell );
129 const amrex::Real (&sx_bx)[depos_order + 1 ] = ((bx_type[0] ==
NODE) ? sx_node : sx_cell );
130 const amrex::Real (&sx_by)[depos_order + 1 - galerkin_interpolation] = ((by_type[0] ==
NODE) ? sx_node_galerkin : sx_cell_galerkin);
131 const amrex::Real (&sx_bz)[depos_order + 1 - galerkin_interpolation] = ((bz_type[0] ==
NODE) ? sx_node_galerkin : sx_cell_galerkin);
132 int const j_ex = ((ex_type[0] ==
NODE) ? j_node_v : j_cell_v);
133 int const j_ey = ((ey_type[0] ==
NODE) ? j_node : j_cell );
134 int const j_ez = ((ez_type[0] ==
NODE) ? j_node : j_cell );
135 int const j_bx = ((bx_type[0] ==
NODE) ? j_node : j_cell );
136 int const j_by = ((by_type[0] ==
NODE) ? j_node_v : j_cell_v);
137 int const j_bz = ((bz_type[0] ==
NODE) ? j_node_v : j_cell_v);
139 #if (AMREX_SPACEDIM == 3) 141 const amrex::Real y = (yp-ymin)*dyi;
142 amrex::Real sy_node[depos_order + 1];
143 amrex::Real sy_cell[depos_order + 1];
144 amrex::Real sy_node_v[depos_order + 1 - galerkin_interpolation];
145 amrex::Real sy_cell_v[depos_order + 1 - galerkin_interpolation];
150 if ((ex_type[1] == NODE) || (ez_type[1] ==
NODE) || (by_type[1] == NODE)) {
151 k_node = compute_shape_factor(sy_node, y);
153 if ((ex_type[1] == CELL) || (ez_type[1] == CELL) || (by_type[1] == CELL)) {
154 k_cell = compute_shape_factor(sy_cell, y - 0.5_rt);
156 if ((ey_type[1] == NODE) || (bx_type[1] == NODE) || (bz_type[1] == NODE)) {
157 k_node_v = compute_shape_factor_galerkin(sy_node_v, y);
159 if ((ey_type[1] == CELL) || (bx_type[1] == CELL) || (bz_type[1] == CELL)) {
160 k_cell_v = compute_shape_factor_galerkin(sy_cell_v, y - 0.5_rt);
162 const amrex::Real (&sy_ex)[depos_order + 1 ] = ((ex_type[1] ==
NODE) ? sy_node : sy_cell );
163 const amrex::Real (&sy_ey)[depos_order + 1 - galerkin_interpolation] = ((ey_type[1] ==
NODE) ? sy_node_v : sy_cell_v);
164 const amrex::Real (&sy_ez)[depos_order + 1 ] = ((ez_type[1] ==
NODE) ? sy_node : sy_cell );
165 const amrex::Real (&sy_bx)[depos_order + 1 - galerkin_interpolation] = ((bx_type[1] ==
NODE) ? sy_node_v : sy_cell_v);
166 const amrex::Real (&sy_by)[depos_order + 1 ] = ((by_type[1] ==
NODE) ? sy_node : sy_cell );
167 const amrex::Real (&sy_bz)[depos_order + 1 - galerkin_interpolation] = ((bz_type[1] ==
NODE) ? sy_node_v : sy_cell_v);
168 int const k_ex = ((ex_type[1] ==
NODE) ? k_node : k_cell );
169 int const k_ey = ((ey_type[1] ==
NODE) ? k_node_v : k_cell_v);
170 int const k_ez = ((ez_type[1] ==
NODE) ? k_node : k_cell );
171 int const k_bx = ((bx_type[1] ==
NODE) ? k_node_v : k_cell_v);
172 int const k_by = ((by_type[1] ==
NODE) ? k_node : k_cell );
173 int const k_bz = ((bz_type[1] ==
NODE) ? k_node_v : k_cell_v);
177 const amrex::Real
z = (zp-zmin)*dzi;
178 amrex::Real sz_node[depos_order + 1];
179 amrex::Real sz_cell[depos_order + 1];
180 amrex::Real sz_node_v[depos_order + 1 - galerkin_interpolation];
181 amrex::Real sz_cell_v[depos_order + 1 - galerkin_interpolation];
186 if ((ex_type[zdir] == NODE) || (ey_type[zdir] ==
NODE) || (bz_type[zdir] == NODE)) {
187 l_node = compute_shape_factor(sz_node, z);
189 if ((ex_type[zdir] == CELL) || (ey_type[zdir] == CELL) || (bz_type[zdir] == CELL)) {
190 l_cell = compute_shape_factor(sz_cell, z - 0.5_rt);
192 if ((ez_type[zdir] == NODE) || (bx_type[zdir] == NODE) || (by_type[zdir] == NODE)) {
193 l_node_v = compute_shape_factor_galerkin(sz_node_v, z);
195 if ((ez_type[zdir] == CELL) || (bx_type[zdir] == CELL) || (by_type[zdir] == CELL)) {
196 l_cell_v = compute_shape_factor_galerkin(sz_cell_v, z - 0.5_rt);
198 const amrex::Real (&sz_ex)[depos_order + 1 ] = ((ex_type[zdir] ==
NODE) ? sz_node : sz_cell );
199 const amrex::Real (&sz_ey)[depos_order + 1 ] = ((ey_type[zdir] ==
NODE) ? sz_node : sz_cell );
200 const amrex::Real (&sz_ez)[depos_order + 1 - galerkin_interpolation] = ((ez_type[zdir] ==
NODE) ? sz_node_v : sz_cell_v);
201 const amrex::Real (&sz_bx)[depos_order + 1 - galerkin_interpolation] = ((bx_type[zdir] ==
NODE) ? sz_node_v : sz_cell_v);
202 const amrex::Real (&sz_by)[depos_order + 1 - galerkin_interpolation] = ((by_type[zdir] ==
NODE) ? sz_node_v : sz_cell_v);
203 const amrex::Real (&sz_bz)[depos_order + 1 ] = ((bz_type[zdir] ==
NODE) ? sz_node : sz_cell );
204 int const l_ex = ((ex_type[zdir] ==
NODE) ? l_node : l_cell );
205 int const l_ey = ((ey_type[zdir] ==
NODE) ? l_node : l_cell );
206 int const l_ez = ((ez_type[zdir] ==
NODE) ? l_node_v : l_cell_v);
207 int const l_bx = ((bx_type[zdir] ==
NODE) ? l_node_v : l_cell_v);
208 int const l_by = ((by_type[zdir] ==
NODE) ? l_node_v : l_cell_v);
209 int const l_bz = ((bz_type[zdir] ==
NODE) ? l_node : l_cell );
216 #if (AMREX_SPACEDIM == 2) 218 for (
int iz=0; iz<=depos_order; iz++){
219 for (
int ix=0; ix<=depos_order; ix++){
220 Eyp += sx_ey[ix]*sz_ey[iz]*
221 ey_arr(lo.x+j_ey+ix, lo.y+l_ey+iz, 0, 0);
226 for (
int iz=0; iz<=depos_order; iz++){
227 for (
int ix=0; ix<=depos_order-galerkin_interpolation; ix++){
228 Exp += sx_ex[ix]*sz_ex[iz]*
229 ex_arr(lo.x+j_ex+ix, lo.y+l_ex+iz, 0, 0);
230 Bzp += sx_bz[ix]*sz_bz[iz]*
231 bz_arr(lo.x+j_bz+ix, lo.y+l_bz+iz, 0, 0);
236 for (
int iz=0; iz<=depos_order-galerkin_interpolation; iz++){
237 for (
int ix=0; ix<=depos_order; ix++){
238 Ezp += sx_ez[ix]*sz_ez[iz]*
239 ez_arr(lo.x+j_ez+ix, lo.y+l_ez+iz, 0, 0);
240 Bxp += sx_bx[ix]*sz_bx[iz]*
241 bx_arr(lo.x+j_bx+ix, lo.y+l_bx+iz, 0, 0);
245 for (
int iz=0; iz<=depos_order-galerkin_interpolation; iz++){
246 for (
int ix=0; ix<=depos_order-galerkin_interpolation; ix++){
247 Byp += sx_by[ix]*sz_by[iz]*
248 by_arr(lo.x+j_by+ix, lo.y+l_by+iz, 0, 0);
254 amrex::Real costheta;
255 amrex::Real sintheta;
266 for (
int imode=1 ; imode < n_rz_azimuthal_modes ; imode++) {
269 for (
int iz=0; iz<=depos_order; iz++){
270 for (
int ix=0; ix<=depos_order; ix++){
271 const amrex::Real dEy = (+ ey_arr(lo.x+j_ey+ix, lo.y+l_ey+iz, 0, 2*imode-1)*xy.real()
272 - ey_arr(lo.x+j_ey+ix, lo.y+l_ey+iz, 0, 2*imode)*xy.imag());
273 Eyp += sx_ey[ix]*sz_ey[iz]*dEy;
278 for (
int iz=0; iz<=depos_order; iz++){
279 for (
int ix=0; ix<=depos_order-galerkin_interpolation; ix++){
280 const amrex::Real dEx = (+ ex_arr(lo.x+j_ex+ix, lo.y+l_ex+iz, 0, 2*imode-1)*xy.real()
281 - ex_arr(lo.x+j_ex+ix, lo.y+l_ex+iz, 0, 2*imode)*xy.imag());
282 Exp += sx_ex[ix]*sz_ex[iz]*dEx;
283 const amrex::Real dBz = (+ bz_arr(lo.x+j_bz+ix, lo.y+l_bz+iz, 0, 2*imode-1)*xy.real()
284 - bz_arr(lo.x+j_bz+ix, lo.y+l_bz+iz, 0, 2*imode)*xy.imag());
285 Bzp += sx_bz[ix]*sz_bz[iz]*dBz;
290 for (
int iz=0; iz<=depos_order-galerkin_interpolation; iz++){
291 for (
int ix=0; ix<=depos_order; ix++){
292 const amrex::Real dEz = (+ ez_arr(lo.x+j_ez+ix, lo.y+l_ez+iz, 0, 2*imode-1)*xy.real()
293 - ez_arr(lo.x+j_ez+ix, lo.y+l_ez+iz, 0, 2*imode)*xy.imag());
294 Ezp += sx_ez[ix]*sz_ez[iz]*dEz;
295 const amrex::Real dBx = (+ bx_arr(lo.x+j_bx+ix, lo.y+l_bx+iz, 0, 2*imode-1)*xy.real()
296 - bx_arr(lo.x+j_bx+ix, lo.y+l_bx+iz, 0, 2*imode)*xy.imag());
297 Bxp += sx_bx[ix]*sz_bx[iz]*dBx;
301 for (
int iz=0; iz<=depos_order-galerkin_interpolation; iz++){
302 for (
int ix=0; ix<=depos_order-galerkin_interpolation; ix++){
303 const amrex::Real dBy = (+ by_arr(lo.x+j_by+ix, lo.y+l_by+iz, 0, 2*imode-1)*xy.real()
304 - by_arr(lo.x+j_by+ix, lo.y+l_by+iz, 0, 2*imode)*xy.imag());
305 Byp += sx_by[ix]*sz_by[iz]*dBy;
312 const amrex::Real Exp_save = Exp;
313 Exp = costheta*Exp - sintheta*Eyp;
314 Eyp = costheta*Eyp + sintheta*Exp_save;
315 const amrex::Real Bxp_save = Bxp;
316 Bxp = costheta*Bxp - sintheta*Byp;
317 Byp = costheta*Byp + sintheta*Bxp_save;
320 #else // (AMREX_SPACEDIM == 3) 322 for (
int iz=0; iz<=depos_order; iz++){
323 for (
int iy=0; iy<=depos_order; iy++){
324 for (
int ix=0; ix<= depos_order - galerkin_interpolation; ix++){
325 Exp += sx_ex[ix]*sy_ex[iy]*sz_ex[iz]*
326 ex_arr(lo.x+j_ex+ix, lo.y+k_ex+iy, lo.z+l_ex+iz);
331 for (
int iz=0; iz<=depos_order; iz++){
332 for (
int iy=0; iy<= depos_order - galerkin_interpolation; iy++){
333 for (
int ix=0; ix<=depos_order; ix++){
334 Eyp += sx_ey[ix]*sy_ey[iy]*sz_ey[iz]*
335 ey_arr(lo.x+j_ey+ix, lo.y+k_ey+iy, lo.z+l_ey+iz);
340 for (
int iz=0; iz<= depos_order - galerkin_interpolation; iz++){
341 for (
int iy=0; iy<=depos_order; iy++){
342 for (
int ix=0; ix<=depos_order; ix++){
343 Ezp += sx_ez[ix]*sy_ez[iy]*sz_ez[iz]*
344 ez_arr(lo.x+j_ez+ix, lo.y+k_ez+iy, lo.z+l_ez+iz);
349 for (
int iz=0; iz<=depos_order; iz++){
350 for (
int iy=0; iy<= depos_order - galerkin_interpolation; iy++){
351 for (
int ix=0; ix<= depos_order - galerkin_interpolation; ix++){
352 Bzp += sx_bz[ix]*sy_bz[iy]*sz_bz[iz]*
353 bz_arr(lo.x+j_bz+ix, lo.y+k_bz+iy, lo.z+l_bz+iz);
358 for (
int iz=0; iz<= depos_order - galerkin_interpolation; iz++){
359 for (
int iy=0; iy<=depos_order; iy++){
360 for (
int ix=0; ix<= depos_order - galerkin_interpolation; ix++){
361 Byp += sx_by[ix]*sy_by[iy]*sz_by[iz]*
362 by_arr(lo.x+j_by+ix, lo.y+k_by+iy, lo.z+l_by+iz);
367 for (
int iz=0; iz<= depos_order - galerkin_interpolation; iz++){
368 for (
int iy=0; iy<= depos_order - galerkin_interpolation; iy++){
369 for (
int ix=0; ix<=depos_order; ix++){
370 Bxp += sx_bx[ix]*sy_bx[iy]*sz_bx[iz]*
371 bx_arr(lo.x+j_bx+ix, lo.y+k_bx+iy, lo.z+l_bx+iz);
394 template <
int depos_order,
int lower_in_v>
397 amrex::ParticleReal *
const Exp, amrex::ParticleReal *
const Eyp,
398 amrex::ParticleReal *
const Ezp, amrex::ParticleReal *
const Bxp,
399 amrex::ParticleReal *
const Byp, amrex::ParticleReal *
const Bzp,
400 amrex::FArrayBox
const *
const exfab,
401 amrex::FArrayBox
const *
const eyfab,
402 amrex::FArrayBox
const *
const ezfab,
403 amrex::FArrayBox
const *
const bxfab,
404 amrex::FArrayBox
const *
const byfab,
405 amrex::FArrayBox
const *
const bzfab,
406 const long np_to_gather,
407 const std::array<amrex::Real, 3>&
dx,
408 const std::array<amrex::Real, 3> xyzmin,
409 const amrex::Dim3 lo,
410 const long n_rz_azimuthal_modes)
413 amrex::GpuArray<amrex::Real, 3> dx_arr = {dx[0], dx[1], dx[2]};
414 amrex::GpuArray<amrex::Real, 3> xyzmin_arr = {xyzmin[0], xyzmin[1], xyzmin[2]};
416 amrex::Array4<const amrex::Real>
const& ex_arr = exfab->array();
417 amrex::Array4<const amrex::Real>
const& ey_arr = eyfab->array();
418 amrex::Array4<const amrex::Real>
const& ez_arr = ezfab->array();
419 amrex::Array4<const amrex::Real>
const& bx_arr = bxfab->array();
420 amrex::Array4<const amrex::Real>
const& by_arr = byfab->array();
421 amrex::Array4<const amrex::Real>
const& bz_arr = bzfab->array();
423 amrex::IndexType
const ex_type = exfab->box().ixType();
424 amrex::IndexType
const ey_type = eyfab->box().ixType();
425 amrex::IndexType
const ez_type = ezfab->box().ixType();
426 amrex::IndexType
const bx_type = bxfab->box().ixType();
427 amrex::IndexType
const by_type = byfab->box().ixType();
428 amrex::IndexType
const bz_type = bzfab->box().ixType();
434 [=] AMREX_GPU_DEVICE (
long ip) {
436 amrex::ParticleReal xp, yp, zp;
437 getPosition(ip, xp, yp, zp);
438 getExternalE(ip, Exp[ip], Eyp[ip], Ezp[ip]);
439 getExternalB(ip, Bxp[ip], Byp[ip], Bzp[ip]);
441 doGatherShapeN<depos_order, lower_in_v>(
442 xp, yp, zp, Exp[ip], Eyp[ip], Ezp[ip], Bxp[ip], Byp[ip], Bzp[ip],
443 ex_arr, ey_arr, ez_arr, bx_arr, by_arr, bz_arr,
444 ex_type, ey_type, ez_type, bx_type, by_type, bz_type,
445 dx_arr, xyzmin_arr, lo, n_rz_azimuthal_modes);
467 AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
469 const amrex::ParticleReal yp,
470 const amrex::ParticleReal zp,
471 amrex::ParticleReal& Exp,
472 amrex::ParticleReal& Eyp,
473 amrex::ParticleReal& Ezp,
474 amrex::ParticleReal& Bxp,
475 amrex::ParticleReal& Byp,
476 amrex::ParticleReal& Bzp,
477 amrex::Array4<amrex::Real const>
const& ex_arr,
478 amrex::Array4<amrex::Real const>
const& ey_arr,
479 amrex::Array4<amrex::Real const>
const& ez_arr,
480 amrex::Array4<amrex::Real const>
const& bx_arr,
481 amrex::Array4<amrex::Real const>
const& by_arr,
482 amrex::Array4<amrex::Real const>
const& bz_arr,
483 const amrex::IndexType ex_type,
484 const amrex::IndexType ey_type,
485 const amrex::IndexType ez_type,
486 const amrex::IndexType bx_type,
487 const amrex::IndexType by_type,
488 const amrex::IndexType bz_type,
489 const amrex::GpuArray<amrex::Real, 3>& dx_arr,
490 const amrex::GpuArray<amrex::Real, 3>& xyzmin_arr,
491 const amrex::Dim3& lo,
492 const long n_rz_azimuthal_modes,
494 const bool galerkin_interpolation)
496 if (galerkin_interpolation) {
498 doGatherShapeN<1,1>(xp, yp, zp, Exp, Eyp, Ezp, Bxp, Byp, Bzp,
499 ex_arr, ey_arr, ez_arr, bx_arr, by_arr, bz_arr,
500 ex_type, ey_type, ez_type, bx_type, by_type, bz_type,
501 dx_arr, xyzmin_arr, lo, n_rz_azimuthal_modes);
502 }
else if (nox == 2) {
503 doGatherShapeN<2,1>(xp, yp, zp, Exp, Eyp, Ezp, Bxp, Byp, Bzp,
504 ex_arr, ey_arr, ez_arr, bx_arr, by_arr, bz_arr,
505 ex_type, ey_type, ez_type, bx_type, by_type, bz_type,
506 dx_arr, xyzmin_arr, lo, n_rz_azimuthal_modes);
507 }
else if (nox == 3) {
508 doGatherShapeN<3,1>(xp, yp, zp, Exp, Eyp, Ezp, Bxp, Byp, Bzp,
509 ex_arr, ey_arr, ez_arr, bx_arr, by_arr, bz_arr,
510 ex_type, ey_type, ez_type, bx_type, by_type, bz_type,
511 dx_arr, xyzmin_arr, lo, n_rz_azimuthal_modes);
515 doGatherShapeN<1,0>(xp, yp, zp, Exp, Eyp, Ezp, Bxp, Byp, Bzp,
516 ex_arr, ey_arr, ez_arr, bx_arr, by_arr, bz_arr,
517 ex_type, ey_type, ez_type, bx_type, by_type, bz_type,
518 dx_arr, xyzmin_arr, lo, n_rz_azimuthal_modes);
519 }
else if (nox == 2) {
520 doGatherShapeN<2,0>(xp, yp, zp, Exp, Eyp, Ezp, Bxp, Byp, Bzp,
521 ex_arr, ey_arr, ez_arr, bx_arr, by_arr, bz_arr,
522 ex_type, ey_type, ez_type, bx_type, by_type, bz_type,
523 dx_arr, xyzmin_arr, lo, n_rz_azimuthal_modes);
524 }
else if (nox == 3) {
525 doGatherShapeN<3,0>(xp, yp, zp, Exp, Eyp, Ezp, Bxp, Byp, Bzp,
526 ex_arr, ey_arr, ez_arr, bx_arr, by_arr, bz_arr,
527 ex_type, ey_type, ez_type, bx_type, by_type, bz_type,
528 dx_arr, xyzmin_arr, lo, n_rz_azimuthal_modes);
533 #endif // FIELDGATHER_H_ amrex::GpuComplex< amrex::Real > Complex
Definition: WarpX_Complex.H:22
Definition: wp_parser.tab.cpp:115
def x
Definition: read_lab_particles.py:25
int dx
Definition: compute_domain.py:35
def z
Definition: read_lab_particles.py:26
Functor that can be used to assign the external B field to a particle inside a ParallelFor kernel...
Definition: GetExternalFields.H:68
Definition: ShapeFactors.H:22
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void doGatherShapeN(const amrex::ParticleReal xp, const amrex::ParticleReal yp, const amrex::ParticleReal zp, amrex::ParticleReal &Exp, amrex::ParticleReal &Eyp, amrex::ParticleReal &Ezp, amrex::ParticleReal &Bxp, amrex::ParticleReal &Byp, amrex::ParticleReal &Bzp, amrex::Array4< amrex::Real const > const &ex_arr, amrex::Array4< amrex::Real const > const &ey_arr, amrex::Array4< amrex::Real const > const &ez_arr, amrex::Array4< amrex::Real const > const &bx_arr, amrex::Array4< amrex::Real const > const &by_arr, amrex::Array4< amrex::Real const > const &bz_arr, const amrex::IndexType ex_type, const amrex::IndexType ey_type, const amrex::IndexType ez_type, const amrex::IndexType bx_type, const amrex::IndexType by_type, const amrex::IndexType bz_type, const amrex::GpuArray< amrex::Real, 3 > &dx, const amrex::GpuArray< amrex::Real, 3 > &xyzmin, const amrex::Dim3 &lo, const long n_rz_azimuthal_modes)
Field gather for a single particle.
Definition: FieldGather.H:38
Functor that can be used to extract the positions of the macroparticles inside a ParallelFor kernel...
Definition: GetAndSetPosition.H:25
Functor that can be used to assign the external E field to a particle inside a ParallelFor kernel...
Definition: GetExternalFields.H:58