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)
29 raise RuntimeError, "Expect MPI_COMM_WORLD size == 5"
32 procs_source = range( nproc_source )
33 procs_target = range( size - nproc_source + 1, size)
35 interface = CommInterface()
37 target_group = MPIProcessorGroup(interface, procs_target)
38 source_group = MPIProcessorGroup(interface, procs_source)
50 dec = NonCoincidentDEC(source_group, target_group)
52 data_dir = os.environ['MED_ROOT_DIR']
53 tmp_dir = os.environ['TMP']
58 filename_xml1 = data_dir + "/share/salome/resources/med/square1_split"
59 filename_xml2 = data_dir + "/share/salome/resources/med/square2_split"
61 MPI_Barrier(MPI_COMM_WORLD)
63 if source_group.containsMyRank():
65 filename = filename_xml1 + str(rank+1) + ".med"
66 meshname = "Mesh_2_" + str(rank+1)
68 mesh = MESH(MED_DRIVER, filename, meshname)
69 support = SUPPORT(mesh, "all elements", MED_CELL)
70 paramesh = ParaMESH(mesh, source_group, "source mesh")
72 parasupport = UnstructuredParaSUPPORT( support, source_group)
73 comptopo = ComponentTopology()
75 parafield = ParaFIELD(parasupport, comptopo)
77 nb_local = support.getNumberOfElements(MED_ALL_ELEMENTS);
79 value = [1.0]*nb_local
81 parafield.getField().setValue(value)
82 icocofield = ICoCo_MEDField(paramesh,parafield)
83 dec.attachLocalField(icocofield,'P0')
86 if target_group.containsMyRank():
88 filename = filename_xml2 + str(rank - nproc_source + 1) + ".med"
89 meshname = "Mesh_3_" + str(rank - nproc_source + 1)
91 mesh = MESH(MED_DRIVER, filename, meshname)
92 support = SUPPORT(mesh, "all elements", MED_CELL)
93 paramesh = ParaMESH(mesh, target_group, "target mesh")
95 parasupport = UnstructuredParaSUPPORT( support, target_group)
96 comptopo = ComponentTopology()
97 parafield = ParaFIELD(parasupport, comptopo)
99 nb_local = support.getNumberOfElements(MED_ALL_ELEMENTS)
100 value = [0.0]*nb_local
102 parafield.getField().setValue(value)
103 icocofield = ICoCo_MEDField(paramesh,parafield)
105 dec.attachLocalField(icocofield, 'P0')
108 field_before_int = [0.0]
109 field_after_int = [0.0]
111 if source_group.containsMyRank():
113 field_before_int = [parafield.getVolumeIntegral(1)]
114 MPI_Bcast(field_before_int, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD);
117 dec.setForcedRenormalization(False)
122 if target_group.containsMyRank():
124 MPI_Bcast(field_before_int, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD)
126 dec.setForcedRenormalization(False)
128 field_after_int = [parafield.getVolumeIntegral(1)]
131 MPI_Bcast(field_before_int, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD)
132 MPI_Bcast(field_after_int , 1, MPI_DOUBLE, size-1, MPI_COMM_WORLD)
135 if abs(field_before_int[0] - field_after_int[0]) > epsilon:
136 print "Field before is not equal field after: %s != %s"%\
137 (field_before_int[0],field_after_int[0])
141 MPI_Barrier(MPI_COMM_WORLD)
143 print "# End of testNonCoincidentDEC"