2 # -*- coding: iso-8859-1 -*-
3 # Copyright (C) 2007-2016 CEA/DEN, EDF R&D
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, or (at your option) any later version.
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
22 from ParaMEDMEM import *
27 size = MPI_Comm_size(MPI_COMM_WORLD)
28 rank = MPI_Comm_rank(MPI_COMM_WORLD)
30 raise RuntimeError("Expect MPI_COMM_WORLD size == 5")
33 procs_source = list(range(nproc_source))
34 procs_target = list(range(size - nproc_source + 1, size))
36 interface = CommInterface()
38 target_group = MPIProcessorGroup(interface, procs_target)
39 source_group = MPIProcessorGroup(interface, procs_source)
51 dec = NonCoincidentDEC(source_group, target_group)
53 data_dir = os.environ['MEDCOUPLING_ROOT_DIR']
54 tmp_dir = os.environ['TMP']
59 filename_xml1 = data_dir + "/share/resources/med/square1_split"
60 filename_xml2 = data_dir + "/share/resources/med/square2_split"
62 MPI_Barrier(MPI_COMM_WORLD)
64 if source_group.containsMyRank():
66 filename = filename_xml1 + str(rank+1) + ".med"
67 meshname = "Mesh_2_" + str(rank+1)
69 mesh = MESH(MED_DRIVER, filename, meshname)
70 support = SUPPORT(mesh, "all elements", MED_CELL)
71 paramesh = ParaMESH(mesh, source_group, "source mesh")
73 parasupport = UnstructuredParaSUPPORT( support, source_group)
74 comptopo = ComponentTopology()
76 parafield = ParaFIELD(parasupport, comptopo)
78 nb_local = support.getNumberOfElements(MED_ALL_ELEMENTS);
80 value = [1.0]*nb_local
82 parafield.getField().setValue(value)
83 icocofield = ICoCo_MEDField(paramesh,parafield)
84 dec.attachLocalField(icocofield,'P0')
87 if target_group.containsMyRank():
89 filename = filename_xml2 + str(rank - nproc_source + 1) + ".med"
90 meshname = "Mesh_3_" + str(rank - nproc_source + 1)
92 mesh = MESH(MED_DRIVER, filename, meshname)
93 support = SUPPORT(mesh, "all elements", MED_CELL)
94 paramesh = ParaMESH(mesh, target_group, "target mesh")
96 parasupport = UnstructuredParaSUPPORT( support, target_group)
97 comptopo = ComponentTopology()
98 parafield = ParaFIELD(parasupport, comptopo)
100 nb_local = support.getNumberOfElements(MED_ALL_ELEMENTS)
101 value = [0.0]*nb_local
103 parafield.getField().setValue(value)
104 icocofield = ICoCo_MEDField(paramesh,parafield)
106 dec.attachLocalField(icocofield, 'P0')
109 field_before_int = [0.0]
110 field_after_int = [0.0]
112 if source_group.containsMyRank():
114 field_before_int = [parafield.getVolumeIntegral(1)]
115 MPI_Bcast(field_before_int, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD);
118 dec.setForcedRenormalization(False)
123 if target_group.containsMyRank():
125 MPI_Bcast(field_before_int, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD)
127 dec.setForcedRenormalization(False)
129 field_after_int = [parafield.getVolumeIntegral(1)]
132 MPI_Bcast(field_before_int, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD)
133 MPI_Bcast(field_after_int , 1, MPI_DOUBLE, size-1, MPI_COMM_WORLD)
136 if abs(field_before_int[0] - field_after_int[0]) > epsilon:
137 print("Field before is not equal field after: %s != %s"%\
138 (field_before_int[0],field_after_int[0]))
142 MPI_Barrier(MPI_COMM_WORLD)
144 print("# End of testNonCoincidentDEC")