Salome HOME
Move medtool folder to MED base repository
[modules/med.git] / medtool / src / ParaMEDMEM_Swig / ParaMEDMEM.i
1 // Copyright (C) 2007-2015  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 %module ParaMEDMEM
21
22 %include "ParaMEDMEM.typemap"
23 %include "MEDLoaderCommon.i"
24
25 %{
26 #include "CommInterface.hxx"
27 #include "ProcessorGroup.hxx"
28 #include "Topology.hxx"
29 #include "MPIProcessorGroup.hxx"
30 #include "DEC.hxx"
31 #include "InterpKernelDEC.hxx"
32 #include "NonCoincidentDEC.hxx"
33 #include "StructuredCoincidentDEC.hxx"
34 #include "ParaMESH.hxx"
35 #include "ParaFIELD.hxx"
36 #include "ICoCoMEDField.hxx"
37 #include "ComponentTopology.hxx"
38
39 #include <mpi.h>
40
41 using namespace ParaMEDMEM;
42 using namespace ICoCo;
43       
44 enum mpi_constants { mpi_comm_world, mpi_comm_self, mpi_double, mpi_int };
45 %}
46
47 %include "CommInterface.hxx"
48 %include "ProcessorGroup.hxx"
49 %include "DECOptions.hxx"
50 %include "ParaMESH.hxx"
51 %include "ParaFIELD.hxx"
52 %include "MPIProcessorGroup.hxx"
53 %include "ComponentTopology.hxx"
54 %include "DEC.hxx"
55 %include "InterpKernelDEC.hxx"
56 %include "StructuredCoincidentDEC.hxx"
57
58 %rename(ICoCoMEDField) ICoCo::MEDField;
59 %include "ICoCoMEDField.hxx"
60
61 %nodefaultctor;
62
63 /* This object can be used only if MED_ENABLE_FVM is defined*/
64 #ifdef MED_ENABLE_FVM
65 class NonCoincidentDEC : public DEC
66 {
67 public:
68   NonCoincidentDEC(ProcessorGroup& source, ProcessorGroup& target);
69 };
70 #endif
71
72 %extend ParaMEDMEM::ParaMESH
73 {
74   PyObject *getGlobalNumberingCell2() const
75   {
76     const int *tmp=self->getGlobalNumberingCell();
77     int size=self->getCellMesh()->getNumberOfCells();
78     PyObject *ret=PyList_New(size);
79     for(int i=0;i<size;i++)
80       PyList_SetItem(ret,i,PyInt_FromLong(tmp[i])); 
81     return ret;
82   }
83
84   PyObject *getGlobalNumberingFace2() const
85   {
86     const int *tmp=self->getGlobalNumberingFace();
87     int size=self->getFaceMesh()->getNumberOfCells();
88     PyObject *ret=PyList_New(size);
89     for(int i=0;i<size;i++)
90       PyList_SetItem(ret,i,PyInt_FromLong(tmp[i])); 
91     return ret;
92   }
93
94   PyObject *getGlobalNumberingNode2() const
95   {
96     const int *tmp=self->getGlobalNumberingNode();
97     int size=self->getCellMesh()->getNumberOfNodes();
98     PyObject *ret=PyList_New(size);
99     for(int i=0;i<size;i++)
100       PyList_SetItem(ret,i,PyInt_FromLong(tmp[i])); 
101     return ret;
102   }
103 }
104
105 //=============================================================================================
106 // Interface for MPI-realization-specific constants like MPI_COMM_WORLD.
107 //
108 // Type and values of constants like MPI_COMM_WORLD depends on MPI realization
109 // and usually such constants actually are macros. To have such symbols in python
110 // and translate them into correct values we use the following technique.
111 // We define some constants (enum mpi_constants) and map them into real MPI values
112 // using typemaps, and we create needed python symbols equal to 'mpi_constants'
113 // via %pythoncode directive.
114
115 // Constants corresponding to similar MPI definitions
116 enum mpi_constants { mpi_comm_world, mpi_comm_self, mpi_double, mpi_int };
117
118 // Map mpi_comm_world and mpi_comm_self -> MPI_COMM_WORLD and MPI_COMM_SELF
119 %typemap(in) MPI_Comm
120
121   switch (PyInt_AsLong($input))
122     {
123     case mpi_comm_world: $1 = MPI_COMM_WORLD; break;
124     case mpi_comm_self:  $1 = MPI_COMM_SELF;  break;
125     default:
126       PyErr_SetString(PyExc_TypeError,"unexpected value of MPI_Comm");
127       return NULL;
128     }
129 }
130 // Map mpi_double and mpi_int -> MPI_DOUBLE and MPI_INT
131 %typemap(in) MPI_Datatype
132
133   switch (PyInt_AsLong($input))
134     {
135     case mpi_double:     $1 = MPI_DOUBLE;     break;
136     case mpi_int:        $1 = MPI_INT;        break;
137     default:
138       PyErr_SetString(PyExc_TypeError,"unexpected value of MPI_Datatype");
139       return NULL;
140     }
141 }
142 // The following code gets inserted into the result python file:
143 // create needed python symbols
144 %pythoncode %{
145 MPI_COMM_WORLD = mpi_comm_world
146 MPI_COMM_SELF  = mpi_comm_self
147 MPI_DOUBLE     = mpi_double
148 MPI_INT        = mpi_int
149 %}
150 //=============================================================================================
151
152 // ==============
153 // MPI_Comm_size
154 // ==============
155 %inline %{ PyObject* MPI_Comm_size(MPI_Comm comm)
156   {
157     int res = 0;
158     int err = MPI_Comm_size(comm, &res);
159     if ( err != MPI_SUCCESS )
160       {
161         PyErr_SetString(PyExc_RuntimeError,"Erorr in MPI_Comm_size()");
162         return NULL;
163       }
164     return PyInt_FromLong( res );
165   } %}
166
167 // ==============
168 // MPI_Comm_rank
169 // ==============
170 %inline %{ PyObject* MPI_Comm_rank(MPI_Comm comm)
171   {
172     int res = 0;
173     int err = MPI_Comm_rank(comm, &res);
174     if ( err != MPI_SUCCESS )
175       {
176         PyErr_SetString(PyExc_RuntimeError,"Erorr in MPI_Comm_rank()");
177         return NULL;
178       }
179     return PyInt_FromLong( res );
180   } 
181   %}
182
183 int MPI_Init(int *argc, char ***argv );
184 int MPI_Barrier(MPI_Comm comm);
185 int MPI_Finalize();
186
187 // ==========
188 // MPI_Bcast
189 // ==========
190
191 %inline %{ PyObject* MPI_Bcast(PyObject* buffer, int nb, MPI_Datatype type, int root, MPI_Comm c)
192   {
193     // buffer must be a list
194     if (!PyList_Check(buffer))
195       {
196         PyErr_SetString(PyExc_TypeError, "buffer is expected to be a list");
197         return NULL;
198       }
199     // check list size
200     int aSize = PyList_Size(buffer);
201     if ( aSize != nb )
202       {
203         std::ostringstream stream; stream << "buffer is expected to be of size " << nb;
204         PyErr_SetString(PyExc_ValueError, stream.str().c_str());
205         return NULL;
206       }
207     // allocate and fill a buffer
208     void* aBuf = 0;
209     int* intBuf = 0;
210     double* dblBuf = 0;
211     if ( type == MPI_DOUBLE )
212       {
213         aBuf = (void*) ( dblBuf = new double[ nb ] );
214         for ( int i = 0; i < aSize; ++i )
215           dblBuf[i] = PyFloat_AS_DOUBLE( PyList_GetItem( buffer, i ));
216       }
217     else if ( type == MPI_INT )
218       {
219         aBuf = (void*) ( intBuf = new int[ nb ] );
220         for ( int i = 0; i < aSize; ++i )
221           intBuf[i] = int( PyInt_AS_LONG( PyList_GetItem( buffer, i )));
222       }
223     else
224       {
225         PyErr_SetString(PyExc_TypeError, "Only MPI_DOUBLE and MPI_INT supported");
226         return NULL;
227       }
228     // call MPI_Bcast
229     int err = MPI_Bcast(aBuf, nb, type, root, c);
230     // treat error
231     if ( err != MPI_SUCCESS )
232       {
233         PyErr_SetString(PyExc_RuntimeError,"Erorr in MPI_Bcast()");
234         delete [] intBuf; delete [] dblBuf;
235         return NULL;
236       }
237     // put recieved data into the list
238     int pyerr = 0;
239     if ( type == MPI_DOUBLE )
240       {
241         for ( int i = 0; i < aSize && !pyerr; ++i )
242           pyerr = PyList_SetItem(buffer, i, PyFloat_FromDouble( dblBuf[i] ));
243         delete [] dblBuf;
244       }
245     else
246       {
247         for ( int i = 0; i < aSize && !pyerr; ++i )
248           pyerr = PyList_SetItem(buffer, i, PyInt_FromLong( intBuf[i] ));
249         delete [] intBuf;
250       }
251     if ( pyerr )
252       {
253         PyErr_SetString(PyExc_RuntimeError, "Error of PyList_SetItem()");
254         return NULL;
255       }
256     return PyInt_FromLong( err );
257
258   }
259   %}
260
261 %pythoncode %{
262 def ParaMEDMEMDataArrayDoublenew(cls,*args):
263     import _ParaMEDMEM
264     return _ParaMEDMEM.DataArrayDouble____new___(cls,args)
265 def ParaMEDMEMDataArrayDoubleIadd(self,*args):
266     import _ParaMEDMEM
267     return _ParaMEDMEM.DataArrayDouble____iadd___(self, self, *args)
268 def ParaMEDMEMDataArrayDoubleIsub(self,*args):
269     import _ParaMEDMEM
270     return _ParaMEDMEM.DataArrayDouble____isub___(self, self, *args)
271 def ParaMEDMEMDataArrayDoubleImul(self,*args):
272     import _ParaMEDMEM
273     return _ParaMEDMEM.DataArrayDouble____imul___(self, self, *args)
274 def ParaMEDMEMDataArrayDoubleIdiv(self,*args):
275     import _ParaMEDMEM
276     return _ParaMEDMEM.DataArrayDouble____idiv___(self, self, *args)
277 def ParaMEDMEMDataArrayDoubleIpow(self,*args):
278     import _ParaMEDMEM
279     return _ParaMEDMEM.DataArrayDouble____ipow___(self, self, *args)
280 def ParaMEDMEMDataArrayDoubleTupleIadd(self,*args):
281     import _ParaMEDMEM
282     return _ParaMEDMEM.DataArrayDoubleTuple____iadd___(self, self, *args)
283 def ParaMEDMEMDataArrayDoubleTupleIsub(self,*args):
284     import _ParaMEDMEM
285     return _ParaMEDMEM.DataArrayDoubleTuple____isub___(self, self, *args)
286 def ParaMEDMEMDataArrayDoubleTupleImul(self,*args):
287     import _ParaMEDMEM
288     return _ParaMEDMEM.DataArrayDoubleTuple____imul___(self, self, *args)
289 def ParaMEDMEMDataArrayDoubleTupleIdiv(self,*args):
290     import _ParaMEDMEM
291     return _ParaMEDMEM.DataArrayDoubleTuple____idiv___(self, self, *args)
292 def ParaMEDMEMMEDCouplingFieldDoublenew(cls,*args):
293     import _ParaMEDMEM
294     return _ParaMEDMEM.MEDCouplingFieldDouble____new___(cls,args)
295 def ParaMEDMEMMEDCouplingFieldDoubleIadd(self,*args):
296     import _ParaMEDMEM
297     return _ParaMEDMEM.MEDCouplingFieldDouble____iadd___(self, self, *args)
298 def ParaMEDMEMMEDCouplingFieldDoubleIsub(self,*args):
299     import _ParaMEDMEM
300     return _ParaMEDMEM.MEDCouplingFieldDouble____isub___(self, self, *args)
301 def ParaMEDMEMMEDCouplingFieldDoubleImul(self,*args):
302     import _ParaMEDMEM
303     return _ParaMEDMEM.MEDCouplingFieldDouble____imul___(self, self, *args)
304 def ParaMEDMEMMEDCouplingFieldDoubleIdiv(self,*args):
305     import _ParaMEDMEM
306     return _ParaMEDMEM.MEDCouplingFieldDouble____idiv___(self, self, *args)
307 def ParaMEDMEMMEDCouplingFieldDoubleIpow(self,*args):
308     import _ParaMEDMEM
309     return _ParaMEDMEM.MEDCouplingFieldDouble____ipow___(self, self, *args)
310 def ParaMEDMEMDataArrayIntnew(cls,*args):
311     import _ParaMEDMEM
312     return _ParaMEDMEM.DataArrayInt____new___(cls,args)
313 def ParaMEDMEMDataArrayIntIadd(self,*args):
314     import _ParaMEDMEM
315     return _ParaMEDMEM.DataArrayInt____iadd___(self, self, *args)
316 def ParaMEDMEMDataArrayIntIsub(self,*args):
317     import _ParaMEDMEM
318     return _ParaMEDMEM.DataArrayInt____isub___(self, self, *args)
319 def ParaMEDMEMDataArrayIntImul(self,*args):
320     import _ParaMEDMEM
321     return _ParaMEDMEM.DataArrayInt____imul___(self, self, *args)
322 def ParaMEDMEMDataArrayIntIdiv(self,*args):
323     import _ParaMEDMEM
324     return _ParaMEDMEM.DataArrayInt____idiv___(self, self, *args)
325 def ParaMEDMEMDataArrayIntImod(self,*args):
326     import _ParaMEDMEM
327     return _ParaMEDMEM.DataArrayInt____imod___(self, self, *args)
328 def ParaMEDMEMDataArrayIntIpow(self,*args):
329     import _ParaMEDMEM
330     return _ParaMEDMEM.DataArrayInt____ipow___(self, self, *args)
331 def ParaMEDMEMDataArrayIntTupleIadd(self,*args):
332     import _ParaMEDMEM
333     return _ParaMEDMEM.DataArrayIntTuple____iadd___(self, self, *args)
334 def ParaMEDMEMDataArrayIntTupleIsub(self,*args):
335     import _ParaMEDMEM
336     return _ParaMEDMEM.DataArrayIntTuple____isub___(self, self, *args)
337 def ParaMEDMEMDataArrayIntTupleImul(self,*args):
338     import _ParaMEDMEM
339     return _ParaMEDMEM.DataArrayIntTuple____imul___(self, self, *args)
340 def ParaMEDMEMDataArrayIntTupleIdiv(self,*args):
341     import _ParaMEDMEM
342     return _ParaMEDMEM.DataArrayIntTuple____idiv___(self, self, *args)
343 def ParaMEDMEMDataArrayIntTupleImod(self,*args):
344     import _ParaMEDMEM
345     return _ParaMEDMEM.DataArrayIntTuple____imod___(self, self, *args)
346 %}
347
348 %include "MEDCouplingFinalize.i"