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