Skip to content

Commit

Permalink
Merge branch 'development'
Browse files Browse the repository at this point in the history
  • Loading branch information
WeiqunZhang committed Mar 2, 2019
2 parents 749721f + 6e1acca commit 91323f2
Show file tree
Hide file tree
Showing 7 changed files with 44 additions and 21 deletions.
7 changes: 7 additions & 0 deletions Src/Base/AMReX_BLassert.H
Expand Up @@ -48,4 +48,11 @@
#define AMREX_ALWAYS_ASSERT_WITH_MESSAGE(EX,MSG) (EX)?((void)0):amrex::Assert( # EX , __FILE__, __LINE__ , # MSG)
#define AMREX_ALWAYS_ASSERT(EX) (EX)?((void)0):amrex::Assert( # EX , __FILE__, __LINE__)


#if !defined(AMREX_USE_GPU) || !defined(AMREX_DEBUG) && !defined(AMREX_USE_ASSERTION)
#define AMREX_GPU_ASSERT(EX) ((void)0)
#else
#define AMREX_GPU_ASSERT(EX) (EX)?((void)0):amrex::Assert( # EX , __FILE__, __LINE__)
#endif

#endif /*BL_BL_ASSERT_H*/
11 changes: 10 additions & 1 deletion Src/Base/AMReX_CArena.H
Expand Up @@ -5,6 +5,8 @@
#include <set>
#include <vector>
#include <mutex>
#include <unordered_set>
#include <functional>

#include <AMReX_Arena.H>

Expand Down Expand Up @@ -89,6 +91,12 @@ protected:
return m_owner == rhs.m_owner;
}

struct hash {
std::size_t operator() (const Node& n) const noexcept {
return std::hash<void*>{}(n.m_block);
}
};

private:
//! The block of memory we reference.
void* m_block;
Expand Down Expand Up @@ -117,7 +125,8 @@ protected:
* \brief The list of busy blocks.
* A block is either on the freelist or on the blocklist, but not on both.
*/
NL m_busylist;
// NL m_busylist;
std::unordered_set<Node, Node::hash> m_busylist;
//! The minimal size of hunks to request via ::operator new().
std::size_t m_hunk;
//! The amount of heap space currently allocated.
Expand Down
2 changes: 1 addition & 1 deletion Src/Base/AMReX_CArena.cpp
Expand Up @@ -116,7 +116,7 @@ CArena::free (void* vp)
//
// `vp' had better be in the busy list.
//
NL::iterator busy_it = m_busylist.find(Node(vp,0,0));
auto busy_it = m_busylist.find(Node(vp,0,0));

BL_ASSERT(!(busy_it == m_busylist.end()));
BL_ASSERT(m_freelist.find(*busy_it) == m_freelist.end());
Expand Down
23 changes: 17 additions & 6 deletions Src/Base/AMReX_FabArray.H
Expand Up @@ -782,7 +782,9 @@ FabArray<FAB>::fabDevicePtr (const MFIter& mfi)
{
BL_ASSERT(mfi.LocalIndex() < indexArray.size());
BL_ASSERT(DistributionMap() == mfi.DistributionMap());
return m_fabs_v[mfi.LocalIndex()];
int li = mfi.LocalIndex();
AMREX_GPU_ASSERT(m_fabs_v[li]->dataPtr() == m_host_fabs_v[li]->dataPtr());
return m_fabs_v[li];
}

template <class FAB>
Expand All @@ -791,7 +793,9 @@ FabArray<FAB>::fabDevicePtr (const MFIter& mfi) const
{
BL_ASSERT(mfi.LocalIndex() < indexArray.size());
BL_ASSERT(DistributionMap() == mfi.DistributionMap());
return m_fabs_v[mfi.LocalIndex()];
int li = mfi.LocalIndex();
AMREX_GPU_ASSERT(m_fabs_v[li]->dataPtr() == m_host_fabs_v[li]->dataPtr());
return m_fabs_v[li];
}

template <class FAB>
Expand All @@ -800,6 +804,7 @@ FabArray<FAB>::fabDevicePtr (int K)
{
int li = localindex(K);
BL_ASSERT(li >=0 && li < indexArray.size());
AMREX_GPU_ASSERT(m_fabs_v[li]->dataPtr() == m_host_fabs_v[li]->dataPtr());
return m_fabs_v[li];
}

