Salome HOME
Still some imp into ParaUMesh
[tools/medcoupling.git] / src / ParaMEDMEM_Swig / ParaMEDMEMCommon.i
1 // Copyright (C) 2017-2020  CEA/DEN, EDF R&D
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 "DEC.hxx"
30 #include "InterpKernelDEC.hxx"
31 #include "NonCoincidentDEC.hxx"
32 #include "StructuredCoincidentDEC.hxx"
33 #include "ParaMESH.hxx"
34 #include "ParaFIELD.hxx"
35 #include "ICoCoMEDField.hxx"
36 #include "ComponentTopology.hxx"
37 #include "ParaUMesh.hxx"
38 #include "ParaSkyLineArray.hxx"
39
40 using namespace INTERP_KERNEL;
41 using namespace MEDCoupling;
42 using namespace ICoCo;
43 %}
44
45 %include "InterpolationOptions.hxx"
46 %include "ProcessorGroup.hxx"
47 %include "DECOptions.hxx"
48 %include "ParaMESH.hxx"
49 %include "ParaFIELD.hxx"
50 %include "MPIProcessorGroup.hxx"
51 %include "ComponentTopology.hxx"
52 %include "DEC.hxx"
53 %include "DisjointDEC.hxx"
54 %include "InterpKernelDEC.hxx"
55 %include "StructuredCoincidentDEC.hxx"
56
57 %include "ICoCoField.hxx"
58 %rename(ICoCoMEDField) ICoCo::MEDField;
59 %include "ICoCoMEDField.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::ParaSkyLineArray::New;
70 %newobject MEDCoupling::ParaSkyLineArray::equiRedistribute;
71 %newobject MEDCoupling::ParaSkyLineArray::getSkyLineArray;
72 %newobject MEDCoupling::ParaSkyLineArray::getGlobalIdsArray;
73
74 %feature("unref") ParaSkyLineArray "$this->decrRef();"
75 %feature("unref") ParaUMesh "$this->decrRef();"
76
77 %nodefaultctor;
78
79 namespace MEDCoupling
80 {
81   class CommInterface
82   {
83   public:
84     CommInterface();
85     virtual ~CommInterface();
86     int worldSize() const;
87     int commSize(MPI_Comm comm, int* size) const;
88     int commRank(MPI_Comm comm, int* rank) const;
89     int commGroup(MPI_Comm comm, MPI_Group* group) const;
90     int groupIncl(MPI_Group group, int size, int* ranks, MPI_Group* group_output) const;
91     int commCreate(MPI_Comm comm, MPI_Group group, MPI_Comm* comm_output) const;
92     int groupFree(MPI_Group* group) const;
93     int commFree(MPI_Comm* comm) const;
94
95     int send(void* buffer, int count, MPI_Datatype datatype, int target, int tag, MPI_Comm comm) const;
96     int recv(void* buffer, int count, MPI_Datatype datatype, int source, int tag, MPI_Comm comm, MPI_Status* status) const;
97     int sendRecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, 
98                  int dest, int sendtag, void* recvbuf, int recvcount, 
99                  MPI_Datatype recvtype, int source, int recvtag, MPI_Comm comm,
100                  MPI_Status* status);
101
102     int Isend(void* buffer, int count, MPI_Datatype datatype, int target,
103               int tag, MPI_Comm comm, MPI_Request *request) const;
104     int Irecv(void* buffer, int count, MPI_Datatype datatype, int source,
105               int tag, MPI_Comm comm, MPI_Request* request) const;
106
107     int wait(MPI_Request *request, MPI_Status *status) const;
108     int test(MPI_Request *request, int *flag, MPI_Status *status) const;
109     int requestFree(MPI_Request *request) const;
110     int waitany(int count, MPI_Request *array_of_requests, int *index, MPI_Status *status) const;
111     int testany(int count, MPI_Request *array_of_requests, int *index, int *flag, MPI_Status *status) const;
112     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); }
113     int testall(int count, MPI_Request *array_of_requests, int *flag, MPI_Status *array_of_status) const;
114     int waitsome(int incount, MPI_Request *array_of_requests,int *outcount, int *array_of_indices, MPI_Status *array_of_status) const;
115     int testsome(int incount, MPI_Request *array_of_requests, int *outcount,
116                  int *array_of_indices, MPI_Status *array_of_status) const;
117     int probe(int source, int tag, MPI_Comm comm, MPI_Status *status) const;
118     int Iprobe(int source, int tag, MPI_Comm comm, int *flag, MPI_Status *status) const;
119     int cancel(MPI_Request *request) const;
120     int testCancelled(MPI_Status *status, int *flag) const;
121     int barrier(MPI_Comm comm) const;
122     int errorString(int errorcode, char *string, int *resultlen) const;
123     int getCount(MPI_Status *status, MPI_Datatype datatype, int *count) const;
124
125     int broadcast(void* buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm) const;
126     int allGather(void* sendbuf, int sendcount, MPI_Datatype sendtype,
127                   void* recvbuf, int recvcount, MPI_Datatype recvtype,
128                   MPI_Comm comm) const;
129     int allToAll(void* sendbuf, int sendcount, MPI_Datatype sendtype,
130                  void* recvbuf, int recvcount, MPI_Datatype recvtype,
131                  MPI_Comm comm) const;
132     int allToAllV(void* sendbuf, int* sendcounts, int* senddispls,
133                   MPI_Datatype sendtype, void* recvbuf, int* recvcounts,
134                   int* recvdispls, MPI_Datatype recvtype, 
135                   MPI_Comm comm) const;
136
137     int reduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm) const;
138     int allReduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) const;
139     %extend
140     {
141       PyObject *allGatherArrays(const DataArrayIdType *array) const
142       {
143         std::vector< MCAuto<DataArrayIdType> > ret;
144         self->allGatherArrays(MPI_COMM_WORLD,array,ret);
145         return convertFromVectorAutoObjToPyObj<DataArrayIdType>(ret,SWIGTITraits<mcIdType>::TI);
146       }
147
148       PyObject *allToAllArrays(PyObject *arrays) const
149       {
150         std::vector< DataArrayIdType * > arraysIn;
151         std::vector< MCAuto<DataArrayIdType> > arrayOut;
152         convertFromPyObjVectorOfObj<MEDCoupling::DataArrayIdType*>(arrays,SWIGTITraits<mcIdType>::TI,"DataArrayIdType",arraysIn);
153         std::vector< MCAuto<DataArrayIdType> > arraysIn2(FromVecToVecAuto<DataArrayIdType>(arraysIn));
154         self->allToAllArrays(MPI_COMM_WORLD,arraysIn2,arrayOut);
155         return convertFromVectorAutoObjToPyObj<DataArrayIdType>(arrayOut,SWIGTITraits<mcIdType>::TI);
156       }
157     }
158   };
159
160   class ParaUMesh : public RefCountObject
161   {
162   public:
163     static ParaUMesh *New(MEDCouplingUMesh *mesh, DataArrayIdType *globalCellIds, DataArrayIdType *globalNodeIds);
164     ParaUMesh *redistributeCells(const DataArrayIdType *globalCellIds) const;
165     DataArrayIdType *redistributeCellField(const DataArrayIdType *globalCellIds, const DataArrayIdType *fieldValueToRed) const;
166     DataArrayIdType *redistributeNodeField(const DataArrayIdType *globalCellIds, const DataArrayIdType *fieldValueToRed) const;
167     %extend
168     {
169       ParaUMesh(MEDCouplingUMesh *mesh, DataArrayIdType *globalCellIds, DataArrayIdType *globalNodeIds)
170       {
171         return ParaUMesh::New(mesh,globalCellIds,globalNodeIds);
172       }
173
174       MEDCouplingUMesh *getMesh()
175       {
176         MEDCouplingUMesh *ret(self->getMesh());
177         if(ret) ret->incrRef();
178         return ret;
179       }
180
181       DataArrayIdType *getGlobalCellIds()
182       {
183         DataArrayIdType *ret(self->getGlobalCellIds());
184         if(ret) ret->incrRef();
185         return ret;
186       }
187
188       DataArrayIdType *getGlobalNodeIds()
189       {
190         DataArrayIdType *ret(self->getGlobalNodeIds());
191         if(ret) ret->incrRef();
192         return ret;
193       }
194
195       DataArrayIdType *getCellIdsLyingOnNodes(const DataArrayIdType *globalNodeIds, bool fullyIn) const
196       { 
197         MCAuto<DataArrayIdType> ret(self->getCellIdsLyingOnNodes(globalNodeIds,fullyIn));
198         return ret.retn();
199       }
200     }
201   };
202
203   class ParaSkyLineArray : public RefCountObject
204   {
205   public:
206     static ParaSkyLineArray *New(MEDCouplingSkyLineArray *ska, DataArrayIdType *globalIds);
207     %extend
208     {
209       ParaSkyLineArray(MEDCouplingSkyLineArray *ska, DataArrayIdType *globalIds)
210       {
211         return ParaSkyLineArray::New(ska,globalIds);
212       }
213
214       ParaSkyLineArray *equiRedistribute(mcIdType nbOfEntities) const
215       {
216         MCAuto<ParaSkyLineArray> ret(self->equiRedistribute(nbOfEntities));
217         return ret.retn();
218       }
219
220       MEDCouplingSkyLineArray *getSkyLineArray() const
221       {
222         MEDCouplingSkyLineArray *ret(self->getSkyLineArray());
223         if(ret)
224           ret->incrRef();
225         return ret;
226       }
227       
228       DataArrayIdType *getGlobalIdsArray() const
229       {
230         DataArrayIdType *ret(self->getGlobalIdsArray());
231         if(ret)
232           ret->incrRef();
233         return ret;
234       }
235     }
236   };
237 }
238
239 /* This object can be used only if MED_ENABLE_FVM is defined*/
240 #ifdef MED_ENABLE_FVM
241 class NonCoincidentDEC : public DEC
242 {
243 public:
244   NonCoincidentDEC(ProcessorGroup& source, ProcessorGroup& target);
245 };
246 #endif
247
248 %extend MEDCoupling::ParaMESH
249 {
250   PyObject *getGlobalNumberingCell2() const
251   {
252     const mcIdType *tmp=self->getGlobalNumberingCell();
253     mcIdType size=self->getCellMesh()->getNumberOfCells();
254     PyObject *ret=PyList_New(size);
255     for(mcIdType i=0;i<size;i++)
256       PyList_SetItem(ret,i,PyInt_FromLong(tmp[i])); 
257     return ret;
258   }
259
260   PyObject *getGlobalNumberingFace2() const
261   {
262     const mcIdType *tmp=self->getGlobalNumberingFace();
263     mcIdType size=self->getFaceMesh()->getNumberOfCells();
264     PyObject *ret=PyList_New(size);
265     for(mcIdType i=0;i<size;i++)
266       PyList_SetItem(ret,i,PyInt_FromLong(tmp[i])); 
267     return ret;
268   }
269
270   PyObject *getGlobalNumberingNode2() const
271   {
272     const mcIdType *tmp=self->getGlobalNumberingNode();
273     mcIdType size=self->getCellMesh()->getNumberOfNodes();
274     PyObject *ret=PyList_New(size);
275     for(mcIdType i=0;i<size;i++)
276       PyList_SetItem(ret,i,PyInt_FromLong(tmp[i])); 
277     return ret;
278   }
279 }