Salome HOME
48cfb7c2116339896184eba93b7690f743eb9a64
[tools/medcoupling.git] / src / ParaMEDMEM_Swig / test_StructuredCoincidentDEC.py
1 #!/usr/bin/env python
2 #  -*- coding: iso-8859-1 -*-
3 # Copyright (C) 2007-2016  CEA/DEN, EDF R&D
4 #
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.
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
22 from ParaMEDMEM import *
23 import sys, os
24 import unittest
25 import math
26
27 class ParaMEDMEMBasicsTest2(unittest.TestCase):
28     def testStructuredCoincidentDEC(self):
29         MPI_Init(sys.argv)
30         #
31         size = MPI_Comm_size(MPI_COMM_WORLD)
32         rank = MPI_Comm_rank(MPI_COMM_WORLD)
33         #
34         if size < 4:
35             raise RuntimeError("Expect MPI_COMM_WORLD size >= 4")
36         #
37         interface = CommInterface()
38         #
39         self_group   = MPIProcessorGroup(interface, rank, rank)
40         target_group = MPIProcessorGroup(interface, 3, size-1)
41         source_group = MPIProcessorGroup(interface, 0, 2)
42         #
43         mesh      = 0
44         support   = 0
45         paramesh  = 0
46         parafield = 0
47         comptopo  = 0
48         icocofield= 0
49         #
50         data_dir = os.environ['MEDCOUPLING_ROOT_DIR']
51         tmp_dir  = os.environ['TMP']
52         if tmp_dir == '':
53             tmp_dir = "/tmp"
54             pass
55
56         filename_xml1    = data_dir + "/share/resources/med/square1_split"
57         filename_2       = data_dir + "/share/resources/med/square1.med"
58         filename_seq_wr  = tmp_dir + "/"
59         filename_seq_med = tmp_dir + "/myWrField_seq_pointe221.med"
60
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=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(IntensiveMaximum)
71             nb_local=mesh.getNumberOfCells()
72             global_numbering=paramesh.getGlobalNumberingCell2()
73             value = []
74             for ielem in range(nb_local):
75                 for icomp in range(6):
76                     value.append(global_numbering[ielem]*6.0+icomp);
77                     pass
78                 pass
79             parafield.getField().setValues(value)
80             icocofield = ICoCoMEDField(mesh,parafield.getField())
81             dec.setMethod("P0")
82             dec.attachLocalField(parafield)
83             dec.synchronize()
84             dec.sendData()
85             pass
86
87         if target_group.containsMyRank():
88             meshname2 = "Mesh_2"
89             mesh=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(IntensiveMaximum)
94             nb_local=mesh.getNumberOfCells()
95             value = [0.0]*(nb_local*comptopo.nbLocalComponents())
96             parafield.getField().setValues(value)
97             icocofield = ICoCoMEDField(mesh,parafield.getField())
98             dec.setMethod("P0")
99             dec.attachLocalField(parafield)
100             dec.synchronize()
101             dec.recvData()
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.assertTrue(math.fabs(recv_value[i*comptopo.nbLocalComponents()+icomp]-
107                                               (float)(i*6+icomp+first))<1e-12)
108                     pass
109                 pass
110             pass
111         comptopo=0
112         interface = 0
113         mesh       =0
114         support    =0
115         paramesh   =0
116         parafield  =0
117         icocofield =0
118         dec=0
119         self_group =0
120         target_group = 0
121         source_group = 0
122         MPI_Barrier(MPI_COMM_WORLD)
123         MPI_Finalize()
124         print("End of test StructuredCoincidentDEC")
125         pass
126
127
128 unittest.main()