//----------------------------------------------------------------------------- void MeshSmoothing::move_interior_vertices(Mesh& mesh, BoundaryMesh& boundary, bool harmonic_smoothing) { // Select smoothing of interior vertices if (harmonic_smoothing) { std::shared_ptr<Mesh> _mesh(&mesh, [](Mesh*){}); ALE::move(_mesh, boundary); } else { if (mesh.geometry().degree() != 1) { dolfin_error("MeshSmoothing.cpp", "move interior vertices", "This function does not support higher-order mesh geometry"); } // Use vertex map to update boundary coordinates of original mesh const MeshFunction<std::size_t>& vertex_map = boundary.entity_map(0); for (VertexIterator v(boundary); !v.end(); ++v) mesh.geometry().set(vertex_map[*v], v->x()); } }
//----------------------------------------------------------------------------- void MeshSmoothing::move_interior_vertices(Mesh& mesh, BoundaryMesh& boundary, bool harmonic_smoothing) { // Select smoothing of interior vertices if (harmonic_smoothing) ALE::move(mesh, boundary); else { // Use vertex map to update boundary coordinates of original mesh const MeshFunction<std::size_t>& vertex_map = boundary.entity_map(0); for (VertexIterator v(boundary); !v.end(); ++v) mesh.geometry().set(vertex_map[*v], v->x()); } }
//----------------------------------------------------------------------------- void BoundaryComputation::compute_boundary(const Mesh& mesh, const std::string type, BoundaryMesh& boundary) { // We iterate over all facets in the mesh and check if they are on // the boundary. A facet is on the boundary if it is connected to // exactly one cell. log(TRACE, "Computing boundary mesh."); bool exterior = true; bool interior = true; if (type == "exterior") interior = false; else if (type == "interior") exterior = false; else if (type != "local") { dolfin_error("BoundaryComputation.cpp", "determine boundary mesh type", "Unknown boundary type (%d)", type.c_str()); } // Get my MPI process rank and number of MPI processes const std::size_t my_rank = MPI::rank(mesh.mpi_comm()); const std::size_t num_processes = MPI::size(mesh.mpi_comm()); // Open boundary mesh for editing const std::size_t D = mesh.topology().dim(); MeshEditor editor; editor.open(boundary, mesh.type().facet_type(), D - 1, mesh.geometry().dim()); // Generate facet - cell connectivity if not generated mesh.init(D - 1, D); // Temporary arrays for assignment of indices to vertices on the boundary std::map<std::size_t, std::size_t> boundary_vertices; // Map of index "owners" (process responsible for assigning global index) std::map< std::size_t, std::size_t > global_index_owner; // Shared vertices for full mesh // FIXME: const_cast const std::map<unsigned int, std::set<unsigned int>> & shared_vertices = const_cast<Mesh&>(mesh).topology().shared_entities(0); // Shared vertices for boundary mesh std::map<unsigned int, std::set<unsigned int>> shared_boundary_vertices; if (exterior) { // Extract shared vertices if vertex is identified as part of globally // exterior facet. std::vector<std::size_t> boundary_global_indices; for (std::map<unsigned int, std::set<unsigned int>>::const_iterator sv_it=shared_vertices.begin(); sv_it != shared_vertices.end(); ++sv_it) { std::size_t local_mesh_index = sv_it->first; Vertex v(mesh, local_mesh_index); for (FacetIterator f(v); !f.end(); ++f) { if (f->num_global_entities(D) == 1) { const std::size_t global_mesh_index = mesh.topology().global_indices(0)[local_mesh_index]; shared_boundary_vertices[local_mesh_index] = sv_it->second; boundary_global_indices.push_back(global_mesh_index); break; } } } // Distribute all shared boundary vertices std::vector<std::vector<std::size_t>> boundary_global_indices_all; MPI::all_gather(mesh.mpi_comm(), boundary_global_indices, boundary_global_indices_all); // Identify and clean up discrepancies between shared vertices of full mesh // and shared vertices of boundary mesh for (auto sbv_it = shared_boundary_vertices.begin(); sbv_it != shared_boundary_vertices.end(); ) { std::size_t local_mesh_index = sbv_it->first; const std::size_t global_mesh_index = mesh.topology().global_indices(0)[local_mesh_index]; // Check if this vertex is identified as boundary vertex on // other processes sharing this vertex std::set<unsigned int> &other_processes = sbv_it->second; for (auto op_it=other_processes.begin(); op_it != other_processes.end(); ) { // Check if vertex is identified as boundary vertex on process *op_it bool is_boundary_vertex = (std::find(boundary_global_indices_all[*op_it].begin(), boundary_global_indices_all[*op_it].end(), global_mesh_index) != boundary_global_indices_all[*op_it].end()); // Erase item if this is not identified as a boundary vertex // on process *op_it, and increment iterator if (!is_boundary_vertex) { // Erase item while carefully avoiding invalidating the // iterator: First increment it to get the next, valid // iterator, and then erase what it pointed to from // other_processes other_processes.erase(op_it++); } else ++op_it; } // Erase item from map if no other processes identify this // vertex as a boundary vertex, and increment iterator if (other_processes.size() == 0) { // Erase carefully as above shared_boundary_vertices.erase(sbv_it++); } else ++sbv_it; } } else { // If interior boundary, shared vertices are the same shared_boundary_vertices = shared_vertices; } // Determine boundary facet, count boundary vertices and facets, and // assign vertex indices std::size_t num_boundary_vertices = 0; std::size_t num_owned_vertices = 0; std::size_t num_boundary_cells = 0; MeshFunction<bool> boundary_facet(reference_to_no_delete_pointer(mesh), D - 1, false); for (FacetIterator f(mesh); !f.end(); ++f) { // Boundary facets are connected to exactly one cell if (f->num_entities(D) == 1) { const bool global_exterior_facet = (f->num_global_entities(D) == 1); if (global_exterior_facet && exterior) boundary_facet[*f] = true; else if (!global_exterior_facet && interior) boundary_facet[*f] = true; if (boundary_facet[*f]) { // Count boundary vertices and assign indices for (VertexIterator v(*f); !v.end(); ++v) { const std::size_t local_mesh_index = v->index(); if (boundary_vertices.find(local_mesh_index) == boundary_vertices.end()) { const std::size_t local_boundary_index = num_boundary_vertices; boundary_vertices[local_mesh_index] = local_boundary_index; // Determine "owner" of global_mesh_index std::size_t owner = my_rank; std::map<unsigned int, std::set<unsigned int>>::const_iterator other_processes_it = shared_boundary_vertices.find(local_mesh_index); if (other_processes_it != shared_boundary_vertices.end() && D > 1) { const std::set<unsigned int>& other_processes = other_processes_it->second; const std::size_t min_process = *std::min_element(other_processes.begin(), other_processes.end()); boundary.topology().shared_entities(0)[local_boundary_index] = other_processes; // FIXME: More sophisticated ownership determination if (min_process < owner) owner = min_process; } const std::size_t global_mesh_index = mesh.topology().global_indices(0)[local_mesh_index]; global_index_owner[global_mesh_index] = owner; // Update counts if (owner == my_rank) num_owned_vertices++; num_boundary_vertices++; } } // Count boundary cells (facets of the mesh) num_boundary_cells++; } } } // Initiate boundary topology /* boundary.topology().init(0, num_boundary_vertices, MPI::sum(mesh.mpi_comm(), num_owned_vertices)); boundary.topology().init(D - 1, num_boundary_cells, MPI::sum(mesh.mpi_comm(), num_boundary_cells)); */ // Specify number of vertices and cells editor.init_vertices_global(num_boundary_vertices, MPI::sum(mesh.mpi_comm(), num_owned_vertices)); editor.init_cells_global(num_boundary_cells, MPI::sum(mesh.mpi_comm(), num_boundary_cells)); // Write vertex map MeshFunction<std::size_t>& vertex_map = boundary.entity_map(0); if (num_boundary_vertices > 0) { vertex_map.init(reference_to_no_delete_pointer(boundary), 0, num_boundary_vertices); } std::map<std::size_t, std::size_t>::const_iterator it; for (it = boundary_vertices.begin(); it != boundary_vertices.end(); ++it) vertex_map[it->second] = it->first; // Get vertex ownership distribution, and find index to start global // numbering from std::vector<std::size_t> ownership_distribution(num_processes); MPI::all_gather(mesh.mpi_comm(), num_owned_vertices, ownership_distribution); std::size_t start_index = 0; for (std::size_t j = 0; j < my_rank; j++) start_index += ownership_distribution[j]; // Set global indices of owned vertices, request global indices for // vertices owned elsewhere std::map<std::size_t, std::size_t> global_indices; std::vector<std::vector<std::size_t>> request_global_indices(num_processes); std::size_t current_index = start_index; for (std::size_t local_boundary_index = 0; local_boundary_index<num_boundary_vertices; local_boundary_index++) { const std::size_t local_mesh_index = vertex_map[local_boundary_index]; const std::size_t global_mesh_index = mesh.topology().global_indices(0)[local_mesh_index]; const std::size_t owner = global_index_owner[global_mesh_index]; if (owner != my_rank) request_global_indices[owner].push_back(global_mesh_index); else global_indices[global_mesh_index] = current_index++; } // Send and receive requests from other processes std::vector<std::vector<std::size_t>> global_index_requests(num_processes); MPI::all_to_all(mesh.mpi_comm(), request_global_indices, global_index_requests); // Find response to requests of global indices std::vector<std::vector<std::size_t>> respond_global_indices(num_processes); for (std::size_t i = 0; i < num_processes; i++) { const std::size_t N = global_index_requests[i].size(); respond_global_indices[i].resize(N); for (std::size_t j = 0; j < N; j++) respond_global_indices[i][j] = global_indices[global_index_requests[i][j]]; } // Scatter responses back to requesting processes std::vector<std::vector<std::size_t>> global_index_responses(num_processes); MPI::all_to_all(mesh.mpi_comm(), respond_global_indices, global_index_responses); // Update global_indices for (std::size_t i = 0; i < num_processes; i++) { const std::size_t N = global_index_responses[i].size(); // Check that responses are the same size as the requests made dolfin_assert(global_index_responses[i].size() == request_global_indices[i].size()); for (std::size_t j = 0; j < N; j++) { const std::size_t global_mesh_index = request_global_indices[i][j]; const std::size_t global_boundary_index = global_index_responses[i][j]; global_indices[global_mesh_index] = global_boundary_index; } } // Create vertices for (std::size_t local_boundary_index = 0; local_boundary_index < num_boundary_vertices; local_boundary_index++) { const std::size_t local_mesh_index = vertex_map[local_boundary_index]; const std::size_t global_mesh_index = mesh.topology().global_indices(0)[local_mesh_index]; const std::size_t global_boundary_index = global_indices[global_mesh_index]; Vertex v(mesh, local_mesh_index); editor.add_vertex_global(local_boundary_index, global_boundary_index, v.point()); } // Find global index to start cell numbering from for current process std::vector<std::size_t> cell_distribution(num_processes); MPI::all_gather(mesh.mpi_comm(), num_boundary_cells, cell_distribution); std::size_t start_cell_index = 0; for (std::size_t i = 0; i < my_rank; i++) start_cell_index += cell_distribution[i]; // Create cells (facets) and map between boundary mesh cells and facets parent MeshFunction<std::size_t>& cell_map = boundary.entity_map(D - 1); if (num_boundary_cells > 0) { cell_map.init(reference_to_no_delete_pointer(boundary), D - 1, num_boundary_cells); } std::vector<std::size_t> cell(boundary.type().num_vertices(boundary.topology().dim())); std::size_t current_cell = 0; for (FacetIterator f(mesh); !f.end(); ++f) { if (boundary_facet[*f]) { // Compute new vertex numbers for cell const unsigned int* vertices = f->entities(0); for (std::size_t i = 0; i < cell.size(); i++) cell[i] = boundary_vertices[vertices[i]]; // Reorder vertices so facet is right-oriented w.r.t. facet // normal reorder(cell, *f); // Create mapping from boundary cell to mesh facet if requested if (!cell_map.empty()) cell_map[current_cell] = f->index(); // Add cell editor.add_cell(current_cell, start_cell_index+current_cell, cell); current_cell++; } } // Close mesh editor. Note the argument order=false to prevent // ordering from destroying the orientation of facets accomplished // by calling reorder() below. editor.close(false); }
//----------------------------------------------------------------------------- std::shared_ptr<MeshDisplacement> HarmonicSmoothing::move(std::shared_ptr<Mesh> mesh, const BoundaryMesh& new_boundary) { dolfin_assert(mesh); // Now this works regardless of reorder_dofs_serial value const bool reorder_dofs_serial = parameters["reorder_dofs_serial"]; if (!reorder_dofs_serial) { warning("The function HarmonicSmoothing::move no longer needs " "parameters[\"reorder_dofs_serial\"] = false"); } const std::size_t D = mesh->topology().dim(); const std::size_t d = mesh->geometry().dim(); // Choose form and function space std::shared_ptr<FunctionSpace> V; std::shared_ptr<Form> form; switch (D) { case 1: V.reset(new Poisson1D::FunctionSpace(mesh)); form.reset(new Poisson1D::BilinearForm(V, V)); break; case 2: V.reset(new Poisson2D::FunctionSpace(mesh)); form.reset(new Poisson2D::BilinearForm(V, V)); break; case 3: V.reset(new Poisson3D::FunctionSpace(mesh)); form.reset(new Poisson3D::BilinearForm(V, V)); break; default: dolfin_error("HarmonicSmoothing.cpp", "move mesh using harmonic smoothing", "Illegal mesh dimension (%d)", D); } // Assemble matrix auto A = std::make_shared<Matrix>(); assemble(*A, *form); // Number of mesh vertices (local) const std::size_t num_vertices = mesh->num_vertices(); // Dof range const dolfin::la_index n0 = V->dofmap()->ownership_range().first; const dolfin::la_index n1 = V->dofmap()->ownership_range().second; const dolfin::la_index num_owned_dofs = n1 - n0; // Mapping of new_boundary vertex numbers to mesh vertex numbers const MeshFunction<std::size_t>& vertex_map_mesh_func = new_boundary.entity_map(0); const std::size_t num_boundary_vertices = vertex_map_mesh_func.size(); const std::vector<std::size_t> vertex_map(vertex_map_mesh_func.values(), vertex_map_mesh_func.values() + num_boundary_vertices); // Mapping of mesh vertex numbers to dofs (including ghost dofs) const std::vector<dolfin::la_index> vertex_to_dofs = vertex_to_dof_map(*V); // Array of all dofs (including ghosts) with local numbering std::vector<dolfin::la_index> all_global_dofs(num_vertices); for (std::size_t i = 0; i < num_vertices; i++) all_global_dofs[i] = vertex_to_dofs[i]; // Create arrays for setting bcs. Their indexing does not matter - // same ordering does. std::size_t num_boundary_dofs = 0; std::vector<dolfin::la_index> boundary_dofs; boundary_dofs.reserve(num_boundary_vertices); std::vector<std::size_t> boundary_vertices; boundary_vertices.reserve(num_boundary_vertices); for (std::size_t vert = 0; vert < num_boundary_vertices; vert++) { // Skip ghosts const dolfin::la_index dof = vertex_to_dofs[vertex_map[vert]]; if (dof < num_owned_dofs) { // Global dof numbers boundary_dofs.push_back(dof + n0); // new_boundary vertex indices boundary_vertices.push_back(vert); num_boundary_dofs++; } } // Modify matrix (insert 1 on diagonal) A->ident(num_boundary_dofs, boundary_dofs.data()); A->apply("insert"); // Arrays for storing Dirichlet condition and solution std::vector<double> boundary_values(num_boundary_dofs); std::vector<double> displacement; displacement.reserve(d*num_vertices); // Pick amg as preconditioner if available const std::string prec(has_krylov_solver_preconditioner("amg") ? "amg" : "default"); // Displacement solution wrapped in Expression subclass // MeshDisplacement std::shared_ptr<MeshDisplacement> u(new MeshDisplacement(mesh)); // RHS vector Vector b(*(*u)[0].vector()); // Prepare solver // NOTE: GMRES needs to be used until Eigen a4b7b6e or 8dcc4ed is widespread; // afterwards CG can be used again KrylovSolver solver("bicgstab", prec); solver.parameters["nonzero_initial_guess"] = true; solver.set_operator(A); // Solve system for each dimension for (std::size_t dim = 0; dim < d; dim++) { // Get solution vector std::shared_ptr<GenericVector> x = (*u)[dim].vector(); if (dim > 0) b.zero(); // Store bc into RHS and solution so that CG solver can be used for (std::size_t i = 0; i < num_boundary_dofs; i++) { boundary_values[i] = new_boundary.geometry().x(boundary_vertices[i], dim) - mesh->geometry().x(vertex_map[boundary_vertices[i]], dim); } b.set(boundary_values.data(), num_boundary_dofs, boundary_dofs.data()); b.apply("insert"); *x = b; // Solve the system solver.solve(*x, b); // Get displacement std::vector<double> _displacement(num_vertices); x->get_local(_displacement.data(), num_vertices, all_global_dofs.data()); displacement.insert(displacement.end(), _displacement.begin(), _displacement.end()); } // Modify mesh coordinates MeshGeometry& geometry = mesh->geometry(); std::vector<double> coord(d); for (std::size_t i = 0; i < num_vertices; i++) { for (std::size_t dim = 0; dim < d; dim++) coord[dim] = displacement[dim*num_vertices + i] + geometry.x(i, dim); geometry.set(i, coord.data()); } // Return calculated displacement return u; }