Back to index

salome-med  6.5.0
test_InterpKernelDEC.py
Go to the documentation of this file.
00001 #!/usr/bin/env python
00002 #  -*- coding: iso-8859-1 -*-
00003 # Copyright (C) 2007-2012  CEA/DEN, EDF R&D
00004 #
00005 # This library is free software; you can redistribute it and/or
00006 # modify it under the terms of the GNU Lesser General Public
00007 # License as published by the Free Software Foundation; either
00008 # version 2.1 of the License.
00009 #
00010 # This library is distributed in the hope that it will be useful,
00011 # but WITHOUT ANY WARRANTY; without even the implied warranty of
00012 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00013 # Lesser General Public License for more details.
00014 #
00015 # You should have received a copy of the GNU Lesser General Public
00016 # License along with this library; if not, write to the Free Software
00017 # Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
00018 #
00019 # See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
00020 #
00021 
00022 from ParaMEDMEM import *
00023 import sys, os
00024 import unittest
00025 import math
00026 
00027 class ParaMEDMEMBasicsTest(unittest.TestCase):
00028     def testInterpKernelDEC_2D(self):
00029         MPI_Init(sys.argv)
00030         size = MPI_Comm_size(MPI_COMM_WORLD)
00031         rank = MPI_Comm_rank(MPI_COMM_WORLD)
00032         if size != 5:
00033             raise RuntimeError, "Expect MPI_COMM_WORLD size == 5"
00034         print rank
00035         nproc_source = 3
00036         procs_source = range( nproc_source )
00037         procs_target = range( size - nproc_source + 1, size)
00038         
00039         interface = CommInterface()
00040         target_group = MPIProcessorGroup(interface, procs_target)
00041         source_group = MPIProcessorGroup(interface, procs_source)
00042         dec = InterpKernelDEC(source_group, target_group)
00043         
00044         mesh       =0
00045         support    =0
00046         paramesh   =0
00047         parafield  =0
00048         icocofield =0
00049         data_dir = os.environ['MED_ROOT_DIR']
00050         tmp_dir  = os.environ['TMP']
00051         
00052         if not tmp_dir or len(tmp_dir)==0:
00053             tmp_dir = "/tmp"
00054             pass
00055         
00056         filename_xml1 = os.path.join(data_dir, "share/salome/resources/med/square1_split")
00057         filename_xml2 = os.path.join(data_dir, "share/salome/resources/med/square2_split")
00058         
00059         MPI_Barrier(MPI_COMM_WORLD)
00060         if source_group.containsMyRank():
00061             filename = filename_xml1 + str(rank+1) + ".med"
00062             meshname = "Mesh_2_" + str(rank+1)
00063             mesh=MEDLoader.ReadUMeshFromFile(filename,meshname,0)
00064             paramesh=ParaMESH(mesh,source_group,"source mesh")
00065             comptopo = ComponentTopology()
00066             parafield = ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo)
00067             parafield.getField().setNature(ConservativeVolumic)
00068             nb_local=mesh.getNumberOfCells()
00069             value = [1.0]*nb_local
00070             parafield.getField().setValues(value)
00071             icocofield = ICoCoMEDField(mesh,parafield.getField())
00072             dec.setMethod("P0")
00073             dec.attachLocalField(icocofield)
00074             pass
00075         else:
00076             filename = filename_xml2 + str(rank - nproc_source + 1) + ".med"
00077             meshname = "Mesh_3_" + str(rank - nproc_source + 1)
00078             mesh=MEDLoader.ReadUMeshFromFile(filename,meshname,0)
00079             paramesh=ParaMESH(mesh,target_group,"target mesh")
00080             comptopo = ComponentTopology()
00081             parafield = ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo)
00082             parafield.getField().setNature(ConservativeVolumic)
00083             nb_local=mesh.getNumberOfCells()
00084             value = [0.0]*nb_local
00085             parafield.getField().setValues(value)
00086             icocofield = ICoCoMEDField(mesh,parafield.getField())
00087             dec.setMethod("P0")
00088             dec.attachLocalField(icocofield)
00089             pass
00090         
00091         if source_group.containsMyRank():
00092             field_before_int = parafield.getVolumeIntegral(0,True)
00093             dec.synchronize()
00094             dec.setForcedRenormalization(False)
00095             dec.sendData()
00096             dec.recvData()
00097             field_after_int=parafield.getVolumeIntegral(0,True);
00098             self.failUnless(math.fabs(field_after_int-field_before_int)<1e-8)
00099             pass
00100         else:
00101             dec.synchronize()
00102             dec.setForcedRenormalization(False)
00103             dec.recvData()
00104             dec.sendData()
00105             pass
00106         ## end
00107         interface = 0
00108         target_group = 0
00109         source_group = 0
00110         dec = 0
00111         mesh       =0
00112         support    =0
00113         paramesh   =0
00114         parafield  =0
00115         icocofield =0
00116         MPI_Barrier(MPI_COMM_WORLD)
00117         MPI_Finalize()
00118         pass
00119     pass
00120 
00121 unittest.main()