Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_AmrParticles.H
Go to the documentation of this file.
1#ifndef AMREX_AmrParticles_H_
2#define AMREX_AmrParticles_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_Particles.H>
7#include <AMReX_AmrParGDB.H>
10
16namespace amrex {
17
18// This is already documented in AMReX_ParticleContainer.H. So we need Doxygen to
19// ignore it to avoid warnings.
21
32template <typename ParticleType, int NArrayReal, int NArrayInt,
33 template<class> class Allocator, class CellAssignor>
34void
35ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
36::AssignDensity (int rho_index,
37 Vector<std::unique_ptr<MultiFab> >& mf_to_be_filled,
38 int lev_min, int ncomp, int finest_level, int ngrow) const
39{
40
41 BL_PROFILE("ParticleContainer::AssignDensity()");
42
43 if (rho_index != 0) { amrex::Abort("AssignDensity only works if rho_index = 0"); }
44
45 BL_ASSERT(NStructReal >= 1);
46 BL_ASSERT(NStructReal >= ncomp);
47 BL_ASSERT(ncomp == 1 || ncomp == AMREX_SPACEDIM+1);
48
49 if (finest_level == -1) {
50 finest_level = this->finestLevel();
51 }
52 while (!this->m_gdb->LevelDefined(finest_level)) {
53 finest_level--;
54 }
55
56 ngrow = std::max(ngrow, 2);
57
58 // Create the space for mf_to_be_filled, regardless of whether we'll need a temporary mf
59 mf_to_be_filled.resize(finest_level+1);
60 for (int lev = lev_min; lev <= finest_level; lev++)
61 {
62 auto ng = lev == lev_min ? IntVect(AMREX_D_DECL(ngrow,ngrow,ngrow)) : this->m_gdb->refRatio(lev-1);
63 mf_to_be_filled[lev] = std::make_unique<MultiFab>(this->m_gdb->boxArray(lev),
64 this->m_gdb->DistributionMap(lev),
65 ncomp, ng);
66 mf_to_be_filled[lev]->setVal(0.0);
67 }
68
69 // Test whether the grid structure of the boxArray is the same as the ParticleBoxArray at all levels
70 bool all_grids_the_same = true;
71 for (int lev = lev_min; lev <= finest_level; lev++)
72 {
73 if (!this->OnSameGrids(lev, *mf_to_be_filled[lev]))
74 {
75 all_grids_the_same = false;
76 break;
77 }
78 }
79
81 if (!all_grids_the_same)
82 {
83 // Create the space for the temporary, mf_part
84 mf_part.resize(finest_level+1);
85 for (int lev = lev_min; lev <= finest_level; lev++)
86 {
87 auto ng = lev == lev_min ? IntVect(AMREX_D_DECL(ngrow,ngrow,ngrow)) : this->m_gdb->refRatio(lev-1);
88 mf_part[lev] = std::make_unique<MultiFab>(this->ParticleBoxArray(lev),
89 this->ParticleDistributionMap(lev),
90 ncomp, ng);
91 mf_part[lev]->setVal(0.0);
92 }
93 }
94
95 auto & mf = (all_grids_the_same) ? mf_to_be_filled : mf_part;
96
97 if (finest_level == 0)
98 {
99 //
100 // Just use the far simpler single-level version.
101 //
102 this->AssignCellDensitySingleLevel(rho_index, *mf[0], 0, ncomp);
103 if ( ! all_grids_the_same) {
104 mf_to_be_filled[0]->ParallelCopy(*mf[0],0,0,ncomp,0,0);
105 }
106 return;
107 }
108
109 // configure this to do a no-op at the physical boundaries.
112 Vector<BCRec> bcs(ncomp, BCRec(lo_bc, hi_bc));
113 PCInterp mapper;
114
115 Vector<std::unique_ptr<MultiFab> > tmp(finest_level+1);
116 for (int lev = lev_min; lev <= finest_level; ++lev) {
117 const BoxArray& ba = mf[lev]->boxArray();
118 const DistributionMapping& dm = mf[lev]->DistributionMap();
119 tmp[lev] = std::make_unique<MultiFab>(ba, dm, ncomp, 0);
120 tmp[lev]->setVal(0.0);
121 }
122
123 for (int lev = lev_min; lev <= finest_level; ++lev) {
124
125 this->AssignCellDensitySingleLevel(rho_index, *mf[lev], lev, ncomp, 0);
126
127 if (lev < finest_level) {
128 PhysBCFunctNoOp cphysbc, fphysbc;
129 amrex::InterpFromCoarseLevel(*tmp[lev+1], 0.0, *mf[lev],
130 rho_index, rho_index, ncomp,
131 this->m_gdb->Geom(lev),
132 this->m_gdb->Geom(lev+1),
133 cphysbc, 0, fphysbc, 0,
134 this->m_gdb->refRatio(lev), &mapper,
135 bcs, 0);
136 }
137
138 if (lev > lev_min) {
139 // Note - this will double count the mass on the coarse level in
140 // regions covered by the fine level, but this will be corrected
141 // below in the call to average_down.
142 amrex::sum_fine_to_coarse(*mf[lev], *mf[lev-1], rho_index, ncomp,
143 this->m_gdb->refRatio(lev-1),
144 this->m_gdb->Geom(lev-1),
145 this->m_gdb->Geom(lev));
146 }
147
148 mf[lev]->plus(*tmp[lev], rho_index, ncomp, 0);
149 }
150
151 for (int lev = finest_level - 1; lev >= lev_min; --lev) {
152 amrex::average_down(*mf[lev+1], *mf[lev], rho_index, ncomp, this->m_gdb->refRatio(lev));
153 }
154
155 if (!all_grids_the_same) {
156 for (int lev = lev_min; lev <= finest_level; lev++) {
157 //
158 // We haven't interpolated the ghost cells so we can't copy them
159 //
160 mf_to_be_filled[lev]->ParallelCopy(*mf_part[lev],0,0,ncomp,0,0);
161 }
162 }
163 if (lev_min > 0) {
164 int nlevels = finest_level - lev_min + 1;
165 for (int i = 0; i < nlevels; i++)
166 {
167 mf_to_be_filled[i] = std::move(mf_to_be_filled[i+lev_min]);
168 }
169 mf_to_be_filled.resize(nlevels);
170 }
171}
172
186template <class PC, class F, std::enable_if_t<IsParticleContainer<PC>::value, int> foo = 0>
187void
188ParticleToMesh (PC const& pc, const Vector<MultiFab*>& mf,
189 int lev_min, int lev_max, F&& f,
190 bool zero_out_input=true, bool vol_weight=true)
191{
192 BL_PROFILE("amrex::ParticleToMesh");
193
194 if (lev_max == -1) { lev_max = pc.finestLevel(); }
195 while (!pc.GetParGDB()->LevelDefined(lev_max)) { lev_max--; }
196
197 if (lev_max == 0)
198 {
199 ParticleToMesh(pc, *mf[0], 0, std::forward<F>(f), zero_out_input);
200 if (vol_weight) {
201 const Real* dx = pc.Geom(0).CellSize();
202 const Real vol = AMREX_D_TERM(dx[0], *dx[1], *dx[2]);
203 mf[0]->mult(Real(1.0)/vol, 0, mf[0]->nComp(), mf[0]->nGrow());
204 }
205 return;
206 }
207
208 int ngrow = amrex::max(AMREX_D_DECL(mf[0]->nGrow(0),
209 mf[0]->nGrow(1),
210 mf[0]->nGrow(2)), 2);
211
212 if (zero_out_input)
213 {
214 for (int lev = lev_min; lev <= lev_max; ++lev)
215 {
216 mf[lev]->setVal(Real(0.0));
217 }
218 }
219
220 Vector<MultiFab> mf_part(lev_max+1);
221 Vector<MultiFab> mf_tmp( lev_max+1);
222 for (int lev = lev_min; lev <= lev_max; ++lev)
223 {
224 mf_part[lev].define(pc.ParticleBoxArray(lev),
225 pc.ParticleDistributionMap(lev),
226 mf[lev]->nComp(), ngrow);
227 mf_tmp[lev].define( pc.ParticleBoxArray(lev),
228 pc.ParticleDistributionMap(lev),
229 mf[lev]->nComp(), 0);
230 mf_part[lev].setVal(0.0);
231 mf_tmp[lev].setVal( 0.0);
232 }
233
234 // configure this to do a no-op at the physical boundaries.
237 Vector<BCRec> bcs(mf_part[lev_min].nComp(), BCRec(lo_bc, hi_bc));
238 PCInterp mapper;
239
240 for (int lev = lev_min; lev <= lev_max; ++lev)
241 {
242 ParticleToMesh(pc, mf_part[lev], lev, f, true); // always zero out input on temp level
243 if (vol_weight) {
244 const Real* dx = pc.Geom(lev).CellSize();
245 const Real vol = AMREX_D_TERM(dx[0], *dx[1], *dx[2]);
246 mf_part[lev].mult(Real(1.0)/vol, 0, mf_part[lev].nComp(), mf_part[lev].nGrow());
247 }
248
249 if (lev < lev_max) {
250 PhysBCFunctNoOp cphysbc, fphysbc;
251 amrex::InterpFromCoarseLevel(mf_tmp[lev+1], 0.0, mf_part[lev],
252 0, 0, mf_part[lev].nComp(),
253 pc.GetParGDB()->Geom(lev),
254 pc.GetParGDB()->Geom(lev+1),
255 cphysbc, 0, fphysbc, 0,
256 pc.GetParGDB()->refRatio(lev), &mapper,
257 bcs, 0);
258 }
259
260 if (lev > lev_min) {
261 // Note - this will double count the mass on the coarse level in
262 // regions covered by the fine level, but this will be corrected
263 // below in the call to average_down.
264 amrex::sum_fine_to_coarse(mf_part[lev], mf_part[lev-1], 0, mf_part[lev].nComp(),
265 pc.GetParGDB()->refRatio(lev-1),
266 pc.GetParGDB()->Geom(lev-1),
267 pc.GetParGDB()->Geom(lev));
268 }
269
270 mf_part[lev].plus(mf_tmp[lev], 0, mf_part[lev].nComp(), 0);
271 }
272
273 for (int lev = lev_max-1; lev >= lev_min; --lev) {
274 amrex::average_down(mf_part[lev+1], mf_part[lev], 0, mf_part[lev].nComp(),
275 pc.GetParGDB()->refRatio(lev));
276 }
277
278 for (int lev = lev_min; lev <= lev_max; lev++)
279 {
280 if (zero_out_input) {
281 mf[lev]->ParallelCopy(mf_part[lev],0,0,mf_part[lev].nComp(),0,0);
282 } else { // do not zero out input
283 mf[lev]->ParallelAdd(mf_part[lev],0,0,mf_part[lev].nComp(),0,0);
284 }
285 }
286}
287
288template <typename T_ParticleType, int NArrayReal=0, int NArrayInt=0,
289 template<class> class Allocator=DefaultAllocator, class CellAssignor=DefaultAssignor>
293class AmrParticleContainer_impl // NOLINT(cppcoreguidelines-virtual-class-destructor)
294 : public ParticleContainer_impl<T_ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
295{
296
297public:
298
299 using ParticleType = T_ParticleType;
300
306
312
322 const Vector<DistributionMapping> & dmap,
323 const Vector<BoxArray> & ba,
324 const Vector<int> & rr)
325 : ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>(geom, dmap, ba, rr)
326 {
327 }
328
329 ~AmrParticleContainer_impl () override = default;
330
333
335 AmrParticleContainer_impl& operator= ( AmrParticleContainer_impl && ) noexcept = default;
336};
337
339template <int T_NStructReal, int T_NStructInt=0, int T_NArrayReal=0, int T_NArrayInt=0,
340 template<class> class Allocator=DefaultAllocator, class CellAssignor=DefaultAssignor>
341using AmrParticleContainer = AmrParticleContainer_impl<Particle<T_NStructReal, T_NStructInt>, T_NArrayReal, T_NArrayInt, Allocator, CellAssignor>;
342
365
366}
367
368#endif
Particle Geometry, DistributionMapping and BoxArray helper shared by AmrCore-based apps.
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
High-level FillPatch helpers for AMR coarse-to-fine synchronization.
Collection of spatial interpolaters used by FillPatch and flux-register logic.
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
Provide basic functionalities to set up an AMR hierarchy.
Definition AMReX_AmrCore.H:31
Particle container that automatically pulls grid metadata from AmrCore when available.
Definition AMReX_AmrParticles.H:295
AmrParticleContainer_impl(const AmrParticleContainer_impl &)=delete
AmrParticleContainer_impl & operator=(const AmrParticleContainer_impl &)=delete
AmrParticleContainer_impl(const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< int > &rr)
Construct from explicit Geometry, DistributionMapping and BoxArray.
Definition AMReX_AmrParticles.H:321
~AmrParticleContainer_impl() override=default
AmrParticleContainer_impl(AmrCore *amr_core)
Construct from amr_core.
Definition AMReX_AmrParticles.H:308
AmrParticleContainer_impl(AmrParticleContainer_impl &&) noexcept=default
T_ParticleType ParticleType
Definition AMReX_AmrParticles.H:299
AmrParticleContainer_impl()
Construct an empty container (call Define later).
Definition AMReX_AmrParticles.H:302
Simple tracer container that wires TracerParticleContainer to an AmrCore.
Definition AMReX_AmrParticles.H:348
AmrTracerParticleContainer(AmrTracerParticleContainer &&) noexcept=default
~AmrTracerParticleContainer() override=default
AmrTracerParticleContainer(const AmrTracerParticleContainer &)=delete
AmrTracerParticleContainer(AmrCore *amr_core)
Construct from amr_core.
Definition AMReX_AmrParticles.H:352
Definition AMReX_GpuAllocators.H:115
Boundary Condition Records. Necessary information and functions for computing boundary conditions.
Definition AMReX_BCRec.H:17
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
Piecewise Constant interpolation on cell centered data.
Definition AMReX_Interpolater.H:418
const ParGDBBase * GetParGDB() const
Get the ParGDB object used to define this container (const version)
Definition AMReX_ParticleContainerBase.H:225
A distributed container for Particles sorted onto the levels, grids, and tiles of a block-structured ...
Definition AMReX_ParticleContainer.H:148
static constexpr int NArrayInt
Number of extra integer components stored in struct-of-array form.
Definition AMReX_ParticleContainer.H:161
static constexpr int NArrayReal
Number of extra Real components stored in struct-of-array form.
Definition AMReX_ParticleContainer.H:159
T_CellAssignor CellAssignor
Definition AMReX_ParticleContainer.H:152
Definition AMReX_PhysBCFunct.H:120
Definition AMReX_TracerParticles.H:11
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
@ int_dir
Definition AMReX_BC_TYPES.H:70
Definition AMReX_Amr.cpp:49
int nComp(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2851
amrex::ArenaAllocator< T > DefaultAllocator
Definition AMReX_GpuAllocators.H:199
void average_down(const MultiFab &S_fine, MultiFab &S_crse, const Geometry &fgeom, const Geometry &cgeom, int scomp, int ncomp, int rr)
Definition AMReX_MultiFabUtil.cpp:336
std::enable_if_t< IsFabArray< MF >::value > InterpFromCoarseLevel(MF &mf, Real time, const MF &cmf, int scomp, int dcomp, int ncomp, const Geometry &cgeom, const Geometry &fgeom, BC &cbc, int cbccomp, BC &fbc, int fbccomp, const IntVect &ratio, Interp *mapper, const Vector< BCRec > &bcs, int bcscomp, const PreInterpHook &pre_interp={}, const PostInterpHook &post_interp={})
Fill with interpolation of coarse level data.
Definition AMReX_FillPatchUtil_I.H:1005
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
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
__host__ __device__ constexpr const T & max(const T &a, const T &b) noexcept
Definition AMReX_Algorithm.H:44
void sum_fine_to_coarse(const MultiFab &S_fine, MultiFab &S_crse, int scomp, int ncomp, const IntVect &ratio, const Geometry &cgeom, const Geometry &)
Definition AMReX_MultiFabUtil.cpp:415
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:240
Definition AMReX_ParticleUtil.H:394
The struct used to store particles.
Definition AMReX_Particle.H:404