1 // Copyright (C) 2007-2016 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 #define MEDCOUPLING_EXPORT
23 #define INTERPKERNEL_EXPORT
25 %include "ParaMEDMEM.typemap"
26 %include "MEDCouplingCommon.i"
29 #include "CommInterface.hxx"
30 #include "ProcessorGroup.hxx"
31 #include "Topology.hxx"
32 #include "MPIProcessorGroup.hxx"
34 #include "InterpKernelDEC.hxx"
35 #include "NonCoincidentDEC.hxx"
36 #include "StructuredCoincidentDEC.hxx"
37 #include "ParaMESH.hxx"
38 #include "ParaFIELD.hxx"
39 #include "ICoCoMEDField.hxx"
40 #include "ComponentTopology.hxx"
44 using namespace INTERP_KERNEL;
45 using namespace MEDCoupling;
46 using namespace ICoCo;
48 enum mpi_constants { mpi_comm_world, mpi_comm_self, mpi_double, mpi_int };
51 %include "InterpolationOptions.hxx"
52 %include "CommInterface.hxx"
53 %include "ProcessorGroup.hxx"
54 %include "DECOptions.hxx"
55 %include "ParaMESH.hxx"
56 %include "ParaFIELD.hxx"
57 %include "MPIProcessorGroup.hxx"
58 %include "ComponentTopology.hxx"
60 %include "DisjointDEC.hxx"
61 %include "InterpKernelDEC.hxx"
62 %include "StructuredCoincidentDEC.hxx"
64 %include "ICoCoField.hxx"
65 %rename(ICoCoMEDField) ICoCo::MEDField;
66 %include "ICoCoMEDField.hxx"
70 /* This object can be used only if MED_ENABLE_FVM is defined*/
72 class NonCoincidentDEC : public DEC
75 NonCoincidentDEC(ProcessorGroup& source, ProcessorGroup& target);
79 %extend MEDCoupling::ParaMESH
81 PyObject *getGlobalNumberingCell2() const
83 const int *tmp=self->getGlobalNumberingCell();
84 int size=self->getCellMesh()->getNumberOfCells();
85 PyObject *ret=PyList_New(size);
86 for(int i=0;i<size;i++)
87 PyList_SetItem(ret,i,PyInt_FromLong(tmp[i]));
91 PyObject *getGlobalNumberingFace2() const
93 const int *tmp=self->getGlobalNumberingFace();
94 int size=self->getFaceMesh()->getNumberOfCells();
95 PyObject *ret=PyList_New(size);
96 for(int i=0;i<size;i++)
97 PyList_SetItem(ret,i,PyInt_FromLong(tmp[i]));
101 PyObject *getGlobalNumberingNode2() const
103 const int *tmp=self->getGlobalNumberingNode();
104 int size=self->getCellMesh()->getNumberOfNodes();
105 PyObject *ret=PyList_New(size);
106 for(int i=0;i<size;i++)
107 PyList_SetItem(ret,i,PyInt_FromLong(tmp[i]));
112 //=============================================================================================
113 // Interface for MPI-realization-specific constants like MPI_COMM_WORLD.
115 // Type and values of constants like MPI_COMM_WORLD depends on MPI realization
116 // and usually such constants actually are macros. To have such symbols in python
117 // and translate them into correct values we use the following technique.
118 // We define some constants (enum mpi_constants) and map them into real MPI values
119 // using typemaps, and we create needed python symbols equal to 'mpi_constants'
120 // via %pythoncode directive.
122 // Constants corresponding to similar MPI definitions
123 enum mpi_constants { mpi_comm_world, mpi_comm_self, mpi_double, mpi_int };
125 // Map mpi_comm_world and mpi_comm_self -> MPI_COMM_WORLD and MPI_COMM_SELF
126 %typemap(in) MPI_Comm
128 switch (PyInt_AsLong($input))
130 case mpi_comm_world: $1 = MPI_COMM_WORLD; break;
131 case mpi_comm_self: $1 = MPI_COMM_SELF; break;
133 PyErr_SetString(PyExc_TypeError,"unexpected value of MPI_Comm");
137 // Map mpi_double and mpi_int -> MPI_DOUBLE and MPI_INT
138 %typemap(in) MPI_Datatype
140 switch (PyInt_AsLong($input))
142 case mpi_double: $1 = MPI_DOUBLE; break;
143 case mpi_int: $1 = MPI_INT; break;
145 PyErr_SetString(PyExc_TypeError,"unexpected value of MPI_Datatype");
149 // The following code gets inserted into the result python file:
150 // create needed python symbols
152 MPI_COMM_WORLD = mpi_comm_world
153 MPI_COMM_SELF = mpi_comm_self
154 MPI_DOUBLE = mpi_double
157 //=============================================================================================
162 %inline %{ PyObject* MPI_Comm_size(MPI_Comm comm)
165 int err = MPI_Comm_size(comm, &res);
166 if ( err != MPI_SUCCESS )
168 PyErr_SetString(PyExc_RuntimeError,"Erorr in MPI_Comm_size()");
171 return PyInt_FromLong( res );
177 %inline %{ PyObject* MPI_Comm_rank(MPI_Comm comm)
180 int err = MPI_Comm_rank(comm, &res);
181 if ( err != MPI_SUCCESS )
183 PyErr_SetString(PyExc_RuntimeError,"Erorr in MPI_Comm_rank()");
186 return PyInt_FromLong( res );
190 int MPI_Init(int *argc, char ***argv );
191 int MPI_Barrier(MPI_Comm comm);
198 %inline %{ PyObject* MPI_Bcast(PyObject* buffer, int nb, MPI_Datatype type, int root, MPI_Comm c)
200 // buffer must be a list
201 if (!PyList_Check(buffer))
203 PyErr_SetString(PyExc_TypeError, "buffer is expected to be a list");
207 int aSize = PyList_Size(buffer);
210 std::ostringstream stream; stream << "buffer is expected to be of size " << nb;
211 PyErr_SetString(PyExc_ValueError, stream.str().c_str());
214 // allocate and fill a buffer
218 if ( type == MPI_DOUBLE )
220 aBuf = (void*) ( dblBuf = new double[ nb ] );
221 for ( int i = 0; i < aSize; ++i )
222 dblBuf[i] = PyFloat_AS_DOUBLE( PyList_GetItem( buffer, i ));
224 else if ( type == MPI_INT )
226 aBuf = (void*) ( intBuf = new int[ nb ] );
227 for ( int i = 0; i < aSize; ++i )
228 intBuf[i] = int( PyInt_AS_LONG( PyList_GetItem( buffer, i )));
232 PyErr_SetString(PyExc_TypeError, "Only MPI_DOUBLE and MPI_INT supported");
236 int err = MPI_Bcast(aBuf, nb, type, root, c);
238 if ( err != MPI_SUCCESS )
240 PyErr_SetString(PyExc_RuntimeError,"Erorr in MPI_Bcast()");
241 delete [] intBuf; delete [] dblBuf;
244 // put received data into the list
246 if ( type == MPI_DOUBLE )
248 for ( int i = 0; i < aSize && !pyerr; ++i )
249 pyerr = PyList_SetItem(buffer, i, PyFloat_FromDouble( dblBuf[i] ));
254 for ( int i = 0; i < aSize && !pyerr; ++i )
255 pyerr = PyList_SetItem(buffer, i, PyInt_FromLong( intBuf[i] ));
260 PyErr_SetString(PyExc_RuntimeError, "Error of PyList_SetItem()");
263 return PyInt_FromLong( err );
269 def MEDCouplingDataArrayDoubleIadd(self,*args):
271 return _ParaMEDMEM.DataArrayDouble____iadd___(self, self, *args)
272 def MEDCouplingDataArrayDoubleIsub(self,*args):
274 return _ParaMEDMEM.DataArrayDouble____isub___(self, self, *args)
275 def MEDCouplingDataArrayDoubleImul(self,*args):
277 return _ParaMEDMEM.DataArrayDouble____imul___(self, self, *args)
278 def MEDCouplingDataArrayDoubleIdiv(self,*args):
280 return _ParaMEDMEM.DataArrayDouble____idiv___(self, self, *args)
281 def MEDCouplingDataArrayDoubleIpow(self,*args):
283 return _ParaMEDMEM.DataArrayDouble____ipow___(self, self, *args)
284 def MEDCouplingFieldDoubleIadd(self,*args):
286 return _ParaMEDMEM.MEDCouplingFieldDouble____iadd___(self, self, *args)
287 def MEDCouplingFieldDoubleIsub(self,*args):
289 return _ParaMEDMEM.MEDCouplingFieldDouble____isub___(self, self, *args)
290 def MEDCouplingFieldDoubleImul(self,*args):
292 return _ParaMEDMEM.MEDCouplingFieldDouble____imul___(self, self, *args)
293 def MEDCouplingFieldDoubleIdiv(self,*args):
295 return _ParaMEDMEM.MEDCouplingFieldDouble____idiv___(self, self, *args)
296 def MEDCouplingFieldDoubleIpow(self,*args):
298 return _ParaMEDMEM.MEDCouplingFieldDouble____ipow___(self, self, *args)
299 def MEDCouplingDataArrayIntIadd(self,*args):
301 return _ParaMEDMEM.DataArrayInt____iadd___(self, self, *args)
302 def MEDCouplingDataArrayIntIsub(self,*args):
304 return _ParaMEDMEM.DataArrayInt____isub___(self, self, *args)
305 def MEDCouplingDataArrayIntImul(self,*args):
307 return _ParaMEDMEM.DataArrayInt____imul___(self, self, *args)
308 def MEDCouplingDataArrayIntIdiv(self,*args):
310 return _ParaMEDMEM.DataArrayInt____idiv___(self, self, *args)
311 def MEDCouplingDataArrayIntImod(self,*args):
313 return _ParaMEDMEM.DataArrayInt____imod___(self, self, *args)
314 def MEDCouplingDataArrayIntIpow(self,*args):
316 return _ParaMEDMEM.DataArrayInt____ipow___(self, self, *args)
317 def MEDCouplingDataArrayFloatIadd(self,*args):
319 return _ParaMEDMEM.DataArrayFloat____iadd___(self, self, *args)
320 def MEDCouplingDataArrayFloatIsub(self,*args):
322 return _ParaMEDMEM.DataArrayFloat____isub___(self, self, *args)
323 def MEDCouplingDataArrayFloatImul(self,*args):
325 return _ParaMEDMEM.DataArrayFloat____imul___(self, self, *args)
326 def MEDCouplingDataArrayFloatIdiv(self,*args):
328 return _ParaMEDMEM.DataArrayFloat____idiv___(self, self, *args)
329 def MEDCouplingDataArrayDoubleTupleIadd(self,*args):
331 return _ParaMEDMEM.DataArrayDoubleTuple____iadd___(self, self, *args)
332 def MEDCouplingDataArrayDoubleTupleIsub(self,*args):
334 return _ParaMEDMEM.DataArrayDoubleTuple____isub___(self, self, *args)
335 def MEDCouplingDataArrayDoubleTupleImul(self,*args):
337 return _ParaMEDMEM.DataArrayDoubleTuple____imul___(self, self, *args)
338 def MEDCouplingDataArrayDoubleTupleIdiv(self,*args):
340 return _ParaMEDMEM.DataArrayDoubleTuple____idiv___(self, self, *args)
341 def MEDCouplingDataArrayIntTupleIadd(self,*args):
343 return _ParaMEDMEM.DataArrayIntTuple____iadd___(self, self, *args)
344 def MEDCouplingDataArrayIntTupleIsub(self,*args):
346 return _ParaMEDMEM.DataArrayIntTuple____isub___(self, self, *args)
347 def MEDCouplingDataArrayIntTupleImul(self,*args):
349 return _ParaMEDMEM.DataArrayIntTuple____imul___(self, self, *args)
350 def MEDCouplingDataArrayIntTupleIdiv(self,*args):
352 return _ParaMEDMEM.DataArrayIntTuple____idiv___(self, self, *args)
353 def MEDCouplingDataArrayIntTupleImod(self,*args):
355 return _ParaMEDMEM.DataArrayIntTuple____imod___(self, self, *args)
356 def ParaMEDMEMDenseMatrixIadd(self,*args):
358 return _ParaMEDMEM.DenseMatrix____iadd___(self, self, *args)
359 def ParaMEDMEMDenseMatrixIsub(self,*args):
361 return _ParaMEDMEM.DenseMatrix____isub___(self, self, *args)
364 %include "MEDCouplingFinalize.i"