2 # -*- coding: iso-8859-1 -*-
3 # Copyright (C) 2007-2013 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.
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 class ParaMEDMEMBasicsTest2(unittest.TestCase):
28 def testStructuredCoincidentDEC(self):
31 size = MPI_Comm_size(MPI_COMM_WORLD)
32 rank = MPI_Comm_rank(MPI_COMM_WORLD)
35 raise RuntimeError, "Expect MPI_COMM_WORLD size >= 4"
37 interface = CommInterface()
39 self_group = MPIProcessorGroup(interface, rank, rank)
40 target_group = MPIProcessorGroup(interface, 3, size-1)
41 source_group = MPIProcessorGroup(interface, 0, 2)
50 data_dir = os.environ['MED_ROOT_DIR']
51 tmp_dir = os.environ['TMP']
56 filename_xml1 = data_dir + "/share/salome/resources/med/square1_split"
57 filename_2 = data_dir + "/share/salome/resources/med/square1.med"
58 filename_seq_wr = tmp_dir + "/"
59 filename_seq_med = tmp_dir + "/myWrField_seq_pointe221.med"
61 dec = StructuredCoincidentDEC(source_group, target_group)
62 MPI_Barrier(MPI_COMM_WORLD)
63 if source_group.containsMyRank():
64 filename = filename_xml1 + str(rank+1) + ".med"
65 meshname = "Mesh_2_" + str(rank+1)
66 mesh=MEDLoader.ReadUMeshFromFile(filename,meshname,0)
67 paramesh=ParaMESH(mesh,source_group,"source mesh")
68 comptopo=ComponentTopology(6)
69 parafield=ParaFIELD(ON_CELLS,NO_TIME,paramesh,comptopo)
70 parafield.getField().setNature(ConservativeVolumic)
71 nb_local=mesh.getNumberOfCells()
72 global_numbering=paramesh.getGlobalNumberingCell2()
74 for ielem in range(nb_local):
75 for icomp in range(6):
76 value.append(global_numbering[ielem]*6.0+icomp);
79 parafield.getField().setValues(value)
80 icocofield = ICoCoMEDField(mesh,parafield.getField())
82 dec.attachLocalField(parafield)
87 if target_group.containsMyRank():
89 mesh=MEDLoader.ReadUMeshFromFile(filename_2, meshname2,0)
90 paramesh=ParaMESH(mesh, self_group, "target mesh")
91 comptopo=ComponentTopology(6,target_group)
92 parafield=ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo)
93 parafield.getField().setNature(ConservativeVolumic)
94 nb_local=mesh.getNumberOfCells()
95 value = [0.0]*(nb_local*comptopo.nbLocalComponents())
96 parafield.getField().setValues(value)
97 icocofield = ICoCoMEDField(mesh,parafield.getField())
99 dec.attachLocalField(parafield)
102 recv_value = parafield.getField().getArray().getValues()
103 for i in range(nb_local):
104 first=comptopo.firstLocalComponent()
105 for icomp in range(comptopo.nbLocalComponents()):
106 self.failUnless(math.fabs(recv_value[i*comptopo.nbLocalComponents()+icomp]-
107 (float)(i*6+icomp+first))<1e-12)
122 MPI_Barrier(MPI_COMM_WORLD)
124 print "End of test StructuredCoincidentDEC"