1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
|
/* Copyright 2019 Axel Huebl, David Grote, Maxence Thevenet
* Revathi Jambunathan, Weiqun Zhang
*
* This file is part of WarpX.
*
* License: BSD-3-Clause-LBNL
*/
#ifndef FIELDGATHER_H_
#define FIELDGATHER_H_
#include "Particles/Gather/GetExternalFields.H"
#include "Particles/Pusher/GetAndSetPosition.H"
#include "Particles/ShapeFactors.H"
#include "Utils/WarpX_Complex.H"
#include <AMReX.H>
/**
* \brief Field gather for a single particle
*
* \tparam depos_order Particle shape order
* \tparam galerkin_interpolation Lower the order of the particle shape by
* this value (0/1) for the parallel field component
* \param xp,yp,zp Particle position coordinates
* \param Exp,Eyp,Ezp Electric field on particles.
* \param Bxp,Byp,Bzp Magnetic field on particles.
* \param ex_arr,ey_arr,ez_arr Array4 of the electric field, either full array or tile.
* \param bx_arr,by_arr,bz_arr Array4 of the magnetic field, either full array or tile.
* \param ex_type,ey_type,ez_type IndexType of the electric field
* \param bx_type,by_type,bz_type IndexType of the magnetic field
* \param dx 3D cell spacing
* \param xyzmin Physical lower bounds of domain in x, y, z.
* \param lo Index lower bounds of domain.
* \param n_rz_azimuthal_modes Number of azimuthal modes when using RZ geometry
*/
template <int depos_order, int galerkin_interpolation>
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 int n_rz_azimuthal_modes)
{
using namespace amrex;
#if defined(WARPX_DIM_XZ)
amrex::ignore_unused(yp);
#endif
#if defined(WARPX_DIM_1D_Z)
amrex::ignore_unused(xp,yp);
#endif
#ifndef WARPX_DIM_RZ
amrex::ignore_unused(n_rz_azimuthal_modes);
#endif
#if (AMREX_SPACEDIM >= 2)
const amrex::Real dxi = 1.0_rt/dx[0];
#endif
const amrex::Real dzi = 1.0_rt/dx[2];
#if defined(WARPX_DIM_3D)
const amrex::Real dyi = 1.0_rt/dx[1];
#endif
#if (AMREX_SPACEDIM >= 2)
const amrex::Real xmin = xyzmin[0];
#endif
#if defined(WARPX_DIM_3D)
const amrex::Real ymin = xyzmin[1];
#endif
const amrex::Real zmin = xyzmin[2];
constexpr int zdir = WARPX_ZINDEX;
constexpr int NODE = amrex::IndexType::NODE;
constexpr int CELL = amrex::IndexType::CELL;
// --- Compute shape factors
Compute_shape_factor< depos_order > const compute_shape_factor;
Compute_shape_factor<depos_order - galerkin_interpolation > const compute_shape_factor_galerkin;
#if (AMREX_SPACEDIM >= 2)
// x direction
// Get particle position
#ifdef WARPX_DIM_RZ
const amrex::Real rp = std::sqrt(xp*xp + yp*yp);
const amrex::Real x = (rp - xmin)*dxi;
#else
const amrex::Real x = (xp-xmin)*dxi;
#endif
// j_[eb][xyz] leftmost grid point in x that the particle touches for the centering of each current
// sx_[eb][xyz] shape factor along x for the centering of each current
// There are only two possible centerings, node or cell centered, so at most only two shape factor
// arrays will be needed.
amrex::Real sx_node[depos_order + 1];
amrex::Real sx_cell[depos_order + 1];
amrex::Real sx_node_galerkin[depos_order + 1 - galerkin_interpolation] = {0._rt};
amrex::Real sx_cell_galerkin[depos_order + 1 - galerkin_interpolation] = {0._rt};
int j_node = 0;
int j_cell = 0;
int j_node_v = 0;
int j_cell_v = 0;
if ((ey_type[0] == NODE) || (ez_type[0] == NODE) || (bx_type[0] == NODE)) {
j_node = compute_shape_factor(sx_node, x);
}
if ((ey_type[0] == CELL) || (ez_type[0] == CELL) || (bx_type[0] == CELL)) {
j_cell = compute_shape_factor(sx_cell, x - 0.5_rt);
}
if ((ex_type[0] == NODE) || (by_type[0] == NODE) || (bz_type[0] == NODE)) {
j_node_v = compute_shape_factor_galerkin(sx_node_galerkin, x);
}
if ((ex_type[0] == CELL) || (by_type[0] == CELL) || (bz_type[0] == CELL)) {
j_cell_v = compute_shape_factor_galerkin(sx_cell_galerkin, x - 0.5_rt);
}
const amrex::Real (&sx_ex)[depos_order + 1 - galerkin_interpolation] = ((ex_type[0] == NODE) ? sx_node_galerkin : sx_cell_galerkin);
const amrex::Real (&sx_ey)[depos_order + 1 ] = ((ey_type[0] == NODE) ? sx_node : sx_cell );
const amrex::Real (&sx_ez)[depos_order + 1 ] = ((ez_type[0] == NODE) ? sx_node : sx_cell );
const amrex::Real (&sx_bx)[depos_order + 1 ] = ((bx_type[0] == NODE) ? sx_node : sx_cell );
const amrex::Real (&sx_by)[depos_order + 1 - galerkin_interpolation] = ((by_type[0] == NODE) ? sx_node_galerkin : sx_cell_galerkin);
const amrex::Real (&sx_bz)[depos_order + 1 - galerkin_interpolation] = ((bz_type[0] == NODE) ? sx_node_galerkin : sx_cell_galerkin);
int const j_ex = ((ex_type[0] == NODE) ? j_node_v : j_cell_v);
int const j_ey = ((ey_type[0] == NODE) ? j_node : j_cell );
int const j_ez = ((ez_type[0] == NODE) ? j_node : j_cell );
int const j_bx = ((bx_type[0] == NODE) ? j_node : j_cell );
int const j_by = ((by_type[0] == NODE) ? j_node_v : j_cell_v);
int const j_bz = ((bz_type[0] == NODE) ? j_node_v : j_cell_v);
#endif
#if defined(WARPX_DIM_3D)
// y direction
const amrex::Real y = (yp-ymin)*dyi;
amrex::Real sy_node[depos_order + 1];
amrex::Real sy_cell[depos_order + 1];
amrex::Real sy_node_v[depos_order + 1 - galerkin_interpolation];
amrex::Real sy_cell_v[depos_order + 1 - galerkin_interpolation];
int k_node = 0;
int k_cell = 0;
int k_node_v = 0;
int k_cell_v = 0;
if ((ex_type[1] == NODE) || (ez_type[1] == NODE) || (by_type[1] == NODE)) {
k_node = compute_shape_factor(sy_node, y);
}
if ((ex_type[1] == CELL) || (ez_type[1] == CELL) || (by_type[1] == CELL)) {
k_cell = compute_shape_factor(sy_cell, y - 0.5_rt);
}
if ((ey_type[1] == NODE) || (bx_type[1] == NODE) || (bz_type[1] == NODE)) {
k_node_v = compute_shape_factor_galerkin(sy_node_v, y);
}
if ((ey_type[1] == CELL) || (bx_type[1] == CELL) || (bz_type[1] == CELL)) {
k_cell_v = compute_shape_factor_galerkin(sy_cell_v, y - 0.5_rt);
}
const amrex::Real (&sy_ex)[depos_order + 1 ] = ((ex_type[1] == NODE) ? sy_node : sy_cell );
const amrex::Real (&sy_ey)[depos_order + 1 - galerkin_interpolation] = ((ey_type[1] == NODE) ? sy_node_v : sy_cell_v);
const amrex::Real (&sy_ez)[depos_order + 1 ] = ((ez_type[1] == NODE) ? sy_node : sy_cell );
const amrex::Real (&sy_bx)[depos_order + 1 - galerkin_interpolation] = ((bx_type[1] == NODE) ? sy_node_v : sy_cell_v);
const amrex::Real (&sy_by)[depos_order + 1 ] = ((by_type[1] == NODE) ? sy_node : sy_cell );
const amrex::Real (&sy_bz)[depos_order + 1 - galerkin_interpolation] = ((bz_type[1] == NODE) ? sy_node_v : sy_cell_v);
int const k_ex = ((ex_type[1] == NODE) ? k_node : k_cell );
int const k_ey = ((ey_type[1] == NODE) ? k_node_v : k_cell_v);
int const k_ez = ((ez_type[1] == NODE) ? k_node : k_cell );
int const k_bx = ((bx_type[1] == NODE) ? k_node_v : k_cell_v);
int const k_by = ((by_type[1] == NODE) ? k_node : k_cell );
int const k_bz = ((bz_type[1] == NODE) ? k_node_v : k_cell_v);
#endif
// z direction
const amrex::Real z = (zp-zmin)*dzi;
amrex::Real sz_node[depos_order + 1];
amrex::Real sz_cell[depos_order + 1];
amrex::Real sz_node_v[depos_order + 1 - galerkin_interpolation];
amrex::Real sz_cell_v[depos_order + 1 - galerkin_interpolation];
int l_node = 0;
int l_cell = 0;
int l_node_v = 0;
int l_cell_v = 0;
if ((ex_type[zdir] == NODE) || (ey_type[zdir] == NODE) || (bz_type[zdir] == NODE)) {
l_node = compute_shape_factor(sz_node, z);
}
if ((ex_type[zdir] == CELL) || (ey_type[zdir] == CELL) || (bz_type[zdir] == CELL)) {
l_cell = compute_shape_factor(sz_cell, z - 0.5_rt);
}
if ((ez_type[zdir] == NODE) || (bx_type[zdir] == NODE) || (by_type[zdir] == NODE)) {
l_node_v = compute_shape_factor_galerkin(sz_node_v, z);
}
if ((ez_type[zdir] == CELL) || (bx_type[zdir] == CELL) || (by_type[zdir] == CELL)) {
l_cell_v = compute_shape_factor_galerkin(sz_cell_v, z - 0.5_rt);
}
const amrex::Real (&sz_ex)[depos_order + 1 ] = ((ex_type[zdir] == NODE) ? sz_node : sz_cell );
const amrex::Real (&sz_ey)[depos_order + 1 ] = ((ey_type[zdir] == NODE) ? sz_node : sz_cell );
const amrex::Real (&sz_ez)[depos_order + 1 - galerkin_interpolation] = ((ez_type[zdir] == NODE) ? sz_node_v : sz_cell_v);
const amrex::Real (&sz_bx)[depos_order + 1 - galerkin_interpolation] = ((bx_type[zdir] == NODE) ? sz_node_v : sz_cell_v);
const amrex::Real (&sz_by)[depos_order + 1 - galerkin_interpolation] = ((by_type[zdir] == NODE) ? sz_node_v : sz_cell_v);
const amrex::Real (&sz_bz)[depos_order + 1 ] = ((bz_type[zdir] == NODE) ? sz_node : sz_cell );
int const l_ex = ((ex_type[zdir] == NODE) ? l_node : l_cell );
int const l_ey = ((ey_type[zdir] == NODE) ? l_node : l_cell );
int const l_ez = ((ez_type[zdir] == NODE) ? l_node_v : l_cell_v);
int const l_bx = ((bx_type[zdir] == NODE) ? l_node_v : l_cell_v);
int const l_by = ((by_type[zdir] == NODE) ? l_node_v : l_cell_v);
int const l_bz = ((bz_type[zdir] == NODE) ? l_node : l_cell );
// Each field is gathered in a separate block of
// AMREX_SPACEDIM nested loops because the deposition
// order can differ for each component of each field
// when galerkin_interpolation is set to 1
#if defined(WARPX_DIM_1D_Z)
// Gather field on particle Eyp from field on grid ey_arr
// Gather field on particle Exp from field on grid ex_arr
// Gather field on particle Bzp from field on grid bz_arr
for (int iz=0; iz<=depos_order; iz++){
Eyp += sz_ey[iz]*
ey_arr(lo.x+l_ey+iz, 0, 0, 0);
Exp += sz_ex[iz]*
ex_arr(lo.x+l_ex+iz, 0, 0, 0);
Bzp += sz_bz[iz]*
bz_arr(lo.x+l_bz+iz, 0, 0, 0);
}
// Gather field on particle Byp from field on grid by_arr
// Gather field on particle Ezp from field on grid ez_arr
// Gather field on particle Bxp from field on grid bx_arr
for (int iz=0; iz<=depos_order-galerkin_interpolation; iz++){
Ezp += sz_ez[iz]*
ez_arr(lo.x+l_ez+iz, 0, 0, 0);
Bxp += sz_bx[iz]*
bx_arr(lo.x+l_bx+iz, 0, 0, 0);
Byp += sz_by[iz]*
by_arr(lo.x+l_by+iz, 0, 0, 0);
}
#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
// Gather field on particle Eyp from field on grid ey_arr
for (int iz=0; iz<=depos_order; iz++){
for (int ix=0; ix<=depos_order; ix++){
Eyp += sx_ey[ix]*sz_ey[iz]*
ey_arr(lo.x+j_ey+ix, lo.y+l_ey+iz, 0, 0);
}
}
// Gather field on particle Exp from field on grid ex_arr
// Gather field on particle Bzp from field on grid bz_arr
for (int iz=0; iz<=depos_order; iz++){
for (int ix=0; ix<=depos_order-galerkin_interpolation; ix++){
Exp += sx_ex[ix]*sz_ex[iz]*
ex_arr(lo.x+j_ex+ix, lo.y+l_ex+iz, 0, 0);
Bzp += sx_bz[ix]*sz_bz[iz]*
bz_arr(lo.x+j_bz+ix, lo.y+l_bz+iz, 0, 0);
}
}
// Gather field on particle Ezp from field on grid ez_arr
// Gather field on particle Bxp from field on grid bx_arr
for (int iz=0; iz<=depos_order-galerkin_interpolation; iz++){
for (int ix=0; ix<=depos_order; ix++){
Ezp += sx_ez[ix]*sz_ez[iz]*
ez_arr(lo.x+j_ez+ix, lo.y+l_ez+iz, 0, 0);
Bxp += sx_bx[ix]*sz_bx[iz]*
bx_arr(lo.x+j_bx+ix, lo.y+l_bx+iz, 0, 0);
}
}
// Gather field on particle Byp from field on grid by_arr
for (int iz=0; iz<=depos_order-galerkin_interpolation; iz++){
for (int ix=0; ix<=depos_order-galerkin_interpolation; ix++){
Byp += sx_by[ix]*sz_by[iz]*
by_arr(lo.x+j_by+ix, lo.y+l_by+iz, 0, 0);
}
}
#ifdef WARPX_DIM_RZ
amrex::Real costheta;
amrex::Real sintheta;
if (rp > 0.) {
costheta = xp/rp;
sintheta = yp/rp;
} else {
costheta = 1.;
sintheta = 0.;
}
const Complex xy0 = Complex{costheta, -sintheta};
Complex xy = xy0;
for (int imode=1 ; imode < n_rz_azimuthal_modes ; imode++) {
// Gather field on particle Eyp from field on grid ey_arr
for (int iz=0; iz<=depos_order; iz++){
for (int ix=0; ix<=depos_order; ix++){
const amrex::Real dEy = (+ ey_arr(lo.x+j_ey+ix, lo.y+l_ey+iz, 0, 2*imode-1)*xy.real()
- ey_arr(lo.x+j_ey+ix, lo.y+l_ey+iz, 0, 2*imode)*xy.imag());
Eyp += sx_ey[ix]*sz_ey[iz]*dEy;
}
}
// Gather field on particle Exp from field on grid ex_arr
// Gather field on particle Bzp from field on grid bz_arr
for (int iz=0; iz<=depos_order; iz++){
for (int ix=0; ix<=depos_order-galerkin_interpolation; ix++){
const amrex::Real dEx = (+ ex_arr(lo.x+j_ex+ix, lo.y+l_ex+iz, 0, 2*imode-1)*xy.real()
- ex_arr(lo.x+j_ex+ix, lo.y+l_ex+iz, 0, 2*imode)*xy.imag());
Exp += sx_ex[ix]*sz_ex[iz]*dEx;
const amrex::Real dBz = (+ bz_arr(lo.x+j_bz+ix, lo.y+l_bz+iz, 0, 2*imode-1)*xy.real()
- bz_arr(lo.x+j_bz+ix, lo.y+l_bz+iz, 0, 2*imode)*xy.imag());
Bzp += sx_bz[ix]*sz_bz[iz]*dBz;
}
}
// Gather field on particle Ezp from field on grid ez_arr
// Gather field on particle Bxp from field on grid bx_arr
for (int iz=0; iz<=depos_order-galerkin_interpolation; iz++){
for (int ix=0; ix<=depos_order; ix++){
const amrex::Real dEz = (+ ez_arr(lo.x+j_ez+ix, lo.y+l_ez+iz, 0, 2*imode-1)*xy.real()
- ez_arr(lo.x+j_ez+ix, lo.y+l_ez+iz, 0, 2*imode)*xy.imag());
Ezp += sx_ez[ix]*sz_ez[iz]*dEz;
const amrex::Real dBx = (+ bx_arr(lo.x+j_bx+ix, lo.y+l_bx+iz, 0, 2*imode-1)*xy.real()
- bx_arr(lo.x+j_bx+ix, lo.y+l_bx+iz, 0, 2*imode)*xy.imag());
Bxp += sx_bx[ix]*sz_bx[iz]*dBx;
}
}
// Gather field on particle Byp from field on grid by_arr
for (int iz=0; iz<=depos_order-galerkin_interpolation; iz++){
for (int ix=0; ix<=depos_order-galerkin_interpolation; ix++){
const amrex::Real dBy = (+ by_arr(lo.x+j_by+ix, lo.y+l_by+iz, 0, 2*imode-1)*xy.real()
- by_arr(lo.x+j_by+ix, lo.y+l_by+iz, 0, 2*imode)*xy.imag());
Byp += sx_by[ix]*sz_by[iz]*dBy;
}
}
xy = xy*xy0;
}
// Convert Exp and Eyp (which are actually Er and Etheta) to Ex and Ey
const amrex::Real Exp_save = Exp;
Exp = costheta*Exp - sintheta*Eyp;
Eyp = costheta*Eyp + sintheta*Exp_save;
const amrex::Real Bxp_save = Bxp;
Bxp = costheta*Bxp - sintheta*Byp;
Byp = costheta*Byp + sintheta*Bxp_save;
#endif
#else // defined(WARPX_DIM_3D)
// Gather field on particle Exp from field on grid ex_arr
for (int iz=0; iz<=depos_order; iz++){
for (int iy=0; iy<=depos_order; iy++){
for (int ix=0; ix<= depos_order - galerkin_interpolation; ix++){
Exp += sx_ex[ix]*sy_ex[iy]*sz_ex[iz]*
ex_arr(lo.x+j_ex+ix, lo.y+k_ex+iy, lo.z+l_ex+iz);
}
}
}
// Gather field on particle Eyp from field on grid ey_arr
for (int iz=0; iz<=depos_order; iz++){
for (int iy=0; iy<= depos_order - galerkin_interpolation; iy++){
for (int ix=0; ix<=depos_order; ix++){
Eyp += sx_ey[ix]*sy_ey[iy]*sz_ey[iz]*
ey_arr(lo.x+j_ey+ix, lo.y+k_ey+iy, lo.z+l_ey+iz);
}
}
}
// Gather field on particle Ezp from field on grid ez_arr
for (int iz=0; iz<= depos_order - galerkin_interpolation; iz++){
for (int iy=0; iy<=depos_order; iy++){
for (int ix=0; ix<=depos_order; ix++){
Ezp += sx_ez[ix]*sy_ez[iy]*sz_ez[iz]*
ez_arr(lo.x+j_ez+ix, lo.y+k_ez+iy, lo.z+l_ez+iz);
}
}
}
// Gather field on particle Bzp from field on grid bz_arr
for (int iz=0; iz<=depos_order; iz++){
for (int iy=0; iy<= depos_order - galerkin_interpolation; iy++){
for (int ix=0; ix<= depos_order - galerkin_interpolation; ix++){
Bzp += sx_bz[ix]*sy_bz[iy]*sz_bz[iz]*
bz_arr(lo.x+j_bz+ix, lo.y+k_bz+iy, lo.z+l_bz+iz);
}
}
}
// Gather field on particle Byp from field on grid by_arr
for (int iz=0; iz<= depos_order - galerkin_interpolation; iz++){
for (int iy=0; iy<=depos_order; iy++){
for (int ix=0; ix<= depos_order - galerkin_interpolation; ix++){
Byp += sx_by[ix]*sy_by[iy]*sz_by[iz]*
by_arr(lo.x+j_by+ix, lo.y+k_by+iy, lo.z+l_by+iz);
}
}
}
// Gather field on particle Bxp from field on grid bx_arr
for (int iz=0; iz<= depos_order - galerkin_interpolation; iz++){
for (int iy=0; iy<= depos_order - galerkin_interpolation; iy++){
for (int ix=0; ix<=depos_order; ix++){
Bxp += sx_bx[ix]*sy_bx[iy]*sz_bx[iz]*
bx_arr(lo.x+j_bx+ix, lo.y+k_bx+iy, lo.z+l_bx+iz);
}
}
}
#endif
}
/**
* \brief Field gather for particles
*
* \tparam depos_order deposition order
* \tparam lower_in_v lower shape order in parallel direction (Galerkin)
* \param getPosition A functor for returning the particle position.
* \param getExternalEB A functor for assigning the external E and B fields.
* \param Exp,Eyp,Ezp Pointer to array of electric field on particles.
* \param Bxp,Byp,Bzp Pointer to array of magnetic field on particles.
* \param exfab,eyfab,ezfab Array4 of the electric field, either full array or tile.
* \param bxfab,byfab,bzfab Array4 of the magnetic field, either full array or tile.
* \param np_to_gather Number of particles for which field is gathered.
* \param dx 3D cell size
* \param xyzmin Physical lower bounds of domain.
* \param lo Index lower bounds of domain.
* \param n_rz_azimuthal_modes Number of azimuthal modes when using RZ geometry
*/
template <int depos_order, int lower_in_v>
void doGatherShapeN(const GetParticlePosition& getPosition,
const GetExternalEBField& getExternalEB,
amrex::ParticleReal * const Exp, amrex::ParticleReal * const Eyp,
amrex::ParticleReal * const Ezp, amrex::ParticleReal * const Bxp,
amrex::ParticleReal * const Byp, amrex::ParticleReal * const Bzp,
amrex::FArrayBox const * const exfab,
amrex::FArrayBox const * const eyfab,
amrex::FArrayBox const * const ezfab,
amrex::FArrayBox const * const bxfab,
amrex::FArrayBox const * const byfab,
amrex::FArrayBox const * const bzfab,
const long np_to_gather,
const std::array<amrex::Real, 3>& dx,
const std::array<amrex::Real, 3> xyzmin,
const amrex::Dim3 lo,
const int n_rz_azimuthal_modes)
{
amrex::GpuArray<amrex::Real, 3> dx_arr = {dx[0], dx[1], dx[2]};
amrex::GpuArray<amrex::Real, 3> xyzmin_arr = {xyzmin[0], xyzmin[1], xyzmin[2]};
amrex::Array4<const amrex::Real> const& ex_arr = exfab->array();
amrex::Array4<const amrex::Real> const& ey_arr = eyfab->array();
amrex::Array4<const amrex::Real> const& ez_arr = ezfab->array();
amrex::Array4<const amrex::Real> const& bx_arr = bxfab->array();
amrex::Array4<const amrex::Real> const& by_arr = byfab->array();
amrex::Array4<const amrex::Real> const& bz_arr = bzfab->array();
amrex::IndexType const ex_type = exfab->box().ixType();
amrex::IndexType const ey_type = eyfab->box().ixType();
amrex::IndexType const ez_type = ezfab->box().ixType();
amrex::IndexType const bx_type = bxfab->box().ixType();
amrex::IndexType const by_type = byfab->box().ixType();
amrex::IndexType const bz_type = bzfab->box().ixType();
// Loop over particles and gather fields from
// {e,b}{x,y,z}_arr to {E,B}{xyz}p.
amrex::ParallelFor(
np_to_gather,
[=] AMREX_GPU_DEVICE (long ip) {
amrex::ParticleReal xp, yp, zp;
getPosition(ip, xp, yp, zp);
getExternalEB(ip, Exp[ip], Eyp[ip], Ezp[ip], Bxp[ip], Byp[ip], Bzp[ip]);
doGatherShapeN<depos_order, lower_in_v>(
xp, yp, zp, Exp[ip], Eyp[ip], Ezp[ip], Bxp[ip], Byp[ip], Bzp[ip],
ex_arr, ey_arr, ez_arr, bx_arr, by_arr, bz_arr,
ex_type, ey_type, ez_type, bx_type, by_type, bz_type,
dx_arr, xyzmin_arr, lo, n_rz_azimuthal_modes);
}
);
}
/**
* \brief Field gather for a single particle
*
* \param xp,yp,zp Particle position coordinates
* \param Exp,Eyp,Ezp Electric field on particles.
* \param Bxp,Byp,Bzp Magnetic field on particles.
* \param ex_arr,ey_arr,ez_arr Array4 of the electric field, either full array or tile.
* \param bx_arr,by_arr,bz_arr Array4 of the magnetic field, either full array or tile.
* \param ex_type,ey_type,ez_type IndexType of the electric field
* \param bx_type,by_type,bz_type IndexType of the magnetic field
* \param dx_arr 3D cell spacing
* \param xyzmin_arr Physical lower bounds of domain in x, y, z.
* \param lo Index lower bounds of domain.
* \param n_rz_azimuthal_modes Number of azimuthal modes when using RZ geometry
* \param nox order of the particle shape function
* \param galerkin_interpolation whether to use lower order in v
*/
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_arr,
const amrex::GpuArray<amrex::Real, 3>& xyzmin_arr,
const amrex::Dim3& lo,
const int n_rz_azimuthal_modes,
const int nox,
const bool galerkin_interpolation)
{
if (galerkin_interpolation) {
if (nox == 1) {
doGatherShapeN<1,1>(xp, yp, zp, Exp, Eyp, Ezp, Bxp, Byp, Bzp,
ex_arr, ey_arr, ez_arr, bx_arr, by_arr, bz_arr,
ex_type, ey_type, ez_type, bx_type, by_type, bz_type,
dx_arr, xyzmin_arr, lo, n_rz_azimuthal_modes);
} else if (nox == 2) {
doGatherShapeN<2,1>(xp, yp, zp, Exp, Eyp, Ezp, Bxp, Byp, Bzp,
ex_arr, ey_arr, ez_arr, bx_arr, by_arr, bz_arr,
ex_type, ey_type, ez_type, bx_type, by_type, bz_type,
dx_arr, xyzmin_arr, lo, n_rz_azimuthal_modes);
} else if (nox == 3) {
doGatherShapeN<3,1>(xp, yp, zp, Exp, Eyp, Ezp, Bxp, Byp, Bzp,
ex_arr, ey_arr, ez_arr, bx_arr, by_arr, bz_arr,
ex_type, ey_type, ez_type, bx_type, by_type, bz_type,
dx_arr, xyzmin_arr, lo, n_rz_azimuthal_modes);
}
} else {
if (nox == 1) {
doGatherShapeN<1,0>(xp, yp, zp, Exp, Eyp, Ezp, Bxp, Byp, Bzp,
ex_arr, ey_arr, ez_arr, bx_arr, by_arr, bz_arr,
ex_type, ey_type, ez_type, bx_type, by_type, bz_type,
dx_arr, xyzmin_arr, lo, n_rz_azimuthal_modes);
} else if (nox == 2) {
doGatherShapeN<2,0>(xp, yp, zp, Exp, Eyp, Ezp, Bxp, Byp, Bzp,
ex_arr, ey_arr, ez_arr, bx_arr, by_arr, bz_arr,
ex_type, ey_type, ez_type, bx_type, by_type, bz_type,
dx_arr, xyzmin_arr, lo, n_rz_azimuthal_modes);
} else if (nox == 3) {
doGatherShapeN<3,0>(xp, yp, zp, Exp, Eyp, Ezp, Bxp, Byp, Bzp,
ex_arr, ey_arr, ez_arr, bx_arr, by_arr, bz_arr,
ex_type, ey_type, ez_type, bx_type, by_type, bz_type,
dx_arr, xyzmin_arr, lo, n_rz_azimuthal_modes);
}
}
}
#endif // FIELDGATHER_H_
|