ARTEMIS
DepositCharge.H
Go to the documentation of this file.
1 /* Copyright 2019-2021 Axel Huebl, Andrew Myers
2  *
3  * This file is part of WarpX.
4  *
5  * License: BSD-3-Clause-LBNL
6  */
7 #ifndef ABLASTR_DEPOSIT_CHARGE_H_
8 #define ABLASTR_DEPOSIT_CHARGE_H_
9 
12 #include "Particles/ShapeFactors.H"
14 #include "ablastr/utils/TextMsg.H"
15 
16 #include <AMReX.H>
17 
18 #include <optional>
19 
20 
22 {
23 
53 template< typename T_PC >
54 static void
55 deposit_charge (typename T_PC::ParIterType& pti,
56  typename T_PC::RealVector const& wp,
57  amrex::Real const charge,
58  int const * const ion_lev,
59  amrex::MultiFab* rho,
60  amrex::FArrayBox& local_rho,
61  int const particle_shape,
62  std::array<amrex::Real, 3> const & dx,
63  std::array<amrex::Real, 3> const & xyzmin,
64  int const n_rz_azimuthal_modes = 0,
65  std::optional<amrex::IntVect> num_rho_deposition_guards = std::nullopt,
66  std::optional<int> depos_lev = std::nullopt,
67  std::optional<amrex::IntVect> rel_ref_ratio = std::nullopt,
68  long const offset = 0,
69  std::optional<long> np_to_depose = std::nullopt,
70  int const icomp = 0, int const nc = 1,
71  amrex::Real * const AMREX_RESTRICT cost = nullptr,
72  long const load_balance_costs_update_algo = 0,
73  bool const do_device_synchronize = true)
74 {
75  // deposition guards
76  amrex::IntVect ng_rho = rho->nGrowVect();
77  if (num_rho_deposition_guards.has_value())
78  ng_rho = num_rho_deposition_guards.value();
79  ABLASTR_ALWAYS_ASSERT_WITH_MESSAGE(ng_rho <= rho->nGrowVect(),
80  "num_rho_deposition_guards are larger than allocated!");
81  // particle shape
82  auto const[nox, noy, noz] = std::array<int, 3>{particle_shape, particle_shape, particle_shape};
83 
84  // used for MR when we want to deposit for a subset of the particles on the level in the
85  // current box; with offset, we start at a later particle index
86  if (!np_to_depose.has_value())
87  np_to_depose = pti.numParticles();
88  ABLASTR_ALWAYS_ASSERT_WITH_MESSAGE(np_to_depose.value() + offset <= pti.numParticles(),
89  "np_to_depose + offset are out-of-bounds for particle iterator");
90 
91  int const lev = pti.GetLevel();
92  if (!depos_lev.has_value())
93  depos_lev = lev;
94  ABLASTR_ALWAYS_ASSERT_WITH_MESSAGE((depos_lev.value() == (lev-1)) ||
95  (depos_lev.value() == (lev )),
96  "Deposition buffers only work for lev or lev-1");
97  if (!rel_ref_ratio.has_value()) {
98  ABLASTR_ALWAYS_ASSERT_WITH_MESSAGE(lev == depos_lev,
99  "rel_ref_ratio must be set if lev != depos_lev");
100  rel_ref_ratio = amrex::IntVect(AMREX_D_DECL(1, 1, 1));
101  }
102 
103  // If no particles, do not do anything
104  if (np_to_depose == 0) return;
105 
106  // Extract deposition order and check that particles shape fits within the guard cells.
107  // NOTE: In specific situations where the staggering of rho and the charge deposition algorithm
108  // are not trivial, this check might be too strict and we might need to relax it, as currently
109  // done for the current deposition.
110 
111 #if defined(WARPX_DIM_1D_Z)
114  const amrex::IntVect shape_extent = amrex::IntVect(static_cast<int>(noz/2+1));
115 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
117  const amrex::IntVect shape_extent = amrex::IntVect(static_cast<int>(nox/2+1),
118  static_cast<int>(noz/2+1));
119 #elif defined(WARPX_DIM_3D)
120  const amrex::IntVect shape_extent = amrex::IntVect(static_cast<int>(nox/2+1),
121  static_cast<int>(noy/2+1),
122  static_cast<int>(noz/2+1));
123 #endif
124 
125  // On CPU: particles deposit on tile arrays, which have a small number of guard cells ng_rho
126  // On GPU: particles deposit directly on the rho array, which usually have a larger number of guard cells
127 #ifndef AMREX_USE_GPU
128  const amrex::IntVect range = ng_rho - shape_extent;
129 #else
130  const amrex::IntVect range = rho->nGrowVect() - shape_extent;
131 #endif
132  amrex::ignore_unused(range); // for release builds
134  amrex::numParticlesOutOfRange(pti, range) == 0,
135  "Particles shape does not fit within tile (CPU) or guard cells (GPU) used for charge deposition");
136 
137  ABLASTR_PROFILE_VAR_NS("ablastr::particles::deposit_charge::ChargeDeposition", blp_ppc_chd, do_device_synchronize);
138  ABLASTR_PROFILE_VAR_NS("ablastr::particles::deposit_charge::Accumulate", blp_accumulate, do_device_synchronize);
139 
140  // Get tile box where charge is deposited.
141  // The tile box is different when depositing in the buffers (depos_lev<lev)
142  // or when depositing inside the level (depos_lev=lev)
143  amrex::Box tilebox;
144  if (lev == depos_lev) {
145  tilebox = pti.tilebox();
146  } else {
147  tilebox = amrex::coarsen(pti.tilebox(), rel_ref_ratio.value());
148  }
149 
150 #ifndef AMREX_USE_GPU
151  // Staggered tile box
152  amrex::Box tb = amrex::convert( tilebox, rho->ixType().toIntVect() );
153 #endif
154 
155  tilebox.grow(ng_rho);
156 
157 #ifdef AMREX_USE_GPU
158  amrex::ignore_unused(local_rho);
159  // GPU, no tiling: rho_fab points to the full rho array
160  amrex::MultiFab rhoi(*rho, amrex::make_alias, icomp*nc, nc);
161  auto & rho_fab = rhoi.get(pti);
162 #else
163  tb.grow(ng_rho);
164 
165  // CPU, tiling: rho_fab points to local_rho
166  local_rho.resize(tb, nc);
167 
168  // local_rho is set to zero
169  local_rho.setVal(0.0);
170 
171  auto & rho_fab = local_rho;
172 #endif
173 
174  const auto GetPosition = GetParticlePosition(pti, offset);
175 
176  // Indices of the lower bound
177  const amrex::Dim3 lo = lbound(tilebox);
178 
179  ABLASTR_PROFILE_VAR_START(blp_ppc_chd, do_device_synchronize);
180 
181  if (nox == 1){
182  doChargeDepositionShapeN<1>(GetPosition, wp.dataPtr()+offset, ion_lev,
183  rho_fab, np_to_depose.value(), dx, xyzmin, lo, charge,
184  n_rz_azimuthal_modes, cost,
185  load_balance_costs_update_algo);
186  } else if (nox == 2){
187  doChargeDepositionShapeN<2>(GetPosition, wp.dataPtr()+offset, ion_lev,
188  rho_fab, np_to_depose.value(), dx, xyzmin, lo, charge,
189  n_rz_azimuthal_modes, cost,
190  load_balance_costs_update_algo);
191  } else if (nox == 3){
192  doChargeDepositionShapeN<3>(GetPosition, wp.dataPtr()+offset, ion_lev,
193  rho_fab, np_to_depose.value(), dx, xyzmin, lo, charge,
194  n_rz_azimuthal_modes, cost,
195  load_balance_costs_update_algo);
196  }
197  ABLASTR_PROFILE_VAR_STOP(blp_ppc_chd, do_device_synchronize);
198 
199 #ifndef AMREX_USE_GPU
200  // CPU, tiling: atomicAdd local_rho into rho
201  ABLASTR_PROFILE_VAR_START(blp_accumulate, do_device_synchronize);
202  (*rho)[pti].atomicAdd(local_rho, tb, tb, 0, icomp*nc, nc);
203  ABLASTR_PROFILE_VAR_STOP(blp_accumulate, do_device_synchronize);
204 #endif
205 }
206 
207 } // namespace ablastr::particles
208 
209 #endif // ABLASTR_DEPOSIT_CHARGE_H_
Array4< int const > offset
#define AMREX_D_DECL(a, b, c)
#define ABLASTR_PROFILE_VAR_START(vname, sync)
Definition: ProfilerWrapper.H:52
#define ABLASTR_PROFILE_VAR_NS(fname, vname, sync)
Definition: ProfilerWrapper.H:51
#define ABLASTR_PROFILE_VAR_STOP(vname, sync)
Definition: ProfilerWrapper.H:53
#define ABLASTR_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition: TextMsg.H:77
void setVal(T const &x, const Box &bx, int dcomp, int ncomp) noexcept
AMREX_GPU_HOST_DEVICE BoxND & grow(int i) noexcept
void resize(const Box &b, int N=1, Arena *ar=nullptr)
IntVect nGrowVect() const noexcept
IndexType ixType() const noexcept
const FAB & get(const MFIter &mfi) const noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE IntVectND< dim > toIntVect() const noexcept
Definition: DepositCharge.H:22
static void deposit_charge(typename T_PC::ParIterType &pti, typename T_PC::RealVector const &wp, amrex::Real const charge, int const *const ion_lev, amrex::MultiFab *rho, amrex::FArrayBox &local_rho, int const particle_shape, std::array< amrex::Real, 3 > const &dx, std::array< amrex::Real, 3 > const &xyzmin, int const n_rz_azimuthal_modes=0, std::optional< amrex::IntVect > num_rho_deposition_guards=std::nullopt, std::optional< int > depos_lev=std::nullopt, std::optional< amrex::IntVect > rel_ref_ratio=std::nullopt, long const offset=0, std::optional< long > np_to_depose=std::nullopt, int const icomp=0, int const nc=1, amrex::Real *const AMREX_RESTRICT cost=nullptr, long const load_balance_costs_update_algo=0, bool const do_device_synchronize=true)
Definition: DepositCharge.H:55
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
IntVectND< AMREX_SPACEDIM > IntVect
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
int numParticlesOutOfRange(Iterator const &pti, int nGrow)
int nc
Definition: check_interp_points_and_weights.py:130
int nox
Definition: stencil.py:442
int dx
Definition: stencil.py:436
int noy
Definition: stencil.py:443
int noz
Definition: stencil.py:444
Functor that can be used to extract the positions of the macroparticles inside a ParallelFor kernel.
Definition: GetAndSetPosition.H:53