1 #ifndef AMREX_EB2_LEVEL_H_
2 #define AMREX_EB2_LEVEL_H_
3 #include <AMReX_Config.H>
22 #include <unordered_map>
25 #include <type_traits>
30 template <
typename G>
class GShopLevel;
112 template <
typename G>
134 template <
typename G>
139 template <
typename G>
144 if (std::is_same<typename G::FunctionType, AllRegularIF>::value) {
153 template <
typename G>
159 amrex::Print() <<
"AMReX WARNING: extend_domain_face=false is not recommended!\n";
164 #ifdef AMREX_USE_FLOAT
165 Real small_volfrac = 1.e-5_rt;
167 Real small_volfrac = 1.e-14;
169 bool cover_multiple_cuts =
false;
174 pp.
queryAdd(
"cover_multiple_cuts", cover_multiple_cuts);
177 maxiter =
std::min(100000, maxiter);
180 m_ngrow =
IntVect{
static_cast<int>(std::ceil(ngrow/16.)) * 16};
183 Box domain_grown = domain;
184 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
188 m_ngrow[idim] =
std::min(m_ngrow[idim], domain_grown.
length(idim));
191 domain_grown.
grow(m_ngrow);
202 for (
int clev = num_crse_opt; clev >= 0; --clev) {
209 covered_boxes.
clear();
213 test_boxes.
swap(cut_boxes);
214 test_boxes.
coarsen(crse_ratio);
218 const Long nboxes = test_boxes.
size();
219 const auto& boxes = test_boxes.
data();
220 for (Long i = iproc; i < nboxes; i += nprocs) {
221 const Box& vbx = boxes[i];
223 auto box_type = gshop.getBoxType(gbx&crse_bounding_box,crse_geom,
RunOn::Gpu);
224 if (box_type == gshop.allcovered) {
226 }
else if (box_type == gshop.mixedcells) {
238 auto grow_at_domain_boundary = [&] (
BoxList& bl)
243 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
244 if (m_ngrow[idim] != 0) {
245 if (
b.smallEnd(idim) == domlo[idim]) {
246 b.growLo(idim,m_ngrow[idim]);
248 if (
b.bigEnd(idim) == domhi[idim]) {
249 b.growHi(idim,m_ngrow[idim]);
255 grow_at_domain_boundary(covered_boxes);
256 grow_at_domain_boundary(cut_boxes);
262 amrex::Abort(
"AMReX_EB2_Level.H: Domain is completely covered");
265 if (!covered_boxes.
isEmpty()) {
266 if (num_crse_opt > 2) {
269 m_covered_grids =
BoxArray(std::move(covered_boxes));
280 m_grids =
BoxArray(std::move(cut_boxes));
283 m_mgf.define(m_grids, m_dmap);
286 mf_info.
SetTag(
"EB2::Level");
287 m_cellflag.define(m_grids, m_dmap, 1, ng, mf_info);
288 m_volfrac.define(m_grids, m_dmap, 1, ng, mf_info);
289 m_centroid.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
290 m_bndryarea.define(m_grids, m_dmap, 1, ng, mf_info);
291 m_bndrycent.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
292 m_bndrynorm.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
293 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
295 m_dmap, 1, ng, mf_info);
297 m_dmap, AMREX_SPACEDIM-1, ng, mf_info);
298 IntVect edge_type{1}; edge_type[idim] = 0;
299 m_edgecent[idim].define(
amrex::convert(m_grids, edge_type), m_dmap, 1, ng, mf_info);
305 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
318 for (; iter < maxiter; ++iter)
323 #pragma omp parallel if (Gpu::notInLaunchRegion()) reduction(+:nsmallcells,nmulticuts)
326 #if (AMREX_SPACEDIM == 3)
332 auto& gfab = m_mgf[mfi];
333 const Box& vbx = gfab.validbox();
335 auto&
levelset = gfab.getLevelSet();
337 gshop.fillFab(
levelset, geom, gshop_run_on, bounding_box);
364 auto& facetype = gfab.getFaceType();
372 #if (AMREX_SPACEDIM == 3)
373 auto& edgetype = gfab.getEdgeType();
386 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
387 edgetype[idim].prefetchToHost();
388 m_edgecent[idim][mfi].prefetchToHost();
393 gshop.getIntercept({xip,yip,zip}, {xdg,ydg,zdg}, clst,
394 geom, gshop_run_on, bounding_box);
398 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
399 edgetype[idim].prefetchToDevice();
400 m_edgecent[idim][mfi].prefetchToDevice();
405 gshop.updateIntercept({xip,yip,zip}, {xdg,ydg,zdg}, clst, geom);
408 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
409 const Box&
b = facetype[idim].box();
410 M2[idim].resize(
b,3);
416 nmc =
build_faces(vbx, cfg, ftx, fty, ftz, xdg, ydg, zdg, lst,
417 xip, yip, zip, apx, apy, apz, fcx, fcy, fcz,
418 xm2, ym2, zm2, dx, problo, cover_multiple_cuts);
420 cellflagtmp.
resize(m_cellflag[mfi].box());
424 build_cells(vbx, cfg, ftx, fty, ftz, apx, apy, apz,
425 fcx, fcy, fcz, xm2, ym2, zm2, dx, vfr, ctr,
426 bar, bct, bnm, cfgtmp, lst,
435 #elif (AMREX_SPACEDIM == 2)
443 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
444 facetype[idim].prefetchToHost();
445 m_edgecent[idim][mfi].prefetchToHost();
451 gshop.getIntercept({xip,yip},
452 {facetype[1].const_array(), facetype[0].const_array()},
453 clst, geom, gshop_run_on, bounding_box);
457 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
458 facetype[idim].prefetchToDevice();
459 m_edgecent[idim][mfi].prefetchToDevice();
464 gshop.updateIntercept({xip,yip},
465 {facetype[1].const_array(), facetype[0].const_array()},
469 nmc =
build_faces(vbx, cfg, ftx, fty, lst, xip, yip, apx, apy, fcx, fcy,
470 dx, problo, cover_multiple_cuts, nsm);
472 build_cells(vbx, cfg, ftx, fty, apx, apy, dx, vfr, ctr,
482 if (nsmallcells == 0 && nmulticuts == 0) {
485 auto ls = m_mgf.getLevelSet();
490 amrex::Print() <<
"AMReX EB: Iter. " << iter+1 <<
" fixed " << nsmallcells
491 <<
" small cells" <<
'\n';
494 amrex::Print() <<
"AMReX EB: Iter. " << iter+1 <<
" fixed " << nmulticuts
495 <<
" multicuts" <<
'\n';
504 #pragma omp parallel if (Gpu::notInLaunchRegion())
508 auto& gfab = m_mgf[mfi];
509 auto const&
levelset = gfab.getLevelSet();
515 #
if (AMREX_SPACEDIM == 3)
516 auto const& edgetype = gfab.getEdgeType();
523 #elif (AMREX_SPACEDIM == 2)
524 auto& facetype = gfab.getFaceType();
532 m_levelset = m_mgf.getLevelSet();
538 template <
typename G>
554 const int coarse_ratio = 2;
555 const int min_width = 8;
556 bool coarsenable = fine_grids.
coarsenable(coarse_ratio, min_width)
557 && (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:129
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:2098
AMREX_FORCE_INLINE Array4< T const > array() const noexcept
Definition: AMReX_BaseFab.H:379
Elixir elixir() noexcept
Definition: AMReX_BaseFab.H:2140
A collection of Boxes stored in an Array.
Definition: AMReX_BoxArray.H:550
bool coarsenable(int refinement_ratio, int min_width=1) const
Coarsen each Box in the BoxArray to the specified ratio.
bool empty() const noexcept
Return whether the BoxArray is empty.
Definition: AMReX_BoxArray.H:603
A class for managing a List of Boxes that share a common IndexType. This class implements operations ...
Definition: AMReX_BoxList.H:52
void clear()
Remove all Boxes from this BoxList.
Vector< Box > & data() noexcept
Returns a reference to the Vector<Box>.
Definition: AMReX_BoxList.H:215
BoxList & coarsen(int ratio)
Coarsen each Box in the BoxList by the ratio.
bool isEmpty() const noexcept
Is this BoxList empty?
Definition: AMReX_BoxList.H:135
BoxList & maxSize(int chunk)
Forces each Box in the BoxList to have sides of length <= chunk.
void swap(BoxList &rhs) noexcept
Definition: AMReX_BoxList.H:219
void push_back(const Box &bn)
Append a Box to this BoxList.
Definition: AMReX_BoxList.H:93
Long size() const noexcept
The number of Boxes in this BoxList.
Definition: AMReX_BoxList.H:113
AMREX_GPU_HOST_DEVICE BoxND & surroundingNodes() noexcept
Convert to NODE type in all directions.
Definition: AMReX_Box.H:946
AMREX_GPU_HOST_DEVICE const IntVectND< dim > & smallEnd() const &noexcept
Get the smallend of the BoxND.
Definition: AMReX_Box.H:105
AMREX_GPU_HOST_DEVICE BoxND & grow(int i) noexcept
Definition: AMReX_Box.H:627
AMREX_GPU_HOST_DEVICE const IntVectND< dim > & bigEnd() const &noexcept
Get the bigend.
Definition: AMReX_Box.H:116
AMREX_GPU_HOST_DEVICE IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition: AMReX_Box.H:146
AMREX_GPU_HOST_DEVICE bool coarsenable(const IntVectND< dim > &refrat, const IntVectND< dim > &min_width) const noexcept
Definition: AMReX_Box.H:745
GpuArray< Real, AMREX_SPACEDIM > CellSizeArray() const noexcept
Definition: AMReX_CoordSys.H:76
Calculates the distribution of FABs to MPI processes.
Definition: AMReX_DistributionMapping.H:41
static constexpr int ng
Definition: AMReX_EB2_MultiGFab.H:19
Definition: AMReX_EB2_Level.H:115
GShopLevel(IndexSpace const *is, const Geometry &geom)
Definition: AMReX_EB2_Level.H:135
void define_fine(G const &gshop, const Geometry &geom, int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
Definition: AMReX_EB2_Level.H:155
static GShopLevel< G > makeAllRegular(IndexSpace const *is, const Geometry &geom)
Definition: AMReX_EB2_Level.H:126
GShopLevel(IndexSpace const *is, int ilev, int max_grid_size, int ngrow, const Geometry &geom, GShopLevel< G > &fineLevel)
Definition: AMReX_EB2_Level.H:539
GShopLevel(IndexSpace const *is, G const &gshop, const Geometry &geom, int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
Definition: AMReX_EB2_Level.H:140
Definition: AMReX_EB2.H:26
Definition: AMReX_EB2_Level.H:33
MultiFab m_bndryarea
Definition: AMReX_EB2_Level.H:89
bool isOK() const noexcept
Definition: AMReX_EB2_Level.H:37
MultiFab m_volfrac
Definition: AMReX_EB2_Level.H:87
MultiGFab m_mgf
Definition: AMReX_EB2_Level.H:84
bool m_allregular
Definition: AMReX_EB2_Level.H:96
iMultiFab m_cutcellmask
Definition: AMReX_EB2_Level.H:95
IndexSpace const * getEBIndexSpace() const noexcept
Definition: AMReX_EB2_Level.H:63
void prepareForCoarsening(const Level &rhs, int max_grid_size, IntVect const &ngrow)
Definition: AMReX_EB2_Level.cpp:10
void fillCutCellMask(iMultiFab &cutcellmask, const Geometry &geom) const
Definition: AMReX_EB2_Level.cpp:915
Array< MultiFab, AMREX_SPACEDIM > m_areafrac
Definition: AMReX_EB2_Level.H:92
MultiFab m_centroid
Definition: AMReX_EB2_Level.H:88
IntVect const & nGrowVect() const noexcept
Definition: AMReX_EB2_Level.H:65
FabArray< EBCellFlagFab > m_cellflag
Definition: AMReX_EB2_Level.H:86
Array< MultiFab, AMREX_SPACEDIM > m_facecent
Definition: AMReX_EB2_Level.H:93
BoxArray m_covered_grids
Definition: AMReX_EB2_Level.H:82
void fillLevelSet(MultiFab &levelset, const Geometry &geom) const
Definition: AMReX_EB2_Level.cpp:880
Geometry m_geom
Definition: AMReX_EB2_Level.H:79
void fillFaceCent(Array< MultiCutFab *, AMREX_SPACEDIM > const &facecent, const Geometry &geom) const
Definition: AMReX_EB2_Level.cpp:775
IntVect m_ngrow
Definition: AMReX_EB2_Level.H:80
IndexSpace const * m_parent
Definition: AMReX_EB2_Level.H:99
BoxArray m_grids
Definition: AMReX_EB2_Level.H:81
Array< MultiFab, AMREX_SPACEDIM > m_edgecent
Definition: AMReX_EB2_Level.H:94
bool hasEBInfo() const noexcept
Definition: AMReX_EB2_Level.H:69
void write_to_chkpt_file(const std::string &fname, bool extend_domain_face, int max_grid_size) const
Definition: AMReX_EB2_Level.cpp:924
void fillEBCellFlag(FabArray< EBCellFlagFab > &cellflag, const Geometry &geom) const
Definition: AMReX_EB2_Level.cpp:426
int coarsenFromFine(Level &fineLevel, bool fill_boundary)
Definition: AMReX_EB2_Level.cpp:84
MultiFab m_bndrynorm
Definition: AMReX_EB2_Level.H:91
bool m_has_eb_info
Definition: AMReX_EB2_Level.H:98
MultiFab m_bndrycent
Definition: AMReX_EB2_Level.H:90
MultiFab m_levelset
Definition: AMReX_EB2_Level.H:85
void setRegularLevel()
Definition: AMReX_EB2_Level.H:105
bool isAllRegular() const noexcept
Definition: AMReX_EB2_Level.H:36
void fillEdgeCent(Array< MultiCutFab *, AMREX_SPACEDIM > const &edgecent, const Geometry &geom) const
Definition: AMReX_EB2_Level.cpp:813
void fillBndryCent(MultiCutFab &bndrycent, const Geometry &geom) const
Definition: AMReX_EB2_Level.cpp:581
void fillBndryNorm(MultiCutFab &bndrynorm, const Geometry &geom) const
Definition: AMReX_EB2_Level.cpp:605
const Geometry & Geom() const noexcept
Definition: AMReX_EB2_Level.H:62
const BoxArray & boxArray() const noexcept
Definition: AMReX_EB2_Level.H:56
void fillVolFrac(MultiFab &vfrac, const Geometry &geom) const
Definition: AMReX_EB2_Level.cpp:477
Level(IndexSpace const *is, const Geometry &geom)
Definition: AMReX_EB2_Level.H:59
const DistributionMapping & DistributionMap() const noexcept
Definition: AMReX_EB2_Level.H:57
friend class GShopLevel
Definition: AMReX_EB2_Level.H:102
bool m_ok
Definition: AMReX_EB2_Level.H:97
void fillCentroid(MultiCutFab ¢roid, const Geometry &geom) const
Definition: AMReX_EB2_Level.cpp:534
DistributionMapping m_dmap
Definition: AMReX_EB2_Level.H:83
void fillBndryArea(MultiCutFab &bndryarea, const Geometry &geom) const
Definition: AMReX_EB2_Level.cpp:558
void buildCellFlag()
Definition: AMReX_EB2_Level.cpp:400
void buildCutCellMask(Level const &fine_level)
Definition: AMReX_EB2_Level.cpp:934
void fillAreaFrac(Array< MultiCutFab *, AMREX_SPACEDIM > const &areafrac, const Geometry &geom) const
Definition: AMReX_EB2_Level.cpp:629
Definition: AMReX_EB2_MultiGFab.H:65
Definition: AMReX_EBCellFlag.H:287
An Array of FortranArrayBox(FAB)-like Objects.
Definition: AMReX_FabArray.H:344
Rectangular problem domain geometry.
Definition: AMReX_Geometry.H:73
GpuArray< Real, AMREX_SPACEDIM > ProbLoArray() const noexcept
Definition: AMReX_Geometry.H:186
Periodicity periodicity() const noexcept
Definition: AMReX_Geometry.H:355
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition: AMReX_Geometry.H:210
bool isPeriodic(int dir) const noexcept
Is the domain periodic in the specified direction?
Definition: AMReX_Geometry.H:331
Definition: AMReX_GpuElixir.H:13
AMREX_GPU_HOST_DEVICE static constexpr AMREX_FORCE_INLINE 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:672
AMREX_GPU_HOST_DEVICE static constexpr AMREX_FORCE_INLINE 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:691
Definition: AMReX_MFIter.H:57
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition: AMReX_MFIter.H:141
Definition: AMReX_MultiCutFab.H:81
A collection (stored as an array) of FArrayBox objects.
Definition: AMReX_MultiFab.H:38
Parse Parameters From Command Line and Input Files.
Definition: AMReX_ParmParse.H:320
int queryAdd(const char *name, T &ref)
If name is found, the value in the ParmParse database will be stored in the ref argument....
Definition: AMReX_ParmParse.H:993
This class provides the user with a few print options.
Definition: AMReX_Print.H:35
Definition: AMReX_iMultiFab.H:32
Definition: AMReX_FabArrayBase.H:32
AMREX_EXPORT int max_grid_size
Definition: AMReX_EB2.cpp:23
int build_faces(Box const &bx, Array4< EBCellFlag > const &cell, Array4< Type_t > const &fx, Array4< Type_t > const &fy, Array4< Real > const &levset, Array4< Real const > const &interx, Array4< Real const > const &intery, Array4< Real > const &apx, Array4< Real > const &apy, Array4< Real > const &fcx, Array4< Real > const &fcy, GpuArray< Real, AMREX_SPACEDIM > const &dx, GpuArray< Real, AMREX_SPACEDIM > const &problo, bool cover_multiple_cuts, int &nsmallfaces) noexcept
Definition: AMReX_EB2_2D_C.cpp:198
AMREX_EXPORT bool extend_domain_face
Definition: AMReX_EB2.cpp:24
void build_cells(Box const &bx, Array4< EBCellFlag > const &cell, Array4< Type_t > const &fx, Array4< Type_t > const &fy, Array4< Real > const &apx, Array4< Real > const &apy, GpuArray< Real, AMREX_SPACEDIM > const &dx, Array4< Real > const &vfrac, Array4< Real > const &vcent, Array4< Real > const &barea, Array4< Real > const &bcent, Array4< Real > const &bnorm, Array4< Real > const &levset, Real small_volfrac, Geometry const &geom, bool extend_domain_face, int &nsmallcells, int const nmulticuts) noexcept
Definition: AMReX_EB2_2D_C.cpp:352
void intercept_to_edge_centroid(AMREX_D_DECL(Array4< Real > const &excent, Array4< Real > const &eycent, Array4< Real > const &ezcent), AMREX_D_DECL(Array4< Type_t const > const &fx, Array4< Type_t const > const &fy, Array4< Type_t const > const &fz), Array4< Real const > const &levset, GpuArray< Real, AMREX_SPACEDIM > const &dx, GpuArray< Real, AMREX_SPACEDIM > const &problo) noexcept
Definition: AMReX_EB2_ND_C.cpp:5
void streamSynchronize() noexcept
Definition: AMReX_GpuDevice.H:237
bool inLaunchRegion() noexcept
Definition: AMReX_GpuControl.H:86
MPI_Comm CommunicatorSub() noexcept
sub-communicator for current frame
Definition: AMReX_ParallelContext.H:70
int MyProc() noexcept
return the rank number local to the current Parallel Context
Definition: AMReX_ParallelDescriptor.H:125
int NProcs() noexcept
return the number of MPI ranks local to the current Parallel Context
Definition: AMReX_ParallelDescriptor.H:243
@ min
Definition: AMReX_ParallelReduce.H:18
@ max
Definition: AMReX_ParallelReduce.H:17
RunOn
Definition: AMReX_GpuControl.H:69
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition: AMReX_Box.H:1211
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Returns a BoxND with different type.
Definition: AMReX_Box.H:1435
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition: AMReX.H:111
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > surroundingNodes(const BoxND< dim > &b, int dir) noexcept
Returns a BoxND with NODE based coordinates in direction dir that encloses BoxND b....
Definition: AMReX_Box.H:1399
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE 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:1006
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > refine(const BoxND< dim > &b, int ref_ratio) noexcept
Definition: AMReX_Box.H:1342
int Verbose() noexcept
Definition: AMReX.cpp:164
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
Coarsen BoxND by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo/rati...
Definition: AMReX_Box.H:1304
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition: AMReX.cpp:225
void AllGatherBoxes(Vector< Box > &bxs, int n_extra_reserve)
Definition: AMReX_Box.cpp:120
std::array< T, N > Array
Definition: AMReX_Array.H:24
Definition: AMReX_Array4.H:61
FabArray memory allocation information.
Definition: AMReX_FabArray.H:66
MFInfo & SetTag() noexcept
Definition: AMReX_FabArray.H:79