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

New constructor for MPMD::Copier #3806

Merged
merged 34 commits into from
Mar 19, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
Show all changes
34 commits
Select commit Hold shift + click to select a range
d89e4d6
First MPMD changes to call MPI split on python side
Feb 6, 2024
8975760
Comments added to new functions in AMReX_MPMD
Feb 7, 2024
4375052
Merge branch 'AMReX-Codes:development' into development
siddanib Feb 7, 2024
77e24f8
Reducing code duplication in AMReX_MPMD
Feb 7, 2024
514cf42
Attempting to remove trailing whitespaces.
Feb 8, 2024
7259543
Merge branch 'AMReX-Codes:development' into development
siddanib Feb 8, 2024
2bdf597
Merge branch 'AMReX-Codes:development' into development
siddanib Feb 8, 2024
5d6b4e8
Merge branch 'AMReX-Codes:development' into development
siddanib Feb 12, 2024
68048ac
Merge branch 'AMReX-Codes:development' into development
siddanib Feb 13, 2024
025a7b1
Merge branch 'AMReX-Codes:development' into development
siddanib Feb 14, 2024
31f8537
Merge branch 'AMReX-Codes:development' into development
siddanib Feb 14, 2024
34fb942
Merge branch 'AMReX-Codes:development' into development
siddanib Feb 15, 2024
b72e91d
Merge branch 'AMReX-Codes:development' into development
siddanib Feb 17, 2024
0435785
Merge branch 'AMReX-Codes:development' into development
siddanib Feb 23, 2024
7d21faf
Merge branch 'AMReX-Codes:development' into development
siddanib Feb 23, 2024
c86feab
Merge branch 'AMReX-Codes:development' into development
siddanib Mar 4, 2024
b200869
Merge branch 'AMReX-Codes:development' into development
siddanib Mar 11, 2024
bc9ef4f
Merge branch 'AMReX-Codes:development' into development
siddanib Mar 13, 2024
83d8e41
IndexType check in MPMD
Mar 13, 2024
4003d88
Merge branch 'AMReX-Codes:development' into development
siddanib Mar 15, 2024
e598d1b
New MPMD::Copier constructor
Mar 15, 2024
55b194a
Merge branch 'AMReX-Codes:development' into development
siddanib Mar 15, 2024
8e18973
Merge branch 'AMReX-Codes:development' into development
siddanib Mar 16, 2024
0afb301
Merge branch 'AMReX-Codes:development' into development
siddanib Mar 17, 2024
2a32c74
Merge branch 'AMReX-Codes:development' into development
siddanib Mar 18, 2024
283f250
Cleaner new Copier ctor
Mar 18, 2024
4983a78
Update Src/Base/AMReX_MPMD.cpp
siddanib Mar 18, 2024
4f50dc4
Update Src/Base/AMReX_MPMD.H
siddanib Mar 18, 2024
d61ac38
Update Src/Base/AMReX_MPMD.H
siddanib Mar 18, 2024
23d09f0
Update Src/Base/AMReX_MPMD.cpp
siddanib Mar 18, 2024
6d644be
Update Src/Base/AMReX_MPMD.cpp
siddanib Mar 18, 2024
6a510a9
minor updates
WeiqunZhang Mar 18, 2024
162ae84
Removed isects_o variable
Mar 18, 2024
7bf84c8
Merge branch 'AMReX-Codes:development' into development
siddanib Mar 19, 2024
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
12 changes: 11 additions & 1 deletion Src/Base/AMReX_MPMD.H
Original file line number Diff line number Diff line change
Expand Up @@ -26,18 +26,27 @@ int MyProgId (); //! Program ID
class Copier
{
public:
Copier (BoxArray const& ba, DistributionMapping const& dm);
Copier (bool);

Copier (BoxArray const& ba, DistributionMapping const& dm,
bool send_ba = false);

template <typename FAB>
void send (FabArray<FAB> const& mf, int icomp, int ncomp) const;

template <typename FAB>
void recv (FabArray<FAB>& mf, int icomp, int ncomp) const;

[[nodiscard]] BoxArray const& boxArray () const;

[[nodiscard]] DistributionMapping const& DistributionMap () const;

private:
std::map<int,FabArrayBase::CopyComTagsContainer> m_SndTags;
std::map<int,FabArrayBase::CopyComTagsContainer> m_RcvTags;
bool m_is_thread_safe;
BoxArray m_ba;
DistributionMapping m_dm;
};

template <typename FAB>
Expand Down Expand Up @@ -176,6 +185,7 @@ void Copier::recv (FabArray<FAB>& mf, int icomp, int ncomp) const
}
}


}

