Salome HOME
[bos #38048] [EDF] (2023-T3) PARAMEDMEM Ergonomy.
[tools/medcoupling.git] / src / ParaMEDMEM_Swig / ParaMEDMEMCommon.i
1 // Copyright (C) 2017-2024  CEA, EDF
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 %include std_set.i
21
22 %template() std::set<int>;
23
24 %{
25 #include "CommInterface.hxx"
26 #include "ProcessorGroup.hxx"
27 #include "Topology.hxx"
28 #include "MPIProcessorGroup.hxx"
29 #include "ByStringMPIProcessorGroup.hxx"
30 #include "DEC.hxx"
31 #include "InterpKernelDEC.hxx"
32 #include "NonCoincidentDEC.hxx"
33 #include "StructuredCoincidentDEC.hxx"
34 #include "OverlapDEC.hxx"
35 #include "ParaMESH.hxx"
36 #include "ParaFIELD.hxx"
37 #include "ICoCoMEDDoubleField.hxx"
38 #include "ICoCoMEDIntField.hxx"
39 #include "ComponentTopology.hxx"
40 #include "ParaUMesh.hxx"
41 #include "ParaSkyLineArray.hxx"
42 #include "ParaDataArray.hxx"
43
44 using namespace INTERP_KERNEL;
45 using namespace MEDCoupling;
46 using namespace ICoCo;
47 %}
48
49 %include "InterpolationOptions.hxx"
50 %include "ProcessorGroup.hxx"
51 %include "DECOptions.hxx"
52 %include "ParaMESH.hxx"
53 %include "ParaFIELD.hxx"
54 %include "MPIProcessorGroup.hxx"
55 %include "ByStringMPIProcessorGroup.hxx"
56 %include "ComponentTopology.hxx"
57 %include "DEC.hxx"
58 %include "DisjointDEC.hxx"
59 %include "StructuredCoincidentDEC.hxx"
60
61 %newobject MEDCoupling::ParaUMesh::New;
62 %newobject MEDCoupling::ParaUMesh::getMesh;
63 %newobject MEDCoupling::ParaUMesh::getGlobalCellIds;
64 %newobject MEDCoupling::ParaUMesh::getGlobalNodeIds;
65 %newobject MEDCoupling::ParaUMesh::getCellIdsLyingOnNodes;
66 %newobject MEDCoupling::ParaUMesh::redistributeCells;
67 %newobject MEDCoupling::ParaUMesh::redistributeCellField;
68 %newobject MEDCoupling::ParaUMesh::redistributeNodeField;
69 %newobject MEDCoupling::ParaDataArrayInt32::New;
70 %newobject MEDCoupling::ParaDataArrayInt32::buildComplement;
71 %newobject MEDCoupling::ParaDataArrayInt64::New;
72 %newobject MEDCoupling::ParaDataArrayInt64::buildComplement;
73 %newobject MEDCoupling::ParaSkyLineArray::New;
74 %newobject MEDCoupling::ParaSkyLineArray::equiRedistribute;
75 %newobject MEDCoupling::ParaSkyLineArray::getSkyLineArray;
76 %newobject MEDCoupling::ParaSkyLineArray::getGlobalIdsArray;
77
78 %newobject MEDCoupling::InterpKernelDEC::_NewWithPG_internal;
79 %newobject MEDCoupling::InterpKernelDEC::_NewWithComm_internal;
80 %newobject MEDCoupling::InterpKernelDEC::retrieveNonFetchedIds;
81 %newobject MEDCoupling::OverlapDEC::_NewWithComm_internal;
82
83 %feature("unref") ParaSkyLineArray "$this->decrRef();"
84 %feature("unref") ParaUMesh "$this->decrRef();"
85 %feature("unref") ParaDataArrayInt32 "$this->decrRef();"
86 %feature("unref") ParaDataArrayInt64 "$this->decrRef();"
87
88 %nodefaultctor;
89
90 namespace MEDCoupling
91 {
92   class CommInterface
93   {
94   public:
95     CommInterface();
96     virtual ~CommInterface();
97     int worldSize() const;
98     int commSize(MPI_Comm comm, int* size) const;
99     int commRank(MPI_Comm comm, int* rank) const;
100     int commGroup(MPI_Comm comm, MPI_Group* group) const;
101     int groupIncl(MPI_Group group, int size, int* ranks, MPI_Group* group_output) const;
102     int commCreate(MPI_Comm comm, MPI_Group group, MPI_Comm* comm_output) const;
103     int groupFree(MPI_Group* group) const;
104     int commFree(MPI_Comm* comm) const;
105
106     int send(void* buffer, int count, MPI_Datatype datatype, int target, int tag, MPI_Comm comm) const;
107     int recv(void* buffer, int count, MPI_Datatype datatype, int source, int tag, MPI_Comm comm, MPI_Status* status) const;
108     int sendRecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, 
109                  int dest, int sendtag, void* recvbuf, int recvcount, 
110                  MPI_Datatype recvtype, int source, int recvtag, MPI_Comm comm,
111                  MPI_Status* status);
112
113     int Isend(void* buffer, int count, MPI_Datatype datatype, int target,
114               int tag, MPI_Comm comm, MPI_Request *request) const;
115     int Irecv(void* buffer, int count, MPI_Datatype datatype, int source,
116               int tag, MPI_Comm comm, MPI_Request* request) const;
117
118     int wait(MPI_Request *request, MPI_Status *status) const;
119     int test(MPI_Request *request, int *flag, MPI_Status *status) const;
120     int requestFree(MPI_Request *request) const;
121     int waitany(int count, MPI_Request *array_of_requests, int *index, MPI_Status *status) const;
122     int testany(int count, MPI_Request *array_of_requests, int *index, int *flag, MPI_Status *status) const;
123     int waitall(int count, MPI_Request *array_of_requests, MPI_Status *array_of_status) const { return MPI_Waitall(count, array_of_requests, array_of_status); }
124     int testall(int count, MPI_Request *array_of_requests, int *flag, MPI_Status *array_of_status) const;
125     int waitsome(int incount, MPI_Request *array_of_requests,int *outcount, int *array_of_indices, MPI_Status *array_of_status) const;
126     int testsome(int incount, MPI_Request *array_of_requests, int *outcount,
127                  int *array_of_indices, MPI_Status *array_of_status) const;
128     int probe(int source, int tag, MPI_Comm comm, MPI_Status *status) const;
129     int Iprobe(int source, int tag, MPI_Comm comm, int *flag, MPI_Status *status) const;
130     int cancel(MPI_Request *request) const;
131     int testCancelled(MPI_Status *status, int *flag) const;
132     int barrier(MPI_Comm comm) const;
133     int errorString(int errorcode, char *string, int *resultlen) const;
134     int getCount(MPI_Status *status, MPI_Datatype datatype, int *count) const;
135
136     int broadcast(void* buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm) const;
137     int allGather(void* sendbuf, int sendcount, MPI_Datatype sendtype,
138                   void* recvbuf, int recvcount, MPI_Datatype recvtype,
139                   MPI_Comm comm) const;
140     int allToAll(void* sendbuf, int sendcount, MPI_Datatype sendtype,
141                  void* recvbuf, int recvcount, MPI_Datatype recvtype,
142                  MPI_Comm comm) const;
143     int allToAllV(void* sendbuf, int* sendcounts, int* senddispls,
144                   MPI_Datatype sendtype, void* recvbuf, int* recvcounts,
145                   int* recvdispls, MPI_Datatype recvtype, 
146                   MPI_Comm comm) const;
147
148     int reduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm) const;
149     int allReduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) const;
150     %extend
151     {
152       PyObject *allGatherArrays(const DataArrayIdType *array) const
153       {
154         std::vector< MCAuto<DataArrayIdType> > ret;
155         self->allGatherArrays(MPI_COMM_WORLD,array,ret);
156         return convertFromVectorAutoObjToPyObj<DataArrayIdType>(ret,SWIGTITraits<mcIdType>::TI);
157       }
158
159       PyObject *allToAllArrays(PyObject *arrays) const
160       {
161         std::vector< DataArrayIdType * > arraysIn;
162         std::vector< MCAuto<DataArrayIdType> > arrayOut;
163         convertFromPyObjVectorOfObj<MEDCoupling::DataArrayIdType*>(arrays,SWIGTITraits<mcIdType>::TI,"DataArrayIdType",arraysIn);
164         std::vector< MCAuto<DataArrayIdType> > arraysIn2(FromVecToVecAuto<DataArrayIdType>(arraysIn));
165         self->allToAllArrays(MPI_COMM_WORLD,arraysIn2,arrayOut);
166         return convertFromVectorAutoObjToPyObj<DataArrayIdType>(arrayOut,SWIGTITraits<mcIdType>::TI);
167       }
168     }
169   };
170
171   class ParaUMesh : public RefCountObject
172   {
173   public:
174     static ParaUMesh *New(MEDCouplingUMesh *mesh, DataArrayIdType *globalCellIds, DataArrayIdType *globalNodeIds);
175     ParaUMesh *redistributeCells(const DataArrayIdType *globalCellIds) const;
176     DataArrayIdType *redistributeCellField(const DataArrayIdType *globalCellIds, const DataArrayIdType *fieldValueToRed) const;
177     DataArrayDouble *redistributeCellField(const DataArrayIdType *globalCellIds, const DataArrayDouble *fieldValueToRed) const;
178     DataArrayIdType *redistributeNodeField(const DataArrayIdType *globalCellIds, const DataArrayIdType *fieldValueToRed) const;
179     DataArrayDouble *redistributeNodeField(const DataArrayIdType *globalCellIds, const DataArrayDouble *fieldValueToRed) const;
180     %extend
181     {
182       ParaUMesh(MEDCouplingUMesh *mesh, DataArrayIdType *globalCellIds, DataArrayIdType *globalNodeIds)
183       {
184         return ParaUMesh::New(mesh,globalCellIds,globalNodeIds);
185       }
186
187       MEDCouplingUMesh *getMesh()
188       {
189         MEDCouplingUMesh *ret(self->getMesh());
190         if(ret) ret->incrRef();
191         return ret;
192       }
193
194       DataArrayIdType *getGlobalCellIds()
195       {
196         DataArrayIdType *ret(self->getGlobalCellIds());
197         if(ret) ret->incrRef();
198         return ret;
199       }
200
201       DataArrayIdType *getGlobalNodeIds()
202       {
203         DataArrayIdType *ret(self->getGlobalNodeIds());
204         if(ret) ret->incrRef();
205         return ret;
206       }
207
208       DataArrayIdType *getCellIdsLyingOnNodes(const DataArrayIdType *globalNodeIds, bool fullyIn) const
209       { 
210         MCAuto<DataArrayIdType> ret(self->getCellIdsLyingOnNodes(globalNodeIds,fullyIn));
211         return ret.retn();
212       }
213     }
214   };
215
216   class ParaDataArray : public RefCountObject
217   {
218   };
219
220   class ParaDataArrayInt32 : public ParaDataArray
221   {
222   public:
223     static ParaDataArrayInt32 *New(DataArrayInt32 *seqDa);
224     DataArrayIdType *buildComplement(int nbOfElems) const;
225     %extend
226     {
227       ParaDataArrayInt32(DataArrayInt32 *seqDa)
228       {
229         return ParaDataArrayInt32::New(seqDa);
230       }
231     }
232   };
233
234   class ParaDataArrayInt64 : public ParaDataArray
235   {
236   public:
237     static ParaDataArrayInt64 *New(DataArrayInt64 *seqDa);
238     DataArrayIdType *buildComplement(long nbOfElems) const;
239     %extend
240     {
241       ParaDataArrayInt64(DataArrayInt64 *seqDa)
242       {
243         return ParaDataArrayInt64::New(seqDa);
244       }
245     }
246   };
247
248   class ParaSkyLineArray : public RefCountObject
249   {
250   public:
251     static ParaSkyLineArray *New(MEDCouplingSkyLineArray *ska, DataArrayIdType *globalIds);
252     %extend
253     {
254       ParaSkyLineArray(MEDCouplingSkyLineArray *ska, DataArrayIdType *globalIds)
255       {
256         return ParaSkyLineArray::New(ska,globalIds);
257       }
258
259       ParaSkyLineArray *equiRedistribute(mcIdType nbOfEntities) const
260       {
261         MCAuto<ParaSkyLineArray> ret(self->equiRedistribute(nbOfEntities));
262         return ret.retn();
263       }
264
265       MEDCouplingSkyLineArray *getSkyLineArray() const
266       {
267         MEDCouplingSkyLineArray *ret(self->getSkyLineArray());
268         if(ret)
269           ret->incrRef();
270         return ret;
271       }
272       
273       DataArrayIdType *getGlobalIdsArray() const
274       {
275         DataArrayIdType *ret(self->getGlobalIdsArray());
276         if(ret)
277           ret->incrRef();
278         return ret;
279       }
280     }
281   };
282
283   /* This object can be used only if MED_ENABLE_FVM is defined*/
284   #ifdef MED_ENABLE_FVM
285   class NonCoincidentDEC : public DEC
286   {
287   public:
288     NonCoincidentDEC(ProcessorGroup& source, ProcessorGroup& target);
289   };
290   #endif
291
292   class InterpKernelDEC : public DisjointDEC, public INTERP_KERNEL::InterpolationOptions
293   {
294     public:
295       InterpKernelDEC();
296       InterpKernelDEC(ProcessorGroup& source_group, ProcessorGroup& target_group);
297       InterpKernelDEC(const std::set<int>& src_ids, const std::set<int>& trg_ids); // hide last optional parameter!
298       InterpKernelDEC(ProcessorGroup& generic_group, const std::string& source_group, const std::string& target_group);
299       InterpKernelDEC(ProcessorGroup& generic_group, const std::string& interaction_group);
300       virtual ~InterpKernelDEC();
301       void release();
302
303       void synchronize();
304       void synchronizeWithDefaultValue(double val);
305       void recvData();
306       void recvData(double time);
307       void sendData();
308       void sendData(double time , double deltatime);
309       void prepareSourceDE();
310       void prepareTargetDE();
311
312       %extend {
313         // Provides a direct ctor for which the communicator can be passed with "MPI._addressof(the_com)":
314         InterpKernelDEC(const std::set<int>& src_ids, const std::set<int>& trg_ids, long long comm_ptr)
315         {
316             return new InterpKernelDEC(src_ids, trg_ids, *((MPI_Comm*)comm_ptr));
317         }
318
319         // This one should really not be called directly by the user since it still has an interface with a pointer to MPI_Comm 
320         // which Swig doesn't handle nicely.
321         // It is just here to provide a constructor taking a **pointer** to a comm - See pythoncode below.
322         static InterpKernelDEC* _NewWithPG_internal(ProcessorGroup& source_group, ProcessorGroup& target_group)
323         {
324           return new InterpKernelDEC(source_group,target_group);
325         }
326
327         static InterpKernelDEC* _NewWithComm_internal(const std::set<int>& src_ids, const std::set<int>& trg_ids, long long another_comm)
328         {
329           return new InterpKernelDEC(src_ids,trg_ids, *(MPI_Comm*)another_comm); // I know, ugly cast ...
330         }
331
332         DataArrayIdType *retrieveNonFetchedIds() const
333         {
334           MCAuto<DataArrayIdType> ret = self->retrieveNonFetchedIds();
335           return ret.retn();
336         }
337       }
338   };
339
340   class OverlapDEC : public DEC, public INTERP_KERNEL::InterpolationOptions
341   {
342       public:
343         OverlapDEC(const std::set<int>& procIds);  // hide optional param comm
344         virtual ~OverlapDEC();
345         void release();
346
347         void sendRecvData(bool way=true);
348         void sendData();
349         void recvData();
350         void synchronize();
351         void attachSourceLocalField(ParaFIELD *field, bool ownPt=false);
352         void attachTargetLocalField(ParaFIELD *field, bool ownPt=false);
353         void attachSourceLocalField(MEDCouplingFieldDouble *field);
354         void attachTargetLocalField(MEDCouplingFieldDouble *field);
355         void attachSourceLocalField(ICoCo::MEDDoubleField *field);
356         void attachTargetLocalField(ICoCo::MEDDoubleField *field);
357         ProcessorGroup *getGroup();
358         bool isInGroup() const;
359
360         void setDefaultValue(double val);
361         void setWorkSharingAlgo(int method);
362
363         void debugPrintWorkSharing(std::ostream & ostr) const;
364
365         %extend {
366           OverlapDEC(const std::set<int>& ids, long long comm_ptr)
367           {
368              return new OverlapDEC(ids, *((MPI_Comm*)comm_ptr));
369           }
370
371           // This one should really not be called directly by the user since it still has an interface with a pointer to MPI_Comm 
372           // which Swig doesn't handle nicely.
373           // It is just here to provide a constructor taking a **pointer** to a comm - See pythoncode below.
374           static OverlapDEC* _NewWithComm_internal(const std::set<int>& ids, long long another_comm)
375           {
376             return new OverlapDEC(ids, *(MPI_Comm*)another_comm); // I know, ugly cast ...
377           }
378         }
379    };
380
381 } // end namespace MEDCoupling
382
383 %extend MEDCoupling::ParaMESH
384 {
385   PyObject *getGlobalNumberingCell2() const
386   {
387     const mcIdType *tmp=self->getGlobalNumberingCell();
388     mcIdType size=self->getCellMesh()->getNumberOfCells();
389     PyObject *ret=PyList_New(size);
390     for(mcIdType i=0;i<size;i++)
391       PyList_SetItem(ret,i,PyInt_FromLong(tmp[i])); 
392     return ret;
393   }
394
395   PyObject *getGlobalNumberingFace2() const
396   {
397     const mcIdType *tmp=self->getGlobalNumberingFace();
398     mcIdType size=self->getFaceMesh()->getNumberOfCells();
399     PyObject *ret=PyList_New(size);
400     for(mcIdType i=0;i<size;i++)
401       PyList_SetItem(ret,i,PyInt_FromLong(tmp[i])); 
402     return ret;
403   }
404
405   PyObject *getGlobalNumberingNode2() const
406   {
407     const mcIdType *tmp=self->getGlobalNumberingNode();
408     mcIdType size=self->getCellMesh()->getNumberOfNodes();
409     PyObject *ret=PyList_New(size);
410     for(mcIdType i=0;i<size;i++)
411       PyList_SetItem(ret,i,PyInt_FromLong(tmp[i])); 
412     return ret;
413   }
414 }
415
416 %pythoncode %{
417 if MEDCouplingUse64BitIDs():
418   ParaDataArrayInt = ParaDataArrayInt64
419 else:
420   ParaDataArrayInt = ParaDataArrayInt32
421 %}
422
423 %pythoncode %{
424
425 # And here we use mpi4py ability to provide its internal (C++) pointer to the communicator:
426 # NB: doing a proper typemap from MPI_Comm from Python to C++ requires the inclusion of mpi4py headers and .i file ... an extra dependency ...
427 def _IKDEC_WithComm_internal(src_procs, tgt_procs, mpicomm=None):
428     from mpi4py import MPI
429     # Check iterable:
430     try:
431         s, t = [el for el in src_procs], [el for el in tgt_procs]
432     except:
433         s, t = None, None
434     msg =  "InterpKernelDEC: invalid type in ctor arguments! Possible signatures are:\n"
435     msg += "   - InterpKernelDEC(ProcessorGroup, ProcessorGroup)\n"
436     msg += "   - InterpKernelDEC(<iterable>, <iterable>)\n"
437     msg += "   - InterpKernelDEC(<iterable>, <iterable>, MPI_Comm*) : WARNING here the address of the communicator should be passed with MPI._addressof(the_com)\n"
438     msg += "   - InterpKernelDEC.New(ProcessorGroup, ProcessorGroup)\n"
439     msg += "   - InterpKernelDEC.New(<iterable>, <iterable>)\n"
440     msg += "   - InterpKernelDEC.New(<iterable>, <iterable>, MPI_Comm)\n"
441     if mpicomm is None:
442         if isinstance(src_procs, ProcessorGroup) and isinstance(tgt_procs, ProcessorGroup):
443             return InterpKernelDEC._NewWithPG_internal(src_procs, tgt_procs)
444         elif not s is None:  # iterable
445             return InterpKernelDEC._NewWithComm_internal(s, t, MPI._addressof(MPI.COMM_WORLD))
446         else:
447             raise InterpKernelException(msg)
448     else:
449         if s is None: raise InterpKernelException(msg)  # must be iterable
450         return InterpKernelDEC._NewWithComm_internal(s, t, MPI._addressof(mpicomm))
451
452 def _ODEC_WithComm_internal(procs, mpicomm=None):
453     from mpi4py import MPI
454     # Check iterable:
455     try:
456         g = [el for el in procs]
457     except:
458         msg =  "OverlapDEC: invalid type in ctor arguments! Possible signatures are:\n"
459         msg += "   - OverlapDEC.New(<iterable>)\n"
460         msg += "   - OverlapDEC.New(<iterable>, MPI_Comm)\n"
461         msg += "   - OverlapDEC(<iterable>)\n"
462         msg += "   - OverlapDEC(<iterable>, MPI_Comm*) : WARNING here the address of the communicator should be passed with MPI._addressof(the_com)\n"
463         raise InterpKernelException(msg)
464     if mpicomm is None:
465         return OverlapDEC(g)
466     else:
467         return OverlapDEC._NewWithComm_internal(g, MPI._addressof(mpicomm))
468
469 InterpKernelDEC.New = _IKDEC_WithComm_internal
470 OverlapDEC.New = _ODEC_WithComm_internal
471
472 %}