Reference documentation for deal.II version 9.1.0-pre
Classes | Enumerations | Functions
GridTools Namespace Reference

Classes

class  Cache
 
struct  CellDataTransferBuffer
 
struct  PeriodicFacePair
 

Enumerations

Functions

template<typename DataType , typename MeshType >
void exchange_cell_data_to_ghosts (const MeshType &mesh, const std::function< boost::optional< DataType >(const typename MeshType::active_cell_iterator &)> &pack, const std::function< void(const typename MeshType::active_cell_iterator &, const DataType &)> &unpack)
 
template<class StreamType >
StreamType & operator<< (StreamType &s, const CacheUpdateFlags u)
 
CacheUpdateFlags operator| (const CacheUpdateFlags f1, const CacheUpdateFlags f2)
 
CacheUpdateFlags operator~ (const CacheUpdateFlags f1)
 
CacheUpdateFlagsoperator|= (CacheUpdateFlags &f1, const CacheUpdateFlags f2)
 
CacheUpdateFlags operator& (const CacheUpdateFlags f1, const CacheUpdateFlags f2)
 
CacheUpdateFlagsoperator&= (CacheUpdateFlags &f1, const CacheUpdateFlags f2)
 
Information about meshes and cells
template<int dim, int spacedim>
double diameter (const Triangulation< dim, spacedim > &tria)
 
template<int dim, int spacedim>
double volume (const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping=(StaticMappingQ1< dim, spacedim >::mapping))
 
template<int dim, int spacedim>
double minimal_cell_diameter (const Triangulation< dim, spacedim > &triangulation, const Mapping< dim, spacedim > &mapping=(StaticMappingQ1< dim, spacedim >::mapping))
 
template<int dim, int spacedim>
double maximal_cell_diameter (const Triangulation< dim, spacedim > &triangulation, const Mapping< dim, spacedim > &mapping=(StaticMappingQ1< dim, spacedim >::mapping))
 
template<int dim>
double cell_measure (const std::vector< Point< dim >> &all_vertices, const unsigned int(&vertex_indices)[GeometryInfo< dim >::vertices_per_cell])
 
template<int dim, typename T >
double cell_measure (const T &,...)
 
template<int dim, int spacedim>
BoundingBox< spacedim > compute_bounding_box (const Triangulation< dim, spacedim > &triangulation)
 
template<typename Iterator >
Point< Iterator::AccessorType::space_dimension > project_to_object (const Iterator &object, const Point< Iterator::AccessorType::space_dimension > &trial_point)
 
Functions supporting the creation of meshes
template<int dim, int spacedim>
void delete_unused_vertices (std::vector< Point< spacedim >> &vertices, std::vector< CellData< dim >> &cells, SubCellData &subcelldata)
 
template<int dim, int spacedim>
void delete_duplicated_vertices (std::vector< Point< spacedim >> &all_vertices, std::vector< CellData< dim >> &cells, SubCellData &subcelldata, std::vector< unsigned int > &considered_vertices, const double tol=1e-12)
 
Rotating, stretching and otherwise transforming meshes
template<int dim, typename Transformation , int spacedim>
void transform (const Transformation &transformation, Triangulation< dim, spacedim > &triangulation)
 
template<int dim, int spacedim>
void shift (const Tensor< 1, spacedim > &shift_vector, Triangulation< dim, spacedim > &triangulation)
 
void rotate (const double angle, Triangulation< 2 > &triangulation)
 
template<int dim>
void rotate (const double angle, const unsigned int axis, Triangulation< dim, 3 > &triangulation)
 
template<int dim>
void laplace_transform (const std::map< unsigned int, Point< dim >> &new_points, Triangulation< dim > &tria, const Function< dim, double > *coefficient=nullptr, const bool solve_for_absolute_positions=false)
 
template<int dim, int spacedim>
std::map< unsigned int, Point< spacedim > > get_all_vertices_at_boundary (const Triangulation< dim, spacedim > &tria)
 
template<int dim, int spacedim>
void scale (const double scaling_factor, Triangulation< dim, spacedim > &triangulation)
 
template<int dim, int spacedim>
void distort_random (const double factor, Triangulation< dim, spacedim > &triangulation, const bool keep_boundary=true)
 
template<int dim, int spacedim>
void remove_hanging_nodes (Triangulation< dim, spacedim > &tria, const bool isotropic=false, const unsigned int max_iterations=100)
 
template<int dim, int spacedim>
void remove_anisotropy (Triangulation< dim, spacedim > &tria, const double max_ratio=1.6180339887, const unsigned int max_iterations=5)
 
template<int dim, int spacedim>
void regularize_corner_cells (Triangulation< dim, spacedim > &tria, const double limit_angle_fraction=.75)
 
Finding cells and vertices of a triangulation
template<int dim, int spacedim>
return_type compute_point_locations (const Cache< dim, spacedim > &cache, const std::vector< Point< spacedim >> &points, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell_hint=typename Triangulation< dim, spacedim >::active_cell_iterator())
 
template<int dim, int spacedim>
return_type distributed_compute_point_locations (const GridTools::Cache< dim, spacedim > &cache, const std::vector< Point< spacedim >> &local_points, const std::vector< std::vector< BoundingBox< spacedim >>> &global_bboxes)
 
template<int dim, int spacedim>
std::map< unsigned int, Point< spacedim > > extract_used_vertices (const Triangulation< dim, spacedim > &container, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim, spacedim >::mapping)
 
template<int spacedim>
unsigned int find_closest_vertex (const std::map< unsigned int, Point< spacedim >> &vertices, const Point< spacedim > &p)
 
template<int dim, template< int, int > class MeshType, int spacedim>
unsigned int find_closest_vertex (const MeshType< dim, spacedim > &mesh, const Point< spacedim > &p, const std::vector< bool > &marked_vertices=std::vector< bool >())
 
template<int dim, template< int, int > class MeshType, int spacedim>
unsigned int find_closest_vertex (const Mapping< dim, spacedim > &mapping, const MeshType< dim, spacedim > &mesh, const Point< spacedim > &p, const std::vector< bool > &marked_vertices=std::vector< bool >())
 
template<int dim, template< int, int > class MeshType, int spacedim>
std::vector< typename MeshType< dim, spacedim >::active_cell_iterator > find_cells_adjacent_to_vertex (const MeshType< dim, spacedim > &container, const unsigned int vertex_index)
 
template<int dim, template< int, int > class MeshType, int spacedim>
MeshType< dim, spacedim >::active_cell_iterator find_active_cell_around_point (const MeshType< dim, spacedim > &mesh, const Point< spacedim > &p, const std::vector< bool > &marked_vertices=std::vector< bool >())
 
template<int dim, template< int, int > class MeshType, int spacedim>
std::pair< typename MeshType< dim, spacedim >::active_cell_iterator, Point< dim > > find_active_cell_around_point (const Mapping< dim, spacedim > &mapping, const MeshType< dim, spacedim > &mesh, const Point< spacedim > &p, const std::vector< bool > &marked_vertices=std::vector< bool >())
 
template<int dim, template< int, int > class MeshType, int spacedim>
std::pair< typename MeshType< dim, spacedim >::active_cell_iterator, Point< dim > > find_active_cell_around_point (const Mapping< dim, spacedim > &mapping, const MeshType< dim, spacedim > &mesh, const Point< spacedim > &p, const std::vector< std::set< typename MeshType< dim, spacedim >::active_cell_iterator >> &vertex_to_cell_map, const std::vector< std::vector< Tensor< 1, spacedim >>> &vertex_to_cell_centers, const typename MeshType< dim, spacedim >::active_cell_iterator &cell_hint=typename MeshType< dim, spacedim >::active_cell_iterator(), const std::vector< bool > &marked_vertices=std::vector< bool >())
 
template<int dim, int spacedim>
std::pair< typename hp::DoFHandler< dim, spacedim >::active_cell_iterator, Point< dim > > find_active_cell_around_point (const hp::MappingCollection< dim, spacedim > &mapping, const hp::DoFHandler< dim, spacedim > &mesh, const Point< spacedim > &p)
 
template<int dim, int spacedim>
std::pair< typename Triangulation< dim, spacedim >::active_cell_iterator, Point< dim > > find_active_cell_around_point (const Cache< dim, spacedim > &cache, const Point< spacedim > &p, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell_hint=typename Triangulation< dim, spacedim >::active_cell_iterator(), const std::vector< bool > &marked_vertices=std::vector< bool >())
 
template<int dim, template< int, int > class MeshType, int spacedim>
std::vector< std::pair< typename MeshType< dim, spacedim >::active_cell_iterator, Point< dim > > > find_all_active_cells_around_point (const Mapping< dim, spacedim > &mapping, const MeshType< dim, spacedim > &mesh, const Point< spacedim > &p, const double tolerance=1e-12, const std::vector< bool > &marked_vertices=std::vector< bool >())
 
template<class MeshType >
std::vector< typename MeshType::active_cell_iterator > get_active_child_cells (const typename MeshType::cell_iterator &cell)
 
template<class MeshType >
void get_active_neighbors (const typename MeshType::active_cell_iterator &cell, std::vector< typename MeshType::active_cell_iterator > &active_neighbors)
 
template<class MeshType >
std::vector< typename MeshType::active_cell_iterator > compute_active_cell_halo_layer (const MeshType &mesh, const std::function< bool(const typename MeshType::active_cell_iterator &)> &predicate)
 
template<class MeshType >
std::vector< typename MeshType::cell_iterator > compute_cell_halo_layer_on_level (const MeshType &mesh, const std::function< bool(const typename MeshType::cell_iterator &)> &predicate, const unsigned int level)
 
template<class MeshType >
std::vector< typename MeshType::active_cell_iterator > compute_ghost_cell_halo_layer (const MeshType &mesh)
 
template<class MeshType >
std::vector< typename MeshType::active_cell_iterator > compute_active_cell_layer_within_distance (const MeshType &mesh, const std::function< bool(const typename MeshType::active_cell_iterator &)> &predicate, const double layer_thickness)
 
template<class MeshType >
std::vector< typename MeshType::active_cell_iterator > compute_ghost_cell_layer_within_distance (const MeshType &mesh, const double layer_thickness)
 
template<class MeshType >
std::pair< Point< MeshType::space_dimension >, Point< MeshType::space_dimension > > compute_bounding_box (const MeshType &mesh, const std::function< bool(const typename MeshType::active_cell_iterator &)> &predicate)
 
template<class MeshType >
std::vector< BoundingBox< MeshType::space_dimension > > compute_mesh_predicate_bounding_box (const MeshType &mesh, const std::function< bool(const typename MeshType::active_cell_iterator &)> &predicate, const unsigned int refinement_level=0, const bool allow_merge=false, const unsigned int max_boxes=numbers::invalid_unsigned_int)
 
template<int spacedim>
return_type guess_point_owner (const std::vector< std::vector< BoundingBox< spacedim >>> &global_bboxes, const std::vector< Point< spacedim >> &points)
 
template<int dim, int spacedim>
std::vector< std::set< typename Triangulation< dim, spacedim >::active_cell_iterator > > vertex_to_cell_map (const Triangulation< dim, spacedim > &triangulation)
 
template<int dim, int spacedim>
std::vector< std::vector< Tensor< 1, spacedim > > > vertex_to_cell_centers_directions (const Triangulation< dim, spacedim > &mesh, const std::vector< std::set< typename Triangulation< dim, spacedim >::active_cell_iterator >> &vertex_to_cells)
 
template<int dim, int spacedim>
unsigned int find_closest_vertex_of_cell (const typename Triangulation< dim, spacedim >::active_cell_iterator &cell, const Point< spacedim > &position)
 
template<int dim, int spacedim>
std::map< unsigned int, types::global_vertex_indexcompute_local_to_global_vertex_index_map (const parallel::distributed::Triangulation< dim, spacedim > &triangulation)
 
template<int dim, int spacedim>
std::pair< unsigned int, double > get_longest_direction (typename Triangulation< dim, spacedim >::active_cell_iterator cell)
 
Partitions and subdomains of triangulations
template<int dim, int spacedim>
void get_face_connectivity_of_cells (const Triangulation< dim, spacedim > &triangulation, DynamicSparsityPattern &connectivity)
 
template<int dim, int spacedim>
void get_vertex_connectivity_of_cells (const Triangulation< dim, spacedim > &triangulation, DynamicSparsityPattern &connectivity)
 
template<int dim, int spacedim>
void get_vertex_connectivity_of_cells_on_level (const Triangulation< dim, spacedim > &triangulation, const unsigned int level, DynamicSparsityPattern &connectivity)
 
template<int dim, int spacedim>
void partition_triangulation (const unsigned int n_partitions, Triangulation< dim, spacedim > &triangulation, const SparsityTools::Partitioner partitioner=SparsityTools::Partitioner::metis)
 