Expand All @@ -818,10 +823,12 @@ FabArray<FAB>::fabHostPtr (const MFIter& mfi)
{
BL_ASSERT(mfi.LocalIndex() < indexArray.size());
BL_ASSERT(DistributionMap() == mfi.DistributionMap());
int li = mfi.LocalIndex();
AMREX_GPU_ASSERT(m_fabs_v[li]->dataPtr() == m_host_fabs_v[li]->dataPtr());
#if AMREX_USE_GPU
return m_host_fabs_v[mfi.LocalIndex()];
return m_host_fabs_v[li];
#else
return m_fabs_v[mfi.LocalIndex()];
return m_fabs_v[li];
#endif
}

Expand All @@ -831,10 +838,12 @@ FabArray<FAB>::fabHostPtr (const MFIter& mfi) const
{
BL_ASSERT(mfi.LocalIndex() < indexArray.size());
BL_ASSERT(DistributionMap() == mfi.DistributionMap());
int li = mfi.LocalIndex();
AMREX_GPU_ASSERT(m_fabs_v[li]->dataPtr() == m_host_fabs_v[li]->dataPtr());
#if AMREX_USE_GPU
return m_host_fabs_v[mfi.LocalIndex()];
return m_host_fabs_v[li];
#else
return m_fabs_v[mfi.LocalIndex()];
return m_fabs_v[li];
#endif
}

