Skip to content

Commit

Permalink
Support for multiple periods in FillBoundary and ParallelCopy (#4106)
Browse files Browse the repository at this point in the history
This is a feature requested by a user who wants to run simulations with
the number of ghost cells more than the domain length in periodic
directions.
  • Loading branch information
WeiqunZhang committed Aug 22, 2024
1 parent b05f69f commit 12002e7
Show file tree
Hide file tree
Showing 8 changed files with 119 additions and 8 deletions.
10 changes: 5 additions & 5 deletions Src/Base/AMReX_FabArrayBase.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -345,7 +345,7 @@ FabArrayBase::CPC::define (const BoxArray& ba_dst, const DistributionMapping& dm

std::vector< std::pair<int,Box> > isects;

const std::vector<IntVect>& pshifts = m_period.shiftIntVect();
const std::vector<IntVect>& pshifts = m_period.shiftIntVect(ng_dst);

auto& send_tags = *m_SndTags;

Expand Down Expand Up @@ -672,7 +672,7 @@ FabArrayBase::define_fb_metadata (CommMetaData& cmd, const IntVect& nghost,
const IntVect ng_ng =nghost - 1;
std::vector< std::pair<int,Box> > isects;

const std::vector<IntVect>& pshifts = period.shiftIntVect();
const std::vector<IntVect>& pshifts = period.shiftIntVect(nghost);

auto& send_tags = *cmd.m_SndTags;

Expand Down Expand Up @@ -901,7 +901,7 @@ FabArrayBase::FB::define_epo (const FabArrayBase& fa)
const IndexType& typ = ba.ixType();
std::vector< std::pair<int,Box> > isects;

const std::vector<IntVect>& pshifts = m_period.shiftIntVect();
const std::vector<IntVect>& pshifts = m_period.shiftIntVect(ng);

auto& send_tags = *m_SndTags;

Expand Down Expand Up @@ -1053,7 +1053,7 @@ void FabArrayBase::FB::tag_one_box (int krcv, BoxArray const& ba, DistributionMa

std::vector<std::pair<int,Box> > isects2;
std::vector<std::tuple<int,Box,IntVect> > isects3;
auto const& pshifts = m_period.shiftIntVect();
auto const& pshifts = m_period.shiftIntVect(m_ngrow);
for (auto const& shft: pshifts) {
ba.intersections(gbx+shft, isects2);
for (auto const& is2 : isects2) {
Expand Down Expand Up @@ -1144,7 +1144,7 @@ FabArrayBase::FB::define_os (const FabArrayBase& fa)

#ifdef AMREX_USE_MPI
if (ParallelDescriptor::NProcs() > 1) {
const std::vector<IntVect>& pshifts = m_period.shiftIntVect();
const std::vector<IntVect>& pshifts = m_period.shiftIntVect(m_ngrow);
std::vector< std::pair<int,Box> > isects;

std::set<int> my_receiver;
Expand Down
2 changes: 1 addition & 1 deletion Src/Base/AMReX_Periodicity.H
Original file line number Diff line number Diff line change
Expand Up @@ -32,7 +32,7 @@ public:
//! Cell-centered domain Box "infinitely" long in non-periodic directions.
[[nodiscard]] Box Domain () const noexcept;

[[nodiscard]] std::vector<IntVect> shiftIntVect () const;
[[nodiscard]] std::vector<IntVect> shiftIntVect (IntVect const& nghost = IntVect(0)) const;

static const Periodicity& NonPeriodic () noexcept;

Expand Down
5 changes: 4 additions & 1 deletion Src/Base/AMReX_Periodicity.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -5,7 +5,7 @@
namespace amrex {

std::vector<IntVect>
Periodicity::shiftIntVect () const
Periodicity::shiftIntVect (IntVect const& nghost) const
{
std::vector<IntVect> r;

Expand All @@ -15,6 +15,9 @@ Periodicity::shiftIntVect () const
for (int i = 0; i < AMREX_SPACEDIM; ++i) {
if (isPeriodic(i)) {
per[i] = jmp[i] = period[i];
while (per[i] < nghost[i]) {
per[i] += period[i];
}
}
}

Expand Down
3 changes: 2 additions & 1 deletion Tests/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -122,7 +122,8 @@ else()
# List of subdirectories to search for CMakeLists.
#
set( AMREX_TESTS_SUBDIRS Amr AsyncOut CLZ CTOParFor DeviceGlobal Enum
MultiBlock Parser Parser2 Reinit RoundoffDomain)
MultiBlock MultiPeriod Parser Parser2 Reinit
RoundoffDomain)

if (AMReX_PARTICLES)
list(APPEND AMREX_TESTS_SUBDIRS Particles)
Expand Down
9 changes: 9 additions & 0 deletions Tests/MultiPeriod/CMakeLists.txt
Original file line number Diff line number Diff line change
@@ -0,0 +1,9 @@
foreach(D IN LISTS AMReX_SPACEDIM)
set(_sources main.cpp)
set(_input_files )

setup_test(${D} _sources _input_files)

unset(_sources)
unset(_input_files)
endforeach()
24 changes: 24 additions & 0 deletions Tests/MultiPeriod/GNUmakefile
Original file line number Diff line number Diff line change
@@ -0,0 +1,24 @@
AMREX_HOME := ../..

DEBUG = FALSE

DIM = 3

COMP = gcc

USE_MPI = FALSE
USE_OMP = FALSE
USE_CUDA = FALSE
USE_HIP = FALSE
USE_SYCL = FALSE

BL_NO_FORT = TRUE

TINY_PROFILE = FALSE

include $(AMREX_HOME)/Tools/GNUMake/Make.defs

include ./Make.package
include $(AMREX_HOME)/Src/Base/Make.package

include $(AMREX_HOME)/Tools/GNUMake/Make.rules
1 change: 1 addition & 0 deletions Tests/MultiPeriod/Make.package
Original file line number Diff line number Diff line change
@@ -0,0 +1 @@
CEXE_sources += main.cpp
73 changes: 73 additions & 0 deletions Tests/MultiPeriod/main.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,73 @@
#include <AMReX.H>
#include <AMReX_Geometry.H>
#include <AMReX_FabArray.H>
#include <AMReX_MFParallelFor.H>
#include <AMReX_ParReduce.H>

using namespace amrex;

int main (int argc, char* argv[])
{
amrex::Initialize(argc, argv);
{
// Domain size: 2 x 128 x 4
Box box(IntVect(0), IntVect(AMREX_D_DECL(1, 127, 3)));
Array<int,AMREX_SPACEDIM> is_periodic{AMREX_D_DECL(1,1,1)};
Geometry geom(box, RealBox(AMREX_D_DECL(Real(0),Real(0),Real(0)),
AMREX_D_DECL(Real(1),Real(1),Real(1))),
CoordSys::cartesian, is_periodic);
BoxArray ba(box);
ba.maxSize(32);
ba.convert(IntVect(AMREX_D_DECL(1,0,0))); // nodal in x-direction
DistributionMapping dm(ba);

FabArray<BaseFab<Long>> mf1(ba,dm,1,IntVect(4));
FabArray<BaseFab<Long>> mf2(ba,dm,1,IntVect(5));

mf1.setVal(-1);
mf2.setVal(-2);

auto const& len = geom.Domain().length3d();
auto expected = [=] AMREX_GPU_DEVICE (int i, int j, int k)
{
while (i < 0 ) { i += len[0]; }
while (i >= len[0]) { i -= len[0]; }
while (j < 0 ) { j += len[1]; }
while (j >= len[1]) { j -= len[1]; }
while (k < 0 ) { k += len[2]; }
while (k >= len[2]) { k -= len[2]; }
return Long(i) + Long(j)*Long(len[0]) + Long(k)*Long(len[0])*Long(len[1]);
};

auto const& ma1 = mf1.arrays();
auto const& ma2 = mf2.arrays();

// Initialize valid region
ParallelFor(mf1, IntVect(0), [=] AMREX_GPU_DEVICE (int b, int i, int j, int k)
{
ma1[b](i,j,k) = expected(i,j,k);
});

mf1.FillBoundary(geom.periodicity());
mf2.ParallelCopy(mf1, 0, 0, 1, IntVect(0), mf2.nGrowVect(), geom.periodicity());

auto r1 = ParReduce(TypeList<ReduceOpSum>{}, TypeList<Long>{}, mf1, mf1.nGrowVect(),
[=] AMREX_GPU_DEVICE (int b, int i, int j, int k) -> GpuTuple<Long>
{
return { Long(expected(i,j,k) != ma1[b](i,j,k)) };
});
auto r2 = ParReduce(TypeList<ReduceOpSum>{}, TypeList<Long>{}, mf2, mf2.nGrowVect(),
[=] AMREX_GPU_DEVICE (int b, int i, int j, int k) -> GpuTuple<Long>
{
return { Long(expected(i,j,k) != ma2[b](i,j,k)) };
});

AMREX_ALWAYS_ASSERT(r1 == 0);
AMREX_ALWAYS_ASSERT(r2 == 0);

if (r1 == 0 && r2 == 0) {
amrex::Print() << "SUCCESS\n";
}
}
amrex::Finalize();
}

0 comments on commit 12002e7

Please sign in to comment.