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

Add option to derefine to AMRErrorTag #2875

Merged
merged 1 commit into from
Jul 18, 2022
Merged
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
5 changes: 5 additions & 0 deletions Src/AmrCore/AMReX_ErrorList.H
Original file line number Diff line number Diff line change
Expand Up @@ -383,6 +383,7 @@ std::ostream& operator << (std::ostream& os, const ErrorList& elst);
Real m_min_time = std::numeric_limits<Real>::lowest();
Real m_max_time = std::numeric_limits<Real>::max();
int m_volume_weighting = 0;
int m_derefine = 0;
RealBox m_realbox;

AMRErrorTagInfo& SetMaxLevel (int max_level) noexcept {
Expand All @@ -405,6 +406,10 @@ std::ostream& operator << (std::ostream& os, const ErrorList& elst);
m_volume_weighting = volume_weighting;
return *this;
}
AMRErrorTagInfo& SetDerefine (int derefine) noexcept {
m_derefine = derefine;
return *this;
}
};

class AMRErrorTag
Expand Down
19 changes: 12 additions & 7 deletions Src/AmrCore/AMReX_ErrorList.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -293,6 +293,11 @@ AMRErrorTag::operator() (TagBoxArray& tba,
auto threshold = m_value[level];
auto const volume_weighting = m_info.m_volume_weighting;
auto geomdata = geom.data();
auto tag_update = tagval;
if (m_info.m_derefine) {
tag_update = clearval;
}

if (m_test == GRAD)
{
ParallelFor(tba, [=] AMREX_GPU_DEVICE (int bi, int i, int j, int k) noexcept
Expand All @@ -301,7 +306,7 @@ AMRErrorTag::operator() (TagBoxArray& tba,
auto ax = amrex::Math::abs(dat(i+1,j,k) - dat(i,j,k));
ax = amrex::max(ax,amrex::Math::abs(dat(i,j,k) - dat(i-1,j,k)));
#if AMREX_SPACEDIM == 1
if (ax >= threshold) { tagma[bi](i,j,k) = tagval;}
if (ax >= threshold) { tagma[bi](i,j,k) = tag_update;}
#else
auto ay = amrex::Math::abs(dat(i,j+1,k) - dat(i,j,k));
ay = amrex::max(ay,amrex::Math::abs(dat(i,j,k) - dat(i,j-1,k)));
Expand All @@ -310,7 +315,7 @@ AMRErrorTag::operator() (TagBoxArray& tba,
az = amrex::max(az,amrex::Math::abs(dat(i,j,k) - dat(i,j,k-1)));
#endif
if (amrex::max(AMREX_D_DECL(ax,ay,az)) >= threshold) {
tagma[bi](i,j,k) = tagval;
tagma[bi](i,j,k) = tag_update;
}
#endif
});
Expand All @@ -323,7 +328,7 @@ AMRErrorTag::operator() (TagBoxArray& tba,
auto ax = amrex::Math::abs(dat(i+1,j,k) - dat(i,j,k));
ax = amrex::max(ax,amrex::Math::abs(dat(i,j,k) - dat(i-1,j,k)));
#if AMREX_SPACEDIM == 1
if (ax >= threshold * amrex::Math::abs(dat(i,j,k))) { tagma[bi](i,j,k) = tagval;}
if (ax >= threshold * amrex::Math::abs(dat(i,j,k))) { tagma[bi](i,j,k) = tag_update;}
#else
auto ay = amrex::Math::abs(dat(i,j+1,k) - dat(i,j,k));
ay = amrex::max(ay,amrex::Math::abs(dat(i,j,k) - dat(i,j-1,k)));
Expand All @@ -333,7 +338,7 @@ AMRErrorTag::operator() (TagBoxArray& tba,
#endif
if (amrex::max(AMREX_D_DECL(ax,ay,az))
>= threshold * amrex::Math::abs(dat(i,j,k))) {
tagma[bi](i,j,k) = tagval;
tagma[bi](i,j,k) = tag_update;
}
#endif
});
Expand All @@ -344,7 +349,7 @@ AMRErrorTag::operator() (TagBoxArray& tba,
{
Real vol = volume_weighting ? Geometry::Volume(IntVect{AMREX_D_DECL(i,j,k)}, geomdata) : 1.0_rt;
if (datma[bi](i,j,k) * vol <= threshold) {
tagma[bi](i,j,k) = tagval;
tagma[bi](i,j,k) = tag_update;
}
});
}
Expand All @@ -354,7 +359,7 @@ AMRErrorTag::operator() (TagBoxArray& tba,
{
Real vol = volume_weighting ? Geometry::Volume(IntVect{AMREX_D_DECL(i,j,k)}, geomdata) : 1.0_rt;
if (datma[bi](i,j,k) * vol >= threshold) {
tagma[bi](i,j,k) = tagval;
tagma[bi](i,j,k) = tag_update;
}
});
}
Expand All @@ -364,7 +369,7 @@ AMRErrorTag::operator() (TagBoxArray& tba,
ParallelFor(tba, [=] AMREX_GPU_DEVICE (int bi, int i, int j, int k) noexcept
{
if (datma[bi](i,j,k) >= fac) {
tagma[bi](i,j,k) = tagval;
tagma[bi](i,j,k) = tag_update;
}
});
}
Expand Down