/[escript]/trunk/finley/test/python/run_inputOutput.py
ViewVC logotype

Contents of /trunk/finley/test/python/run_inputOutput.py

Parent Directory Parent Directory | Revision Log Revision Log


Revision 4154 - (show annotations)
Tue Jan 22 09:30:23 2013 UTC (6 years, 10 months ago) by jfenwick
File MIME type: text/x-python
File size: 9282 byte(s)
Round 1 of copyright fixes
1
2 ##############################################################################
3 #
4 # Copyright (c) 2003-2013 by University of Queensland
5 # http://www.uq.edu.au
6 #
7 # Primary Business: Queensland, Australia
8 # Licensed under the Open Software License version 3.0
9 # http://www.opensource.org/licenses/osl-3.0.php
10 #
11 # Development until 2012 by Earth Systems Science Computational Center (ESSCC)
12 # Development since 2012 by School of Earth Sciences
13 #
14 ##############################################################################
15
16 from __future__ import print_function
17
18 __copyright__="""Copyright (c) 2003-2013 by University of Queensland
19 http://www.uq.edu.au
20 Primary Business: Queensland, Australia"""
21 __license__="""Licensed under the Open Software License version 3.0
22 http://www.opensource.org/licenses/osl-3.0.php"""
23 __url__="https://launchpad.net/escript-finley"
24
25 """
26 Test suite for input and output of meshes and data objects
27
28 :remark:
29
30 :var __author__: name of author
31 :var __licence__: licence agreement
32 :var __url__: url entry point on documentation
33 :var __version__: version
34 :var __date__: date of the version
35 """
36
37 import unittest, sys
38
39 from esys.escript import *
40 from esys.finley import Rectangle, Brick, LoadMesh, ReadMesh, GetMeshFromFile, ReadGmsh
41
42 try:
43 FINLEY_WORKDIR=os.environ['FINLEY_WORKDIR']
44 except KeyError:
45 FINLEY_WORKDIR='.'
46
47 try:
48 FINLEY_TEST_DATA=os.environ['FINLEY_TEST_DATA']
49 except KeyError:
50 FINLEY_TEST_DATA='.'
51
52 FINLEY_TEST_MESH_PATH=os.path.join(FINLEY_TEST_DATA,"data_meshes")
53
54 REL_TOL=1.e-6
55
56 # Number of elements scales up with number of MPI processes
57 NE0 = 7 * getMPISizeWorld()
58 NE1 = 11
59 NE2 = 5
60
61 class InputOutput(unittest.TestCase):
62
63 # Check that two domains are equal using Fourier integrals
64 # We cannot compare the X coordinates since they are on different domains
65 def domainsEqual(self, m1, m2, nft=100):
66 self.assertTrue(m1.getDim() == m2.getDim(), "Dimensions differ")
67 self.assertTrue(m1.getNumDataPointsGlobal() == m2.getNumDataPointsGlobal(), "Num data points differ")
68 for tagName in m1.showTagNames().split(", "):
69 self.assertTrue(m2.isValidTagName(tagName), "m1 has a tag '%s' not present in m2" % tagName)
70 for tagName in m2.showTagNames().split(", "):
71 self.assertTrue(m1.isValidTagName(tagName), "m2 has a tag '%s' not present in m1" % tagName)
72 self.assertTrue(m1.getTag(tagName) == m2.getTag(tagName), "values of tag '%s' differ" % tagName)
73 for fs in ["Solution", "ReducedSolution", "Function", "ReducedFunction", "ContinuousFunction", "ReducedContinuousFunction"]:
74 fs1 = eval("%s(m1)" % fs)
75 fs2 = eval("%s(m2)" % fs)
76 x1 = fs1.getX()
77 x2 = fs2.getX()
78 for n in range(1, nft+1):
79 integ1 = integrate(sin(n*x1))
80 integ2 = integrate(sin(n*x2))
81 self.assertTrue(Lsup(abs(integ1-integ2)) <= REL_TOL, "integrals for n=%d differ" % n)
82 return True
83
84 # Does optimize=True change Rectangle for order=1?
85 def test_Rectangle_optimize_order1(self):
86 mydomain1 = Rectangle(n0=NE0, n1=NE1, order=1, l0=1., l1=1., optimize=False,useElementsOnFace=0)
87 mydomain2 = Rectangle(n0=NE0, n1=NE1, order=1, l0=1., l1=1., optimize=True,useElementsOnFace=0)
88 self.domainsEqual(mydomain1, mydomain2)
89
90 # Does optimize=True change Rectangle for order=2?
91 def test_Rectangle_optimize_order2(self):
92 mydomain1 = Rectangle(n0=NE0, n1=NE1, order=2, l0=1., l1=1., optimize=False,useElementsOnFace=0)
93 mydomain2 = Rectangle(n0=NE0, n1=NE1, order=2, l0=1., l1=1., optimize=True,useElementsOnFace=0)
94 self.domainsEqual(mydomain1, mydomain2)
95
96 # Does optimize=True change Rectangle for order=-1?
97 def test_Rectangle_optimize_macro(self):
98 mydomain1 = Rectangle(n0=NE0, n1=NE1, order=-1, l0=1., l1=1., optimize=False, useElementsOnFace=0)
99 mydomain2 = Rectangle(n0=NE0, n1=NE1, order=-1, l0=1., l1=1., optimize=True,useElementsOnFace=0)
100 self.domainsEqual(mydomain1, mydomain2)
101
102 # Does optimize=True change Brick for order=1?
103 def test_Brick_optimize_order1(self):
104 mydomain1 = Brick(n0=NE0, n1=NE1, n2=NE2, order=1, l0=1., l1=1., l2=1., optimize=False,useElementsOnFace=0)
105 mydomain2 = Brick(n0=NE0, n1=NE1, n2=NE2, order=1, l0=1., l1=1., l2=1., optimize=True,useElementsOnFace=0)
106 self.domainsEqual(mydomain1, mydomain2)
107
108 # Does optimize=True change Brick for order=2?
109 def test_Brick_optimize_order2(self):
110 mydomain1 = Brick(n0=NE0, n1=NE1, n2=NE2, order=2, l0=1., l1=1., l2=1., optimize=False,useElementsOnFace=0)
111 mydomain2 = Brick(n0=NE0, n1=NE1, n2=NE2, order=2, l0=1., l1=1., l2=1., optimize=True,useElementsOnFace=0)
112 self.domainsEqual(mydomain1, mydomain2)
113 # Does optimize=True change Brick for order=-1?
114 def test_Brick_optimize_macro(self):
115 mydomain1 = Brick(n0=NE0, n1=NE1, n2=NE2, order=-1, l0=1., l1=1., l2=1., optimize=False,useElementsOnFace=0)
116 mydomain2 = Brick(n0=NE0, n1=NE1, n2=NE2, order=-1, l0=1., l1=1., l2=1., optimize=True,useElementsOnFace=0)
117 self.domainsEqual(mydomain1, mydomain2)
118
119 def test_data_dump_to_NetCDF_rectangle(self):
120 if loadIsConfigured():
121 mydomain1 = Rectangle(n0=NE0, n1=NE1, order=1, l0=1., l1=1., optimize=False,useElementsOnFace=0)
122 d1=Data(mydomain1.getMPIRank(), Function(mydomain1))
123 d1.expand()
124 dumpfile=os.path.join(FINLEY_WORKDIR, "tempfile.dump.nc")
125 d1.dump(dumpfile)
126 d2=load(dumpfile, mydomain1)
127 self.assertTrue(Lsup(abs(d1-d2)) <= REL_TOL, "data objects differ")
128
129 def test_data_dump_to_NetCDF_brick(self):
130 if loadIsConfigured():
131 mydomain1 = Brick(n0=NE0, n1=NE1, n2=NE2, order=2, l0=1., l1=1., l2=1., optimize=False)
132 d1=Data(mydomain1.getMPIRank(), Function(mydomain1))
133 d1.expand()
134 dumpfile=os.path.join(FINLEY_WORKDIR, "tempfile.dump.nc")
135 d1.dump(dumpfile)
136 d2=load(dumpfile, mydomain1)
137 self.assertTrue(Lsup(abs(d1-d2)) <= REL_TOL, "data objects differ")
138
139 def test_mesh_dump_to_NetCDF_rectangle(self):
140 if loadIsConfigured():
141 mydomain1 = Rectangle(n0=NE0, n1=NE1, order=1, l0=1., l1=1., optimize=False)
142 dumpfile=os.path.join(FINLEY_WORKDIR, "tempfile.mesh.nc")
143 mydomain1.dump(dumpfile)
144 mydomain2=LoadMesh(dumpfile)
145 self.domainsEqual(mydomain1, mydomain2)
146
147 def test_gmshTags(self):
148 if getEscriptParamInt('MPIBUILD',0)==0:
149 dom=ReadGmsh(os.path.join(FINLEY_TEST_MESH_PATH, "tagtest.msh"),2)
150 tags=dom.showTagNames().split(', ')
151 self.assertEqual(tags,['tag1', 'tag2', 'tag3'],'error with tags')
152 self.assertEqual(dom.getTag('tag1'),1,'error with tag1')
153 self.assertEqual(dom.getTag('tag2'),2,'error with tag2')
154 self.assertEqual(dom.getTag('tag3'),3,'error with tag3')
155 self.assertRaises(RuntimeError, dom.getTag, 'tag4')
156
157 else:
158 print("Test supressed due to MPI build")
159 def test_flyTags(self):
160 dom=ReadMesh(os.path.join(FINLEY_TEST_MESH_PATH, "rectangle_8x10.fly"))
161 tags=dom.showTagNames().split(', ')
162 self.assertEqual(tags,['top', 'bottom', 'left', 'right'])
163 self.assertEqual(dom.getTag('top'),20,'error with top')
164 self.assertEqual(dom.getTag('bottom'),10,'error with bottom,')
165 self.assertEqual(dom.getTag('left'),1,'error with left')
166 self.assertEqual(dom.getTag('right'),2,'error with reight')
167 self.assertRaises(RuntimeError, dom.getTag, 'tag4')
168
169
170 def test_mesh_dump_to_NetCDF_brick(self):
171 if loadIsConfigured():
172 mydomain1 = Brick(n0=NE0, n1=NE1, n2=NE2, order=2, l0=1., l1=1., l2=1., optimize=False)
173 dumpfile=os.path.join(FINLEY_WORKDIR, "tempfile.mesh.nc")
174 mydomain1.dump(dumpfile)
175 mydomain2=LoadMesh(dumpfile)
176 self.domainsEqual(mydomain1, mydomain2)
177
178 def test_mesh_read_rectangle_from_finley_file(self):
179 if getMPISizeWorld() < 16:
180 mydomain1 = Rectangle(n0=8, n1=10, order=1, l0=1., l1=1., optimize=False)
181 mydomain2 = ReadMesh(os.path.join(FINLEY_TEST_MESH_PATH,"rectangle_8x10.fly"))
182 self.domainsEqual(mydomain1, mydomain2)
183
184 def test_mesh_read_brick_from_finley_file(self):
185 if getMPISizeWorld() < 16:
186 mydomain1 = Brick(n0=8, n1=10, n2=12, order=1, l0=1., l1=1., l2=1., optimize=False)
187 mydomain2 = ReadMesh(os.path.join(FINLEY_TEST_MESH_PATH,"brick_8x10x12.fly"))
188 self.domainsEqual(mydomain1, mydomain2)
189
190 def test_GetMeshFromFile(self):
191 if getMPISizeWorld() <2:
192 m=GetMeshFromFile(os.path.join(FINLEY_TEST_MESH_PATH,'tet10_gmsh.msh'), numDim=3)
193 del m
194 m=GetMeshFromFile(os.path.join(FINLEY_TEST_MESH_PATH, 'tet10.fly'))
195 # now we try some params
196 m=GetMeshFromFile(os.path.join(FINLEY_TEST_MESH_PATH,'tet10_gmsh.msh'), numDim=3, integrationOrder=2)
197
198 if __name__ == '__main__':
199 suite = unittest.TestSuite()
200 suite.addTest(unittest.makeSuite(InputOutput))
201 s=unittest.TextTestRunner(verbosity=2).run(suite)
202 if not s.wasSuccessful(): sys.exit(1)
203

  ViewVC Help
Powered by ViewVC 1.1.26