Skip to content

Commit

Permalink
Fix assertions involving IntVect (#3919)
Browse files Browse the repository at this point in the history
A MultiFab used to have the same number of ghost cells in all
directions. This has no longer been the case for many years. But many
assertions have not been updated for this change. This commit fixes
various assertions to make it more appropriate.

Also added are overloaded functions for element-wise comparison between
an IntVect and an integer. This makes the code cleaner.
  • Loading branch information
WeiqunZhang committed May 6, 2024
1 parent 07352f7 commit a72fb51
Show file tree
Hide file tree
Showing 22 changed files with 109 additions and 79 deletions.
6 changes: 3 additions & 3 deletions Src/Amr/AMReX_AmrLevel.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -1531,7 +1531,7 @@ AmrLevel::FillCoarsePatch (MultiFab& mf,
//
BL_ASSERT(level != 0);
BL_ASSERT(ncomp <= (mf.nComp()-dcomp));
BL_ASSERT(nghost <= mf.nGrow());
BL_ASSERT(mf.nGrowVect().allGE(nghost));
BL_ASSERT(0 <= idx && idx < desc_lst.size());

int DComp = dcomp;
Expand Down Expand Up @@ -2186,7 +2186,7 @@ AmrLevel::FillPatch (AmrLevel& amrlevel,
{
BL_PROFILE("AmrLevel::FillPatch()");
BL_ASSERT(dcomp+ncomp-1 <= leveldata.nComp());
BL_ASSERT(boxGrow <= leveldata.nGrow());
BL_ASSERT(leveldata.nGrowVect().allGE(boxGrow));
FillPatchIterator fpi(amrlevel, leveldata, boxGrow, time, index, scomp, ncomp);
const MultiFab& mf_fillpatched = fpi.get_mf();
MultiFab::Copy(leveldata, mf_fillpatched, 0, dcomp, ncomp, boxGrow);
Expand All @@ -2204,7 +2204,7 @@ AmrLevel::FillPatchAdd (AmrLevel& amrlevel,
{
BL_PROFILE("AmrLevel::FillPatchAdd()");
BL_ASSERT(dcomp+ncomp-1 <= leveldata.nComp());
BL_ASSERT(boxGrow <= leveldata.nGrow());
BL_ASSERT(leveldata.nGrowVect().allGE(boxGrow));
FillPatchIterator fpi(amrlevel, leveldata, boxGrow, time, index, scomp, ncomp);
const MultiFab& mf_fillpatched = fpi.get_mf();
MultiFab::Add(leveldata, mf_fillpatched, 0, dcomp, ncomp, boxGrow);
Expand Down
2 changes: 1 addition & 1 deletion Src/Amr/AMReX_Extrapolater.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -17,7 +17,7 @@ namespace amrex::Extrapolater

void FirstOrderExtrap (MultiFab& mf, const Geometry& geom, int scomp, int ncomp, int ngrow)
{
BL_ASSERT(mf.nGrow() >= ngrow);
BL_ASSERT(mf.nGrowVect().allGE(ngrow));
BL_ASSERT(scomp >= 0);
BL_ASSERT((scomp+ncomp) <= mf.nComp());

Expand Down
8 changes: 4 additions & 4 deletions Src/Amr/AMReX_StateData.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -140,10 +140,10 @@ StateData::copyOld (const StateData& state)
const MultiFab& MF = state.oldData();

int nc = MF.nComp();
int ng = MF.nGrow();
auto ng = MF.nGrowVect();

BL_ASSERT(nc == (*old_data).nComp());
BL_ASSERT(ng == (*old_data).nGrow());
BL_ASSERT(ng == (*old_data).nGrowVect());

MultiFab::Copy(*old_data, state.oldData(), 0, 0, nc, ng);

Expand All @@ -156,10 +156,10 @@ StateData::copyNew (const StateData& state)
const MultiFab& MF = state.newData();

int nc = MF.nComp();
int ng = MF.nGrow();
auto ng = MF.nGrowVect();

BL_ASSERT(nc == (*new_data).nComp());
BL_ASSERT(ng == (*new_data).nGrow());
BL_ASSERT(ng == (*new_data).nGrowVect());

MultiFab::Copy(*new_data, state.newData(), 0, 0, nc, ng);

Expand Down
2 changes: 1 addition & 1 deletion Src/AmrCore/AMReX_Interpolater.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -1110,7 +1110,7 @@ CellConservativeProtected::protect (const FArrayBox& /*crse*/,
Vector<BCRec>& /*bcr*/,
RunOn runon)
{
AMREX_ALWAYS_ASSERT(ratio.allGT(IntVect(1)));
AMREX_ALWAYS_ASSERT(ratio.allGT(1));

#if (AMREX_SPACEDIM == 1)
amrex::ignore_unused(fine,fine_state,
Expand Down
2 changes: 1 addition & 1 deletion Src/Base/AMReX_BCUtil.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -36,7 +36,7 @@ struct dummy_gpu_fill_extdir
void FillDomainBoundary (MultiFab& phi, const Geometry& geom, const Vector<BCRec>& bc)
{
if (geom.isAllPeriodic()) { return; }
if (phi.nGrow() == 0) { return; }
if (phi.nGrowVect() == 0) { return; }

AMREX_ALWAYS_ASSERT(phi.ixType().cellCentered());

Expand Down
4 changes: 2 additions & 2 deletions Src/Base/AMReX_Box.H
Original file line number Diff line number Diff line change
Expand Up @@ -78,7 +78,7 @@ public:
bigend(big),
btype(typ)
{
BL_ASSERT(typ.allGE(IntVect::TheZeroVector()) && typ.allLE(IntVect::TheUnitVector()));
BL_ASSERT(typ.allGE(0) && typ.allLE(1));
}

//! Construct dimension specific Boxes.
Expand Down Expand Up @@ -829,7 +829,7 @@ AMREX_FORCE_INLINE
Box&
Box::convert (const IntVect& typ) noexcept
{
BL_ASSERT(typ.allGE(IntVect::TheZeroVector()) && typ.allLE(IntVect::TheUnitVector()));
BL_ASSERT(typ.allGE(0) && typ.allLE(1));
IntVect shft(typ - btype.ixType());
bigend += shft;
btype = IndexType(typ);
Expand Down
8 changes: 4 additions & 4 deletions Src/Base/AMReX_FabArray.H
Original file line number Diff line number Diff line change
Expand Up @@ -2035,7 +2035,7 @@ FabArray<FAB>::define (const BoxArray& bxs,

define_function_called = true;

AMREX_ASSERT(ngrow.allGE(IntVect::TheZeroVector()));
AMREX_ASSERT(ngrow.allGE(0));
AMREX_ASSERT(boxarray.empty());
FabArrayBase::define(bxs, dm, nvar, ngrow);

Expand Down Expand Up @@ -2509,7 +2509,7 @@ FabArray<FAB>::setVal (value_type val,
int ncomp,
const IntVect& nghost)
{
AMREX_ASSERT(nghost.allGE(IntVect::TheZeroVector()) && nghost.allLE(n_grow));
AMREX_ASSERT(nghost.allGE(0) && nghost.allLE(n_grow));
AMREX_ALWAYS_ASSERT(comp+ncomp <= n_comp);

BL_PROFILE("FabArray::setVal()");
Expand Down Expand Up @@ -2564,7 +2564,7 @@ FabArray<FAB>::setVal (value_type val,
int ncomp,
const IntVect& nghost)
{
AMREX_ASSERT(nghost.allGE(IntVect::TheZeroVector()) && nghost.allLE(n_grow));
AMREX_ASSERT(nghost.allGE(0) && nghost.allLE(n_grow));
AMREX_ALWAYS_ASSERT(comp+ncomp <= n_comp);

BL_PROFILE("FabArray::setVal(val,region,comp,ncomp,nghost)");
Expand Down Expand Up @@ -2617,7 +2617,7 @@ template <class F, std::enable_if_t<IsBaseFab<F>::value,int>Z>
void
FabArray<FAB>::abs (int comp, int ncomp, const IntVect& nghost)
{
AMREX_ASSERT(nghost.allGE(IntVect::TheZeroVector()) && nghost.allLE(n_grow));
AMREX_ASSERT(nghost.allGE(0) && nghost.allLE(n_grow));
AMREX_ALWAYS_ASSERT(comp+ncomp <= n_comp);
BL_PROFILE("FabArray::abs()");

Expand Down
4 changes: 2 additions & 2 deletions Src/Base/AMReX_FabArrayBase.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -184,7 +184,7 @@ FabArrayBase::define (const BoxArray& bxs,
int nvar,
const IntVect& ngrow)
{
BL_ASSERT(ngrow.allGE(IntVect::TheZeroVector()));
BL_ASSERT(ngrow.allGE(0));
BL_ASSERT(boxarray.empty());
indexArray.clear();
ownership.clear();
Expand Down Expand Up @@ -875,7 +875,7 @@ FabArrayBase::define_fb_metadata (CommMetaData& cmd, const IntVect& nghost,
void
FabArrayBase::FB::define_fb (const FabArrayBase& fa)
{
AMREX_ASSERT(m_multi_ghost ? fa.nGrow() >= 2 : true); // must have >= 2 ghost nodes
AMREX_ASSERT(m_multi_ghost ? fa.nGrowVect().allGE(2) : true); // must have >= 2 ghost nodes
AMREX_ASSERT(m_multi_ghost ? !m_period.isAnyPeriodic() : true); // this only works for non-periodic

fa.define_fb_metadata(*this, m_ngrow, m_cross, m_period, m_multi_ghost);
Expand Down
34 changes: 33 additions & 1 deletion Src/Base/AMReX_IntVect.H
Original file line number Diff line number Diff line change
Expand Up @@ -299,6 +299,14 @@ public:
return AMREX_D_TERM(vect[0] < rhs[0], && vect[1] < rhs[1], && vect[2] < rhs[2]);
}
/**
* \brief Returns true if this is less than argument for all components.
*/
[[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
bool allLT (int rhs) const noexcept
{
return AMREX_D_TERM(vect[0] < rhs, && vect[1] < rhs, && vect[2] < rhs);
}
/**
* \brief Returns true if this is less than or equal to argument for all components.
* NOTE: This is NOT a strict weak ordering usable by STL sorting algorithms.
*/
Expand All @@ -308,6 +316,14 @@ public:
return AMREX_D_TERM(vect[0] <= rhs[0], && vect[1] <= rhs[1], && vect[2] <= rhs[2]);
}
/**
* \brief Returns true if this is less than or equal to argument for all components.
*/
[[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
bool allLE (int rhs) const noexcept
{
return AMREX_D_TERM(vect[0] <= rhs, && vect[1] <= rhs, && vect[2] <= rhs);
}
/**
* \brief Returns true if this is greater than argument for all components.
* NOTE: This is NOT a strict weak ordering usable by STL sorting algorithms.
*/
Expand All @@ -317,6 +333,14 @@ public:
return AMREX_D_TERM(vect[0] > rhs[0], && vect[1] > rhs[1], && vect[2] > rhs[2]);
}
/**
* \brief Returns true if this is greater than argument for all components.
*/
[[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
bool allGT (int rhs) const noexcept
{
return AMREX_D_TERM(vect[0] > rhs, && vect[1] > rhs, && vect[2] > rhs);
}
/**
* \brief Returns true if this is greater than or equal to argument for all components.
* NOTE: This is NOT a strict weak ordering usable by STL sorting algorithms.
*/
Expand All @@ -325,6 +349,14 @@ public:
{
return AMREX_D_TERM(vect[0] >= rhs[0], && vect[1] >= rhs[1], && vect[2] >= rhs[2]);
}
/**
* \brief Returns true if this is greater than or equal to argument for all components.
*/
[[nodiscard]] AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
bool allGE (int rhs) const noexcept
{
return AMREX_D_TERM(vect[0] >= rhs, && vect[1] >= rhs, && vect[2] >= rhs);
}
//! Unary plus -- for completeness.
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
IntVect operator+ () const noexcept { return *this; }
Expand Down Expand Up @@ -580,7 +612,7 @@ AMREX_FORCE_INLINE
IntVect&
IntVect::coarsen (const IntVect& p) noexcept
{
BL_ASSERT(p.allGT(IntVect::TheZeroVector()));
BL_ASSERT(p.allGT(0));
AMREX_D_TERM(vect[0] = amrex::coarsen(vect[0], p.vect[0]);,
vect[1] = amrex::coarsen(vect[1], p.vect[1]);,
vect[2] = amrex::coarsen(vect[2], p.vect[2]););
Expand Down
34 changes: 17 additions & 17 deletions Src/Base/AMReX_MultiFab.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -44,7 +44,7 @@ MultiFab::Dot (const MultiFab& x, int xcomp,
Real
MultiFab::Dot (const MultiFab& x, int xcomp, int numcomp, int nghost, bool local)
{
BL_ASSERT(x.nGrow() >= nghost);
BL_ASSERT(x.nGrowVect().allGE(nghost));

BL_PROFILE("MultiFab::Dot()");

Expand Down Expand Up @@ -97,8 +97,8 @@ MultiFab::Dot (const iMultiFab& mask,
BL_ASSERT(x.boxArray() == mask.boxArray());
BL_ASSERT(x.DistributionMap() == y.DistributionMap());
BL_ASSERT(x.DistributionMap() == mask.DistributionMap());
BL_ASSERT(x.nGrow() >= nghost && y.nGrow() >= nghost);
BL_ASSERT(mask.nGrow() >= nghost);
BL_ASSERT(x.nGrowVect().allGE(nghost) && y.nGrowVect().allGE(nghost));
BL_ASSERT(mask.nGrowVect().allGE(nghost));

Real sm = Real(0.0);
#ifdef AMREX_USE_GPU
Expand Down Expand Up @@ -712,15 +712,15 @@ MultiFab::contains_inf (int scomp, int ncomp, int ngrow, bool local) const
bool
MultiFab::contains_inf (bool local) const
{
return contains_inf(0,nComp(),nGrow(),local);
return contains_inf(0,nComp(),nGrowVect(),local);
}

Real
MultiFab::min (int comp, int nghost, bool local) const
{
BL_PROFILE("MultiFab::min()");

BL_ASSERT(nghost >= 0 && nghost <= n_grow.min());
BL_ASSERT(nghost >= 0 && n_grow.allGE(nghost));

Real mn = std::numeric_limits<Real>::max();

Expand Down Expand Up @@ -801,7 +801,7 @@ MultiFab::min (int comp, int nghost, bool local) const
Real
MultiFab::min (const Box& region, int comp, int nghost, bool local) const
{
BL_ASSERT(nghost >= 0 && nghost <= n_grow.min());
BL_ASSERT(nghost >= 0 && n_grow.allGE(nghost));

BL_PROFILE("MultiFab::min(region)");

Expand Down Expand Up @@ -847,7 +847,7 @@ MultiFab::min (const Box& region, int comp, int nghost, bool local) const
Real
MultiFab::max (int comp, int nghost, bool local) const
{
BL_ASSERT(nghost >= 0 && nghost <= n_grow.min());
BL_ASSERT(nghost >= 0 && n_grow.allGE(nghost));

BL_PROFILE("MultiFab::max()");

Expand Down Expand Up @@ -1006,15 +1006,15 @@ indexFromValue (MultiFab const& mf, int comp, int nghost, Real value, MPI_Op mml
IntVect
MultiFab::minIndex (int comp, int nghost) const
{
BL_ASSERT(nghost >= 0 && nghost <= n_grow.min());
BL_ASSERT(nghost >= 0 && n_grow.allGE(nghost));
Real mn = this->min(comp, nghost, true);
return indexFromValue(*this, comp, nghost, mn, MPI_MINLOC);
}

IntVect
MultiFab::maxIndex (int comp, int nghost) const
{
BL_ASSERT(nghost >= 0 && nghost <= n_grow.min());
BL_ASSERT(nghost >= 0 && n_grow.allGE(nghost));
Real mx = this->max(comp, nghost, true);
return indexFromValue(*this, comp, nghost, mx, MPI_MAXLOC);
}
Expand Down Expand Up @@ -1374,7 +1374,7 @@ void
MultiFab::plus (Real val, int comp, int num_comp, int nghost)
{

BL_ASSERT(nghost >= 0 && nghost <= n_grow.min());
BL_ASSERT(nghost >= 0 && n_grow.allGE(nghost));
BL_ASSERT(comp+num_comp <= n_comp);
BL_ASSERT(num_comp > 0);

Expand All @@ -1384,7 +1384,7 @@ MultiFab::plus (Real val, int comp, int num_comp, int nghost)
void
MultiFab::plus (Real val, const Box& region, int comp, int num_comp, int nghost)
{
BL_ASSERT(nghost >= 0 && nghost <= n_grow.min());
BL_ASSERT(nghost >= 0 && n_grow.allGE(nghost));
BL_ASSERT(comp+num_comp <= n_comp);
BL_ASSERT(num_comp > 0);

Expand All @@ -1400,7 +1400,7 @@ MultiFab::plus (const MultiFab& mf, int strt_comp, int num_comp, int nghost)
void
MultiFab::mult (Real val, int comp, int num_comp, int nghost)
{
BL_ASSERT(nghost >= 0 && nghost <= n_grow.min());
BL_ASSERT(nghost >= 0 && n_grow.allGE(nghost));
BL_ASSERT(comp+num_comp <= n_comp);
BL_ASSERT(num_comp > 0);

Expand All @@ -1410,7 +1410,7 @@ MultiFab::mult (Real val, int comp, int num_comp, int nghost)
void
MultiFab::mult (Real val, const Box& region, int comp, int num_comp, int nghost)
{
BL_ASSERT(nghost >= 0 && nghost <= n_grow.min());
BL_ASSERT(nghost >= 0 && n_grow.allGE(nghost));
BL_ASSERT(comp+num_comp <= n_comp);
BL_ASSERT(num_comp > 0);

Expand All @@ -1420,7 +1420,7 @@ MultiFab::mult (Real val, const Box& region, int comp, int num_comp, int nghost)
void
MultiFab::invert (Real numerator, int comp, int num_comp, int nghost)
{
BL_ASSERT(nghost >= 0 && nghost <= n_grow.min());
BL_ASSERT(nghost >= 0 && n_grow.allGE(nghost));
BL_ASSERT(comp+num_comp <= n_comp);
BL_ASSERT(num_comp > 0);

Expand All @@ -1430,7 +1430,7 @@ MultiFab::invert (Real numerator, int comp, int num_comp, int nghost)
void
MultiFab::invert (Real numerator, const Box& region, int comp, int num_comp, int nghost)
{
BL_ASSERT(nghost >= 0 && nghost <= n_grow.min());
BL_ASSERT(nghost >= 0 && n_grow.allGE(nghost));
BL_ASSERT(comp+num_comp <= n_comp);
BL_ASSERT(num_comp > 0);

Expand All @@ -1440,7 +1440,7 @@ MultiFab::invert (Real numerator, const Box& region, int comp, int num_comp, int
void
MultiFab::negate (int comp, int num_comp, int nghost)
{
BL_ASSERT(nghost >= 0 && nghost <= n_grow.min());
BL_ASSERT(nghost >= 0 && n_grow.allGE(nghost));
BL_ASSERT(comp+num_comp <= n_comp);

FabArray<FArrayBox>::mult(-1., comp, num_comp, nghost);
Expand All @@ -1449,7 +1449,7 @@ MultiFab::negate (int comp, int num_comp, int nghost)
void
MultiFab::negate (const Box& region, int comp, int num_comp, int nghost)
{
BL_ASSERT(nghost >= 0 && nghost <= n_grow.min());
BL_ASSERT(nghost >= 0 && n_grow.allGE(nghost));
BL_ASSERT(comp+num_comp <= n_comp);

FabArray<FArrayBox>::mult(-1.,region,comp,num_comp,nghost);
Expand Down

0 comments on commit a72fb51

Please sign in to comment.