Skip to content

Commit

Permalink
Fix for small cells (#2781)
Browse files Browse the repository at this point in the history
Move small cells outside domain extent into iterative approach.
  • Loading branch information
drangara authored May 19, 2022
1 parent 994073b commit 0540fed
Show file tree
Hide file tree
Showing 2 changed files with 88 additions and 86 deletions.
51 changes: 26 additions & 25 deletions Src/EB/AMReX_EB2_2D_C.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -322,31 +322,6 @@ void build_cells (Box const& bx, Array4<EBCellFlag> const& cell,
}
});

smc.copyToHost();
nsmallcells += *hp;

if (nsmallcells > 0 || nmulticuts > 0) {
Box const& nbxg1 = amrex::surroundingNodes(bxg1);
AMREX_HOST_DEVICE_FOR_3D(nbxg1, i, j, k,
{
if (levset(i,j,k) < Real(0.0)) {
if (bxg1.contains(i-1,j-1,k)
&& cell(i-1,j-1,k).isCovered()) {
levset(i,j,k) = Real(0.0);
} else if (bxg1.contains(i ,j-1,k)
&& cell(i ,j-1,k).isCovered()) {
levset(i,j,k) = Real(0.0);
} else if (bxg1.contains(i-1,j ,k)
&& cell(i-1,j ,k).isCovered()) {
levset(i,j,k) = Real(0.0);
} else if (bxg1.contains(i ,j ,k)
&& cell(i ,j ,k).isCovered()) {
levset(i,j,k) = Real(0.0);
}
}
});
}

// set cells in the extended region to covered if the
// corresponding cell on the domain face is covered
if(extend_domain_face) {
Expand Down Expand Up @@ -393,12 +368,38 @@ void build_cells (Box const& bx, Array4<EBCellFlag> const& cell,
if( in_extended_domain && (! cell(i,j,k).isCovered())
&& cell(ii,jj,kk).isCovered() )
{
Gpu::Atomic::Add(dp, 1);
set_covered(i,j,cell,vfrac,vcent,barea,bcent,bnorm);
}
});
}
}

smc.copyToHost();
nsmallcells += *hp;

if (nsmallcells > 0 || nmulticuts > 0) {
Box const& nbxg1 = amrex::surroundingNodes(bxg1);
AMREX_HOST_DEVICE_FOR_3D(nbxg1, i, j, k,
{
if (levset(i,j,k) < Real(0.0)) {
if (bxg1.contains(i-1,j-1,k)
&& cell(i-1,j-1,k).isCovered()) {
levset(i,j,k) = Real(0.0);
} else if (bxg1.contains(i ,j-1,k)
&& cell(i ,j-1,k).isCovered()) {
levset(i,j,k) = Real(0.0);
} else if (bxg1.contains(i-1,j ,k)
&& cell(i-1,j ,k).isCovered()) {
levset(i,j,k) = Real(0.0);
} else if (bxg1.contains(i ,j ,k)
&& cell(i ,j ,k).isCovered()) {
levset(i,j,k) = Real(0.0);
}
}
});
}

// Build neighbors. By default, all neighbors are already set.
AMREX_HOST_DEVICE_FOR_3D ( bxg1, i, j, k,
{
Expand Down
123 changes: 62 additions & 61 deletions Src/EB/AMReX_EB2_3D_C.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -787,6 +787,68 @@ void build_cells (Box const& bx, Array4<EBCellFlag> const& cell,
}
});

// set cells in the extended region to covered if the
// corresponding cell on the domain face is covered
if(extend_domain_face) {

Box gdomain = geom.Domain();
for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
if (geom.isPeriodic(idim)) {
gdomain.setSmall(idim, std::min(gdomain.smallEnd(idim), bxg1.smallEnd(idim)));
gdomain.setBig(idim, std::max(gdomain.bigEnd(idim), bxg1.bigEnd(idim)));
}
}

if (! gdomain.contains(bxg1)) {
AMREX_HOST_DEVICE_FOR_3D ( bxg1, i, j, k,
{
const auto & dlo = gdomain.loVect();
const auto & dhi = gdomain.hiVect();

// find the cell(ii,jj,kk) on the corr. domain face
// this would have already been set to correct value
bool in_extended_domain = false;
int ii = i;
int jj = j;
int kk = k;
if(i < dlo[0]) {
in_extended_domain = true;
ii = dlo[0];
}
else if(i > dhi[0]) {
in_extended_domain = true;
ii = dhi[0];
}

if(j < dlo[1]) {
in_extended_domain = true;
jj = dlo[1];
}
else if(j > dhi[1]) {
in_extended_domain = true;
jj = dhi[1];
}

if(k < dlo[2]) {
in_extended_domain = true;
kk = dlo[2];
}
else if(k > dhi[2]) {
in_extended_domain = true;
kk = dhi[2];
}

// set cell in extendable region to covered if necessary
if( in_extended_domain && (! cell(i,j,k).isCovered())
&& cell(ii,jj,kk).isCovered() )
{
Gpu::Atomic::Add(dp, 1);
set_covered(i,j,k,cell,vfrac,vcent,barea,bcent,bnorm);
}
});
}
}

n_smallcell_multicuts.copyToHost();
nsmallcells += hp[0];
nmulticuts += hp[1];
Expand Down Expand Up @@ -874,67 +936,6 @@ void build_cells (Box const& bx, Array4<EBCellFlag> const& cell,
return;
}

// set cells in the extended region to covered if the
// corresponding cell on the domain face is covered
if(extend_domain_face) {

Box gdomain = geom.Domain();
for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
if (geom.isPeriodic(idim)) {
gdomain.setSmall(idim, std::min(gdomain.smallEnd(idim), bxg1.smallEnd(idim)));
gdomain.setBig(idim, std::max(gdomain.bigEnd(idim), bxg1.bigEnd(idim)));
}
}

if (! gdomain.contains(bxg1)) {
AMREX_HOST_DEVICE_FOR_3D ( bxg1, i, j, k,
{
const auto & dlo = gdomain.loVect();
const auto & dhi = gdomain.hiVect();

// find the cell(ii,jj,kk) on the corr. domain face
// this would have already been set to correct value
bool in_extended_domain = false;
int ii = i;
int jj = j;
int kk = k;
if(i < dlo[0]) {
in_extended_domain = true;
ii = dlo[0];
}
else if(i > dhi[0]) {
in_extended_domain = true;
ii = dhi[0];
}

if(j < dlo[1]) {
in_extended_domain = true;
jj = dlo[1];
}
else if(j > dhi[1]) {
in_extended_domain = true;
jj = dhi[1];
}

if(k < dlo[2]) {
in_extended_domain = true;
kk = dlo[2];
}
else if(k > dhi[2]) {
in_extended_domain = true;
kk = dhi[2];
}

// set cell in extendable region to covered if necessary
if( in_extended_domain && (! cell(i,j,k).isCovered())
&& cell(ii,jj,kk).isCovered() )
{
set_covered(i,j,k,cell,vfrac,vcent,barea,bcent,bnorm);
}
});
}
}

// Build neighbors. By default all 26 neighbors are already set.
AMREX_HOST_DEVICE_FOR_3D ( bxg1, i, j, k,
{
Expand Down

0 comments on commit 0540fed

Please sign in to comment.