1 // Copyright (C) 2017-2020 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 "ParaMESH.hxx"
34 #include "ParaFIELD.hxx"
35 #include "ICoCoMEDField.hxx"
36 #include "ComponentTopology.hxx"
37 #include "ParaUMesh.hxx"
38 #include "ParaSkyLineArray.hxx"
40 using namespace INTERP_KERNEL;
41 using namespace MEDCoupling;
42 using namespace ICoCo;
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"
53 %include "DisjointDEC.hxx"
54 %include "InterpKernelDEC.hxx"
55 %include "StructuredCoincidentDEC.hxx"
57 %include "ICoCoField.hxx"
58 %rename(ICoCoMEDField) ICoCo::MEDField;
59 %include "ICoCoMEDField.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::ParaSkyLineArray::New;
70 %newobject MEDCoupling::ParaSkyLineArray::equiRedistribute;
71 %newobject MEDCoupling::ParaSkyLineArray::getSkyLineArray;
72 %newobject MEDCoupling::ParaSkyLineArray::getGlobalIdsArray;
74 %feature("unref") ParaSkyLineArray "$this->decrRef();"
75 %feature("unref") ParaUMesh "$this->decrRef();"
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;
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,
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;
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;
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;
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;
141 PyObject *allGatherArrays(const DataArrayIdType *array) const
143 std::vector< MCAuto<DataArrayIdType> > ret;
144 self->allGatherArrays(MPI_COMM_WORLD,array,ret);
145 return convertFromVectorAutoObjToPyObj<DataArrayIdType>(ret,SWIGTITraits<mcIdType>::TI);
148 PyObject *allToAllArrays(PyObject *arrays) const
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);
160 class ParaUMesh : public RefCountObject
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 DataArrayDouble *redistributeCellField(const DataArrayIdType *globalCellIds, const DataArrayDouble *fieldValueToRed) const;
167 DataArrayIdType *redistributeNodeField(const DataArrayIdType *globalCellIds, const DataArrayIdType *fieldValueToRed) const;
168 DataArrayDouble *redistributeNodeField(const DataArrayIdType *globalCellIds, const DataArrayDouble *fieldValueToRed) const;
171 ParaUMesh(MEDCouplingUMesh *mesh, DataArrayIdType *globalCellIds, DataArrayIdType *globalNodeIds)
173 return ParaUMesh::New(mesh,globalCellIds,globalNodeIds);
176 MEDCouplingUMesh *getMesh()
178 MEDCouplingUMesh *ret(self->getMesh());
179 if(ret) ret->incrRef();
183 DataArrayIdType *getGlobalCellIds()
185 DataArrayIdType *ret(self->getGlobalCellIds());
186 if(ret) ret->incrRef();
190 DataArrayIdType *getGlobalNodeIds()
192 DataArrayIdType *ret(self->getGlobalNodeIds());
193 if(ret) ret->incrRef();
197 DataArrayIdType *getCellIdsLyingOnNodes(const DataArrayIdType *globalNodeIds, bool fullyIn) const
199 MCAuto<DataArrayIdType> ret(self->getCellIdsLyingOnNodes(globalNodeIds,fullyIn));
205 class ParaSkyLineArray : public RefCountObject
208 static ParaSkyLineArray *New(MEDCouplingSkyLineArray *ska, DataArrayIdType *globalIds);
211 ParaSkyLineArray(MEDCouplingSkyLineArray *ska, DataArrayIdType *globalIds)
213 return ParaSkyLineArray::New(ska,globalIds);
216 ParaSkyLineArray *equiRedistribute(mcIdType nbOfEntities) const
218 MCAuto<ParaSkyLineArray> ret(self->equiRedistribute(nbOfEntities));
222 MEDCouplingSkyLineArray *getSkyLineArray() const
224 MEDCouplingSkyLineArray *ret(self->getSkyLineArray());
230 DataArrayIdType *getGlobalIdsArray() const
232 DataArrayIdType *ret(self->getGlobalIdsArray());
241 /* This object can be used only if MED_ENABLE_FVM is defined*/
242 #ifdef MED_ENABLE_FVM
243 class NonCoincidentDEC : public DEC
246 NonCoincidentDEC(ProcessorGroup& source, ProcessorGroup& target);
250 %extend MEDCoupling::ParaMESH
252 PyObject *getGlobalNumberingCell2() const
254 const mcIdType *tmp=self->getGlobalNumberingCell();
255 mcIdType size=self->getCellMesh()->getNumberOfCells();
256 PyObject *ret=PyList_New(size);
257 for(mcIdType i=0;i<size;i++)
258 PyList_SetItem(ret,i,PyInt_FromLong(tmp[i]));
262 PyObject *getGlobalNumberingFace2() const
264 const mcIdType *tmp=self->getGlobalNumberingFace();
265 mcIdType size=self->getFaceMesh()->getNumberOfCells();
266 PyObject *ret=PyList_New(size);
267 for(mcIdType i=0;i<size;i++)
268 PyList_SetItem(ret,i,PyInt_FromLong(tmp[i]));
272 PyObject *getGlobalNumberingNode2() const
274 const mcIdType *tmp=self->getGlobalNumberingNode();
275 mcIdType size=self->getCellMesh()->getNumberOfNodes();
276 PyObject *ret=PyList_New(size);
277 for(mcIdType i=0;i<size;i++)
278 PyList_SetItem(ret,i,PyInt_FromLong(tmp[i]));