PageRenderTime 60ms CodeModel.GetById 21ms RepoModel.GetById 0ms app.codeStats 0ms

/test/unit/fem/python/Form.py

https://bitbucket.org/pefarrell/dolfin
Python | 439 lines | 307 code | 88 blank | 44 comment | 13 complexity | e2dccbf16bf9b4fc6a1e73aa39dcfb41 MD5 | raw file
  1. """Unit tests for the fem interface"""
  2. # Copyright (C) 2011 Johan Hake
  3. #
  4. # This file is part of DOLFIN.
  5. #
  6. # DOLFIN is free software: you can redistribute it and/or modify
  7. # it under the terms of the GNU Lesser General Public License as published by
  8. # the Free Software Foundation, either version 3 of the License, or
  9. # (at your option) any later version.
  10. #
  11. # DOLFIN is distributed in the hope that it will be useful,
  12. # but WITHOUT ANY WARRANTY; without even the implied warranty of
  13. # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  14. # GNU Lesser General Public License for more details.
  15. #
  16. # You should have received a copy of the GNU Lesser General Public License
  17. # along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
  18. #
  19. # First added: 2011-12-02
  20. # Last changed: 2012-12-11
  21. #
  22. # Modified by Marie E. Rognes (meg@simula.no), 2012
  23. import unittest
  24. import numpy
  25. import ufl
  26. from dolfin import *
  27. class FormTest(unittest.TestCase):
  28. def setUp(self):
  29. self.mesh = UnitSquareMesh(10, 10)
  30. self.V = FunctionSpace(self.mesh, "CG", 1)
  31. self.f = Expression("sin(pi*x[0]*x[1])")
  32. self.v = TestFunction(self.V)
  33. self.u = TrialFunction(self.V)
  34. def test_assemble(self):
  35. ufl_form = self.f*self.u*self.v*dx
  36. dolfin_form = Form(ufl_form)
  37. ufc_form = dolfin_form._compiled_form
  38. A_ufl_norm = assemble(ufl_form).norm("frobenius")
  39. A_dolfin_norm = assemble(dolfin_form).norm("frobenius")
  40. A_ufc_norm = assemble(ufc_form, coefficients=[self.f],
  41. function_spaces=[self.V, self.V]).norm("frobenius")
  42. self.assertAlmostEqual(A_ufl_norm, A_dolfin_norm)
  43. self.assertAlmostEqual(A_ufl_norm, A_ufc_norm)
  44. class FormTestsOverManifolds(unittest.TestCase):
  45. def setUp(self):
  46. # 1D in 2D spaces
  47. self.square = UnitSquareMesh(2, 2)
  48. self.mesh1 = BoundaryMesh(self.square, "exterior")
  49. self.V1 = FunctionSpace(self.mesh1, "CG", 1)
  50. self.VV1 = VectorFunctionSpace(self.mesh1, "CG", 1)
  51. self.Q1 = FunctionSpace(self.mesh1, "DG", 0)
  52. # 2D in 3D spaces
  53. self.cube = UnitCubeMesh(2, 2, 2)
  54. self.mesh2 = BoundaryMesh(self.cube, "exterior")
  55. self.V2 = FunctionSpace(self.mesh2, "CG", 1)
  56. self.VV2 = VectorFunctionSpace(self.mesh2, "CG", 1)
  57. self.Q2 = FunctionSpace(self.mesh2, "DG", 0)
  58. def test_assemble_functional(self):
  59. u = Function(self.V1)
  60. u.vector()[:] = 1.0
  61. surfacearea = assemble(u*dx)
  62. self.assertAlmostEqual(surfacearea, 4.0)
  63. u = Function(self.V2)
  64. u.vector()[:] = 1.0
  65. surfacearea = assemble(u*dx)
  66. self.assertAlmostEqual(surfacearea, 6.0)
  67. f = Expression("1.0")
  68. u = interpolate(f, self.V1)
  69. surfacearea = assemble(u*dx)
  70. self.assertAlmostEqual(surfacearea, 4.0)
  71. f = Expression("1.0")
  72. u = interpolate(f, self.V2)
  73. surfacearea = assemble(u*dx)
  74. self.assertAlmostEqual(surfacearea, 6.0)
  75. def test_assemble_linear(self):
  76. u = Function(self.V1)
  77. w = TestFunction(self.Q1)
  78. u.vector()[:] = 0.5
  79. facetareas = MPI.sum(self.mesh1.mpi_comm(),
  80. assemble(u*w*dx).array().sum())
  81. self.assertAlmostEqual(facetareas, 2.0)
  82. u = Function(self.V2)
  83. w = TestFunction(self.Q2)
  84. u.vector()[:] = 0.5
  85. a = u*w*dx
  86. b = assemble(a)
  87. facetareas = MPI.sum(self.mesh2.mpi_comm(),
  88. assemble(u*w*dx).array().sum())
  89. self.assertAlmostEqual(facetareas, 3.0)
  90. mesh = UnitSquareMesh(8, 8)
  91. bdry = BoundaryMesh(mesh, "exterior")
  92. V = FunctionSpace(mesh, "CG", 1)
  93. u = TrialFunction(V)
  94. v = TestFunction(V)
  95. BV = FunctionSpace(bdry, "CG", 1)
  96. bu = TrialFunction(BV)
  97. bv = TestFunction(BV)
  98. a = MPI.sum(mesh.mpi_comm(),
  99. assemble(inner(u, v)*ds).array().sum())
  100. b = MPI.sum(bdry.mpi_comm(),
  101. assemble(inner(bu, bv)*dx).array().sum())
  102. self.assertAlmostEqual(a, b)
  103. def test_assemble_bilinear_1D_2D(self):
  104. # SubMesh not running in parallel
  105. if MPI.size(self.cube.mpi_comm()) > 1:
  106. return
  107. V = FunctionSpace(self.square, 'CG', 1)
  108. u = TrialFunction(V)
  109. v = TestFunction(V)
  110. bu = TrialFunction(self.V1)
  111. bv = TestFunction(self.V1)
  112. a = MPI.sum(self.square.mpi_comm(),
  113. assemble(inner(u, v)*ds).array().sum())
  114. b = MPI.sum(self.mesh1.mpi_comm(),
  115. assemble(inner(bu, bv)*dx).array().sum())
  116. self.assertAlmostEqual(a, b)
  117. bottom = CompiledSubDomain("near(x[1], 0.0)")
  118. foo = MPI.sum(self.square.mpi_comm(),
  119. abs(assemble(inner(grad(u)[0], grad(v)[0])*ds(0),
  120. exterior_facet_domains=bottom).array()).sum())
  121. BV = FunctionSpace(SubMesh(self.mesh1, bottom), "CG", 1)
  122. bu = TrialFunction(BV)
  123. bv = TestFunction(BV)
  124. bar = MPI.sum(self.mesh1.mpi_comm(),
  125. abs(assemble(inner(grad(bu), grad(bv))*dx).array()).sum())
  126. self.assertAlmostEqual(bar, foo)
  127. def test_assemble_bilinear_2D_3D(self):
  128. # SubMesh not running in parallel
  129. if MPI.size(self.cube.mpi_comm()) > 1:
  130. return
  131. V = FunctionSpace(self.cube, 'CG', 1)
  132. u = TrialFunction(V)
  133. v = TestFunction(V)
  134. bu = TrialFunction(self.V2)
  135. bv = TestFunction(self.V2)
  136. a = MPI.sum(self.mesh1.mpi_comm(),
  137. assemble(inner(u, v)*ds).array().sum())
  138. b = MPI.sum(self.mesh2.mpi_comm(),
  139. assemble(inner(bu, bv)*dx).array().sum())
  140. self.assertAlmostEqual(a, b)
  141. bottom = CompiledSubDomain("near(x[1], 0.0)")
  142. foo = MPI.sum(self.mesh1.mpi_comm(),
  143. abs(assemble(inner(grad(u)[0], grad(v)[0])*ds(0),
  144. exterior_facet_domains=bottom).array()).sum())
  145. BV = FunctionSpace(SubMesh(self.mesh1, bottom), "CG", 1)
  146. bu = TrialFunction(BV)
  147. bv = TestFunction(BV)
  148. bar = MPI.sum(self.mesh1.mpi_comm(),
  149. abs(assemble(inner(grad(bu), grad(bv))*dx).array()).sum())
  150. self.assertAlmostEqual(bar, foo)
  151. class FormTestsOverFunnySpaces(unittest.TestCase):
  152. def setUp(self):
  153. # Set-up meshes
  154. n = 16
  155. plane = CompiledSubDomain("near(x[1], 1.0)")
  156. self.square = UnitSquareMesh(n, n)
  157. # SubMesh not running in parallel
  158. if MPI.size(self.square.mpi_comm()) > 1:
  159. return
  160. self.square3d = SubMesh(BoundaryMesh(UnitCubeMesh(n, n, n), "exterior"), plane)
  161. # Define global normal and create orientation map
  162. global_normal = Expression(("0.0", "1.0", "0.0"))
  163. self.square3d.init_cell_orientations(global_normal)
  164. self.CG2 = VectorFunctionSpace(self.square, "CG", 1)
  165. self.CG3 = VectorFunctionSpace(self.square3d, "CG", 1)
  166. self.RT2 = FunctionSpace(self.square, "RT", 1)
  167. self.RT3 = FunctionSpace(self.square3d, "RT", 1)
  168. self.DG2 = FunctionSpace(self.square, "DG", 0)
  169. self.DG3 = FunctionSpace(self.square3d, "DG", 0)
  170. self.W2 = self.RT2*self.DG2
  171. self.W3 = self.RT3*self.DG3
  172. def test_basic_rt(self):
  173. # SubMesh not running in parallel
  174. if MPI.size(self.square.mpi_comm()) > 1:
  175. return
  176. f2 = Expression(("2.0", "1.0"))
  177. f3 = Expression(("1.0", "0.0", "2.0"))
  178. u2 = TrialFunction(self.RT2)
  179. u3 = TrialFunction(self.RT3)
  180. v2 = TestFunction(self.RT2)
  181. v3 = TestFunction(self.RT3)
  182. # Project
  183. pw2 = project(f2, self.RT2)
  184. pw3 = project(f3, self.RT3)
  185. pa2 = assemble(inner(pw2, pw2)*dx)
  186. pa3 = assemble(inner(pw3, pw3)*dx)
  187. # Project explicitly
  188. a2 = inner(u2, v2)*dx
  189. a3 = inner(u3, v3)*dx
  190. L2 = inner(f2, v2)*dx
  191. L3 = inner(f3, v3)*dx
  192. w2 = Function(self.RT2)
  193. w3 = Function(self.RT3)
  194. A2 = assemble(a2)
  195. b2 = assemble(L2)
  196. A3 = assemble(a3)
  197. b3 = assemble(L3)
  198. solve(A2, w2.vector(), b2)
  199. solve(A3, w3.vector(), b3)
  200. a2 = assemble(inner(w2, w2)*dx)
  201. a3 = assemble(inner(w3, w3)*dx)
  202. # Compare various results
  203. self.assertAlmostEqual((w2.vector() - pw2.vector()).norm("l2"), 0.0, \
  204. places=6)
  205. self.assertAlmostEqual(a3, 5.0)
  206. self.assertAlmostEqual(a2, a3)
  207. self.assertAlmostEqual(pa2, a2)
  208. self.assertAlmostEqual(pa2, pa3)
  209. def test_mixed_poisson_solve(self):
  210. # SubMesh not running in parallel
  211. if MPI.size(self.square.mpi_comm()) > 1:
  212. return
  213. f = Constant(1.0)
  214. # Solve mixed Poisson on standard unit square
  215. (sigma2, u2) = TrialFunctions(self.W2)
  216. (tau2, v2) = TestFunctions(self.W2)
  217. a = (inner(sigma2, tau2) + div(tau2)*u2 + div(sigma2)*v2)*dx
  218. L = f*v2*dx
  219. w2 = Function(self.W2)
  220. solve(a == L, w2)
  221. # Solve mixed Poisson on unit square in 3D
  222. (sigma3, u3) = TrialFunctions(self.W3)
  223. (tau3, v3) = TestFunctions(self.W3)
  224. a = (inner(sigma3, tau3) + div(tau3)*u3 + div(sigma3)*v3)*dx
  225. L = f*v3*dx
  226. w3 = Function(self.W3)
  227. solve(a == L, w3)
  228. # Check that results are about the same
  229. self.assertAlmostEqual(assemble(inner(w2, w2)*dx),
  230. assemble(inner(w3, w3)*dx))
  231. class TestGeometricQuantitiesOverManifolds(unittest.TestCase):
  232. def setUp(self):
  233. m = 3
  234. self.m = m
  235. self.mesh2 = BoundaryMesh(UnitCubeMesh(m, m, m), "exterior")
  236. # SubMesh not running in parallel
  237. if MPI.size(self.mesh2.mpi_comm()) > 1:
  238. return
  239. plane = CompiledSubDomain("near(x[1], 0.0)")
  240. self.mesh1 = BoundaryMesh(UnitSquareMesh(m, m), "exterior")
  241. self.bottom1 = SubMesh(self.mesh1, plane)
  242. self.bottom2 = SubMesh(self.mesh2, plane)
  243. line = CompiledSubDomain("near(x[0], 0.0)")
  244. self.mesh3 = BoundaryMesh(SubMesh(self.mesh2, plane), "exterior")
  245. self.bottom3 = SubMesh(self.mesh3, line)
  246. def test_normals_2D_1D(self):
  247. # SubMesh not running in parallel
  248. if MPI.size(self.mesh2.mpi_comm()) > 1:
  249. return
  250. "Testing assembly of normals for 1D meshes embedded in 2D"
  251. n = ufl.Cell("interval", geometric_dimension=2).n
  252. a = inner(n, n)*ds
  253. value_bottom1 = assemble(a, mesh=self.bottom1)
  254. self.assertAlmostEqual(value_bottom1, 2.0)
  255. b = inner(n('+'), n('+'))*dS
  256. b1 = assemble(b, mesh=self.bottom1)
  257. c = inner(n('+'), n('-'))*dS
  258. c1 = assemble(c, mesh=self.bottom1)
  259. self.assertAlmostEqual(b1, self.m-1)
  260. self.assertAlmostEqual(c1, - b1)
  261. def test_normals_3D_1D(self):
  262. "Testing assembly of normals for 1D meshes embedded in 3D"
  263. # SubMesh not running in parallel
  264. if MPI.size(self.mesh2.mpi_comm()) > 1:
  265. return
  266. n = ufl.Cell("interval", geometric_dimension=3).n
  267. a = inner(n, n)*ds
  268. v1 = assemble(a, mesh=self.bottom3)
  269. self.assertAlmostEqual(v1, 2.0)
  270. b = inner(n('+'), n('+'))*dS
  271. b1 = assemble(b, mesh=self.bottom3)
  272. c = inner(n('+'), n('-'))*dS
  273. c1 = assemble(c, mesh=self.bottom3)
  274. self.assertAlmostEqual(b1, self.m-1)
  275. self.assertAlmostEqual(c1, - b1)
  276. def test_normals_3D_2D(self):
  277. "Testing assembly of normals for 2D meshes embedded in 3D"
  278. # SubMesh not running in parallel
  279. if MPI.size(self.mesh2.mpi_comm()) > 1:
  280. return
  281. n = ufl.Cell("triangle", geometric_dimension=3).n
  282. a = inner(n, n)*ds
  283. v1 = assemble(a, mesh=self.bottom2)
  284. self.assertAlmostEqual(v1, 4.0)
  285. b = inner(n('+'), n('+'))*dS
  286. b1 = assemble(b, mesh=self.bottom2)
  287. c = inner(n('+'), n('-'))*dS
  288. c1 = assemble(c, mesh=self.bottom2)
  289. self.assertAlmostEqual(c1, - b1)
  290. def test_cell_volume(self):
  291. "Testing assembly of volume for embedded meshes"
  292. # SubMesh not running in parallel
  293. if MPI.size(self.mesh2.mpi_comm()) > 1:
  294. return
  295. volume = ufl.Cell("interval", geometric_dimension=2).volume
  296. a = volume*dx
  297. b = assemble(a, mesh=self.bottom1)
  298. self.assertAlmostEqual(b, 1.0/self.m)
  299. volume = ufl.Cell("interval", geometric_dimension=3).volume
  300. a = volume*dx
  301. b = assemble(a, mesh=self.bottom3)
  302. self.assertAlmostEqual(b, 1.0/self.m)
  303. volume = ufl.Cell("triangle", geometric_dimension=3).volume
  304. a = volume*dx
  305. b = assemble(a, mesh=self.bottom2)
  306. self.assertAlmostEqual(b, 1.0/(2*self.m*self.m))
  307. def test_circumradius(self):
  308. "Testing assembly of circumradius for embedded meshes"
  309. # SubMesh not running in parallel
  310. if MPI.size(self.mesh2.mpi_comm()) > 1:
  311. return
  312. r = ufl.Cell("interval", geometric_dimension=2).circumradius
  313. a = r*dx
  314. b = assemble(a, mesh=self.bottom1)
  315. self.assertAlmostEqual(b, 0.5*(1.0/self.m))
  316. r = ufl.Cell("interval", geometric_dimension=3).circumradius
  317. a = r*dx
  318. b = assemble(a, mesh=self.bottom3)
  319. self.assertAlmostEqual(b, 0.5*(1.0/self.m))
  320. r = ufl.Cell("triangle", geometric_dimension=2).circumradius
  321. a = r*dx
  322. b0 = assemble(a, mesh=UnitSquareMesh(self.m, self.m))
  323. r = ufl.Cell("triangle", geometric_dimension=3).circumradius
  324. a = r*dx
  325. b1 = assemble(a, mesh=self.bottom2)
  326. self.assertAlmostEqual(b0, b1)
  327. def test_facetarea(self):
  328. "Testing assembly of facet area for embedded meshes"
  329. # SubMesh not running in parallel
  330. if MPI.size(self.mesh2.mpi_comm()) > 1:
  331. return
  332. area = ufl.Cell("interval", geometric_dimension=2).facet_area
  333. a = area*ds
  334. b = assemble(a, mesh=self.bottom1)
  335. self.assertAlmostEqual(b, 2.0)
  336. area = ufl.Cell("interval", geometric_dimension=3).facet_area
  337. a = area*ds
  338. b = assemble(a, mesh=self.bottom3)
  339. self.assertAlmostEqual(b, 2.0)
  340. area = ufl.Cell("triangle", geometric_dimension=2).facet_area
  341. a = area*ds
  342. b0 = assemble(a, mesh=UnitSquareMesh(self.m, self.m))
  343. area = ufl.Cell("triangle", geometric_dimension=3).facet_area
  344. a = area*ds
  345. b1 = assemble(a, mesh=self.bottom2)
  346. self.assertAlmostEqual(b0, b1)
  347. if __name__ == "__main__":
  348. print ""
  349. print "Testing PyDOLFIN Form operations"
  350. print "------------------------------------------------"
  351. unittest.main()