1 // Copyright (C) 2017-2021 CEA/DEN, EDF R&D
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.
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.
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
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
22 %template() std::set<int>;
25 #include "CommInterface.hxx"
26 #include "ProcessorGroup.hxx"
27 #include "Topology.hxx"
28 #include "MPIProcessorGroup.hxx"
30 #include "InterpKernelDEC.hxx"
31 #include "NonCoincidentDEC.hxx"
32 #include "StructuredCoincidentDEC.hxx"
33 #include "OverlapDEC.hxx"
34 #include "ParaMESH.hxx"
35 #include "ParaFIELD.hxx"
36 #include "ICoCoMEDDoubleField.hxx"
37 #include "ICoCoMEDIntField.hxx"
38 #include "ComponentTopology.hxx"
39 #include "ParaUMesh.hxx"
40 #include "ParaSkyLineArray.hxx"
41 #include "ParaDataArray.hxx"
43 using namespace INTERP_KERNEL;
44 using namespace MEDCoupling;
45 using namespace ICoCo;
48 %include "InterpolationOptions.hxx"
49 %include "ProcessorGroup.hxx"
50 %include "DECOptions.hxx"
51 %include "ParaMESH.hxx"
52 %include "ParaFIELD.hxx"
53 %include "MPIProcessorGroup.hxx"
54 %include "ComponentTopology.hxx"
56 %include "DisjointDEC.hxx"
57 %include "InterpKernelDEC.hxx"
58 %include "StructuredCoincidentDEC.hxx"
59 %include "OverlapDEC.hxx"
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;
78 %feature("unref") ParaSkyLineArray "$this->decrRef();"
79 %feature("unref") ParaUMesh "$this->decrRef();"
80 %feature("unref") ParaDataArrayInt32 "$this->decrRef();"
81 %feature("unref") ParaDataArrayInt64 "$this->decrRef();"
91 virtual ~CommInterface();
92 int worldSize() const;
93 int commSize(MPI_Comm comm, int* size) const;
94 int commRank(MPI_Comm comm, int* rank) const;
95 int commGroup(MPI_Comm comm, MPI_Group* group) const;
96 int groupIncl(MPI_Group group, int size, int* ranks, MPI_Group* group_output) const;
97 int commCreate(MPI_Comm comm, MPI_Group group, MPI_Comm* comm_output) const;
98 int groupFree(MPI_Group* group) const;
99 int commFree(MPI_Comm* comm) const;
101 int send(void* buffer, int count, MPI_Datatype datatype, int target, int tag, MPI_Comm comm) const;
102 int recv(void* buffer, int count, MPI_Datatype datatype, int source, int tag, MPI_Comm comm, MPI_Status* status) const;
103 int sendRecv(void* sendbuf, int sendcount, MPI_Datatype sendtype,
104 int dest, int sendtag, void* recvbuf, int recvcount,
105 MPI_Datatype recvtype, int source, int recvtag, MPI_Comm comm,
108 int Isend(void* buffer, int count, MPI_Datatype datatype, int target,
109 int tag, MPI_Comm comm, MPI_Request *request) const;
110 int Irecv(void* buffer, int count, MPI_Datatype datatype, int source,
111 int tag, MPI_Comm comm, MPI_Request* request) const;
113 int wait(MPI_Request *request, MPI_Status *status) const;
114 int test(MPI_Request *request, int *flag, MPI_Status *status) const;
115 int requestFree(MPI_Request *request) const;
116 int waitany(int count, MPI_Request *array_of_requests, int *index, MPI_Status *status) const;
117 int testany(int count, MPI_Request *array_of_requests, int *index, int *flag, MPI_Status *status) const;
118 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); }
119 int testall(int count, MPI_Request *array_of_requests, int *flag, MPI_Status *array_of_status) const;
120 int waitsome(int incount, MPI_Request *array_of_requests,int *outcount, int *array_of_indices, MPI_Status *array_of_status) const;
121 int testsome(int incount, MPI_Request *array_of_requests, int *outcount,
122 int *array_of_indices, MPI_Status *array_of_status) const;
123 int probe(int source, int tag, MPI_Comm comm, MPI_Status *status) const;
124 int Iprobe(int source, int tag, MPI_Comm comm, int *flag, MPI_Status *status) const;
125 int cancel(MPI_Request *request) const;
126 int testCancelled(MPI_Status *status, int *flag) const;
127 int barrier(MPI_Comm comm) const;
128 int errorString(int errorcode, char *string, int *resultlen) const;
129 int getCount(MPI_Status *status, MPI_Datatype datatype, int *count) const;
131 int broadcast(void* buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm) const;
132 int allGather(void* sendbuf, int sendcount, MPI_Datatype sendtype,
133 void* recvbuf, int recvcount, MPI_Datatype recvtype,
134 MPI_Comm comm) const;
135 int allToAll(void* sendbuf, int sendcount, MPI_Datatype sendtype,
136 void* recvbuf, int recvcount, MPI_Datatype recvtype,
137 MPI_Comm comm) const;
138 int allToAllV(void* sendbuf, int* sendcounts, int* senddispls,
139 MPI_Datatype sendtype, void* recvbuf, int* recvcounts,
140 int* recvdispls, MPI_Datatype recvtype,
141 MPI_Comm comm) const;
143 int reduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm) const;
144 int allReduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) const;
147 PyObject *allGatherArrays(const DataArrayIdType *array) const
149 std::vector< MCAuto<DataArrayIdType> > ret;
150 self->allGatherArrays(MPI_COMM_WORLD,array,ret);
151 return convertFromVectorAutoObjToPyObj<DataArrayIdType>(ret,SWIGTITraits<mcIdType>::TI);
154 PyObject *allToAllArrays(PyObject *arrays) const
156 std::vector< DataArrayIdType * > arraysIn;
157 std::vector< MCAuto<DataArrayIdType> > arrayOut;
158 convertFromPyObjVectorOfObj<MEDCoupling::DataArrayIdType*>(arrays,SWIGTITraits<mcIdType>::TI,"DataArrayIdType",arraysIn);
159 std::vector< MCAuto<DataArrayIdType> > arraysIn2(FromVecToVecAuto<DataArrayIdType>(arraysIn));
160 self->allToAllArrays(MPI_COMM_WORLD,arraysIn2,arrayOut);
161 return convertFromVectorAutoObjToPyObj<DataArrayIdType>(arrayOut,SWIGTITraits<mcIdType>::TI);
166 class ParaUMesh : public RefCountObject
169 static ParaUMesh *New(MEDCouplingUMesh *mesh, DataArrayIdType *globalCellIds, DataArrayIdType *globalNodeIds);
170 ParaUMesh *redistributeCells(const DataArrayIdType *globalCellIds) const;
171 DataArrayIdType *redistributeCellField(const DataArrayIdType *globalCellIds, const DataArrayIdType *fieldValueToRed) const;
172 DataArrayDouble *redistributeCellField(const DataArrayIdType *globalCellIds, const DataArrayDouble *fieldValueToRed) const;
173 DataArrayIdType *redistributeNodeField(const DataArrayIdType *globalCellIds, const DataArrayIdType *fieldValueToRed) const;
174 DataArrayDouble *redistributeNodeField(const DataArrayIdType *globalCellIds, const DataArrayDouble *fieldValueToRed) const;
177 ParaUMesh(MEDCouplingUMesh *mesh, DataArrayIdType *globalCellIds, DataArrayIdType *globalNodeIds)
179 return ParaUMesh::New(mesh,globalCellIds,globalNodeIds);
182 MEDCouplingUMesh *getMesh()
184 MEDCouplingUMesh *ret(self->getMesh());
185 if(ret) ret->incrRef();
189 DataArrayIdType *getGlobalCellIds()
191 DataArrayIdType *ret(self->getGlobalCellIds());
192 if(ret) ret->incrRef();
196 DataArrayIdType *getGlobalNodeIds()
198 DataArrayIdType *ret(self->getGlobalNodeIds());
199 if(ret) ret->incrRef();
203 DataArrayIdType *getCellIdsLyingOnNodes(const DataArrayIdType *globalNodeIds, bool fullyIn) const
205 MCAuto<DataArrayIdType> ret(self->getCellIdsLyingOnNodes(globalNodeIds,fullyIn));
211 class ParaDataArray : public RefCountObject
215 class ParaDataArrayInt32 : public ParaDataArray
218 static ParaDataArrayInt32 *New(DataArrayInt32 *seqDa);
219 DataArrayIdType *buildComplement(int nbOfElems) const;
222 ParaDataArrayInt32(DataArrayInt32 *seqDa)
224 return ParaDataArrayInt32::New(seqDa);
229 class ParaDataArrayInt64 : public ParaDataArray
232 static ParaDataArrayInt64 *New(DataArrayInt64 *seqDa);
233 DataArrayIdType *buildComplement(long nbOfElems) const;
236 ParaDataArrayInt64(DataArrayInt64 *seqDa)
238 return ParaDataArrayInt64::New(seqDa);
243 class ParaSkyLineArray : public RefCountObject
246 static ParaSkyLineArray *New(MEDCouplingSkyLineArray *ska, DataArrayIdType *globalIds);
249 ParaSkyLineArray(MEDCouplingSkyLineArray *ska, DataArrayIdType *globalIds)
251 return ParaSkyLineArray::New(ska,globalIds);
254 ParaSkyLineArray *equiRedistribute(mcIdType nbOfEntities) const
256 MCAuto<ParaSkyLineArray> ret(self->equiRedistribute(nbOfEntities));
260 MEDCouplingSkyLineArray *getSkyLineArray() const
262 MEDCouplingSkyLineArray *ret(self->getSkyLineArray());
268 DataArrayIdType *getGlobalIdsArray() const
270 DataArrayIdType *ret(self->getGlobalIdsArray());
279 /* This object can be used only if MED_ENABLE_FVM is defined*/
280 #ifdef MED_ENABLE_FVM
281 class NonCoincidentDEC : public DEC
284 NonCoincidentDEC(ProcessorGroup& source, ProcessorGroup& target);
288 %extend MEDCoupling::ParaMESH
290 PyObject *getGlobalNumberingCell2() const
292 const mcIdType *tmp=self->getGlobalNumberingCell();
293 mcIdType size=self->getCellMesh()->getNumberOfCells();
294 PyObject *ret=PyList_New(size);
295 for(mcIdType i=0;i<size;i++)
296 PyList_SetItem(ret,i,PyInt_FromLong(tmp[i]));
300 PyObject *getGlobalNumberingFace2() const
302 const mcIdType *tmp=self->getGlobalNumberingFace();
303 mcIdType size=self->getFaceMesh()->getNumberOfCells();
304 PyObject *ret=PyList_New(size);
305 for(mcIdType i=0;i<size;i++)
306 PyList_SetItem(ret,i,PyInt_FromLong(tmp[i]));
310 PyObject *getGlobalNumberingNode2() const
312 const mcIdType *tmp=self->getGlobalNumberingNode();
313 mcIdType size=self->getCellMesh()->getNumberOfNodes();
314 PyObject *ret=PyList_New(size);
315 for(mcIdType i=0;i<size;i++)
316 PyList_SetItem(ret,i,PyInt_FromLong(tmp[i]));
322 if MEDCouplingUse64BitIDs():
323 ParaDataArrayInt = ParaDataArrayInt64
325 ParaDataArrayInt = ParaDataArrayInt32