1 #ifndef AMREX_AmrParticles_H_
2 #define AMREX_AmrParticles_H_
3 #include <AMReX_Config.H>
13 template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
14 template<
class>
class Allocator,
class CellAssignor>
18 Vector<std::unique_ptr<MultiFab> >& mf_to_be_filled,
19 int lev_min,
int ncomp,
int finest_level,
int ngrow)
const
22 BL_PROFILE(
"ParticleContainer::AssignDensity()");
24 if (rho_index != 0) {
amrex::Abort(
"AssignDensity only works if rho_index = 0"); }
28 BL_ASSERT(ncomp == 1 || ncomp == AMREX_SPACEDIM+1);
30 if (finest_level == -1) {
31 finest_level = this->finestLevel();
33 while (!this->m_gdb->LevelDefined(finest_level)) {
40 mf_to_be_filled.resize(finest_level+1);
41 for (
int lev = lev_min; lev <= finest_level; lev++)
43 auto ng = lev == lev_min ?
IntVect(
AMREX_D_DECL(ngrow,ngrow,ngrow)) : this->m_gdb->refRatio(lev-1);
44 mf_to_be_filled[lev] = std::make_unique<MultiFab>(this->m_gdb->boxArray(lev),
45 this->m_gdb->DistributionMap(lev),
47 mf_to_be_filled[lev]->setVal(0.0);
51 bool all_grids_the_same =
true;
52 for (
int lev = lev_min; lev <= finest_level; lev++)
54 if (!this->OnSameGrids(lev, *mf_to_be_filled[lev]))
56 all_grids_the_same =
false;
62 if (!all_grids_the_same)
65 mf_part.resize(finest_level+1);
66 for (
int lev = lev_min; lev <= finest_level; lev++)
68 auto ng = lev == lev_min ?
IntVect(
AMREX_D_DECL(ngrow,ngrow,ngrow)) : this->m_gdb->refRatio(lev-1);
69 mf_part[lev] = std::make_unique<MultiFab>(this->ParticleBoxArray(lev),
70 this->ParticleDistributionMap(lev),
72 mf_part[lev]->setVal(0.0);
76 auto & mf = (all_grids_the_same) ? mf_to_be_filled : mf_part;
78 if (finest_level == 0)
83 this->AssignCellDensitySingleLevel(rho_index, *mf[0], 0, ncomp);
84 if ( ! all_grids_the_same) {
85 mf_to_be_filled[0]->ParallelCopy(*mf[0],0,0,ncomp,0,0);
91 int lo_bc[] = {BCType::int_dir, BCType::int_dir, BCType::int_dir};
92 int hi_bc[] = {BCType::int_dir, BCType::int_dir, BCType::int_dir};
97 for (
int lev = lev_min; lev <= finest_level; ++lev) {
98 const BoxArray& ba = mf[lev]->boxArray();
100 tmp[lev] = std::make_unique<MultiFab>(ba, dm, ncomp, 0);
101 tmp[lev]->setVal(0.0);
104 for (
int lev = lev_min; lev <= finest_level; ++lev) {
106 this->AssignCellDensitySingleLevel(rho_index, *mf[lev], lev, ncomp, 0);
108 if (lev < finest_level) {
111 rho_index, rho_index, ncomp,
112 this->m_gdb->Geom(lev),
113 this->m_gdb->Geom(lev+1),
114 cphysbc, 0, fphysbc, 0,
115 this->m_gdb->refRatio(lev), &mapper,
124 this->m_gdb->refRatio(lev-1),
125 this->m_gdb->Geom(lev-1),
126 this->m_gdb->Geom(lev));
129 mf[lev]->plus(*tmp[lev], rho_index, ncomp, 0);
132 for (
int lev = finest_level - 1; lev >= lev_min; --lev) {
136 if (!all_grids_the_same) {
137 for (
int lev = lev_min; lev <= finest_level; lev++) {
141 mf_to_be_filled[lev]->ParallelCopy(*mf_part[lev],0,0,ncomp,0,0);
145 int nlevels = finest_level - lev_min + 1;
146 for (
int i = 0; i < nlevels; i++)
148 mf_to_be_filled[i] = std::move(mf_to_be_filled[i+lev_min]);
150 mf_to_be_filled.resize(nlevels);
154 template <class PC, class F, std::enable_if_t<IsParticleContainer<PC>::value,
int> foo = 0>
157 int lev_min,
int lev_max, F&&
f,
158 bool zero_out_input=
true,
bool vol_weight=
true)
162 if (lev_max == -1) { lev_max = pc.finestLevel(); }
163 while (!pc.GetParGDB()->LevelDefined(lev_max)) { lev_max--; }
169 const Real* dx = pc.Geom(0).CellSize();
171 mf[0]->mult(Real(1.0)/vol, 0, mf[0]->
nComp(), mf[0]->nGrow());
178 mf[0]->nGrow(2)), 2);
182 for (
int lev = lev_min; lev <= lev_max; ++lev)
184 mf[lev]->setVal(Real(0.0));
190 for (
int lev = lev_min; lev <= lev_max; ++lev)
192 mf_part[lev].define(pc.ParticleBoxArray(lev),
193 pc.ParticleDistributionMap(lev),
194 mf[lev]->nComp(), ngrow);
195 mf_tmp[lev].define( pc.ParticleBoxArray(lev),
196 pc.ParticleDistributionMap(lev),
197 mf[lev]->nComp(), 0);
198 mf_part[lev].setVal(0.0);
199 mf_tmp[lev].setVal( 0.0);
203 int lo_bc[] = {BCType::int_dir, BCType::int_dir, BCType::int_dir};
204 int hi_bc[] = {BCType::int_dir, BCType::int_dir, BCType::int_dir};
208 for (
int lev = lev_min; lev <= lev_max; ++lev)
212 const Real* dx = pc.Geom(lev).CellSize();
214 mf_part[lev].mult(Real(1.0)/vol, 0, mf_part[lev].
nComp(), mf_part[lev].nGrow());
220 0, 0, mf_part[lev].
nComp(),
221 pc.GetParGDB()->Geom(lev),
222 pc.GetParGDB()->Geom(lev+1),
223 cphysbc, 0, fphysbc, 0,
224 pc.GetParGDB()->refRatio(lev), &mapper,
233 pc.GetParGDB()->refRatio(lev-1),
234 pc.GetParGDB()->Geom(lev-1),
235 pc.GetParGDB()->Geom(lev));
238 mf_part[lev].plus(mf_tmp[lev], 0, mf_part[lev].
nComp(), 0);
241 for (
int lev = lev_max-1; lev >= lev_min; --lev) {
243 pc.GetParGDB()->refRatio(lev));
246 for (
int lev = lev_min; lev <= lev_max; lev++)
248 mf[lev]->ParallelCopy(mf_part[lev],0,0,mf_part[lev].
nComp(),0,0);
252 template <
typename T_ParticleType,
int NArrayReal=0,
int NArrayInt=0,
253 template<
class>
class Allocator=
DefaultAllocator,
class CellAssignor=DefaultAssignor>
289 template <
int T_NStructReal,
int T_NStructInt=0,
int T_NArrayReal=0,
int T_NArrayInt=0,
#define BL_PROFILE(a)
Definition: AMReX_BLProfiler.H:551
#define BL_ASSERT(EX)
Definition: AMReX_BLassert.H:39
#define AMREX_D_TERM(a, b, c)
Definition: AMReX_SPACE.H:129
#define AMREX_D_DECL(a, b, c)
Definition: AMReX_SPACE.H:104
Provide basic functionalities to set up an AMR hierarchy.
Definition: AMReX_AmrCore.H:25
Definition: AMReX_AmrParticles.H:256
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)
Definition: AMReX_AmrParticles.H:272
~AmrParticleContainer_impl() override=default
AmrParticleContainer_impl(AmrCore *amr_core)
Definition: AMReX_AmrParticles.H:267
AmrParticleContainer_impl(AmrParticleContainer_impl &&) noexcept=default
T_ParticleType ParticleType
Definition: AMReX_AmrParticles.H:260
AmrParticleContainer_impl()
Definition: AMReX_AmrParticles.H:262
Definition: AMReX_AmrParticles.H:295
AmrTracerParticleContainer(AmrTracerParticleContainer &&) noexcept=default
~AmrTracerParticleContainer() override=default
AmrTracerParticleContainer(const AmrTracerParticleContainer &)=delete
AmrTracerParticleContainer(AmrCore *amr_core)
Definition: AMReX_AmrParticles.H:298
Definition: AMReX_GpuAllocators.H:114
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:550
Calculates the distribution of FABs to MPI processes.
Definition: AMReX_DistributionMapping.H:41
Piecewise Constant interpolation on cell centered data.
Definition: AMReX_Interpolater.H:451
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:145
void AssignDensity(int rho_index, Vector< std::unique_ptr< MultiFab > > &mf_to_be_filled, int lev_min, int ncomp, int finest_level, int ngrow=2) const
Functions depending the layout of the data. Use with caution.
Definition: AMReX_AmrParticles.H:17
static constexpr int NArrayInt
Number of extra integer components stored in struct-of-array form.
Definition: AMReX_ParticleContainer.H:158
static constexpr int NArrayReal
Number of extra Real components stored in struct-of-array form.
Definition: AMReX_ParticleContainer.H:156
DefaultAssignor CellAssignor
Definition: AMReX_ParticleContainer.H:149
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:27
static int f(amrex::Real t, N_Vector y_data, N_Vector y_rhs, void *user_data)
Definition: AMReX_SundialsIntegrator.H:44
@ max
Definition: AMReX_ParallelReduce.H:17
Definition: AMReX_Amr.cpp:49
int nComp(FabArrayBase const &fa)
amrex::ArenaAllocator< T > DefaultAllocator
Definition: AMReX_GpuAllocators.H:194
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:314
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE const T & max(const T &a, const T &b) noexcept
Definition: AMReX_Algorithm.H:35
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:984
IntVectND< AMREX_SPACEDIM > IntVect
Definition: AMReX_BaseFwd.H:30
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)
Definition: AMReX_AmrParticles.H:156
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:393
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition: AMReX.cpp:225
Definition: AMReX_ParticleUtil.H:432
The struct used to store particles.
Definition: AMReX_Particle.H:295