1#ifndef AMREX_EB2_LEVEL_H_
2#define AMREX_EB2_LEVEL_H_
3#include <AMReX_Config.H>
14#if (AMREX_SPACEDIM == 3)
18#include <AMReX_EB2_C.H>
29#include <unordered_map>
39template <
typename G>
class GShopLevel;
134 void write_to_chkpt_file (
const std::string& fname,
bool extend_domain_face,
int max_grid_size)
const;
144 void setShift (
int direction,
int ncells);
168#if (AMREX_SPACEDIM == 3)
193 void setRegularLevel () {
209 int ngrow,
bool extend_domain_face,
int num_crse_opt);
219 int max_grid_size,
int ngrow,
bool extend_domain_face,
int num_crse_opt);
221#if (AMREX_SPACEDIM == 3)
222 template <
typename GS = G, std::enable_if_t<std::is_same_v<GS,STLtools>,
int> FOO = 0>
227 int max_grid_size,
int ngrow,
bool extend_domain_face,
int num_crse_opt);
245 void prepare_grids (G
const& gshop,
Geometry const& geom,
int max_grid_size,
int ngrow,
246 bool extend_domain_face,
int num_crse_opt);
257 int max_grid_size,
int ngrow,
bool extend_domain_face,
int num_crse_opt)
260 if (std::is_same_v<typename G::FunctionType, AllRegularIF>) {
266 define_fine(gshop, geom, max_grid_size, ngrow, extend_domain_face, num_crse_opt);
272 int max_grid_size,
int ngrow,
bool extend_domain_face,
int num_crse_opt)
275 amrex::Print() <<
"AMReX WARNING: extend_domain_face=false is not recommended!\n";
278 BL_PROFILE(
"EB2::GShopLevel()-prepare_grids");
281 m_ngrow =
IntVect{
static_cast<int>(std::ceil(ngrow/16.)) * 16};
284 Box domain_grown = domain;
285 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
289 m_ngrow[idim] = std::min(m_ngrow[idim], domain_grown.length(idim));
292 domain_grown.grow(m_ngrow);
293 m_bounding_box = (extend_domain_face) ? domain : domain_grown;
294 m_bounding_box.surroundingNodes();
297 BoxList covered_boxes;
302 num_crse_opt = std::max(0,std::min(8,num_crse_opt));
303 for (
int clev = num_crse_opt; clev >= 0; --clev) {
305 if (domain.coarsenable(crse_ratio)) {
309 if (cut_boxes.isEmpty()) {
310 covered_boxes.clear();
311 test_boxes = BoxList(crse_geom.Domain());
312 test_boxes.maxSize(max_grid_size);
314 test_boxes.swap(cut_boxes);
315 test_boxes.coarsen(crse_ratio);
316 test_boxes.maxSize(max_grid_size);
319 const Long nboxes = test_boxes.size();
320 const auto& boxes = test_boxes.data();
321 for (
Long i = iproc; i < nboxes; i += nprocs) {
322 const Box& vbx = boxes[i];
324 auto box_type = gshop.getBoxType(gbx&crse_bounding_box,crse_geom,
RunOn::Gpu);
325 if (box_type == gshop.allcovered) {
327 }
else if (box_type == gshop.mixedcells) {
339 auto grow_at_domain_boundary = [&] (BoxList& bl)
341 const IntVect& domlo = domain.smallEnd();
342 const IntVect& domhi = domain.bigEnd();
344 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
345 if (m_ngrow[idim] != 0) {
346 if (b.smallEnd(idim) == domlo[idim]) {
347 b.growLo(idim,m_ngrow[idim]);
349 if (b.bigEnd(idim) == domhi[idim]) {
350 b.growHi(idim,m_ngrow[idim]);
356 grow_at_domain_boundary(covered_boxes);
357 grow_at_domain_boundary(cut_boxes);
360 if ( cut_boxes.isEmpty() &&
361 !covered_boxes.isEmpty())
363 amrex::Abort(
"AMReX_EB2_Level.H: Domain is completely covered");
366 if (!covered_boxes.isEmpty()) {
367 if (num_crse_opt > 2) {
368 covered_boxes.maxSize(max_grid_size*4);
370 m_covered_grids = BoxArray(std::move(covered_boxes));
373 if (cut_boxes.isEmpty()) {
374 m_grids = BoxArray();
375 m_dmap = DistributionMapping();
379 m_grids = BoxArray(std::move(cut_boxes));
380 m_dmap = DistributionMapping(m_grids);
387 int max_grid_size,
int ngrow,
bool extend_domain_face,
int num_crse_opt)
391#ifdef AMREX_USE_FLOAT
392 Real small_volfrac = 1.e-5_rt;
394 Real small_volfrac = 1.e-14;
396 bool cover_multiple_cuts =
false;
401 pp.
queryAdd(
"cover_multiple_cuts", cover_multiple_cuts);
404 maxiter = std::min(100000, maxiter);
406 prepare_grids(gshop, geom, max_grid_size, ngrow, extend_domain_face, num_crse_opt);
408 if (m_allregular) {
return; }
410 m_mgf.define(m_grids, m_dmap);
413 mf_info.
SetTag(
"EB2::Level");
414 m_cellflag.define(m_grids, m_dmap, 1, ng, mf_info);
415 m_volfrac.define(m_grids, m_dmap, 1, ng, mf_info);
416 m_centroid.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
417 m_bndryarea.define(m_grids, m_dmap, 1, ng, mf_info);
418 m_bndrycent.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
419 m_bndrynorm.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
420 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
422 m_dmap, 1, ng, mf_info);
424 m_dmap, AMREX_SPACEDIM-1, ng, mf_info);
425 IntVect edge_type{1}; edge_type[idim] = 0;
426 m_edgecent[idim].define(
amrex::convert(m_grids, edge_type), m_dmap, 1, ng, mf_info);
432 auto bounding_box = m_bounding_box;
433 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
434 if (!extend_domain_face || geom.
isPeriodic(idim)) {
446 for (; iter < maxiter; ++iter)
451#pragma omp parallel if (Gpu::notInLaunchRegion()) reduction(+:nsmallcells,nmulticuts)
454#if (AMREX_SPACEDIM == 3)
460 auto& gfab = m_mgf[mfi];
461 const Box& vbx = gfab.validbox();
463 auto&
levelset = gfab.getLevelSet();
465 gshop.fillFab(
levelset, geom, gshop_run_on, bounding_box);
492 auto& facetype = gfab.getFaceType();
500#if (AMREX_SPACEDIM == 3)
501 auto& edgetype = gfab.getEdgeType();
514 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
515 edgetype[idim].prefetchToHost();
516 m_edgecent[idim][mfi].prefetchToHost();
521 gshop.getIntercept({xip,yip,zip}, {xdg,ydg,zdg}, clst,
522 geom, gshop_run_on, bounding_box);
526 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
527 edgetype[idim].prefetchToDevice();
528 m_edgecent[idim][mfi].prefetchToDevice();
533 gshop.updateIntercept({xip,yip,zip}, {xdg,ydg,zdg}, clst, geom);
536 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
537 const Box& b = facetype[idim].box();
544 nmc = build_faces(vbx, cfg, ftx, fty, ftz, xdg, ydg, zdg, lst,
546 xm2, ym2, zm2, dx, problo, cover_multiple_cuts);
548 cellflagtmp.
resize(m_cellflag[mfi].box());
552 build_cells(vbx, cfg, ftx, fty, ftz,
apx,
apy,
apz,
553 fcx,
fcy,
fcz, xm2, ym2, zm2, dx, vfr, ctr,
554 bar, bct, bnm, cfgtmp, lst,
555 small_volfrac, geom, extend_domain_face, cover_multiple_cuts,
563#elif (AMREX_SPACEDIM == 2)
571 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
572 facetype[idim].prefetchToHost();
573 m_edgecent[idim][mfi].prefetchToHost();
579 gshop.getIntercept({xip,yip},
580 {facetype[1].const_array(), facetype[0].const_array()},
581 clst, geom, gshop_run_on, bounding_box);
585 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
586 facetype[idim].prefetchToDevice();
587 m_edgecent[idim][mfi].prefetchToDevice();
592 gshop.updateIntercept({xip,yip},
593 {facetype[1].const_array(), facetype[0].const_array()},
597 nmc = build_faces(vbx, cfg, ftx, fty, lst, xip, yip,
apx,
apy,
fcx,
fcy,
598 dx, problo, cover_multiple_cuts, nsm);
600 build_cells(vbx, cfg, ftx, fty,
apx,
apy, dx, vfr, ctr,
601 bar, bct, bnm, lst, small_volfrac, geom, extend_domain_face,
610 if (nsmallcells == 0 && nmulticuts == 0) {
613 auto ls = m_mgf.getLevelSet();
618 amrex::Print() <<
"AMReX EB: Iter. " << iter+1 <<
" fixed " << nsmallcells
619 <<
" small cells" <<
'\n';
622 amrex::Print() <<
"AMReX EB: Iter. " << iter+1 <<
" fixed " << nmulticuts
623 <<
" multicuts" <<
'\n';
632#pragma omp parallel if (Gpu::notInLaunchRegion())
636 auto& gfab = m_mgf[mfi];
637 auto const&
levelset = gfab.getLevelSet();
643#
if (AMREX_SPACEDIM == 3)
644 auto const& edgetype = gfab.getEdgeType();
649 intercept_to_edge_centroid(xip, yip, zip, xdg, ydg, zdg, clst, dx, problo);
651#elif (AMREX_SPACEDIM == 2)
652 auto& facetype = gfab.getFaceType();
656 intercept_to_edge_centroid(xip, yip, fty, ftx, clst, dx, problo);
660 m_levelset = m_mgf.getLevelSet();
665#if (AMREX_SPACEDIM == 3)
667template <
typename GS, std::enable_if_t<std::is_same_v<GS,STLtools>,
int> FOO>
670 int max_grid_size,
int ngrow,
bool extend_domain_face,
int num_crse_opt)
674 prepare_grids(gshop, geom, max_grid_size, ngrow, extend_domain_face, num_crse_opt);
676 if (m_allregular) {
return; }
681 gshop.fillSignedDistance(m_sdf, m_sdf.nGrowVect(), geom);
684 m_marching_cubes[mfi.index()] = std::make_unique<MC::MCFab>();
689#if defined(AMREX_USE_OMP) && !defined(AMREX_USE_GPU)
719 const int coarse_ratio = 2;
720 const int min_width = 8;
721 bool coarsenable = fine_grids.
coarsenable(coarse_ratio, min_width)
722 && (fine_covered_grids.
empty() || fine_covered_grids.
coarsenable(coarse_ratio));
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition AMReX_BLassert.H:49
amrex::ParmParse pp
Input file parser instance for the given namespace.
Definition AMReX_HypreIJIface.cpp:15
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
void resize(const Box &b, int N=1, Arena *ar=nullptr)
This function resizes a BaseFab so it covers the Box b with N components.
Definition AMReX_BaseFab.H:1629
Array4< T const > array() const noexcept
Definition AMReX_BaseFab.H:382
Elixir elixir() noexcept
Definition AMReX_BaseFab.H:1671
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
bool coarsenable(int refinement_ratio, int min_width=1) const
Coarsen each Box in the BoxArray to the specified ratio.
Definition AMReX_BoxArray.cpp:615
bool empty() const noexcept
Return whether the BoxArray is empty.
Definition AMReX_BoxArray.H:621
__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
GpuArray< Real, 3 > CellSizeArray() const noexcept
Definition AMReX_CoordSys.H:76
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
static constexpr int ng
Definition AMReX_EB2_MultiGFab.H:26
Definition AMReX_EB2_Level.H:203
GShopLevel(IndexSpace const *is, const Geometry &geom)
Construct an empty level (usually regular) bound to geom.
Definition AMReX_EB2_Level.H:251
void define_fine(G const &gshop, const Geometry &geom, int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
Define data for the finest level using the supplied GeometryShop.
Definition AMReX_EB2_Level.H:386
GShopLevel(IndexSpace const *is, int ilev, int max_grid_size, int ngrow, const Geometry &geom, GShopLevel< G > &fineLevel)
Construct a coarse level by coarsening fineLevel.
Definition AMReX_EB2_Level.H:704
void define_fine_mvmc(GS const &gshop, const Geometry &geom, int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
Define a fine level for multi-valued multi-cut STL geometries (3-D only).
Definition AMReX_EB2_Level.H:669
GShopLevel(IndexSpace const *is, G const &gshop, const Geometry &geom, int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
Build a level directly from a GeometryShop object.
Definition AMReX_EB2_Level.H:256
static GShopLevel< G > makeAllRegular(IndexSpace const *is, const Geometry &geom)
Build a regular (no EB) level.
Definition AMReX_EB2_Level.H:234
Definition AMReX_EB2.H:28
Definition AMReX_EB2_Level.H:42
MultiFab m_bndryarea
Definition AMReX_EB2_Level.H:176
bool isOK() const noexcept
Definition AMReX_EB2_Level.H:48
void fillFaceCent(Array< MultiCutFab *, 3 > const &facecent, const Geometry &geom) const
Populate face centroids for each direction.
Definition AMReX_EB2_Level.cpp:794
MultiFab m_volfrac
Definition AMReX_EB2_Level.H:174
MultiGFab m_mgf
Definition AMReX_EB2_Level.H:171
bool m_allregular
Definition AMReX_EB2_Level.H:183
IntVect m_shift
Definition AMReX_EB2_Level.H:186
iMultiFab m_cutcellmask
Definition AMReX_EB2_Level.H:182
Array< MultiFab, 3 > m_facecent
Definition AMReX_EB2_Level.H:180
Array< MultiFab, 3 > m_edgecent
Definition AMReX_EB2_Level.H:181
void prepareForCoarsening(const Level &rhs, int max_grid_size, IntVect const &ngrow)
Prepare for coarsening by copying from rhs, respecting max grid size and grow vector.
Definition AMReX_EB2_Level.cpp:10
void fillCutCellMask(iMultiFab &cutcellmask, const Geometry &geom) const
Populate the cut-cell mask into cutcellmask.
Definition AMReX_EB2_Level.cpp:946
MultiFab m_centroid
Definition AMReX_EB2_Level.H:175
FabArray< EBCellFlagFab > m_cellflag
Definition AMReX_EB2_Level.H:173
const DistributionMapping & DistributionMap() const noexcept
Definition AMReX_EB2_Level.H:110
BoxArray m_covered_grids
Definition AMReX_EB2_Level.H:165
void fillLevelSet(MultiFab &levelset, const Geometry &geom) const
Write the implicit function values into levelset.
Definition AMReX_EB2_Level.cpp:907
Geometry m_geom
Definition AMReX_EB2_Level.H:162
Array< MultiFab, 3 > m_areafrac
Definition AMReX_EB2_Level.H:179
IntVect m_ngrow
Definition AMReX_EB2_Level.H:163
IndexSpace const * m_parent
Definition AMReX_EB2_Level.H:187
BoxArray m_grids
Definition AMReX_EB2_Level.H:164
const Geometry & Geom() const noexcept
Definition AMReX_EB2_Level.H:120
std::map< int, std::unique_ptr< MC::MCFab > > m_marching_cubes
Definition AMReX_EB2_Level.H:169
bool hasEBInfo() const noexcept
Definition AMReX_EB2_Level.H:137
void write_to_chkpt_file(const std::string &fname, bool extend_domain_face, int max_grid_size) const
Write this level’s EB data to a checkpoint file fname.
Definition AMReX_EB2_Level.cpp:956
void fillEBCellFlag(FabArray< EBCellFlagFab > &cellflag, const Geometry &geom) const
Fill cellflag with EBCellFlag data using geometry geom.
Definition AMReX_EB2_Level.cpp:429
void fillEdgeCent(Array< MultiCutFab *, 3 > const &edgecent, const Geometry &geom) const
Populate edge centroids for each direction.
Definition AMReX_EB2_Level.cpp:835
int coarsenFromFine(Level &fineLevel, bool fill_boundary)
Coarsen EB data from fineLevel into this level.
Definition AMReX_EB2_Level.cpp:84
MultiFab m_bndrynorm
Definition AMReX_EB2_Level.H:178
bool m_has_eb_info
Definition AMReX_EB2_Level.H:185
void fillAreaFrac(Array< MultiCutFab *, 3 > const &areafrac, const Geometry &geom) const
Populate face area fractions for each direction.
Definition AMReX_EB2_Level.cpp:640
MultiFab m_bndrycent
Definition AMReX_EB2_Level.H:177
MultiFab m_levelset
Definition AMReX_EB2_Level.H:172
void setShift(int direction, int ncells)
Shift this level by ncells cells along direction direction (for periodic tiling).
Definition AMReX_EB2_Level.cpp:1047
bool isAllRegular() const noexcept
Definition AMReX_EB2_Level.H:46
void fillBndryCent(MultiCutFab &bndrycent, const Geometry &geom) const
Populate boundary centroids into bndrycent.
Definition AMReX_EB2_Level.cpp:592
void fillBndryNorm(MultiCutFab &bndrynorm, const Geometry &geom) const
Populate boundary normals into bndrynorm.
Definition AMReX_EB2_Level.cpp:616
MultiFab m_sdf
Definition AMReX_EB2_Level.H:167
IntVect const & nGrowVect() const noexcept
Definition AMReX_EB2_Level.H:125
void fillVolFrac(MultiFab &vfrac, const Geometry &geom) const
Populate volume fractions into vfrac using geom.
Definition AMReX_EB2_Level.cpp:483
const BoxArray & boxArray() const noexcept
Definition AMReX_EB2_Level.H:108
Level(IndexSpace const *is, const Geometry &geom)
Construct a Level bound to EB index space is and Geometry geom.
Definition AMReX_EB2_Level.H:113
friend class GShopLevel
Definition AMReX_EB2_Level.H:190
bool m_ok
Definition AMReX_EB2_Level.H:184
void fillCentroid(MultiCutFab ¢roid, const Geometry &geom) const
Populate cut-cell centroids into a MultiCutFab.
Definition AMReX_EB2_Level.cpp:544
DistributionMapping m_dmap
Definition AMReX_EB2_Level.H:166
void fillBndryArea(MultiCutFab &bndryarea, const Geometry &geom) const
Populate boundary-face areas into bndryarea.
Definition AMReX_EB2_Level.cpp:568
void buildCellFlag()
Build EBCellFlag data from the current implicit function field.
Definition AMReX_EB2_Level.cpp:403
void buildCutCellMask(Level const &fine_level)
Build the cut-cell mask by coarsening from fine_level.
Definition AMReX_EB2_Level.cpp:966
IndexSpace const * getEBIndexSpace() const noexcept
Definition AMReX_EB2_Level.H:122
LayoutData wrapper that allocates a GFab per FAB/Box.
Definition AMReX_EB2_MultiGFab.H:86
FAB storing EBCellFlag data per cell.
Definition AMReX_EBCellFlag.H:324
An Array of FortranArrayBox(FAB)-like Objects.
Definition AMReX_FabArray.H:350
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:211
Periodicity periodicity() const noexcept
Definition AMReX_Geometry.H:356
GpuArray< Real, 3 > ProbLoArray() const noexcept
Definition AMReX_Geometry.H:187
bool isPeriodic(int dir) const noexcept
Is the domain periodic in the specified direction?
Definition AMReX_Geometry.H:332
Definition AMReX_GpuElixir.H:13
__host__ static __device__ constexpr IntVectND< dim > TheZeroVector() noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition AMReX_IntVect.H:771
__host__ static __device__ constexpr IntVectND< dim > TheDimensionVector(int d) noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition AMReX_IntVect.H:790
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:172
FabArray façade that only allocates cut-cell Fabs (skipping regular cells).
Definition AMReX_MultiCutFab.H:87
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
Parse Parameters From Command Line and Input Files.
Definition AMReX_ParmParse.H:349
int queryAdd(std::string_view name, T &ref)
If name is found, the value in the ParmParse database will be stored in the ref argument....
Definition AMReX_ParmParse.H:1045
This class provides the user with a few print options.
Definition AMReX_Print.H:35
A Collection of IArrayBoxes.
Definition AMReX_iMultiFab.H:34
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
amrex_long Long
Definition AMReX_INT.H:30
__host__ __device__ BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
Coarsen BoxND by given (positive) coarsening ratio.
Definition AMReX_Box.H:1409
__host__ __device__ BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition AMReX_Box.H:1280
__host__ __device__ BoxND< dim > refine(const BoxND< dim > &b, int ref_ratio) noexcept
Definition AMReX_Box.H:1459
std::array< T, N > Array
Definition AMReX_Array.H:26
int MyProc() noexcept
Definition AMReX_ParallelDescriptor.H:128
int NProcs() noexcept
Definition AMReX_ParallelDescriptor.H:255
Definition AMReX_FabArrayBase.H:33
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
bool inLaunchRegion() noexcept
Definition AMReX_GpuControl.H:92
void write_stl(std::string const &filename, std::map< int, std::unique_ptr< MCFab > > const &mc_fabs)
Write the collected marching-cubes output to an STL file.
Definition AMReX_MarchingCubes.cpp:924
void marching_cubes(Geometry const &geom, FArrayBox &sdf_fab, MCFab &mc_fab)
Run marching cubes on signed-distance field sdf_fab.
Definition AMReX_MarchingCubes.cpp:693
void Initialize()
Initialize internal lookup tables and device buffers for marching cubes.
Definition AMReX_MarchingCubes.cpp:651
MPI_Comm CommunicatorSub() noexcept
sub-communicator for current frame
Definition AMReX_ParallelContext.H:70
Definition AMReX_Amr.cpp:49
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Return a BoxND with different type.
Definition AMReX_Box.H:1558
__host__ __device__ BoxND< dim > surroundingNodes(const BoxND< dim > &b, int dir) noexcept
Return a BoxND with NODE based coordinates in direction dir that encloses BoxND b....
Definition AMReX_Box.H:1522
RunOn
Definition AMReX_GpuControl.H:69
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
__host__ __device__ constexpr IntVectND< dim > scale(const IntVectND< dim > &p, int s) noexcept
Returns a IntVectND obtained by multiplying each of the components of this IntVectND by s.
Definition AMReX_IntVect.H:1105
int Verbose() noexcept
Definition AMReX.cpp:179
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:240
void AllGatherBoxes(Vector< Box > &bxs, int n_extra_reserve)
Definition AMReX_Box.cpp:124
A multidimensional array accessor.
Definition AMReX_Array4.H:283
FabArray memory allocation information.
Definition AMReX_FabArray.H:66
MFInfo & SetTag() noexcept
Definition AMReX_FabArray.H:79
Definition AMReX_MFIter.H:20
MFItInfo & SetDynamic(bool f) noexcept
Definition AMReX_MFIter.H:43