parallel pickling with mesh-merging

This commit is contained in:
Joachim Schöberl 2020-08-19 14:50:11 +02:00
parent 9e105c48ea
commit 5e09626466
10 changed files with 294 additions and 33 deletions

View File

@ -567,6 +567,10 @@ namespace ngcore
virtual void FlushBuffer() {}
bool parallel = false;
bool IsParallel() const { return parallel; }
void SetParallel (bool _parallel) { parallel = _parallel; }
private:
template<typename T, typename ... Bases>
friend class RegisterClassForArchive;

View File

@ -858,7 +858,7 @@ namespace ngcore
size++;
}
NETGEN_INLINE Array<T> & operator += (const T & el)
NETGEN_INLINE Array & operator += (const T & el)
{
Append (el);
return *this;

View File

@ -26,6 +26,9 @@ namespace ngcore
template <> struct MPI_typetrait<char> {
static MPI_Datatype MPIType () { return MPI_CHAR; } };
template <> struct MPI_typetrait<signed char> {
static MPI_Datatype MPIType () { return MPI_CHAR; } };
template <> struct MPI_typetrait<unsigned char> {
static MPI_Datatype MPIType () { return MPI_CHAR; } };
@ -44,6 +47,10 @@ namespace ngcore
return MPI_typetrait<T>::MPIType();
}
template <class T>
inline MPI_Datatype GetMPIType (T &) {
return GetMPIType<T>();
}
class NgMPI_Comm
{
@ -139,8 +146,8 @@ namespace ngcore
MPI_Send( const_cast<char*> (&s[0]), s.length(), MPI_CHAR, dest, tag, comm);
}
template<typename T, typename T2 = decltype(GetMPIType<T>())>
void Send(FlatArray<T> s, int dest, int tag) const {
template<typename T, typename TI, typename T2 = decltype(GetMPIType<T>())>
void Send(FlatArray<T,TI> s, int dest, int tag) const {
MPI_Send (s.Data(), s.Size(), GetMPIType<T>(), dest, tag, comm);
}
@ -160,13 +167,13 @@ namespace ngcore
}
template <typename T, typename T2 = decltype(GetMPIType<T>())>
void Recv (FlatArray <T> s, int src, int tag) const {
template <typename T, typename TI, typename T2 = decltype(GetMPIType<T>())>
void Recv (FlatArray <T,TI> s, int src, int tag) const {
MPI_Recv (s.Data(), s.Size(), GetMPIType<T> (), src, tag, comm, MPI_STATUS_IGNORE);
}
template <typename T, typename T2 = decltype(GetMPIType<T>())>
void Recv (Array <T> & s, int src, int tag) const
template <typename T, typename TI, typename T2 = decltype(GetMPIType<T>())>
void Recv (Array <T,TI> & s, int src, int tag) const
{
MPI_Status status;
int len;

View File

@ -8,6 +8,8 @@ using std::string;
namespace ngcore
{
bool ngcore_have_numpy = false;
bool parallel_pickling = false;
void SetFlag(Flags &flags, string s, py::object value)
{
if (py::isinstance<py::dict>(value))

View File

@ -17,7 +17,8 @@ namespace py = pybind11;
namespace ngcore
{
NGCORE_API extern bool ngcore_have_numpy;
NGCORE_API extern bool parallel_pickling;
// Python class name type traits
template <typename T>
struct PyNameTraits {
@ -297,26 +298,14 @@ namespace ngcore
return pybind11::pickle([](T* self)
{
PyArchive<T_ARCHIVE_OUT> ar;
ar.SetParallel(parallel_pickling);
ar & self;
static Timer t("ngspickle 2");
t.Start();
auto output = pybind11::make_tuple(ar.WriteOut());
t.Stop();
/*
GetLogger("Archive")->trace("Pickling output for object of type {} = {}",
Demangle(typeid(T).name()),
std::string(pybind11::str(output)));
*/
return output;
},
[](const pybind11::tuple & state)
{
T* val = nullptr;
/*
GetLogger("Archive")->trace("State for unpickling of object of type {} = {}",
Demangle(typeid(T).name()),
std::string(pybind11::str(state[0])));
*/
PyArchive<T_ARCHIVE_IN> ar(state[0]);
ar & val;
return val;

View File

@ -1281,13 +1281,132 @@ namespace netgen
void Mesh :: DoArchive (Archive & archive)
{
static Timer t("Mesh::Archive"); RegionTimer r(t);
static Timer tvol("Mesh::Archive vol elements");
#ifdef PARALLEL
auto comm = GetCommunicator();
if (archive.IsParallel() && comm.Size() > 1)
{ // parallel pickling supported only for output archives
if (comm.Rank() == 0)
archive & dimension;
// merge points
auto & partop = GetParallelTopology();
Array<int, PointIndex> globnum(points.Size());
int maxglob = 0;
for (auto pi : Range(points))
{
globnum[pi] = partop.GetGlobalPNum(pi);
maxglob = max(globnum[pi], maxglob);
}
maxglob = comm.AllReduce (maxglob, MPI_MAX);
int numglob = maxglob+1-PointIndex::BASE;
if (comm.Rank() > 0)
{
comm.Send (globnum, 0, 200);
comm.Send (points, 0, 200);
}
else
{
Array<int, PointIndex> globnumi;
Array<MeshPoint, PointIndex> pointsi;
Array<MeshPoint, PointIndex> globpoints(numglob);
for (int j = 1; j < comm.Size(); j++)
{
comm.Recv (globnumi, j, 200);
comm.Recv (pointsi, j, 200);
for (auto i : Range(globnumi))
globpoints[globnumi[i]] = pointsi[i];
}
archive & globpoints;
}
// sending surface elements
auto copy_el2d (surfelements);
for (auto & el : copy_el2d)
for (auto & pi : el.PNums())
pi = globnum[pi];
if (comm.Rank() > 0)
comm.Send(copy_el2d, 0, 200);
else
{
Array<Element2d, SurfaceElementIndex> el2di;
for (int j = 1; j < comm.Size(); j++)
{
comm.Recv(el2di, j, 200);
for (auto & el : el2di)
copy_el2d += el;
}
archive & copy_el2d;
}
// sending volume elements
auto copy_el3d (volelements);
for (auto & el : copy_el3d)
for (auto & pi : el.PNums())
pi = globnum[pi];
if (comm.Rank() > 0)
comm.Send(copy_el3d, 0, 200);
else
{
Array<Element, ElementIndex> el3di;
for (int j = 1; j < comm.Size(); j++)
{
comm.Recv(el3di, j, 200);
for (auto & el : el3di)
copy_el3d += el;
}
archive & copy_el3d;
}
// sending 1D elements
auto copy_el1d (segments);
for (auto & el : copy_el1d)
for (auto & pi : el.pnums)
if (pi != PointIndex(PointIndex::INVALID))
pi = globnum[pi];
if (comm.Rank() > 0)
comm.Send(copy_el1d, 0, 200);
else
{
Array<Segment, SegmentIndex> el1di;
for (int j = 1; j < comm.Size(); j++)
{
comm.Recv(el1di, j, 200);
for (auto & el : el1di)
copy_el1d += el;
}
archive & copy_el1d;
}
if (comm.Rank() == 0)
{
archive & facedecoding;
archive & materials & bcnames & cd2names & cd3names;
auto mynv = numglob;
archive & mynv; // numvertices;
archive & *ident;
archive.Shallow(geometry);
archive & *curvedelems;
}
if (comm.Rank() == 0)
return;
}
#endif
archive & dimension;
archive & points;
archive & surfelements;
tvol.Start();
archive & volelements;
tvol.Stop();
archive & segments;
archive & facedecoding;
archive & materials & bcnames & cd2names & cd3names;

View File

@ -45,6 +45,102 @@ namespace netgen
}
return type;
}
MPI_Datatype Element2d :: MyGetMPIType ( )
{
static MPI_Datatype type = MPI_DATATYPE_NULL;
static MPI_Datatype htype = MPI_DATATYPE_NULL;
if (type == MPI_DATATYPE_NULL)
{
Element2d hel;
int blocklen[] = { ELEMENT2D_MAXPOINTS, 1, 1, 1 };
MPI_Aint displ[] =
{ (char*)&hel.pnum[0] - (char*)&hel,
(char*)&hel.index - (char*)&hel,
(char*)&hel.typ - (char*)&hel,
(char*)&hel.np - (char*)&hel
};
MPI_Datatype types[] = { GetMPIType<PointIndex>(), GetMPIType(hel.index),
GetMPIType(hel.typ), GetMPIType(hel.np) };
// *testout << "displ = " << displ[0] << ", " << displ[1] << ", " << displ[2] << endl;
// *testout << "sizeof = " << sizeof (MeshPoint) << endl;
MPI_Type_create_struct (4, blocklen, displ, types, &htype);
MPI_Type_commit ( &htype );
MPI_Aint lb, ext;
MPI_Type_get_extent (htype, &lb, &ext);
// *testout << "lb = " << lb << endl;
// *testout << "ext = " << ext << endl;
ext = sizeof (Element2d);
MPI_Type_create_resized (htype, lb, ext, &type);
MPI_Type_commit ( &type );
}
return type;
}
MPI_Datatype Element :: MyGetMPIType ( )
{
static MPI_Datatype type = MPI_DATATYPE_NULL;
static MPI_Datatype htype = MPI_DATATYPE_NULL;
if (type == MPI_DATATYPE_NULL)
{
Element hel;
int blocklen[] = { ELEMENT_MAXPOINTS, 1, 1, 1 };
MPI_Aint displ[] =
{ (char*)&hel.pnum[0] - (char*)&hel,
(char*)&hel.index - (char*)&hel,
(char*)&hel.typ - (char*)&hel,
(char*)&hel.np - (char*)&hel
};
MPI_Datatype types[] = { GetMPIType<PointIndex>(), GetMPIType(hel.index),
GetMPIType(hel.typ), GetMPIType(hel.np) };
// *testout << "displ = " << displ[0] << ", " << displ[1] << ", " << displ[2] << endl;
// *testout << "sizeof = " << sizeof (MeshPoint) << endl;
MPI_Type_create_struct (4, blocklen, displ, types, &htype);
MPI_Type_commit ( &htype );
MPI_Aint lb, ext;
MPI_Type_get_extent (htype, &lb, &ext);
// *testout << "lb = " << lb << endl;
// *testout << "ext = " << ext << endl;
ext = sizeof (Element);
MPI_Type_create_resized (htype, lb, ext, &type);
MPI_Type_commit ( &type );
}
return type;
}
MPI_Datatype Segment :: MyGetMPIType ( )
{
static MPI_Datatype type = MPI_DATATYPE_NULL;
static MPI_Datatype htype = MPI_DATATYPE_NULL;
if (type == MPI_DATATYPE_NULL)
{
Segment hel;
int blocklen[] = { 3, 1, 1, 1 };
MPI_Aint displ[] =
{ (char*)&hel.pnums[0] - (char*)&hel,
(char*)&hel.edgenr - (char*)&hel,
(char*)&hel.cd2i - (char*)&hel,
(char*)&hel.si - (char*)&hel
};
MPI_Datatype types[] = {
GetMPIType<PointIndex>(), GetMPIType(hel.edgenr), GetMPIType(hel.cd2i), GetMPIType(hel.si)
};
// *testout << "displ = " << displ[0] << ", " << displ[1] << ", " << displ[2] << endl;
// *testout << "sizeof = " << sizeof (MeshPoint) << endl;
MPI_Type_create_struct (4, blocklen, displ, types, &htype);
MPI_Type_commit ( &htype );
MPI_Aint lb, ext;
MPI_Type_get_extent (htype, &lb, &ext);
// *testout << "lb = " << lb << endl;
// *testout << "ext = " << ext << endl;
ext = sizeof (Segment);
MPI_Type_create_resized (htype, lb, ext, &type);
MPI_Type_commit ( &type );
}
return type;
}
#endif
@ -158,7 +254,8 @@ namespace netgen
<< " si = " << seg.si << ", edgenr = " << seg.edgenr;
return s;
}
/*
// needed, e.g. for MPI communication
Element2d :: Element2d ()
{
for (int i = 0; i < ELEMENT2D_MAXPOINTS; i++)
@ -177,7 +274,7 @@ namespace netgen
strongrefflag = false;
is_curved = false;
}
*/
Element2d :: Element2d (int anp)
{
for (int i = 0; i < ELEMENT2D_MAXPOINTS; i++)

View File

@ -228,8 +228,8 @@ namespace netgen
{
int i;
public:
ElementIndex () { ; }
ElementIndex (int ai) : i(ai) { ; }
ElementIndex () = default;
constexpr ElementIndex (int ai) : i(ai) { ; }
ElementIndex & operator= (const ElementIndex & ai) { i = ai.i; return *this; }
ElementIndex & operator= (int ai) { i = ai; return *this; }
operator int () const { return i; }
@ -288,8 +288,8 @@ namespace netgen
{
int i;
public:
SegmentIndex () { ; }
SegmentIndex (int ai) : i(ai) { ; }
SegmentIndex () = default;
constexpr SegmentIndex (int ai) : i(ai) { ; }
SegmentIndex & operator= (const SegmentIndex & ai)
{ i = ai.i; return *this; }
SegmentIndex & operator= (int ai) { i = ai; return *this; }
@ -393,7 +393,7 @@ namespace netgen
///
ELEMENT_TYPE typ;
/// number of points
unsigned int np:4;
int8_t np;
bool badel:1;
bool refflag:1; // marked for refinement
bool strongrefflag:1;
@ -417,7 +417,7 @@ namespace netgen
public:
///
Element2d () = default;
Element2d ();
Element2d (const Element2d &) = default;
Element2d (Element2d &&) = default;
Element2d & operator= (const Element2d &) = default;
@ -548,6 +548,11 @@ namespace netgen
ar & pnum[i];
}
#ifdef PARALLEL
static MPI_Datatype MyGetMPIType();
#endif
void SetIndex (int si) { index = si; }
///
int GetIndex () const { return index; }
@ -699,7 +704,7 @@ namespace netgen
///
ELEMENT_TYPE typ;
/// number of points (4..tet, 5..pyramid, 6..prism, 8..hex, 10..quad tet, 12..quad prism)
int np:6;
int8_t np;
///
class flagstruct {
public:
@ -826,6 +831,10 @@ namespace netgen
for (size_t i = 0; i < np; i++)
ar & pnum[i];
}
#ifdef PARALLEL
static MPI_Datatype MyGetMPIType();
#endif
///
void SetIndex (int si) { index = si; }
@ -1074,6 +1083,10 @@ namespace netgen
*/
void DoArchive (Archive & ar);
#ifdef PARALLEL
static MPI_Datatype MyGetMPIType();
#endif
};
ostream & operator<<(ostream & s, const Segment & seg);
@ -1547,6 +1560,33 @@ namespace netgen
}
#ifdef PARALLEL
namespace ngcore
{
template <> struct MPI_typetrait<netgen::PointIndex> {
static MPI_Datatype MPIType () { return MPI_INT; }
};
template <> struct MPI_typetrait<netgen::ELEMENT_TYPE> {
static MPI_Datatype MPIType () { return MPI_CHAR; }
};
template <> struct MPI_typetrait<netgen::MeshPoint> {
static MPI_Datatype MPIType () { return netgen::MeshPoint::MyGetMPIType(); }
};
template <> struct MPI_typetrait<netgen::Element> {
static MPI_Datatype MPIType () { return netgen::Element::MyGetMPIType(); }
};
template <> struct MPI_typetrait<netgen::Element2d> {
static MPI_Datatype MPIType () { return netgen::Element2d::MyGetMPIType(); }
};
template <> struct MPI_typetrait<netgen::Segment> {
static MPI_Datatype MPIType () { return netgen::Segment::MyGetMPIType(); }
};
}
#endif
#endif

View File

@ -25,10 +25,12 @@ namespace metis {
using namespace metis;
#endif
/*
namespace ngcore {
template <> struct MPI_typetrait<netgen::PointIndex> {
static MPI_Datatype MPIType () { return MPI_INT; } };
}
*/
namespace netgen
{

View File

@ -1251,6 +1251,7 @@ grow_edges : bool = False
py::class_<ClearSolutionClass> (m, "ClearSolutionClass")
.def(py::init<>())
;
m.def("SetParallelPickling", [](bool par) { parallel_pickling = par; });
}
PYBIND11_MODULE(libmesh, m) {