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.
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.
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
19 # See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 from ParaMEDMEM import *
26 size = MPI_Comm_size(MPI_COMM_WORLD)
27 rank = MPI_Comm_rank(MPI_COMM_WORLD)
30 raise RuntimeError, "Expect MPI_COMM_WORLD size >= 4"
32 interface = CommInterface()
34 self_group = MPIProcessorGroup(interface, rank, rank)
35 target_group = MPIProcessorGroup(interface, 3, size-1)
36 source_group = MPIProcessorGroup(interface, 0, 2)
43 data_dir = os.environ['MED_ROOT_DIR']
44 tmp_dir = os.environ['TMP']
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"
54 dec = StructuredCoincidentDEC(source_group, target_group)
55 MPI_Barrier(MPI_COMM_WORLD)
57 if source_group.containsMyRank():
58 filename = filename_xml1 + str(rank+1) + ".med"
59 meshname = "Mesh_2_" + str(rank+1)
61 mesh = MESH(MED_DRIVER, filename, meshname)
62 support = SUPPORT(mesh, "all elements", MED_CELL)
63 paramesh = ParaMESH(mesh, source_group, "source mesh")
65 parasupport = UnstructuredParaSUPPORT( support, source_group);
66 comptopo = ComponentTopology(6)
67 parafield = ParaFIELD(parasupport, comptopo)
69 nb_local = support.getNumberOfElements(MED_ALL_ELEMENTS);
70 global_numbering = parasupport.getGlobalNumbering()
73 for ielem in range(nb_local):
74 for icomp in range(6):
75 value.append(global_numbering[ielem]*6.0+icomp);
79 parafield.getField().setValue(value)
80 icocofield = ICoCo_MEDField(paramesh,parafield)
81 dec.attachLocalField(icocofield)
86 if target_group.containsMyRank():
89 mesh = MESH(MED_DRIVER, filename_2, meshname2)
90 support = SUPPORT(mesh, "all elements", MED_CELL)
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)
97 nb_local = support.getNumberOfElements(MED_ALL_ELEMENTS)
98 value = [0.0]*(nb_local*comptopo.nbLocalComponents())
100 parafield.getField().setValue(value)
101 icocofield = ICoCo_MEDField(paramesh,parafield)
103 dec.attachLocalField(icocofield)
107 recv_value = parafield.getField().getValue()
112 print "End of test StructuredCoincidentDEC"