def test_assemble_derivatives(): """This test checks the original_coefficient_positions, which may change under differentiation (some coefficients and constants are eliminated)""" mesh = create_unit_square(MPI.COMM_WORLD, 12, 12) Q = FunctionSpace(mesh, ("Lagrange", 1)) u = Function(Q) v = ufl.TestFunction(Q) du = ufl.TrialFunction(Q) b = Function(Q) c1 = Constant(mesh, np.array([[1.0, 0.0], [3.0, 4.0]], PETSc.ScalarType)) c2 = Constant(mesh, PETSc.ScalarType(2.0)) b.x.array[:] = 2.0 # derivative eliminates 'u' and 'c1' L = ufl.inner(c1, c1) * v * dx + c2 * b * inner(u, v) * dx a = form(derivative(L, u, du)) A1 = assemble_matrix(a) A1.assemble() a = form(c2 * b * inner(du, v) * dx) A2 = assemble_matrix(a) A2.assemble() assert (A1 - A2).norm() == pytest.approx(0.0, rel=1e-12, abs=1e-12)
def test_basic_assembly_constant(mode): """Tests assembly with Constant The following test should be sensitive to order of flattening the matrix-valued constant. """ mesh = create_unit_square(MPI.COMM_WORLD, 5, 5, ghost_mode=mode) V = FunctionSpace(mesh, ("Lagrange", 1)) u, v = ufl.TrialFunction(V), ufl.TestFunction(V) c = Constant(mesh, np.array([[1.0, 2.0], [5.0, 3.0]], PETSc.ScalarType)) a = inner(c[1, 0] * u, v) * dx + inner(c[1, 0] * u, v) * ds L = inner(c[1, 0], v) * dx + inner(c[1, 0], v) * ds a, L = form(a), form(L) # Initial assembly A1 = assemble_matrix(a) A1.assemble() b1 = assemble_vector(L) b1.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) c.value = [[1.0, 2.0], [3.0, 4.0]] A2 = assemble_matrix(a) A2.assemble() b2 = assemble_vector(L) b2.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) assert (A1 * 3.0 - A2 * 5.0).norm() == pytest.approx(0.0) assert (b1 * 3.0 - b2 * 5.0).norm() == pytest.approx(0.0)
def J_mono(self, snes, x, J, P): J.zeroEntries() assemble_matrix(J, self.a, bcs=self.bcs, diagonal=1.0) J.assemble() if self.a_precon is not None: P.zeroEntries() assemble_matrix(P, self.a_precon, bcs=self.bcs, diagonal=1.0) P.assemble()
def test_assembly_dx_domains(mode, meshtags_factory): mesh = create_unit_square(MPI.COMM_WORLD, 10, 10, ghost_mode=mode) V = FunctionSpace(mesh, ("Lagrange", 1)) u, v = ufl.TrialFunction(V), ufl.TestFunction(V) # Prepare a marking structures # indices cover all cells # values are [1, 2, 3, 3, ...] cell_map = mesh.topology.index_map(mesh.topology.dim) num_cells = cell_map.size_local + cell_map.num_ghosts indices = np.arange(0, num_cells) values = np.full(indices.shape, 3, dtype=np.intc) values[0] = 1 values[1] = 2 marker = meshtags_factory(mesh, mesh.topology.dim, indices, values) dx = ufl.Measure('dx', subdomain_data=marker, domain=mesh) w = Function(V) w.x.array[:] = 0.5 # Assemble matrix a = form(w * ufl.inner(u, v) * (dx(1) + dx(2) + dx(3))) A = assemble_matrix(a) A.assemble() a2 = form(w * ufl.inner(u, v) * dx) A2 = assemble_matrix(a2) A2.assemble() assert (A - A2).norm() < 1.0e-12 bc = dirichletbc(Function(V), range(30)) # Assemble vector L = form(ufl.inner(w, v) * (dx(1) + dx(2) + dx(3))) b = assemble_vector(L) apply_lifting(b, [a], [[bc]]) b.ghostUpdate(addv=PETSc.InsertMode.ADD_VALUES, mode=PETSc.ScatterMode.REVERSE) set_bc(b, [bc]) L2 = form(ufl.inner(w, v) * dx) b2 = assemble_vector(L2) apply_lifting(b2, [a], [[bc]]) b2.ghostUpdate(addv=PETSc.InsertMode.ADD_VALUES, mode=PETSc.ScatterMode.REVERSE) set_bc(b2, [bc]) assert (b - b2).norm() < 1.0e-12 # Assemble scalar L = form(w * (dx(1) + dx(2) + dx(3))) s = assemble_scalar(L) s = mesh.comm.allreduce(s, op=MPI.SUM) assert s == pytest.approx(0.5, 1.0e-12) L2 = form(w * dx) s2 = assemble_scalar(L2) s2 = mesh.comm.allreduce(s2, op=MPI.SUM) assert s == pytest.approx(s2, 1.0e-12)
def test_refine_create_form(): """Check that forms can be assembled on refined mesh""" mesh = create_unit_cube(MPI.COMM_WORLD, 3, 3, 3) mesh.topology.create_entities(1) mesh = refine(mesh, redistribute=True) V = FunctionSpace(mesh, ("Lagrange", 1)) # Define variational problem u = ufl.TrialFunction(V) v = ufl.TestFunction(V) a = form(ufl.inner(ufl.grad(u), ufl.grad(v)) * ufl.dx) assemble_matrix(a)
def test_complex_assembly(): """Test assembly of complex matrices and vectors""" mesh = create_unit_square(MPI.COMM_WORLD, 10, 10) P2 = ufl.FiniteElement("Lagrange", mesh.ufl_cell(), 2) V = FunctionSpace(mesh, P2) u = ufl.TrialFunction(V) v = ufl.TestFunction(V) g = -2 + 3.0j j = 1.0j a_real = form(inner(u, v) * dx) L1 = form(inner(g, v) * dx) b = assemble_vector(L1) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) bnorm = b.norm(PETSc.NormType.N1) b_norm_ref = abs(-2 + 3.0j) assert bnorm == pytest.approx(b_norm_ref) A = assemble_matrix(a_real) A.assemble() A0_norm = A.norm(PETSc.NormType.FROBENIUS) x = ufl.SpatialCoordinate(mesh) a_imag = form(j * inner(u, v) * dx) f = 1j * ufl.sin(2 * np.pi * x[0]) L0 = form(inner(f, v) * dx) A = assemble_matrix(a_imag) A.assemble() A1_norm = A.norm(PETSc.NormType.FROBENIUS) assert A0_norm == pytest.approx(A1_norm) b = assemble_vector(L0) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) b1_norm = b.norm(PETSc.NormType.N2) a_complex = form((1 + j) * inner(u, v) * dx) f = ufl.sin(2 * np.pi * x[0]) L2 = form(inner(f, v) * dx) A = assemble_matrix(a_complex) A.assemble() A2_norm = A.norm(PETSc.NormType.FROBENIUS) assert A1_norm == pytest.approx(A2_norm / np.sqrt(2)) b = assemble_vector(L2) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) b2_norm = b.norm(PETSc.NormType.N2) assert b2_norm == pytest.approx(b1_norm)
def test_assemble_manifold(): """Test assembly of poisson problem on a mesh with topological dimension 1 but embedded in 2D (gdim=2)""" points = np.array([[0.0, 0.0], [0.2, 0.0], [0.4, 0.0], [0.6, 0.0], [0.8, 0.0], [1.0, 0.0]], dtype=np.float64) cells = np.array([[0, 1], [1, 2], [2, 3], [3, 4], [4, 5]], dtype=np.int32) cell = ufl.Cell("interval", geometric_dimension=points.shape[1]) domain = ufl.Mesh(ufl.VectorElement("Lagrange", cell, 1)) mesh = create_mesh(MPI.COMM_WORLD, cells, points, domain) assert mesh.geometry.dim == 2 assert mesh.topology.dim == 1 U = FunctionSpace(mesh, ("P", 1)) u, v = ufl.TrialFunction(U), ufl.TestFunction(U) a = ufl.inner(ufl.grad(u), ufl.grad(v)) * ufl.dx(mesh) L = ufl.inner(1.0, v) * ufl.dx(mesh) a, L = form(a), form(L) bcdofs = locate_dofs_geometrical(U, lambda x: np.isclose(x[0], 0.0)) bcs = [dirichletbc(PETSc.ScalarType(0), bcdofs, U)] A = assemble_matrix(a, bcs=bcs) A.assemble() b = assemble_vector(L) apply_lifting(b, [a], bcs=[bcs]) set_bc(b, bcs) assert np.isclose(b.norm(), 0.41231) assert np.isclose(A.norm(), 25.0199)
def test_assembly_bcs(mode): mesh = create_unit_square(MPI.COMM_WORLD, 12, 12, ghost_mode=mode) V = FunctionSpace(mesh, ("Lagrange", 1)) u, v = ufl.TrialFunction(V), ufl.TestFunction(V) a = inner(u, v) * dx + inner(u, v) * ds L = inner(1.0, v) * dx a, L = form(a), form(L) bdofsV = locate_dofs_geometrical( V, lambda x: np.logical_or(np.isclose(x[0], 0.0), np.isclose(x[0], 1.0))) bc = dirichletbc(PETSc.ScalarType(1), bdofsV, V) # Assemble and apply 'global' lifting of bcs A = assemble_matrix(a) A.assemble() b = assemble_vector(L) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) g = b.duplicate() with g.localForm() as g_local: g_local.set(0.0) set_bc(g, [bc]) f = b - A * g set_bc(f, [bc]) # Assemble vector and apply lifting of bcs during assembly b_bc = assemble_vector(L) apply_lifting(b_bc, [a], [[bc]]) b_bc.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) set_bc(b_bc, [bc]) assert (f - b_bc).norm() == pytest.approx(0.0, rel=1e-12, abs=1e-12)
def test_vector_assemble_matrix_interior(): mesh = create_unit_square(MPI.COMM_WORLD, 3, 3) V = VectorFunctionSpace(mesh, ("Lagrange", 1)) u, v = ufl.TrialFunction(V), ufl.TestFunction(V) a = form(ufl.inner(ufl.jump(u), ufl.jump(v)) * ufl.dS) A = assemble_matrix(a) A.assemble()
def test_nullspace_check(mesh, degree): V = VectorFunctionSpace(mesh, ('Lagrange', degree)) u, v = TrialFunction(V), TestFunction(V) E, nu = 2.0e2, 0.3 mu = E / (2.0 * (1.0 + nu)) lmbda = E * nu / ((1.0 + nu) * (1.0 - 2.0 * nu)) def sigma(w, gdim): return 2.0 * mu * ufl.sym(grad(w)) + lmbda * ufl.tr( grad(w)) * ufl.Identity(gdim) a = form(inner(sigma(u, mesh.geometry.dim), grad(v)) * dx) # Assemble matrix and create compatible vector A = assemble_matrix(a) A.assemble() # Create null space basis and test nullspace = build_elastic_nullspace(V) la.orthonormalize(nullspace) ns = PETSc.NullSpace().create(vectors=nullspace) assert ns.test(A) # Create incorrect null space basis and test nullspace = build_broken_elastic_nullspace(V) la.orthonormalize(nullspace) ns = PETSc.NullSpace().create(vectors=nullspace) assert not ns.test(A)
def test_krylov_solver_lu(): mesh = create_unit_square(MPI.COMM_WORLD, 12, 12) V = FunctionSpace(mesh, ("Lagrange", 1)) u, v = TrialFunction(V), TestFunction(V) a = form(inner(u, v) * dx) L = form(inner(1.0, v) * dx) A = assemble_matrix(a) A.assemble() b = assemble_vector(L) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) norm = 13.0 solver = PETSc.KSP().create(mesh.comm) solver.setOptionsPrefix("test_lu_") opts = PETSc.Options("test_lu_") opts["ksp_type"] = "preonly" opts["pc_type"] = "lu" solver.setFromOptions() x = A.createVecRight() solver.setOperators(A) solver.solve(b, x) # *Tight* tolerance for LU solves assert x.norm(PETSc.NormType.N2) == pytest.approx(norm, abs=1.0e-12)
def amg_solve(N, method): # Elasticity parameters E = 1.0e9 nu = 0.3 mu = E / (2.0 * (1.0 + nu)) lmbda = E * nu / ((1.0 + nu) * (1.0 - 2.0 * nu)) # Stress computation def sigma(v): return 2.0 * mu * sym(grad(v)) + lmbda * tr(sym( grad(v))) * Identity(2) # Define problem mesh = create_unit_square(MPI.COMM_WORLD, N, N) V = VectorFunctionSpace(mesh, 'Lagrange', 1) u = TrialFunction(V) v = TestFunction(V) facetdim = mesh.topology.dim - 1 bndry_facets = locate_entities_boundary( mesh, facetdim, lambda x: np.full(x.shape[1], True)) bdofs = locate_dofs_topological(V.sub(0), V, facetdim, bndry_facets) bc = dirichletbc(PETSc.ScalarType(0), bdofs, V.sub(0)) # Forms a, L = inner(sigma(u), grad(v)) * dx, dot(ufl.as_vector( (1.0, 1.0)), v) * dx # Assemble linear algebra objects A = assemble_matrix(a, [bc]) A.assemble() b = assemble_vector(L) apply_lifting(b, [a], [[bc]]) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) set_bc(b, [bc]) # Create solution function u = Function(V) # Create near null space basis and orthonormalize null_space = build_nullspace(V, u.vector) # Attached near-null space to matrix A.set_near_nullspace(null_space) # Test that basis is orthonormal assert null_space.is_orthonormal() # Create PETSC smoothed aggregation AMG preconditioner, and # create CG solver solver = PETSc.KSP().create(mesh.comm) solver.setType("cg") # Set matrix operator solver.setOperators(A) # Compute solution and return number of iterations return solver.solve(b, u.vector)
def run_scalar_test(mesh, V, degree): """ Manufactured Poisson problem, solving u = x[1]**p, where p is the degree of the Lagrange function space. """ u, v = TrialFunction(V), TestFunction(V) a = inner(grad(u), grad(v)) * dx # Get quadrature degree for bilinear form integrand (ignores effect of non-affine map) a = inner(grad(u), grad(v)) * dx(metadata={"quadrature_degree": -1}) a.integrals()[0].metadata( )["quadrature_degree"] = ufl.algorithms.estimate_total_polynomial_degree(a) a = form(a) # Source term x = SpatialCoordinate(mesh) u_exact = x[1]**degree f = -div(grad(u_exact)) # Set quadrature degree for linear form integrand (ignores effect of non-affine map) L = inner(f, v) * dx(metadata={"quadrature_degree": -1}) L.integrals()[0].metadata( )["quadrature_degree"] = ufl.algorithms.estimate_total_polynomial_degree(L) L = form(L) u_bc = Function(V) u_bc.interpolate(lambda x: x[1]**degree) # Create Dirichlet boundary condition facetdim = mesh.topology.dim - 1 mesh.topology.create_connectivity(facetdim, mesh.topology.dim) bndry_facets = np.where( np.array(compute_boundary_facets(mesh.topology)) == 1)[0] bdofs = locate_dofs_topological(V, facetdim, bndry_facets) bc = dirichletbc(u_bc, bdofs) b = assemble_vector(L) apply_lifting(b, [a], bcs=[[bc]]) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) set_bc(b, [bc]) a = form(a) A = assemble_matrix(a, bcs=[bc]) A.assemble() # Create LU linear solver solver = PETSc.KSP().create(MPI.COMM_WORLD) solver.setType(PETSc.KSP.Type.PREONLY) solver.getPC().setType(PETSc.PC.Type.LU) solver.setOperators(A) uh = Function(V) solver.solve(b, uh.vector) uh.x.scatter_forward() M = (u_exact - uh)**2 * dx M = form(M) error = mesh.comm.allreduce(assemble_scalar(M), op=MPI.SUM) assert np.absolute(error) < 1.0e-14
def test_basic_assembly(mode): mesh = create_unit_square(MPI.COMM_WORLD, 12, 12, ghost_mode=mode) V = FunctionSpace(mesh, ("Lagrange", 1)) u, v = ufl.TrialFunction(V), ufl.TestFunction(V) f = Function(V) f.x.array[:] = 10.0 a = inner(f * u, v) * dx + inner(u, v) * ds L = inner(f, v) * dx + inner(2.0, v) * ds a, L = form(a), form(L) # Initial assembly A = assemble_matrix(a) A.assemble() assert isinstance(A, PETSc.Mat) b = assemble_vector(L) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) assert isinstance(b, PETSc.Vec) # Second assembly normA = A.norm() A.zeroEntries() A = assemble_matrix(A, a) A.assemble() assert isinstance(A, PETSc.Mat) assert normA == pytest.approx(A.norm()) normb = b.norm() with b.localForm() as b_local: b_local.set(0.0) b = assemble_vector(b, L) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) assert isinstance(b, PETSc.Vec) assert normb == pytest.approx(b.norm()) # Vector re-assembly - no zeroing (but need to zero ghost entries) with b.localForm() as b_local: b_local.array[b.local_size:] = 0.0 assemble_vector(b, L) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) assert 2.0 * normb == pytest.approx(b.norm()) # Matrix re-assembly (no zeroing) assemble_matrix(A, a) A.assemble() assert 2.0 * normA == pytest.approx(A.norm())
def assemble_div_matrix(k, offset): mesh = create_quad_mesh(offset) V = FunctionSpace(mesh, ("DQ", k)) W = FunctionSpace(mesh, ("RTCF", k + 1)) u, w = ufl.TrialFunction(V), ufl.TestFunction(W) a = form(ufl.inner(u, ufl.div(w)) * ufl.dx) A = assemble_matrix(a) A.assemble() return A[:, :]
def test_overlapping_bcs(): """Test that, when boundaries condition overlap, the last provided boundary condition is applied""" n = 23 mesh = create_unit_square(MPI.COMM_WORLD, n, n) V = FunctionSpace(mesh, ("Lagrange", 1)) u, v = ufl.TrialFunction(V), ufl.TestFunction(V) a = form(inner(u, v) * dx) L = form(inner(1, v) * dx) dofs_left = locate_dofs_geometrical(V, lambda x: x[0] < 1.0 / (2.0 * n)) dofs_top = locate_dofs_geometrical(V, lambda x: x[1] > 1.0 - 1.0 / (2.0 * n)) dof_corner = np.array(list(set(dofs_left).intersection(set(dofs_top))), dtype=np.int64) # Check only one dof pair is found globally assert len(set(np.concatenate(MPI.COMM_WORLD.allgather(dof_corner)))) == 1 bcs = [ dirichletbc(PETSc.ScalarType(0), dofs_left, V), dirichletbc(PETSc.ScalarType(123.456), dofs_top, V) ] A, b = create_matrix(a), create_vector(L) assemble_matrix(A, a, bcs=bcs) A.assemble() # Check the diagonal (only on the rank that owns the row) d = A.getDiagonal() if len(dof_corner) > 0 and dof_corner[0] < V.dofmap.index_map.size_local: assert np.isclose(d.array_r[dof_corner[0]], 1.0) with b.localForm() as b_loc: b_loc.set(0) assemble_vector(b, L) apply_lifting(b, [a], [bcs]) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) set_bc(b, bcs) b.ghostUpdate(addv=PETSc.InsertMode.INSERT, mode=PETSc.ScatterMode.FORWARD) if len(dof_corner) > 0: with b.localForm() as b_loc: assert b_loc[dof_corner[0]] == 123.456
def test_assemble_empty_rank_mesh(): """Assembly on mesh where some ranks are empty""" comm = MPI.COMM_WORLD cell_type = CellType.triangle domain = ufl.Mesh( ufl.VectorElement("Lagrange", ufl.Cell(cell_type.name), 1)) def partitioner(comm, nparts, local_graph, num_ghost_nodes, ghosting): """Leave cells on the curent rank""" dest = np.full(len(cells), comm.rank, dtype=np.int32) return graph.create_adjacencylist(dest) if comm.rank == 0: # Put cells on rank 0 cells = np.array([[0, 1, 2], [0, 2, 3]], dtype=np.int64) cells = graph.create_adjacencylist(cells) x = np.array([[0., 0.], [1., 0.], [1., 1.], [0., 1.]]) else: # No cells onm other ranks cells = graph.create_adjacencylist(np.empty((0, 3), dtype=np.int64)) x = np.empty((0, 2), dtype=np.float64) mesh = create_mesh(comm, cells, x, domain, GhostMode.none, partitioner) V = FunctionSpace(mesh, ("Lagrange", 2)) u, v = ufl.TrialFunction(V), ufl.TestFunction(V) f, k, zero = Function(V), Function(V), Function(V) f.x.array[:] = 10.0 k.x.array[:] = 1.0 zero.x.array[:] = 0.0 a = form(inner(k * u, v) * dx + inner(zero * u, v) * ds) L = form(inner(f, v) * dx + inner(zero, v) * ds) M = form(2 * k * dx + k * ds) sum = comm.allreduce(assemble_scalar(M), op=MPI.SUM) assert sum == pytest.approx(6.0) # Assemble A = assemble_matrix(a) A.assemble() b = assemble_vector(L) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) # Solve ksp = PETSc.KSP() ksp.create(mesh.comm) ksp.setOperators(A) ksp.setTolerances(rtol=1.0e-9, max_it=50) ksp.setFromOptions() x = b.copy() ksp.solve(b, x) assert np.allclose(x.array, 10.0)
def test_plus_minus_matrix(cell_type, pm1, pm2): """Test that ('+') and ('-') match up with the correct DOFs for DG functions""" results = [] spaces = [] orders = [] for count in range(3): for agree in [True, False]: # Two cell mesh with randomly numbered points mesh, order = two_unit_cells(cell_type, agree, return_order=True) V = FunctionSpace(mesh, ("DG", 1)) u, v = ufl.TrialFunction(V), ufl.TestFunction(V) # Assemble matrices with combinations of + and - for a few # different numberings a = form(ufl.inner(u(pm1), v(pm2)) * ufl.dS) result = assemble_matrix(a, []) result.assemble() spaces.append(V) results.append(result) orders.append(order) # Check that the above matrices all have the same values, but # permuted due to differently ordered dofs dofmap0 = spaces[0].mesh.geometry.dofmap for result, space in zip(results[1:], spaces[1:]): # Get the data relating to two results dofmap1 = space.mesh.geometry.dofmap dof_order = [] # For each cell for cell in range(2): # For each point in cell 0 in the first mesh for dof0, point0 in zip(spaces[0].dofmap.cell_dofs(cell), dofmap0.links(cell)): # Find the point in the cell 0 in the second mesh for dof1, point1 in zip(space.dofmap.cell_dofs(cell), dofmap1.links(cell)): if np.allclose(spaces[0].mesh.geometry.x[point0], space.mesh.geometry.x[point1]): break else: # If no matching point found, fail assert False dof_order.append((dof0, dof1)) # For all dof pairs, check that entries in the matrix agree for a, b in dof_order: for c, d in dof_order: assert np.isclose(results[0][a, c], result[b, d])
def test_custom_mesh_loop_cffi_rank2(set_vals): """Test numba assembler for bilinear form""" mesh = create_unit_square(MPI.COMM_WORLD, 64, 64) V = FunctionSpace(mesh, ("Lagrange", 1)) # Test against generated code and general assembler u, v = ufl.TrialFunction(V), ufl.TestFunction(V) a = form(inner(u, v) * dx) A0 = assemble_matrix(a) A0.assemble() A0.zeroEntries() start = time.time() assemble_matrix(A0, a) end = time.time() print("Time (C++, pass 2):", end - start) A0.assemble() # Unpack mesh and dofmap data num_owned_cells = mesh.topology.index_map(mesh.topology.dim).size_local num_cells = num_owned_cells + mesh.topology.index_map( mesh.topology.dim).num_ghosts x_dofs = mesh.geometry.dofmap.array.reshape(num_cells, 3) x = mesh.geometry.x dofmap = V.dofmap.list.array.reshape(num_cells, 3).astype(np.dtype(PETSc.IntType)) A1 = A0.copy() for i in range(2): A1.zeroEntries() start = time.time() assemble_matrix_cffi(A1.handle, (x_dofs, x), dofmap, num_owned_cells, set_vals, PETSc.InsertMode.ADD_VALUES) end = time.time() print("Time (Numba, pass {}): {}".format(i, end - start)) A1.assemble() assert (A1 - A0).norm() == pytest.approx(0.0)
def test_complex_assembly_solve(): """Solve a positive definite helmholtz problem and verify solution with the method of manufactured solutions""" degree = 3 mesh = create_unit_square(MPI.COMM_WORLD, 20, 20) P = ufl.FiniteElement("Lagrange", mesh.ufl_cell(), degree) V = FunctionSpace(mesh, P) x = ufl.SpatialCoordinate(mesh) # Define source term A = 1.0 + 2.0 * (2.0 * np.pi)**2 f = (1. + 1j) * A * ufl.cos(2 * np.pi * x[0]) * ufl.cos(2 * np.pi * x[1]) # Variational problem u, v = ufl.TrialFunction(V), ufl.TestFunction(V) C = 1.0 + 1.0j a = form(C * inner(grad(u), grad(v)) * dx + C * inner(u, v) * dx) L = form(inner(f, v) * dx) # Assemble A = assemble_matrix(a) A.assemble() b = assemble_vector(L) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) # Create solver solver = PETSc.KSP().create(mesh.comm) solver.setOptionsPrefix("test_lu_") opts = PETSc.Options("test_lu_") opts["ksp_type"] = "preonly" opts["pc_type"] = "lu" solver.setFromOptions() x = A.createVecRight() solver.setOperators(A) solver.solve(b, x) # Reference Solution def ref_eval(x): return np.cos(2 * np.pi * x[0]) * np.cos(2 * np.pi * x[1]) u_ref = Function(V) u_ref.interpolate(ref_eval) diff = (x - u_ref.vector).norm(PETSc.NormType.N2) assert diff == pytest.approx(0.0, abs=1e-1)
def test_custom_mesh_loop_ctypes_rank2(): """Test numba assembler for bilinear form""" # Create mesh and function space mesh = create_unit_square(MPI.COMM_WORLD, 64, 64) V = FunctionSpace(mesh, ("Lagrange", 1)) # Extract mesh and dofmap data num_owned_cells = mesh.topology.index_map(mesh.topology.dim).size_local num_cells = num_owned_cells + mesh.topology.index_map( mesh.topology.dim).num_ghosts x_dofs = mesh.geometry.dofmap.array.reshape(num_cells, 3) x = mesh.geometry.x dofmap = V.dofmap.list.array.reshape(num_cells, 3).astype(np.dtype(PETSc.IntType)) # Generated case with general assembler u, v = ufl.TrialFunction(V), ufl.TestFunction(V) a = form(inner(u, v) * dx) A0 = assemble_matrix(a) A0.assemble() A0.zeroEntries() start = time.time() dolfinx.fem.petsc.assemble_matrix(A0, a) end = time.time() print("Time (C++, pass 2):", end - start) A0.assemble() # Custom case A1 = A0.copy() for i in range(2): A1.zeroEntries() mat = A1.handle start = time.time() assemble_matrix_ctypes(mat, (x_dofs, x), dofmap, num_owned_cells, MatSetValues_ctypes, PETSc.InsertMode.ADD_VALUES) end = time.time() print("Time (numba, pass {}): {}".format(i, end - start)) A1.assemble() assert (A0 - A1).norm() == pytest.approx(0.0, abs=1.0e-9)
def test_basic_interior_facet_assembly(): mesh = create_rectangle( MPI.COMM_WORLD, [np.array([0.0, 0.0]), np.array([1.0, 1.0])], [5, 5], cell_type=CellType.triangle, ghost_mode=GhostMode.shared_facet) V = FunctionSpace(mesh, ("DG", 1)) u, v = ufl.TrialFunction(V), ufl.TestFunction(V) a = ufl.inner(ufl.avg(u), ufl.avg(v)) * ufl.dS a = form(a) A = assemble_matrix(a) A.assemble() assert isinstance(A, PETSc.Mat) L = ufl.conj(ufl.avg(v)) * ufl.dS L = form(L) b = assemble_vector(L) b.assemble() assert isinstance(b, PETSc.Vec)
def run_vector_test(mesh, V, degree): """Projection into H(div/curl) spaces""" u, v = ufl.TrialFunction(V), ufl.TestFunction(V) a = form(inner(u, v) * dx) # Source term x = SpatialCoordinate(mesh) u_exact = x[0]**degree L = form(inner(u_exact, v[0]) * dx) b = assemble_vector(L) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) A = assemble_matrix(a) A.assemble() # Create LU linear solver (Note: need to use a solver that # re-orders to handle pivots, e.g. not the PETSc built-in LU solver) solver = PETSc.KSP().create(MPI.COMM_WORLD) solver.setType("preonly") solver.getPC().setType('lu') solver.setOperators(A) # Solve uh = Function(V) solver.solve(b, uh.vector) uh.x.scatter_forward() # Calculate error M = (u_exact - uh[0])**2 * dx for i in range(1, mesh.topology.dim): M += uh[i]**2 * dx M = form(M) error = mesh.comm.allreduce(assemble_scalar(M), op=MPI.SUM) assert np.absolute(error) < 1.0e-14
def test_matrix_assembly_block_nl(): """Test assembly of block matrices and vectors into (a) monolithic blocked structures, PETSc Nest structures, and monolithic structures in the nonlinear setting """ mesh = create_unit_square(MPI.COMM_WORLD, 4, 8) p0, p1 = 1, 2 P0 = ufl.FiniteElement("Lagrange", mesh.ufl_cell(), p0) P1 = ufl.FiniteElement("Lagrange", mesh.ufl_cell(), p1) V0 = FunctionSpace(mesh, P0) V1 = FunctionSpace(mesh, P1) def initial_guess_u(x): return np.sin(x[0]) * np.sin(x[1]) def initial_guess_p(x): return -x[0]**2 - x[1]**3 def bc_value(x): return np.cos(x[0]) * np.cos(x[1]) facetdim = mesh.topology.dim - 1 bndry_facets = locate_entities_boundary( mesh, facetdim, lambda x: np.logical_or(np.isclose(x[0], 0.0), np.isclose(x[0], 1.0))) u_bc = Function(V1) u_bc.interpolate(bc_value) bdofs = locate_dofs_topological(V1, facetdim, bndry_facets) bc = dirichletbc(u_bc, bdofs) # Define variational problem du, dp = ufl.TrialFunction(V0), ufl.TrialFunction(V1) u, p = Function(V0), Function(V1) v, q = ufl.TestFunction(V0), ufl.TestFunction(V1) u.interpolate(initial_guess_u) p.interpolate(initial_guess_p) f = 1.0 g = -3.0 F0 = inner(u, v) * dx + inner(p, v) * dx - inner(f, v) * dx F1 = inner(u, q) * dx + inner(p, q) * dx - inner(g, q) * dx a_block = form([[derivative(F0, u, du), derivative(F0, p, dp)], [derivative(F1, u, du), derivative(F1, p, dp)]]) L_block = form([F0, F1]) # Monolithic blocked x0 = create_vector_block(L_block) scatter_local_vectors(x0, [u.vector.array_r, p.vector.array_r], [(u.function_space.dofmap.index_map, u.function_space.dofmap.index_map_bs), (p.function_space.dofmap.index_map, p.function_space.dofmap.index_map_bs)]) x0.ghostUpdate(addv=PETSc.InsertMode.INSERT, mode=PETSc.ScatterMode.FORWARD) # Ghosts are updated inside assemble_vector_block A0 = assemble_matrix_block(a_block, bcs=[bc]) b0 = assemble_vector_block(L_block, a_block, bcs=[bc], x0=x0, scale=-1.0) A0.assemble() assert A0.getType() != "nest" Anorm0 = A0.norm() bnorm0 = b0.norm() # Nested (MatNest) x1 = create_vector_nest(L_block) for x1_soln_pair in zip(x1.getNestSubVecs(), (u, p)): x1_sub, soln_sub = x1_soln_pair soln_sub.vector.ghostUpdate(addv=PETSc.InsertMode.INSERT, mode=PETSc.ScatterMode.FORWARD) soln_sub.vector.copy(result=x1_sub) x1_sub.ghostUpdate(addv=PETSc.InsertMode.INSERT, mode=PETSc.ScatterMode.FORWARD) A1 = assemble_matrix_nest(a_block, bcs=[bc]) b1 = assemble_vector_nest(L_block) apply_lifting_nest(b1, a_block, bcs=[bc], x0=x1, scale=-1.0) for b_sub in b1.getNestSubVecs(): b_sub.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) bcs0 = bcs_by_block([L.function_spaces[0] for L in L_block], [bc]) set_bc_nest(b1, bcs0, x1, scale=-1.0) A1.assemble() assert A1.getType() == "nest" assert nest_matrix_norm(A1) == pytest.approx(Anorm0, 1.0e-12) assert b1.norm() == pytest.approx(bnorm0, 1.0e-12) # Monolithic version E = P0 * P1 W = FunctionSpace(mesh, E) dU = ufl.TrialFunction(W) U = Function(W) u0, u1 = ufl.split(U) v0, v1 = ufl.TestFunctions(W) U.sub(0).interpolate(initial_guess_u) U.sub(1).interpolate(initial_guess_p) F = inner(u0, v0) * dx + inner(u1, v0) * dx + inner(u0, v1) * dx + inner(u1, v1) * dx \ - inner(f, v0) * ufl.dx - inner(g, v1) * dx J = derivative(F, U, dU) F, J = form(F), form(J) bdofsW_V1 = locate_dofs_topological((W.sub(1), V1), facetdim, bndry_facets) bc = dirichletbc(u_bc, bdofsW_V1, W.sub(1)) A2 = assemble_matrix(J, bcs=[bc]) A2.assemble() b2 = assemble_vector(F) apply_lifting(b2, [J], bcs=[[bc]], x0=[U.vector], scale=-1.0) b2.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) set_bc(b2, [bc], x0=U.vector, scale=-1.0) assert A2.getType() != "nest" assert A2.norm() == pytest.approx(Anorm0, 1.0e-12) assert b2.norm() == pytest.approx(bnorm0, 1.0e-12)
def J(self, x, A): """Assemble Jacobian matrix.""" A.zeroEntries() assemble_matrix(A, self.a, bcs=[self.bc]) A.assemble()
def test_eigen_assembly(tempdir): # noqa: F811 """Compare assembly into scipy.CSR matrix with PETSc assembly""" def compile_eigen_csr_assembler_module(): dolfinx_pc = dolfinx.pkgconfig.parse("dolfinx") eigen_dir = dolfinx.pkgconfig.parse("eigen3")["include_dirs"] cpp_code_header = f""" <% setup_pybind11(cfg) cfg['include_dirs'] = {dolfinx_pc["include_dirs"] + [petsc4py.get_include()] + [pybind11.get_include()] + [str(pybind_inc())] + eigen_dir} cfg['compiler_args'] = ["-std=c++17", "-Wno-comment"] cfg['libraries'] = {dolfinx_pc["libraries"]} cfg['library_dirs'] = {dolfinx_pc["library_dirs"]} %> """ cpp_code = """ #include <pybind11/pybind11.h> #include <pybind11/eigen.h> #include <pybind11/stl.h> #include <vector> #include <Eigen/Sparse> #include <petscsys.h> #include <dolfinx/fem/assembler.h> #include <dolfinx/fem/DirichletBC.h> #include <dolfinx/fem/Form.h> template<typename T> Eigen::SparseMatrix<T, Eigen::RowMajor> assemble_csr(const dolfinx::fem::Form<T>& a, const std::vector<std::shared_ptr<const dolfinx::fem::DirichletBC<T>>>& bcs) { std::vector<Eigen::Triplet<T>> triplets; auto mat_add = [&triplets](const xtl::span<const std::int32_t>& rows, const xtl::span<const std::int32_t>& cols, const xtl::span<const T>& v) { for (std::size_t i = 0; i < rows.size(); ++i) for (std::size_t j = 0; j < cols.size(); ++j) triplets.emplace_back(rows[i], cols[j], v[i * cols.size() + j]); return 0; }; dolfinx::fem::assemble_matrix(mat_add, a, bcs); auto map0 = a.function_spaces().at(0)->dofmap()->index_map; int bs0 = a.function_spaces().at(0)->dofmap()->index_map_bs(); auto map1 = a.function_spaces().at(1)->dofmap()->index_map; int bs1 = a.function_spaces().at(1)->dofmap()->index_map_bs(); Eigen::SparseMatrix<T, Eigen::RowMajor> mat( bs0 * (map0->size_local() + map0->num_ghosts()), bs1 * (map1->size_local() + map1->num_ghosts())); mat.setFromTriplets(triplets.begin(), triplets.end()); return mat; } PYBIND11_MODULE(eigen_csr, m) { m.def("assemble_matrix", &assemble_csr<PetscScalar>); } """ path = pathlib.Path(tempdir) open(pathlib.Path(tempdir, "eigen_csr.cpp"), "w").write(cpp_code + cpp_code_header) rel_path = path.relative_to(pathlib.Path(__file__).parent) p = str(rel_path).replace("/", ".") + ".eigen_csr" return cppimport.imp(p) def assemble_csr_matrix(a, bcs): """Assemble bilinear form into an SciPy CSR matrix, in serial.""" module = compile_eigen_csr_assembler_module() A = module.assemble_matrix(a, bcs) if a.function_spaces[0] is a.function_spaces[1]: for bc in bcs: if a.function_spaces[0].contains(bc.function_space): bc_dofs, _ = bc.dof_indices() # See https://github.com/numpy/numpy/issues/14132 # for why we copy bc_dofs as a work-around dofs = bc_dofs.copy() A[dofs, dofs] = 1.0 return A mesh = create_unit_square(MPI.COMM_SELF, 12, 12) Q = FunctionSpace(mesh, ("Lagrange", 1)) u = ufl.TrialFunction(Q) v = ufl.TestFunction(Q) a = form(ufl.inner(ufl.grad(u), ufl.grad(v)) * ufl.dx) bdofsQ = locate_dofs_geometrical( Q, lambda x: np.logical_or(np.isclose(x[0], 0.0), np.isclose(x[0], 1.0))) bc = dirichletbc(PETSc.ScalarType(1), bdofsQ, Q) A1 = assemble_matrix(a, [bc]) A1.assemble() A2 = assemble_csr_matrix(a, [bc]) assert np.isclose(A1.norm(), scipy.sparse.linalg.norm(A2))
def J(self, snes, x, J, P): """Assemble Jacobian matrix.""" J.zeroEntries() assemble_matrix(J, self.a, bcs=[self.bc]) J.assemble()
def run_dg_test(mesh, V, degree): """Manufactured Poisson problem, solving u = x[component]**n, where n is the degree of the Lagrange function space. """ u, v = TrialFunction(V), TestFunction(V) # Exact solution x = SpatialCoordinate(mesh) u_exact = x[1]**degree # Coefficient k = Function(V) k.x.array[:] = 2.0 # Source term f = -div(k * grad(u_exact)) # Mesh normals and element size n = FacetNormal(mesh) h = CellDiameter(mesh) h_avg = (h("+") + h("-")) / 2.0 # Penalty parameter alpha = 32 dx_ = dx(metadata={"quadrature_degree": -1}) ds_ = ds(metadata={"quadrature_degree": -1}) dS_ = dS(metadata={"quadrature_degree": -1}) a = inner(k * grad(u), grad(v)) * dx_ \ - k("+") * inner(avg(grad(u)), jump(v, n)) * dS_ \ - k("+") * inner(jump(u, n), avg(grad(v))) * dS_ \ + k("+") * (alpha / h_avg) * inner(jump(u, n), jump(v, n)) * dS_ \ - inner(k * grad(u), v * n) * ds_ \ - inner(u * n, k * grad(v)) * ds_ \ + (alpha / h) * inner(k * u, v) * ds_ L = inner(f, v) * dx_ - inner(k * u_exact * n, grad(v)) * ds_ \ + (alpha / h) * inner(k * u_exact, v) * ds_ for integral in a.integrals(): integral.metadata( )["quadrature_degree"] = ufl.algorithms.estimate_total_polynomial_degree( a) for integral in L.integrals(): integral.metadata( )["quadrature_degree"] = ufl.algorithms.estimate_total_polynomial_degree( L) a, L = form(a), form(L) b = assemble_vector(L) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) A = assemble_matrix(a, []) A.assemble() # Create LU linear solver solver = PETSc.KSP().create(MPI.COMM_WORLD) solver.setType(PETSc.KSP.Type.PREONLY) solver.getPC().setType(PETSc.PC.Type.LU) solver.setOperators(A) # Solve uh = Function(V) solver.solve(b, uh.vector) uh.x.scatter_forward() # Calculate error M = (u_exact - uh)**2 * dx M = form(M) error = mesh.comm.allreduce(assemble_scalar(M), op=MPI.SUM) assert np.absolute(error) < 1.0e-14
def test_biharmonic(): """Manufactured biharmonic problem. Solved using rotated Regge mixed finite element method. This is equivalent to the Hellan-Herrmann-Johnson (HHJ) finite element method in two-dimensions.""" mesh = create_rectangle( MPI.COMM_WORLD, [np.array([0.0, 0.0]), np.array([1.0, 1.0])], [32, 32], CellType.triangle) element = ufl.MixedElement([ ufl.FiniteElement("Regge", ufl.triangle, 1), ufl.FiniteElement("Lagrange", ufl.triangle, 2) ]) V = FunctionSpace(mesh, element) sigma, u = ufl.TrialFunctions(V) tau, v = ufl.TestFunctions(V) x = ufl.SpatialCoordinate(mesh) u_exact = ufl.sin(ufl.pi * x[0]) * ufl.sin(ufl.pi * x[0]) * ufl.sin( ufl.pi * x[1]) * ufl.sin(ufl.pi * x[1]) f_exact = div(grad(div(grad(u_exact)))) sigma_exact = grad(grad(u_exact)) # sigma and tau are tangential-tangential continuous according to the # H(curl curl) continuity of the Regge space. However, for the biharmonic # problem we require normal-normal continuity H (div div). Theorem 4.2 of # Lizao Li's PhD thesis shows that the latter space can be constructed by # the former through the action of the operator S: def S(tau): return tau - ufl.Identity(2) * ufl.tr(tau) sigma_S = S(sigma) tau_S = S(tau) # Discrete duality inner product eq. 4.5 Lizao Li's PhD thesis def b(tau_S, v): n = FacetNormal(mesh) return inner(tau_S, grad(grad(v))) * dx \ - ufl.dot(ufl.dot(tau_S('+'), n('+')), n('+')) * jump(grad(v), n) * dS \ - ufl.dot(ufl.dot(tau_S, n), n) * ufl.dot(grad(v), n) * ds # Non-symmetric formulation a = form(inner(sigma_S, tau_S) * dx - b(tau_S, u) + b(sigma_S, v)) L = form(inner(f_exact, v) * dx) V_1 = V.sub(1).collapse()[0] zero_u = Function(V_1) zero_u.x.array[:] = 0.0 # Strong (Dirichlet) boundary condition boundary_facets = locate_entities_boundary( mesh, mesh.topology.dim - 1, lambda x: np.full(x.shape[1], True, dtype=bool)) boundary_dofs = locate_dofs_topological( (V.sub(1), V_1), mesh.topology.dim - 1, boundary_facets) bcs = [dirichletbc(zero_u, boundary_dofs, V.sub(1))] A = assemble_matrix(a, bcs=bcs) A.assemble() b = assemble_vector(L) apply_lifting(b, [a], bcs=[bcs]) b.ghostUpdate(addv=PETSc.InsertMode.ADD, mode=PETSc.ScatterMode.REVERSE) set_bc(b, bcs) # Solve solver = PETSc.KSP().create(MPI.COMM_WORLD) PETSc.Options()["ksp_type"] = "preonly" PETSc.Options()["pc_type"] = "lu" # PETSc.Options()["pc_factor_mat_solver_type"] = "mumps" solver.setFromOptions() solver.setOperators(A) x_h = Function(V) solver.solve(b, x_h.vector) x_h.x.scatter_forward() # Recall that x_h has flattened indices. u_error_numerator = np.sqrt( mesh.comm.allreduce(assemble_scalar( form( inner(u_exact - x_h[4], u_exact - x_h[4]) * dx(mesh, metadata={"quadrature_degree": 5}))), op=MPI.SUM)) u_error_denominator = np.sqrt( mesh.comm.allreduce(assemble_scalar( form( inner(u_exact, u_exact) * dx(mesh, metadata={"quadrature_degree": 5}))), op=MPI.SUM)) assert np.absolute(u_error_numerator / u_error_denominator) < 0.05 # Reconstruct tensor from flattened indices. # Apply inverse transform. In 2D we have S^{-1} = S. sigma_h = S(ufl.as_tensor([[x_h[0], x_h[1]], [x_h[2], x_h[3]]])) sigma_error_numerator = np.sqrt( mesh.comm.allreduce(assemble_scalar( form( inner(sigma_exact - sigma_h, sigma_exact - sigma_h) * dx(mesh, metadata={"quadrature_degree": 5}))), op=MPI.SUM)) sigma_error_denominator = np.sqrt( mesh.comm.allreduce(assemble_scalar( form( inner(sigma_exact, sigma_exact) * dx(mesh, metadata={"quadrature_degree": 5}))), op=MPI.SUM)) assert np.absolute(sigma_error_numerator / sigma_error_denominator) < 0.005
def test_curl_curl_eigenvalue(family, order): """curl curl eigenvalue problem. Solved using H(curl)-conforming finite element method. See https://www-users.cse.umn.edu/~arnold/papers/icm2002.pdf for details. """ slepc4py = pytest.importorskip("slepc4py") # noqa: F841 from slepc4py import SLEPc mesh = create_rectangle( MPI.COMM_WORLD, [np.array([0.0, 0.0]), np.array([np.pi, np.pi])], [24, 24], CellType.triangle) element = ufl.FiniteElement(family, ufl.triangle, order) V = FunctionSpace(mesh, element) u = ufl.TrialFunction(V) v = ufl.TestFunction(V) a = inner(ufl.curl(u), ufl.curl(v)) * dx b = inner(u, v) * dx boundary_facets = locate_entities_boundary( mesh, mesh.topology.dim - 1, lambda x: np.full(x.shape[1], True, dtype=bool)) boundary_dofs = locate_dofs_topological(V, mesh.topology.dim - 1, boundary_facets) zero_u = Function(V) zero_u.x.array[:] = 0.0 bcs = [dirichletbc(zero_u, boundary_dofs)] a, b = form(a), form(b) A = assemble_matrix(a, bcs=bcs) A.assemble() B = assemble_matrix(b, bcs=bcs, diagonal=0.01) B.assemble() eps = SLEPc.EPS().create() eps.setOperators(A, B) PETSc.Options()["eps_type"] = "krylovschur" PETSc.Options()["eps_gen_hermitian"] = "" PETSc.Options()["eps_target_magnitude"] = "" PETSc.Options()["eps_target"] = 5.0 PETSc.Options()["eps_view"] = "" PETSc.Options()["eps_nev"] = 12 eps.setFromOptions() eps.solve() num_converged = eps.getConverged() eigenvalues_unsorted = np.zeros(num_converged, dtype=np.complex128) for i in range(0, num_converged): eigenvalues_unsorted[i] = eps.getEigenvalue(i) assert np.isclose(np.imag(eigenvalues_unsorted), 0.0).all() eigenvalues_sorted = np.sort(np.real(eigenvalues_unsorted))[:-1] eigenvalues_sorted = eigenvalues_sorted[np.logical_not( eigenvalues_sorted < 1E-8)] eigenvalues_exact = np.array([1.0, 1.0, 2.0, 4.0, 4.0, 5.0, 5.0, 8.0, 9.0]) assert np.isclose(eigenvalues_sorted[0:eigenvalues_exact.shape[0]], eigenvalues_exact, rtol=1E-2).all()