1#ifndef AMREX_PARTICLEMESH_H_
2#define AMREX_PARTICLEMESH_H_
3#include <AMReX_Config.H>
13namespace particle_detail {
15template <
typename F,
typename T,
typename T_ParticleType,
template<
class,
int,
int>
class PTDType,
int NAR,
int NAI>
17auto call_f (F
const& f,
18 const PTDType<T_ParticleType, NAR, NAI>& p,
19 const int i, Array4<T>
const& fabarr,
20 GpuArray<Real,AMREX_SPACEDIM>
const& plo,
21 GpuArray<Real,AMREX_SPACEDIM>
const& dxi)
noexcept
23 using PTDTypeT = std::remove_const_t<std::remove_reference_t<
decltype(p)>>;
24 if constexpr ( ! T_ParticleType::is_soa_particle &&
25 IsCallable<
F,
typename PTDTypeT::ParticleRefType,
decltype(fabarr),
decltype(plo),
decltype(dxi)>::value) {
26 return f(p.m_aos[i], fabarr, plo, dxi);
27 }
else if constexpr ( ! T_ParticleType::is_soa_particle &&
28 IsCallable<
F,
typename PTDTypeT::ParticleRefType,
decltype(fabarr)>::value) {
29 return f(p.m_aos[i], fabarr);
30 }
else if constexpr (IsCallable<
F,
decltype(p.getSuperParticle(i)),
decltype(fabarr),
decltype(plo),
decltype(dxi)>::value) {
31 return f(p.getSuperParticle(i), fabarr, plo, dxi);
32 }
else if constexpr (IsCallable<
F,
decltype(p.getSuperParticle(i)),
decltype(fabarr)>::value) {
33 return f(p.getSuperParticle(i), fabarr);
34 }
else if constexpr (IsCallable<
F,
decltype(p),
int,
decltype(fabarr),
decltype(plo),
decltype(dxi)>::value) {
35 return f(p, i, fabarr, plo, dxi);
37 return f(p, i, fabarr);
41template <
typename F,
typename T,
template<
class,
class>
class PTDType,
class RType,
class IType>
43auto call_f (F
const& f,
44 const PTDType<RType, IType>& p,
45 const int i, Array4<T>
const& fabarr,
46 GpuArray<Real,AMREX_SPACEDIM>
const& plo,
47 GpuArray<Real,AMREX_SPACEDIM>
const& dxi)
noexcept
49 if constexpr (IsCallable<
F,
decltype(p),
int,
decltype(fabarr),
decltype(plo),
decltype(dxi)>::value) {
50 return f(p, i, fabarr, plo, dxi);
52 return f(p, i, fabarr);
58template <class PC, class MF, class F, std::enable_if_t<IsParticleContainer<PC>::value,
int> foo = 0>
60ParticleToMesh (PC
const& pc, MF& mf,
int lev,
F const& f,
bool zero_out_input=
true)
64 if (zero_out_input) { mf.setVal(0.0); }
68 if (pc.OnSameGrids(lev, mf) && zero_out_input)
72 mf_pointer =
new MF(pc.ParticleBoxArray(lev),
73 pc.ParticleDistributionMap(lev),
74 mf.nComp(), mf.nGrowVect());
75 mf_pointer->setVal(0.0);
78 const auto plo = pc.Geom(lev).ProbLoArray();
79 const auto dxi = pc.Geom(lev).InvCellSizeArray();
81 using ParIter =
typename PC::ParConstIterType;
88 const auto np = tile.numParticles();
89 const auto& ptd = tile.getConstParticleTileData();
91 auto& fab = (*mf_pointer)[pti];
92 auto fabarr = fab.array();
96 particle_detail::call_f(f, ptd, i, fabarr, plo, dxi);
104#pragma omp parallel if (Gpu::notInLaunchRegion())
107 typename MF::FABType::value_type local_fab;
110 const auto& tile = pti.GetParticleTile();
111 const auto np = tile.numParticles();
112 const auto& ptd = tile.getConstParticleTileData();
114 auto& fab = (*mf_pointer)[pti];
116 Box tile_box = pti.tilebox();
117 tile_box.
grow(mf_pointer->nGrowVect());
118 local_fab.
resize(tile_box,mf_pointer->nComp());
119 local_fab.template setVal<RunOn::Host>(0.0);
120 auto fabarr = local_fab.array();
124 particle_detail::call_f(f, ptd, i, fabarr, plo, dxi);
127 fab.template atomicAdd<RunOn::Host>(local_fab, tile_box, tile_box,
128 0, 0, mf_pointer->nComp());
133 if (mf_pointer != &mf)
135 mf.ParallelAdd(*mf_pointer, 0, 0, mf_pointer->nComp(),
136 mf_pointer->nGrowVect(),
IntVect(0), pc.Geom(lev).periodicity());
139 mf_pointer->SumBoundary(pc.Geom(lev).periodicity());
143template <class PC, class MF, class F, std::enable_if_t<IsParticleContainer<PC>::value,
int> foo = 0>
149 MF* mf_pointer = pc.OnSameGrids(lev, mf) ?
150 const_cast<MF*
>(&mf) :
new MF(pc.ParticleBoxArray(lev),
151 pc.ParticleDistributionMap(lev),
152 mf.nComp(), mf.nGrowVect());
154 if (mf_pointer != &mf) {mf_pointer->ParallelCopy(mf,0,0,mf.nComp(),mf.nGrowVect(),mf.nGrowVect()); }
156 const auto plo = pc.Geom(lev).ProbLoArray();
157 const auto dxi = pc.Geom(lev).InvCellSizeArray();
159 using ParIter =
typename PC::ParIterType;
161#pragma omp parallel if (Gpu::notInLaunchRegion())
166 const auto np = tile.numParticles();
167 const auto& ptd = tile.getParticleTileData();
169 const auto& fab = (*mf_pointer)[pti];
170 auto fabarr = fab.array();
174 particle_detail::call_f(f, ptd, i, fabarr, plo, dxi);
178 if (mf_pointer != &mf) {
delete mf_pointer; }
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_FORCE_INLINE
Definition AMReX_Extension.H:119
#define AMREX_FOR_1D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:97
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
__host__ __device__ BoxND & grow(int i) noexcept
Definition AMReX_Box.H:641
__host__ __device__ BoxND< new_dim > resize() const noexcept
Return a new BoxND of size new_dim by either shrinking or expanding this BoxND.
Definition AMReX_Box.H:893
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:169
ParticleTileRef GetParticleTile() const
Definition AMReX_ParIter.H:87
Definition AMReX_ParIter.H:118
bool inLaunchRegion() noexcept
Definition AMReX_GpuControl.H:92
Definition AMReX_Amr.cpp:49
void MeshToParticle(PC &pc, MF const &mf, int lev, F const &f)
Definition AMReX_ParticleMesh.H:145
void ParticleToMesh(PC const &pc, const Vector< MultiFab * > &mf, int lev_min, int lev_max, F &&f, bool zero_out_input=true, bool vol_weight=true)
Deposit particles onto a hierarchy of MultiFabs.
Definition AMReX_AmrParticles.H:188