Skip to content

Commit 0cbbb6d

Browse files
committed
Add: MultiFab::sum_unique
This provides a new method to sum values in a `MultiFab`. For non-cell-centered data, `MultiFab::sum` double counts box boundary values that are owned by multiple boxes. This provides a function that does not double count these and provides a quick way to get only the sum of physically unique values.
1 parent 3f715d2 commit 0cbbb6d

File tree

2 files changed

+54
-0
lines changed

2 files changed

+54
-0
lines changed

Src/Base/AMReX_MultiFab.H

Lines changed: 6 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -232,6 +232,12 @@ public:
232232
*/
233233
Real sum (int comp = 0, bool local = false) const;
234234
/**
235+
* \brief Same as sum with local=false, but for non-cell-centered data, this
236+
* skips non-unique points that are owned by multiple boxes.
237+
*/
238+
Real sum_unique (int comp = 0,
239+
const Periodicity& period = Periodicity::NonPeriodic()) const;
240+
/**
235241
* \brief Adds the scalar value val to the value of each cell in the
236242
* specified subregion of the MultiFab. The subregion consists
237243
* of the num_comp components starting at component comp.

Src/Base/AMReX_MultiFab.cpp

Lines changed: 48 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -5,6 +5,7 @@
55
#include <AMReX_BLProfiler.H>
66
#include <AMReX_iMultiFab.H>
77
#include <AMReX_FabArrayUtility.H>
8+
#include <AMReX_REAL.H>
89

910
#ifdef AMREX_MEM_PROFILING
1011
#include <AMReX_MemProfiler.H>
@@ -1586,6 +1587,53 @@ MultiFab::sum (int comp, bool local) const
15861587
return sm;
15871588
}
15881589

1590+
Real
1591+
MultiFab::sum_unique (int comp,
1592+
const Periodicity& period) const
1593+
{
1594+
BL_PROFILE("MultiFab::sum_unique()");
1595+
1596+
using namespace amrex::literals;
1597+
1598+
// Owner is the grid with the lowest grid number containing the data
1599+
std::unique_ptr<iMultiFab> owner_mask = OwnerMask(period);
1600+
1601+
Real sm = Real(0.0);
1602+
#ifdef AMREX_USE_GPU
1603+
if (Gpu::inLaunchRegion()) {
1604+
auto const& ma = this->const_arrays();
1605+
auto const& msk = owner_mask->const_array();
1606+
sm = ParReduce(TypeList<ReduceOpSum>{}, TypeList<Real>{}, *this, IntVect(0),
1607+
[=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept
1608+
-> GpuTuple<Real>
1609+
{
1610+
return msk(i,j,k,comp) ? ma[box_no](i,j,k,comp) : 0.0_rt;
1611+
});
1612+
} else
1613+
#endif
1614+
{
1615+
#ifdef AMREX_USE_OMP
1616+
#pragma omp parallel if (!system::regtest_reduction) reduction(+:sm)
1617+
#endif
1618+
for (MFIter mfi(*this,true); mfi.isValid(); ++mfi)
1619+
{
1620+
Array4<int const> const& msk = owner_mask->const_array(mfi);
1621+
Box const& bx = mfi.tilebox();
1622+
auto const& a = this->const_array(mfi);
1623+
Real tmp = 0.0_rt;
1624+
AMREX_LOOP_3D(bx, i, j, k,
1625+
{
1626+
tmp += msk(i,j,k,comp) ? a(i,j,k,comp) : 0.0_rt;
1627+
});
1628+
sm += tmp; // Do it this way so that it does not break regression tests.
1629+
}
1630+
}
1631+
1632+
ParallelAllReduce::Sum(sm, ParallelContext::CommunicatorSub());
1633+
1634+
return sm;
1635+
}
1636+
15891637
void
15901638
MultiFab::minus (const MultiFab& mf, int strt_comp, int num_comp, int nghost)
15911639
{

0 commit comments

Comments
 (0)