1#ifndef AMREX_EB2_LEVEL_H_
2#define AMREX_EB2_LEVEL_H_
3#include <AMReX_Config.H>
22#include <unordered_map>
30template <
typename G>
class GShopLevel;
72 void setShift (
int direction,
int ncells);
147 if (std::is_same<typename G::FunctionType, AllRegularIF>::value) {
162 amrex::Print() <<
"AMReX WARNING: extend_domain_face=false is not recommended!\n";
167#ifdef AMREX_USE_FLOAT
168 Real small_volfrac = 1.e-5_rt;
170 Real small_volfrac = 1.e-14;
172 bool cover_multiple_cuts =
false;
177 pp.
queryAdd(
"cover_multiple_cuts", cover_multiple_cuts);
180 maxiter = std::min(100000, maxiter);
183 m_ngrow =
IntVect{
static_cast<int>(std::ceil(ngrow/16.)) * 16};
186 Box domain_grown = domain;
187 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
191 m_ngrow[idim] = std::min(m_ngrow[idim], domain_grown.
length(idim));
194 domain_grown.
grow(m_ngrow);
204 num_crse_opt = std::max(0,std::min(8,num_crse_opt));
205 for (
int clev = num_crse_opt; clev >= 0; --clev) {
212 covered_boxes.
clear();
216 test_boxes.
swap(cut_boxes);
217 test_boxes.
coarsen(crse_ratio);
221 const Long nboxes = test_boxes.
size();
222 const auto& boxes = test_boxes.
data();
223 for (Long i = iproc; i < nboxes; i += nprocs) {
224 const Box& vbx = boxes[i];
226 auto box_type = gshop.getBoxType(gbx&crse_bounding_box,crse_geom,
RunOn::Gpu);
227 if (box_type == gshop.allcovered) {
229 }
else if (box_type == gshop.mixedcells) {
241 auto grow_at_domain_boundary = [&] (
BoxList& bl)
246 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
247 if (m_ngrow[idim] != 0) {
248 if (
b.smallEnd(idim) == domlo[idim]) {
249 b.growLo(idim,m_ngrow[idim]);
251 if (
b.bigEnd(idim) == domhi[idim]) {
252 b.growHi(idim,m_ngrow[idim]);
258 grow_at_domain_boundary(covered_boxes);
259 grow_at_domain_boundary(cut_boxes);
265 amrex::Abort(
"AMReX_EB2_Level.H: Domain is completely covered");
268 if (!covered_boxes.
isEmpty()) {
269 if (num_crse_opt > 2) {
272 m_covered_grids =
BoxArray(std::move(covered_boxes));
283 m_grids =
BoxArray(std::move(cut_boxes));
286 m_mgf.define(m_grids, m_dmap);
289 mf_info.
SetTag(
"EB2::Level");
290 m_cellflag.define(m_grids, m_dmap, 1, ng, mf_info);
291 m_volfrac.define(m_grids, m_dmap, 1, ng, mf_info);
292 m_centroid.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
293 m_bndryarea.define(m_grids, m_dmap, 1, ng, mf_info);
294 m_bndrycent.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
295 m_bndrynorm.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
296 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
298 m_dmap, 1, ng, mf_info);
300 m_dmap, AMREX_SPACEDIM-1, ng, mf_info);
301 IntVect edge_type{1}; edge_type[idim] = 0;
302 m_edgecent[idim].define(
amrex::convert(m_grids, edge_type), m_dmap, 1, ng, mf_info);
308 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
321 for (; iter < maxiter; ++iter)
326#pragma omp parallel if (Gpu::notInLaunchRegion()) reduction(+:nsmallcells,nmulticuts)
329#if (AMREX_SPACEDIM == 3)
335 auto& gfab = m_mgf[mfi];
336 const Box& vbx = gfab.validbox();
338 auto&
levelset = gfab.getLevelSet();
340 gshop.fillFab(
levelset, geom, gshop_run_on, bounding_box);
367 auto& facetype = gfab.getFaceType();
375#if (AMREX_SPACEDIM == 3)
376 auto& edgetype = gfab.getEdgeType();
389 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
390 edgetype[idim].prefetchToHost();
391 m_edgecent[idim][mfi].prefetchToHost();
396 gshop.getIntercept({xip,yip,zip}, {xdg,ydg,zdg}, clst,
397 geom, gshop_run_on, bounding_box);
401 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
402 edgetype[idim].prefetchToDevice();
403 m_edgecent[idim][mfi].prefetchToDevice();
408 gshop.updateIntercept({xip,yip,zip}, {xdg,ydg,zdg}, clst, geom);
411 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
412 const Box&
b = facetype[idim].box();
419 nmc =
build_faces(vbx, cfg, ftx, fty, ftz, xdg, ydg, zdg, lst,
421 xm2, ym2, zm2, dx, problo, cover_multiple_cuts);
423 cellflagtmp.
resize(m_cellflag[mfi].box());
428 fcx,
fcy,
fcz, xm2, ym2, zm2, dx, vfr, ctr,
429 bar, bct, bnm, cfgtmp, lst,
438#elif (AMREX_SPACEDIM == 2)
446 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
447 facetype[idim].prefetchToHost();
448 m_edgecent[idim][mfi].prefetchToHost();
454 gshop.getIntercept({xip,yip},
455 {facetype[1].const_array(), facetype[0].const_array()},
456 clst, geom, gshop_run_on, bounding_box);
460 for (
int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
461 facetype[idim].prefetchToDevice();
462 m_edgecent[idim][mfi].prefetchToDevice();
467 gshop.updateIntercept({xip,yip},
468 {facetype[1].const_array(), facetype[0].const_array()},
472 nmc =
build_faces(vbx, cfg, ftx, fty, lst, xip, yip,
apx,
apy,
fcx,
fcy,
473 dx, problo, cover_multiple_cuts, nsm);
485 if (nsmallcells == 0 && nmulticuts == 0) {
488 auto ls = m_mgf.getLevelSet();
493 amrex::Print() <<
"AMReX EB: Iter. " << iter+1 <<
" fixed " << nsmallcells
494 <<
" small cells" <<
'\n';
497 amrex::Print() <<
"AMReX EB: Iter. " << iter+1 <<
" fixed " << nmulticuts
498 <<
" multicuts" <<
'\n';
507#pragma omp parallel if (Gpu::notInLaunchRegion())
511 auto& gfab = m_mgf[mfi];
512 auto const&
levelset = gfab.getLevelSet();
518#
if (AMREX_SPACEDIM == 3)
519 auto const& edgetype = gfab.getEdgeType();
526#elif (AMREX_SPACEDIM == 2)
527 auto& facetype = gfab.getFaceType();
535 m_levelset = m_mgf.getLevelSet();
557 const int coarse_ratio = 2;
558 const int min_width = 8;
559 bool coarsenable = fine_grids.
coarsenable(coarse_ratio, min_width)
560 && (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:1622
Array4< T const > array() const noexcept
Definition AMReX_BaseFab.H:375
Elixir elixir() noexcept
Definition AMReX_BaseFab.H:1664
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:551
bool coarsenable(int refinement_ratio, int min_width=1) const
Coarsen each Box in the BoxArray to the specified ratio.
Definition AMReX_BoxArray.cpp:608
bool empty() const noexcept
Return whether the BoxArray is empty.
Definition AMReX_BoxArray.H:604
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.
Definition AMReX_BoxList.cpp:65
BoxList & coarsen(int ratio)
Coarsen each Box in the BoxList by the ratio.
Definition AMReX_BoxList.cpp:525
Vector< Box > & data() noexcept
Returns a reference to the Vector<Box>.
Definition AMReX_BoxList.H:215
bool isEmpty() const noexcept
Is this BoxList empty?
Definition AMReX_BoxList.H:135
void swap(BoxList &rhs) noexcept
Definition AMReX_BoxList.H:219
BoxList & maxSize(int chunk)
Forces each Box in the BoxList to have sides of length <= chunk.
Definition AMReX_BoxList.cpp:821
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
__host__ __device__ BoxND & grow(int i) noexcept
Definition AMReX_Box.H:630
__host__ __device__ const IntVectND< dim > & bigEnd() const &noexcept
Get the bigend.
Definition AMReX_Box.H:119
__host__ __device__ BoxND< new_dim > resize() const noexcept
Returns a new BoxND of size new_dim by either shrinking or expanding this BoxND.
Definition AMReX_Box.H:849
__host__ __device__ IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:149
__host__ __device__ bool coarsenable(const IntVectND< dim > &refrat, const IntVectND< dim > &min_width) const noexcept
Definition AMReX_Box.H:748
__host__ __device__ BoxND & surroundingNodes() noexcept
Convert to NODE type in all directions.
Definition AMReX_Box.H:964
__host__ __device__ const IntVectND< dim > & smallEnd() const &noexcept
Get the smallend of the BoxND.
Definition AMReX_Box.H:108
GpuArray< Real, 3 > 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:118
GShopLevel(IndexSpace const *is, const Geometry &geom)
Definition AMReX_EB2_Level.H:138
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:158
GShopLevel(IndexSpace const *is, int ilev, int max_grid_size, int ngrow, const Geometry &geom, GShopLevel< G > &fineLevel)
Definition AMReX_EB2_Level.H:542
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:143
static GShopLevel< G > makeAllRegular(IndexSpace const *is, const Geometry &geom)
Definition AMReX_EB2_Level.H:129
Definition AMReX_EB2.H:26
Definition AMReX_EB2_Level.H:33
MultiFab m_bndryarea
Definition AMReX_EB2_Level.H:91
bool isOK() const noexcept
Definition AMReX_EB2_Level.H:37
void fillFaceCent(Array< MultiCutFab *, 3 > const &facecent, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:791
MultiFab m_volfrac
Definition AMReX_EB2_Level.H:89
MultiGFab m_mgf
Definition AMReX_EB2_Level.H:86
bool m_allregular
Definition AMReX_EB2_Level.H:98
IntVect m_shift
Definition AMReX_EB2_Level.H:101
iMultiFab m_cutcellmask
Definition AMReX_EB2_Level.H:97
Array< MultiFab, 3 > m_facecent
Definition AMReX_EB2_Level.H:95
Array< MultiFab, 3 > m_edgecent
Definition AMReX_EB2_Level.H:96
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:943
MultiFab m_centroid
Definition AMReX_EB2_Level.H:90
FabArray< EBCellFlagFab > m_cellflag
Definition AMReX_EB2_Level.H:88
const DistributionMapping & DistributionMap() const noexcept
Definition AMReX_EB2_Level.H:57
BoxArray m_covered_grids
Definition AMReX_EB2_Level.H:84
void fillLevelSet(MultiFab &levelset, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:904
Geometry m_geom
Definition AMReX_EB2_Level.H:81
Array< MultiFab, 3 > m_areafrac
Definition AMReX_EB2_Level.H:94
IntVect m_ngrow
Definition AMReX_EB2_Level.H:82
IndexSpace const * m_parent
Definition AMReX_EB2_Level.H:102
BoxArray m_grids
Definition AMReX_EB2_Level.H:83
const Geometry & Geom() const noexcept
Definition AMReX_EB2_Level.H:62
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:953
void fillEBCellFlag(FabArray< EBCellFlagFab > &cellflag, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:426
void fillEdgeCent(Array< MultiCutFab *, 3 > const &edgecent, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:832
int coarsenFromFine(Level &fineLevel, bool fill_boundary)
Definition AMReX_EB2_Level.cpp:84
MultiFab m_bndrynorm
Definition AMReX_EB2_Level.H:93
bool m_has_eb_info
Definition AMReX_EB2_Level.H:100
void fillAreaFrac(Array< MultiCutFab *, 3 > const &areafrac, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:637
MultiFab m_bndrycent
Definition AMReX_EB2_Level.H:92
MultiFab m_levelset
Definition AMReX_EB2_Level.H:87
void setRegularLevel()
Definition AMReX_EB2_Level.H:108
void setShift(int direction, int ncells)
Definition AMReX_EB2_Level.cpp:1044
bool isAllRegular() const noexcept
Definition AMReX_EB2_Level.H:36
void fillBndryCent(MultiCutFab &bndrycent, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:589
void fillBndryNorm(MultiCutFab &bndrynorm, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:613
IntVect const & nGrowVect() const noexcept
Definition AMReX_EB2_Level.H:65
void fillVolFrac(MultiFab &vfrac, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:480
const BoxArray & boxArray() const noexcept
Definition AMReX_EB2_Level.H:56
Level(IndexSpace const *is, const Geometry &geom)
Definition AMReX_EB2_Level.H:59
friend class GShopLevel
Definition AMReX_EB2_Level.H:105
bool m_ok
Definition AMReX_EB2_Level.H:99
void fillCentroid(MultiCutFab ¢roid, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:541
DistributionMapping m_dmap
Definition AMReX_EB2_Level.H:85
void fillBndryArea(MultiCutFab &bndryarea, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:565
void buildCellFlag()
Definition AMReX_EB2_Level.cpp:400
void buildCutCellMask(Level const &fine_level)
Definition AMReX_EB2_Level.cpp:963
IndexSpace const * getEBIndexSpace() const noexcept
Definition AMReX_EB2_Level.H:63
Definition AMReX_EB2_MultiGFab.H:65
Definition AMReX_EBCellFlag.H:287
An Array of FortranArrayBox(FAB)-like Objects.
Definition AMReX_FabArray.H:345
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:73
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:210
Periodicity periodicity() const noexcept
Definition AMReX_Geometry.H:355
GpuArray< Real, 3 > ProbLoArray() const noexcept
Definition AMReX_Geometry.H:186
bool isPeriodic(int dir) const noexcept
Is the domain periodic in the specified direction?
Definition AMReX_Geometry.H:331
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:677
__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:696
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:343
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:1037
This class provides the user with a few print options.
Definition AMReX_Print.H:35
Definition AMReX_iMultiFab.H:32
Definition AMReX_FabArrayBase.H:33
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, 3 > 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
int max_grid_size
Definition AMReX_EB2.cpp:23
bool extend_domain_face
Definition AMReX_EB2.cpp:24
void intercept_to_edge_centroid(Array4< Real > const &excent, Array4< Real > const &eycent, Array4< Real > const &ezcent, 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, 3 > const &dx, GpuArray< Real, 3 > const &problo) noexcept
Definition AMReX_EB2_ND_C.cpp:5
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, 3 > const &dx, GpuArray< Real, 3 > const &problo, bool cover_multiple_cuts, int &nsmallfaces) noexcept
Definition AMReX_EB2_2D_C.cpp:198
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:260
bool inLaunchRegion() noexcept
Definition AMReX_GpuControl.H:92
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
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:138
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Returns a BoxND with different type.
Definition AMReX_Box.H:1453
__host__ __device__ 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:1417
RunOn
Definition AMReX_GpuControl.H:69
__host__ __device__ 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:1322
__host__ __device__ BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition AMReX_Box.H:1229
int Verbose() noexcept
Definition AMReX.cpp:169
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:230
void AllGatherBoxes(Vector< Box > &bxs, int n_extra_reserve)
Definition AMReX_Box.cpp:120
__host__ __device__ BoxND< dim > refine(const BoxND< dim > &b, int ref_ratio) noexcept
Definition AMReX_Box.H:1360
__host__ __device__ 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:1011
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