From: abn Date: Fri, 12 May 2023 13:10:52 +0000 (+0200) Subject: [DEC] Enhance ctors on Python side to accept mpi4py communicators X-Git-Tag: V9_11_0b1~1 X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=c04fb8eb0dfbaac9c3bb1145701665cb5c9b27c3;p=tools%2Fmedcoupling.git [DEC] Enhance ctors on Python side to accept mpi4py communicators + factory pattern with 'New' static method + can use 'MPI._addressof()' in direct constructors --- diff --git a/CMakeLists.txt b/CMakeLists.txt index f130ad150..952609497 100644 --- a/CMakeLists.txt +++ b/CMakeLists.txt @@ -231,6 +231,7 @@ IF(MEDCOUPLING_BUILD_DOC) SALOME_LOG_OPTIONAL_PACKAGE(Sphinx MEDCOUPLING_BUILD_DOC) ENDIF(MEDCOUPLING_BUILD_DOC) + # Detection report SALOME_PACKAGE_REPORT_AND_CHECK() diff --git a/doc/developer/doxygen/doxfiles/reference/misc/swig_mpi.dox b/doc/developer/doxygen/doxfiles/reference/misc/swig_mpi.dox index 84004333e..a11cc962b 100644 --- a/doc/developer/doxygen/doxfiles/reference/misc/swig_mpi.dox +++ b/doc/developer/doxygen/doxfiles/reference/misc/swig_mpi.dox @@ -21,4 +21,11 @@ void method_taking_a_MPI_Comm(MPI_Comm mpicomm); where the path to mpi4py.i needs to be adapted to your configuration. +In MEDCoupling DEC APIs (InterpKernelDEC and OverlapDEC constructors), the Python interface deals with that using the following trick to get the C++ address of the MPI communicator: + +\code{.py} +from mpi4py import MPI +MPI._addressof(mpicomm) # returns the C++ address of the MPI communicator +\endcode + */ diff --git a/resources/dev/mc_suppr_valgrind b/resources/dev/mc_suppr_valgrind index 7f536fdc9..9609c3887 100644 --- a/resources/dev/mc_suppr_valgrind +++ b/resources/dev/mc_suppr_valgrind @@ -62,5 +62,11 @@ fun:*PyInit* } - +{ + + Memcheck:Leak + fun:*alloc + ... + fun:_dl_catch_exception +} diff --git a/src/MEDCoupling_Swig/CMakeLists.txt b/src/MEDCoupling_Swig/CMakeLists.txt index 125dda076..387a34a08 100644 --- a/src/MEDCoupling_Swig/CMakeLists.txt +++ b/src/MEDCoupling_Swig/CMakeLists.txt @@ -37,7 +37,7 @@ ENDIF() SET(SWIG_MODULE_MEDCoupling_EXTRA_FLAGS "${NUMPY_DEFINITIONS};${SCIPY_DEFINITIONS}") IF(MEDCOUPLING_USE_64BIT_IDS) STRING(APPEND SWIG_MODULE_MEDCoupling_EXTRA_FLAGS ";-DMEDCOUPLING_USE_64BIT_IDS") -ENDIF(MEDCOUPLING_USE_64BIT_IDS) +ENDIF() SET (MEDCoupling_SWIG_DPYS_FILES MEDCouplingCommon.i diff --git a/src/ParaMEDMEM/InterpKernelDEC.hxx b/src/ParaMEDMEM/InterpKernelDEC.hxx index 95c40f705..f4b156329 100644 --- a/src/ParaMEDMEM/InterpKernelDEC.hxx +++ b/src/ParaMEDMEM/InterpKernelDEC.hxx @@ -130,8 +130,7 @@ namespace MEDCoupling public: InterpKernelDEC(); InterpKernelDEC(ProcessorGroup& source_group, ProcessorGroup& target_group); - InterpKernelDEC(const std::set& src_ids, const std::set& trg_ids, - const MPI_Comm& world_comm=MPI_COMM_WORLD); + InterpKernelDEC(const std::set& src_ids, const std::set& trg_ids, const MPI_Comm& world_comm=MPI_COMM_WORLD); virtual ~InterpKernelDEC(); void release(); diff --git a/src/ParaMEDMEM_Swig/ParaMEDMEMCommon.i b/src/ParaMEDMEM_Swig/ParaMEDMEMCommon.i index 290da2a8a..d47b8bf48 100644 --- a/src/ParaMEDMEM_Swig/ParaMEDMEMCommon.i +++ b/src/ParaMEDMEM_Swig/ParaMEDMEMCommon.i @@ -54,9 +54,7 @@ using namespace ICoCo; %include "ComponentTopology.hxx" %include "DEC.hxx" %include "DisjointDEC.hxx" -%include "InterpKernelDEC.hxx" %include "StructuredCoincidentDEC.hxx" -%include "OverlapDEC.hxx" %newobject MEDCoupling::ParaUMesh::New; %newobject MEDCoupling::ParaUMesh::getMesh; @@ -75,6 +73,10 @@ using namespace ICoCo; %newobject MEDCoupling::ParaSkyLineArray::getSkyLineArray; %newobject MEDCoupling::ParaSkyLineArray::getGlobalIdsArray; +%newobject MEDCoupling::InterpKernelDEC::_NewWithPG_internal; +%newobject MEDCoupling::InterpKernelDEC::_NewWithComm_internal; +%newobject MEDCoupling::OverlapDEC::_NewWithComm_internal; + %feature("unref") ParaSkyLineArray "$this->decrRef();" %feature("unref") ParaUMesh "$this->decrRef();" %feature("unref") ParaDataArrayInt32 "$this->decrRef();" @@ -274,16 +276,97 @@ namespace MEDCoupling } } }; -} -/* This object can be used only if MED_ENABLE_FVM is defined*/ -#ifdef MED_ENABLE_FVM -class NonCoincidentDEC : public DEC -{ -public: - NonCoincidentDEC(ProcessorGroup& source, ProcessorGroup& target); -}; -#endif + /* This object can be used only if MED_ENABLE_FVM is defined*/ + #ifdef MED_ENABLE_FVM + class NonCoincidentDEC : public DEC + { + public: + NonCoincidentDEC(ProcessorGroup& source, ProcessorGroup& target); + }; + #endif + + class InterpKernelDEC : public DisjointDEC, public INTERP_KERNEL::InterpolationOptions + { + public: + InterpKernelDEC(); + InterpKernelDEC(ProcessorGroup& source_group, ProcessorGroup& target_group); + InterpKernelDEC(const std::set& src_ids, const std::set& trg_ids); // hide last optional parameter! + virtual ~InterpKernelDEC(); + void release(); + + void synchronize(); + void recvData(); + void recvData(double time); + void sendData(); + void sendData(double time , double deltatime); + void prepareSourceDE(); + void prepareTargetDE(); + + %extend { + // Provides a direct ctor for which the communicator can be passed with "MPI._addressof(the_com)": + InterpKernelDEC(const std::set& src_ids, const std::set& trg_ids, long long comm_ptr) + { + return new InterpKernelDEC(src_ids, trg_ids, *((MPI_Comm*)comm_ptr)); + } + + // This one should really not be called directly by the user since it still has an interface with a pointer to MPI_Comm + // which Swig doesn't handle nicely. + // It is just here to provide a constructor taking a **pointer** to a comm - See pythoncode below. + static InterpKernelDEC* _NewWithPG_internal(ProcessorGroup& source_group, ProcessorGroup& target_group) + { + return new InterpKernelDEC(source_group,target_group); + } + + static InterpKernelDEC* _NewWithComm_internal(const std::set& src_ids, const std::set& trg_ids, long long another_comm) + { + return new InterpKernelDEC(src_ids,trg_ids, *(MPI_Comm*)another_comm); // I know, ugly cast ... + } + } + }; + + class OverlapDEC : public DEC, public INTERP_KERNEL::InterpolationOptions + { + public: + OverlapDEC(const std::set& procIds); // hide optional param comm + virtual ~OverlapDEC(); + void release(); + + void sendRecvData(bool way=true); + void sendData(); + void recvData(); + void synchronize(); + void attachSourceLocalField(ParaFIELD *field, bool ownPt=false); + void attachTargetLocalField(ParaFIELD *field, bool ownPt=false); + void attachSourceLocalField(MEDCouplingFieldDouble *field); + void attachTargetLocalField(MEDCouplingFieldDouble *field); + void attachSourceLocalField(ICoCo::MEDDoubleField *field); + void attachTargetLocalField(ICoCo::MEDDoubleField *field); + ProcessorGroup *getGroup(); + bool isInGroup() const; + + void setDefaultValue(double val); + void setWorkSharingAlgo(int method); + + void debugPrintWorkSharing(std::ostream & ostr) const; + + %extend { + OverlapDEC(const std::set& ids, long long comm_ptr) + { + return new OverlapDEC(ids, *((MPI_Comm*)comm_ptr)); + } + + // This one should really not be called directly by the user since it still has an interface with a pointer to MPI_Comm + // which Swig doesn't handle nicely. + // It is just here to provide a constructor taking a **pointer** to a comm - See pythoncode below. + static OverlapDEC* _NewWithComm_internal(const std::set& ids, long long another_comm) + { + return new OverlapDEC(ids, *(MPI_Comm*)another_comm); // I know, ugly cast ... + } + } + }; + +} // end namespace MEDCoupling %extend MEDCoupling::ParaMESH { @@ -324,3 +407,54 @@ if MEDCouplingUse64BitIDs(): else: ParaDataArrayInt = ParaDataArrayInt32 %} + +%pythoncode %{ + +# And here we use mpi4py ability to provide its internal (C++) pointer to the communicator: +# NB: doing a proper typemap from MPI_Comm from Python to C++ requires the inclusion of mpi4py headers and .i file ... an extra dependency ... +def _IKDEC_WithComm_internal(src_procs, tgt_procs, mpicomm=None): + from mpi4py import MPI + # Check iterable: + try: + s, t = [el for el in src_procs], [el for el in tgt_procs] + except: + s, t = None, None + msg = "InterpKernelDEC: invalid type in ctor arguments! Possible signatures are:\n" + msg += " - InterpKernelDEC(ProcessorGroup, ProcessorGroup)\n" + msg += " - InterpKernelDEC(, )\n" + msg += " - InterpKernelDEC(, , MPI_Comm*) : WARNING here the address of the communicator should be passed with MPI._addressof(the_com)\n" + msg += " - InterpKernelDEC.New(ProcessorGroup, ProcessorGroup)\n" + msg += " - InterpKernelDEC.New(, )\n" + msg += " - InterpKernelDEC.New(, , MPI_Comm)\n" + if mpicomm is None: + if isinstance(src_procs, ProcessorGroup) and isinstance(tgt_procs, ProcessorGroup): + return InterpKernelDEC._NewWithPG_internal(src_procs, tgt_procs) + elif not s is None: # iterable + return InterpKernelDEC._NewWithComm_internal(s, t, MPI._addressof(MPI.COMM_WORLD)) + else: + raise InterpKernelException(msg) + else: + if s is None: raise InterpKernelException(msg) # must be iterable + return InterpKernelDEC._NewWithComm_internal(s, t, MPI._addressof(mpicomm)) + +def _ODEC_WithComm_internal(procs, mpicomm=None): + from mpi4py import MPI + # Check iterable: + try: + g = [el for el in procs] + except: + msg = "OverlapDEC: invalid type in ctor arguments! Possible signatures are:\n" + msg += " - OverlapDEC.New()\n" + msg += " - OverlapDEC.New(, MPI_Comm)\n" + msg += " - OverlapDEC()\n" + msg += " - OverlapDEC(, MPI_Comm*) : WARNING here the address of the communicator should be passed with MPI._addressof(the_com)\n" + raise InterpKernelException(msg) + if mpicomm is None: + return OverlapDEC(g) + else: + return OverlapDEC._NewWithComm_internal(g, MPI._addressof(mpicomm)) + +InterpKernelDEC.New = _IKDEC_WithComm_internal +OverlapDEC.New = _ODEC_WithComm_internal + +%} diff --git a/src/ParaMEDMEM_Swig/test_InterpKernelDEC.py b/src/ParaMEDMEM_Swig/test_InterpKernelDEC.py index ef68d76a0..fe30d14cb 100755 --- a/src/ParaMEDMEM_Swig/test_InterpKernelDEC.py +++ b/src/ParaMEDMEM_Swig/test_InterpKernelDEC.py @@ -82,6 +82,37 @@ class ParaMEDMEM_IK_DEC_Tests(unittest.TestCase): fld.setMesh(sub_m) return sub_m, fld + def testInterpKernelDEC_ctor(self): + """ Test the various Python ctors """ + size = MPI.COMM_WORLD.size + if size != 4: + print("Should be run on 4 procs!") + return + # Define two processor groups + nproc_source = 2 + l1, l2 = range(nproc_source), range(size - nproc_source, size) + # With 2 iterables: + i1 = InterpKernelDEC.New(l1, l2) + # Should also work directly: + i2 = InterpKernelDEC(l1, l2) + # With 2 proc groups: + interface = CommInterface() + source_group = MPIProcessorGroup(interface, list(l1)) + target_group = MPIProcessorGroup(interface, list(l2)) + i3 = InterpKernelDEC.New(source_group, target_group) + # Should also work directly: + i4 = InterpKernelDEC(source_group, target_group) + # With 2 iterables and a custom comm: + i5 = InterpKernelDEC.New(l1, l2, MPI.COMM_WORLD) + # Work directly with the **hack** + i6 = InterpKernelDEC(l1, l2, MPI._addressof(MPI.COMM_WORLD)) + # Should fail with 2 proc groups **and** a communicator + self.assertRaises(InterpKernelException, InterpKernelDEC.New, source_group, target_group, MPI.COMM_WORLD) + self.assertRaises(NotImplementedError, InterpKernelDEC, source_group, target_group, MPI.COMM_WORLD) + i6.release(); i5.release(); i4.release(); i3.release(); i2.release(); i1.release() + source_group.release() + target_group.release() + @WriteInTmpDir def testInterpKernelDEC_2D_py_1(self): """ This test illustrates a basic use of the InterpKernelDEC. diff --git a/src/ParaMEDMEM_Swig/test_OverlapDEC.py b/src/ParaMEDMEM_Swig/test_OverlapDEC.py index 48df8f4dc..778f40b03 100755 --- a/src/ParaMEDMEM_Swig/test_OverlapDEC.py +++ b/src/ParaMEDMEM_Swig/test_OverlapDEC.py @@ -88,6 +88,25 @@ class ParaMEDMEM_O_DEC_Tests(unittest.TestCase): fld.setMesh(sub_m) return sub_m, fld + def testOverlapDEC_ctor(self): + """ Test the various Python ctors """ + size = MPI.COMM_WORLD.size + if size != 4: + print("Should be run on 4 procs!") + return + # Define processor group + proc_group = list(range(size)) + # With 2 iterables: + o1 = OverlapDEC.New(proc_group) + # Should also work directly: + o2 = OverlapDEC(proc_group) + # With an iterable and a custom comm: + o3 = OverlapDEC.New(proc_group, MPI.COMM_WORLD) + # Also work directly with the **hack** on the comm: + o4 = OverlapDEC(proc_group, MPI._addressof(MPI.COMM_WORLD)) + self.assertRaises(NotImplementedError, OverlapDEC, proc_group, MPI.COMM_WORLD) + o4.release(); o3.release(); o2.release(); o1.release() + @WriteInTmpDir def testOverlapDEC_2D_py_1(self): """ The main method of the test """