def test_higher_order_coordinate_map(points, celltype): """ Computes physical coordinates of a cell, based on the coordinate map. """ cells = np.array([range(len(points))]) mesh = Mesh(MPI.comm_world, celltype, points, cells, [], GhostMode.none) V = FunctionSpace(mesh, ("Lagrange", mesh.degree())) X = V.element.dof_reference_coordinates() coord_dofs = mesh.coordinate_dofs().entity_points() x_g = mesh.geometry.points cmap = fem.create_coordinate_map(mesh.ufl_domain()) x_coord_new = np.zeros([len(points), mesh.geometry.dim]) i = 0 for node in range(len(points)): x_coord_new[i] = x_g[coord_dofs[0, node], :mesh.geometry.dim] i += 1 x = np.zeros(X.shape) cmap.push_forward(x, X, x_coord_new) assert (np.allclose(x[:, 0], X[:, 0])) assert (np.allclose(x[:, 1], 2 * X[:, 1])) if mesh.geometry.dim == 3: assert (np.allclose(x[:, 2], 3 * X[:, 2]))
def test_third_order_tri(): # *---*---*---* 3--11--10--2 # | \ | | \ | # * * * * 8 7 15 13 # | \ | | \ | # * * * * 9 14 6 12 # | \ | | \ | # *---*---*---* 0--4---5---1 for H in (1.0, 2.0): for Z in (0.0, 0.5): L = 1 points = np.array([ [0, 0, 0], [L, 0, 0], [L, H, Z], [0, H, Z], # 0, 1, 2, 3 [L / 3, 0, 0], [2 * L / 3, 0, 0], # 4, 5 [2 * L / 3, H / 3, 0], [L / 3, 2 * H / 3, 0], # 6, 7 [0, 2 * H / 3, 0], [0, H / 3, 0], # 8, 9 [2 * L / 3, H, Z], [L / 3, H, Z], # 10, 11 [L, H / 3, 0], [L, 2 * H / 3, 0], # 12, 13 [L / 3, H / 3, 0], # 14 [2 * L / 3, 2 * H / 3, 0] ]) # 15 cells = np.array([[0, 1, 3, 4, 5, 6, 7, 8, 9, 14], [1, 2, 3, 12, 13, 10, 11, 7, 6, 15]]) cells = permute_cell_ordering( cells, permutation_vtk_to_dolfin(CellType.triangle, cells.shape[1])) mesh = Mesh(MPI.comm_world, CellType.triangle, points, cells, [], GhostMode.none) def e2(x): return x[2] + x[0] * x[1] degree = mesh.degree() # Interpolate function V = FunctionSpace(mesh, ("CG", degree)) u = Function(V) cmap = fem.create_coordinate_map(mesh.ufl_domain()) mesh.geometry.coord_mapping = cmap u.interpolate(e2) intu = assemble_scalar(u * dx(metadata={"quadrature_degree": 40})) intu = MPI.sum(mesh.mpi_comm(), intu) nodes = [0, 9, 8, 3] ref = sympy_scipy(points, nodes, L, H) assert ref == pytest.approx(intu, rel=1e-6)
def test_fourth_order_tri(): L = 1 # *--*--*--*--* 3-21-20-19--2 # | \ | | \ | # * * * * * 10 9 24 23 18 # | \ | | \ | # * * * * * 11 15 8 22 17 # | \ | | \ | # * * * * * 12 13 14 7 16 # | \ | | \ | # *--*--*--*--* 0--4--5--6--1 for H in (1.0, 2.0): for Z in (0.0, 0.5): points = np.array( [[0, 0, 0], [L, 0, 0], [L, H, Z], [0, H, Z], # 0, 1, 2, 3 [L / 4, 0, 0], [L / 2, 0, 0], [3 * L / 4, 0, 0], # 4, 5, 6 [3 / 4 * L, H / 4, Z / 2], [L / 2, H / 2, 0], # 7, 8 [L / 4, 3 * H / 4, 0], [0, 3 * H / 4, 0], # 9, 10 [0, H / 2, 0], [0, H / 4, Z / 2], # 11, 12 [L / 4, H / 4, Z / 2], [L / 2, H / 4, Z / 2], [L / 4, H / 2, 0], # 13, 14, 15 [L, H / 4, Z / 2], [L, H / 2, 0], [L, 3 * H / 4, 0], # 16, 17, 18 [3 * L / 4, H, Z], [L / 2, H, Z], [L / 4, H, Z], # 19, 20, 21 [3 * L / 4, H / 2, 0], [3 * L / 4, 3 * H / 4, 0], # 22, 23 [L / 2, 3 * H / 4, 0]] # 24 ) cells = np.array([[0, 1, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15], [1, 2, 3, 16, 17, 18, 19, 20, 21, 9, 8, 7, 22, 23, 24]]) cells = permute_cell_ordering(cells, permutation_vtk_to_dolfin(CellType.triangle, cells.shape[1])) mesh = Mesh(MPI.comm_world, CellType.triangle, points, cells, [], GhostMode.none) def e2(x): return x[2] + x[0] * x[1] degree = mesh.degree() # Interpolate function V = FunctionSpace(mesh, ("CG", degree)) u = Function(V) cmap = fem.create_coordinate_map(mesh.ufl_domain()) mesh.geometry.coord_mapping = cmap u.interpolate(e2) intu = assemble_scalar(u * dx(metadata={"quadrature_degree": 50})) intu = MPI.sum(mesh.mpi_comm(), intu) nodes = [0, 3, 10, 11, 12] ref = sympy_scipy(points, nodes, L, H) assert ref == pytest.approx(intu, rel=1e-4)
def test_second_order_tri(): # Test second order mesh by computing volume of two cells # *-----*-----* 3----6-----2 # | \ | | \ | # | \ | | \ | # * * * 7 8 5 # | \ | | \ | # | \ | | \ | # *-----*-----* 0----4-----1 for H in (1.0, 2.0): for Z in (0.0, 0.5): L = 1 points = np.array([[0, 0, 0], [L, 0, 0], [L, H, Z], [0, H, Z], [L / 2, 0, 0], [L, H / 2, 0], [L / 2, H, Z], [0, H / 2, 0], [L / 2, H / 2, 0]]) cells = np.array([[0, 1, 3, 4, 8, 7], [1, 2, 3, 5, 6, 8]]) cells = permute_cell_ordering( cells, permutation_vtk_to_dolfin(CellType.triangle, cells.shape[1])) mesh = Mesh(MPI.comm_world, CellType.triangle, points, cells, [], GhostMode.none) def e2(x): return x[2] + x[0] * x[1] degree = mesh.degree() # Interpolate function V = FunctionSpace(mesh, ("CG", degree)) u = Function(V) cmap = fem.create_coordinate_map(mesh.ufl_domain()) mesh.geometry.coord_mapping = cmap u.interpolate(e2) intu = assemble_scalar( u * dx(mesh, metadata={"quadrature_degree": 20})) intu = MPI.sum(mesh.mpi_comm(), intu) nodes = [0, 3, 7] ref = sympy_scipy(points, nodes, L, H) assert ref == pytest.approx(intu, rel=1e-6)
def test_mesh_construction_pygmsh(): pygmsh = pytest.importorskip("pygmsh") if MPI.rank(MPI.comm_world) == 0: geom = pygmsh.opencascade.Geometry() geom.add_ball([0.0, 0.0, 0.0], 1.0, char_length=0.2) pygmsh_mesh = pygmsh.generate_mesh(geom) points, cells = pygmsh_mesh.points, pygmsh_mesh.cells else: points = np.zeros([0, 3]) cells = { "tetra": np.zeros([0, 4], dtype=np.int64), "triangle": np.zeros([0, 3], dtype=np.int64), "line": np.zeros([0, 2], dtype=np.int64) } mesh = Mesh(MPI.comm_world, dolfinx.cpp.mesh.CellType.tetrahedron, points, cells['tetra'], [], cpp.mesh.GhostMode.none) assert mesh.degree() == 1 assert mesh.geometry.dim == 3 assert mesh.topology.dim == 3 mesh = Mesh(MPI.comm_world, dolfinx.cpp.mesh.CellType.triangle, points, cells['triangle'], [], cpp.mesh.GhostMode.none) assert mesh.degree() == 1 assert mesh.geometry.dim == 3 assert mesh.topology.dim == 2 mesh = Mesh(MPI.comm_world, dolfinx.cpp.mesh.CellType.interval, points, cells['line'], [], cpp.mesh.GhostMode.none) assert mesh.degree() == 1 assert mesh.geometry.dim == 3 assert mesh.topology.dim == 1 if MPI.rank(MPI.comm_world) == 0: print("Generate mesh") geom = pygmsh.opencascade.Geometry() geom.add_ball([0.0, 0.0, 0.0], 1.0, char_length=0.2) pygmsh_mesh = pygmsh.generate_mesh( geom, extra_gmsh_arguments=['-order', '2']) points, cells = pygmsh_mesh.points, pygmsh_mesh.cells print("End Generate mesh", cells.keys()) else: points = np.zeros([0, 3]) cells = { "tetra10": np.zeros([0, 10], dtype=np.int64), "triangle6": np.zeros([0, 6], dtype=np.int64), "line3": np.zeros([0, 3], dtype=np.int64) } mesh = Mesh(MPI.comm_world, dolfinx.cpp.mesh.CellType.tetrahedron, points, cells['tetra10'], [], cpp.mesh.GhostMode.none) assert mesh.degree() == 2 assert mesh.geometry.dim == 3 assert mesh.topology.dim == 3 mesh = Mesh(MPI.comm_world, dolfinx.cpp.mesh.CellType.triangle, points, cells['triangle6'], [], cpp.mesh.GhostMode.none) assert mesh.degree() == 2 assert mesh.geometry.dim == 3 assert mesh.topology.dim == 2