Salome HOME
Merge from BR_V5_DEV 16Feb09
[modules/med.git] / src / ParaMEDMEM_Swig / test_StructuredCoincodentDEC.py
1 #!/usr/bin/env python
2
3 # Copyright (C) 2005  OPEN CASCADE, CEA, EDF R&D, LEG
4 #           PRINCIPIA R&D, EADS CCR, Lip6, BV, CEDRAT
5 # This library is free software; you can redistribute it and/or
6 # modify it under the terms of the GNU Lesser General Public
7 # License as published by the Free Software Foundation; either 
8 # version 2.1 of the License.
9
10 # This library is distributed in the hope that it will be useful 
11 # but WITHOUT ANY WARRANTY; without even the implied warranty of 
12 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU 
13 # Lesser General Public License for more details.
14
15 # You should have received a copy of the GNU Lesser General Public  
16 # License along with this library; if not, write to the Free Software 
17 # Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
18
19 # See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
20
21 from ParaMEDMEM import *
22 import sys, os
23
24 MPI_Init(sys.argv)
25
26 size = MPI_Comm_size(MPI_COMM_WORLD)
27 rank = MPI_Comm_rank(MPI_COMM_WORLD)
28     
29 if size < 4:
30     raise RuntimeError, "Expect MPI_COMM_WORLD size >= 4"
31     
32 interface = CommInterface()
33     
34 self_group   = MPIProcessorGroup(interface, rank, rank)
35 target_group = MPIProcessorGroup(interface, 3, size-1)
36 source_group = MPIProcessorGroup(interface, 0, 2)
37
38 mesh      = 0
39 support   = 0
40 paramesh  = 0
41 parafield = 0
42
43 data_dir = os.environ['MED_ROOT_DIR']
44 tmp_dir  = os.environ['TMP']
45 if tmp_dir == '':
46     tmp_dir = "/tmp"
47     pass
48
49 filename_xml1    = data_dir + "/share/salome/resources/med/square1_split"
50 filename_2       = data_dir + "/share/salome/resources/med/square1.med"
51 filename_seq_wr  = tmp_dir + "/"
52 filename_seq_med = tmp_dir + "/myWrField_seq_pointe221.med"
53
54 dec = StructuredCoincidentDEC(source_group, target_group)
55 MPI_Barrier(MPI_COMM_WORLD)
56     
57 if source_group.containsMyRank():
58     filename = filename_xml1 + str(rank+1) + ".med"
59     meshname = "Mesh_2_" + str(rank+1)
60
61     mesh     = MESH(MED_DRIVER, filename, meshname)
62     support  = SUPPORT(mesh, "all elements", MED_CELL)
63     paramesh = ParaMESH(mesh, source_group, "source mesh")
64
65     parasupport = UnstructuredParaSUPPORT( support, source_group);
66     comptopo    = ComponentTopology(6)
67     parafield   = ParaFIELD(parasupport, comptopo)
68
69     nb_local = support.getNumberOfElements(MED_ALL_ELEMENTS);
70     global_numbering = parasupport.getGlobalNumbering()
71
72     value = []
73     for ielem in range(nb_local):
74         for icomp in range(6):
75             value.append(global_numbering[ielem]*6.0+icomp);
76             pass
77         pass
78
79     parafield.getField().setValue(value)
80     icocofield = ICoCo_MEDField(paramesh,parafield)
81     dec.attachLocalField(icocofield)
82     dec.synchronize()
83     dec.sendData()
84     pass
85
86 if target_group.containsMyRank():
87
88     meshname2 = "Mesh_2"
89     mesh      = MESH(MED_DRIVER, filename_2, meshname2)
90     support   = SUPPORT(mesh, "all elements", MED_CELL)
91
92     paramesh    = ParaMESH(mesh, self_group, "target mesh")
93     parasupport = UnstructuredParaSUPPORT( support, self_group)
94     comptopo    = ComponentTopology(6, target_group)
95     parafield   = ParaFIELD(parasupport, comptopo)
96
97     nb_local  = support.getNumberOfElements(MED_ALL_ELEMENTS)
98     value = [0.0]*(nb_local*comptopo.nbLocalComponents())
99
100     parafield.getField().setValue(value)
101     icocofield = ICoCo_MEDField(paramesh,parafield)
102
103     dec.attachLocalField(icocofield)
104     dec.synchronize()
105     dec.recvData()
106
107     recv_value = parafield.getField().getValue()
108     pass
109
110 MPI_Finalize()
111
112 print "End of test StructuredCoincidentDEC"