From b81cc0ede23ee68b972c177e1ecd37b551484a19 Mon Sep 17 00:00:00 2001 From: Weiqun Zhang Date: Sun, 1 Oct 2023 11:28:49 -0700 Subject: [PATCH] Fix Boundary Centroid in a Corner Case in 2D (#3568) Due to roundoff errors, It's possible to have a cell with area fraction of 1-epsilon, 1, 1, and 1 at its four faces. In this corner case, the previous way of computing boundary centroid was not correct. In this PR, it's fixed by using the level set. --- Src/EB/AMReX_EB2_2D_C.cpp | 39 +++++++++++++++++++++++++++++++++------ 1 file changed, 33 insertions(+), 6 deletions(-) diff --git a/Src/EB/AMReX_EB2_2D_C.cpp b/Src/EB/AMReX_EB2_2D_C.cpp index b77b2ebd00e..b2bbde200c5 100644 --- a/Src/EB/AMReX_EB2_2D_C.cpp +++ b/Src/EB/AMReX_EB2_2D_C.cpp @@ -9,7 +9,7 @@ void set_eb_data (const int i, const int j, GpuArray const& dx, Array4 const& vfrac, Array4 const& vcent, Array4 const& barea, Array4 const& bcent, - Array4 const& bnorm) noexcept + Array4 const& bnorm, Array4 const& levset) noexcept { #ifdef AMREX_USE_FLOAT constexpr Real almostone = 1.0_rt-1.e-6_rt; @@ -37,7 +37,20 @@ void set_eb_data (const int i, const int j, const Real nyabs = std::abs(ny); Real x_ym, x_yp, y_xm, y_xp; - if (nx > 0.0_rt) { + if (nx == 0.0_rt) { + if (apx(i,j,0) == 1.0_rt && apx(i+1,j,0) == 1.0_rt) { + if (levset(i,j,0) > 0.0_rt || levset(i,j+1,0) > 0.0_rt) { + x_ym = 0.5_rt*dx[0] - aym; + x_yp = 0.5_rt*dx[0] - ayp; + } else { + x_ym = -0.5_rt*dx[0] + aym; + x_yp = -0.5_rt*dx[0] + ayp; + } + } else { + x_ym = 0.0_rt; + x_yp = 0.0_rt; + } + } else if (nx > 0.0_rt) { x_ym = -0.5_rt*dx[0] + aym; x_yp = -0.5_rt*dx[0] + ayp; } else { @@ -45,7 +58,20 @@ void set_eb_data (const int i, const int j, x_yp = 0.5_rt*dx[0] - ayp; } - if (ny > 0.0_rt) { + if (ny == 0.0_rt) { + if (apy(i,j,0) == 1.0_rt && apy(i,j+1,0) == 1.0_rt) { + if (levset(i,j,0) > 0.0_rt || levset(i+1,j,0) > 0.0_rt) { + y_xm = 0.5_rt*dx[1] - axm; + y_xp = 0.5_rt*dx[1] - axp; + } else { + y_xm = -0.5_rt*dx[1] + axm; + y_xp = -0.5_rt*dx[1] + axp; + } + } else { + y_xm = 0.0_rt; + y_xp = 0.0_rt; + } + } else if (ny > 0.0_rt) { y_xm = -0.5_rt*dx[1] + axm; y_xp = -0.5_rt*dx[1] + axp; } else { @@ -135,7 +161,8 @@ bool set_eb_cell (int i, int j, Array4 const& cell, GpuArray const& dx, Array4 const& vfrac, Array4 const& vcent, Array4 const& barea, Array4 const& bcent, - Array4 const& bnorm, Real small_volfrac) noexcept + Array4 const& bnorm, Array4 const& levset, + Real small_volfrac) noexcept { bool is_small_cell = false; if (cell(i,j,0).isRegular()) { @@ -157,7 +184,7 @@ bool set_eb_cell (int i, int j, Array4 const& cell, bnorm(i,j,0,0) = 0.0_rt; bnorm(i,j,0,1) = 0.0_rt; } else { - set_eb_data(i,j,apx,apy,dx,vfrac,vcent,barea,bcent,bnorm); + set_eb_data(i,j,apx,apy,dx,vfrac,vcent,barea,bcent,bnorm,levset); // remove small cells if (vfrac(i,j,0) < small_volfrac) { set_covered(i,j,cell,vfrac,vcent,barea,bcent,bnorm); @@ -341,7 +368,7 @@ void build_cells (Box const& bx, Array4 const& cell, { amrex::ignore_unused(k); bool is_small = set_eb_cell(i, j, cell, apx, apy, dx, vfrac, vcent, barea, bcent, - bnorm, small_volfrac); + bnorm, levset, small_volfrac); if (is_small) { Gpu::Atomic::Add(dp, 1); }