Salome HOME
5df6ca3637c1f32d548b20e19f7c939a509868b0
[tools/medcoupling.git] / src / ParaMEDMEM_Swig / test_NonCoincidentDEC.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
25 MPI_Init(sys.argv)
26
27 size = MPI_Comm_size(MPI_COMM_WORLD)
28 rank = MPI_Comm_rank(MPI_COMM_WORLD)
29 if size != 5:
30     raise RuntimeError("Expect MPI_COMM_WORLD size == 5")
31
32 nproc_source = 3
33 procs_source = list(range(nproc_source))
34 procs_target = list(range(size - nproc_source + 1, size))
35
36 interface = CommInterface()
37
38 target_group = MPIProcessorGroup(interface, procs_target)
39 source_group = MPIProcessorGroup(interface, procs_source)
40
41 source_mesh= 0
42 target_mesh= 0
43 parasupport= 0
44 mesh       = 0
45 support    = 0
46 field      = 0
47 paramesh   = 0
48 parafield  = 0
49 icocofield = 0
50
51 dec = NonCoincidentDEC(source_group, target_group)
52
53 data_dir = os.environ['MEDCOUPLING_ROOT_DIR']
54 tmp_dir  = os.environ['TMP']
55 if tmp_dir == '':
56     tmp_dir = "/tmp"
57     pass
58
59 filename_xml1 = data_dir + "/share/resources/med/square1_split"
60 filename_xml2 = data_dir + "/share/resources/med/square2_split"
61
62 MPI_Barrier(MPI_COMM_WORLD)
63
64 if source_group.containsMyRank():
65
66     filename = filename_xml1 + str(rank+1) + ".med"
67     meshname = "Mesh_2_" + str(rank+1)
68
69     mesh = MESH(MED_DRIVER, filename, meshname)
70     support = SUPPORT(mesh, "all elements", MED_CELL)
71     paramesh = ParaMESH(mesh, source_group, "source mesh")
72
73     parasupport = UnstructuredParaSUPPORT( support, source_group)
74     comptopo = ComponentTopology()
75
76     parafield = ParaFIELD(parasupport, comptopo)
77
78     nb_local = support.getNumberOfElements(MED_ALL_ELEMENTS);
79
80     value = [1.0]*nb_local
81
82     parafield.getField().setValue(value)
83     icocofield = ICoCo_MEDField(paramesh,parafield)
84     dec.attachLocalField(icocofield,'P0')
85     pass
86
87 if target_group.containsMyRank():
88
89     filename = filename_xml2 + str(rank - nproc_source + 1) + ".med"
90     meshname = "Mesh_3_" + str(rank - nproc_source + 1)
91
92     mesh = MESH(MED_DRIVER, filename, meshname)
93     support = SUPPORT(mesh, "all elements", MED_CELL)
94     paramesh = ParaMESH(mesh, target_group, "target mesh")
95
96     parasupport = UnstructuredParaSUPPORT( support, target_group)
97     comptopo = ComponentTopology()
98     parafield = ParaFIELD(parasupport, comptopo)
99
100     nb_local = support.getNumberOfElements(MED_ALL_ELEMENTS)
101     value = [0.0]*nb_local
102
103     parafield.getField().setValue(value)
104     icocofield = ICoCo_MEDField(paramesh,parafield)
105
106     dec.attachLocalField(icocofield, 'P0')
107     pass
108
109 field_before_int = [0.0]
110 field_after_int = [0.0]
111
112 if source_group.containsMyRank():
113
114     field_before_int = [parafield.getVolumeIntegral(1)]
115     MPI_Bcast(field_before_int, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD);
116     dec.synchronize()
117     print("DEC usage")
118     dec.setForcedRenormalization(False)
119
120     dec.sendData()
121     pass
122
123 if target_group.containsMyRank():
124
125     MPI_Bcast(field_before_int, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD)
126     dec.synchronize()
127     dec.setForcedRenormalization(False)
128     dec.recvData()
129     field_after_int = [parafield.getVolumeIntegral(1)]
130     pass
131
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)
134
135 epsilon = 1e-6
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]))
139     pass
140
141
142 MPI_Barrier(MPI_COMM_WORLD)
143 MPI_Finalize()
144 print("# End of testNonCoincidentDEC")