Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Vis multicuts #4146

Open
wants to merge 4 commits into
base: development
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
9 changes: 9 additions & 0 deletions Docs/sphinx_documentation/source/RuntimeParameters.rst
Original file line number Diff line number Diff line change
Expand Up @@ -787,6 +787,15 @@ Embedded Boundary
.. tip:: Because AMReX currently does not support multi-cut cells, it
would be a runtime error if multi-cut cells are left unfixed.

.. py:data:: eb2.plt_multiple_cuts
:type: bool
:value: false

If this parameter is set to true and multi-cut cells are present,
we produce plot files to see the location of the multi-cut issue.
In the plot file, a cell value of 2 indicates regular, 0 means covered,
and values greater than 2 indicate multi-cut cells.

.. py:data:: eb2.maxiter
:type: int
:value: 32
Expand Down
12 changes: 10 additions & 2 deletions Src/EB/AMReX_EB2_2D_C.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -203,7 +203,8 @@ int build_faces (Box const& bx, Array4<EBCellFlag> const& cell,
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
bool cover_multiple_cuts, int& nsmallfaces,
bool plt_multiple_cuts, Array4<Real> const& mcx) noexcept
{
#ifdef AMREX_USE_FLOAT
constexpr Real small = 1.e-5_rt;
Expand Down Expand Up @@ -311,6 +312,10 @@ int build_faces (Box const& bx, Array4<EBCellFlag> const& cell,
if (ncuts > 2) {
Gpu::Atomic::Add(dp,1);
}
if (plt_multiple_cuts)
{
mcx(i,j,k) = ncuts;
}
}
}
});
Expand Down Expand Up @@ -342,7 +347,10 @@ int build_faces (Box const& bx, Array4<EBCellFlag> const& cell,
nsmallfaces += *(hp+1);

if (*hp > 0 && !cover_multiple_cuts) {
amrex::Abort("amrex::EB2::build_faces: more than 2 cuts not supported");
if (!plt_multiple_cuts)
{
amrex::Abort("amrex::EB2::build_faces: more than 2 cuts not supported");
}
}

return *hp;
Expand Down
28 changes: 24 additions & 4 deletions Src/EB/AMReX_EB2_3D_C.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -370,7 +370,8 @@ int build_faces (Box const& bx, Array4<EBCellFlag> const& cell,
Array4<Real> const& m2z,
GpuArray<Real,AMREX_SPACEDIM> const& dx,
GpuArray<Real,AMREX_SPACEDIM> const& problo,
bool cover_multiple_cuts) noexcept
bool cover_multiple_cuts, bool plt_multiple_cuts,
Array4<Real> const& mcx, Array4<Real> const& mcy, Array4<Real> const& mcz) noexcept
{
Gpu::Buffer<int> nmulticuts = {0};
int* hp = nmulticuts.hostData();
Expand Down Expand Up @@ -491,6 +492,11 @@ int build_faces (Box const& bx, Array4<EBCellFlag> const& cell,
} else if (apx(i,j,k) == 1.0_rt) {
fx(i,j,k) = Type::regular;
}

if (plt_multiple_cuts){
mcx(i,j,k) = ncuts;
}

}
});

Expand Down Expand Up @@ -599,6 +605,10 @@ int build_faces (Box const& bx, Array4<EBCellFlag> const& cell,
} else if (apy(i,j,k) == 1.0_rt) {
fy(i,j,k) = Type::regular;
}

if (plt_multiple_cuts){
mcy(i,j,k) = ncuts;
}
}
});

Expand Down Expand Up @@ -707,6 +717,10 @@ int build_faces (Box const& bx, Array4<EBCellFlag> const& cell,
} else if (apz(i,j,k) == 1.0_rt) {
fz(i,j,k) = Type::regular;
}

if (plt_multiple_cuts){
mcz(i,j,k) = ncuts;
}
}
});

Expand Down Expand Up @@ -768,7 +782,10 @@ int build_faces (Box const& bx, Array4<EBCellFlag> const& cell,
}
});
} else {
amrex::Abort("amrex::EB2::build_faces: more than 2 cuts not supported");
if (!plt_multiple_cuts)
{
amrex::Abort("amrex::EB2::build_faces: more than 2 cuts not supported");
}
}
}

