ARTEMIS
ChargeDeposition.H
Go to the documentation of this file.
1 /* Copyright 2019 Axel Huebl, Andrew Myers, David Grote, Maxence Thevenet
2  * Weiqun Zhang
3  *
4  * This file is part of WarpX.
5  *
6  * License: BSD-3-Clause-LBNL
7  */
8 #ifndef CHARGEDEPOSITION_H_
9 #define CHARGEDEPOSITION_H_
10 
14 #include "Particles/ShapeFactors.H"
17 #ifdef WARPX_DIM_RZ
18 # include "Utils/WarpX_Complex.H"
19 #endif
20 
21 #include <AMReX.H>
22 
23 /* \brief Perform charge deposition on a tile
24  * \param GetPosition : A functor for returning the particle position.
25  * \param wp : Pointer to array of particle weights.
26  * \param ion_lev : Pointer to array of particle ionization level. This is
27  required to have the charge of each macroparticle
28  since q is a scalar. For non-ionizable species,
29  ion_lev is a null pointer.
30  * \param rho_fab : FArrayBox of charge density, either full array or tile.
31  * \param np_to_depose : Number of particles for which current is deposited.
32  * \param dx : 3D cell size
33  * \param xyzmin : Physical lower bounds of domain.
34  * \param lo : Index lower bounds of domain.
35  * \param q : species charge.
36  * \param n_rz_azimuthal_modes: Number of azimuthal modes when using RZ geometry.
37  * \param cost: Pointer to (load balancing) cost corresponding to box where present particles deposit current.
38  * \param load_balance_costs_update_algo: Selected method for updating load balance costs.
39  */
40 template <int depos_order>
42  const amrex::ParticleReal * const wp,
43  const int * const ion_lev,
44  amrex::FArrayBox& rho_fab,
45  const long np_to_depose,
46  const std::array<amrex::Real,3>& dx,
47  const std::array<amrex::Real, 3> xyzmin,
48  const amrex::Dim3 lo,
49  const amrex::Real q,
50  const int n_rz_azimuthal_modes,
51  amrex::Real* cost,
52  const long load_balance_costs_update_algo)
53 {
54  using namespace amrex;
55 
56 #if !defined(AMREX_USE_GPU)
57  amrex::ignore_unused(cost, load_balance_costs_update_algo);
58 #endif
59 
60  // Whether ion_lev is a null pointer (do_ionization=0) or a real pointer
61  // (do_ionization=1)
62  const bool do_ionization = ion_lev;
63  const amrex::Real dzi = 1.0_rt/dx[2];
64 #if defined(WARPX_DIM_1D_Z)
65  const amrex::Real invvol = dzi;
66 #endif
67 #if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
68  const amrex::Real dxi = 1.0_rt/dx[0];
69  const amrex::Real invvol = dxi*dzi;
70 #elif defined(WARPX_DIM_3D)
71  const amrex::Real dxi = 1.0_rt/dx[0];
72  const amrex::Real dyi = 1.0_rt/dx[1];
73  const amrex::Real invvol = dxi*dyi*dzi;
74 #endif
75 
76 #if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) || defined(WARPX_DIM_3D)
77  const amrex::Real xmin = xyzmin[0];
78 #endif
79 #if defined(WARPX_DIM_3D)
80  const amrex::Real ymin = xyzmin[1];
81 #endif
82  const amrex::Real zmin = xyzmin[2];
83 
84  amrex::Array4<amrex::Real> const& rho_arr = rho_fab.array();
85  amrex::IntVect const rho_type = rho_fab.box().type();
86 
87  constexpr int NODE = amrex::IndexType::NODE;
88  constexpr int CELL = amrex::IndexType::CELL;
89 
90  // Loop over particles and deposit into rho_fab
91 #if defined(WARPX_USE_GPUCLOCK)
92  amrex::Real* cost_real = nullptr;
93  if( load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::GpuClock) {
94  cost_real = (amrex::Real *) amrex::The_Managed_Arena()->alloc(sizeof(amrex::Real));
95  *cost_real = 0.;
96  }
97 #endif
99  np_to_depose,
100  [=] AMREX_GPU_DEVICE (long ip) {
101 #if defined(WARPX_USE_GPUCLOCK)
102  const auto KernelTimer = ablastr::parallelization::KernelTimer(
103  cost && (load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::GpuClock),
104  cost_real);
105 #endif
106  // --- Get particle quantities
107  amrex::Real wq = q*wp[ip]*invvol;
108  if (do_ionization){
109  wq *= ion_lev[ip];
110  }
111 
112  amrex::ParticleReal xp, yp, zp;
113  GetPosition(ip, xp, yp, zp);
114 
115  // --- Compute shape factors
116  Compute_shape_factor< depos_order > const compute_shape_factor;
117 #if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) || defined(WARPX_DIM_3D)
118  // x direction
119  // Get particle position in grid coordinates
120 #if defined(WARPX_DIM_RZ)
121  const amrex::Real rp = std::sqrt(xp*xp + yp*yp);
122  amrex::Real costheta;
123  amrex::Real sintheta;
124  if (rp > 0.) {
125  costheta = xp/rp;
126  sintheta = yp/rp;
127  } else {
128  costheta = 1._rt;
129  sintheta = 0._rt;
130  }
131  const Complex xy0 = Complex{costheta, sintheta};
132  const amrex::Real x = (rp - xmin)*dxi;
133 #else
134  const amrex::Real x = (xp - xmin)*dxi;
135 #endif
136 
137  // Compute shape factor along x
138  // i: leftmost grid point that the particle touches
139  amrex::Real sx[depos_order + 1] = {0._rt};
140  int i = 0;
141  if (rho_type[0] == NODE) {
142  i = compute_shape_factor(sx, x);
143  } else if (rho_type[0] == CELL) {
144  i = compute_shape_factor(sx, x - 0.5_rt);
145  }
146 #endif //defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) || defined(WARPX_DIM_3D)
147 #if defined(WARPX_DIM_3D)
148  // y direction
149  const amrex::Real y = (yp - ymin)*dyi;
150  amrex::Real sy[depos_order + 1] = {0._rt};
151  int j = 0;
152  if (rho_type[1] == NODE) {
153  j = compute_shape_factor(sy, y);
154  } else if (rho_type[1] == CELL) {
155  j = compute_shape_factor(sy, y - 0.5_rt);
156  }
157 #endif
158  // z direction
159  const amrex::Real z = (zp - zmin)*dzi;
160  amrex::Real sz[depos_order + 1] = {0._rt};
161  int k = 0;
162  if (rho_type[WARPX_ZINDEX] == NODE) {
163  k = compute_shape_factor(sz, z);
164  } else if (rho_type[WARPX_ZINDEX] == CELL) {
165  k = compute_shape_factor(sz, z - 0.5_rt);
166  }
167 
168  // Deposit charge into rho_arr
169 #if defined(WARPX_DIM_1D_Z)
170  for (int iz=0; iz<=depos_order; iz++){
172  &rho_arr(lo.x+k+iz, 0, 0, 0),
173  sz[iz]*wq);
174  }
175 #endif
176 #if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
177  for (int iz=0; iz<=depos_order; iz++){
178  for (int ix=0; ix<=depos_order; ix++){
180  &rho_arr(lo.x+i+ix, lo.y+k+iz, 0, 0),
181  sx[ix]*sz[iz]*wq);
182 #if defined(WARPX_DIM_RZ)
183  Complex xy = xy0; // Throughout the following loop, xy takes the value e^{i m theta}
184  for (int imode=1 ; imode < n_rz_azimuthal_modes ; imode++) {
185  // The factor 2 on the weighting comes from the normalization of the modes
186  amrex::Gpu::Atomic::AddNoRet( &rho_arr(lo.x+i+ix, lo.y+k+iz, 0, 2*imode-1), 2._rt*sx[ix]*sz[iz]*wq*xy.real());
187  amrex::Gpu::Atomic::AddNoRet( &rho_arr(lo.x+i+ix, lo.y+k+iz, 0, 2*imode ), 2._rt*sx[ix]*sz[iz]*wq*xy.imag());
188  xy = xy*xy0;
189  }
190 #endif
191  }
192  }
193 #elif defined(WARPX_DIM_3D)
194  for (int iz=0; iz<=depos_order; iz++){
195  for (int iy=0; iy<=depos_order; iy++){
196  for (int ix=0; ix<=depos_order; ix++){
198  &rho_arr(lo.x+i+ix, lo.y+j+iy, lo.z+k+iz),
199  sx[ix]*sy[iy]*sz[iz]*wq);
200  }
201  }
202  }
203 #endif
204  }
205  );
206 #if defined(WARPX_USE_GPUCLOCK)
207  if (cost && load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::GpuClock) {
209  *cost += *cost_real;
210  amrex::The_Managed_Arena()->free(cost_real);
211  }
212 #endif
213 
214 #ifndef WARPX_DIM_RZ
215  amrex::ignore_unused(n_rz_azimuthal_modes);
216 #endif
217 }
218 
219 /* \brief Perform charge deposition on a tile using shared memory
220  * \param GetPosition : A functor for returning the particle position.
221  * \param wp : Pointer to array of particle weights.
222  * \param ion_lev : Pointer to array of particle ionization level. This is
223  required to have the charge of each macroparticle
224  since q is a scalar. For non-ionizable species,
225  ion_lev is a null pointer.
226  * \param rho_fab : FArrayBox of charge density, either full array or tile.
227  * \param np_to_deposit : Number of particles for which charge is deposited.
228  * \param dx : 3D cell size
229  * \param xyzmin : Physical lower bounds of domain.
230  * \param lo : Index lower bounds of domain.
231  * \param q : species charge.
232  * \param n_rz_azimuthal_modes: Number of azimuthal modes when using RZ geometry.
233  * \param cost: Pointer to (load balancing) cost corresponding to box where present particles deposit current.
234  * \param load_balance_costs_update_algo: Selected method for updating load balance costs.
235  */
236 template <int depos_order>
238  const amrex::ParticleReal * const wp,
239  const int * const ion_lev,
240  amrex::FArrayBox& rho_fab,
241  const amrex::IntVect& ix_type,
242  const long np_to_deposit,
243  const std::array<amrex::Real,3>& dx,
244  const std::array<amrex::Real, 3> xyzmin,
245  const amrex::Dim3 lo,
246  const amrex::Real q,
247  const int n_rz_azimuthal_modes,
248  amrex::Real* cost,
249  const long load_balance_costs_update_algo,
251  const amrex::Box& box,
252  const amrex::Geometry& geom,
253  const amrex::IntVect& a_tbox_max_size)
254 {
255  using namespace amrex;
256 
257  auto permutation = a_bins.permutationPtr();
258 
259 #if !defined(AMREX_USE_GPU)
260  amrex::ignore_unused(ix_type, cost, load_balance_costs_update_algo, a_bins, box, geom, a_tbox_max_size);
261 #endif
262 
263  // Whether ion_lev is a null pointer (do_ionization=0) or a real pointer
264  // (do_ionization=1)
265  const bool do_ionization = ion_lev;
266  const amrex::Real dzi = 1.0_rt/dx[2];
267 #if defined(WARPX_DIM_1D_Z)
268  const amrex::Real invvol = dzi;
269 #endif
270 #if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
271  const amrex::Real dxi = 1.0_rt/dx[0];
272  const amrex::Real invvol = dxi*dzi;
273 #elif defined(WARPX_DIM_3D)
274  const amrex::Real dxi = 1.0_rt/dx[0];
275  const amrex::Real dyi = 1.0_rt/dx[1];
276  const amrex::Real invvol = dxi*dyi*dzi;
277 #endif
278 
279 #if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) || defined(WARPX_DIM_3D)
280  const amrex::Real xmin = xyzmin[0];
281 #endif
282 #if defined(WARPX_DIM_3D)
283  const amrex::Real ymin = xyzmin[1];
284 #endif
285  const amrex::Real zmin = xyzmin[2];
286 
287  amrex::Array4<amrex::Real> const& rho_arr = rho_fab.array();
288  auto rho_box = rho_fab.box();
289  amrex::IntVect const rho_type = rho_box.type();
290 
291  constexpr int NODE = amrex::IndexType::NODE;
292  constexpr int CELL = amrex::IndexType::CELL;
293 
294  // Loop over particles and deposit into rho_fab
295 #if defined(WARPX_USE_GPUCLOCK)
296  amrex::Real* cost_real = nullptr;
297  if( load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::GpuClock) {
298  cost_real = (amrex::Real *) amrex::The_Managed_Arena()->alloc(sizeof(amrex::Real));
299  *cost_real = 0.;
300  }
301 #endif
302 
303 #if defined(AMREX_USE_CUDA) || defined(AMREX_USE_HIP)
304  const auto dxiarr = geom.InvCellSizeArray();
305  const auto plo = geom.ProbLoArray();
306  const auto domain = geom.Domain();
307 
308  const amrex::Box sample_tbox(IntVect(AMREX_D_DECL(0, 0, 0)), a_tbox_max_size - 1);
309  amrex::Box sample_tbox_x = convert(sample_tbox, ix_type);
310 
311  sample_tbox_x.grow(depos_order);
312 
313  const auto npts = sample_tbox_x.numPts();
314 
315  const int nblocks = a_bins.numBins();
316  const auto offsets_ptr = a_bins.offsetsPtr();
317  const int threads_per_block = 256;
318 
319  std::size_t shared_mem_bytes = npts*sizeof(amrex::Real);
321 
322  const std::size_t max_shared_mem_bytes = amrex::Gpu::Device::sharedMemPerBlock();
323 
324  WARPX_ALWAYS_ASSERT_WITH_MESSAGE(shared_mem_bytes <= max_shared_mem_bytes,
325  "Tile size too big for GPU shared memory charge deposition");
326 #endif
327 
328 #if defined(AMREX_USE_CUDA) || defined(AMREX_USE_HIP)
329  amrex::ignore_unused(np_to_deposit);
330  // Loop with one block per tile (the shared memory is allocated on a per-block basis)
331  // The threads within each block loop over the particles of its tile
332  // (Each threads processes a different set of particles.)
334  nblocks, threads_per_block, shared_mem_bytes, amrex::Gpu::gpuStream(),
335  [=] AMREX_GPU_DEVICE () noexcept
336 #else // defined(AMREX_USE_CUDA) || defined(AMREX_USE_HIP)
337  amrex::ParallelFor(np_to_deposit, [=] AMREX_GPU_DEVICE (long ip_orig) noexcept
338 #endif
339  {
340 #if defined(AMREX_USE_CUDA) || defined(AMREX_USE_HIP)
341  const int bin_id = blockIdx.x;
342  const unsigned int bin_start = offsets_ptr[bin_id];
343  const unsigned int bin_stop = offsets_ptr[bin_id+1];
344 
345  if (bin_start == bin_stop) { return; }
346 
347  amrex::Box buffer_box;
348  {
349  ParticleReal xp, yp, zp;
350  GetPosition(permutation[bin_start], xp, yp, zp);
351 #if defined(WARPX_DIM_3D)
352  IntVect iv = IntVect(int(amrex::Math::floor((xp-plo[0])*dxiarr[0])),
353  int(amrex::Math::floor((yp-plo[1])*dxiarr[1])),
354  int(amrex::Math::floor((zp-plo[2])*dxiarr[2])));
355 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
356  IntVect iv = IntVect(
357  int(amrex::Math::floor((xp-plo[0])*dxiarr[0])),
358  int(amrex::Math::floor((zp-plo[1])*dxiarr[1])));
359 #elif defined(WARPX_DIM_1D_Z)
360  IntVect iv = IntVect(int(amrex::Math::floor((zp-plo[0])*dxiarr[0])));
361 #endif
362  iv += domain.smallEnd();
363  getTileIndex(iv, box, true, bin_size, buffer_box);
364  }
365 
366  Box tbx = convert( buffer_box, ix_type);
367  tbx.grow(depos_order);
368 
370  amrex::Real* const shared = gsm.dataPtr();
371 
372  amrex::Array4<amrex::Real> buf(shared, amrex::begin(tbx), amrex::end(tbx), 1);
373 
374  // Zero-initialize the temporary array in shared memory
375  volatile amrex::Real* vs = shared;
376  for (int i = threadIdx.x; i < tbx.numPts(); i += blockDim.x) {
377  vs[i] = 0.0;
378  }
379  __syncthreads();
380 #else
381  amrex::Array4<amrex::Real> const &buf = rho_arr;
382 #endif // defined(AMREX_USE_CUDA) || defined(AMREX_USE_HIP)
383 
384 #if defined(AMREX_USE_CUDA) || defined(AMREX_USE_HIP)
385  // Loop over macroparticles: each threads loops over particles with a stride of `blockDim.x`
386  for (unsigned int ip_orig = bin_start + threadIdx.x; ip_orig < bin_stop; ip_orig += blockDim.x)
387 #endif
388  {
389  unsigned int ip = permutation[ip_orig];
390 
391 #if defined(WARPX_USE_GPUCLOCK)
392  const auto KernelTimer = ablastr::parallelization::KernelTimer(
393  cost && (load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::GpuClock),
394  cost_real);
395 #endif
396  // --- Get particle quantities
397  amrex::Real wq = q*wp[ip]*invvol;
398  if (do_ionization){
399  wq *= ion_lev[ip];
400  }
401 
402  amrex::ParticleReal xp, yp, zp;
403  GetPosition(ip, xp, yp, zp);
404 
405  // --- Compute shape factors
406  Compute_shape_factor< depos_order > const compute_shape_factor;
407 #if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) || defined(WARPX_DIM_3D)
408  // x direction
409  // Get particle position in grid coordinates
410 #if defined(WARPX_DIM_RZ)
411  const amrex::Real rp = std::sqrt(xp*xp + yp*yp);
412  amrex::Real costheta;
413  amrex::Real sintheta;
414  if (rp > 0.) {
415  costheta = xp/rp;
416  sintheta = yp/rp;
417  } else {
418  costheta = 1._rt;
419  sintheta = 0._rt;
420  }
421  const Complex xy0 = Complex{costheta, sintheta};
422  const amrex::Real x = (rp - xmin)*dxi;
423 #else
424  const amrex::Real x = (xp - xmin)*dxi;
425 #endif
426 
427  // Compute shape factor along x
428  // i: leftmost grid point that the particle touches
429  amrex::Real sx[depos_order + 1] = {0._rt};
430  int i = 0;
431  if (rho_type[0] == NODE) {
432  i = compute_shape_factor(sx, x);
433  } else if (rho_type[0] == CELL) {
434  i = compute_shape_factor(sx, x - 0.5_rt);
435  }
436 #endif //defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ) || defined(WARPX_DIM_3D)
437 #if defined(WARPX_DIM_3D)
438  // y direction
439  const amrex::Real y = (yp - ymin)*dyi;
440  amrex::Real sy[depos_order + 1] = {0._rt};
441  int j = 0;
442  if (rho_type[1] == NODE) {
443  j = compute_shape_factor(sy, y);
444  } else if (rho_type[1] == CELL) {
445  j = compute_shape_factor(sy, y - 0.5_rt);
446  }
447 #endif
448  // z direction
449  const amrex::Real z = (zp - zmin)*dzi;
450  amrex::Real sz[depos_order + 1] = {0._rt};
451  int k = 0;
452  if (rho_type[WARPX_ZINDEX] == NODE) {
453  k = compute_shape_factor(sz, z);
454  } else if (rho_type[WARPX_ZINDEX] == CELL) {
455  k = compute_shape_factor(sz, z - 0.5_rt);
456  }
457 
458  // Deposit charge into buf
459 #if defined(WARPX_DIM_1D_Z)
460  for (int iz=0; iz<=depos_order; iz++){
462  &buf(lo.x+k+iz, 0, 0, 0),
463  sz[iz]*wq);
464  }
465 #endif
466 #if defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
467  for (int iz=0; iz<=depos_order; iz++){
468  for (int ix=0; ix<=depos_order; ix++){
470  &buf(lo.x+i+ix, lo.y+k+iz, 0, 0),
471  sx[ix]*sz[iz]*wq);
472 #if defined(WARPX_DIM_RZ)
473  Complex xy = xy0; // Throughout the following loop, xy takes the value e^{i m theta}
474  for (int imode=1 ; imode < n_rz_azimuthal_modes ; imode++) {
475  // The factor 2 on the weighting comes from the normalization of the modes
476  amrex::Gpu::Atomic::AddNoRet( &buf(lo.x+i+ix, lo.y+k+iz, 0, 2*imode-1), 2._rt*sx[ix]*sz[iz]*wq*xy.real());
477  amrex::Gpu::Atomic::AddNoRet( &buf(lo.x+i+ix, lo.y+k+iz, 0, 2*imode ), 2._rt*sx[ix]*sz[iz]*wq*xy.imag());
478  xy = xy*xy0;
479  }
480 #endif
481  }
482  }
483 #elif defined(WARPX_DIM_3D)
484  for (int iz=0; iz<=depos_order; iz++){
485  for (int iy=0; iy<=depos_order; iy++){
486  for (int ix=0; ix<=depos_order; ix++){
488  &buf(lo.x+i+ix, lo.y+j+iy, lo.z+k+iz),
489  sx[ix]*sy[iy]*sz[iz]*wq);
490  }
491  }
492  }
493 #endif
494  }
495 
496 #if defined(AMREX_USE_CUDA) || defined(AMREX_USE_HIP)
497  __syncthreads();
498 
499  addLocalToGlobal(tbx, rho_arr, buf);
500 #endif // defined(AMREX_USE_CUDA) || defined(AMREX_USE_HIP)
501  }
502  );
503 #if defined(WARPX_USE_GPUCLOCK)
504  if(cost && load_balance_costs_update_algo == LoadBalanceCostsUpdateAlgo::GpuClock) {
506  *cost += *cost_real;
507  amrex::The_Managed_Arena()->free(cost_real);
508  }
509 #endif
510 
511 #ifndef WARPX_DIM_RZ
512  amrex::ignore_unused(n_rz_azimuthal_modes);
513 #endif
514 }
515 
516 #endif // CHARGEDEPOSITION_H_
#define AMREX_GPU_DEVICE
#define AMREX_D_DECL(a, b, c)
void doChargeDepositionSharedShapeN(const GetParticlePosition &GetPosition, const amrex::ParticleReal *const wp, const int *const ion_lev, amrex::FArrayBox &rho_fab, const amrex::IntVect &ix_type, const long np_to_deposit, const std::array< amrex::Real, 3 > &dx, const std::array< amrex::Real, 3 > xyzmin, const amrex::Dim3 lo, const amrex::Real q, const int n_rz_azimuthal_modes, amrex::Real *cost, const long load_balance_costs_update_algo, const amrex::DenseBins< WarpXParticleContainer::ParticleType > &a_bins, const amrex::Box &box, const amrex::Geometry &geom, const amrex::IntVect &a_tbox_max_size)
Definition: ChargeDeposition.H:237
void doChargeDepositionShapeN(const GetParticlePosition &GetPosition, const amrex::ParticleReal *const wp, const int *const ion_lev, amrex::FArrayBox &rho_fab, const long np_to_depose, const std::array< amrex::Real, 3 > &dx, const std::array< amrex::Real, 3 > xyzmin, const amrex::Dim3 lo, const amrex::Real q, const int n_rz_azimuthal_modes, amrex::Real *cost, const long load_balance_costs_update_algo)
Definition: ChargeDeposition.H:41
#define WARPX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition: TextMsg.H:13
NODE
static amrex::IntVect shared_tilesize
tileSize to use for shared current deposition operations
Definition: WarpX.H:370
Defines a timer object to be used on GPU; measures summed thread cycles.
Definition: KernelTimer.H:27
virtual void free(void *pt)=0
virtual void * alloc(std::size_t sz)=0
const Box & box() const noexcept
AMREX_FORCE_INLINE Array4< T const > array() const noexcept
AMREX_GPU_HOST_DEVICE BoxND & grow(int i) noexcept
AMREX_GPU_HOST_DEVICE IntVectND< dim > type() const noexcept
AMREX_GPU_HOST_DEVICE Long numPts() const noexcept
GpuArray< Real, AMREX_SPACEDIM > InvCellSizeArray() const noexcept
index_type * offsetsPtr() noexcept
index_type * permutationPtr() noexcept
Long numBins() const noexcept
const Box & Domain() const noexcept
GpuArray< Real, AMREX_SPACEDIM > ProbLoArray() const noexcept
static std::size_t sharedMemPerBlock() noexcept
def y
Definition: Excitation_Flag_Generator.py:76
def z
Definition: Excitation_Flag_Generator.py:77
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void AddNoRet(T *sum, T value) noexcept
void streamSynchronize() noexcept
gpuStream_t gpuStream() noexcept
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
void launch(T const &n, L &&f) noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 end(BoxND< dim > const &box) noexcept
Arena * The_Managed_Arena()
IntVectND< AMREX_SPACEDIM > IntVect
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 begin(BoxND< dim > const &box) noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int getTileIndex(const IntVect &iv, const Box &box, const bool a_do_tiling, const IntVect &a_tile_size, Box &tbx)
i
Definition: check_interp_points_and_weights.py:174
int dx
Definition: stencil.py:436
Definition: ShapeFactors.H:27
Functor that can be used to extract the positions of the macroparticles inside a ParallelFor kernel.
Definition: GetAndSetPosition.H:53
@ GpuClock
Definition: WarpXAlgorithmSelection.H:138
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE constexpr T real() const noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE constexpr T imag() const noexcept
AMREX_GPU_DEVICE T * dataPtr() noexcept