template<int dim, int spacedim>
void partition_triangulation (const unsigned int n_partitions, const std::vector< unsigned int > &cell_weights, Triangulation< dim, spacedim > &triangulation, const SparsityTools::Partitioner partitioner=SparsityTools::Partitioner::metis)
 
template<int dim, int spacedim>
void partition_triangulation (const unsigned int n_partitions, const SparsityPattern &cell_connection_graph, Triangulation< dim, spacedim > &triangulation, const SparsityTools::Partitioner partitioner=SparsityTools::Partitioner::metis)
 
template<int dim, int spacedim>
void partition_triangulation (const unsigned int n_partitions, const std::vector< unsigned int > &cell_weights, const SparsityPattern &cell_connection_graph, Triangulation< dim, spacedim > &triangulation, const SparsityTools::Partitioner partitioner=SparsityTools::Partitioner::metis)
 
template<int dim, int spacedim>
void partition_triangulation_zorder (const unsigned int n_partitions, Triangulation< dim, spacedim > &triangulation)
 
template<int dim, int spacedim>
void partition_multigrid_levels (Triangulation< dim, spacedim > &triangulation)
 
template<int dim, int spacedim>
void get_subdomain_association (const Triangulation< dim, spacedim > &triangulation, std::vector< types::subdomain_id > &subdomain)
 
template<int dim, int spacedim>
unsigned int count_cells_with_subdomain_association (const Triangulation< dim, spacedim > &triangulation, const types::subdomain_id subdomain)
 
template<int dim, int spacedim>
std::vector< bool > get_locally_owned_vertices (const Triangulation< dim, spacedim > &triangulation)
 
Comparing different meshes
template<typename MeshType >
std::list< std::pair< typename MeshType::cell_iterator, typename MeshType::cell_iterator > > get_finest_common_cells (const MeshType &mesh_1, const MeshType &mesh_2)
 
template<int dim, int spacedim>
bool have_same_coarse_mesh (const Triangulation< dim, spacedim > &mesh_1, const Triangulation< dim, spacedim > &mesh_2)
 
template<typename MeshType >
bool have_same_coarse_mesh (const MeshType &mesh_1, const MeshType &mesh_2)
 
Dealing with distorted cells
template<int dim, int spacedim>
Triangulation< dim, spacedim >::DistortedCellList fix_up_distorted_child_cells (const typename Triangulation< dim, spacedim >::DistortedCellList &distorted_cells, Triangulation< dim, spacedim > &triangulation)
 
Extracting and creating patches of cells surrounding a single cell,

and creating triangulation out of them

template<class MeshType >
std::vector< typename MeshType::active_cell_iterator > get_patch_around_cell (const typename MeshType::active_cell_iterator &cell)
 
template<class Container >
std::vector< typename Container::cell_iterator > get_cells_at_coarsest_common_level (const std::vector< typename Container::active_cell_iterator > &patch_cells)
 
template<class Container >
void build_triangulation_from_patch (const std::vector< typename Container::active_cell_iterator > &patch, Triangulation< Container::dimension, Container::space_dimension > &local_triangulation, std::map< typename Triangulation< Container::dimension, Container::space_dimension >::active_cell_iterator, typename Container::active_cell_iterator > &patch_to_global_tria_map)
 
template<class DoFHandlerType >
std::map< types::global_dof_index, std::vector< typename DoFHandlerType::active_cell_iterator > > get_dof_to_support_patch_map (DoFHandlerType &dof_handler)
 
Dealing with periodic domains
template<typename FaceIterator >
bool orthogonal_equality (std::bitset< 3 > &orientation, const FaceIterator &face1, const FaceIterator &face2, const int direction, const Tensor< 1, FaceIterator::AccessorType::space_dimension > &offset=Tensor< 1, FaceIterator::AccessorType::space_dimension >(), const FullMatrix< double > &matrix=FullMatrix< double >())
 
template<typename FaceIterator >
bool orthogonal_equality (const FaceIterator &face1, const FaceIterator &face2, const int direction, const Tensor< 1, FaceIterator::AccessorType::space_dimension > &offset=Tensor< 1, FaceIterator::AccessorType::space_dimension >(), const FullMatrix< double > &matrix=FullMatrix< double >())
 
template<typename MeshType >
void collect_periodic_faces (const MeshType &mesh, const types::boundary_id b_id1, const types::boundary_id b_id2, const int direction, std::vector< PeriodicFacePair< typename MeshType::cell_iterator >> &matched_pairs, const Tensor< 1, MeshType::space_dimension > &offset=::Tensor< 1, MeshType::space_dimension >(), const FullMatrix< double > &matrix=FullMatrix< double >())
 
template<typename MeshType >
void collect_periodic_faces (const MeshType &mesh, const types::boundary_id b_id, const int direction, std::vector< PeriodicFacePair< typename MeshType::cell_iterator >> &matched_pairs, const ::Tensor< 1, MeshType::space_dimension > &offset=::Tensor< 1, MeshType::space_dimension >(), const FullMatrix< double > &matrix=FullMatrix< double >())
 
Dealing with boundary and manifold ids
template<int dim, int spacedim>
void copy_boundary_to_manifold_id (Triangulation< dim, spacedim > &tria, const bool reset_boundary_ids=false)
 
template<int dim, int spacedim>
void map_boundary_to_manifold_ids (const std::vector< types::boundary_id > &src_boundary_ids, const std::vector< types::manifold_id > &dst_manifold_ids, Triangulation< dim, spacedim > &tria, const std::vector< types::boundary_id > &reset_boundary_ids={})
 
template<int dim, int spacedim>
void copy_material_to_manifold_id (Triangulation< dim, spacedim > &tria, const bool compute_face_ids=false)
 
Exceptions
static::ExceptionBase & ExcInvalidNumberOfPartitions (int arg1)
 
static::ExceptionBase & ExcNonExistentSubdomain (int arg1)
 
static::ExceptionBase & ExcTriangulationHasBeenRefined ()
 
static::ExceptionBase & ExcScalingFactorNotPositive (double arg1)
 
template<int N>
static::ExceptionBase & ExcPointNotFoundInCoarseGrid (Point< N > arg1)
 
template<int N>
static::ExceptionBase & ExcPointNotFound (Point< N > arg1)
 
static::ExceptionBase & ExcVertexNotUsed (unsigned int arg1)
 

Detailed Description

This namespace is a collection of algorithms working on triangulations, such as shifting or rotating triangulations, but also finding a cell that contains a given point. See the descriptions of the individual functions for more information.

Enumeration Type Documentation

The enum type given to the Cache class to select what information to update.

You can select more than one flag by concatenation using the bitwise or operator|(CacheUpdateFlags,CacheUpdateFlags).

Author
Luca Heltai, 2017.
Enumerator
update_nothing 

Update Nothing.

update_vertex_to_cell_map 

Update vertex_to_cell_map, as returned by GridTools::vertex_to_cell_map().

update_vertex_to_cell_centers_directions 

Update vertex_to_cell_centers_directions, as returned by GridTools::vertex_to_cell_centers_directions()

update_vertex_kdtree 

Update a KDTree object, initialized with the vertices of the Triangulation.

update_used_vertices 

Update a mapping of used vertices.

update_all 

Update all objects.

Definition at line 35 of file grid_tools_cache_update_flags.h.

Function Documentation

template<int dim, int spacedim>
double GridTools::diameter ( const Triangulation< dim, spacedim > &  tria)

Return the diameter of a triangulation. The diameter is computed using only the vertices, i.e. if the diameter should be larger than the maximal distance between boundary vertices due to a higher order mapping, then this function will not catch this.

Definition at line 75 of file grid_tools.cc.

template<int dim, int spacedim>
double GridTools::volume ( const Triangulation< dim, spacedim > &  tria,
const Mapping< dim, spacedim > &  mapping = (StaticMappingQ1<dim, spacedim>::mapping) 
)

Compute the volume (i.e. the dim-dimensional measure) of the triangulation. We compute the measure using the integral \(\sum_K \int_K 1 \; dx\) where \(K\) are the cells of the given triangulation. The integral is approximated via quadrature for which we need the mapping argument.