Expand All @@ -787,7 +804,7 @@ void build_cells (Box const& bx, Array4<EBCellFlag> const& cell,
Array4<Real> const& barea, Array4<Real> const& bcent,
Array4<Real> const& bnorm, Array4<EBCellFlag> const& ctmp,
Array4<Real> const& levset, Real small_volfrac, Geometry const& geom,
bool extend_domain_face, bool cover_multiple_cuts,
bool extend_domain_face, bool cover_multiple_cuts, bool plt_multiple_cuts,
int& nsmallcells, int& nmulticuts) noexcept
{
Gpu::Buffer<int> n_smallcell_multicuts = {0,0};
Expand Down Expand Up @@ -932,7 +949,10 @@ void build_cells (Box const& bx, Array4<EBCellFlag> const& cell,

if (nsmallcells > 0 || nmulticuts > 0) {
if (!cover_multiple_cuts && nmulticuts > 0) {
amrex::Abort("amrex::EB2::build_cells: multi-cuts not supported");
if (!plt_multiple_cuts)
{
amrex::Abort("amrex::EB2::build_cells: multi-cuts not supported");
}
}
return;
} else {
Expand Down
8 changes: 5 additions & 3 deletions Src/EB/AMReX_EB2_C.H
Original file line number Diff line number Diff line change
Expand Up @@ -25,7 +25,8 @@ int build_faces (Box const& bx, Array4<EBCellFlag> const& cell,
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;
bool cover_multiple_cuts, int& nsmallfaces,
bool plt_multiple_cuts, Array4<Real> const& mcx) noexcept;

void build_cells (Box const& bx, Array4<EBCellFlag> const& cell,
Array4<Type_t> const& fx, Array4<Type_t> const& fy,
Expand Down Expand Up @@ -55,7 +56,8 @@ int build_faces (Box const& bx, Array4<EBCellFlag> const& cell,
Array4<Real> const& m2z,
GpuArray<Real,AMREX_SPACEDIM> const& dx,
GpuArray<Real,AMREX_SPACEDIM> const& problo,
bool cover_multiple_cuts) noexcept;
bool cover_multiple_cuts, bool plt_multiple_cuts,
Array4<Real> const& mcx, Array4<Real> const& mcy, Array4<Real> const& mcz) noexcept;

void build_cells (Box const& bx, Array4<EBCellFlag> const& cell,
Array4<Type_t> const& fx, Array4<Type_t> const& fy,
Expand All @@ -69,7 +71,7 @@ void build_cells (Box const& bx, Array4<EBCellFlag> const& cell,
Array4<Real> const& barea, Array4<Real> const& bcent,
Array4<Real> const& bnorm, Array4<EBCellFlag> const& ctmp,
Array4<Real> const& levset, Real small_volfrac, Geometry const& geom,
bool extend_domain_face, bool cover_multiple_cuts,
bool extend_domain_face, bool cover_multiple_cuts, bool plt_multiple_cuts,
int& nsmallcells, int& nmulticuts) noexcept;

void set_connection_flags(Box const& bx, Box const& bxg1,
Expand Down
52 changes: 49 additions & 3 deletions Src/EB/AMReX_EB2_Level.H
Original file line number Diff line number Diff line change
Expand Up @@ -14,6 +14,7 @@
#include <AMReX_EB2_MultiGFab.H>
#include <AMReX_EB2_C.H>
#include <AMReX_EB2_IF_AllRegular.H>
#include <AMReX_PlotFileUtil.H>

#ifdef AMREX_USE_OMP
#include <omp.h>
Expand Down Expand Up @@ -92,6 +93,9 @@ protected:
Array<MultiFab,AMREX_SPACEDIM> m_areafrac;
Array<MultiFab,AMREX_SPACEDIM> m_facecent;
Array<MultiFab,AMREX_SPACEDIM> m_edgecent;
MultiFab m_multicut_fcx;
MultiFab m_multicut_fcy;
MultiFab m_multicut_fcz;
iMultiFab m_cutcellmask;
bool m_allregular = false;
bool m_ok = false;
Expand Down Expand Up @@ -167,12 +171,14 @@ GShopLevel<G>::define_fine (G const& gshop, const Geometry& geom,
Real small_volfrac = 1.e-14;
#endif
bool cover_multiple_cuts = false;
bool plt_multiple_cuts = false;
int maxiter = 32;
{
ParmParse pp("eb2");
pp.queryAdd("small_volfrac", small_volfrac);
pp.queryAdd("cover_multiple_cuts", cover_multiple_cuts);
pp.queryAdd("maxiter", maxiter);
pp.queryAdd("plt_multiple_cuts", plt_multiple_cuts);
}
maxiter = std::min(100000, maxiter);

Expand Down Expand Up @@ -290,6 +296,9 @@ GShopLevel<G>::define_fine (G const& gshop, const Geometry& geom,
m_bndryarea.define(m_grids, m_dmap, 1, ng, mf_info);
m_bndrycent.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
m_bndrynorm.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
m_multicut_fcx.define(m_grids, m_dmap, 1, ng, mf_info);
m_multicut_fcy.define(m_grids, m_dmap, 1, ng, mf_info);
m_multicut_fcz.define(m_grids, m_dmap, 1, ng, mf_info);
for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
m_areafrac[idim].define(amrex::convert(m_grids, IntVect::TheDimensionVector(idim)),
m_dmap, 1, ng, mf_info);
Expand Down Expand Up @@ -413,9 +422,14 @@ GShopLevel<G>::define_fine (G const& gshop, const Geometry& geom,
Array4<Real> const& ym2 = M2[1].array();
Array4<Real> const& zm2 = M2[2].array();

Array4<Real> const& mcx = m_multicut_fcx.array(mfi);
Array4<Real> const& mcy = m_multicut_fcy.array(mfi);
Array4<Real> const& mcz = m_multicut_fcz.array(mfi);

nmc = build_faces(vbx, cfg, ftx, fty, ftz, xdg, ydg, zdg, lst,
xip, yip, zip, apx, apy, apz, fcx, fcy, fcz,
xm2, ym2, zm2, dx, problo, cover_multiple_cuts);
xm2, ym2, zm2, dx, problo, cover_multiple_cuts,
plt_multiple_cuts, mcx, mcy, mcz);

cellflagtmp.resize(m_cellflag[mfi].box());
Elixir cellflagtmp_eli = cellflagtmp.elixir();
Expand All @@ -424,7 +438,8 @@ GShopLevel<G>::define_fine (G const& gshop, const Geometry& geom,
build_cells(vbx, cfg, ftx, fty, ftz, apx, apy, apz,
fcx, fcy, fcz, xm2, ym2, zm2, dx, vfr, ctr,
bar, bct, bnm, cfgtmp, lst,
small_volfrac, geom, extend_domain_face, cover_multiple_cuts,
small_volfrac, geom, extend_domain_face,
cover_multiple_cuts, plt_multiple_cuts,
nsm, nmc);

// Because it is used in a synchronous reduction kernel in
Expand Down Expand Up @@ -466,8 +481,11 @@ GShopLevel<G>::define_fine (G const& gshop, const Geometry& geom,
clst, geom);
}

Array4<Real> const& mcx = m_multicut_fcx.array(mfi);

nmc = build_faces(vbx, cfg, ftx, fty, lst, xip, yip, apx, apy, fcx, fcy,
dx, problo, cover_multiple_cuts, nsm);
dx, problo, cover_multiple_cuts, nsm,
plt_multiple_cuts, mcx);

build_cells(vbx, cfg, ftx, fty, apx, apy, dx, vfr, ctr,
bar, bct, bnm, lst, small_volfrac, geom, extend_domain_face,
Expand All @@ -479,6 +497,34 @@ GShopLevel<G>::define_fine (G const& gshop, const Geometry& geom,
}

ParallelAllReduce::Sum<int>({nsmallcells,nmulticuts}, ParallelContext::CommunicatorSub());

if (plt_multiple_cuts && nmulticuts > 0)
{
amrex::Print() << "Total number of multicuts = " << nmulticuts << "\n";
amrex::Print() << "Plotting multicut locations..." << "\n";

#if (AMREX_SPACEDIM == 3)
if (m_multicut_fcx.max(0) > 2)
{
WriteSingleLevelPlotfile("plt.multicut.x", m_multicut_fcx, {"multicut_fcx"}, geom, 0.0, 0);
}
if (m_multicut_fcy.max(0) > 2)
{
WriteSingleLevelPlotfile("plt.multicut.y", m_multicut_fcy, {"multicut_fcy"}, geom, 0.0, 0);
}
if (m_multicut_fcz.max(0) > 2)
{
WriteSingleLevelPlotfile("plt.multicut.z", m_multicut_fcz, {"multicut_fcz"}, geom, 0.0, 0);
}
#elif (AMREX_SPACEDIM == 2)
if (m_multicut_fcx.max(0) > 2)
{
WriteSingleLevelPlotfile("plt.multicut", m_multicut_fcx, {"multicut_fcx"}, geom, 0.0, 0);
}
#endif
amrex::Abort("amrex::EB2:: more than 2 cuts not supported");
}

if (nsmallcells == 0 && nmulticuts == 0) {
break;
} else {
Expand Down