Expand All @@ -844,6 +853,7 @@ FabArray<FAB>::fabHostPtr (int K)
{
int li = localindex(K);
BL_ASSERT(li >=0 && li < indexArray.size());
AMREX_GPU_ASSERT(m_fabs_v[li]->dataPtr() == m_host_fabs_v[li]->dataPtr());
#if AMREX_USE_GPU
return m_host_fabs_v[li];
#else
Expand All @@ -857,6 +867,7 @@ FabArray<FAB>::fabHostPtr (int K) const
{
int li = localindex(K);
BL_ASSERT(li >=0 && li < indexArray.size());
AMREX_GPU_ASSERT(m_fabs_v[li]->dataPtr() == m_host_fabs_v[li]->dataPtr());
#if AMREX_USE_GPU
return m_host_fabs_v[li];
#else
Expand Down
8 changes: 2 additions & 6 deletions Src/Base/AMReX_FabArrayCommI.H
Expand Up @@ -301,8 +301,6 @@ FabArray<FAB>::FBEP_nowait (int scomp, int ncomp, const IntVect& nghost,

if (Gpu::inLaunchRegion())
{
LayoutData<Array<Vector<FabCopyTag<FAB> >,AMREX_SPACEDIM> >
dtags(loc_copy_tags.boxArray(), loc_copy_tags.DistributionMap());
for (MFIter mfi(*this); mfi.isValid(); ++mfi)
{
FAB* dfab = this->fabPtr(mfi);
Expand All @@ -318,7 +316,7 @@ FabArray<FAB>::FBEP_nowait (int scomp, int ncomp, const IntVect& nghost,
gbx_lo.grow(jdim,nghost[jdim]);
}
const int ncells_to_hi = vbx.length(idim)+nghost[idim];
const Box gbx_hi = gbx_lo + IntVect{ncells_to_hi};
const Box gbx_hi = amrex::shift(gbx_lo, idim, ncells_to_hi);

Vector<FabCopyTag<FAB> > tags_lo, tags_hi;
for (const auto& lc_tag : lc_tags) {
Expand Down Expand Up @@ -1503,8 +1501,6 @@ FillBoundary (Vector<FabArray<FAB>*> const& mf, const Periodicity& period)

if (Gpu::inLaunchRegion())
{
LayoutData<Array<Vector<FabCopyTag<FAB> >,AMREX_SPACEDIM> >
dtags(loc_copy_tags.boxArray(), loc_copy_tags.DistributionMap());
for (MFIter mfi(fa); mfi.isValid(); ++mfi)
{
FAB* dfab = fa.fabPtr(mfi);
Expand All @@ -1520,7 +1516,7 @@ FillBoundary (Vector<FabArray<FAB>*> const& mf, const Periodicity& period)
gbx_lo.grow(jdim,ng[jdim]);
}
const int ncells_to_hi = vbx.length(idim)+ng[idim];
const Box gbx_hi = gbx_lo + IntVect{ncells_to_hi};
const Box gbx_hi = amrex::shift(gbx_lo, idim, ncells_to_hi);

Vector<FabCopyTag<FAB> > vtags_lo, vtags_hi;
for (const auto& lc_tag : lc_tags) {
Expand Down
10 changes: 5 additions & 5 deletions Src/LinearSolvers/MLMG/AMReX_MLEBABecLap.cpp
Expand Up @@ -603,7 +603,7 @@ MLEBABecLap::FFlux (int amrlev, const MFIter& mfi, const Array<FArrayBox*,AMREX_
const FArrayBox& sol, Location loc, const int face_only) const
{
// todo: gpu
Gpu::LaunchSafeGuard(false);
Gpu::LaunchSafeGuard lg(false);
BL_PROFILE("MLEBABecLap::FFlux()");
const int at_centroid = (Location::FaceCentroid == loc) ? 1 : 0;
const int mglev = 0;
Expand Down Expand Up @@ -665,7 +665,7 @@ MLEBABecLap::compGrad (int amrlev, const Array<MultiFab*,AMREX_SPACEDIM>& grad,
MultiFab& sol, Location loc) const
{
// todo: gpu
Gpu::LaunchSafeGuard(false);
Gpu::LaunchSafeGuard lg(false);
BL_PROFILE("MLEBABecLap::compGrad()");

const int at_centroid = (Location::FaceCentroid == loc) ? 1 : 0;
Expand Down Expand Up @@ -753,7 +753,7 @@ void
MLEBABecLap::normalize (int amrlev, int mglev, MultiFab& mf) const
{
// todo: gpu
Gpu::LaunchSafeGuard(false);
Gpu::LaunchSafeGuard lg(false);
const MultiFab& acoef = m_a_coeffs[amrlev][mglev];
AMREX_D_TERM(const MultiFab& bxcoef = m_b_coeffs[amrlev][mglev][0];,
const MultiFab& bycoef = m_b_coeffs[amrlev][mglev][1];,
Expand Down Expand Up @@ -840,7 +840,7 @@ void
MLEBABecLap::interpolation (int amrlev, int fmglev, MultiFab& fine, const MultiFab& crse) const
{
// todo: gpu
Gpu::LaunchSafeGuard(false);
Gpu::LaunchSafeGuard lg(false);
auto factory = dynamic_cast<EBFArrayBoxFactory const*>(m_factory[amrlev][fmglev].get());
const FabArray<EBCellFlagFab>* flags = (factory) ? &(factory->getMultiEBCellFlagFab()) : nullptr;

Expand Down Expand Up @@ -891,7 +891,7 @@ MLEBABecLap::applyBC (int amrlev, int mglev, MultiFab& in, BCMode bc_mode, State
const MLMGBndry* bndry, bool skip_fillboundary) const
{
// todo: gpu
Gpu::LaunchSafeGuard(false);
Gpu::LaunchSafeGuard lg(false);
BL_PROFILE("MLEBABecLap::applyBC()");

// No coarsened boundary values, cannot apply inhomog at mglev>0.
Expand Down
4 changes: 2 additions & 2 deletions Src/LinearSolvers/MLMG/AMReX_MLMG.cpp
Expand Up @@ -561,7 +561,7 @@ MLMG::interpCorrection (int alev)

if (linop.isCellCentered())
{
Gpu::LaunchSafeGuard(!isEB); // turn off gpu for eb for now TODO
Gpu::LaunchSafeGuard lg(!isEB); // turn off gpu for eb for now TODO
MFItInfo mfi_info;
if (Gpu::notInLaunchRegion()) mfi_info.EnableTiling().SetDynamic(true);
#ifdef _OPENMP
Expand Down Expand Up @@ -692,7 +692,7 @@ MLMG::interpCorrection (int alev, int mglev)

if (linop.isCellCentered())
{
Gpu::LaunchSafeGuard(!isEB); // turn off gpu for eb for now TODO
Gpu::LaunchSafeGuard lg(!isEB); // turn off gpu for eb for now TODO
MFItInfo mfi_info;
if (Gpu::notInLaunchRegion()) mfi_info.EnableTiling().SetDynamic(true);
#ifdef _OPENMP
Expand Down

0 comments on commit 91323f2

Please sign in to comment.