template <int dim>
void LaplaceBeltrami<dim>::make_mesh ()
{
- std::map<typename Triangulation<dim-1,dim>::cell_iterator,
- typename Triangulation<dim,dim>::face_iterator>
- surface_to_volume_mapping;
-
HyperBallBoundary<dim> boundary_description;
Triangulation<dim> volume_mesh;
GridGenerator::half_hyper_ball(volume_mesh);
boundary_ids.insert(0);
GridTools::extract_boundary_mesh (volume_mesh, triangulation,
- surface_to_volume_mapping,
boundary_ids);
triangulation.refine_global(3);
* mesh"). The boundary to be extracted
* is specified by a list of
* boundary_ids. If none is specified
- * the whole boundary will be extracted.
+ * the whole boundary will be
+ * extracted. The function is used in
+ * step-38.
*
* It also builds a mapping linking the
* cells on the surface mesh to the
- * corresponding faces on the volume one.
+ * corresponding faces on the volume
+ * one. This mapping is the return value
+ * of the function.
*
* @note The function builds the surface
* mesh by creating a coarse mesh from
* curved boundaries.
*/
template <template <int,int> class Container, int dim, int spacedim>
- static void
+ static
+ std::map<typename Container<dim-1,spacedim>::cell_iterator,
+ typename Container<dim,spacedim>::face_iterator>
extract_boundary_mesh (const Container<dim,spacedim> &volume_mesh,
Container<dim-1,spacedim> &surface_mesh,
- std::map<typename Container<dim-1,spacedim>::cell_iterator,
- typename Container<dim,spacedim>::face_iterator>
- &surface_to_volume_mapping,
const std::set<unsigned char> &boundary_ids
= std::set<unsigned char>());
template <template <int,int> class Container, int dim, int spacedim>
-void
+std::map<typename Container<dim-1,spacedim>::cell_iterator,
+ typename Container<dim,spacedim>::face_iterator>
GridTools::extract_boundary_mesh (const Container<dim,spacedim> &volume_mesh,
Container<dim-1,spacedim> &surface_mesh,
- std::map<typename Container<dim-1,spacedim>::cell_iterator,
- typename Container<dim,spacedim>::face_iterator>
- &surface_to_volume_mapping,
const std::set<unsigned char> &boundary_ids)
{
// Assumption:
// We are relying below on the fact that Triangulation::create_triangulation(...) will keep the order
// pass by CellData and that it will not reorder the vertices.
- surface_to_volume_mapping.clear();
+ std::map<typename Container<dim-1,spacedim>::cell_iterator,
+ typename Container<dim,spacedim>::face_iterator>
+ surface_to_volume_mapping;
const unsigned int boundary_dim = dim-1; //dimension of the boundary mesh
break;
}
while (true);
+
+ return surface_to_volume_mapping;
}
// explicit instantiations
{
#if deal_II_dimension != 1
template
- void
+ std::map< Triangulation<deal_II_dimension-1,deal_II_dimension>::cell_iterator,
+ Triangulation<deal_II_dimension>::face_iterator>
GridTools::extract_boundary_mesh (const Triangulation<deal_II_dimension> &volume_mesh,
Triangulation<deal_II_dimension-1,deal_II_dimension> &surface_mesh,
- std::map< Triangulation<deal_II_dimension-1,deal_II_dimension>::cell_iterator,
- Triangulation<deal_II_dimension>::face_iterator>
- &surface_to_volume_mapping,
const std::set<unsigned char> &boundary_ids);
#endif
}
#if deal_II_dimension != 1
template
- void
+ std::map< Container<deal_II_dimension-1,deal_II_dimension>::cell_iterator,
+ Container<deal_II_dimension>::face_iterator>
GridTools::extract_boundary_mesh (const Container<deal_II_dimension> &volume_mesh,
Container<deal_II_dimension-1,deal_II_dimension> &surface_mesh,
- std::map< Container<deal_II_dimension-1,deal_II_dimension>::cell_iterator,
- Container<deal_II_dimension>::face_iterator>
- &surface_to_volume_mapping,
const std::set<unsigned char> &boundary_ids);
#endif
- }
\ No newline at end of file
+ }
const unsigned int dim = spacedim-1;
Triangulation<dim,spacedim> boundary_mesh;
- map<Triangulation<dim,spacedim>::cell_iterator,Triangulation<spacedim,spacedim>::face_iterator >
- surface_to_volume_mapping;
Triangulation<spacedim> volume_mesh;
GridGenerator::hyper_cube(volume_mesh);
- GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh,
- surface_to_volume_mapping);
+ GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh);
for (Triangulation<dim,spacedim>::active_cell_iterator
cell = boundary_mesh.begin_active();
cell != boundary_mesh.end(); ++cell)
const unsigned int dim = spacedim-1;
Triangulation<dim,spacedim> boundary_mesh;
- map<Triangulation<dim,spacedim>::cell_iterator,Triangulation<spacedim,spacedim>::face_iterator >
- surface_to_volume_mapping;
Triangulation<spacedim> volume_mesh;
GridGenerator::hyper_cube(volume_mesh);
- GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh,
- surface_to_volume_mapping);
+ GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh);
for (Triangulation<dim,spacedim>::active_cell_iterator
cell = boundary_mesh.begin_active();
cell != boundary_mesh.end(); ++cell)
Triangulation<dim-1,dim> boundary_mesh;
-
- GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh,
- surface_to_volume_mapping);
+ surface_to_volume_mapping
+ = GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh);
if (test_vertices_orientation(boundary_mesh, surface_to_volume_mapping))
deallog << "Passed.";
Triangulation<dim-1,dim> boundary_mesh;
-
- GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh,
- surface_to_volume_mapping);
+ surface_to_volume_mapping
+ = GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh);
if (test_vertices_orientation(boundary_mesh, surface_to_volume_mapping))
deallog << "Passed.";
set<unsigned char> boundary_ids;
boundary_ids.insert(0);
- GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh,
- surface_to_volume_mapping,
- boundary_ids);
+ surface_to_volume_mapping
+ = GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh,
+ boundary_ids);
if (test_vertices_orientation(boundary_mesh, surface_to_volume_mapping))
deallog << "Passed.";
Triangulation<dim-1,dim> boundary_mesh;
- GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh,
- surface_to_volume_mapping);
+ surface_to_volume_mapping
+ = GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh);
Assert (test_vertices_orientation(boundary_mesh, surface_to_volume_mapping,2),
ExcInternalError());
Triangulation<dim-1,dim> boundary_mesh;
boundary_mesh.set_boundary (0, surface_description);
- GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh,
- surface_to_volume_mapping);
+ surface_to_volume_mapping
+ = GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh);
deallog << volume_mesh.n_active_cells () << std::endl;
deallog << boundary_mesh.n_active_cells () << std::endl;
save_mesh(boundary_mesh);
Triangulation<dim-1,dim> boundary_mesh;
boundary_mesh.set_boundary (1, surface_description);
- GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh,
- surface_to_volume_mapping);
+ surface_to_volume_mapping
+ = GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh);
deallog << volume_mesh.n_active_cells () << std::endl;
deallog << boundary_mesh.n_active_cells () << std::endl;
save_mesh(boundary_mesh);
set<unsigned char> boundary_ids;
boundary_ids.insert(1);
- GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh,
- surface_to_volume_mapping,
- boundary_ids);
+ surface_to_volume_mapping
+ = GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh,
+ boundary_ids);
deallog << volume_mesh.n_active_cells () << std::endl;
deallog << boundary_mesh.n_active_cells () << std::endl;
Triangulation<dim> volume_mesh;
GridGenerator::hyper_cube(volume_mesh);
- GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh,
- surface_to_volume_mapping);
+ surface_to_volume_mapping
+ = GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh);
FE_Q <dim-1,dim> boundary_fe (1);
DoFHandler<dim-1,dim> boundary_dh(boundary_mesh);
Triangulation<spacedim> volume_mesh;
GridGenerator::hyper_cube(volume_mesh);
volume_mesh.refine_global(1);
- GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh,
- surface_to_volume_mapping);
+ surface_to_volume_mapping
+ = GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh);
boundary_mesh.refine_global(1);
for (Triangulation<dim,spacedim>::active_cell_iterator
Triangulation<spacedim> volume_mesh;
GridGenerator::hyper_cube(volume_mesh);
- GridTools::extract_boundary_mesh (volume_mesh, tria,
- surface_to_volume_mapping);
+ surface_to_volume_mapping
+ = GridTools::extract_boundary_mesh (volume_mesh, tria);
FE_Q<dim,spacedim> fe(2);
DoFHandler<dim,spacedim> dof_handler (tria);
std::set<unsigned char> boundary_ids;
boundary_ids.insert(0);
- GridTools::extract_boundary_mesh (volume_mesh, tria,
- surface_to_volume_mapping,
- boundary_ids);
+ surface_to_volume_mapping
+ = GridTools::extract_boundary_mesh (volume_mesh, tria,
+ boundary_ids);
FE_Q<dim-1,dim> fe (1);
DoFHandler<dim-1,dim> dh(tria);