Back to index

salome-med  6.5.0
test_NonCoincidentDEC.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 
00025 MPI_Init(sys.argv)
00026 
00027 size = MPI_Comm_size(MPI_COMM_WORLD)
00028 rank = MPI_Comm_rank(MPI_COMM_WORLD)
00029 if size != 5:
00030     raise RuntimeError, "Expect MPI_COMM_WORLD size == 5"
00031 
00032 nproc_source = 3
00033 procs_source = range( nproc_source )
00034 procs_target = range( size - nproc_source + 1, size)
00035 
00036 interface = CommInterface()
00037 
00038 target_group = MPIProcessorGroup(interface, procs_target)
00039 source_group = MPIProcessorGroup(interface, procs_source)
00040 
00041 source_mesh= 0
00042 target_mesh= 0
00043 parasupport= 0
00044 mesh       = 0
00045 support    = 0
00046 field      = 0
00047 paramesh   = 0
00048 parafield  = 0
00049 icocofield = 0
00050 
00051 dec = NonCoincidentDEC(source_group, target_group)
00052 
00053 data_dir = os.environ['MED_ROOT_DIR']
00054 tmp_dir  = os.environ['TMP']
00055 if tmp_dir == '':
00056     tmp_dir = "/tmp"
00057     pass
00058 
00059 filename_xml1 = data_dir + "/share/salome/resources/med/square1_split"
00060 filename_xml2 = data_dir + "/share/salome/resources/med/square2_split"
00061 
00062 MPI_Barrier(MPI_COMM_WORLD)
00063     
00064 if source_group.containsMyRank():
00065 
00066     filename = filename_xml1 + str(rank+1) + ".med"
00067     meshname = "Mesh_2_" + str(rank+1)
00068 
00069     mesh = MESH(MED_DRIVER, filename, meshname)
00070     support = SUPPORT(mesh, "all elements", MED_CELL)
00071     paramesh = ParaMESH(mesh, source_group, "source mesh")
00072 
00073     parasupport = UnstructuredParaSUPPORT( support, source_group)
00074     comptopo = ComponentTopology()
00075 
00076     parafield = ParaFIELD(parasupport, comptopo)
00077 
00078     nb_local = support.getNumberOfElements(MED_ALL_ELEMENTS);
00079 
00080     value = [1.0]*nb_local
00081 
00082     parafield.getField().setValue(value)
00083     icocofield = ICoCo_MEDField(paramesh,parafield)
00084     dec.attachLocalField(icocofield,'P0')
00085     pass
00086 
00087 if target_group.containsMyRank():
00088 
00089     filename = filename_xml2 + str(rank - nproc_source + 1) + ".med"
00090     meshname = "Mesh_3_" + str(rank - nproc_source + 1)
00091 
00092     mesh = MESH(MED_DRIVER, filename, meshname)
00093     support = SUPPORT(mesh, "all elements", MED_CELL)
00094     paramesh = ParaMESH(mesh, target_group, "target mesh")
00095 
00096     parasupport = UnstructuredParaSUPPORT( support, target_group)
00097     comptopo = ComponentTopology()
00098     parafield = ParaFIELD(parasupport, comptopo)
00099 
00100     nb_local = support.getNumberOfElements(MED_ALL_ELEMENTS)
00101     value = [0.0]*nb_local
00102 
00103     parafield.getField().setValue(value)
00104     icocofield = ICoCo_MEDField(paramesh,parafield)
00105 
00106     dec.attachLocalField(icocofield, 'P0')
00107     pass
00108 
00109 field_before_int = [0.0]
00110 field_after_int = [0.0]
00111 
00112 if source_group.containsMyRank():
00113 
00114     field_before_int = [parafield.getVolumeIntegral(1)]
00115     MPI_Bcast(field_before_int, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD);
00116     dec.synchronize()
00117     print "DEC usage"
00118     dec.setForcedRenormalization(False)
00119 
00120     dec.sendData()
00121     pass
00122     
00123 if target_group.containsMyRank():
00124 
00125     MPI_Bcast(field_before_int, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD)
00126     dec.synchronize()
00127     dec.setForcedRenormalization(False)
00128     dec.recvData()
00129     field_after_int = [parafield.getVolumeIntegral(1)]
00130     pass
00131 
00132 MPI_Bcast(field_before_int, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD)
00133 MPI_Bcast(field_after_int , 1, MPI_DOUBLE, size-1, MPI_COMM_WORLD)
00134 
00135 epsilon = 1e-6
00136 if abs(field_before_int[0] - field_after_int[0]) > epsilon:
00137     print "Field before is not equal field after: %s != %s"%\
00138           (field_before_int[0],field_after_int[0])
00139     pass
00140 
00141 
00142 MPI_Barrier(MPI_COMM_WORLD)
00143 MPI_Finalize()
00144 print "# End of testNonCoincidentDEC"