Salome HOME
539d452c7972044994d963326b50d8eb1df0f810
[tools/medcoupling.git] / src / ParaMEDMEM_Swig / test_InterpKernelDEC.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 from MEDLoader import ReadUMeshFromFile
24 import sys, os
25 import unittest
26 import math
27 from mpi4py import MPI
28
29
30 class ParaMEDMEMBasicsTest(unittest.TestCase):
31     def testInterpKernelDEC_2D(self):
32         size = MPI.COMM_WORLD.size
33         rank = MPI.COMM_WORLD.rank
34         if size != 5:
35             raise RuntimeError("Expect MPI_COMM_WORLD size == 5")
36         print(rank)
37         nproc_source = 3
38         procs_source = list(range(nproc_source))
39         procs_target = list(range(size - nproc_source + 1, size))
40
41         interface = CommInterface()
42         target_group = MPIProcessorGroup(interface, procs_target)
43         source_group = MPIProcessorGroup(interface, procs_source)
44         dec = InterpKernelDEC(source_group, target_group)
45
46         mesh       =0
47         support    =0
48         paramesh   =0
49         parafield  =0
50         icocofield =0
51         data_dir = os.environ['MEDCOUPLING_ROOT_DIR']
52         tmp_dir  = os.environ['TMP']
53
54         if not tmp_dir or len(tmp_dir)==0:
55             tmp_dir = "/tmp"
56             pass
57
58         filename_xml1 = os.path.join(data_dir, "share/resources/med/square1_split")
59         filename_xml2 = os.path.join(data_dir, "share/resources/med/square2_split")
60
61         MPI.COMM_WORLD.Barrier()
62         if source_group.containsMyRank():
63             filename = filename_xml1 + str(rank+1) + ".med"
64             meshname = "Mesh_2_" + str(rank+1)
65             mesh=ReadUMeshFromFile(filename,meshname,0)
66             paramesh=ParaMESH(mesh,source_group,"source mesh")
67             comptopo = ComponentTopology()
68             parafield = ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo)
69             parafield.getField().setNature(IntensiveMaximum)
70             nb_local=mesh.getNumberOfCells()
71             value = [1.0]*nb_local
72             parafield.getField().setValues(value)
73             icocofield = ICoCoMEDField(parafield.getField())
74             dec.attachLocalField(icocofield)
75             pass
76         else:
77             filename = filename_xml2 + str(rank - nproc_source + 1) + ".med"
78             meshname = "Mesh_3_" + str(rank - nproc_source + 1)
79             mesh=ReadUMeshFromFile(filename,meshname,0)
80             paramesh=ParaMESH(mesh,target_group,"target mesh")
81             comptopo = ComponentTopology()
82             parafield = ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo)
83             parafield.getField().setNature(IntensiveMaximum)
84             nb_local=mesh.getNumberOfCells()
85             value = [0.0]*nb_local
86             parafield.getField().setValues(value)
87             icocofield = ICoCoMEDField(parafield.getField())
88             dec.attachLocalField(icocofield)
89             pass
90
91         if source_group.containsMyRank():
92             field_before_int = parafield.getVolumeIntegral(0,True)
93             dec.synchronize()
94             dec.setForcedRenormalization(False)
95             dec.sendData()
96             dec.recvData()
97             field_after_int=parafield.getVolumeIntegral(0,True);
98             self.assertTrue(math.fabs(field_after_int-field_before_int)<1e-8)
99             pass
100         else:
101             dec.synchronize()
102             dec.setForcedRenormalization(False)
103             dec.recvData()
104             dec.sendData()
105             pass
106         ## end
107         interface = 0
108         target_group = 0
109         source_group = 0
110         dec = 0
111         mesh       =0
112         support    =0
113         paramesh   =0
114         parafield  =0
115         icocofield =0
116         MPI.COMM_WORLD.Barrier()
117         MPI.Finalize()
118         pass
119     pass
120
121 unittest.main()