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