]> SALOME platform Git repositories - tools/medcoupling.git/blob - src/ParaMEDMEM_Swig/libParaMEDMEM_Swig.i
Salome HOME
Merge from BR_V5_DEV 16Feb09
[tools/medcoupling.git] / src / ParaMEDMEM_Swig / libParaMEDMEM_Swig.i
1 // Copyright (C) 2005  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
2 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
3 // 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either 
7 // version 2.1 of the License.
8 // 
9 // This library is distributed in the hope that it will be useful 
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of 
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU 
12 // Lesser General Public License for more details.
13 //
14 // You should have received a copy of the GNU Lesser General Public  
15 // License along with this library; if not, write to the Free Software 
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
17 //
18 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
19 //
20 %module libParaMEDMEM_Swig
21
22 %include "libParaMEDMEM_Swig.typemap"
23 %include "../MEDMEM_SWIG/libMEDMEM_Swig.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 <IntersectionDEC.hxx>
32 #include <NonCoincidentDEC.hxx>
33 #include <StructuredCoincidentDEC.hxx>
34 #include <ParaMESH.hxx>
35 #include <UnstructuredParaSUPPORT.hxx>
36 #include <StructuredParaSUPPORT.hxx>
37 #include <ParaFIELD.hxx>
38 #include <ICoCoMEDField.hxx>
39
40 #include <mpi.h>
41
42 using namespace ParaMEDMEM;
43 using namespace ICoCo;
44
45  enum mpi_constants { mpi_comm_world, mpi_comm_self, mpi_double, mpi_int };
46
47 %}
48
49
50 class CommInterface {
51 public:
52   CommInterface();
53 };
54
55 class ProcessorGroup {
56 public:
57   %extend {
58   // avoid error when SWIG creates a default constructor of an abstract class
59   ProcessorGroup() { return NULL; }
60   }
61 };
62
63 class MPIProcessorGroup: public ProcessorGroup {
64 public:
65   MPIProcessorGroup(const CommInterface& interface);
66   MPIProcessorGroup(const CommInterface& interface, set<int> proc_ids);
67   MPIProcessorGroup(const CommInterface& interface,int pstart, int pend);
68
69   int translateRank(const ProcessorGroup* group, int rank) const;
70   ProcessorGroup* createComplementProcGroup() const;
71   ProcessorGroup* fuse (const ProcessorGroup&) const;
72
73   bool containsMyRank() const;
74   int myRank() const;
75 };
76
77 %rename(ICoCo_MEDField) ICoCo::MEDField;
78 namespace ICoCo {
79 class Field {
80 }; 
81
82 class MEDField: public Field {
83 public:
84  MEDField(ParaMESH* mesh, ParaFIELD* field);
85 };
86 }
87
88 class DEC {
89 public:
90   %extend {
91   // avoid error when SWIG creates a default constructor of an abstract class
92   DEC() { return NULL; };
93   }
94   void recvData();
95   void sendData();
96   void synchronize();
97
98   void attachLocalField(FIELD<double>*      field, const char* method="P0");
99   void attachLocalField(const ParaFIELD*    field, const char* method="P0");
100   void attachLocalField(const ICoCo::Field* field, const char* method="P0");
101
102   void renormalizeTargetField();
103 };
104
105
106 typedef enum{WithoutTimeInterp,LinearTimeInterp} TimeInterpolationMethod;
107
108 typedef enum{Native,PointToPoint} AllToAllMethod;
109
110 class IntersectionDEC: public DEC { 
111 public:
112
113   IntersectionDEC(ProcessorGroup& local_group, ProcessorGroup&  distant_group);
114   
115   void recvData();
116   void recvData( double time );
117
118   void sendData();
119   void sendData( double time, double deltatime );
120
121   void setTimeInterpolationMethod(TimeInterpolationMethod it);
122   TimeInterpolationMethod getTimeInterpolationMethod();
123
124   void setAsynchronous( bool dr);
125   bool getAsynchronous();
126
127   AllToAllMethod getAllToAllMethod();
128   void setAllToAllMethod(AllToAllMethod sp);
129
130   bool getForcedRenormalization();
131   void setForcedRenormalization( bool dr);
132 };
133
134
135
136 /* This object can be used only if MED_ENABLE_FVM is defined*/
137 #ifdef MED_ENABLE_FVM
138 class NonCoincidentDEC: public DEC {
139 public:
140   NonCoincidentDEC(ProcessorGroup& source, ProcessorGroup& target);
141 };
142 #endif
143
144
145 class StructuredCoincidentDEC: public DEC {
146 public:
147   StructuredCoincidentDEC(ProcessorGroup& source, ProcessorGroup& target);
148   void prepareSourceDE();
149   void prepareTargetDE();
150 };
151
152
153 class ParaMESH {
154 public:
155   ParaMESH(driverTypes driver_type, const char* file_name, const ProcessorGroup& group);
156   ParaMESH(MESH& subdomain_mesh, const ProcessorGroup& proc_group, const char* name);
157
158   void write(driverTypes driverType, const char* fileName="");
159   MESH* getMesh() const;
160 };
161
162
163
164 class ParaSUPPORT {
165 public:
166   %extend {
167   // avoid error when SWIG creates a default constructor of an abstract class
168   ParaSUPPORT() { return NULL; }
169
170   PyObject* getGlobalNumbering() {
171     const int * numIndex = self->getGlobalNumbering();
172     int aSize = self->getSupport()->getNumberOfElements(999);
173     TYPEMAP_OUTPUT_ARRAY(numIndex, aSize, PyInt_FromLong, ParaSUPPORT::getGlobalNumbering);
174   }
175   }
176 };
177
178
179 class UnstructuredParaSUPPORT: public ParaSUPPORT {
180 public:
181   UnstructuredParaSUPPORT(const ParaMESH* const mesh, const SUPPORT* support );
182   UnstructuredParaSUPPORT(const ParaMESH* const mesh, const medEntityMesh entity);
183   UnstructuredParaSUPPORT(const SUPPORT* support, const ProcessorGroup& group);
184 };
185
186 class StructuredParaSUPPORT: public ParaSUPPORT {
187 public:
188   StructuredParaSUPPORT(const ParaGRID* const grid, const medEntityMesh entity);
189   StructuredParaSUPPORT(const ParaMESH* const mesh, const medEntityMesh entity);
190 };
191
192 class ComponentTopology {
193 public:
194   ComponentTopology();
195   ComponentTopology(int nb_comp);
196   ComponentTopology(int nb_comp, int nb_blocks);
197   ComponentTopology(int nb_comp, ProcessorGroup* group);
198
199   int nbComponents();
200   int nbLocalComponents();
201 };
202
203
204 class ParaFIELD {
205 public:
206   ParaFIELD(const ParaSUPPORT* support,
207             const ComponentTopology& component_topology);
208
209   ParaFIELD(driverTypes driver_type,
210             const char* file_name, 
211             const char* driver_name,
212             const ComponentTopology& component_topology);
213
214   ParaFIELD(FIELDDOUBLE* field, const ProcessorGroup& group);
215
216   void write(driverTypes driverType, 
217              const char* fileName="", 
218              const char* meshName="");
219
220   FIELDDOUBLE* getField() const;
221
222   double getVolumeIntegral(int icomp) const;
223 };
224
225 //=============================================================================================
226 // Interface for MPI-realization-specific constants like MPI_COMM_WORLD.
227 //
228 // Type and values of constants like MPI_COMM_WORLD depends on MPI realization
229 // and usually such constants actually are macros. To have such symbols in python
230 // and translate them into correct values we use the following technique.
231 // We define some constants (enum mpi_constants) and map them into real MPI values
232 // using typemaps, and we create needed python symbols equal to 'mpi_constants'
233 // via %pythoncode directive.
234
235 // Constants corresponding to similar MPI definitions
236 enum mpi_constants { mpi_comm_world, mpi_comm_self, mpi_double, mpi_int };
237
238 // Map mpi_comm_world and mpi_comm_self -> MPI_COMM_WORLD and MPI_COMM_SELF
239 %typemap(in) MPI_Comm
240
241   switch (PyInt_AsLong($input)) {
242   case mpi_comm_world: $1 = MPI_COMM_WORLD; break;
243   case mpi_comm_self:  $1 = MPI_COMM_SELF;  break;
244   default:
245     PyErr_SetString(PyExc_TypeError,"unexpected value of MPI_Comm");
246     return NULL;
247   }
248 }
249 // Map mpi_double and mpi_int -> MPI_DOUBLE and MPI_INT
250 %typemap(in) MPI_Datatype
251
252   switch (PyInt_AsLong($input)) {
253   case mpi_double:     $1 = MPI_DOUBLE;     break;
254   case mpi_int:        $1 = MPI_INT;        break;
255   default:
256     PyErr_SetString(PyExc_TypeError,"unexpected value of MPI_Datatype");
257     return NULL;
258   }
259 }
260 // The following code gets inserted into the result python file:
261 // create needed python symbols
262 %pythoncode %{
263   MPI_COMM_WORLD = mpi_comm_world
264   MPI_COMM_SELF  = mpi_comm_self
265   MPI_DOUBLE     = mpi_double
266   MPI_INT        = mpi_int
267 %}
268 //=============================================================================================
269
270 // ==============
271 // MPI_Comm_size
272 // ==============
273 %inline %{ PyObject* MPI_Comm_size(MPI_Comm comm) {
274   int res = 0;
275   int err = MPI_Comm_size(comm, &res);
276   if ( err != MPI_SUCCESS ) {
277     PyErr_SetString(PyExc_RuntimeError,"Erorr in MPI_Comm_size()");
278     return NULL;
279   }
280   return PyInt_FromLong( res );
281 } %}
282
283 // ==============
284 // MPI_Comm_rank
285 // ==============
286 %inline %{ PyObject* MPI_Comm_rank(MPI_Comm comm) {
287   int res = 0;
288   int err = MPI_Comm_rank(comm, &res);
289   if ( err != MPI_SUCCESS ) {
290     PyErr_SetString(PyExc_RuntimeError,"Erorr in MPI_Comm_rank()");
291     return NULL;
292   }
293   return PyInt_FromLong( res );
294 } %}
295
296 int MPI_Init(int *argc, char ***argv );
297 int MPI_Barrier(MPI_Comm comm);
298 int MPI_Finalize();
299
300 // ==========
301 // MPI_Bcast
302 // ==========
303
304 %inline %{ PyObject* MPI_Bcast(PyObject* buffer, int nb, MPI_Datatype type, int root, MPI_Comm c) {
305   // buffer must be a list
306   if (!PyList_Check(buffer)) {
307     PyErr_SetString(PyExc_TypeError, "buffer is expected to be a list");
308     return NULL;
309   }
310   // check list size
311   int aSize = PyList_Size(buffer);
312   if ( aSize != nb ) {
313     PyErr_SetString(PyExc_ValueError, MEDMEM::STRING("buffer is expected to be of size ")<<nb);
314     return NULL;
315   }
316   // allocate and fill a buffer
317   void* aBuf = 0;
318   int* intBuf = 0;
319   double* dblBuf = 0;
320   if ( type == MPI_DOUBLE ) {
321     aBuf = (void*) ( dblBuf = new double[ nb ] );
322     for ( int i = 0; i < aSize; ++i )
323       dblBuf[i] = PyFloat_AS_DOUBLE( PyList_GetItem( buffer, i ));
324   }
325   else if ( type == MPI_INT ) {
326     aBuf = (void*) ( intBuf = new int[ nb ] );
327     for ( int i = 0; i < aSize; ++i )
328       intBuf[i] = int( PyInt_AS_LONG( PyList_GetItem( buffer, i )));
329   }
330   else {
331     PyErr_SetString(PyExc_TypeError, "Only MPI_DOUBLE and MPI_INT supported");
332     return NULL;
333   }
334   // call MPI_Bcast
335   int err = MPI_Bcast(aBuf, nb, type, root, c);
336   // treat error
337   if ( err != MPI_SUCCESS ) {
338     PyErr_SetString(PyExc_RuntimeError,"Erorr in MPI_Bcast()");
339     delete [] intBuf; delete [] dblBuf;
340     return NULL;
341   }
342   // put recieved data into the list
343   int pyerr = 0;
344   if ( type == MPI_DOUBLE ) {
345     for ( int i = 0; i < aSize && !pyerr; ++i )
346       pyerr = PyList_SetItem(buffer, i, PyFloat_FromDouble( dblBuf[i] ));
347     delete [] dblBuf;
348   }
349   else {
350     for ( int i = 0; i < aSize && !pyerr; ++i )
351       pyerr = PyList_SetItem(buffer, i, PyInt_FromLong( intBuf[i] ));
352     delete [] intBuf;
353   }
354   if ( pyerr ) {
355     PyErr_SetString(PyExc_RuntimeError, "Error of PyList_SetItem()");
356     return NULL;
357   }
358   return PyInt_FromLong( err );
359
360 } %}
361