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