If the triangulation is a dim-dimensional one embedded in a higher dimensional space of dimension spacedim, then the value returned is the dim-dimensional measure. For example, for a two-dimensional triangulation in three-dimensional space, the value returned is the area of the surface so described. (This obviously makes sense since the spacedim-dimensional measure of a dim-dimensional triangulation would always be zero if dim < spacedim.

This function also works for objects of type parallel::distributed::Triangulation, in which case the function is a collective operation.

Parameters
triaThe triangulation.
mappingAn optional argument used to denote the mapping that should be used when describing whether cells are bounded by straight or curved faces. The default is to use a \(Q_1\) mapping, which corresponds to straight lines bounding the cells.
Returns
The dim-dimensional measure of the domain described by the triangulation, as discussed above.

Definition at line 133 of file grid_tools.cc.

template<int dim, int spacedim>
double GridTools::minimal_cell_diameter ( const Triangulation< dim, spacedim > &  triangulation,
const Mapping< dim, spacedim > &  mapping = (StaticMappingQ1<dim, spacedim>::mapping) 
)

Return an approximation of the diameter of the smallest active cell of a triangulation. See step-24 for an example of use of this function.

Notice that, even if you pass a non-trivial mapping, the returned value is computed only using information on the vertices of the triangulation, possibly transformed by the mapping. While this is accurate most of the times, it may fail to give the correct result when the triangulation contains very distorted cells.

Definition at line 2888 of file grid_tools.cc.

template<int dim, int spacedim>
double GridTools::maximal_cell_diameter ( const Triangulation< dim, spacedim > &  triangulation,
const Mapping< dim, spacedim > &  mapping = (StaticMappingQ1<dim, spacedim>::mapping) 
)

Return an approximation of the diameter of the largest active cell of a triangulation.

Notice that, even if you pass a non-trivial mapping to this function, the returned value is computed only using information on the vertices of the triangulation, possibly transformed by the mapping. While this is accurate most of the times, it may fail to give the correct result when the triangulation contains very distorted cells.

Definition at line 2916 of file grid_tools.cc.

template<int dim>
double GridTools::cell_measure ( const std::vector< Point< dim >> &  all_vertices,
const unsigned int(&)  vertex_indices[GeometryInfo< dim >::vertices_per_cell] 
)

Given a list of vertices (typically obtained using Triangulation::get_vertices) as the first, and a list of vertex indices that characterize a single cell as the second argument, return the measure (area, volume) of this cell. If this is a real cell, then you can get the same result using cell->measure(), but this function also works for cells that do not exist except that you make it up by naming its vertices from the list.

template<int dim, typename T >
double GridTools::cell_measure ( const T &  ,
  ... 
)

A version of the last function that can accept input for nonzero codimension cases. This function only exists to aid generic programming and calling it will just raise an exception.

template<int dim, int spacedim>
BoundingBox< spacedim > GridTools::compute_bounding_box ( const Triangulation< dim, spacedim > &  triangulation)

Compute the smallest box containing the entire triangulation.

If the input triangulation is a parallel::distributed::Triangulation, then each processor will compute a bounding box enclosing all locally owned, ghost, and artificial cells. In the case of a domain without curved boundaries, these bounding boxes will all agree between processors because the union of the areas occupied by artificial and ghost cells equals the union of the areas occupied by the cells that other processors own. However, if the domain has curved boundaries, this is no longer the case. The bounding box returned may be appropriate for the current processor, but different from the bounding boxes computed on other processors.

Definition at line 420 of file grid_tools.cc.

template<typename Iterator >
Point<Iterator::AccessorType::space_dimension> GridTools::project_to_object ( const Iterator &  object,
const Point< Iterator::AccessorType::space_dimension > &  trial_point 
)

Return the point on the geometrical object object closest to the given point trial_point. For example, if object is a one-dimensional line or edge, then the returned point will be a point on the geodesic that connects the vertices as the manifold associated with the object sees it (i.e., the geometric line may be curved if it lives in a higher dimensional space). If the iterator points to a quadrilateral in a higher dimensional space, then the returned point lies within the convex hull of the vertices of the quad as seen by the associated manifold.

Note
This projection is usually not well-posed since there may be multiple points on the object that minimize the distance. The algorithm used in this function is robust (and the output is guaranteed to be on the given object) but may only provide a few correct digits if the object has high curvature. If your manifold supports it then the specialized function Manifold::project_to_manifold() may perform better.
Author
Luca Heltai, David Wells, 2017.
template<int dim, int spacedim>
void GridTools::delete_unused_vertices ( std::vector< Point< spacedim >> &  vertices,
std::vector< CellData< dim >> &  cells,
SubCellData subcelldata 
)

Remove vertices that are not referenced by any of the cells. This function is called by all GridIn::read_* functions to eliminate vertices that are listed in the input files but are not used by the cells in the input file. While these vertices should not be in the input from the beginning, they sometimes are, most often when some cells have been removed by hand without wanting to update the vertex lists, as they might be lengthy.

This function is called by all GridIn::read_* functions as the triangulation class requires them to be called with used vertices only. This is so, since the vertices are copied verbatim by that class, so we have to eliminate unused vertices beforehand.

Not implemented for the codimension one case.

Definition at line 434 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::delete_duplicated_vertices ( std::vector< Point< spacedim >> &  all_vertices,
std::vector< CellData< dim >> &  cells,
SubCellData subcelldata,
std::vector< unsigned int > &  considered_vertices,
const double  tol = 1e-12 
)

Remove vertices that are duplicated, due to the input of a structured grid, for example. If these vertices are not removed, the faces bounded by these vertices become part of the boundary, even if they are in the interior of the mesh.

This function is called by some GridIn::read_* functions. Only the vertices with indices in considered_vertices are tested for equality. This speeds up the algorithm, which is quadratic and thus quite slow to begin with. However, if you wish to consider all vertices, simply pass an empty vector. In that case, the function fills considered_vertices with all vertices.

Two vertices are considered equal if their difference in each coordinate direction is less than tol.

Definition at line 533 of file grid_tools.cc.

template<int dim, typename Transformation , int spacedim>
void GridTools::transform ( const Transformation &  transformation,
Triangulation< dim, spacedim > &  triangulation 
)

Transform the vertices of the given triangulation by applying the function object provided as first argument to all its vertices.

The transformation given as argument is used to transform each vertex. Its respective type has to offer a function-like syntax, i.e. the predicate is either an object of a type that has an operator(), or it is a pointer to the function. In either case, argument and return value have to be of type Point<spacedim>.

Note
If you are using a parallel::distributed::Triangulation you will have hanging nodes in your local Triangulation even if your "global" mesh has no hanging nodes. This will cause issues with wrong positioning of hanging nodes in ghost cells if you call the current functions: The vertices of all locally owned cells will be correct, but the vertices of some ghost cells may not. This means that computations like KellyErrorEstimator may give wrong answers. A safe approach is to use this function prior to any refinement in parallel, if that is possible, but not after you refine the mesh.

This function is used in the "Possibilities for extensions" section of step-38. It is also used in step-49 and step-53.

template<int dim, int spacedim>
void GridTools::shift ( const Tensor< 1, spacedim > &  shift_vector,
Triangulation< dim, spacedim > &  triangulation 
)

Shift each vertex of the triangulation by the given shift vector. This function uses the transform() function above, so the requirements on the triangulation stated there hold for this function as well.

Definition at line 698 of file grid_tools.cc.

void GridTools::rotate ( const double  angle,
Triangulation< 2 > &  triangulation 
)

Rotate all vertices of the given two-dimensional triangulation in counter-clockwise sense around the origin of the coordinate system by the given angle (given in radians, rather than degrees). This function uses the transform() function above, so the requirements on the triangulation stated there hold for this function as well.

Definition at line 707 of file grid_tools.cc.

template<int dim>
void GridTools::rotate ( const double  angle,
const unsigned int  axis,
Triangulation< dim, 3 > &  triangulation 
)

Rotate all vertices of the given triangulation in counter-clockwise direction around the axis with the given index. Otherwise like the function above.

Parameters
[in]angleAngle in radians to rotate the Triangulation by.
[in]axisIndex of the coordinate axis to rotate around, keeping that coordinate fixed (0=x axis, 1=y axis, 2=z axis).
[in,out]triangulationThe Triangulation object to rotate.
Note
Implemented for dim=1, 2, and 3.

Definition at line 714 of file grid_tools.cc.

template<int dim>
void GridTools::laplace_transform ( const std::map< unsigned int, Point< dim >> &  new_points,
Triangulation< dim > &  tria,
const Function< dim, double > *  coefficient = nullptr,
const bool  solve_for_absolute_positions = false 
)

Transform the given triangulation smoothly to a different domain where, typically, each of the vertices at the boundary of the triangulation is mapped to the corresponding points in the new_points map.

The unknown displacement field \(u_d(\mathbf x)\) in direction \(d\) is obtained from the minimization problem

\[ \min\, \int \frac{1}{2} c(\mathbf x) \mathbf \nabla u_d(\mathbf x) \cdot \mathbf \nabla u_d(\mathbf x) \,\rm d x \]

subject to prescribed constraints. The minimizer is obtained by solving the Laplace equation of the dim components of a displacement field that maps the current domain into one described by new_points . Linear finite elements with four Gaussian quadrature points in each direction are used. The difference between the vertex positions specified in new_points and their current value in tria therefore represents the prescribed values of this displacement field at the boundary of the domain, or more precisely at all of those locations for which new_points provides values (which may be at part of the boundary, or even in the interior of the domain). The function then evaluates this displacement field at each unconstrained vertex and uses it to place the mapped vertex where the displacement field locates it. Because the solution of the Laplace equation is smooth, this guarantees a smooth mapping from the old domain to the new one.

Parameters
[in]new_pointsThe locations where a subset of the existing vertices are to be placed. Typically, this would be a map from the vertex indices of all nodes on the boundary to their new locations, thus completely specifying the geometry of the mapped domain. However, it may also include interior points if necessary and it does not need to include all boundary vertices (although you then lose control over the exact shape of the mapped domain).
[in,out]triaThe Triangulation object. This object is changed in- place, i.e., the previous locations of vertices are overwritten.
[in]coefficientAn optional coefficient for the Laplace problem. Larger values make cells less prone to deformation (effectively increasing their stiffness). The coefficient is evaluated in the coordinate system of the old, undeformed configuration of the triangulation as input, i.e., before the transformation is applied. Should this function be provided, sensible results can only be expected if all coefficients are positive.
[in]solve_for_absolute_positionsIf set to true, the minimization problem is formulated with respect to the final vertex positions as opposed to their displacement. The two formulations are equivalent for the homogeneous problem (default value of coefficient), but they result in very different mesh motion otherwise. Since in most cases one will be using a non-constant coefficient in displacement formulation, the default value of this parameter is false.
Note
This function is not currently implemented for the 1d case.
template<int dim, int spacedim>
std::map< unsigned int, Point< spacedim > > GridTools::get_all_vertices_at_boundary ( const Triangulation< dim, spacedim > &  tria)

Return a std::map with all vertices of faces located in the boundary

Parameters
[in]triaThe Triangulation object.

Definition at line 875 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::scale ( const double  scaling_factor,
Triangulation< dim, spacedim > &  triangulation 
)

Scale the entire triangulation by the given factor. To preserve the orientation of the triangulation, the factor must be positive.

This function uses the transform() function above, so the requirements on the triangulation stated there hold for this function as well.

Definition at line 725 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::distort_random ( const double  factor,
Triangulation< dim, spacedim > &  triangulation,
const bool  keep_boundary = true 
)

Distort the given triangulation by randomly moving around all the vertices of the grid. The direction of movement of each vertex is random, while the length of the shift vector has a value of factor times the minimal length of the active edges adjacent to this vertex. Note that factor should obviously be well below 0.5.

If keep_boundary is set to true (which is the default), then boundary vertices are not moved.

Distort a triangulation in some random way.

Definition at line 908 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::remove_hanging_nodes ( Triangulation< dim, spacedim > &  tria,
const bool  isotropic = false,
const unsigned int  max_iterations = 100 
)

Remove hanging nodes from a grid. If the isotropic parameter is set to false (default) this function detects cells with hanging nodes and refines the neighbours in the direction that removes hanging nodes. If the isotropic parameter is set to true, the neighbours refinement is made in each directions. In order to remove all hanging nodes this procedure has to be repeated: this could require a large number of iterations. To avoid this a max number (max_iterations) of iteration is provided.

Consider the following grid:

remove_hanging_nodes-hanging.png

isotropic == false would return:

remove_hanging_nodes-aniso.png

isotropic == true would return:

remove_hanging_nodes-isotro.png
Parameters
[in,out]triaTriangulation to refine.
[in]isotropicIf true refine cells in each directions, otherwise (default value) refine the cell in the direction that removes hanging node.
[in]max_iterationsAt each step only closest cells to hanging nodes are refined. The code may require a lot of iterations to remove all hanging nodes. max_iterations is the maximum number of iteration allowed. If max_iterations == numbers::invalid_unsigned_int this function continues refining until there are no hanging nodes.
Note
In the case of parallel codes, this function should be combined with GridGenerator::flatten_triangulation.
Author
Mauro Bardelloni, Luca Heltai, Andrea Mola, 2016

Definition at line 3608 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::remove_anisotropy ( Triangulation< dim, spacedim > &  tria,
const double  max_ratio = 1.6180339887,
const unsigned int  max_iterations = 5 
)

Refine a mesh anisotropically such that the resulting mesh is composed by cells with maximum ratio between dimensions less than max_ratio. This procedure requires an algorithm that may not terminate. Consequently, it is possible to set a maximum number of iterations through the max_iterations parameter.

Starting from a cell like this:

remove_anisotropy-coarse.png

This function would return:

remove_anisotropy-refined.png
Parameters
[in,out]triaTriangulation to refine.
[in]max_ratioMaximum value allowed among the ratio between the dimensions of each cell.
[in]max_iterationsMaximum number of iterations allowed.
Note
In the case of parallel codes, this function should be combined with GridGenerator::flatten_triangulation and GridTools::remove_hanging_nodes.
Author
Mauro Bardelloni, Luca Heltai, Andrea Mola, 2016

Definition at line 3645 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::regularize_corner_cells ( Triangulation< dim, spacedim > &  tria,
const double  limit_angle_fraction = .75 
)

Analyze the boundary cells of a mesh, and if one cell is found at a corner position (with dim adjacent faces on the boundary), and its dim-dimensional angle fraction exceeds limit_angle_fraction, refine globally once, and replace the children of such cell with children where the corner is no longer offending the given angle fraction.

If no boundary cells exist with two adjacent faces on the boundary, then the triangulation is left untouched. If instead we do have cells with dim adjacent faces on the boundary, then the fraction between the dim-dimensional solid angle and dim*pi/2 is checked against the parameter limit_angle_fraction. If it is higher, the grid is refined once, and the children of the offending cell are replaced with some cells that instead respect the limit. After this process the triangulation is flattened, and all Manifold objects are restored as they were in the original triangulation.

An example is given by the following mesh, obtained by attaching a SphericalManifold to a mesh generated using GridGenerator::hyper_cube:

regularize_mesh_01.png

The four cells that were originally the corners of a square will give you some troubles during computations, as the jacobian of the transformation from the reference cell to those cells will go to zero, affecting the error constants of the finite element estimates.

Those cells have a corner with an angle that is very close to 180 degrees, i.e., an angle fraction very close to one.

The same code, adding a call to regularize_corner_cells:

generates a mesh that has a much better behaviour w.r.t. the jacobian of the Mapping:

regularize_mesh_02.png

This mesh is very similar to the one obtained by GridGenerator::hyper_ball. However, using GridTools::regularize_corner_cells one has the freedom to choose when to apply the regularization, i.e., one could in principle first refine a few times, and then call the regularize_corner_cells function:

This generates the following mesh:

regularize_mesh_03.png

The function is currently implemented only for dim = 2 and will throw an exception if called with dim = 3.

Parameters
[in,out]triaTriangulation to regularize.
[in]limit_angle_fractionMaximum ratio of angle or solid angle that is allowed for a corner element in the mesh.
Author
Luca Heltai, Martin Kronbichler, 2017

Definition at line 3678 of file grid_tools.cc.

template<int dim, int spacedim>
return_type GridTools::compute_point_locations ( const Cache< dim, spacedim > &  cache,
const std::vector< Point< spacedim >> &  points,
const typename Triangulation< dim, spacedim >::active_cell_iterator &  cell_hint = typename Triangulation<dim, spacedim>::active_cell_iterator() 
)

Given a Triangulation's cache and a list of points create the quadrature rules.

Parameters
[in]cacheThe triangulation's GridTools::Cache .
[in]pointsThe point's vector.
Returns
A tuple containing the following information:
  • Cells, is a vector of a vector cells of the all cells containing at least one of the points .
  • A vector qpoints of vector of points, containing in qpoints[i] the reference positions of all points that fall within the cell cells[i] .
  • A vector indices of vector of integers, containing the mapping between local numbering in qpoints, and global index in points

If points[a] and points[b] are the only two points that fall in cells[c], then qpoints[c][0] and qpoints[c][1] are the reference positions of points[a] and points[b] in cells[c], and indices[c][0] = a, indices[c][1] = b. The function Mapping::transform_unit_to_real(qpoints[c][0]) returns points[a].

The algorithm assumes it's easier to look for a point in the cell that was used previously. For this reason random points are, computationally speaking, the worst case scenario while points grouped by the cell to which they belong are the best case. Pre-sorting points, trying to minimize distances between them, might make the function extremely faster.

Note
The actual return type of this function, i.e., the type referenced above as return_type, is
std::tuple<
std::vector<
std::vector<std::vector<Point<dim>>>,
std::vector<std::vector<unsigned int>>>
The type is abbreviated in the online documentation to improve readability of this page.
Author
Giovanni Alzetta, 2017

Definition at line 3999 of file grid_tools.cc.

template<int dim, int spacedim>
return_type GridTools::distributed_compute_point_locations ( const GridTools::Cache< dim, spacedim > &  cache,
const std::vector< Point< spacedim >> &  local_points,
const std::vector< std::vector< BoundingBox< spacedim >>> &  global_bboxes 
)

Given a cache and a list of local_points for each process, find the points lying on the locally owned part of the mesh and compute the quadrature rules for them. Distributed compute point locations is a function similar to GridTools::compute_point_locations but working for parallel::Triangulation objects and, unlike its serial version, also for a distributed triangulation (see parallel::distributed::Triangulation).

Parameters
[in]cachea GridTools::Cache object
[in]local_pointsthe array of points owned by the current process. Every process can have a different array of points which can be empty and not contained within the locally owned part of the triangulation
[in]global_bboxesa vector of vectors of bounding boxes; it describes the locally owned part of the mesh for each process. The bounding boxes describing which part of the mesh is locally owned by process with rank rk are contained in global_bboxes[rk]. The local description can be obtained from GridTools::compute_mesh_predicate_bounding_box; then the global one can be obtained using either GridTools::exchange_local_bounding_boxes or Utilities::MPI::all_gather
Returns
A tuple containing the quadrature information

The elements of the output tuple are:

  • cells : a vector of cells of the all cells containing at least a point.
  • qpoints : a vector of vector of points; containing in qpoints[i] the reference positions of all points that fall within the cell cells[i] .
  • maps : a vector of vector of integers, containing the mapping between the numbering in qpoints (previous element of the tuple), and the vector of local points of the process owning the points.
  • points : a vector of vector of points. points[i][j] is the point in the real space corresponding. to qpoints[i][j] . Notice points are the points lying on the locally owned part of the mesh; thus these can be either copies of local_points or points received from other processes i.e. local_points for other processes
  • owners : a vector of vectors; owners[i][j] contains the rank of the process owning the point[i][j] (previous element of the tuple).

The function uses the triangulation's mpi communicator: for this reason it throws an assert error if the Triangulation is not derived from parallel::Triangulation .

In a serial execution the first three elements of the tuple are the same as in GridTools::compute_point_locations .

Note: this function is a collective operation.

Note
The actual return type of this function, i.e., the type referenced above as return_type, is
std::tuple<
std::vector<
typename Triangulation<dim, spacedim>::active_cell_iterator>,
std::vector<std::vector<Point<dim>>>,
std::vector<std::vector<unsigned int>>,
std::vector<std::vector<Point<spacedim>>>,
std::vector<std::vector<unsigned int>>>
The type is abbreviated in the online documentation to improve readability of this page.
Author
Giovanni Alzetta, 2017-2018

Definition at line 4506 of file grid_tools.cc.

template<int dim, int spacedim>
std::map< unsigned int, Point< spacedim > > GridTools::extract_used_vertices ( const Triangulation< dim, spacedim > &  container,
const Mapping< dim, spacedim > &  mapping = StaticMappingQ1<dim, spacedim>::mapping 
)

Return a map vertex index -> Point<spacedim> containing the used vertices of the given container. The key of the returned map (i.e., the first element of the pair above) is the global index in the triangulation, whereas the value of each pair is the physical location of the corresponding vertex. The used vertices are obtained by looping over all cells, and querying for each cell where its vertices are through the (optional) mapping argument.

In serial Triangulation objects and parallel::shared::Triangulation objects, the size of the returned map equals Triangulation::n_used_vertices() (not Triangulation::n_vertices()). Note that in parallel::distributed::Triangulation objects, only vertices in locally owned cells and ghost cells are returned, as for all other vertices their real location might not be known (e.g. for distributed computations using MappingQEulerian).

If you use the default mapping, the returned map satisfies the following equality:

const auto used_vertices = extract_used_vertices(tria);
auto all_vertices = tria.get_vertices();
for(const auto &id_and_v : used_vertices)
all_vertices[id_and_v.first] == id_and_v.second; // true

Notice that the above is not satisfied for mappings that change the location of vertices, like MappingQEulerian.

MeshType concept.

Parameters
containerThe container to extract vertices from.
mappingThe mapping to use to compute the points locations.
Author
Luca Heltai, 2017.

Definition at line 4740 of file grid_tools.cc.

template<int spacedim>
unsigned int GridTools::find_closest_vertex ( const std::map< unsigned int, Point< spacedim >> &  vertices,
const Point< spacedim > &  p 
)

Find and return the index of the closest vertex to a given point in the map of vertices passed as the first argument.

Parameters
verticesA map of index->vertex, as returned by GridTools::extract_used_vertices().
pThe target point.
Returns
The index of the vertex that is closest to the target point p.
Author
Luca Heltai, 2017.

Definition at line 4759 of file grid_tools.cc.

template<int dim, template< int, int > class MeshType, int spacedim>
unsigned int GridTools::find_closest_vertex ( const MeshType< dim, spacedim > &  mesh,
const Point< spacedim > &  p,
const std::vector< bool > &  marked_vertices = std::vector<bool>() 
)

Find and return the index of the used vertex (or marked vertex) in a given mesh that is located closest to a given point.

This function uses the locations of vertices as stored in the triangulation. This is usually sufficient, unless you are using a Mapping that moves the vertices around (for example, MappingQEulerian). In this case, you should call the function with the same name and with an additional Mapping argument.

Parameters
meshA variable of a type that satisfies the requirements of the MeshType concept.
pThe point for which we want to find the closest vertex.
marked_verticesAn array of bools indicating which vertices of mesh will be considered within the search as the potentially closest vertex. On receiving a non-empty marked_vertices, the function will only search among marked_vertices for the closest vertex. The size of this array should be equal to the value returned by Triangulation::n_vertices() for the triangulation underlying the given mesh (as opposed to the value returned by Triangulation::n_used_vertices()).
Returns
The index of the closest vertex found.
Author
Ralf B. Schulz, 2006

Definition at line 1156 of file grid_tools.cc.

template<int dim, template< int, int > class MeshType, int spacedim>
unsigned int GridTools::find_closest_vertex ( const Mapping< dim, spacedim > &  mapping,
const MeshType< dim, spacedim > &  mesh,
const Point< spacedim > &  p,
const std::vector< bool > &  marked_vertices = std::vector<bool>() 
)

Find and return the index of the used vertex (or marked vertex) in a given mesh that is located closest to a given point. Use the given mapping to compute the actual location of the vertices.

If the Mapping does not modify the position of the mesh vertices (like, for example, MappingQEulerian does), then this function is equivalent to the one with the same name, and without the mapping argument.

Parameters
mappingA mapping used to compute the vertex locations
meshA variable of a type that satisfies the requirements of the MeshType concept.
pThe point for which we want to find the closest vertex.
marked_verticesAn array of bools indicating which vertices of mesh will be considered within the search as the potentially closest vertex. On receiving a non-empty marked_vertices, the function will only search among marked_vertices for the closest vertex. The size of this array should be equal to the value returned by Triangulation::n_vertices() for the triangulation underlying the given mesh (as opposed to the value returned by Triangulation::n_used_vertices()).
Returns
The index of the closest vertex found.
Author
Luca Heltai, 2017

Definition at line 1232 of file grid_tools.cc.

template<int dim, template< int, int > class MeshType, int spacedim>
std::vector< typename MeshType< dim, spacedim >::active_cell_iterator > GridTools::find_cells_adjacent_to_vertex ( const MeshType< dim, spacedim > &  container,
const unsigned int  vertex_index 
)

Find and return a vector of iterators to active cells that surround a given vertex with index vertex_index.

For locally refined grids, the vertex itself might not be a vertex of all adjacent cells that are returned. However, it will always be either a vertex of a cell or be a hanging node located on a face or an edge of it.

Parameters
containerA variable of a type that satisfies the requirements of the MeshType concept.
vertex_indexThe index of the vertex for which we try to find adjacent cells.
Returns
A vector of cells that lie adjacent to the given vertex.
Note
If the point requested does not lie in any of the cells of the mesh given, then this function throws an exception of type GridTools::ExcPointNotFound. You can catch this exception and decide what to do in that case.
It isn't entirely clear at this time whether the function does the right thing with anisotropically refined meshes. It needs to be checked for this case.

Definition at line 1298 of file grid_tools.cc.

template<int dim, template< int, int > class MeshType, int spacedim>
MeshType< dim, spacedim >::active_cell_iterator GridTools::find_active_cell_around_point ( const MeshType< dim, spacedim > &  mesh,
const Point< spacedim > &  p,
const std::vector< bool > &  marked_vertices = std::vector<bool>() 
)

Find and return an iterator to the active cell that surrounds a given point. This function simply calls the following one with a MappingQ1 for the mapping argument. See the following function for a more thorough discussion.

Parameters
meshA variable of a type that satisfies the requirements of the MeshType concept.
pThe point for which we want to find the surrounding cell.
marked_verticesAn array of bools indicating whether an entry in the vertex array should be considered (and the others must be ignored) as the potentially closest vertex to the specified point. On specifying a non-default marked_vertices, find_closest_vertex() would only search among marked_vertices for the closest vertex. The size of this array should be equal to n_vertices() of the triangulation (as opposed to n_used_vertices() ).
Returns
An iterator into the mesh that points to the surrounding cell.
Note
If the point requested does not lie in any of the cells of the mesh given, then this function throws an exception of type GridTools::ExcPointNotFound. You can catch this exception and decide what to do in that case.

Definition at line 549 of file grid_tools_dof_handlers.cc.

template<int dim, template< int, int > class MeshType, int spacedim>
std::pair< typename MeshType< dim, spacedim >::active_cell_iterator, Point< dim > > GridTools::find_active_cell_around_point ( const Mapping< dim, spacedim > &  mapping,
const MeshType< dim, spacedim > &  mesh,
const Point< spacedim > &  p,
const std::vector< bool > &  marked_vertices = std::vector<bool>() 
)

Find and return an iterator to the active cell that surrounds a given point p.

The algorithm used in this function proceeds by first looking for the vertex located closest to the given point, see GridTools::find_closest_vertex(). Secondly, all adjacent cells to this vertex are found in the mesh, see GridTools::find_cells_adjacent_to_vertex(). Lastly, for each of these cells, the function tests whether the point is inside. This check is performed using the given mapping argument to determine whether cells have straight or curved boundaries, and if the latter then how exactly they are curved.

If a point lies on the boundary of two or more cells, then the algorithm tries to identify the cell that is of highest refinement level.

Parameters
mappingThe mapping used to determine whether the given point is inside a given cell.
meshA variable of a type that satisfies the requirements of the MeshType concept.
pThe point for which we want to find the surrounding cell.
marked_verticesAn array of bools indicating whether an entry in the vertex array should be considered (and the others must be ignored) as the potentially closest vertex to the specified point. On specifying a non-default marked_vertices, find_closest_vertex() would only search among marked_vertices for the closest vertex. The size of this array should be equal to n_vertices() of the triangulation (as opposed to n_used_vertices() ).
Returns
An pair of an iterators into the mesh that points to the surrounding cell, and of the coordinates of that point inside the cell in the reference coordinates of that cell. This local position might be located slightly outside an actual unit cell, due to numerical roundoff. Therefore, the point returned by this function should be projected onto the unit cell, using GeometryInfo::project_to_unit_cell(). This is not automatically performed by the algorithm.
Note
When marked_vertices is specified the function should always be called inside a try block to catch the exception that the function might throw in the case it couldn't find an active cell surrounding the point. The motivation of using marked_vertices is to cut down the search space of vertices if one has a priori knowledge of a collection of vertices that the point of interest may be close to. For instance, in the case when a parallel::shared::Triangulation is employed and we are looking for a point that we know is inside the locally owned part of the mesh, then it would make sense to pass an array for marked_vertices that flags only the vertices of all locally owned active cells. If, however, the function throws an exception, then that would imply that the point lies outside locally owned active cells.
If the point requested does not lie in any of the cells of the mesh given, then this function throws an exception of type GridTools::ExcPointNotFound. You can catch this exception and decide what to do in that case.
When applied to a triangulation or DoF handler object based on a parallel::distributed::Triangulation object, the cell returned may in fact be a ghost or artificial cell (see GlossArtificialCell and GlossGhostCell). If so, many of the operations one may want to do on this cell (e.g., evaluating the solution) may not be possible and you will have to decide what to do in that case.
Floating point arithmetic implies that a point will, in general, never lie exactly on an edge or a face. It may, however, lie on a vertex of a cell. In either case, it is not predictable which of the cells adjacent to a vertex or an edge/face this function returns when given a point that lies on a vertex or within floating point precision of an edge or face. Consequently, algorithms that call this function need to take into account that the returned cell will only contain the point approximately (to within round-off error) and that these cells may also be ghost cells or artificial cells if the triangulation is a parallel one. The latter may even be true if the given point is in fact a vertex of a locally owned cell: the returned cell may still be a ghost cell that happens to share this vertex with a locally owned one. The reason for this behavior is that it is the only way to guarantee that all processors that participate in a parallel triangulation will agree which cell contains a point. In other words, two processors that own two cells that come together at one vertex will return the same cell when called with this vertex. One of them will then return a locally owned cell and the other one a ghost cell.

Definition at line 568 of file grid_tools_dof_handlers.cc.

template<int dim, template< int, int > class MeshType, int spacedim>
std::pair< typename MeshType< dim, spacedim >::active_cell_iterator, Point< dim > > GridTools::find_active_cell_around_point ( const Mapping< dim, spacedim > &  mapping,
const MeshType< dim, spacedim > &  mesh,
const Point< spacedim > &  p,
const std::vector< std::set< typename MeshType< dim, spacedim >::active_cell_iterator >> &  vertex_to_cell_map,
const std::vector< std::vector< Tensor< 1, spacedim >>> &  vertex_to_cell_centers,
const typename MeshType< dim, spacedim >::active_cell_iterator &  cell_hint = typename MeshType<dim, spacedim>::active_cell_iterator(),
const std::vector< bool > &  marked_vertices = std::vector<bool>() 
)

A version of the previous function that exploits an already existing map between vertices and cells, constructed using the function GridTools::vertex_to_cell_map, a map of vertex_to_cell_centers, obtained through GridTools::vertex_to_cell_centers_directions, and a guess cell_hint.

Author
Luca Heltai, Rene Gassmoeller, 2017

Definition at line 1499 of file grid_tools.cc.

template<int dim, int spacedim>
std::pair< typename hp::DoFHandler< dim, spacedim >::active_cell_iterator, Point< dim > > GridTools::find_active_cell_around_point ( const hp::MappingCollection< dim, spacedim > &  mapping,
const hp::DoFHandler< dim, spacedim > &  mesh,
const Point< spacedim > &  p 
)

A version of the previous function where we use that mapping on a given cell that corresponds to the active finite element index of that cell. This is obviously only useful for hp problems, since the active finite element index for all other DoF handlers is always zero.

Note
If the point requested does not lie in any of the cells of the mesh given, then this function throws an exception of type GridTools::ExcPointNotFound. You can catch this exception and decide what to do in that case.
When applied to a triangulation or DoF handler object based on a parallel::distributed::Triangulation object, the cell returned may in fact be a ghost or artificial cell (see GlossArtificialCell and GlossGhostCell). If so, many of the operations one may want to do on this cell (e.g., evaluating the solution) may not be possible and you will have to decide what to do in that case.

Definition at line 1263 of file grid_tools_dof_handlers.cc.

template<int dim, int spacedim>
std::pair< typename Triangulation< dim, spacedim >::active_cell_iterator, Point< dim > > GridTools::find_active_cell_around_point ( const Cache< dim, spacedim > &  cache,
const Point< spacedim > &  p,
const typename Triangulation< dim, spacedim >::active_cell_iterator &  cell_hint = typename Triangulation<dim, spacedim>::active_cell_iterator(),
const std::vector< bool > &  marked_vertices = std::vector<bool>() 
)

A version of the previous function that exploits an already existing GridTools::Cache<dim,spacedim> object.

Author
Luca Heltai, 2017

Definition at line 4776 of file grid_tools.cc.

template<int dim, template< int, int > class MeshType, int spacedim>
std::vector< std::pair< typename MeshType< dim, spacedim >::active_cell_iterator, Point< dim > > > GridTools::find_all_active_cells_around_point ( const Mapping< dim, spacedim > &  mapping,
const MeshType< dim, spacedim > &  mesh,
const Point< spacedim > &  p,
const double  tolerance = 1e-12,
const std::vector< bool > &  marked_vertices = std::vector<bool>() 
)

A variant of the previous find_active_cell_around_point() function that, instead of returning only the first matching cell, identifies all cells around a point for a given tolerance level tolerance in terms of unit coordinates. More precisely, whenever the point returned by find_active_cell_around_point() is within the given tolerance from the surface of the unit cell, all corresponding neighbors are also identified, including the location of the point in unit coordinates on any of these cells.

This function is useful e.g. for discontinuous function spaces where, for the case the given point p coincides with a vertex or an edge, several cells might hold independent values of the solution that get combined in some way in a user code.

Author
Niklas Fehn, Martin Kronbichler, 2018

Definition at line 589 of file grid_tools_dof_handlers.cc.

template<class MeshType >
std::vector<typename MeshType::active_cell_iterator> GridTools::get_active_child_cells ( const typename MeshType::cell_iterator &  cell)

Return a list of all descendants of the given cell that are active. For example, if the current cell is once refined but none of its children are any further refined, then the returned list will contain all its children.

If the current cell is already active, then the returned list is empty (because the cell has no children that may be active).

Template Parameters
MeshTypeA type that satisfies the requirements of the MeshType concept.
Parameters
cellAn iterator pointing to a cell of the mesh.
Returns
A list of active descendants of the given cell
Note
Since in C++ the MeshType template argument can not be deduced from a function call, you will have to specify it after the function name, as for example in
GridTools::get_active_child_cells<DoFHandler<dim> > (cell)
template<class MeshType >
void GridTools::get_active_neighbors ( const typename MeshType::active_cell_iterator &  cell,
std::vector< typename MeshType::active_cell_iterator > &  active_neighbors 
)

Extract the active cells around a given cell cell and return them in the vector active_neighbors.

Template Parameters
MeshTypeA type that satisfies the requirements of the MeshType concept.
Parameters
[in]cellAn iterator pointing to a cell of the mesh.
[out]active_neighborsA list of active descendants of the given cell
template<class MeshType >
std::vector< typename MeshType::active_cell_iterator > GridTools::compute_active_cell_halo_layer ( const MeshType &  mesh,
const std::function< bool(const typename MeshType::active_cell_iterator &)> &  predicate 
)

Extract and return the active cell layer around a subdomain (set of active cells) in the mesh (i.e. those that share a common set of vertices with the subdomain but are not a part of it). Here, the "subdomain" consists of exactly all of those cells for which the predicate returns true.

An example of a custom predicate is one that checks for a given material id

template <int dim>
bool
pred_mat_id(const typename Triangulation<dim>::active_cell_iterator & cell)
{
return cell->material_id() == 1;
}

and we can then extract the layer of cells around this material with the following call:

Predicates that are frequently useful can be found in namespace IteratorFilters. For example, it is possible to extract a layer of cells around all of those cells with a given material id,

or around all cells with one of a set of active FE indices for an hp::DoFHandler

Note that in the last two examples we ensure that the predicate returns true only for locally owned cells. This means that the halo layer will not contain any artificial cells.

Template Parameters
MeshTypeA type that satisfies the requirements of the MeshType concept.
Parameters
[in]meshA mesh (i.e. objects of type Triangulation, DoFHandler, or hp::DoFHandler).
[in]predicateA function (or object of a type with an operator()) defining the subdomain around which the halo layer is to be extracted. It is a function that takes in an active cell and returns a boolean.
Returns
A list of active cells sharing at least one common vertex with the predicated subdomain.
Author
Jean-Paul Pelteret, Denis Davydov, Wolfgang Bangerth, 2015

Definition at line 737 of file grid_tools_dof_handlers.cc.

template<class MeshType >
std::vector< typename MeshType::cell_iterator > GridTools::compute_cell_halo_layer_on_level ( const MeshType &  mesh,
const std::function< bool(const typename MeshType::cell_iterator &)> &  predicate,
const unsigned int  level 
)

Extract and return the cell layer around a subdomain (set of cells) on a specified level of the mesh (i.e. those cells on that level that share a common set of vertices with the subdomain but are not a part of it). Here, the "subdomain" consists of exactly all of those cells for which the predicate returns true.

Definition at line 782 of file grid_tools_dof_handlers.cc.

template<class MeshType >
std::vector< typename MeshType::active_cell_iterator > GridTools::compute_ghost_cell_halo_layer ( const MeshType &  mesh)

Extract and return ghost cells which are the active cell layer around all locally owned cells. This is most relevant for parallel::shared::Triangulation where it will return a subset of all ghost cells on a processor, but for parallel::distributed::Triangulation this will return all the ghost cells.

Template Parameters
MeshTypeA type that satisfies the requirements of the MeshType concept.
Parameters
[in]meshA mesh (i.e. objects of type Triangulation, DoFHandler, or hp::DoFHandler).
Returns
A list of ghost cells
Author
Jean-Paul Pelteret, Denis Davydov, Wolfgang Bangerth, 2015

Definition at line 867 of file grid_tools_dof_handlers.cc.

template<class MeshType >
std::vector< typename MeshType::active_cell_iterator > GridTools::compute_active_cell_layer_within_distance ( const MeshType &  mesh,
const std::function< bool(const typename MeshType::active_cell_iterator &)> &  predicate,
const double  layer_thickness 
)

Extract and return the set of active cells within a geometric distance of layer_thickness around a subdomain (set of active cells) in the mesh. Here, the "subdomain" consists of exactly all of those cells for which the predicate returns true.

The function first computes the cells that form the 'surface' of the subdomain that consists of all of the active cells for which the predicate is true. Using compute_bounding_box(), a bounding box is computed for this subdomain and extended by layer_thickness. These cells are called interior subdomain boundary cells. The active cells with all of their vertices outside the extended bounding box are ignored. The cells that are inside the extended bounding box are then checked for their proximity to the interior subdomain boundary cells. This implies checking the distance between a pair of arbitrarily oriented cells, which is not trivial in general. To simplify this, the algorithm checks the distance between the two enclosing spheres of the cells. This will definitely result in slightly more cells being marked but also greatly simplifies the arithmetic complexity of the algorithm.

active_cell_layer_within_distance.png

The image shows a mesh generated by subdivided_hyper_rectangle(). The cells are marked using three different colors. If the grey colored cells in the image are the cells for which the predicate is true, then the function compute_active_cell_layer_within_distance() will return a set of cell iterators corresponding to the cells colored in red. The red colored cells are the active cells that are within a given distance to the grey colored cells.

Template Parameters
MeshTypeA type that satisfies the requirements of the MeshType concept.
Parameters
meshA mesh (i.e. objects of type Triangulation, DoFHandler, or hp::DoFHandler).
predicateA function (or object of a type with an operator()) defining the subdomain around which the halo layer is to be extracted. It is a function that takes in an active cell and returns a boolean.
layer_thicknessspecifies the geometric distance within which the function searches for active cells from the predicate domain. If the minimal distance between the enclosing sphere of the an active cell and the enclosing sphere of any of the cells for which the predicate returns true is less than layer_thickness, then the active cell is an active_cell_within_distance.
Returns
A list of active cells within a given geometric distance layer_thickness from the set of active cells for which the predicate returns true.

See compute_active_cell_halo_layer().

Author
Vishal Boddu, Denis Davydov, 2017

Definition at line 889 of file grid_tools_dof_handlers.cc.

template<class MeshType >
std::vector< typename MeshType::active_cell_iterator > GridTools::compute_ghost_cell_layer_within_distance ( const MeshType &  mesh,
const double  layer_thickness 
)

Extract and return a set of ghost cells which are within a layer_thickness around all locally owned cells. This is most relevant for parallel::shared::Triangulation where it will return a subset of all ghost cells on a process, but for parallel::distributed::Triangulation this will return all the ghost cells. All the cells for the parallel::shared::Triangulation class that are not owned by the current processor can be considered as ghost cells; in particular, they do not only form a single layer of cells around the locally owned ones.

Template Parameters
MeshTypeA type that satisfies the requirements of the MeshType concept.
Parameters
meshA mesh (i.e. objects of type Triangulation, DoFHandler, or hp::DoFHandler).
layer_thicknessspecifies the geometric distance within which the function searches for active cells from the locally owned cells.
Returns
A subset of ghost cells within a given geometric distance of layer_thickness from the locally owned cells of a current process.

Also see compute_ghost_cell_halo_layer() and compute_active_cell_layer_within_distance().

Author
Vishal Boddu, Denis Davydov, 2017

Definition at line 1044 of file grid_tools_dof_handlers.cc.

template<class MeshType >
std::pair< Point< MeshType::space_dimension >, Point< MeshType::space_dimension > > GridTools::compute_bounding_box ( const MeshType &  mesh,
const std::function< bool(const typename MeshType::active_cell_iterator &)> &  predicate 
)

Compute and return a bounding box, defined through a pair of points bottom left and top right, that surrounds a subdomain of the mesh. Here, the "subdomain" consists of exactly all of those active cells for which the predicate returns true.

For a description of how predicate works, see compute_active_cell_halo_layer().

Note
This function was written before the BoundingBox class was invented. Consequently, it returns a pair of points, rather than a BoundingBox object as one may expect. However, BoundingBox has a conversion constructor from pairs of points, so the result of this function can still be assigned to a BoundingBox object.

Definition at line 1080 of file grid_tools_dof_handlers.cc.

template<class MeshType >
std::vector< BoundingBox< MeshType::space_dimension > > GridTools::compute_mesh_predicate_bounding_box ( const MeshType &  mesh,
const std::function< bool(const typename MeshType::active_cell_iterator &)> &  predicate,
const unsigned int  refinement_level = 0,
const bool  allow_merge = false,
const unsigned int  max_boxes = numbers::invalid_unsigned_int 
)

Compute a collection of bounding boxes so that all active cells for which the given predicate is true, are completely enclosed in at least one of the bounding boxes. Notice the cover is only guaranteed to contain all these active cells but it's not necessarily exact i.e. it can include a bigger area than their union.

For each cell at a given refinement level containing active cells for which predicate is true, the function creates a bounding box of its children for which predicate is true.

This results in a cover of all active cells for which predicate is true; the parameters allow_merge and max_boxes are used to reduce the number of cells at a computational cost and covering a bigger n-dimensional volume.

The parameters to control the algorithm are:

  • predicate : the property of the cells to enclose e.g. IteratorFilters::LocallyOwnedCell . The predicate is tested only on active cells.
  • refinement_level : it defines the level at which the initial bounding box are created. The refinement should be set to a coarse refinement level. A bounding box is created for each active cell at coarser level than refinement_level; if refinement_level is higher than the number of levels of the triangulation an exception is thrown.
  • allow_merge : This flag allows for box merging and, by default, is false. The algorithm has a cost of O(N^2) where N is the number of the bounding boxes created from the refinement level; for this reason, if the flag is set to true, make sure to choose wisely a coarse enough refinement_level.
  • max_boxes : the maximum number of bounding boxes to compute. If more are created the smaller ones are merged with neighbors. By default after merging the boxes which can be expressed as a single one no more boxes are merged. See the BoundingBox::get_neighbor_type () function for details. Notice only neighboring cells are merged (see the get_neighbor_type function in bounding box class): if the target number of bounding boxes max_boxes can't be reached by merging neighbors an exception is thrown

The following image describes an example of the algorithm with refinement_level = 2, allow_merge = true and max_boxes = 1. The cells with the property predicate are in red, the area of a bounding box is slightly orange.

bounding_box_predicate.png
  • 1. In black we can see the cells of the current level.
  • 2. For each cell containing the red area a bounding box is created: by default these are returned.
  • 3. Because allow_merge = true the number of bounding boxes is reduced while not changing the cover. If max_boxes was left as default or bigger than 1 these two boxes would be returned.
  • 4. Because max_boxes = 1 the smallest bounding box is merged to the bigger one. Notice it is important to choose the parameters wisely. For instance, allow_merge = false and refinement_level = 1 returns the very same bounding box but with a fraction of the computational cost.

This function does not take into account the curvature of cells and thus it is not suited for handling curved geometry: the mapping is assumed to be linear.

Definition at line 1720 of file grid_tools.cc.

template<int spacedim>
return_type GridTools::guess_point_owner ( const std::vector< std::vector< BoundingBox< spacedim >>> &  global_bboxes,
const std::vector< Point< spacedim >> &  points 
)

Given an array of points, use the global bounding box description obtained using GridTools::compute_mesh_predicate_bounding_box to guess, for each of them, which process might own it.

Parameters
[in]global_bboxesVector of bounding boxes describing the portion of mesh with a property for each process.
[in]pointsArray of points to test.
Returns
A tuple containing the following information:
  • A vector indicized with ranks of processes. For each rank it contains a vector of the indices of points it might own.
  • A map from the index unsigned int of the point in points to the rank of the owner.
  • A map from the index unsigned int of the point in points to the ranks of the guessed owners.
Note
The actual return type of this function, i.e., the type referenced above as return_type, is
std::tuple<std::vector<std::vector<unsigned int>>,
std::map< unsigned int, unsigned int>,
std::map< unsigned int, std::vector<unsigned int>>>
The type is abbreviated in the online documentation to improve readability of this page.
Author
Giovanni Alzetta, 2017

Definition at line 1867 of file grid_tools.cc.

template<int dim, int spacedim>
std::vector< std::set< typename Triangulation< dim, spacedim >::active_cell_iterator > > GridTools::vertex_to_cell_map ( const Triangulation< dim, spacedim > &  triangulation)

Return the adjacent cells of all the vertices. If a vertex is also a hanging node, the associated coarse cell is also returned. The vertices are ordered by the vertex index. This is the number returned by the function cell->vertex_index(). Notice that only the indices marked in the array returned by Triangulation<dim,spacedim>::get_used_vertices() are used.

Definition at line 1919 of file grid_tools.cc.

template<int dim, int spacedim>
std::vector< std::vector< Tensor< 1, spacedim > > > GridTools::vertex_to_cell_centers_directions ( const Triangulation< dim, spacedim > &  mesh,
const std::vector< std::set< typename Triangulation< dim, spacedim >::active_cell_iterator >> &  vertex_to_cells 
)

Return a vector of normalized tensors for each vertex-cell combination of the output of GridTools::vertex_to_cell_map() (which is expected as input parameter for this function). Each tensor represents a geometric vector from the vertex to the respective cell center.

An assertion will be thrown if the size of the input vector is not equal to the number of vertices of the triangulation.

result[v][c] is a unit Tensor for vertex index v, indicating the direction of the center of the c-th cell with respect to the vertex v.

Author
Rene Gassmoeller, Luca Heltai, 2017.

Definition at line 1436 of file grid_tools.cc.

template<int dim, int spacedim>
unsigned int GridTools::find_closest_vertex_of_cell ( const typename Triangulation< dim, spacedim >::active_cell_iterator &  cell,
const Point< spacedim > &  position 
)

Return the local vertex index of cell cell that is closest to the given location position.

Author
Rene Gassmoeller, Luca Heltai, 2017.

Definition at line 1637 of file grid_tools.cc.

template<int dim, int spacedim>
std::map< unsigned int, types::global_vertex_index > GridTools::compute_local_to_global_vertex_index_map ( const parallel::distributed::Triangulation< dim, spacedim > &  triangulation)

Compute a globally unique index for each vertex and hanging node associated with a locally owned active cell. The vertices of a ghost cell that are hanging nodes of a locally owned cells have a global index. However, the other vertices of the cells that do not touch an active cell do not have a global index on this processor.

The key of the map is the local index of the vertex and the value is the global index. The indices need to be recomputed after refinement or coarsening and may be different.

Definition at line 1970 of file grid_tools.cc.

template<int dim, int spacedim>
std::pair< unsigned int, double > GridTools::get_longest_direction ( typename Triangulation< dim, spacedim >::active_cell_iterator  cell)

Return the highest value among ratios between extents in each of the coordinate directions of a cell. Moreover, return the dimension relative to the highest elongation.

Parameters
[in]cellan iterator pointing to the cell.
Returns
A std::pair<unsigned int, double> such that the first value is the dimension of the highest elongation and the second value is the ratio among the dimensions of the cell.
Author
Mauro Bardelloni, Luca Heltai, Andrea Mola, 2016

Definition at line 3576 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::get_face_connectivity_of_cells ( const Triangulation< dim, spacedim > &  triangulation,
DynamicSparsityPattern connectivity 
)

Produce a sparsity pattern in which nonzero entries indicate that two cells are connected via a common face. The diagonal entries of the sparsity pattern are also set.

The rows and columns refer to the cells as they are traversed in their natural order using cell iterators.

Definition at line 2327 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::get_vertex_connectivity_of_cells ( const Triangulation< dim, spacedim > &  triangulation,
DynamicSparsityPattern connectivity 
)

Produce a sparsity pattern in which nonzero entries indicate that two cells are connected via a common vertex. The diagonal entries of the sparsity pattern are also set.

The rows and columns refer to the cells as they are traversed in their natural order using cell iterators.

Definition at line 2383 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::get_vertex_connectivity_of_cells_on_level ( const Triangulation< dim, spacedim > &  triangulation,
const unsigned int  level,
DynamicSparsityPattern connectivity 
)

Produce a sparsity pattern for a given level mesh in which nonzero entries indicate that two cells are connected via a common vertex. The diagonal entries of the sparsity pattern are also set.

The rows and columns refer to the cells as they are traversed in their natural order using cell iterators.

Definition at line 2418 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::partition_triangulation ( const unsigned int  n_partitions,
Triangulation< dim, spacedim > &  triangulation,
const SparsityTools::Partitioner  partitioner = SparsityTools::Partitioner::metis 
)

Use graph partitioner to partition the active cells making up the entire domain. After calling this function, the subdomain ids of all active cells will have values between zero and n_partitions-1. You can access the subdomain id of a cell by using cell->subdomain_id().

Use the third argument to select between partitioning algorithms provided by METIS or ZOLTAN. METIS is the default partitioner.

If deal.II was not installed with ZOLTAN or METIS, this function will generate an error when corresponding partition method is chosen, unless n_partitions is one. I.e., you can write a program so that it runs in the single-processor single-partition case without packages installed, and only requires them installed when multiple partitions are required.

Note
If the cell_weight signal has been attached to the triangulation, then this will be used and passed to the partitioner.

Definition at line 2454 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::partition_triangulation ( const unsigned int  n_partitions,
const std::vector< unsigned int > &  cell_weights,
Triangulation< dim, spacedim > &  triangulation,
const SparsityTools::Partitioner  partitioner = SparsityTools::Partitioner::metis 
)

This function performs the same operation as the one above, except that it takes into consideration a specific set of cell_weights, which allow the partitioner to balance the graph while taking into consideration the computational effort expended on each cell.

Note
If the cell_weights vector is empty, then no weighting is taken into consideration. If not then the size of this vector must equal to the number of active cells in the triangulation.

Definition at line 2486 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::partition_triangulation ( const unsigned int  n_partitions,
const SparsityPattern cell_connection_graph,
Triangulation< dim, spacedim > &  triangulation,
const SparsityTools::Partitioner  partitioner = SparsityTools::Partitioner::metis 
)

This function does the same as the previous one, i.e. it partitions a triangulation using a partitioning algorithm into a number of subdomains identified by the cell->subdomain_id() flag.

The difference to the previous function is the second argument, a sparsity pattern that represents the connectivity pattern between cells.

While the function above builds it directly from the triangulation by considering which cells neighbor each other, this function can take a more refined connectivity graph. The sparsity pattern needs to be of size \(N\times N\), where \(N\) is the number of active cells in the triangulation. If the sparsity pattern contains an entry at position \((i,j)\), then this means that cells \(i\) and \(j\) (in the order in which they are traversed by active cell iterators) are to be considered connected; partitioning algorithm will then try to partition the domain in such a way that (i) the subdomains are of roughly equal size, and (ii) a minimal number of connections are broken.

This function is mainly useful in cases where connections between cells exist that are not present in the triangulation alone (otherwise the previous function would be the simpler one to use). Such connections may include that certain parts of the boundary of a domain are coupled through symmetric boundary conditions or integrals (e.g. friction contact between the two sides of a crack in the domain), or if a numerical scheme is used that not only connects immediate neighbors but a larger neighborhood of cells (e.g. when solving integral equations).

In addition, this function may be useful in cases where the default sparsity pattern is not entirely sufficient. This can happen because the default is to just consider face neighbors, not neighboring cells that are connected by edges or vertices. While the latter couple when using continuous finite elements, they are typically still closely connected in the neighborship graph, and partitioning algorithm will not usually cut important connections in this case. However, if there are vertices in the mesh where many cells (many more than the common 4 or 6 in 2d and 3d, respectively) come together, then there will be a significant number of cells that are connected across a vertex, but several degrees removed in the connectivity graph built only using face neighbors. In a case like this, partitioning algorithm may sometimes make bad decisions and you may want to build your own connectivity graph.

Note
If the cell_weight signal has been attached to the triangulation, then this will be used and passed to the partitioner.

Definition at line 2534 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::partition_triangulation ( const unsigned int  n_partitions,
const std::vector< unsigned int > &  cell_weights,
const SparsityPattern cell_connection_graph,
Triangulation< dim, spacedim > &  triangulation,
const SparsityTools::Partitioner  partitioner = SparsityTools::Partitioner::metis 
)

This function performs the same operation as the one above, except that it takes into consideration a specific set of cell_weights, which allow the partitioner to balance the graph while taking into consideration the computational effort expended on each cell.

Note
If the cell_weights vector is empty, then no weighting is taken into consideration. If not then the size of this vector must equal to the number of active cells in the triangulation.

Definition at line 2568 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::partition_triangulation_zorder ( const unsigned int  n_partitions,
Triangulation< dim, spacedim > &  triangulation 
)

Generates a partitioning of the active cells making up the entire domain using the same partitioning scheme as in the p4est library. After calling this function, the subdomain ids of all active cells will have values between zero and n_partitions-1. You can access the subdomain id of a cell by using cell->subdomain_id().

Definition at line 2657 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::partition_multigrid_levels ( Triangulation< dim, spacedim > &  triangulation)

Partitions the cells of a multigrid hierarchy by assigning level subdomain ids using the "youngest child" rule, that is, each cell in the hierarchy is owned by the processor who owns its left most child in the forest, and active cells have the same subdomain id and level subdomain id. You can access the level subdomain id of a cell by using cell->level_subdomain_id().

Note: This function assumes that the active cells have already been partitioned.

Definition at line 2763 of file grid_tools.cc.

template<int dim, int spacedim>
void GridTools::get_subdomain_association ( const Triangulation< dim, spacedim > &  triangulation,
std::vector< types::subdomain_id > &  subdomain 
)

For each active cell, return in the output array to which subdomain (as given by the cell->subdomain_id() function) it belongs. The output array is supposed to have the right size already when calling this function.

This function returns the association of each cell with one subdomain. If you are looking for the association of each DoF with a subdomain, use the DoFTools::get_subdomain_association function.

Definition at line 2790 of file grid_tools.cc.

template<int dim, int spacedim>
unsigned int GridTools::count_cells_with_subdomain_association ( const Triangulation< dim, spacedim > &  triangulation,
const types::subdomain_id  subdomain 
)

Count how many cells are uniquely associated with the given subdomain index.

This function may return zero if there are no cells with the given subdomain index. This can happen, for example, if you try to partition a coarse mesh into more partitions (one for each processor) than there are cells in the mesh.

This function returns the number of cells associated with one subdomain. If you are looking for the association of DoFs with this subdomain, use the DoFTools::count_dofs_with_subdomain_association function.

Definition at line 2807 of file grid_tools.cc.

template<int dim, int spacedim>
std::vector< bool > GridTools::get_locally_owned_vertices ( const Triangulation< dim, spacedim > &  triangulation)

For a triangulation, return a mask that represents which of its vertices are "owned" by the current process in the same way as we talk about locally owned cells or degrees of freedom (see GlossLocallyOwnedCell and GlossLocallyOwnedDof). For the purpose of this function, we define a locally owned vertex as follows: a vertex is owned by that processor with the smallest subdomain id (which equals the MPI rank of that processor) among all owners of cells adjacent to this vertex. In other words, vertices that are in the interior of a partition of the triangulation are owned by the owner of this partition; for vertices that lie on the boundary between two or more partitions, the owner is the processor with the least subdomain_id among all adjacent subdomains.

For sequential triangulations (as opposed to, for example, parallel::distributed::Triangulation), every user vertex is of course owned by the current processor, i.e., the function returns Triangulation::get_used_vertices(). For parallel triangulations, the returned mask is a subset of what Triangulation::get_used_vertices() returns.

Parameters
triangulationThe triangulation of which the function evaluates which vertices are locally owned.
Returns
The subset of vertices, as described above. The length of the returned array equals Triangulation.n_vertices() and may, consequently, be larger than Triangulation::n_used_vertices().

Definition at line 2826 of file grid_tools.cc.

template<typename MeshType >
std::list< std::pair< typename MeshType::cell_iterator, typename MeshType::cell_iterator > > GridTools::get_finest_common_cells ( const MeshType &  mesh_1,
const MeshType &  mesh_2 
)

Given two meshes (i.e. objects of type Triangulation, DoFHandler, or hp::DoFHandler) that are based on the same coarse mesh, this function figures out a set of cells that are matched between the two meshes and where at most one of the meshes is more refined on this cell. In other words, it finds the smallest cells that are common to both meshes, and that together completely cover the domain.

This function is useful, for example, in time-dependent or nonlinear application, where one has to integrate a solution defined on one mesh (e.g., the one from the previous time step or nonlinear iteration) against the shape functions of another mesh (the next time step, the next nonlinear iteration). If, for example, the new mesh is finer, then one has to obtain the solution on the coarse mesh (mesh_1) and interpolate it to the children of the corresponding cell of mesh_2. Conversely, if the new mesh is coarser, one has to express the coarse cell shape function by a linear combination of fine cell shape functions. In either case, one needs to loop over the finest cells that are common to both triangulations. This function returns a list of pairs of matching iterators to cells in the two meshes that can be used to this end.

Note that the list of these iterators is not necessarily ordered, and does also not necessarily coincide with the order in which cells are traversed in one, or both, of the meshes given as arguments.

Template Parameters
MeshTypeA type that satisfies the requirements of the MeshType concept.

Definition at line 1134 of file grid_tools_dof_handlers.cc.

template<int dim, int spacedim>
bool GridTools::have_same_coarse_mesh ( const Triangulation< dim, spacedim > &  mesh_1,
const Triangulation< dim, spacedim > &  mesh_2 
)

Return true if the two triangulations are based on the same coarse mesh. This is determined by checking whether they have the same number of cells on the coarsest level, and then checking that they have the same vertices.

The two meshes may have different refinement histories beyond the coarse mesh.

Definition at line 1220 of file grid_tools_dof_handlers.cc.

template<typename MeshType >
bool GridTools::have_same_coarse_mesh ( const MeshType &  mesh_1,
const MeshType &  mesh_2 
)

The same function as above, but working on arguments of type DoFHandler, or hp::DoFHandler. This function is provided to allow calling have_same_coarse_mesh for all types of containers representing triangulations or the classes built on triangulations.

Template Parameters
MeshTypeA type that satisfies the requirements of the MeshType concept.

Definition at line 1252 of file grid_tools_dof_handlers.cc.

template<int dim, int spacedim>
Triangulation< dim, spacedim >::DistortedCellList GridTools::fix_up_distorted_child_cells ( const typename Triangulation< dim, spacedim >::DistortedCellList &  distorted_cells,
Triangulation< dim, spacedim > &  triangulation 
)

Given a triangulation and a list of cells whose children have become distorted as a result of mesh refinement, try to fix these cells up by moving the center node around.

The function returns a list of cells with distorted children that couldn't be fixed up for whatever reason. The returned list is therefore a subset of the input argument.

For a definition of the concept of distorted cells, see the glossary entry. The first argument passed to the current function is typically the exception thrown by the Triangulation::execute_coarsening_and_refinement function.

Definition at line 3411 of file grid_tools.cc.

template<class MeshType >
std::vector< typename MeshType::active_cell_iterator > GridTools::get_patch_around_cell ( const typename MeshType::active_cell_iterator &  cell)

This function returns a list of all the active neighbor cells of the given, active cell. Here, a neighbor is defined as one having at least part of a face in common with the given cell, but not edge (in 3d) or vertex neighbors (in 2d and 3d).

The first element of the returned list is the cell provided as argument. The remaining ones are neighbors: The function loops over all faces of that given cell and checks if that face is not on the boundary of the domain. Then, if the neighbor cell does not have any children (i.e., it is either at the same refinement level as the current cell, or coarser) then this neighbor cell is added to the list of cells. Otherwise, if the neighbor cell is refined and therefore has children, then this function loops over all subfaces of current face adds the neighbors behind these sub-faces to the list to be returned.

Template Parameters
MeshTypeA type that satisfies the requirements of the MeshType concept. In C++, the compiler can not determine MeshType from the function call. You need to specify it as an explicit template argument following the function name.
Parameters
[in]cellAn iterator pointing to a cell of the mesh.
Returns
A list of active cells that form the patch around the given cell
Note
Patches are often used in defining error estimators that require the solution of a local problem on the patch surrounding each of the cells of the mesh. This also requires manipulating the degrees of freedom associated with the cells of a patch. To this end, there are further functions working on patches in namespace DoFTools.
In the context of a parallel distributed computation, it only makes sense to call this function on locally owned cells. This is because the neighbors of locally owned cells are either locally owned themselves, or ghost cells. For both, we know that these are in fact the real cells of the complete, parallel triangulation. We can also query the degrees of freedom on these.
Author
Arezou Ghesmati, Wolfgang Bangerth, 2014

Definition at line 1393 of file grid_tools_dof_handlers.cc.

template<class Container >
std::vector< typename Container::cell_iterator > GridTools::get_cells_at_coarsest_common_level ( const std::vector< typename Container::active_cell_iterator > &  patch_cells)

This function takes a vector of active cells (hereafter named patch_cells) as input argument, and returns a vector of their parent cells with the coarsest common level of refinement. In other words, find that set of cells living at the same refinement level so that all cells in the input vector are children of the cells in the set, or are in the set itself.

Template Parameters
ContainerIn C++, the compiler can not determine the type of Container from the function call. You need to specify it as an explicit template argument following the function name. This type has to satisfy the requirements of a mesh container (see ConceptMeshType).
Parameters
[in]patch_cellsA vector of active cells for which this function finds the parents at the coarsest common level. This vector of cells typically results from calling the function GridTools::get_patch_around_cell().
Returns
A list of cells with the coarsest common level of refinement of the input cells.
Author
Arezou Ghesmati, Wolfgang Bangerth, 2015

Definition at line 1442 of file grid_tools_dof_handlers.cc.

template<class Container >
void GridTools::build_triangulation_from_patch ( const std::vector< typename Container::active_cell_iterator > &  patch,
Triangulation< Container::dimension, Container::space_dimension > &  local_triangulation,
std::map< typename Triangulation< Container::dimension, Container::space_dimension >::active_cell_iterator, typename Container::active_cell_iterator > &  patch_to_global_tria_map 
)

This function constructs a Triangulation (named local_triangulation) from a given vector of active cells. This vector (which we think of the cells corresponding to a "patch") contains active cells that are part of an existing global Triangulation. The goal of this function is to build a local Triangulation that contains only the active cells given in patch (and potentially a minimum number of additional cells required to form a valid Triangulation). The function also returns a map that allows to identify the cells in the output Triangulation and corresponding cells in the input list.

The function copies the location of vertices of cells from the cells of the source triangulation to the triangulation that is built from the list of patch cells. This adds support for triangulations which have been perturbed or smoothed in some manner which makes the triangulation deviate from the standard deal.ii refinement strategy of placing new vertices at midpoints of faces or edges.

The operation implemented by this function is frequently used in the definition of error estimators that need to solve "local" problems on each cell and its neighbors. A similar construction is necessary in the definition of the Clement interpolation operator in which one needs to solve a local problem on all cells within the support of a shape function. This function then builds a complete Triangulation from a list of cells that make up such a patch; one can then later attach a DoFHandler to such a Triangulation.

If the list of input cells contains only cells at the same refinement level, then the output Triangulation simply consists of a Triangulation containing only exactly these patch cells. On the other hand, if the input cells live on different refinement levels, i.e., the Triangulation of which they are part is adaptively refined, then the construction of the output Triangulation is not so simple because the coarsest level of a Triangulation can not contain hanging nodes. Rather, we first have to find the common refinement level of all input cells, along with their common parents (see GridTools::get_cells_at_coarsest_common_level()), build a Triangulation from those, and then adaptively refine it so that the input cells all also exist in the output Triangulation.

A consequence of this procedure is that that output Triangulation may contain more active cells than the ones that exist in the input vector. On the other hand, one typically wants to solve the local problem not on the entire output Triangulation, but only on those cells of it that correspond to cells in the input list. In this case, a user typically wants to assign degrees of freedom only on cells that are part of the "patch", and somehow ignore those excessive cells. The current function supports this common requirement by setting the user flag for the cells in the output Triangulation that match with cells in the input list. Cells which are not part of the original patch will not have their user_flag set; we can then avoid assigning degrees of freedom using the FE_Nothing<dim> element.

Template Parameters
ContainerIn C++, the compiler can not determine the type of Container from the function call. You need to specify it as an explicit template argument following the function name. This type that satisfies the requirements of a mesh container (see ConceptMeshType).
Parameters
[in]patchA vector of active cells from a common triangulation. These cells may or may not all be at the same refinement level.
[out]local_triangulationA triangulation whose active cells correspond to the given vector of active cells in patch.
[out]patch_to_global_tria_mapA map between the local triangulation which is built as explained above, and the cell iterators in the input list.
Author
Arezou Ghesmati, Wolfgang Bangerth, 2015

Definition at line 1488 of file grid_tools_dof_handlers.cc.

template<class DoFHandlerType >
std::map< types::global_dof_index, std::vector< typename DoFHandlerType::active_cell_iterator > > GridTools::get_dof_to_support_patch_map ( DoFHandlerType &  dof_handler)

This function runs through the degrees of freedom defined by the DoFHandlerType and for each dof constructs a vector of active_cell_iterators representing the cells of support of the associated basis element at that degree of freedom. This function was originally designed for the implementation of local projections, for instance the Clement interpolant, in conjunction with other local patch functions like GridTools::build_triangulation_from_patch.

DoFHandlerType's built on top of Triangulation or parallel:distributed::Triangulation are supported and handled appropriately.

The result is the patch of cells representing the support of the basis element associated to the degree of freedom. For instance using an FE_Q finite element, we obtain the standard patch of cells touching the degree of freedom and then add other cells that take care of possible hanging node constraints. Using a FE_DGQ finite element, the degrees of freedom are logically considered to be "interior" to the cells so the patch would consist exclusively of the single cell on which the degree of freedom is located.

Template Parameters
DoFHandlerTypeThe DoFHandlerType should be a DoFHandler or hp::DoFHandler.
Parameters
[in]dof_handlerThe DoFHandlerType which could be built on a Triangulation or a parallel::distributed::Triangulation with a finite element that has degrees of freedom that are logically associated to a vertex, line, quad, or hex.
Returns
A map from the global_dof_index of degrees of freedom on locally relevant cells to vectors containing DoFHandlerType::active_cell_iterators of cells in the support of the basis function at that degree of freedom.
Author
Spencer Patty, 2016

Definition at line 1712 of file grid_tools_dof_handlers.cc.

template<typename FaceIterator >
bool GridTools::orthogonal_equality ( std::bitset< 3 > &  orientation,
const FaceIterator &  face1,
const FaceIterator &  face2,
const int  direction,
const Tensor< 1, FaceIterator::AccessorType::space_dimension > &  offset = Tensor<1, FaceIterator::AccessorType::space_dimension>(),
const FullMatrix< double > &  matrix = FullMatrix<double>() 
)
inline

An orthogonal equality test for faces.

face1 and face2 are considered equal, if a one to one matching between its vertices can be achieved via an orthogonal equality relation.

Here, two vertices v_1 and v_2 are considered equal, if \(M\cdot v_1 + offset - v_2\) is parallel to the unit vector in unit direction direction. If the parameter matrix is a reference to a spacedim x spacedim matrix, \(M\) is set to matrix, otherwise \(M\) is the identity matrix.

If the matching was successful, the relative orientation of face1 with respect to face2 is returned in the bitset orientation, where

orientation[0] -> face_orientation
orientation[1] -> face_flip
orientation[2] -> face_rotation

In 2D face_orientation is always true, face_rotation is always false, and face_flip has the meaning of line_flip. More precisely in 3d:

face_orientation: true if face1 and face2 have the same orientation. Otherwise, the vertex indices of face1 match the vertex indices of face2 in the following manner:

face1: face2:
1 - 3 2 - 3
| | <--> | |
0 - 2 0 - 1

face_flip: true if the matched vertices are rotated by 180 degrees:

face1: face2:
1 - 0 2 - 3
| | <--> | |
3 - 2 0 - 1

face_rotation: true if the matched vertices are rotated by 90 degrees counterclockwise:

face1: face2:
0 - 2 2 - 3
| | <--> | |
1 - 3 0 - 1

and any combination of that... More information on the topic can be found in the glossary article.

Author
Matthias Maier, 2012

Definition at line 2384 of file grid_tools_dof_handlers.cc.

template<typename FaceIterator >
bool GridTools::orthogonal_equality ( const FaceIterator &  face1,
const FaceIterator &  face2,
const int  direction,
const Tensor< 1, FaceIterator::AccessorType::space_dimension > &  offset = Tensor<1, FaceIterator::AccessorType::space_dimension>(),
const FullMatrix< double > &  matrix = FullMatrix<double>() 
)
inline

Same function as above, but doesn't return the actual orientation

Definition at line 2435 of file grid_tools_dof_handlers.cc.

template<typename MeshType >
void GridTools::collect_periodic_faces ( const MeshType &  mesh,
const types::boundary_id  b_id1,
const types::boundary_id  b_id2,
const int  direction,
std::vector< PeriodicFacePair< typename MeshType::cell_iterator >> &  matched_pairs,
const Tensor< 1, MeshType::space_dimension > &  offset = ::Tensor<1, MeshType::space_dimension>(),
const FullMatrix< double > &  matrix = FullMatrix<double>() 
)

This function will collect periodic face pairs on the coarsest mesh level of the given mesh (a Triangulation or DoFHandler) and add them to the vector matched_pairs leaving the original contents intact.

Define a 'first' boundary as all boundary faces having boundary_id b_id1 and a 'second' boundary consisting of all faces belonging to b_id2.

This function tries to match all faces belonging to the first boundary with faces belonging to the second boundary with the help of orthogonal_equality().

The bitset that is returned inside of PeriodicFacePair encodes the relative orientation of the first face with respect to the second face, see the documentation of orthogonal_equality() for further details.

The direction refers to the space direction in which periodicity is enforced. When matching periodic faces this vector component is ignored.

The offset is a vector tangential to the faces that is added to the location of vertices of the 'first' boundary when attempting to match them to the corresponding vertices of the 'second' boundary. This can be used to implement conditions such as \(u(0,y)=u(1,y+1)\).

Optionally, a \(dim\times dim\) rotation matrix can be specified that describes how vector valued DoFs of the first face should be modified prior to constraining to the DoFs of the second face. The matrix is used in two places. First, matrix will be supplied to orthogonal_equality() and used for matching faces: Two vertices \(v_1\) and \(v_2\) match if \(\text{matrix}\cdot v_1 + \text{offset} - v_2\) is parallel to the unit vector in unit direction direction. (For more details see DoFTools::make_periodicity_constraints(), the glossary glossary entry on periodic conditions and step-45). Second, matrix will be stored in the PeriodicFacePair collection matched_pairs for further use.

Template Parameters
MeshTypeA type that satisfies the requirements of the MeshType concept.
Note
The created std::vector can be used in DoFTools::make_periodicity_constraints() and in parallel::distributed::Triangulation::add_periodicity() to enforce periodicity algebraically.
Because elements will be added to matched_pairs (and existing entries will be preserved), it is possible to call this function several times with different boundary ids to generate a vector with all periodic pairs.
Since the periodic face pairs are found on the coarsest mesh level, it is necessary to ensure that the coarsest level faces have the correct boundary indicators set. In general, this means that one must first set all boundary indicators on the coarse grid before performing any global or local grid refinement.
Author
Daniel Arndt, Matthias Maier, 2013 - 2015

Definition at line 2176 of file grid_tools_dof_handlers.cc.

template<typename MeshType >
void GridTools::collect_periodic_faces ( const MeshType &  mesh,
const types::boundary_id  b_id,
const int  direction,
std::vector< PeriodicFacePair< typename MeshType::cell_iterator >> &  matched_pairs,
const ::Tensor< 1, MeshType::space_dimension > &  offset = ::Tensor< 1, MeshType::space_dimension >(),
const FullMatrix< double > &  matrix = FullMatrix< double >() 
)

This compatibility version of collect_periodic_faces() only works on grids with cells in standard orientation.

Instead of defining a 'first' and 'second' boundary with the help of two boundary_ids this function defines a 'left' boundary as all faces with local face index 2*dimension and boundary indicator b_id and, similarly, a 'right' boundary consisting of all face with local face index 2*dimension+1 and boundary indicator b_id.

This function will collect periodic face pairs on the coarsest mesh level and add them to matched_pairs leaving the original contents intact.

See above function for further details.

Note
This version of collect_periodic_faces() will not work on meshes with cells not in standard orientation.
Author
Daniel Arndt, Matthias Maier, 2013 - 2015
template<typename DataType , typename MeshType >
void GridTools::exchange_cell_data_to_ghosts ( const MeshType &  mesh,
const std::function< boost::optional< DataType >(const typename MeshType::active_cell_iterator &)> &  pack,
const std::function< void(const typename MeshType::active_cell_iterator &, const DataType &)> &  unpack 
)

Exchange arbitrary data of type DataType provided by the function objects from locally owned cells to ghost cells on other processors.

After this call, you typically will have received data from unpack on every ghost cell as it was given by pack on the owning processor. Whether you do or do not receive information to unpack on a given ghost cell depends on whether the pack function decided that something needs to be sent. It does so using the boost::optional mechanism: if the boost::optional return object of the pack function is empty, then this implies that no data has to be sent for the locally owned cell it was called on. In that case, unpack will also not be called on the ghost cell that corresponds to it on the receiving side. On the other hand, if the boost::optional object is not empty, then the data stored within it will be sent to the received and the unpack function called with it.

Template Parameters
DataTypeThe type of the data to be communicated. It is assumed to be serializable by boost::serialization. In many cases, this data type can not be deduced by the compiler, e.g., if you provide lambda functions for the second and third argument to this function. In this case, you have to explicitly specify the DataType as a template argument to the function call.
MeshTypeThe type of mesh.
Parameters
meshA variable of a type that satisfies the requirements of the MeshType concept.
packThe function that will be called on each locally owned cell that is a ghost cell somewhere else. As mentioned above, the function may return a regular data object of type DataType to indicate that data should be sent, or an empty boost::optional<DataType> to indicate that nothing has to be sent for this cell.
unpackThe function that will be called for each ghost cell for which data was sent, i.e., for which the pack function on the sending side returned a non-empty boost::optional object. The unpack function is then called with the data sent by the processor that owns that cell.

An example

Here is an example that shows how this function is to be used in a concrete context. It is taken from the code that makes sure that the active_fe_index (a single unsigned integer) is transported from locally owned cells where one can set it in hp::DoFHandler objects, to the corresponding ghost cells on other processors to ensure that one can query the right value also on those processors:

using active_cell_iterator =
typename ::hp::DoFHandler<dim,spacedim>::active_cell_iterator;
auto pack = [] (const active_cell_iterator &cell) -> unsigned int
{
return cell->active_fe_index();
};
auto unpack = [] (const active_cell_iterator &cell,
const unsigned int &active_fe_index) -> void
{
cell->set_active_fe_index(active_fe_index);
};
unsigned int, ::hp::DoFHandler<dim,spacedim>> (dof_handler,

You will notice that the pack lambda function returns an unsigned int, not a boost::optional<unsigned int>. The former converts automatically to the latter, implying that data will always be transported to the other processor.

(In reality, the unpack function needs to be a bit more complicated because it is not allowed to call DoFAccessor::set_active_fe_index() on ghost cells. Rather, the unpack function directly accesses internal data structures. But you get the idea – the code could, just as well, have exchanged material ids, user indices, boundary indicators, or any kind of other data with similar calls as the ones above.)

template<class StreamType >
StreamType& GridTools::operator<< ( StreamType &  s,
const CacheUpdateFlags  u 
)
inline

Output operator which outputs assemble flags as a set of or'd text values.

CacheUpdateFlags

Definition at line 80 of file grid_tools_cache_update_flags.h.

CacheUpdateFlags GridTools::operator| ( const CacheUpdateFlags  f1,
const CacheUpdateFlags  f2 
)
inline

Global operator which returns an object in which all bits are set which are either set in the first or the second argument. This operator exists since if it did not then the result of the bit-or operator | would be an integer which would in turn trigger a compiler warning when we tried to assign it to an object of type CacheUpdateFlags.

CacheUpdateFlags

Definition at line 105 of file grid_tools_cache_update_flags.h.

CacheUpdateFlags GridTools::operator~ ( const CacheUpdateFlags  f1)
inline

Global operator which returns an object in which all bits are set which are not set in the argument. This operator exists since if it did not then the result of the bit-negation operator ~ would be an integer which would in turn trigger a compiler warning when we tried to assign it to an object of type CacheUpdateFlags.

CacheUpdateFlags

Definition at line 121 of file grid_tools_cache_update_flags.h.

CacheUpdateFlags& GridTools::operator|= ( CacheUpdateFlags f1,
const CacheUpdateFlags  f2 
)
inline

Global operator which sets the bits from the second argument also in the first one.

CacheUpdateFlags

Definition at line 136 of file grid_tools_cache_update_flags.h.

CacheUpdateFlags GridTools::operator& ( const CacheUpdateFlags  f1,
const CacheUpdateFlags  f2 
)
inline

Global operator which returns an object in which all bits are set which are set in the first as well as the second argument. This operator exists since if it did not then the result of the bit-and operator & would be an integer which would in turn trigger a compiler warning when we tried to assign it to an object of type CacheUpdateFlags.

CacheUpdateFlags

Definition at line 152 of file grid_tools_cache_update_flags.h.

CacheUpdateFlags& GridTools::operator&= ( CacheUpdateFlags f1,
const CacheUpdateFlags  f2 
)
inline

Global operator which clears all the bits in the first argument if they are not also set in the second argument.

CacheUpdateFlags

Definition at line 167 of file grid_tools_cache_update_flags.h.