Salome HOME
Initiating medtool
[modules/med.git] / src / ParaMEDMEM_Swig / test_InterpKernelDEC.py
1 #!/usr/bin/env python
2 #  -*- coding: iso-8859-1 -*-
3 # Copyright (C) 2007-2015  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 ParaMEDMEMBasicsTest(unittest.TestCase):
28     def testInterpKernelDEC_2D(self):
29         MPI_Init(sys.argv)
30         size = MPI_Comm_size(MPI_COMM_WORLD)
31         rank = MPI_Comm_rank(MPI_COMM_WORLD)
32         if size != 5:
33             raise RuntimeError, "Expect MPI_COMM_WORLD size == 5"
34         print rank
35         nproc_source = 3
36         procs_source = range( nproc_source )
37         procs_target = range( size - nproc_source + 1, size)
38         
39         interface = CommInterface()
40         target_group = MPIProcessorGroup(interface, procs_target)
41         source_group = MPIProcessorGroup(interface, procs_source)
42         dec = InterpKernelDEC(source_group, target_group)
43         
44         mesh       =0
45         support    =0
46         paramesh   =0
47         parafield  =0
48         icocofield =0
49         data_dir = os.environ['MED_ROOT_DIR']
50         tmp_dir  = os.environ['TMP']
51         
52         if not tmp_dir or len(tmp_dir)==0:
53             tmp_dir = "/tmp"
54             pass
55         
56         filename_xml1 = os.path.join(data_dir, "share/salome/resources/med/square1_split")
57         filename_xml2 = os.path.join(data_dir, "share/salome/resources/med/square2_split")
58         
59         MPI_Barrier(MPI_COMM_WORLD)
60         if source_group.containsMyRank():
61             filename = filename_xml1 + str(rank+1) + ".med"
62             meshname = "Mesh_2_" + str(rank+1)
63             mesh=MEDLoader.ReadUMeshFromFile(filename,meshname,0)
64             paramesh=ParaMESH(mesh,source_group,"source mesh")
65             comptopo = ComponentTopology()
66             parafield = ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo)
67             parafield.getField().setNature(ConservativeVolumic)
68             nb_local=mesh.getNumberOfCells()
69             value = [1.0]*nb_local
70             parafield.getField().setValues(value)
71             icocofield = ICoCoMEDField(mesh,parafield.getField())
72             dec.setMethod("P0")
73             dec.attachLocalField(icocofield)
74             pass
75         else:
76             filename = filename_xml2 + str(rank - nproc_source + 1) + ".med"
77             meshname = "Mesh_3_" + str(rank - nproc_source + 1)
78             mesh=MEDLoader.ReadUMeshFromFile(filename,meshname,0)
79             paramesh=ParaMESH(mesh,target_group,"target mesh")
80             comptopo = ComponentTopology()
81             parafield = ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo)
82             parafield.getField().setNature(ConservativeVolumic)
83             nb_local=mesh.getNumberOfCells()
84             value = [0.0]*nb_local
85             parafield.getField().setValues(value)
86             icocofield = ICoCoMEDField(mesh,parafield.getField())
87             dec.setMethod("P0")
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.failUnless(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_Barrier(MPI_COMM_WORLD)
117         MPI_Finalize()
118         pass
119     pass
120
121 unittest.main()