#endif
Expand Down
173 changes: 151 additions & 22 deletions Src/Base/AMReX_MPMD.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -139,7 +139,9 @@ int MyProgId ()
return (myproc == ParallelDescriptor::MyProc()) ? 0 : 1;
}

Copier::Copier (BoxArray const& ba, DistributionMapping const& dm)
Copier::Copier (BoxArray const& ba, DistributionMapping const& dm,
bool send_ba)
siddanib marked this conversation as resolved.
Show resolved Hide resolved
: m_ba(ba), m_dm(dm)
{
int rank_offset = myproc - ParallelDescriptor::MyProc();
int this_root, other_root;
Expand All @@ -152,7 +154,6 @@ Copier::Copier (BoxArray const& ba, DistributionMapping const& dm)
}

Vector<Box> bv = ba.boxList().data();

int this_nboxes = static_cast<int>(ba.size());
Vector<int> procs = dm.ProcessorMap();
if (rank_offset != 0) {
Expand All @@ -163,34 +164,46 @@ Copier::Copier (BoxArray const& ba, DistributionMapping const& dm)

Vector<Box> obv;
Vector<int> oprocs;
int other_nboxes;
int other_nboxes = this_nboxes;
if (myproc == this_root) {
if (rank_offset == 0) // the first program
{
MPI_Send(&this_nboxes, 1, MPI_INT, other_root, 0, MPI_COMM_WORLD);
MPI_Recv(&other_nboxes, 1, MPI_INT, other_root, 1, MPI_COMM_WORLD,
if (!send_ba)
{
MPI_Recv(&other_nboxes, 1, MPI_INT, other_root, 1, MPI_COMM_WORLD,
MPI_STATUS_IGNORE);
obv.resize(other_nboxes);
obv.resize(other_nboxes);
}
MPI_Send(bv.data(), this_nboxes,
ParallelDescriptor::Mpi_typemap<Box>::type(),
other_root, 2, MPI_COMM_WORLD);
MPI_Recv(obv.data(), other_nboxes,
if (!send_ba)
{
MPI_Recv(obv.data(), other_nboxes,
ParallelDescriptor::Mpi_typemap<Box>::type(),
other_root, 3, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
oprocs.resize(other_nboxes);
}
MPI_Send(procs.data(), this_nboxes, MPI_INT, other_root, 4, MPI_COMM_WORLD);
oprocs.resize(other_nboxes);
MPI_Recv(oprocs.data(), other_nboxes, MPI_INT, other_root, 5, MPI_COMM_WORLD,
MPI_STATUS_IGNORE);
}
else // the second program
{
MPI_Recv(&other_nboxes, 1, MPI_INT, other_root, 0, MPI_COMM_WORLD,
if (!send_ba)
{
MPI_Recv(&other_nboxes, 1, MPI_INT, other_root, 0, MPI_COMM_WORLD,
MPI_STATUS_IGNORE);
obv.resize(other_nboxes);
}
MPI_Send(&this_nboxes, 1, MPI_INT, other_root, 1, MPI_COMM_WORLD);
obv.resize(other_nboxes);
MPI_Recv(obv.data(), other_nboxes,
if (!send_ba)
{
MPI_Recv(obv.data(), other_nboxes,
ParallelDescriptor::Mpi_typemap<Box>::type(),
other_root, 2, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
}
MPI_Send(bv.data(), this_nboxes,
ParallelDescriptor::Mpi_typemap<Box>::type(),
other_root, 3, MPI_COMM_WORLD);
Expand All @@ -201,30 +214,48 @@ Copier::Copier (BoxArray const& ba, DistributionMapping const& dm)
}
}

ParallelDescriptor::Bcast(&other_nboxes, 1);
if (obv.empty()) {
obv.resize(other_nboxes);
if (!send_ba) {
ParallelDescriptor::Bcast(&other_nboxes, 1);
if (obv.empty()){
obv.resize(other_nboxes);
}
ParallelDescriptor::Bcast(obv.data(), obv.size());
}

if (oprocs.empty()) {
oprocs.resize(other_nboxes);
}
ParallelDescriptor::Bcast(obv.data(), obv.size());
ParallelDescriptor::Bcast(oprocs.data(), oprocs.size());

BoxArray oba(BoxList(std::move(obv)));
BoxArray oba;
if (!obv.empty()) {
oba.define(BoxList(std::move(obv)));
}

// At this point, ba and bv hold our boxes, and oba holds the other
// program's boxes. procs holds mpi ranks of our boxes, and oprocs holds
// mpi ranks of the other program's boxes. All mpi ranks are in
// MPI_COMM_WORLD.

// Build communication meta-data
AMREX_ALWAYS_ASSERT(ba.ixType() == oba.ixType());
m_is_thread_safe = ba.ixType().cellCentered();
if (!send_ba){
AMREX_ALWAYS_ASSERT(ba.ixType() == oba.ixType());
m_is_thread_safe = ba.ixType().cellCentered();
}else{
m_is_thread_safe = true;
}

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

for (int i = 0; i < this_nboxes; ++i) {
if (procs[i] == myproc) {
oba.intersections(bv[i], isects);
if (!send_ba){
oba.intersections(bv[i], isects);
}
else{
isects.resize(0);
isects.emplace_back(i,bv[i]);
}
for (auto const& isec : isects) {
const int oi = isec.first;
const Box& bx = isec.second;
Expand All @@ -235,14 +266,112 @@ Copier::Copier (BoxArray const& ba, DistributionMapping const& dm)
}
}

for (auto& kv : m_SndTags) {
std::sort(kv.second.begin(), kv.second.end());
if (!send_ba){
for (auto& kv : m_SndTags) {
std::sort(kv.second.begin(), kv.second.end());
}
for (auto& kv : m_RcvTags) {
std::sort(kv.second.begin(), kv.second.end());
}
}
}

Copier::Copier (bool)
: m_is_thread_safe(true)
{
int rank_offset = myproc - ParallelDescriptor::MyProc();
int this_root, other_root;
if (rank_offset == 0) { // First program
this_root = 0;
other_root = ParallelDescriptor::NProcs();
} else {
this_root = rank_offset;
other_root = 0;
}

Vector<Box> bv;
int this_nboxes;

if (myproc == this_root) {
int tags[2];
if (rank_offset == 0) // the first program
{
tags[0] = 1;
tags[1] = 3;
}
else // the second program
{
tags[0] = 0;
tags[1] = 2;
}

MPI_Recv(&this_nboxes, 1, MPI_INT, other_root, tags[0], MPI_COMM_WORLD,
MPI_STATUS_IGNORE);
bv.resize(this_nboxes);
MPI_Recv(bv.data(), this_nboxes,
ParallelDescriptor::Mpi_typemap<Box>::type(),
other_root, tags[1], MPI_COMM_WORLD, MPI_STATUS_IGNORE);
}

ParallelDescriptor::Bcast(&this_nboxes, 1);
if (bv.empty()) {
bv.resize(this_nboxes);
}

ParallelDescriptor::Bcast(bv.data(), bv.size());
m_ba.define(BoxList(std::move(bv)));
m_dm.define(m_ba);
Vector<int> procs = m_dm.ProcessorMap();
if (rank_offset != 0) {
for (int i = 0; i < this_nboxes; ++i) {
procs[i] += rank_offset;
}
}

Vector<int> oprocs(this_nboxes);
if (myproc == this_root) {
if (rank_offset == 0) // the first program
{
MPI_Send(procs.data(), this_nboxes, MPI_INT, other_root, 4, MPI_COMM_WORLD);
MPI_Recv(oprocs.data(), this_nboxes, MPI_INT, other_root, 5, MPI_COMM_WORLD,
MPI_STATUS_IGNORE);
}
else // the second program
{
MPI_Recv(oprocs.data(), this_nboxes, MPI_INT, other_root, 4, MPI_COMM_WORLD,
MPI_STATUS_IGNORE);
MPI_Send(procs.data(), this_nboxes, MPI_INT, other_root, 5, MPI_COMM_WORLD);
}
}
for (auto& kv : m_RcvTags) {
std::sort(kv.second.begin(), kv.second.end());

ParallelDescriptor::Bcast(oprocs.data(), oprocs.size());

// procs holds mpi ranks of our boxes, and oprocs holds
// mpi ranks of the other program's boxes. All mpi ranks are in
// MPI_COMM_WORLD.

// Build communication meta-data

for (int i = 0; i < this_nboxes; ++i) {
if (procs[i] == myproc) {
const Box& bx = m_ba[i];
const int orank = oprocs[i];
m_SndTags[orank].emplace_back(bx, bx, i, i);
m_RcvTags[orank].emplace_back(bx, bx, i, i);
}
}
}

BoxArray const& Copier::boxArray () const
{
return m_ba;
}

DistributionMapping const& Copier::DistributionMap () const
{
return m_dm;
}

}

#endif
Loading