def test_RefineUnitCubeMesh_repartition(): """Refine mesh of unit cube.""" mesh = UnitCubeMesh(MPI.COMM_WORLD, 5, 7, 9, ghost_mode=GhostMode.none) mesh.topology.create_entities(1) mesh = refine(mesh, redistribute=True) assert mesh.topology.index_map(0).size_global == 3135 assert mesh.topology.index_map(3).size_global == 15120 mesh = UnitCubeMesh(MPI.COMM_WORLD, 5, 7, 9, ghost_mode=GhostMode.shared_facet) mesh.topology.create_entities(1) mesh = refine(mesh, redistribute=True) assert mesh.topology.index_map(0).size_global == 3135 assert mesh.topology.index_map(3).size_global == 15120 Q = FunctionSpace(mesh, ("CG", 1)) assert(Q)
def test_RefineUnitSquareMesh(): """Refine mesh of unit square.""" mesh = UnitSquareMesh(MPI.COMM_WORLD, 5, 7, ghost_mode=GhostMode.none) mesh.topology.create_entities(1) mesh = refine(mesh, redistribute=False) assert mesh.topology.index_map(0).size_global == 165 assert mesh.topology.index_map(2).size_global == 280
def test_refine_from_cells(): """Check user interface for using local cells to define edges""" Nx = 8 Ny = 3 assert Nx % 2 == 0 mesh = create_unit_square(MPI.COMM_WORLD, Nx, Ny, diagonal=DiagonalType.left, ghost_mode=GhostMode.none) mesh.topology.create_entities(1) def left_side(x, tol=1e-16): return x[0] <= 0.5 + tol cells = locate_entities(mesh, mesh.topology.dim, left_side) if MPI.COMM_WORLD.size == 0: assert cells.__len__() == Nx * Ny edges = compute_incident_entities(mesh, cells, 2, 1) if MPI.COMM_WORLD.size == 0: assert edges.__len__() == Nx // 2 * (2 * Ny + 1) + (Nx // 2 + 1) * Ny mesh2 = refine(mesh, edges, redistribute=True) num_cells_global = mesh2.topology.index_map(2).size_global actual_cells = 3 * (Nx * Ny) + 3 * Ny + 2 * Nx * Ny assert num_cells_global == actual_cells
def test_RefineUnitCubeMesh_keep_partition(): """Refine mesh of unit cube.""" mesh = UnitCubeMesh(MPI.COMM_WORLD, 5, 7, 9) mesh.topology.create_entities(1) mesh = refine(mesh, redistribute=False) assert mesh.topology.index_map(0).size_global == 3135 assert mesh.topology.index_map(3).size_global == 15120 Q = FunctionSpace(mesh, ("CG", 1)) assert(Q)
def test_Refinecreate_unit_cube_keep_partition(): """Refine mesh of unit cube.""" mesh = create_unit_cube(MPI.COMM_WORLD, 5, 7, 9, ghost_mode=GhostMode.none) mesh.topology.create_entities(1) mesh = refine(mesh, redistribute=False) assert mesh.topology.index_map(0).size_global == 3135 assert mesh.topology.index_map(3).size_global == 15120 Q = FunctionSpace(mesh, ("Lagrange", 1)) assert Q
def test_refine_create_form(): """Check that forms can be assembled on refined mesh""" mesh = dolfinx.UnitCubeMesh(MPI.COMM_WORLD, 3, 3, 3) mesh.topology.create_entities(1) mesh = refine(mesh, redistribute=True) V = dolfinx.FunctionSpace(mesh, ("CG", 1)) # Define variational problem u = ufl.TrialFunction(V) v = ufl.TestFunction(V) a = ufl.inner(ufl.grad(u), ufl.grad(v)) * ufl.dx dolfinx.fem.assemble_matrix(a)
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_sub_refine(): """Test that refinement of a subset of edges works""" mesh = create_unit_square(MPI.COMM_WORLD, 3, 4, diagonal=DiagonalType.left, ghost_mode=GhostMode.none) mesh.topology.create_entities(1) def left_corner_edge(x, tol=1e-16): return logical_and(isclose(x[0], 0), x[1] < 1 / 4 + tol) edges = locate_entities_boundary(mesh, 1, left_corner_edge) if MPI.COMM_WORLD.size == 0: assert edges == 1 mesh2 = refine(mesh, edges, redistribute=False) assert mesh.topology.index_map(2).size_global + 3 == mesh2.topology.index_map(2).size_global
def xtest_refinement_gdim(): """Test that 2D refinement is still 2D""" mesh = UnitSquareMesh(MPI.COMM_WORLD, 3, 4, ghost_mode=GhostMode.none) mesh2 = refine(mesh, redistribute=True) assert mesh.geometry.dim == mesh2.geometry.dim
def reference_periodic(tetra: bool, r_lvl: int = 0, out_hdf5: h5py.File = None, xdmf: bool = False, boomeramg: bool = False, kspview: bool = False, degree: int = 1): # Create mesh and finite element if tetra: # Tet setup N = 3 mesh = create_unit_cube(MPI.COMM_WORLD, N, N, N) for i in range(r_lvl): mesh.topology.create_entities(mesh.topology.dim - 2) mesh = refine(mesh, redistribute=True) N *= 2 else: # Hex setup N = 3 for i in range(r_lvl): N *= 2 mesh = create_unit_cube(MPI.COMM_WORLD, N, N, N, CellType.hexahedron) V = FunctionSpace(mesh, ("CG", degree)) # Create Dirichlet boundary condition def dirichletboundary(x): return np.logical_or( np.logical_or(np.isclose(x[1], 0), np.isclose(x[1], 1)), np.logical_or(np.isclose(x[2], 0), np.isclose(x[2], 1))) mesh.topology.create_connectivity(2, 1) geometrical_dofs = locate_dofs_geometrical(V, dirichletboundary) bc = dirichletbc(PETSc.ScalarType(0), geometrical_dofs, V) bcs = [bc] # Define variational problem u = TrialFunction(V) v = TestFunction(V) a = inner(grad(u), grad(v)) * dx x = SpatialCoordinate(mesh) dx_ = x[0] - 0.9 dy_ = x[1] - 0.5 dz_ = x[2] - 0.1 f = x[0] * sin(5.0 * pi * x[1]) + 1.0 * exp( -(dx_ * dx_ + dy_ * dy_ + dz_ * dz_) / 0.02) rhs = inner(f, v) * dx # Assemble rhs, RHS and apply lifting bilinear_form = form(a) linear_form = form(rhs) A_org = assemble_matrix(bilinear_form, bcs) A_org.assemble() L_org = assemble_vector(linear_form) apply_lifting(L_org, [bilinear_form], [bcs]) L_org.ghostUpdate(addv=PETSc.InsertMode.ADD_VALUES, mode=PETSc.ScatterMode.REVERSE) set_bc(L_org, bcs) # Create PETSc nullspace nullspace = PETSc.NullSpace().create(constant=True) PETSc.Mat.setNearNullSpace(A_org, nullspace) # Set PETSc options opts = PETSc.Options() if boomeramg: opts["ksp_type"] = "cg" opts["ksp_rtol"] = 1.0e-5 opts["pc_type"] = "hypre" opts['pc_hypre_type'] = 'boomeramg' opts["pc_hypre_boomeramg_max_iter"] = 1 opts["pc_hypre_boomeramg_cycle_type"] = "v" # opts["pc_hypre_boomeramg_print_statistics"] = 1 else: opts["ksp_type"] = "cg" opts["ksp_rtol"] = 1.0e-12 opts["pc_type"] = "gamg" opts["pc_gamg_type"] = "agg" opts["pc_gamg_sym_graph"] = True # Use Chebyshev smoothing for multigrid opts["mg_levels_ksp_type"] = "richardson" opts["mg_levels_pc_type"] = "sor" # opts["help"] = None # List all available options # opts["ksp_view"] = None # List progress of solver # Initialize PETSc solver, set options and operator solver = PETSc.KSP().create(MPI.COMM_WORLD) solver.setFromOptions() solver.setOperators(A_org) # Solve linear problem u_ = Function(V) start = perf_counter() with Timer("Solve"): solver.solve(L_org, u_.vector) end = perf_counter() u_.vector.ghostUpdate(addv=PETSc.InsertMode.INSERT, mode=PETSc.ScatterMode.FORWARD) if kspview: solver.view() it = solver.getIterationNumber() num_dofs = V.dofmap.index_map.size_global * V.dofmap.index_map_bs if out_hdf5 is not None: d_set = out_hdf5.get("its") d_set[r_lvl] = it d_set = out_hdf5.get("num_dofs") d_set[r_lvl] = num_dofs d_set = out_hdf5.get("solve_time") d_set[r_lvl, MPI.COMM_WORLD.rank] = end - start if MPI.COMM_WORLD.rank == 0: print("Rlvl {0:d}, Iterations {1:d}".format(r_lvl, it)) # Output solution to XDMF if xdmf: ext = "tet" if tetra else "hex" fname = "results/reference_periodic_{0:d}_{1:s}.xdmf".format( r_lvl, ext) u_.name = "u_" + ext + "_unconstrained" with XDMFFile(MPI.COMM_WORLD, fname, "w") as out_periodic: out_periodic.write_mesh(mesh) out_periodic.write_function( u_, 0.0, "Xdmf/Domain/" + "Grid[@Name='{0:s}'][1]".format(mesh.name))
def ref_elasticity(tetra: bool = True, r_lvl: int = 0, out_hdf5: h5py.File = None, xdmf: bool = False, boomeramg: bool = False, kspview: bool = False, degree: int = 1): if tetra: N = 3 if degree == 1 else 2 mesh = create_unit_cube(MPI.COMM_WORLD, N, N, N) else: N = 3 mesh = create_unit_cube(MPI.COMM_WORLD, N, N, N, CellType.hexahedron) for i in range(r_lvl): # set_log_level(LogLevel.INFO) N *= 2 if tetra: mesh = refine(mesh, redistribute=True) else: mesh = create_unit_cube(MPI.COMM_WORLD, N, N, N, CellType.hexahedron) # set_log_level(LogLevel.ERROR) N = degree * N fdim = mesh.topology.dim - 1 V = VectorFunctionSpace(mesh, ("Lagrange", int(degree))) # Generate Dirichlet BC on lower boundary (Fixed) u_bc = Function(V) with u_bc.vector.localForm() as u_local: u_local.set(0.0) def boundaries(x): return np.isclose(x[0], np.finfo(float).eps) facets = locate_entities_boundary(mesh, fdim, boundaries) topological_dofs = locate_dofs_topological(V, fdim, facets) bc = dirichletbc(u_bc, topological_dofs) bcs = [bc] # Create traction meshtag def traction_boundary(x): return np.isclose(x[0], 1) t_facets = locate_entities_boundary(mesh, fdim, traction_boundary) facet_values = np.ones(len(t_facets), dtype=np.int32) arg_sort = np.argsort(t_facets) mt = meshtags(mesh, fdim, t_facets[arg_sort], facet_values[arg_sort]) # Elasticity parameters E = PETSc.ScalarType(1.0e4) nu = 0.1 mu = Constant(mesh, E / (2.0 * (1.0 + nu))) lmbda = Constant(mesh, E * nu / ((1.0 + nu) * (1.0 - 2.0 * nu))) g = Constant(mesh, PETSc.ScalarType((0, 0, -1e2))) x = SpatialCoordinate(mesh) f = Constant(mesh, PETSc.ScalarType(1e4)) * \ as_vector((0, -(x[2] - 0.5)**2, (x[1] - 0.5)**2)) # Stress computation def sigma(v): return (2.0 * mu * sym(grad(v)) + lmbda * tr(sym(grad(v))) * Identity(len(v))) # Define variational problem u = TrialFunction(V) v = TestFunction(V) a = inner(sigma(u), grad(v)) * dx rhs = inner(g, v) * ds(domain=mesh, subdomain_data=mt, subdomain_id=1) + inner(f, v) * dx num_dofs = V.dofmap.index_map.size_global * V.dofmap.index_map_bs if MPI.COMM_WORLD.rank == 0: print("Problem size {0:d} ".format(num_dofs)) # Generate reference matrices and unconstrained solution bilinear_form = form(a) A_org = assemble_matrix(bilinear_form, bcs) A_org.assemble() null_space_org = rigid_motions_nullspace(V) A_org.setNearNullSpace(null_space_org) linear_form = form(rhs) L_org = assemble_vector(linear_form) apply_lifting(L_org, [bilinear_form], [bcs]) L_org.ghostUpdate(addv=PETSc.InsertMode.ADD_VALUES, mode=PETSc.ScatterMode.REVERSE) set_bc(L_org, bcs) opts = PETSc.Options() if boomeramg: opts["ksp_type"] = "cg" opts["ksp_rtol"] = 1.0e-5 opts["pc_type"] = "hypre" opts['pc_hypre_type'] = 'boomeramg' opts["pc_hypre_boomeramg_max_iter"] = 1 opts["pc_hypre_boomeramg_cycle_type"] = "v" # opts["pc_hypre_boomeramg_print_statistics"] = 1 else: opts["ksp_rtol"] = 1.0e-8 opts["pc_type"] = "gamg" opts["pc_gamg_type"] = "agg" opts["pc_gamg_coarse_eq_limit"] = 1000 opts["pc_gamg_sym_graph"] = True opts["mg_levels_ksp_type"] = "chebyshev" opts["mg_levels_pc_type"] = "jacobi" opts["mg_levels_esteig_ksp_type"] = "cg" opts["matptap_via"] = "scalable" opts["pc_gamg_square_graph"] = 2 opts["pc_gamg_threshold"] = 0.02 # opts["help"] = None # List all available options # opts["ksp_view"] = None # List progress of solver # Create solver, set operator and options solver = PETSc.KSP().create(MPI.COMM_WORLD) solver.setFromOptions() solver.setOperators(A_org) # Solve linear problem u_ = Function(V) start = perf_counter() with Timer("Ref solve"): solver.solve(L_org, u_.vector) end = perf_counter() u_.x.scatter_forward() if kspview: solver.view() it = solver.getIterationNumber() if out_hdf5 is not None: d_set = out_hdf5.get("its") d_set[r_lvl] = it d_set = out_hdf5.get("num_dofs") d_set[r_lvl] = num_dofs d_set = out_hdf5.get("solve_time") d_set[r_lvl, MPI.COMM_WORLD.rank] = end - start if MPI.COMM_WORLD.rank == 0: print("Refinement level {0:d}, Iterations {1:d}".format(r_lvl, it)) # List memory usage mem = sum(MPI.COMM_WORLD.allgather( resource.getrusage(resource.RUSAGE_SELF).ru_maxrss)) if MPI.COMM_WORLD.rank == 0: print("{1:d}: Max usage after trad. solve {0:d} (kb)" .format(mem, r_lvl)) if xdmf: # Name formatting of functions u_.name = "u_unconstrained" fname = "results/ref_elasticity_{0:d}.xdmf".format(r_lvl) with XDMFFile(MPI.COMM_WORLD, fname, "w") as out_xdmf: out_xdmf.write_mesh(mesh) out_xdmf.write_function(u_, 0.0, "Xdmf/Domain/Grid[@Name='{0:s}'][1]".format(mesh.name))
def mesh_3D_dolfin(theta=0, ct=CellType.tetrahedron, ext="tetrahedron", num_refinements=0, N0=5): timer = Timer("Create mesh") def find_plane_function(p0, p1, p2): """ Find plane function given three points: http://www.nabla.hr/CG-LinesPlanesIn3DA3.htm """ v1 = np.array(p1) - np.array(p0) v2 = np.array(p2) - np.array(p0) n = np.cross(v1, v2) D = -(n[0] * p0[0] + n[1] * p0[1] + n[2] * p0[2]) return lambda x: np.isclose(0, np.dot(n, x) + D) def over_plane(p0, p1, p2): """ Returns function that checks if a point is over a plane defined by the points p0, p1 and p2. """ v1 = np.array(p1) - np.array(p0) v2 = np.array(p2) - np.array(p0) n = np.cross(v1, v2) D = -(n[0] * p0[0] + n[1] * p0[1] + n[2] * p0[2]) return lambda x: n[0] * x[0] + n[1] * x[1] + D > -n[2] * x[2] tmp_mesh_name = "tmp_mesh.xdmf" r_matrix = rotation_matrix([1 / np.sqrt(2), 1 / np.sqrt(2), 0], -theta) if MPI.COMM_WORLD.rank == 0: # Create two coarse meshes and merge them mesh0 = create_unit_cube(MPI.COMM_SELF, N0, N0, N0, ct) mesh0.geometry.x[:, 2] += 1 mesh1 = create_unit_cube(MPI.COMM_SELF, 2 * N0, 2 * N0, 2 * N0, ct) tdim0 = mesh0.topology.dim num_cells0 = mesh0.topology.index_map(tdim0).size_local cells0 = entities_to_geometry( mesh0, tdim0, np.arange(num_cells0, dtype=np.int32).reshape((-1, 1)), False) tdim1 = mesh1.topology.dim num_cells1 = mesh1.topology.index_map(tdim1).size_local cells1 = entities_to_geometry( mesh1, tdim1, np.arange(num_cells1, dtype=np.int32).reshape((-1, 1)), False) cells1 += mesh0.geometry.x.shape[0] # Concatenate points and cells points = np.vstack([mesh0.geometry.x, mesh1.geometry.x]) cells = np.vstack([cells0, cells1]) cell = Cell(ext, geometric_dimension=points.shape[1]) domain = Mesh(VectorElement("Lagrange", cell, 1)) # Rotate mesh points = np.dot(r_matrix, points.T).T mesh = create_mesh(MPI.COMM_SELF, cells, points, domain) with XDMFFile(MPI.COMM_SELF, tmp_mesh_name, "w") as xdmf: xdmf.write_mesh(mesh) MPI.COMM_WORLD.barrier() with XDMFFile(MPI.COMM_WORLD, tmp_mesh_name, "r") as xdmf: mesh = xdmf.read_mesh() # Refine coarse mesh for i in range(num_refinements): mesh.topology.create_entities(mesh.topology.dim - 2) mesh = refine(mesh, redistribute=True) tdim = mesh.topology.dim fdim = tdim - 1 # Find information about facets to be used in meshtags bottom_points = np.dot( r_matrix, np.array([[0, 0, 0], [1, 0, 0], [0, 1, 0], [1, 1, 0]]).T) bottom = find_plane_function(bottom_points[:, 0], bottom_points[:, 1], bottom_points[:, 2]) bottom_facets = locate_entities_boundary(mesh, fdim, bottom) top_points = np.dot( r_matrix, np.array([[0, 0, 2], [1, 0, 2], [0, 1, 2], [1, 1, 2]]).T) top = find_plane_function(top_points[:, 0], top_points[:, 1], top_points[:, 2]) top_facets = locate_entities_boundary(mesh, fdim, top) # Determine interface facets if_points = np.dot( r_matrix, np.array([[0, 0, 1], [1, 0, 1], [0, 1, 1], [1, 1, 1]]).T) interface = find_plane_function(if_points[:, 0], if_points[:, 1], if_points[:, 2]) i_facets = locate_entities_boundary(mesh, fdim, interface) mesh.topology.create_connectivity(fdim, tdim) top_interface = [] bottom_interface = [] facet_to_cell = mesh.topology.connectivity(fdim, tdim) num_cells = mesh.topology.index_map(tdim).size_local # Find top and bottom interface facets cell_midpoints = compute_midpoints(mesh, tdim, range(num_cells)) top_cube = over_plane(if_points[:, 0], if_points[:, 1], if_points[:, 2]) for facet in i_facets: i_cells = facet_to_cell.links(facet) assert (len(i_cells == 1)) i_cell = i_cells[0] if top_cube(cell_midpoints[i_cell]): top_interface.append(facet) else: bottom_interface.append(facet) # Create cell tags num_cells = mesh.topology.index_map(tdim).size_local cell_midpoints = compute_midpoints(mesh, tdim, range(num_cells)) top_cube_marker = 2 indices = [] values = [] for cell_index in range(num_cells): if top_cube(cell_midpoints[cell_index]): indices.append(cell_index) values.append(top_cube_marker) ct = meshtags(mesh, tdim, np.array(indices, dtype=np.intc), np.array(values, dtype=np.intc)) # Create meshtags for facet data markers = { 3: top_facets, 4: bottom_interface, 9: top_interface, 5: bottom_facets } # , 6: left_facets, 7: right_facets} indices = np.array([], dtype=np.intc) values = np.array([], dtype=np.intc) for key in markers.keys(): indices = np.append(indices, markers[key]) values = np.append(values, np.full(len(markers[key]), key, dtype=np.intc)) sorted_indices = np.argsort(indices) mt = meshtags(mesh, fdim, indices[sorted_indices], values[sorted_indices]) mt.name = "facet_tags" fname = f"meshes/mesh_{ext}_{theta:.2f}.xdmf" with XDMFFile(MPI.COMM_WORLD, fname, "w") as o_f: o_f.write_mesh(mesh) o_f.write_meshtags(ct) o_f.write_meshtags(mt) timer.stop()
def test_refinement_gdim(): """Test that 2D refinement is still 2D""" mesh = create_unit_square(MPI.COMM_WORLD, 3, 4, ghost_mode=GhostMode.none) mesh.topology.create_entities(1) mesh2 = refine(mesh, redistribute=True) assert mesh.geometry.dim == mesh2.geometry.dim
def bench_elasticity_one(r_lvl: int = 0, out_hdf5: h5py.File = None, xdmf: bool = False, boomeramg: bool = False, kspview: bool = False): N = 3 mesh = create_unit_cube(MPI.COMM_WORLD, N, N, N) for i in range(r_lvl): mesh.topology.create_entities(mesh.topology.dim - 2) mesh = refine(mesh, redistribute=True) fdim = mesh.topology.dim - 1 V = VectorFunctionSpace(mesh, ("Lagrange", 1)) # Generate Dirichlet BC on lower boundary (Fixed) u_bc = Function(V) with u_bc.vector.localForm() as u_local: u_local.set(0.0) def boundaries(x): return np.isclose(x[0], np.finfo(float).eps) facets = locate_entities_boundary(mesh, fdim, boundaries) topological_dofs = locate_dofs_topological(V, fdim, facets) bc = dirichletbc(u_bc, topological_dofs) bcs = [bc] # Create traction meshtag def traction_boundary(x): return np.isclose(x[0], 1) t_facets = locate_entities_boundary(mesh, fdim, traction_boundary) facet_values = np.ones(len(t_facets), dtype=np.int32) arg_sort = np.argsort(t_facets) mt = meshtags(mesh, fdim, t_facets[arg_sort], facet_values[arg_sort]) # Define variational problem u = TrialFunction(V) v = TestFunction(V) # Elasticity parameters E = 1.0e4 nu = 0.1 mu = Constant(mesh, E / (2.0 * (1.0 + nu))) lmbda = Constant(mesh, E * nu / ((1.0 + nu) * (1.0 - 2.0 * nu))) g = Constant(mesh, (0, 0, -1e2)) # Stress computation def sigma(v): return 2.0 * mu * sym(grad(v)) + lmbda * tr(sym(grad(v))) * Identity( len(v)) # Define variational problem u = TrialFunction(V) v = TestFunction(V) a = inner(sigma(u), grad(v)) * dx rhs = inner(g, v) * ds(domain=mesh, subdomain_data=mt, subdomain_id=1) # Create MPC with Timer("~Elasticity: Init constraint"): def l2b(li): return np.array(li, dtype=np.float64).tobytes() s_m_c = {l2b([1, 0, 0]): {l2b([1, 0, 1]): 0.5}} mpc = MultiPointConstraint(V) mpc.create_general_constraint(s_m_c, 2, 2) mpc.finalize() # Setup MPC system bilinear_form = form(a) linear_form = form(rhs) with Timer("~Elasticity: Assemble LHS and RHS"): A = assemble_matrix(bilinear_form, mpc, bcs=bcs) b = assemble_vector(linear_form, mpc) # Apply boundary conditions apply_lifting(b, [bilinear_form], [bcs], mpc) b.ghostUpdate(addv=PETSc.InsertMode.ADD_VALUES, mode=PETSc.ScatterMode.REVERSE) set_bc(b, bcs) # Create functionspace and function for mpc vector # Solve Linear problem solver = PETSc.KSP().create(MPI.COMM_WORLD) opts = PETSc.Options() if boomeramg: opts["ksp_type"] = "cg" opts["ksp_rtol"] = 1.0e-5 opts["pc_type"] = "hypre" opts['pc_hypre_type'] = 'boomeramg' opts["pc_hypre_boomeramg_max_iter"] = 1 opts["pc_hypre_boomeramg_cycle_type"] = "v" # opts["pc_hypre_boomeramg_print_statistics"] = 1 else: opts["ksp_rtol"] = 1.0e-10 opts["pc_type"] = "gamg" opts["pc_gamg_type"] = "agg" opts["pc_gamg_coarse_eq_limit"] = 1000 opts["pc_gamg_sym_graph"] = True opts["mg_levels_ksp_type"] = "chebyshev" opts["mg_levels_pc_type"] = "jacobi" opts["mg_levels_esteig_ksp_type"] = "cg" opts["matptap_via"] = "scalable" # opts["help"] = None # List all available options # opts["ksp_view"] = None # List progress of solver with Timer("~Elasticity: Solve problem") as timer: null_space = rigid_motions_nullspace(mpc.function_space) A.setNearNullSpace(null_space) solver.setFromOptions() solver.setOperators(A) # Solve linear problem uh = b.copy() uh.set(0) solver.solve(b, uh) uh.ghostUpdate(addv=PETSc.InsertMode.INSERT, mode=PETSc.ScatterMode.FORWARD) mpc.backsubstitution(uh) solver_time = timer.elapsed() it = solver.getIterationNumber() if kspview: solver.view() # Print max usage of summary mem = sum( MPI.COMM_WORLD.allgather( resource.getrusage(resource.RUSAGE_SELF).ru_maxrss)) num_dofs = V.dofmap.index_map.size_global * V.dofmap.index_map_bs if MPI.COMM_WORLD.rank == 0: print(f"Rlvl {r_lvl}, Iterations {it}") print(f"Rlvl {r_lvl}, Max usage {mem} (kb), #dofs {num_dofs}") if out_hdf5 is not None: d_set = out_hdf5.get("its") d_set[r_lvl] = it d_set = out_hdf5.get("num_dofs") d_set[r_lvl] = num_dofs d_set = out_hdf5.get("solve_time") d_set[r_lvl, MPI.COMM_WORLD.rank] = solver_time[0] if xdmf: # Write solution to file u_h = Function(mpc.function_space) u_h.vector.setArray(uh.array) u_h.name = "u_mpc" fname = f"results/bench_elasticity_{r_lvl}.xdmf" with XDMFFile(MPI.COMM_WORLD, fname, "w") as outfile: outfile.write_mesh(mesh) outfile.write_function(u_h)