Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_ParticleMesh.H
Go to the documentation of this file.
1#ifndef AMREX_PARTICLEMESH_H_
2#define AMREX_PARTICLEMESH_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_TypeTraits.H>
6#include <AMReX_MultiFab.H>
8#include <type_traits>
9
10namespace amrex {
11
13namespace particle_detail {
14
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
22{
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);
36 } else {
37 return f(p, i, fabarr);
38 }
39}
40
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
48{
49 if constexpr (IsCallable<F, decltype(p), int, decltype(fabarr), decltype(plo), decltype(dxi)>::value) {
50 return f(p, i, fabarr, plo, dxi);
51 } else {
52 return f(p, i, fabarr);
53 }
54}
55}
57
58template <class PC, class MF, class F>
59requires (IsParticleContainer<PC>::value)
60void
61ParticleToMesh (PC const& pc, MF& mf, int lev, F const& f, bool zero_out_input=true)
62{
63 BL_PROFILE("amrex::ParticleToMesh");
64
65 if (zero_out_input) { mf.setVal(0.0); }
66
67 MF* mf_pointer;
68
69 if (pc.OnSameGrids(lev, mf) && zero_out_input)
70 {
71 mf_pointer = &mf;
72 } else {
73 mf_pointer = new MF(pc.ParticleBoxArray(lev),
74 pc.ParticleDistributionMap(lev),
75 mf.nComp(), mf.nGrowVect());
76 mf_pointer->setVal(0.0);
77 }
78
79 const auto plo = pc.Geom(lev).ProbLoArray();
80 const auto dxi = pc.Geom(lev).InvCellSizeArray();
81
82 using ParIter = typename PC::ParConstIterType;
83#ifdef AMREX_USE_GPU
85 {
86 for(ParIter pti(pc, lev); pti.isValid(); ++pti)
87 {
88 const auto& tile = pti.GetParticleTile();
89 const auto np = tile.numParticles();
90 const auto& ptd = tile.getConstParticleTileData();
91
92 auto& fab = (*mf_pointer)[pti];
93 auto fabarr = fab.array();
94
95 AMREX_FOR_1D( np, i,
96 {
97 particle_detail::call_f(f, ptd, i, fabarr, plo, dxi);
98 });
99 }
100 }
101 else
102#endif
103 {
104#ifdef AMREX_USE_OMP
105#pragma omp parallel if (Gpu::notInLaunchRegion())
106#endif
107 {
108 typename MF::FABType::value_type local_fab;
109 for(ParIter pti(pc, lev); pti.isValid(); ++pti)
110 {
111 const auto& tile = pti.GetParticleTile();
112 const auto np = tile.numParticles();
113 const auto& ptd = tile.getConstParticleTileData();
114
115 auto& fab = (*mf_pointer)[pti];
116
117 Box tile_box = pti.tilebox();
118 tile_box.grow(mf_pointer->nGrowVect());
119 local_fab.resize(tile_box,mf_pointer->nComp());
120 local_fab.template setVal<RunOn::Host>(0.0);
121 auto fabarr = local_fab.array();
122
123 AMREX_FOR_1D( np, i,
124 {
125 particle_detail::call_f(f, ptd, i, fabarr, plo, dxi);
126 });
127
128 fab.template atomicAdd<RunOn::Host>(local_fab, tile_box, tile_box,
129 0, 0, mf_pointer->nComp());
130 }
131 }
132 }
133
134 if (mf_pointer != &mf)
135 {
136 mf.ParallelAdd(*mf_pointer, 0, 0, mf_pointer->nComp(),
137 mf_pointer->nGrowVect(), IntVect(0), pc.Geom(lev).periodicity());
138 delete mf_pointer;
139 } else {
140 mf_pointer->SumBoundary(pc.Geom(lev).periodicity());
141 }
142}
143
144template <class PC, class MF, class F>
145requires (IsParticleContainer<PC>::value)
146void
147MeshToParticle (PC& pc, MF const& mf, int lev, F const& f)
148{
149 BL_PROFILE("amrex::MeshToParticle");
150
151 MF* mf_pointer = pc.OnSameGrids(lev, mf) ?
152 const_cast<MF*>(&mf) : new MF(pc.ParticleBoxArray(lev),
153 pc.ParticleDistributionMap(lev),
154 mf.nComp(), mf.nGrowVect());
155
156 if (mf_pointer != &mf) {mf_pointer->ParallelCopy(mf,0,0,mf.nComp(),mf.nGrowVect(),mf.nGrowVect()); }
157
158 const auto plo = pc.Geom(lev).ProbLoArray();
159 const auto dxi = pc.Geom(lev).InvCellSizeArray();
160
161 using ParIter = typename PC::ParIterType;
162#ifdef AMREX_USE_OMP
163#pragma omp parallel if (Gpu::notInLaunchRegion())
164#endif
165 for(ParIter pti(pc, lev); pti.isValid(); ++pti)
166 {
167 auto& tile = pti.GetParticleTile();
168 const auto np = tile.numParticles();
169 const auto& ptd = tile.getParticleTileData();
170
171 const auto& fab = (*mf_pointer)[pti];
172 auto fabarr = fab.array();
173
174 AMREX_FOR_1D( np, i,
175 {
176 particle_detail::call_f(f, ptd, i, fabarr, plo, dxi);
177 });
178 }
179
180 if (mf_pointer != &mf) { delete mf_pointer; }
181}
182
183}
184#endif
#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:649
__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:901
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:172
ParticleTileRef GetParticleTile() const
Definition AMReX_ParIter.H:87
Definition AMReX_ParIter.H:118
bool inLaunchRegion() noexcept
Definition AMReX_GpuControl.H:88
Definition AMReX_Amr.cpp:50
void MeshToParticle(PC &pc, MF const &mf, int lev, F const &f)
Definition AMReX_ParticleMesh.H:147
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:189