Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Classes | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | Friends | List of all members
Triangulation< dim, spacedim > Class Template Reference

#include <deal.II/grid/tria.h>

Inheritance diagram for Triangulation< dim, spacedim >:
Inheritance graph
[legend]

Classes

struct  CellWeightSum
 
struct  DistortedCellList
 
struct  Signals
 

Public Types

enum  MeshSmoothing {
  none = 0x0 , limit_level_difference_at_vertices = 0x1 , eliminate_unrefined_islands = 0x2 , patch_level_1 = 0x4 ,
  coarsest_level_1 = 0x8 , allow_anisotropic_smoothing = 0x10 , eliminate_refined_inner_islands = 0x100 , eliminate_refined_boundary_islands = 0x200 ,
  do_not_produce_unrefined_islands = 0x400 , smoothing_on_refinement , smoothing_on_coarsening , maximum_smoothing = 0xffff ^ allow_anisotropic_smoothing
}
 
using cell_iterator = TriaIterator< CellAccessor< dim, spacedim > >
 
using level_cell_iterator = cell_iterator
 
using active_cell_iterator = TriaActiveIterator< CellAccessor< dim, spacedim > >
 
using face_iterator = TriaIterator< TriaAccessor< dim - 1, dim, spacedim > >
 
using active_face_iterator = TriaActiveIterator< TriaAccessor< dim - 1, dim, spacedim > >
 
using vertex_iterator = TriaIterator<::TriaAccessor< 0, dim, spacedim > >
 
using active_vertex_iterator = TriaActiveIterator<::TriaAccessor< 0, dim, spacedim > >
 
using line_iterator = typename IteratorSelector::line_iterator
 
using active_line_iterator = typename IteratorSelector::active_line_iterator
 
using quad_iterator = typename IteratorSelector::quad_iterator
 
using active_quad_iterator = typename IteratorSelector::active_quad_iterator
 
using hex_iterator = typename IteratorSelector::hex_iterator
 
using active_hex_iterator = typename IteratorSelector::active_hex_iterator
 
using CellStatus = ::CellStatus
 

Public Member Functions

 Triangulation (const MeshSmoothing smooth_grid=none, const bool check_for_distorted_cells=false)
 
 Triangulation (const Triangulation< dim, spacedim > &)=delete
 
 Triangulation (Triangulation< dim, spacedim > &&tria) noexcept
 
Triangulationoperator= (Triangulation< dim, spacedim > &&tria) noexcept
 
virtual ~Triangulation () override
 
virtual void clear ()
 
virtual MPI_Comm get_communicator () const
 
virtual std::weak_ptr< const Utilities::MPI::Partitionerglobal_active_cell_index_partitioner () const
 
virtual std::weak_ptr< const Utilities::MPI::Partitionerglobal_level_cell_index_partitioner (const unsigned int level) const
 
virtual void set_mesh_smoothing (const MeshSmoothing mesh_smoothing)
 
virtual const MeshSmoothingget_mesh_smoothing () const
 
void set_manifold (const types::manifold_id number, const Manifold< dim, spacedim > &manifold_object)
 
void reset_manifold (const types::manifold_id manifold_number)
 
void reset_all_manifolds ()
 
void set_all_manifold_ids (const types::manifold_id number)
 
void set_all_manifold_ids_on_boundary (const types::manifold_id number)
 
void set_all_manifold_ids_on_boundary (const types::boundary_id b_id, const types::manifold_id number)
 
const Manifold< dim, spacedim > & get_manifold (const types::manifold_id number) const
 
virtual std::vector< types::boundary_idget_boundary_ids () const
 
virtual std::vector< types::manifold_idget_manifold_ids () const
 
virtual void copy_triangulation (const Triangulation< dim, spacedim > &other_tria)
 
virtual void create_triangulation (const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
 
virtual void create_triangulation (const TriangulationDescription::Description< dim, spacedim > &construction_data)
 
void flip_all_direction_flags ()
 
virtual std::size_t memory_consumption () const
 
template<class Archive >
void save (Archive &ar, const unsigned int version) const
 
template<class Archive >
void load (Archive &ar, const unsigned int version)
 
virtual void save (const std::string &filename) const
 
virtual void load (const std::string &filename)
 
virtual void add_periodicity (const std::vector< GridTools::PeriodicFacePair< cell_iterator > > &)
 
const std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, std::bitset< 3 > > > & get_periodic_face_map () const
 
const std::vector< ReferenceCell > & get_reference_cells () const
 
bool all_reference_cells_are_hyper_cube () const
 
bool all_reference_cells_are_simplex () const
 
bool is_mixed_mesh () const
 
template<class Archive >
void serialize (Archive &archive, const unsigned int version)
 
unsigned int n_quads () const
 
unsigned int n_quads (const unsigned int) const
 
unsigned int n_raw_quads (const unsigned int) const
 
unsigned int n_raw_hexs (const unsigned int) const
 
unsigned int n_active_quads (const unsigned int) const
 
unsigned int n_active_quads () const
 
unsigned int n_quads () const
 
unsigned int n_quads (const unsigned int) const
 
unsigned int n_raw_quads (const unsigned int) const
 
unsigned int n_raw_hexs (const unsigned int) const
 
unsigned int n_active_quads (const unsigned int) const
 
unsigned int n_active_quads () const
 
unsigned int n_quads () const
 
unsigned int n_quads (const unsigned int) const
 
unsigned int n_raw_quads (const unsigned int) const
 
unsigned int n_raw_hexs (const unsigned int) const
 
unsigned int n_active_quads (const unsigned int) const
 
unsigned int n_active_quads () const
 
unsigned int n_raw_quads (const unsigned int level) const
 
unsigned int n_raw_quads (const unsigned int level) const
 
unsigned int n_raw_quads (const unsigned int) const
 
unsigned int n_raw_quads () const
 
unsigned int n_hexs () const
 
unsigned int n_hexs (const unsigned int level) const
 
unsigned int n_raw_hexs (const unsigned int level) const
 
unsigned int n_active_hexs () const
 
unsigned int n_active_hexs (const unsigned int level) const
 
unsigned int max_adjacent_cells () const
 
unsigned int max_adjacent_cells () const
 
unsigned int max_adjacent_cells () const
 
bool prepare_coarsening_and_refinement ()
 
bool prepare_coarsening_and_refinement ()
 
bool prepare_coarsening_and_refinement ()
 
Mesh refinement
void set_all_refine_flags ()
 
void refine_global (const unsigned int times=1)
 
void coarsen_global (const unsigned int times=1)
 
virtual void execute_coarsening_and_refinement ()
 
virtual bool prepare_coarsening_and_refinement ()
 
History of a triangulation
void save_refine_flags (std::ostream &out) const
 
void save_refine_flags (std::vector< bool > &v) const
 
void load_refine_flags (std::istream &in)
 
void load_refine_flags (const std::vector< bool > &v)
 
void save_coarsen_flags (std::ostream &out) const
 
void save_coarsen_flags (std::vector< bool > &v) const
 
void load_coarsen_flags (std::istream &out)
 
void load_coarsen_flags (const std::vector< bool > &v)
 
bool get_anisotropic_refinement_flag () const
 
User data
void clear_user_flags ()
 
void save_user_flags (std::ostream &out) const
 
void save_user_flags (std::vector< bool > &v) const
 
void load_user_flags (std::istream &in)
 
void load_user_flags (const std::vector< bool > &v)
 
void clear_user_flags_line ()
 
void save_user_flags_line (std::ostream &out) const
 
void save_user_flags_line (std::vector< bool > &v) const
 
void load_user_flags_line (std::istream &in)
 
void load_user_flags_line (const std::vector< bool > &v)
 
void clear_user_flags_quad ()
 
void save_user_flags_quad (std::ostream &out) const
 
void save_user_flags_quad (std::vector< bool > &v) const
 
void load_user_flags_quad (std::istream &in)
 
void load_user_flags_quad (const std::vector< bool > &v)
 
void clear_user_flags_hex ()
 
void save_user_flags_hex (std::ostream &out) const
 
void save_user_flags_hex (std::vector< bool > &v) const
 
void load_user_flags_hex (std::istream &in)
 
void load_user_flags_hex (const std::vector< bool > &v)
 
void clear_user_data ()
 
void save_user_indices (std::vector< unsigned int > &v) const
 
void load_user_indices (const std::vector< unsigned int > &v)
 
void save_user_pointers (std::vector< void * > &v) const
 
void load_user_pointers (const std::vector< void * > &v)
 
void save_user_indices_line (std::vector< unsigned int > &v) const
 
void load_user_indices_line (const std::vector< unsigned int > &v)
 
void save_user_indices_quad (std::vector< unsigned int > &v) const
 
void load_user_indices_quad (const std::vector< unsigned int > &v)
 
void save_user_indices_hex (std::vector< unsigned int > &v) const
 
void load_user_indices_hex (const std::vector< unsigned int > &v)
 
void save_user_pointers_line (std::vector< void * > &v) const
 
void load_user_pointers_line (const std::vector< void * > &v)
 
void save_user_pointers_quad (std::vector< void * > &v) const
 
void load_user_pointers_quad (const std::vector< void * > &v)
 
void save_user_pointers_hex (std::vector< void * > &v) const
 
void load_user_pointers_hex (const std::vector< void * > &v)
 
Cell iterator functions
cell_iterator begin (const unsigned int level=0) const
 
active_cell_iterator begin_active (const unsigned int level=0) const
 
cell_iterator end () const
 
cell_iterator end (const unsigned int level) const
 
active_cell_iterator end_active (const unsigned int level) const
 
cell_iterator last () const
 
active_cell_iterator last_active () const
 
cell_iterator create_cell_iterator (const CellId &cell_id) const
 
bool contains_cell (const CellId &cell_id) const
 
Cell iterator functions returning ranges of iterators
IteratorRange< cell_iteratorcell_iterators () const
 
IteratorRange< active_cell_iteratoractive_cell_iterators () const
 
IteratorRange< cell_iteratorcell_iterators_on_level (const unsigned int level) const
 
IteratorRange< active_cell_iteratoractive_cell_iterators_on_level (const unsigned int level) const
 
Face iterator functions
face_iterator begin_face () const
 
active_face_iterator begin_active_face () const
 
face_iterator end_face () const
 
IteratorRange< active_face_iteratoractive_face_iterators () const
 
Vertex iterator functions
vertex_iterator begin_vertex () const
 
active_vertex_iterator begin_active_vertex () const
 
vertex_iterator end_vertex () const
 
Information about the triangulation
unsigned int n_lines () const
 
unsigned int n_lines (const unsigned int level) const
 
unsigned int n_active_lines () const
 
unsigned int n_active_lines (const unsigned int level) const
 
unsigned int n_quads () const
 
unsigned int n_quads (const unsigned int level) const
 
unsigned int n_active_quads () const
 
unsigned int n_active_quads (const unsigned int level) const
 
unsigned int n_hexs () const
 
unsigned int n_hexs (const unsigned int level) const
 
unsigned int n_active_hexs () const
 
unsigned int n_active_hexs (const unsigned int level) const
 
unsigned int n_cells () const
 
unsigned int n_cells (const unsigned int level) const
 
unsigned int n_active_cells () const
 
virtual types::global_cell_index n_global_active_cells () const
 
unsigned int n_active_cells (const unsigned int level) const
 
virtual types::coarse_cell_id n_global_coarse_cells () const
 
unsigned int n_faces () const
 
unsigned int n_active_faces () const
 
unsigned int n_levels () const
 
virtual unsigned int n_global_levels () const
 
virtual bool has_hanging_nodes () const
 
unsigned int n_vertices () const
 
const std::vector< Point< spacedim > > & get_vertices () const
 
unsigned int n_used_vertices () const
 
bool vertex_used (const unsigned int index) const
 
const std::vector< bool > & get_used_vertices () const
 
unsigned int max_adjacent_cells () const
 
virtual types::subdomain_id locally_owned_subdomain () const
 
Triangulation< dim, spacedim > & get_triangulation ()
 
const Triangulation< dim, spacedim > & get_triangulation () const
 
Internal information about the number of objects
unsigned int n_raw_lines () const
 
unsigned int n_raw_lines (const unsigned int level) const
 
unsigned int n_raw_quads () const
 
unsigned int n_raw_quads (const unsigned int level) const
 
unsigned int n_raw_hexs (const unsigned int level) const
 
unsigned int n_raw_cells (const unsigned int level) const
 
unsigned int n_raw_faces () const
 
Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 

Static Public Member Functions

static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
Exceptions
static ::ExceptionBaseExcInvalidLevel (int arg1, int arg2)
 
static ::ExceptionBaseExcTriangulationNotEmpty (int arg1, int arg2)
 
static ::ExceptionBaseExcGridReadError ()
 
static ::ExceptionBaseExcFacesHaveNoLevel ()
 
static ::ExceptionBaseExcEmptyLevel (int arg1)
 
static ::ExceptionBaseExcBoundaryIdNotFound (types::boundary_id arg1)
 
static ::ExceptionBaseExcInconsistentCoarseningFlags ()
 

Public Attributes

Keeping up with what happens to a triangulation
Signals signals
 

Static Public Attributes

static constexpr unsigned int dimension = dim
 
static constexpr unsigned int space_dimension = spacedim
 
static constexpr auto CELL_PERSIST
 
static constexpr auto CELL_REFINE
 
static constexpr auto CELL_COARSEN
 
static constexpr auto CELL_INVALID
 

Protected Member Functions

void update_periodic_face_map ()
 
virtual void update_reference_cells ()
 

Static Protected Member Functions

static void write_bool_vector (const unsigned int magic_number1, const std::vector< bool > &v, const unsigned int magic_number2, std::ostream &out)
 
static void read_bool_vector (const unsigned int magic_number1, std::vector< bool > &v, const unsigned int magic_number2, std::istream &in)
 

Protected Attributes

MeshSmoothing smooth_grid
 
std::vector< ReferenceCellreference_cells
 

Private Types

using IteratorSelector = ::internal::TriangulationImplementation::Iterators< dim, spacedim >
 
using raw_cell_iterator = TriaRawIterator< CellAccessor< dim, spacedim > >
 
using raw_face_iterator = TriaRawIterator< TriaAccessor< dim - 1, dim, spacedim > >
 
using raw_vertex_iterator = TriaRawIterator<::TriaAccessor< 0, dim, spacedim > >
 
using raw_line_iterator = typename IteratorSelector::raw_line_iterator
 
using raw_quad_iterator = typename IteratorSelector::raw_quad_iterator
 
using raw_hex_iterator = typename IteratorSelector::raw_hex_iterator
 
using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

void clear_despite_subscriptions ()
 
void reset_policy ()
 
void reset_active_cell_indices ()
 
void reset_global_cell_indices ()
 
void reset_cell_vertex_indices_cache ()
 
DistortedCellList execute_refinement ()
 
void execute_coarsening ()
 
void fix_coarsen_flags ()
 
virtual unsigned int coarse_cell_id_to_coarse_cell_index (const types::coarse_cell_id coarse_cell_id) const
 
virtual types::coarse_cell_id coarse_cell_index_to_coarse_cell_id (const unsigned int coarse_cell_index) const
 
void check_no_subscribers () const noexcept
 
Cell iterator functions for internal use
raw_cell_iterator begin_raw (const unsigned int level=0) const
 
raw_cell_iterator end_raw (const unsigned int level) const
 
Line iterator functions for internal use
raw_line_iterator begin_raw_line (const unsigned int level=0) const
 
line_iterator begin_line (const unsigned int level=0) const
 
active_line_iterator begin_active_line (const unsigned int level=0) const
 
line_iterator end_line () const
 
Quad iterator functions for internal use
raw_quad_iterator begin_raw_quad (const unsigned int level=0) const
 
quad_iterator begin_quad (const unsigned int level=0) const
 
active_quad_iterator begin_active_quad (const unsigned int level=0) const
 
quad_iterator end_quad () const
 
Hex iterator functions for internal use
raw_hex_iterator begin_raw_hex (const unsigned int level=0) const
 
hex_iterator begin_hex (const unsigned int level=0) const
 
active_hex_iterator begin_active_hex (const unsigned int level=0) const
 
hex_iterator end_hex () const
 

Private Attributes

std::unique_ptr< ::internal::TriangulationImplementation::Policy< dim, spacedim > > policy
 
std::vector< GridTools::PeriodicFacePair< cell_iterator > > periodic_face_pairs_level_0
 
std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, std::bitset< 3 > > > periodic_face_map
 
std::vector< std::unique_ptr<::internal::TriangulationImplementation::TriaLevel > > levels
 
std::unique_ptr<::internal::TriangulationImplementation::TriaFacesfaces
 
std::vector< Point< spacedim > > vertices
 
std::vector< boolvertices_used
 
std::map< types::manifold_id, std::unique_ptr< const Manifold< dim, spacedim > > > manifolds
 
bool anisotropic_refinement
 
const bool check_for_distorted_cells
 
::internal::TriangulationImplementation::NumberCache< dim > number_cache
 
std::unique_ptr< std::map< unsigned int, types::boundary_id > > vertex_to_boundary_id_map_1d
 
std::unique_ptr< std::map< unsigned int, types::manifold_id > > vertex_to_manifold_id_map_1d
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Friends

template<int , int , int >
class TriaAccessorBase
 
template<int , int , int >
class TriaAccessor
 
class TriaAccessor< 0, 1, spacedim >
 
class CellAccessor< dim, spacedim >
 
struct ::internal::TriaAccessorImplementation::Implementation
 
struct ::internal::TriangulationImplementation::Implementation
 
struct ::internal::TriangulationImplementation::ImplementationMixedMesh
 
class ::internal::TriangulationImplementation::TriaObjects
 

Serialization facilities.

internal::CellAttachedData< dim, spacedim > cell_attached_data
 
std::vector< typename internal::CellAttachedDataSerializer< dim, spacedim >::cell_relation_t > local_cell_relations
 
internal::CellAttachedDataSerializer< dim, spacedim > data_serializer
 
unsigned int register_data_attach (const std::function< std::vector< char >(const cell_iterator &, const ::CellStatus)> &pack_callback, const bool returns_variable_size_data)
 
void notify_ready_to_unpack (const unsigned int handle, const std::function< void(const cell_iterator &, const ::CellStatus, const boost::iterator_range< std::vector< char >::const_iterator > &)> &unpack_callback)
 
void save_attached_data (const unsigned int global_first_cell, const unsigned int global_num_cells, const std::string &filename) const
 
void load_attached_data (const unsigned int global_first_cell, const unsigned int global_num_cells, const unsigned int local_num_cells, const std::string &filename, const unsigned int n_attached_deserialize_fixed, const unsigned int n_attached_deserialize_variable)
 
virtual void update_cell_relations ()
 

Detailed Description

template<int dim, int spacedim = dim>
class Triangulation< dim, spacedim >

A triangulation is a collection of cells that, jointly, cover the domain on which one typically wants to solve a partial differential equation. This domain, and the mesh that covers it, represents a dim -dimensional manifold and lives in spacedim spatial dimensions, where dim and spacedim are the template arguments of this class. (If spacedim is not specified, it takes the default value spacedim=dim.)

Thus, for example, an object of type Triangulation<1,1> (or simply Triangulation<1> since spacedim==dim by default) is used to represent and handle the usual one-dimensional triangulation used in the finite element method (so, segments on a straight line). On the other hand, objects such as Triangulation<1,2> or Triangulation<2,3> (that are associated with curves in 2d or surfaces in 3d) are the ones one wants to use in the boundary element method.

The name of the class is mostly hierarchical and is not meant to imply that a Triangulation can only consist of triangles. Instead, triangulations consist of line segments in 1d (i.e., if dim==1), and of three-dimensional cells (if dim==3). Moreover, historically, deal.II only supported quadrilaterals (cells with four vertices: deformed rectangles) in 2d and hexahedra (cells with six sides and eight vertices that are deformed boxes), neither of which are triangles. In other words, the term "triangulation" in the deal.II language is synonymous with "mesh" and is to be understood separate from its linguistic origin.

This class is written to be as independent of the dimension as possible (thus the complex construction of the internal::TriangulationImplementation::TriaLevel classes) to allow code-sharing, to allow reducing the need to mirror changes in the code for one dimension to the code for other dimensions. Nonetheless, some of the functions are dependent of the dimension and there only exist specialized versions for distinct dimensions.

This class satisfies the MeshType concept requirements.

Structure and iterators

The actual data structure of a Triangulation object is rather complex and quite inconvenient if one attempted to operate on it directly, since data is spread over quite a lot of arrays and other places. However, there are ways powerful enough to work on these data structures without knowing their exact relations. deal.II uses class local alias (see below) to make things as easy and dimension independent as possible.

The Triangulation class provides iterators which enable looping over all cells without knowing the exact representation used to describe them. For more information see the documentation of TriaIterator. Their names are alias imported from the Iterators class (thus making them local types to this class) and are as follows:

For dim==1, these iterators are mapped as follows:

typename IteratorSelector::line_iterator line_iterator
Definition tria.h:1657
typename IteratorSelector::active_line_iterator active_line_iterator
Definition tria.h:1672

while for dim==2 we have the additional face iterator:

typename IteratorSelector::active_quad_iterator active_quad_iterator
Definition tria.h:1696
typename IteratorSelector::quad_iterator quad_iterator
Definition tria.h:1681

By using the cell iterators, you can write code independent of the spatial dimension. The same applies for substructure iterators, where a substructure is defined as a face of a cell. The face of a cell is a vertex in 1d and a line in 2d; however, vertices are handled in a different way and therefore lines have no faces.

The Triangulation class offers functions like begin_active() which gives you an iterator to the first active cell. There are quite a lot of functions returning iterators. Take a look at the class doc to get an overview.

Usage of these iterators is similar to usage of standard container iterators. Some examples taken from the Triangulation source code follow (notice that in the last two examples the template parameter spacedim has been omitted, so it takes the default value dim).

Usage

Usage of a Triangulation is mainly done through the use of iterators. An example probably shows best how to use it:

int main ()
{
// read in a coarse grid file
// we want to log the refinement history
ofstream history ("mesh.history");
// refine first cell
tria.begin_active()->set_refine_flag();
tria.save_refine_flags (history);
// refine first active cell on coarsest level
tria.begin_active()->set_refine_flag ();
tria.save_refine_flags (history);
for (int i=0; i<17; ++i)
{
// refine the presently second last cell 17 times
cell = tria.last_active(tria.n_levels()-1);
--cell;
cell->set_refine_flag ();
tria.save_refine_flags (history);
};
// output the grid
ofstream out("grid.1");
}
void write_gnuplot(const Triangulation< dim, spacedim > &tria, std::ostream &out, const Mapping< dim, spacedim > *mapping=nullptr) const
Definition grid_out.cc:4598
active_cell_iterator last_active() const
unsigned int n_levels() const
void save_refine_flags(std::ostream &out) const
active_cell_iterator begin_active(const unsigned int level=0) const

Creating a triangulation

There are several possibilities to create a triangulation:

Finally, there is a special function for folks who like bad grids: distort_random(). It moves all the vertices in the grid a bit around by a random value, leaving behind a distorted mesh. Note that you should apply this function to the final mesh, since refinement smoothes the mesh a bit.

The function will make sure that vertices on restricted faces (hanging nodes) will end up in the correct place, i.e. in the middle of the two other vertices of the mother line, and the analogue in higher space dimensions (vertices on the boundary are not corrected, so don't distort boundary vertices in more than two space dimension, i.e. in dimensions where boundary vertices can be hanging nodes). Applying the algorithm has another drawback related to the placement of cells, however: the children of a cell will not occupy the same region of the domain as the mother cell does. While this is the usual behavior with cells at the boundary, here you may get into trouble when using multigrid algorithms or when transferring solutions from coarse to fine grids and back. In general, the use of this function is only safe if you only use the most refined level of the triangulation for computations.

Refinement and coarsening of a triangulation

Refinement of a triangulation may be done through several ways. The most low-level way is directly through iterators: let i be an iterator to an active cell (i.e. the cell pointed to has no children), then the function call i->set_refine_flag() marks the respective cell for refinement. Marking non-active cells results in an error.

After all the cells you wanted to mark for refinement, call execute_coarsening_and_refinement() to actually perform the refinement. This function itself first calls the prepare_coarsening_and_refinement function to regularize the resulting triangulation: since a face between two adjacent cells may only be subdivided once (i.e. the levels of two adjacent cells may differ by one at most; it is not possible to have a cell refined twice while the neighboring one is not refined), some additional cells are flagged for refinement to smooth the grid. This enlarges the number of resulting cells but makes the grid more regular, thus leading to better approximation properties and, above all, making the handling of data structures and algorithms much easier. To be honest, this is mostly an algorithmic step than one needed by the finite element method.

To coarsen a grid, the same way as above is possible by using i->set_coarsen_flag and calling execute_coarsening_and_refinement().

The reason for first coarsening, then refining is that the refinement usually adds some additional cells to keep the triangulation regular and thus satisfies all refinement requests, while the coarsening does not delete cells not requested for; therefore the refinement will often revert some effects of coarsening while the opposite is not true. The stated order of coarsening before refinement will thus normally lead to a result closer to the intended one.

Marking cells for refinement 'by hand' through iterators is one way to produce a new grid, especially if you know what kind of grid you are looking for, e.g. if you want to have a grid successively refined towards the boundary or always at the center (see the example programs, they do exactly these things). There are more advanced functions, however, which are more suitable for automatic generation of hierarchical grids in the context of a posteriori error estimation and adaptive finite elements. These functions can be found in the GridRefinement class.

Smoothing of a triangulation

Some degradation of approximation properties has been observed for grids which are too unstructured. Therefore, prepare_coarsening_and_refinement() which is automatically called by execute_coarsening_and_refinement() can do some smoothing of the triangulation. Note that mesh smoothing is only done for two or more space dimensions, no smoothing is available at present for one spatial dimension. In the following, let execute_* stand for execute_coarsening_and_refinement().

For the purpose of smoothing, the Triangulation constructor takes an argument specifying whether a smoothing step shall be performed on the grid each time execute_* is called. The default is that such a step not be done, since this results in additional cells being produced, which may not be necessary in all cases. If switched on, calling execute_* results in flagging additional cells for refinement to avoid vertices as the ones mentioned. The algorithms for both regularization and smoothing of triangulations are described below in the section on technical issues. The reason why this parameter must be given to the constructor rather than to execute_* is that it would result in algorithmic problems if you called execute_* once without and once with smoothing, since then in some refinement steps would need to be refined twice.

The parameter taken by the constructor is an integer which may be composed bitwise by the constants defined in the enum MeshSmoothing (see there for the possibilities).

Note
While it is possible to pass all of the flags in MeshSmoothing to objects of type parallel::distributed::Triangulation, it is not always possible to honor all of these smoothing options if they would require knowledge of refinement/coarsening flags on cells not locally owned by this processor. As a consequence, for some of these flags, the ultimate number of cells of the parallel triangulation may depend on the number of processors into which it is partitioned.

Material and boundary information

Each cell, face or edge stores information denoting the material or the part of the boundary that an object belongs to. The material id of a cell is typically used to identify which cells belong to a particular part of the domain, e.g., when you have different materials (steel, concrete, wood) that are all part of the same domain. One would then usually query the material id associated with a cell during assembly of the bilinear form, and use it to determine (e.g., by table lookup, or a sequence of if-else statements) what the correct material coefficients would be for that cell. See also this glossary entry.

This material_id may be set upon construction of a triangulation (through the CellData data structure), or later through use of cell iterators. For a typical use of this functionality, see the step-28 tutorial program. The functions of the GridGenerator namespace typically set the material ID of all cells to zero. When reading a triangulation through the GridIn class, different input file formats have different conventions, but typically either explicitly specify the material id, or if they don't, then GridIn simply sets them to zero. Because the material of a cell is intended to pertain to a particular region of the domain, material ids are inherited by child cells from their parent upon mesh refinement.

Boundary indicators on lower dimensional objects (these have no material id) indicate the number of a boundary component. The weak formulation of the partial differential equation may have different boundary conditions on different parts of the boundary. The boundary indicator can be used in creating the matrix or the right hand side vector to indicate these different parts of the model (this use is like the material id of cells). Boundary indicators may be in the range from zero to numbers::internal_face_boundary_id-1. The value numbers::internal_face_boundary_id is reserved to denote interior lines (in 2d) and interior lines and quads (in 3d), which do not have a boundary indicator. This way, a program can easily determine, whether such an object is at the boundary or not. Material indicators may be in the range from zero to numbers::invalid_material_id-1.

Lines in two dimensions and quads in three dimensions inherit their boundary indicator to their children upon refinement. You should therefore make sure that if you have different boundary parts, the different parts are separated by a vertex (in 2d) or a line (in 3d) such that each boundary line or quad has a unique boundary indicator.

By default (unless otherwise specified during creation of a triangulation), all parts of the boundary have boundary indicator zero. As a historical wart, this isn't true for 1d meshes, however: For these, leftmost vertices have boundary indicator zero while rightmost vertices have boundary indicator one. In either case, the boundary indicator of a face can be changed using a call of the kind cell->face(1)->set_boundary_id(42);.

See also
Glossary entry on boundary indicators

History of a triangulation

It is possible to reconstruct a grid from its refinement history, which can be stored and loaded through the save_refine_flags and load_refine_flags functions. Normally, the code will look like this:

// open output file
std::ofstream history("mesh.history");
// do 10 refinement steps
for (unsigned int step=0; step<10; ++step)
{
...;
// flag cells according to some criterion
...;
tria.save_refine_flags (history);
}

If you want to re-create the grid from the stored information, you write:

// open input file
std::ifstream history("mesh.history");
// do 10 refinement steps
for (unsigned int step=0; step<10; ++step)
{
tria.load_refine_flags (history);
}
void load_refine_flags(std::istream &in)

The same scheme is employed for coarsening and the coarsening flags.

You may write other information to the output file between different sets of refinement information, as long as you read it upon re-creation of the grid. You should make sure that the other information in the new triangulation which is to be created from the saved flags, matches that of the old triangulation, for example the smoothing level; if not, the cells actually created from the flags may be other ones, since smoothing adds additional cells, but their number may be depending on the smoothing level.

There actually are two sets of save_*_flags and load_*_flags functions. One takes a stream as argument and reads/writes the information from/to the stream, thus enabling storing flags to files. The other set takes an argument of type vector<bool>. This enables the user to temporarily store some flags, e.g. if another function needs them, and restore them afterwards.

User flags and data

A triangulation offers one bit per subobject for user flags. This field can be accessed as all other data using iterators. Normally, this user flag is used if an algorithm walks over all cells and needs information whether another cell, e.g. a neighbor, has already been processed. See the glossary for more information.

There is another set of user data, which can be either an unsigned int or a void *, for each subobject. You can access these through the functions listed under User data in the accessor classes. Again, see the glossary for more information.

The value of these user indices or pointers is nullptr by default. Note that the pointers are not inherited to children upon refinement. Still, after a remeshing they are available on all cells, where they were set on the previous mesh.

The usual warning about the missing type safety of void pointers are obviously in place here; responsibility for correctness of types etc lies entirely with the user of the pointer.

Note
User pointers and user indices are stored in the same place. In order to avoid unwanted conversions, Triangulation checks which one of them is in use and does not allow access to the other one, until clear_user_data() has been called.

Describing curved geometries

deal.II implements all geometries (curved and otherwise) with classes inheriting from Manifold; see the documentation of Manifold, step-49, or the Manifold description for triangulations module for examples and a complete description of the algorithms. By default, all cells in a Triangulation have a flat geometry, meaning that all lines in the Triangulation are assumed to be straight. If a cell has a manifold_id that is not equal to numbers::flat_manifold_id then the Triangulation uses the associated Manifold object for computations on that cell (e.g., cell refinement). Here is a quick example, taken from the implementation of GridGenerator::hyper_ball(), that sets up a polar grid:

int main ()
{
const std::vector<Point<2>> vertices = {{-1.0,-1.0},
{+1.0,-1.0},
{-0.5,-0.5},
{+0.5,-0.5},
{-0.5,+0.5},
{+1.0,+1.0},
{-1.0,+1.0},
{+1.0,+1.0}};
const std::vector<std::array<int,GeometryInfo<2>::vertices_per_cell>>
cell_vertices = {{0, 1, 2, 3},
{0, 2, 6, 4},
{2, 3, 4, 5},
{1, 7, 3, 5},
{6, 4, 7, 5}};
std::vector<CellData<2>> cells(cell_vertices.size(), CellData<2>());
for (unsigned int i=0; i<cell_vertices.size(); ++i)
for (unsigned int j=0; j<GeometryInfo<2>::vertices_per_cell; ++j)
cells[i].vertices[j] = cell_vertices[i][j];
triangulation.create_triangulation (vertices, cells, SubCellData());
triangulation.set_all_manifold_ids_on_boundary(42);
// set_manifold stores a copy of its second argument,
// so a temporary is okay
triangulation.set_manifold(42, PolarManifold<2>());
for (unsigned int i = 0; i < 4; ++i)
{
// refine all boundary cells
for (const auto &cell : triangulation.active_cell_iterators())
if (cell->at_boundary())
cell->set_refine_flag();
triangulation.execute_coarsening_and_refinement();
}
}
std::vector< Point< spacedim > > vertices
Definition tria.h:4479
if(marked_vertices.size() !=0) for(auto it
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation

This will set up a grid where the boundary lines will be refined by performing calculations in polar coordinates. When the mesh is refined the cells adjacent to the boundary will use this new line midpoint (as well as the other three midpoints and original cell vertices) to calculate the cell midpoint with a transfinite interpolation: this propagates the curved boundary into the interior in a smooth way. It is possible to generate a better grid (which interpolates across all cells between two different Manifold descriptions, instead of just going one cell at a time) by using TransfiniteInterpolationManifold; see the documentation of that class for more information.

You should take note of one caveat: if you have concave boundaries, you must make sure that a new boundary vertex does not lie too much inside the cell which is to be refined. The reason is that the center vertex is placed at the point which is a weighted average of the vertices of the original cell, new face midpoints, and (in 3d) new line midpoints. Therefore if your new boundary vertex is too near the center of the old quadrilateral or hexahedron, the distance to the midpoint vertex will become too small, thus generating distorted cells. This issue is discussed extensively in distorted cells.

Getting notice when a triangulation changes

There are cases where one object would like to know whenever a triangulation is being refined, copied, or modified in a number of other ways. This could of course be achieved if, in your user code, you tell every such object whenever you are about to refine the triangulation, but this will get tedious and is error prone. The Triangulation class implements a more elegant way to achieve this: signals.

In essence, a signal is an object (a member of the Triangulation class) that another object can connect to. A connection is in essence that the connecting object passes a function object taking a certain number and kind of arguments. Whenever the owner of the signal wants to indicate a certain kind of event, it 'triggers' the signal, which in turn means that all connections of the signal are triggered: in other word, the function objects are executed and can take the action that is necessary.

As a simple example, the following code will print something to the output every time the triangulation has just been refined:

void f()
{
std::cout << "Triangulation has been refined." << std::endl;
}
void run ()
{
// fill it somehow
triangulation.refine_global (2);
}
Signals signals
Definition tria.h:2529
boost::signals2::signal< void()> post_refinement
Definition tria.h:2327

This code will produce output twice, once for each refinement cycle.

A more interesting application would be the following, akin to what the FEValues class does. This class stores a pointer to a triangulation and also an iterator to the cell last handled (so that it can compare the current cell with the previous one and, for example, decide that there is no need to re-compute the Jacobian matrix if the new cell is a simple translation of the previous one). However, whenever the triangulation is modified, the iterator to the previously handled cell needs to be invalidated since it now no longer points to any useful cell (or, at the very least, points to something that may not necessarily resemble the cells previously handled). The code would look something like this (the real code has some more error checking and has to handle the case that subsequent cells might actually belong to different triangulation, but that is of no concern to us here):

template <int dim>
class FEValues
{
Triangulation<dim>::active_cell_iterator current_cell, previous_cell;
public:
void invalidate_previous_cell ();
};
template <int dim>
void
{
if (previous_cell.status() != valid)
{
// previous_cell has not been set. set it now, and register with the
// triangulation that we want to be informed about mesh refinement
previous_cell = current_cell;
previous_cell->get_triangulation().signals.post_refinement.connect(
[this]()
{
this->invalidate_previous_cell();
});
}
else
previous_cell = current_cell;
current_cell = cell;
// ... do something with the cell...
}
template <int dim>
void
{
}
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
TriaActiveIterator< CellAccessor< dim, spacedim > > active_cell_iterator
Definition tria.h:1594

Here, whenever the triangulation is refined, it triggers the post-refinement signal which calls the function object attached to it. This function object is the member function FEValues<dim>::invalidate_previous_cell where we have bound the single argument (the this pointer of a member function that otherwise takes no arguments) to the this pointer of the FEValues object. Note how here there is no need for the code that owns the triangulation and the FEValues object to inform the latter if the former is refined. (In practice, the function would want to connect to some of the other signals that the triangulation offers as well, in particular to creation and deletion signals.)

The Triangulation class has a variety of signals that indicate different actions by which the triangulation can modify itself and potentially require follow-up action elsewhere. Please refer to Triangulation::Signals for details.

Serializing (loading or storing) triangulations

Like many other classes in deal.II, the Triangulation class can stream its contents to an archive using BOOST's serialization facilities. The data so stored can later be retrieved again from the archive to restore the contents of this object. This facility is frequently used to save the state of a program to disk for possible later resurrection, often in the context of checkpoint/restart strategies for long running computations or on computers that aren't very reliable (e.g. on very large clusters where individual nodes occasionally fail and then bring down an entire MPI job).

For technical reasons, writing and restoring a Triangulation object is not trivial. The primary reason is that unlike many other objects, triangulations rely on many other objects to which they store pointers or with which they interface; for example, triangulations store pointers to objects describing boundaries and manifolds, and they have signals that store pointers to other objects so they can be notified of changes in the triangulation (see the section on signals in this introduction). Since these objects are owned by the user space (for example the user can create a custom manifold object), they may not be serializable. So in cases like this, boost::serialize can store a reference to an object instead of the pointer, but the reference will never be satisfied at write time because the object pointed to is not serialized. Clearly, at load time, boost::serialize will not know where to let the pointer point to because it never gets to re-create the object originally pointed to.

For these reasons, saving a triangulation to an archive does not store all information, but only certain parts. More specifically, the information that is stored is everything that defines the mesh such as vertex locations, vertex indices, how vertices are connected to cells, boundary indicators, subdomain ids, material ids, etc. On the other hand, the following information is not stored:

On the other hand, since these are objects that are usually set in user code, they can typically easily be set again in that part of your code in which you re-load triangulations.

In a sense, this approach to serialization means that re-loading a triangulation is more akin to calling the Triangulation::create_triangulation() function and filling it with some additional content, as that function also does not touch the signals and Manifold objects that belong to this triangulation. In keeping with this analogy, the Triangulation::load() function also triggers the same kinds of signal as Triangulation::create_triangulation().

Technical details

Algorithms for mesh regularization and smoothing upon refinement

We chose an inductive point of view: since upon creation of the triangulation all cells are on the same level, all regularity assumptions regarding the maximum difference in level of cells sharing a common face, edge or vertex hold. Since we use the regularization and smoothing in each step of the mesh history, when coming to the point of refining it further the assumptions also hold.

The regularization and smoothing is done in the prepare_coarsening_and_refinement function, which is called by execute_coarsening_and_refinement at the very beginning. It decides which additional cells to flag for refinement by looking at the old grid and the refinement flags for each cell.

Regularization and smoothing are a bit complementary in that we check whether we need to set additional refinement flags when being on a cell flagged for refinement (regularization) or on a cell not flagged for refinement. This makes readable programming easier.

All the described algorithms apply only for more than one space dimension, since for one dimension no restrictions apply. It may be necessary to apply some smoothing for multigrid algorithms, but this has to be decided upon later.

Warning

It seems impossible to preserve constness of a triangulation through iterator usage. Thus, if you declare pointers to a const triangulation object, you should be well aware that you might involuntarily alter the data stored in the triangulation.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

Definition at line 1336 of file tria.h.

Member Typedef Documentation

◆ IteratorSelector

template<int dim, int spacedim = dim>
using Triangulation< dim, spacedim >::IteratorSelector = ::internal::TriangulationImplementation::Iterators<dim, spacedim>
private

An internal alias to make the definition of the iterator classes simpler.

Definition at line 1343 of file tria.h.

◆ level_cell_iterator

template<int dim, int spacedim = dim>
using Triangulation< dim, spacedim >::level_cell_iterator = cell_iterator

The same as above to allow the usage of the "MeshType concept" also on the refinement levels.

Definition at line 1576 of file tria.h.

◆ CellStatus

template<int dim, int spacedim = dim>
using Triangulation< dim, spacedim >::CellStatus = ::CellStatus

The elements of this enum are used to inform functions how a specific cell is going to change. This is used in the course of transferring data from one mesh to a refined or coarsened version of the mesh, for example. Note that this may me different than the refine_flag() and coarsen_flag() set on a cell, for example in parallel calculations, because of refinement constraints that an individual machine does not see.

Deprecated:
This is an alias for backward compatibility. Use CellStatus directly.

Definition at line 2244 of file tria.h.

◆ raw_cell_iterator

template<int dim, int spacedim = dim>
using Triangulation< dim, spacedim >::raw_cell_iterator = TriaRawIterator<CellAccessor<dim, spacedim> >
private

Declare a number of iterator types for raw iterators, i.e., iterators that also iterate over holes in the list of cells left by cells that have been coarsened away in previous mesh refinement cycles.

Since users should never have to access these internal properties of how we store data, these iterator types are made private.

Definition at line 4103 of file tria.h.

◆ raw_face_iterator

template<int dim, int spacedim = dim>
using Triangulation< dim, spacedim >::raw_face_iterator = TriaRawIterator<TriaAccessor<dim - 1, dim, spacedim> >
private

Definition at line 4104 of file tria.h.

◆ raw_vertex_iterator

template<int dim, int spacedim = dim>
using Triangulation< dim, spacedim >::raw_vertex_iterator = TriaRawIterator<::TriaAccessor<0, dim, spacedim> >
private

Definition at line 4106 of file tria.h.

◆ raw_line_iterator

template<int dim, int spacedim = dim>
using Triangulation< dim, spacedim >::raw_line_iterator = typename IteratorSelector::raw_line_iterator
private

Definition at line 4108 of file tria.h.

◆ raw_quad_iterator

template<int dim, int spacedim = dim>
using Triangulation< dim, spacedim >::raw_quad_iterator = typename IteratorSelector::raw_quad_iterator
private

Definition at line 4109 of file tria.h.

◆ raw_hex_iterator

template<int dim, int spacedim = dim>
using Triangulation< dim, spacedim >::raw_hex_iterator = typename IteratorSelector::raw_hex_iterator
private

Definition at line 4110 of file tria.h.

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 229 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 234 of file subscriptor.h.

Member Enumeration Documentation

◆ MeshSmoothing

template<int dim, int spacedim = dim>
enum Triangulation::MeshSmoothing

Declare some symbolic names for mesh smoothing algorithms. The meaning of these flags is documented in the Triangulation class.

Enumerator
none 

No mesh smoothing at all, except that meshes have to remain one-irregular.

limit_level_difference_at_vertices 

It can be shown, that degradation of approximation occurs if the triangulation contains vertices which are member of cells with levels differing by more than one. One such example is the following:

It would seem that in two space dimensions, the maximum jump in levels between cells sharing a common vertex is two (as in the example above). However, this is not true if more than four cells meet at a vertex. It is not uncommon that a coarse (initial) mesh contains vertices at which six or even eight cells meet, when small features of the domain have to be resolved even on the coarsest mesh. In that case, the maximum difference in levels is three or four, respectively. The problem gets even worse in three space dimensions.

Looking at an interpolation of the second derivative of the finite element solution (assuming bilinear finite elements), one sees that the numerical solution is almost totally wrong, compared with the true second derivative. Indeed, on regular meshes, there exist sharp estimations that the H2-error is only of order one, so we should not be surprised; however, the numerical solution may show a value for the second derivative which may be a factor of ten away from the true value. These problems are located on the small cell adjacent to the center vertex, where cells of non-subsequent levels meet, as well as on the upper and right neighbor of this cell (but with a less degree of deviation from the true value).

If the smoothing indicator given to the constructor contains the bit for limit_level_difference_at_vertices, situations as the above one are eliminated by also marking the upper right cell for refinement.

In case of anisotropic refinement, the level of a cell is not linked to the refinement of a cell as directly as in case of isotropic refinement. Furthermore, a cell can be strongly refined in one direction and not or at least much less refined in another. Therefore, it is very difficult to decide, which cases should be excluded from the refinement process. As a consequence, when using anisotropic refinement, the limit_level_difference_at_vertices flag must not be set. On the other hand, the implementation of multigrid methods in deal.II requires that this bit be set.

eliminate_unrefined_islands 

Single cells which are not refined and are surrounded by cells which are refined usually also lead to a sharp decline in approximation properties locally. The reason is that the nodes on the faces between unrefined and refined cells are not real degrees of freedom but carry constraints. The patch without additional degrees of freedom is thus significantly larger then the unrefined cell itself. If in the parameter passed to the constructor the bit for eliminate_unrefined_islands is set, all cells which are not flagged for refinement but which are surrounded by more refined cells than unrefined cells are flagged for refinement. Cells which are not yet refined but flagged for that are accounted for the number of refined neighbors. Cells on the boundary are not accounted for at all. An unrefined island is, by this definition also a cell which (in 2d) is surrounded by three refined cells and one unrefined one, or one surrounded by two refined cells, one unrefined one and is at the boundary on one side. It is thus not a true island, as the name of the flag may indicate. However, no better name came to mind to the author by now.

patch_level_1 

A triangulation of patch level 1 consists of patches, i.e. of cells that are refined once. This flag ensures that a mesh of patch level 1 is still of patch level 1 after coarsening and refinement. It is, however, the user's responsibility to ensure that the mesh is of patch level 1 before calling Triangulation::execute_coarsening_and_refinement() the first time. The easiest way to achieve this is by calling global_refine(1) straight after creation of the triangulation. It follows that if at least one of the children of a cell is or will be refined than all children need to be refined. If the patch_level_1 flag is set, than the flags eliminate_unrefined_islands, eliminate_refined_inner_islands and eliminate_refined_boundary_islands will be ignored as they will be fulfilled automatically.

coarsest_level_1 

Each coarse grid cell is refined at least once, i.e., the triangulation might have active cells on level 1 but not on level 0. This flag ensures that a mesh which has coarsest_level_1 has still coarsest_level_1 after coarsening and refinement. It is, however, the user's responsibility to ensure that the mesh has coarsest_level_1 before calling execute_coarsening_and_refinement the first time. The easiest way to achieve this is by calling global_refine(1) straight after creation of the triangulation. It follows that active cells on level 1 may not be coarsened.

The main use of this flag is to ensure that each cell has at least one neighbor in each coordinate direction (i.e. each cell has at least a left or right, and at least an upper or lower neighbor in 2d). This is a necessary precondition for some algorithms that compute finite differences between cells. The DerivativeApproximation class is one of these algorithms that require that a triangulation is coarsest_level_1 unless all cells already have at least one neighbor in each coordinate direction on the coarsest level.

allow_anisotropic_smoothing 

This flag is not included in maximum_smoothing. The flag is concerned with the following case: consider the case that an unrefined and a refined cell share a common face and that one of the children of the refined cell along the common face is flagged for further refinement. In that case, the resulting mesh would have more than one hanging node along one or more of the edges of the triangulation, a situation that is not allowed. Consequently, in order to perform the refinement, the coarser of the two original cells is also going to be refined.

However, in many cases it is sufficient to refine the coarser of the two original cells in an anisotropic way to avoid the case of multiple hanging vertices on a single edge. Doing only the minimal anisotropic refinement can save cells and degrees of freedom. By specifying this flag, the library can produce these anisotropic refinements.

The flag is not included by default since it may lead to anisotropically refined meshes even though no cell has ever been refined anisotropically explicitly by a user command. This surprising fact may lead to programs that do the wrong thing since they are not written for the additional cases that can happen with anisotropic meshes, see the discussion in the introduction to step-30.

eliminate_refined_inner_islands 

This algorithm seeks for isolated cells which are refined or flagged for refinement. This definition is unlike that for eliminate_unrefined_islands, which would mean that an island is defined as a cell which is refined but more of its neighbors are not refined than are refined. For example, in 2d, a cell's refinement would be reverted if at most one of its neighbors is also refined (or refined but flagged for coarsening).

The reason for the change in definition of an island is, that this option would be a bit dangerous, since if you consider a chain of refined cells (e.g. along a kink in the solution), the cells at the two ends would be coarsened, after which the next outermost cells would need to be coarsened. Therefore, only one loop of flagging cells like this could be done to avoid eating up the whole chain of refined cells (‘chain reaction’...).

This algorithm also takes into account cells which are not actually refined but are flagged for refinement. If necessary, it takes away the refinement flag.

Actually there are two versions of this flag, eliminate_refined_inner_islands and eliminate_refined_boundary_islands. The first eliminates islands defined by the definition above which are in the interior of the domain, while the second eliminates only those islands if the cell is at the boundary. The reason for this split of flags is that one often wants to eliminate such islands in the interior while those at the boundary may well be wanted, for example if one refines the mesh according to a criterion associated with a boundary integral or if one has rough boundary data.

eliminate_refined_boundary_islands 

The result of this flag is very similar to eliminate_refined_inner_islands. See the documentation there.

do_not_produce_unrefined_islands 

This flag prevents the occurrence of unrefined islands. In more detail: It prohibits the coarsening of a cell if 'most of the neighbors' will be refined after the step.

smoothing_on_refinement 

This flag sums up all smoothing algorithms which may be performed upon refinement by flagging some more cells for refinement.

smoothing_on_coarsening 

This flag sums up all smoothing algorithms which may be performed upon coarsening by flagging some more cells for coarsening.

maximum_smoothing 

This flag includes all the above ones (therefore combines all smoothing algorithms implemented), with the exception of anisotropic smoothing.

Definition at line 1351 of file tria.h.

Constructor & Destructor Documentation

◆ Triangulation() [1/3]

template<int dim, int spacedim = dim>
Triangulation< dim, spacedim >::Triangulation ( const MeshSmoothing  smooth_grid = none,
const bool  check_for_distorted_cells = false 
)

Create an empty triangulation. Do not create any cells.

Parameters
smooth_gridDetermines the level of smoothness of the mesh size function that should be enforced upon mesh refinement.
check_for_distorted_cellsDetermines whether the triangulation should check whether any of the cells that are created by create_triangulation() or execute_coarsening_and_refinement() are distorted (see distorted cells). If set, these two functions may throw an exception if they encounter distorted cells.

◆ Triangulation() [2/3]

template<int dim, int spacedim = dim>
Triangulation< dim, spacedim >::Triangulation ( const Triangulation< dim, spacedim > &  )
delete

Copy constructor.

You should really use the copy_triangulation function, so this constructor is deleted. The reason for this is that we may want to use triangulation objects in collections. However, C++ containers require that the objects stored in them are copyable, so we need to provide a copy constructor. On the other hand, copying triangulations is so expensive that we do not want such objects copied by accident, for example in compiler-generated temporary objects. By defining a copy constructor but throwing an error, we satisfy the formal requirements of containers, but at the same time disallow actual copies. Finally, through the exception, one easily finds the places where code has to be changed to avoid copies.

◆ Triangulation() [3/3]

template<int dim, int spacedim = dim>
Triangulation< dim, spacedim >::Triangulation ( Triangulation< dim, spacedim > &&  tria)
noexcept

Move constructor.

Create a new triangulation by stealing the internal data of another triangulation.

◆ ~Triangulation()

template<int dim, int spacedim = dim>
virtual Triangulation< dim, spacedim >::~Triangulation ( )
overridevirtual

Member Function Documentation

◆ operator=()

template<int dim, int spacedim = dim>
Triangulation & Triangulation< dim, spacedim >::operator= ( Triangulation< dim, spacedim > &&  tria)
noexcept

Move assignment operator.

◆ clear()

template<int dim, int spacedim = dim>
virtual void Triangulation< dim, spacedim >::clear ( )
virtual

Reset this triangulation into a virgin state by deleting all data.

Note that this operation is only allowed if no subscriptions to this object exist any more, such as DoFHandler objects using it.

Reimplemented in parallel::distributed::Triangulation< dim, spacedim >, parallel::DistributedTriangulationBase< dim, spacedim >, parallel::DistributedTriangulationBase< 1, spacedim >, and parallel::DistributedTriangulationBase< dim, dim >.

◆ get_communicator()

template<int dim, int spacedim = dim>
virtual MPI_Comm Triangulation< dim, spacedim >::get_communicator ( ) const
virtual

Return MPI communicator used by this triangulation. In the case of a serial Triangulation object, MPI_COMM_SELF is returned.

Reimplemented in parallel::TriangulationBase< dim, spacedim >, parallel::TriangulationBase< dim, dim >, and parallel::TriangulationBase< dim, spacedim >.

◆ global_active_cell_index_partitioner()

template<int dim, int spacedim = dim>
virtual std::weak_ptr< const Utilities::MPI::Partitioner > Triangulation< dim, spacedim >::global_active_cell_index_partitioner ( ) const
virtual

Return the partitioner for the global indices of the cells on the active level of the triangulation, which is returned by the function CellAccessor::global_active_cell_index().

Reimplemented in parallel::TriangulationBase< dim, spacedim >, parallel::TriangulationBase< dim, dim >, and parallel::TriangulationBase< dim, spacedim >.

◆ global_level_cell_index_partitioner()

template<int dim, int spacedim = dim>
virtual std::weak_ptr< const Utilities::MPI::Partitioner > Triangulation< dim, spacedim >::global_level_cell_index_partitioner ( const unsigned int  level) const
virtual

Return the partitioner for the global indices of the cells on the given level of the triangulation, which is returned by the function CellAccessor::global_level_cell_index().

Reimplemented in parallel::TriangulationBase< dim, spacedim >, parallel::TriangulationBase< dim, dim >, and parallel::TriangulationBase< dim, spacedim >.

◆ set_mesh_smoothing()

template<int dim, int spacedim = dim>
virtual void Triangulation< dim, spacedim >::set_mesh_smoothing ( const MeshSmoothing  mesh_smoothing)
virtual

Set the mesh smoothing to mesh_smoothing. This overrides the MeshSmoothing given to the constructor.

◆ get_mesh_smoothing()

template<int dim, int spacedim = dim>
virtual const MeshSmoothing & Triangulation< dim, spacedim >::get_mesh_smoothing ( ) const
virtual

Return the mesh smoothing requirements that are obeyed.

◆ get_boundary_ids()

template<int dim, int spacedim = dim>
virtual std::vector< types::boundary_id > Triangulation< dim, spacedim >::get_boundary_ids ( ) const
virtual

Return a vector containing all boundary indicators assigned to boundary faces of active cells of this Triangulation object. Note, that each boundary indicator is reported only once. The size of the return vector will represent the number of different indicators (which is greater or equal one).

See also
Glossary entry on boundary indicators

Reimplemented in parallel::TriangulationBase< dim, spacedim >, parallel::TriangulationBase< dim, dim >, and parallel::TriangulationBase< dim, spacedim >.

◆ copy_triangulation()

template<int dim, int spacedim = dim>
virtual void Triangulation< dim, spacedim >::copy_triangulation ( const Triangulation< dim, spacedim > &  other_tria)
virtual

Copy other_tria to this triangulation. This operation is not cheap, so you should be careful with using this. We do not implement this function as a copy constructor, since it makes it easier to maintain collections of triangulations if you can assign them values later on.

This triangulation must be empty beforehand.

The function is made virtual since some derived classes might want to disable or extend the functionality of this function.

Note
Calling this function triggers the 'copy' signal on other_tria, i.e. the triangulation being copied from. It also triggers the 'create' signal of the current triangulation. See the section on signals in the general documentation for more information.
The list of connections to signals is not copied from the old to the new triangulation since these connections were established to monitor how the old triangulation changes, not how any triangulation it may be copied to changes.

Reimplemented in PersistentTriangulation< dim, spacedim >.

◆ create_triangulation() [1/2]

template<int dim, int spacedim = dim>
virtual void Triangulation< dim, spacedim >::create_triangulation ( const std::vector< Point< spacedim > > &  vertices,
const std::vector< CellData< dim > > &  cells,
const SubCellData subcelldata 
)
virtual

Create a triangulation from a list of vertices and a list of cells, each of the latter being a list of 1<<dim vertex indices. The triangulation must be empty upon calling this function and the cell list should be useful (connected domain, etc.). The result of calling this function is a coarse mesh.

Material data for the cells is given within the cells array, while boundary information is given in the subcelldata field.

The numbering of vertices within the cells array is subject to some constraints; see the general class documentation for this.

For conditions when this function can generate a valid triangulation, see the documentation of this class, and the GridIn and GridTools::consistently_order_cells() function.

If the check_for_distorted_cells flag was specified upon creation of this object, at the very end of its operation, the current function walks over all cells and verifies that none of the cells is deformed (see the entry on distorted cells in the glossary), where we call a cell deformed if the determinant of the Jacobian of the mapping from reference cell to real cell is negative at least at one of the vertices (this computation is done using the GeometryInfo::jacobian_determinants_at_vertices function). If there are deformed cells, this function throws an exception of kind DistortedCellList. Since this happens after all data structures have been set up, you can catch and ignore this exception if you know what you do – for example, it may be that the determinant is zero (indicating that you have collapsed edges in a cell) but that this is ok because you didn't intend to integrate on this cell anyway. On the other hand, deformed cells are often a sign of a mesh that is too coarse to resolve the geometry of the domain, and in this case ignoring the exception is probably unwise.

Note
This function is used in step-14 and step-19.
This function triggers the "create" signal after doing its work. See the section on signals in the general documentation of this class. For example as a consequence of this, all DoFHandler objects connected to this triangulation will be reinitialized via DoFHandler::reinit().
The check for distorted cells is only done if dim==spacedim, as otherwise cells can legitimately be twisted if the manifold they describe is twisted.

Reimplemented in parallel::shared::Triangulation< dim, spacedim >, parallel::distributed::Triangulation< dim, spacedim >, and PersistentTriangulation< dim, spacedim >.

◆ create_triangulation() [2/2]

template<int dim, int spacedim = dim>
virtual void Triangulation< dim, spacedim >::create_triangulation ( const TriangulationDescription::Description< dim, spacedim > &  construction_data)
virtual

Create a triangulation from the provided TriangulationDescription::Description.

Note
Don't forget to attach the manifolds with set_manifold() before calling this function if manifolds are needed.
The namespace TriangulationDescription::Utilities contains functions to create TriangulationDescription::Description.
Parameters
construction_dataThe data needed for this process.

Reimplemented in parallel::fullydistributed::Triangulation< dim, spacedim >, parallel::shared::Triangulation< dim, spacedim >, parallel::distributed::Triangulation< dim, spacedim >, and PersistentTriangulation< dim, spacedim >.

◆ flip_all_direction_flags()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::flip_all_direction_flags ( )

Revert or flip the direction flags of a triangulation with dim==spacedim-1, see GlossDirectionFlag.

This function throws an exception if dim==spacedim or if dim<spacedim-1.

◆ set_all_refine_flags()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::set_all_refine_flags ( )

Flag all active cells for refinement. This will refine all cells of all levels which are not already refined (i.e. only cells are refined which do not yet have children). The cells are only flagged, not refined, thus you have the chance to save the refinement flags.

◆ refine_global()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::refine_global ( const unsigned int  times = 1)

Refine all cells times times. In other words, in each one of the times iterations, loop over all cells and refine each cell uniformly into \(2^\text{dim}\) children. In practice, this function repeats the following operations times times: call set_all_refine_flags() followed by execute_coarsening_and_refinement(). The end result is that the number of cells increases by a factor of \((2^\text{dim})^\text{times}=2^{\text{dim} \times \text{times}}\).

The execute_coarsening_and_refinement() function called in this loop may throw an exception if it creates cells that are distorted (see its documentation for an explanation). This exception will be propagated through this function if that happens, and you may not get the actual number of refinement steps in that case.

Note
This function triggers the pre- and post-refinement signals before and after doing each individual refinement cycle (i.e. more than once if times > 1) . See the section on signals in the general documentation of this class.

◆ coarsen_global()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::coarsen_global ( const unsigned int  times = 1)

Coarsen all cells the given number of times.

In each of one of the times iterations, all cells will be marked for coarsening. If an active cell is already on the coarsest level, it will be ignored.

Note
This function triggers the pre- and post-refinement signals before and after doing each individual coarsening cycle (i.e. more than once if times > 1) . See the section on signals in the general documentation of this class.

◆ execute_coarsening_and_refinement()

template<int dim, int spacedim = dim>
virtual void Triangulation< dim, spacedim >::execute_coarsening_and_refinement ( )
virtual

Execute both refinement and coarsening of the triangulation.

The function resets all refinement and coarsening flags to false. It uses the user flags for internal purposes. They will therefore be overwritten by undefined content.

To allow user programs to fix up these cells if that is desired, this function after completing all other work may throw an exception of type DistortedCellList that contains a list of those cells that have been refined and have at least one child that is distorted. The function does not create such an exception if no cells have created distorted children. Note that for the check for distorted cells to happen, the check_for_distorted_cells flag has to be specified upon creation of a triangulation object.

See the general docs for more information.

Note
This function triggers the pre- and post-refinement signals before and after doing its work. See the section on signals in the general documentation of this class.
If the boundary description is sufficiently irregular, it can happen that some of the children produced by mesh refinement are distorted (see the extensive discussion on distorted cells).
This function is virtual to allow derived classes to insert hooks, such as saving refinement flags and the like (see e.g. the PersistentTriangulation class).

Reimplemented in parallel::fullydistributed::Triangulation< dim, spacedim >, parallel::shared::Triangulation< dim, spacedim >, parallel::distributed::Triangulation< dim, spacedim >, and PersistentTriangulation< dim, spacedim >.

◆ prepare_coarsening_and_refinement() [1/4]

template<int dim, int spacedim = dim>
virtual bool Triangulation< dim, spacedim >::prepare_coarsening_and_refinement ( )
virtual

Do both preparation for refinement and coarsening as well as mesh smoothing.

Regarding the refinement process it fixes the closure of the refinement in dim>=2 (make sure that no two cells are adjacent with a refinement level differing with more than one), etc. It performs some mesh smoothing if the according flag was given to the constructor of this class. The function returns whether additional cells have been flagged for refinement.

See the general doc of this class for more information on smoothing upon refinement.

Regarding the coarsening part, flagging and deflagging cells in preparation of the actual coarsening step are done. This includes deleting coarsen flags from cells which may not be deleted (e.g. because one neighbor is more refined than the cell), doing some smoothing, etc.

The effect is that only those cells are flagged for coarsening which will actually be coarsened. This includes the fact that all flagged cells belong to parent cells of which all children are flagged.

The function returns whether some cells' flagging has been changed in the process.

This function uses the user flags, so store them if you still need them afterwards.

Reimplemented in parallel::fullydistributed::Triangulation< dim, spacedim >, and parallel::distributed::Triangulation< dim, spacedim >.

◆ save_refine_flags() [1/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_refine_flags ( std::ostream &  out) const

Save the addresses of the cells which are flagged for refinement to out. For usage, read the general documentation for this class.

◆ save_refine_flags() [2/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_refine_flags ( std::vector< bool > &  v) const

Same as above, but store the flags to a bitvector rather than to a file.

◆ load_refine_flags() [1/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_refine_flags ( std::istream &  in)

Read the information stored by save_refine_flags.

◆ load_refine_flags() [2/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_refine_flags ( const std::vector< bool > &  v)

Read the information stored by save_refine_flags.

◆ save_coarsen_flags() [1/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_coarsen_flags ( std::ostream &  out) const

Analogue to save_refine_flags.

◆ save_coarsen_flags() [2/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_coarsen_flags ( std::vector< bool > &  v) const

Same as above, but store the flags to a bitvector rather than to a file.

◆ load_coarsen_flags() [1/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_coarsen_flags ( std::istream &  out)

Analogue to load_refine_flags.

◆ load_coarsen_flags() [2/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_coarsen_flags ( const std::vector< bool > &  v)

Analogue to load_refine_flags.

◆ get_anisotropic_refinement_flag()

template<int dim, int spacedim = dim>
bool Triangulation< dim, spacedim >::get_anisotropic_refinement_flag ( ) const

Return whether this triangulation has ever undergone anisotropic (as opposed to only isotropic) refinement.

◆ clear_user_flags()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::clear_user_flags ( )

Clear all user flags. See also GlossUserFlags.

◆ save_user_flags() [1/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_flags ( std::ostream &  out) const

Save all user flags. See the general documentation for this class and the documentation for the save_refine_flags for more details. See also GlossUserFlags.

◆ save_user_flags() [2/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_flags ( std::vector< bool > &  v) const

Same as above, but store the flags to a bitvector rather than to a file. The output vector is resized if necessary. See also GlossUserFlags.

◆ load_user_flags() [1/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_flags ( std::istream &  in)

Read the information stored by save_user_flags. See also GlossUserFlags.

◆ load_user_flags() [2/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_flags ( const std::vector< bool > &  v)

Read the information stored by save_user_flags. See also GlossUserFlags.

◆ clear_user_flags_line()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::clear_user_flags_line ( )

Clear all user flags on lines. See also GlossUserFlags.

◆ save_user_flags_line() [1/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_flags_line ( std::ostream &  out) const

Save the user flags on lines. See also GlossUserFlags.

◆ save_user_flags_line() [2/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_flags_line ( std::vector< bool > &  v) const

Same as above, but store the flags to a bitvector rather than to a file. The output vector is resized if necessary. See also GlossUserFlags.

◆ load_user_flags_line() [1/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_flags_line ( std::istream &  in)

Load the user flags located on lines. See also GlossUserFlags.

◆ load_user_flags_line() [2/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_flags_line ( const std::vector< bool > &  v)

Load the user flags located on lines. See also GlossUserFlags.

◆ clear_user_flags_quad()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::clear_user_flags_quad ( )

Clear all user flags on quads. See also GlossUserFlags.

◆ save_user_flags_quad() [1/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_flags_quad ( std::ostream &  out) const

Save the user flags on quads. See also GlossUserFlags.

◆ save_user_flags_quad() [2/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_flags_quad ( std::vector< bool > &  v) const

Same as above, but store the flags to a bitvector rather than to a file. The output vector is resized if necessary. See also GlossUserFlags.

◆ load_user_flags_quad() [1/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_flags_quad ( std::istream &  in)

Load the user flags located on quads. See also GlossUserFlags.

◆ load_user_flags_quad() [2/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_flags_quad ( const std::vector< bool > &  v)

Load the user flags located on quads. See also GlossUserFlags.

◆ clear_user_flags_hex()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::clear_user_flags_hex ( )

Clear all user flags on quads. See also GlossUserFlags.

◆ save_user_flags_hex() [1/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_flags_hex ( std::ostream &  out) const

Save the user flags on hexs. See also GlossUserFlags.

◆ save_user_flags_hex() [2/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_flags_hex ( std::vector< bool > &  v) const

Same as above, but store the flags to a bitvector rather than to a file. The output vector is resized if necessary. See also GlossUserFlags.

◆ load_user_flags_hex() [1/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_flags_hex ( std::istream &  in)

Load the user flags located on hexs. See also GlossUserFlags.

◆ load_user_flags_hex() [2/2]

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_flags_hex ( const std::vector< bool > &  v)

Load the user flags located on hexs. See also GlossUserFlags.

◆ clear_user_data()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::clear_user_data ( )

Clear all user pointers and indices and allow the use of both for next access. See also GlossUserData.

◆ save_user_indices()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_indices ( std::vector< unsigned int > &  v) const

Save all user indices. The output vector is resized if necessary. See also GlossUserData.

◆ load_user_indices()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_indices ( const std::vector< unsigned int > &  v)

Read the information stored by save_user_indices(). See also GlossUserData.

◆ save_user_pointers()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_pointers ( std::vector< void * > &  v) const

Save all user pointers. The output vector is resized if necessary. See also GlossUserData.

◆ load_user_pointers()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_pointers ( const std::vector< void * > &  v)

Read the information stored by save_user_pointers(). See also GlossUserData.

◆ save_user_indices_line()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_indices_line ( std::vector< unsigned int > &  v) const

Save the user indices on lines. The output vector is resized if necessary. See also GlossUserData.

◆ load_user_indices_line()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_indices_line ( const std::vector< unsigned int > &  v)

Load the user indices located on lines. See also GlossUserData.

◆ save_user_indices_quad()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_indices_quad ( std::vector< unsigned int > &  v) const

Save the user indices on quads. The output vector is resized if necessary. See also GlossUserData.

◆ load_user_indices_quad()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_indices_quad ( const std::vector< unsigned int > &  v)

Load the user indices located on quads. See also GlossUserData.

◆ save_user_indices_hex()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_indices_hex ( std::vector< unsigned int > &  v) const

Save the user indices on hexes. The output vector is resized if necessary. See also GlossUserData.

◆ load_user_indices_hex()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_indices_hex ( const std::vector< unsigned int > &  v)

Load the user indices located on hexs. See also GlossUserData.

◆ save_user_pointers_line()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_pointers_line ( std::vector< void * > &  v) const

Save the user indices on lines. The output vector is resized if necessary. See also GlossUserData.

◆ load_user_pointers_line()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_pointers_line ( const std::vector< void * > &  v)

Load the user pointers located on lines. See also GlossUserData.

◆ save_user_pointers_quad()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_pointers_quad ( std::vector< void * > &  v) const

Save the user pointers on quads. The output vector is resized if necessary. See also GlossUserData.

◆ load_user_pointers_quad()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_pointers_quad ( const std::vector< void * > &  v)

Load the user pointers located on quads. See also GlossUserData.

◆ save_user_pointers_hex()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_user_pointers_hex ( std::vector< void * > &  v) const

Save the user pointers on hexes. The output vector is resized if necessary. See also GlossUserData.

◆ load_user_pointers_hex()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_user_pointers_hex ( const std::vector< void * > &  v)

Load the user pointers located on hexs. See also GlossUserData.

◆ begin()

template<int dim, int spacedim = dim>
cell_iterator Triangulation< dim, spacedim >::begin ( const unsigned int  level = 0) const

Iterator to the first used cell on level level.

Note
The given level argument needs to correspond to a level of the triangulation, i.e., should be less than the value returned by n_levels(). On the other hand, for parallel computations using a parallel::distributed::Triangulation object, it is often convenient to write loops over the cells of all levels of the global mesh, even if the local portion of the triangulation does not actually have cells at one of the higher levels. In those cases, the level argument is accepted if it is less than what the n_global_levels() function returns. If the given level is between the values returned by n_levels() and n_global_levels(), then no cells exist in the local portion of the triangulation at this level, and the function simply returns what end() would return.

◆ begin_active()

template<int dim, int spacedim = dim>
active_cell_iterator Triangulation< dim, spacedim >::begin_active ( const unsigned int  level = 0) const

Iterator to the first active cell on level level. If the given level does not contain any active cells (i.e., all cells on this level are further refined, then this function returns end_active(level) so that loops of the kind

for (const auto cell=tria.begin_active(level);
cell!=tria.end_active(level);
++cell)
{
...
}
active_cell_iterator end_active(const unsigned int level) const

have zero iterations, as may be expected if there are no active cells on this level.

Note
The given level argument needs to correspond to a level of the triangulation, i.e., should be less than the value returned by n_levels(). On the other hand, for parallel computations using a parallel::distributed::Triangulation object, it is often convenient to write loops over the cells of all levels of the global mesh, even if the local portion of the triangulation does not actually have cells at one of the higher levels. In those cases, the level argument is accepted if it is less than what the n_global_levels() function returns. If the given level is between the values returned by n_levels() and n_global_levels(), then no cells exist in the local portion of the triangulation at this level, and the function simply returns what end() would return.

◆ end() [1/2]

template<int dim, int spacedim = dim>
cell_iterator Triangulation< dim, spacedim >::end ( ) const

Iterator past the end; this iterator serves for comparisons of iterators with past-the-end or before-the-beginning states.

◆ end() [2/2]

template<int dim, int spacedim = dim>
cell_iterator Triangulation< dim, spacedim >::end ( const unsigned int  level) const

Return an iterator which is the first iterator not on level. If level is the last level, then this returns end().

Note
The given level argument needs to correspond to a level of the triangulation, i.e., should be less than the value returned by n_levels(). On the other hand, for parallel computations using a parallel::distributed::Triangulation object, it is often convenient to write loops over the cells of all levels of the global mesh, even if the local portion of the triangulation does not actually have cells at one of the higher levels. In those cases, the level argument is accepted if it is less than what the n_global_levels() function returns. If the given level is between the values returned by n_levels() and n_global_levels(), then no cells exist in the local portion of the triangulation at this level, and the function simply returns what end() would return.

◆ end_active()

template<int dim, int spacedim = dim>
active_cell_iterator Triangulation< dim, spacedim >::end_active ( const unsigned int  level) const

Return an active iterator which is the first active iterator not on the given level. If level is the last level, then this returns end().

Note
The given level argument needs to correspond to a level of the triangulation, i.e., should be less than the value returned by n_levels(). On the other hand, for parallel computations using a parallel::distributed::Triangulation object, it is often convenient to write loops over the cells of all levels of the global mesh, even if the local portion of the triangulation does not actually have cells at one of the higher levels. In those cases, the level argument is accepted if it is less than what the n_global_levels() function returns. If the given level is between the values returned by n_levels() and n_global_levels(), then no cells exist in the local portion of the triangulation at this level, and the function simply returns what end() would return.

◆ last()

template<int dim, int spacedim = dim>
cell_iterator Triangulation< dim, spacedim >::last ( ) const

Return an iterator pointing to the last used cell.

◆ last_active()

template<int dim, int spacedim = dim>
active_cell_iterator Triangulation< dim, spacedim >::last_active ( ) const

Return an iterator pointing to the last active cell.

◆ create_cell_iterator()

template<int dim, int spacedim = dim>
cell_iterator Triangulation< dim, spacedim >::create_cell_iterator ( const CellId cell_id) const

Return an iterator to a cell of this Triangulation object constructed from an independent CellId object.

Note
See the documentation of contains_cell() about which CellId objects are valid.

◆ contains_cell()

template<int dim, int spacedim = dim>
bool Triangulation< dim, spacedim >::contains_cell ( const CellId cell_id) const

Check if the triangulation contains a cell with the id cell_id. If the given argument corresponds to a valid cell in this triangulation, this operation will always return true for sequential triangulations where the current processor stores all cells that are part of the triangulation. On the other hand, if this is a parallel triangulation, then the current processor may not actually know about this cell. In this case, this operation will return true for locally relevant cells, but may return false for artificial cells that are less refined on the current processor.

◆ begin_face()

template<int dim, int spacedim = dim>
face_iterator Triangulation< dim, spacedim >::begin_face ( ) const

Iterator to the first used face.

◆ begin_active_face()

template<int dim, int spacedim = dim>
active_face_iterator Triangulation< dim, spacedim >::begin_active_face ( ) const

Iterator to the first active face.

◆ end_face()

template<int dim, int spacedim = dim>
face_iterator Triangulation< dim, spacedim >::end_face ( ) const

Iterator past the end; this iterator serves for comparisons of iterators with past-the-end or before-the-beginning states.

◆ begin_vertex()

template<int dim, int spacedim = dim>
vertex_iterator Triangulation< dim, spacedim >::begin_vertex ( ) const

Iterator to the first used vertex. This function can only be used if dim is not one.

◆ begin_active_vertex()

template<int dim, int spacedim = dim>
active_vertex_iterator Triangulation< dim, spacedim >::begin_active_vertex ( ) const

Iterator to the first active vertex. Because all vertices are active, begin_vertex() and begin_active_vertex() return the same vertex. This function can only be used if dim is not one.

◆ end_vertex()

template<int dim, int spacedim = dim>
vertex_iterator Triangulation< dim, spacedim >::end_vertex ( ) const

Iterator past the end; this iterator serves for comparisons of iterators with past-the-end or before-the-beginning states. This function can only be used if dim is not one.

◆ n_lines() [1/2]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_lines ( ) const

In the following, most functions are provided in two versions, with and without an argument describing the level. The versions with this argument are only applicable for objects describing the cells of the present triangulation. For example: in 2d n_lines(level) cannot be called, only n_lines(), as lines are faces in 2d and therefore have no level. Return the total number of used lines, active or not.

◆ n_lines() [2/2]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_lines ( const unsigned int  level) const

Return the total number of used lines, active or not on level level.

◆ n_active_lines() [1/2]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_active_lines ( ) const

Return the total number of active lines.

◆ n_active_lines() [2/2]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_active_lines ( const unsigned int  level) const

Return the total number of active lines, on level level.

◆ n_quads() [1/8]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_quads ( ) const

Return the total number of used quads, active or not.

◆ n_quads() [2/8]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_quads ( const unsigned int  level) const

Return the total number of used quads, active or not on level level.

◆ n_active_quads() [1/8]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_active_quads ( ) const

Return the total number of active quads, active or not.

◆ n_active_quads() [2/8]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_active_quads ( const unsigned int  level) const

Return the total number of active quads, active or not on level level.

◆ n_hexs() [1/4]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_hexs ( ) const

Return the total number of used hexahedra, active or not.

◆ n_hexs() [2/4]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_hexs ( const unsigned int  level) const

Return the total number of used hexahedra, active or not on level level.

◆ n_active_hexs() [1/4]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_active_hexs ( ) const

Return the total number of active hexahedra, active or not.

◆ n_active_hexs() [2/4]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_active_hexs ( const unsigned int  level) const

Return the total number of active hexahedra, active or not on level level.

◆ n_cells() [1/2]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_cells ( ) const

Return the total number of used cells, active or not. Maps to n_lines() in one space dimension and so on.

◆ n_cells() [2/2]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_cells ( const unsigned int  level) const

Return the total number of used cells, active or not, on level level. Maps to n_lines(level) in one space dimension and so on.

◆ n_active_cells() [1/2]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_active_cells ( ) const

Return the total number of active cells. Maps to n_active_lines() in one space dimension and so on.

◆ n_global_active_cells()

template<int dim, int spacedim = dim>
virtual types::global_cell_index Triangulation< dim, spacedim >::n_global_active_cells ( ) const
virtual

Return the total number of active cells. For the current class, this is the same as n_active_cells(). However, the function may be overloaded in derived classes (e.g., in parallel::distributed::Triangulation) where it may return a value greater than the number of active cells reported by the triangulation object on the current processor.

Reimplemented in parallel::TriangulationBase< dim, spacedim >, parallel::TriangulationBase< dim, dim >, and parallel::TriangulationBase< dim, spacedim >.

◆ n_active_cells() [2/2]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_active_cells ( const unsigned int  level) const

Return the total number of active cells on level level. Maps to n_active_lines(level) in one space dimension and so on.

◆ n_global_coarse_cells()

template<int dim, int spacedim = dim>
virtual types::coarse_cell_id Triangulation< dim, spacedim >::n_global_coarse_cells ( ) const
virtual

Return the total number of coarse cells. If the coarse mesh is replicated on each process, this simply returns n_cells(0).

Reimplemented in parallel::TriangulationBase< dim, spacedim >, parallel::TriangulationBase< dim, dim >, and parallel::TriangulationBase< dim, spacedim >.

◆ n_faces()

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_faces ( ) const

Return the total number of used faces, active or not. In 2d, the result equals n_lines(), in 3d it equals n_quads(), while in 1d it equals the number of used vertices.

◆ n_active_faces()

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_active_faces ( ) const

Return the total number of active faces. In 2d, the result equals n_active_lines(), in 3d it equals n_active_quads(), while in 1d it equals the number of used vertices.

◆ n_levels()

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_levels ( ) const

Return the number of levels in this triangulation.

Note
Internally, triangulations store data in levels, and there may be more levels in this data structure than one may think – for example, imagine a triangulation that we just got by coarsening the highest level so that it was completely depopulated. That level is not removed, since it will most likely be repopulated soon by the next refinement process. As a consequence, if you happened to run through raw cell iterators (which you can't do as a user of this class, but can internally), then the number of objects in the levels hierarchy is larger than the level of the most refined cell plus one. On the other hand, since this is rarely what a user of this class cares about, the function really just returns the level of the most refined active cell plus one. (The plus one is because in a coarse, unrefined mesh, all cells have level zero – making the number of levels equal to one.)

◆ n_global_levels()

template<int dim, int spacedim = dim>
virtual unsigned int Triangulation< dim, spacedim >::n_global_levels ( ) const
virtual

Return the number of levels in use. This function is equivalent to n_levels() for a serial Triangulation, but gives the maximum of n_levels() over all processors for a parallel::distributed::Triangulation and therefore can be larger than n_levels().

Reimplemented in parallel::TriangulationBase< dim, spacedim >, parallel::TriangulationBase< dim, dim >, and parallel::TriangulationBase< dim, spacedim >.

◆ has_hanging_nodes()

template<int dim, int spacedim = dim>
virtual bool Triangulation< dim, spacedim >::has_hanging_nodes ( ) const
virtual

Return true if the triangulation has hanging nodes.

The function is made virtual since the result can be interpreted in different ways, depending on whether the triangulation lives only on a single processor, or may be distributed as done in the parallel::distributed::Triangulation class (see there for a description of what the function is supposed to do in the parallel context).

Reimplemented in parallel::DistributedTriangulationBase< dim, spacedim >, parallel::DistributedTriangulationBase< 1, spacedim >, and parallel::DistributedTriangulationBase< dim, dim >.

◆ n_vertices()

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_vertices ( ) const

Return the total number of vertices. Some of them may not be used, which usually happens upon coarsening of a triangulation when some vertices are discarded, but we do not want to renumber the remaining ones, leading to holes in the numbers of used vertices. You can get the number of used vertices using n_used_vertices function.

◆ get_vertices()

template<int dim, int spacedim = dim>
const std::vector< Point< spacedim > > & Triangulation< dim, spacedim >::get_vertices ( ) const

Return a constant reference to all the vertices present in this triangulation. Note that not necessarily all vertices in this array are actually used; for example, if you coarsen a mesh, then some vertices are deleted, but their positions in this array are unchanged as the indices of vertices are only allocated once. You can find out about which vertices are actually used by the function get_used_vertices().

◆ n_used_vertices()

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_used_vertices ( ) const

Return the number of vertices that are presently in use, i.e. belong to at least one used element.

◆ vertex_used()

template<int dim, int spacedim = dim>
bool Triangulation< dim, spacedim >::vertex_used ( const unsigned int  index) const

Return true if the vertex with this index is used.

◆ get_used_vertices()

template<int dim, int spacedim = dim>
const std::vector< bool > & Triangulation< dim, spacedim >::get_used_vertices ( ) const

Return a constant reference to the array of bools indicating whether an entry in the vertex array is used or not.

◆ max_adjacent_cells() [1/4]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::max_adjacent_cells ( ) const

Return the maximum number of cells meeting at a common vertex. Since this number is an invariant under refinement, only the cells on the coarsest level are considered. The operation is thus reasonably fast. The invariance is only true for sufficiently many cells in the coarsest triangulation (e.g. for a single cell one would be returned), so a minimum of four is returned in two dimensions, 8 in three dimensions, etc, which is how many cells meet if the triangulation is refined.

In one space dimension, two is returned.

◆ locally_owned_subdomain()

template<int dim, int spacedim = dim>
virtual types::subdomain_id Triangulation< dim, spacedim >::locally_owned_subdomain ( ) const
virtual

This function always returns invalid_subdomain_id but is there for compatibility with the derived parallel::distributed::Triangulation class. For distributed parallel triangulations this function returns the subdomain id of those cells that are owned by the current processor.

Reimplemented in parallel::TriangulationBase< dim, spacedim >, parallel::TriangulationBase< dim, dim >, and parallel::TriangulationBase< dim, spacedim >.

◆ get_triangulation() [1/2]

template<int dim, int spacedim = dim>
Triangulation< dim, spacedim > & Triangulation< dim, spacedim >::get_triangulation ( )

Return a reference to the current object.

This doesn't seem to be very useful but allows to write code that can access the underlying triangulation for anything that satisfies the MeshType concept (which may not only be a triangulation, but also a DoFHandler, for example).

◆ get_triangulation() [2/2]

template<int dim, int spacedim = dim>
const Triangulation< dim, spacedim > & Triangulation< dim, spacedim >::get_triangulation ( ) const

Return a reference to the current object. This is the const-version of the previous function.

◆ n_raw_lines() [1/2]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_raw_lines ( ) const

Total number of lines, used or unused.

Note
This function really exports internal information about the triangulation. It shouldn't be used in applications. The function is only part of the public interface of this class because it is used in some of the other classes that build very closely on it (in particular, the DoFHandler class).

◆ n_raw_lines() [2/2]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_raw_lines ( const unsigned int  level) const

Number of lines, used or unused, on the given level.

Note
This function really exports internal information about the triangulation. It shouldn't be used in applications. The function is only part of the public interface of this class because it is used in some of the other classes that build very closely on it (in particular, the DoFHandler class).

◆ n_raw_quads() [1/9]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_raw_quads ( ) const

Total number of quads, used or unused.

Note
This function really exports internal information about the triangulation. It shouldn't be used in applications. The function is only part of the public interface of this class because it is used in some of the other classes that build very closely on it (in particular, the DoFHandler class).

◆ n_raw_quads() [2/9]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_raw_quads ( const unsigned int  level) const

Number of quads, used or unused, on the given level.

Note
This function really exports internal information about the triangulation. It shouldn't be used in applications. The function is only part of the public interface of this class because it is used in some of the other classes that build very closely on it (in particular, the DoFHandler class).

◆ n_raw_hexs() [1/5]

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_raw_hexs ( const unsigned int  level) const

Number of hexs, used or unused, on the given level.

Note
This function really exports internal information about the triangulation. It shouldn't be used in applications. The function is only part of the public interface of this class because it is used in some of the other classes that build very closely on it (in particular, the DoFHandler class).

◆ n_raw_cells()

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_raw_cells ( const unsigned int  level) const

Number of cells, used or unused, on the given level.

Note
This function really exports internal information about the triangulation. It shouldn't be used in applications. The function is only part of the public interface of this class because it is used in some of the other classes that build very closely on it (in particular, the DoFHandler class).

◆ n_raw_faces()

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::n_raw_faces ( ) const

Return the total number of faces, used or not. In 2d, the result equals n_raw_lines(), in 3d it equals n_raw_quads(), while in 1d it equals the number of vertices.

Note
This function really exports internal information about the triangulation. It shouldn't be used in applications. The function is only part of the public interface of this class because it is used in some of the other classes that build very closely on it (in particular, the DoFHandler class).

◆ memory_consumption()

template<int dim, int spacedim = dim>
virtual std::size_t Triangulation< dim, spacedim >::memory_consumption ( ) const
virtual

Determine an estimate for the memory consumption (in bytes) of this object.

This function is made virtual, since a triangulation object might be accessed through a pointer to this base class, even if the actual object is a derived class.

Reimplemented in parallel::fullydistributed::Triangulation< dim, spacedim >, parallel::distributed::Triangulation< dim, spacedim >, parallel::TriangulationBase< dim, spacedim >, parallel::TriangulationBase< dim, dim >, parallel::TriangulationBase< dim, spacedim >, and PersistentTriangulation< dim, spacedim >.

◆ save() [1/2]

template<int dim, int spacedim = dim>
template<class Archive >
void Triangulation< dim, spacedim >::save ( Archive &  ar,
const unsigned int  version 
) const

Write the data of this object to a stream for the purpose of serialization using the BOOST serialization library.

Note
This function does not save all member variables of the current triangulation. Rather, only certain kinds of information are stored. For more information see the general documentation of this class.

◆ load() [1/2]

template<int dim, int spacedim = dim>
template<class Archive >
void Triangulation< dim, spacedim >::load ( Archive &  ar,
const unsigned int  version 
)

Read the data of this object from a stream for the purpose of serialization using the BOOST serialization library. Throw away the previous content.

Note
This function does not reset all member variables of the current triangulation to the ones of the triangulation that was previously stored to an archive. Rather, only certain kinds of information are loaded. For more information see the general documentation of this class.
This function calls the Triangulation::clear() function and consequently triggers the "clear" signal. After loading all data from the archive, it then triggers the "create" signal. For more information on signals, see the general documentation of this class.

◆ save() [2/2]

template<int dim, int spacedim = dim>
virtual void Triangulation< dim, spacedim >::save ( const std::string &  filename) const
virtual

Save the triangulation into the given file. Internally, this function calls the save function which uses BOOST archives. This is a placeholder implementation that, in the near future, will also attach the data associated with the triangulation

Reimplemented in parallel::fullydistributed::Triangulation< dim, spacedim >, parallel::distributed::Triangulation< dim, spacedim >, and parallel::distributed::Triangulation< 1, spacedim >.

◆ load() [2/2]

template<int dim, int spacedim = dim>
virtual void Triangulation< dim, spacedim >::load ( const std::string &  filename)
virtual

◆ add_periodicity()

template<int dim, int spacedim = dim>
virtual void Triangulation< dim, spacedim >::add_periodicity ( const std::vector< GridTools::PeriodicFacePair< cell_iterator > > &  )
virtual

Declare the (coarse) face pairs given in the argument of this function as periodic. This way it is possible to obtain neighbors across periodic boundaries.

The vector can be filled by the function GridTools::collect_periodic_faces.

For more information on periodic boundary conditions see GridTools::collect_periodic_faces, DoFTools::make_periodicity_constraints and step-45.

Note
Before this function can be used the Triangulation has to be initialized and must not be refined.

◆ get_periodic_face_map()

template<int dim, int spacedim = dim>
const std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, std::bitset< 3 > > > & Triangulation< dim, spacedim >::get_periodic_face_map ( ) const

Return the periodic_face_map.

◆ get_reference_cells()

template<int dim, int spacedim = dim>
const std::vector< ReferenceCell > & Triangulation< dim, spacedim >::get_reference_cells ( ) const

Return vector filled with the used reference-cell types of this triangulation.

◆ all_reference_cells_are_hyper_cube()

template<int dim, int spacedim = dim>
bool Triangulation< dim, spacedim >::all_reference_cells_are_hyper_cube ( ) const

Indicate if the triangulation only consists of hypercube-like cells, i.e., lines, quadrilaterals, or hexahedra.

◆ all_reference_cells_are_simplex()

template<int dim, int spacedim = dim>
bool Triangulation< dim, spacedim >::all_reference_cells_are_simplex ( ) const

Indicate if the triangulation only consists of simplex-like cells, i.e., lines, triangles, or tetrahedra.

◆ is_mixed_mesh()

template<int dim, int spacedim = dim>
bool Triangulation< dim, spacedim >::is_mixed_mesh ( ) const

Indicate if the triangulation consists of different cell types (mix of simplices, hypercubes, ...) or different face types, as in the case of pyramids or wedges..

◆ serialize()

template<int dim, int spacedim = dim>
template<class Archive >
void Triangulation< dim, spacedim >::serialize ( Archive &  archive,
const unsigned int  version 
)

Write and read the data of this object from a stream for the purpose of serialization. using the BOOST serialization library.

◆ register_data_attach()

template<int dim, int spacedim = dim>
unsigned int Triangulation< dim, spacedim >::register_data_attach ( const std::function< std::vector< char >(const cell_iterator &, const ::CellStatus)> &  pack_callback,
const bool  returns_variable_size_data 
)

Register a function that can be used to attach data of fixed size to cells. This is useful for two purposes: (i) Upon refinement and coarsening of a triangulation (e.g. in parallel::distributed::Triangulation::execute_coarsening_and_refinement()), one needs to be able to store one or more data vectors per cell that characterizes the solution values on the cell so that this data can then be transferred to the new owning processor of the cell (or its parent/children) when the mesh is re-partitioned; (ii) when serializing a computation to a file, it is necessary to attach data to cells so that it can be saved (e.g. in parallel::distributed::Triangulation::save()) along with the cell's other information and, if necessary, later be reloaded from disk with a different subdivision of cells among the processors.

The way this function works is that it allows any number of interest parties to register their intent to attach data to cells. One example of classes that do this is parallel::distributed::SolutionTransfer where each parallel::distributed::SolutionTransfer object that works on the current Triangulation object then needs to register its intent. Each of these parties registers a callback function (the first argument here, pack_callback) that will be called whenever the triangulation's execute_coarsening_and_refinement() or save() functions are called.

The current function then returns an integer handle that corresponds to the number of data set that the callback provided here will attach. While this number could be given a precise meaning, this is not important: You will never actually have to do anything with this number except return it to the notify_ready_to_unpack() function. In other words, each interested party (i.e., the caller of the current function) needs to store their respective returned handle for later use when unpacking data in the callback provided to notify_ready_to_unpack().

Whenever pack_callback is then called by execute_coarsening_and_refinement() or load() on a given cell, it receives a number of arguments. In particular, the first argument passed to the callback indicates the cell for which it is supposed to attach data. This is always an active cell.

The second, CellStatus, argument provided to the callback function will tell you if the given cell will be coarsened, refined, or will persist as is. (This status may be different than the refinement or coarsening flags set on that cell, to accommodate things such as the "one hanging node per edge" rule.). These flags need to be read in context with the p4est quadrant they belong to, as their relations are gathered in local_cell_relations.

Specifically, the values for this argument mean the following:

  • CellStatus::cell_will_persist: The cell won't be refined/coarsened, but might be moved to a different processor. If this is the case, the callback will want to pack up the data on this cell into an array and store it at the provided address for later unpacking wherever this cell may land.
  • CellStatus::cell_will_be_refined: This cell will be refined into 4 or 8 cells (in 2d and 3d, respectively). However, because these children don't exist yet, you cannot access them at the time when the callback is called. Thus, in local_cell_relations, the corresponding p4est quadrants of the children cells are linked to the deal.II cell which is going to be refined. To be specific, only the very first child is marked with CellStatus::cell_will_be_refined, whereas the others will be marked with CellStatus::cell_invalid, which indicates that these cells will be ignored by default during the packing or unpacking process. This ensures that data is only transferred once onto or from the parent cell. If the callback is called with CellStatus::cell_will_be_refined, the callback will want to pack up the data on this cell into an array and store it at the provided address for later unpacking in a way so that it can then be transferred to the children of the cell that will then be available. In other words, if the data the callback will want to pack up corresponds to a finite element field, then the prolongation from parent to (new) children will have to happen during unpacking.
  • CellStatus::children_will_be_coarsened: The children of this cell will be coarsened into the given cell. These children still exist, so if this is the value given to the callback as second argument, the callback will want to transfer data from the children to the current parent cell and pack it up so that it can later be unpacked again on a cell that then no longer has any children (and may also be located on a different processor). In other words, if the data the callback will want to pack up corresponds to a finite element field, then it will need to do the restriction from children to parent at this point.
  • CellStatus::cell_invalid: See CellStatus::cell_will_be_refined.
Note
If this function is used for serialization of data using save() and load(), then the cell status argument with which the callback is called will always be CellStatus::cell_will_persist.

The callback function is expected to return a memory chunk of the format std::vector<char>, representing the packed data on a certain cell.

The second parameter returns_variable_size_data indicates whether the returned size of the memory region from the callback function varies by cell (=true) or stays constant on each one throughout the whole domain (=false).

Note
The purpose of this function is to register intent to attach data for a single, subsequent call to execute_coarsening_and_refinement() and notify_ready_to_unpack(), save(), load(). Consequently, notify_ready_to_unpack(), save(), and load() all forget the registered callbacks once these callbacks have been called, and you will have to re-register them with a triangulation if you want them to be active for another call to these functions.

◆ notify_ready_to_unpack()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::notify_ready_to_unpack ( const unsigned int  handle,
const std::function< void(const cell_iterator &, const ::CellStatus, const boost::iterator_range< std::vector< char >::const_iterator > &)> &  unpack_callback 
)

This function is the opposite of register_data_attach(). It is called after the execute_coarsening_and_refinement() or save()/load() functions are done when classes and functions that have previously attached data to a triangulation for either transfer to other processors, across mesh refinement, or serialization of data to a file are ready to receive that data back. The important part about this process is that the triangulation cannot do this right away from the end of execute_coarsening_and_refinement() or load() via a previously attached callback function (as the register_data_attach() function does) because the classes that eventually want the data back may need to do some setup between the point in time where the mesh has been recreated and when the data can actually be received. An example is the parallel::distributed::SolutionTransfer class that can really only receive the data once not only the mesh is completely available again on the current processor, but only after a DoFHandler has been reinitialized and distributed degrees of freedom. In other words, there is typically a significant amount of set up that needs to happen in user space before the classes that can receive data attached to cell are ready to actually do so. When they are, they use the current function to tell the triangulation object that now is the time when they are ready by calling the current function.

The supplied callback function is then called for each newly locally owned cell. The first argument to the callback is an iterator that designates the cell; the second argument indicates the status of the cell in question; and the third argument localizes a memory area by two iterators that contains the data that was previously saved from the callback provided to register_data_attach().

The CellStatus will indicate if the cell was refined, coarsened, or persisted unchanged. The cell_iterator argument to the callback will then either be an active, locally owned cell (if the cell was not refined), or the immediate parent if it was refined during execute_coarsening_and_refinement(). Therefore, contrary to during register_data_attach(), you can now access the children if the status is CellStatus::cell_will_be_refined but no longer for callbacks with status CellStatus::children_will_be_coarsened.

The first argument to this function, handle, corresponds to the return value of register_data_attach(). (The precise meaning of what the numeric value of this handle is supposed to represent is neither important, nor should you try to use it for anything other than transmit information between a call to register_data_attach() to the corresponding call to notify_ready_to_unpack().)

◆ save_attached_data()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::save_attached_data ( const unsigned int  global_first_cell,
const unsigned int  global_num_cells,
const std::string &  filename 
) const
protected

Save additional cell-attached data into the given file. The first arguments are used to determine the offsets where to write buffers to.

Called by save.

◆ load_attached_data()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::load_attached_data ( const unsigned int  global_first_cell,
const unsigned int  global_num_cells,
const unsigned int  local_num_cells,
const std::string &  filename,
const unsigned int  n_attached_deserialize_fixed,
const unsigned int  n_attached_deserialize_variable 
)
protected

Load additional cell-attached data from the given file, if any was saved. The first arguments are used to determine the offsets where to read buffers from.

Called by load.

◆ update_cell_relations()

template<int dim, int spacedim = dim>
virtual void Triangulation< dim, spacedim >::update_cell_relations ( )
inlineprotectedvirtual

A function to record the CellStatus of currently active cells that are locally owned. This information is mandatory to transfer data between meshes during adaptation or serialization, e.g., using parallel::distributed::SolutionTransfer.

Relations will be stored in the private member local_cell_relations. For an extensive description of CellStatus, see the documentation for the member function register_data_attach().

Reimplemented in parallel::fullydistributed::Triangulation< dim, spacedim >, parallel::distributed::Triangulation< dim, spacedim >, and parallel::distributed::Triangulation< 1, spacedim >.

Definition at line 3916 of file tria.h.

◆ write_bool_vector()

template<int dim, int spacedim = dim>
static void Triangulation< dim, spacedim >::write_bool_vector ( const unsigned int  magic_number1,
const std::vector< bool > &  v,
const unsigned int  magic_number2,
std::ostream &  out 
)
staticprotected

Write a bool vector to the given stream, writing a pre- and a postfix magic number. The vector is written in an almost binary format, i.e. the bool flags are packed but the data is written as ASCII text.

The flags are stored in a binary format: for each true, a 1 bit is stored, a 0 bit otherwise. The bits are stored as unsigned char, thus avoiding endianness. They are written to out in plain text, thus amounting to 3.6 bits in the output per bits in the input on the average. Other information (magic numbers and number of elements of the input vector) is stored as plain text as well. The format should therefore be interplatform compatible.

◆ read_bool_vector()

template<int dim, int spacedim = dim>
static void Triangulation< dim, spacedim >::read_bool_vector ( const unsigned int  magic_number1,
std::vector< bool > &  v,
const unsigned int  magic_number2,
std::istream &  in 
)
staticprotected

Re-read a vector of bools previously written by write_bool_vector and compare with the magic numbers.

◆ update_periodic_face_map()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::update_periodic_face_map ( )
protected

Recreate information about periodic neighbors from periodic_face_pairs_level_0.

◆ update_reference_cells()

template<int dim, int spacedim = dim>
virtual void Triangulation< dim, spacedim >::update_reference_cells ( )
protectedvirtual

◆ begin_raw()

template<int dim, int spacedim = dim>
raw_cell_iterator Triangulation< dim, spacedim >::begin_raw ( const unsigned int  level = 0) const
private

Iterator to the first cell, used or not, on level level. If a level has no cells, a past-the-end iterator is returned.

◆ end_raw()

template<int dim, int spacedim = dim>
raw_cell_iterator Triangulation< dim, spacedim >::end_raw ( const unsigned int  level) const
private

Return a raw iterator which is the first iterator not on level. If level is the last level, then this returns end().

◆ begin_raw_line()

template<int dim, int spacedim = dim>
raw_line_iterator Triangulation< dim, spacedim >::begin_raw_line ( const unsigned int  level = 0) const
private

Iterator to the first line, used or not, on level level. If a level has no lines, a past-the-end iterator is returned. If lines are no cells, i.e. for dim>1 no level argument must be given. The same applies for all the other functions above, of course.

◆ begin_line()

template<int dim, int spacedim = dim>
line_iterator Triangulation< dim, spacedim >::begin_line ( const unsigned int  level = 0) const
private

Iterator to the first used line on level level.

Note
The given level argument needs to correspond to a level of the triangulation, i.e., should be less than the value returned by n_levels(). On the other hand, for parallel computations using a parallel::distributed::Triangulation object, it is often convenient to write loops over the cells of all levels of the global mesh, even if the local portion of the triangulation does not actually have cells at one of the higher levels. In those cases, the level argument is accepted if it is less than what the n_global_levels() function returns. If the given level is between the values returned by n_levels() and n_global_levels(), then no cells exist in the local portion of the triangulation at this level, and the function simply returns what end() would return.

◆ begin_active_line()

template<int dim, int spacedim = dim>
active_line_iterator Triangulation< dim, spacedim >::begin_active_line ( const unsigned int  level = 0) const
private

Iterator to the first active line on level level.

Note
The given level argument needs to correspond to a level of the triangulation, i.e., should be less than the value returned by n_levels(). On the other hand, for parallel computations using a parallel::distributed::Triangulation object, it is often convenient to write loops over the cells of all levels of the global mesh, even if the local portion of the triangulation does not actually have cells at one of the higher levels. In those cases, the level argument is accepted if it is less than what the n_global_levels() function returns. If the given level is between the values returned by n_levels() and n_global_levels(), then no cells exist in the local portion of the triangulation at this level, and the function simply returns what end() would return.

◆ end_line()

template<int dim, int spacedim = dim>
line_iterator Triangulation< dim, spacedim >::end_line ( ) const
private

Iterator past the end; this iterator serves for comparisons of iterators with past-the-end or before-the-beginning states.

◆ begin_raw_quad()

template<int dim, int spacedim = dim>
raw_quad_iterator Triangulation< dim, spacedim >::begin_raw_quad ( const unsigned int  level = 0) const
private

Iterator to the first quad, used or not, on the given level. If a level has no quads, a past-the-end iterator is returned. If quads are no cells, i.e. for \(dim>2\) no level argument must be given.

Note
The given level argument needs to correspond to a level of the triangulation, i.e., should be less than the value returned by n_levels(). On the other hand, for parallel computations using a parallel::distributed::Triangulation object, it is often convenient to write loops over the cells of all levels of the global mesh, even if the local portion of the triangulation does not actually have cells at one of the higher levels. In those cases, the level argument is accepted if it is less than what the n_global_levels() function returns. If the given level is between the values returned by n_levels() and n_global_levels(), then no cells exist in the local portion of the triangulation at this level, and the function simply returns what end() would return.

◆ begin_quad()

template<int dim, int spacedim = dim>
quad_iterator Triangulation< dim, spacedim >::begin_quad ( const unsigned int  level = 0) const
private

Iterator to the first used quad on level level.

Note
The given level argument needs to correspond to a level of the triangulation, i.e., should be less than the value returned by n_levels(). On the other hand, for parallel computations using a parallel::distributed::Triangulation object, it is often convenient to write loops over the cells of all levels of the global mesh, even if the local portion of the triangulation does not actually have cells at one of the higher levels. In those cases, the level argument is accepted if it is less than what the n_global_levels() function returns. If the given level is between the values returned by n_levels() and n_global_levels(), then no cells exist in the local portion of the triangulation at this level, and the function simply returns what end() would return.

◆ begin_active_quad()

template<int dim, int spacedim = dim>
active_quad_iterator Triangulation< dim, spacedim >::begin_active_quad ( const unsigned int  level = 0) const
private

Iterator to the first active quad on level level.

Note
The given level argument needs to correspond to a level of the triangulation, i.e., should be less than the value returned by n_levels(). On the other hand, for parallel computations using a parallel::distributed::Triangulation object, it is often convenient to write loops over the cells of all levels of the global mesh, even if the local portion of the triangulation does not actually have cells at one of the higher levels. In those cases, the level argument is accepted if it is less than what the n_global_levels() function returns. If the given level is between the values returned by n_levels() and n_global_levels(), then no cells exist in the local portion of the triangulation at this level, and the function simply returns what end() would return.

◆ end_quad()

template<int dim, int spacedim = dim>
quad_iterator Triangulation< dim, spacedim >::end_quad ( ) const
private

Iterator past the end; this iterator serves for comparisons of iterators with past-the-end or before-the-beginning states.

◆ begin_raw_hex()

template<int dim, int spacedim = dim>
raw_hex_iterator Triangulation< dim, spacedim >::begin_raw_hex ( const unsigned int  level = 0) const
private

Iterator to the first hex, used or not, on level level. If a level has no hexes, a past-the-end iterator is returned.

Note
The given level argument needs to correspond to a level of the triangulation, i.e., should be less than the value returned by n_levels(). On the other hand, for parallel computations using a parallel::distributed::Triangulation object, it is often convenient to write loops over the cells of all levels of the global mesh, even if the local portion of the triangulation does not actually have cells at one of the higher levels. In those cases, the level argument is accepted if it is less than what the n_global_levels() function returns. If the given level is between the values returned by n_levels() and n_global_levels(), then no cells exist in the local portion of the triangulation at this level, and the function simply returns what end() would return.

◆ begin_hex()

template<int dim, int spacedim = dim>
hex_iterator Triangulation< dim, spacedim >::begin_hex ( const unsigned int  level = 0) const
private

Iterator to the first used hex on level level.

Note
The given level argument needs to correspond to a level of the triangulation, i.e., should be less than the value returned by n_levels(). On the other hand, for parallel computations using a parallel::distributed::Triangulation object, it is often convenient to write loops over the cells of all levels of the global mesh, even if the local portion of the triangulation does not actually have cells at one of the higher levels. In those cases, the level argument is accepted if it is less than what the n_global_levels() function returns. If the given level is between the values returned by n_levels() and n_global_levels(), then no cells exist in the local portion of the triangulation at this level, and the function simply returns what end() would return.

◆ begin_active_hex()

template<int dim, int spacedim = dim>
active_hex_iterator Triangulation< dim, spacedim >::begin_active_hex ( const unsigned int  level = 0) const
private

Iterator to the first active hex on level level.

Note
The given level argument needs to correspond to a level of the triangulation, i.e., should be less than the value returned by n_levels(). On the other hand, for parallel computations using a parallel::distributed::Triangulation object, it is often convenient to write loops over the cells of all levels of the global mesh, even if the local portion of the triangulation does not actually have cells at one of the higher levels. In those cases, the level argument is accepted if it is less than what the n_global_levels() function returns. If the given level is between the values returned by n_levels() and n_global_levels(), then no cells exist in the local portion of the triangulation at this level, and the function simply returns what end() would return.

◆ end_hex()

template<int dim, int spacedim = dim>
hex_iterator Triangulation< dim, spacedim >::end_hex ( ) const
private

Iterator past the end; this iterator serves for comparisons of iterators with past-the-end or before-the-beginning states.

◆ clear_despite_subscriptions()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::clear_despite_subscriptions ( )
private

The (public) function clear() will only work when the triangulation is not subscribed to by other users. The clear_despite_subscriptions() function now allows the triangulation being cleared even when there are subscriptions.

Make sure, you know what you do, when calling this function, as its use is reasonable in very rare cases, only. For example, when the subscriptions were for the initially empty Triangulation and the Triangulation object wants to release its memory before throwing an assertion due to input errors (e.g. in the create_triangulation() function).

◆ reset_policy()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::reset_policy ( )
private

Reset triangulation policy.

◆ reset_active_cell_indices()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::reset_active_cell_indices ( )
private

For all cells, set the active cell indices so that active cells know the how many-th active cell they are, and all other cells have an invalid value. This function is called after mesh creation, refinement, and serialization.

◆ reset_global_cell_indices()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::reset_global_cell_indices ( )
private

Reset global cell ids and global level cell ids.

◆ reset_cell_vertex_indices_cache()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::reset_cell_vertex_indices_cache ( )
private

Reset cache for the cells' vertex indices.

◆ execute_refinement()

template<int dim, int spacedim = dim>
DistortedCellList Triangulation< dim, spacedim >::execute_refinement ( )
private

Refine all cells on all levels which were previously flagged for refinement.

Note, that this function uses the line->user_flags for dim=2,3 and the quad->user_flags for dim=3.

The function returns a list of cells that have produced children that satisfy the criteria of distorted cells if the check_for_distorted_cells flag was specified upon creation of this object, at

◆ execute_coarsening()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::execute_coarsening ( )
private

Coarsen all cells which were flagged for coarsening, or rather: delete all children of those cells of which all child cells are flagged for coarsening and several other constraints hold (see the general doc of this class).

◆ fix_coarsen_flags()

template<int dim, int spacedim = dim>
void Triangulation< dim, spacedim >::fix_coarsen_flags ( )
private

Make sure that either all or none of the children of a cell are tagged for coarsening.

◆ coarse_cell_id_to_coarse_cell_index()

template<int dim, int spacedim = dim>
virtual unsigned int Triangulation< dim, spacedim >::coarse_cell_id_to_coarse_cell_index ( const types::coarse_cell_id  coarse_cell_id) const
privatevirtual

Translate the unique id of a coarse cell to its index. See the glossary entry on coarse cell IDs for more information.

Note
For serial and shared triangulation both id and index are the same. For distributed triangulations setting both might differ, since the id might correspond to a global id and the index to a local id. If a cell does not exist locally, the returned value is numbers::invalid_unsigned_int.
Parameters
coarse_cell_idUnique id of the coarse cell.
Returns
Index of the coarse cell within the current triangulation.

Reimplemented in parallel::fullydistributed::Triangulation< dim, spacedim >, parallel::distributed::Triangulation< dim, spacedim >, and parallel::distributed::Triangulation< 1, spacedim >.

◆ coarse_cell_index_to_coarse_cell_id()

template<int dim, int spacedim = dim>
virtual types::coarse_cell_id Triangulation< dim, spacedim >::coarse_cell_index_to_coarse_cell_id ( const unsigned int  coarse_cell_index) const
privatevirtual

Translate the index of coarse cell to its unique id. See the glossary entry on coarse cell IDs for more information.

Note
See the note of the method coarse_cell_id_to_coarse_cell_index().
Parameters
coarse_cell_indexIndex of the coarse cell.
Returns
Id of the coarse cell.

Reimplemented in parallel::fullydistributed::Triangulation< dim, spacedim >, parallel::distributed::Triangulation< dim, spacedim >, and parallel::distributed::Triangulation< 1, spacedim >.

◆ n_quads() [3/8]

unsigned int Triangulation< 1, 1 >::n_quads ( ) const

Definition at line 15193 of file tria.cc.

◆ n_quads() [4/8]

unsigned int Triangulation< 1, 1 >::n_quads ( const unsigned int  ) const

Definition at line 15201 of file tria.cc.

◆ n_raw_quads() [3/9]

unsigned int Triangulation< 1, 1 >::n_raw_quads ( const unsigned int  ) const

Definition at line 15209 of file tria.cc.

◆ n_raw_hexs() [2/5]

unsigned int Triangulation< 1, 1 >::n_raw_hexs ( const unsigned int  ) const

Definition at line 15217 of file tria.cc.

◆ n_active_quads() [3/8]

unsigned int Triangulation< 1, 1 >::n_active_quads ( const unsigned int  ) const

Definition at line 15225 of file tria.cc.

◆ n_active_quads() [4/8]

unsigned int Triangulation< 1, 1 >::n_active_quads ( ) const

Definition at line 15233 of file tria.cc.

◆ n_quads() [5/8]

unsigned int Triangulation< 1, 2 >::n_quads ( ) const

Definition at line 15242 of file tria.cc.

◆ n_quads() [6/8]

unsigned int Triangulation< 1, 2 >::n_quads ( const unsigned int  ) const

Definition at line 15250 of file tria.cc.

◆ n_raw_quads() [4/9]

unsigned int Triangulation< 1, 2 >::n_raw_quads ( const unsigned int  ) const

Definition at line 15258 of file tria.cc.

◆ n_raw_hexs() [3/5]

unsigned int Triangulation< 1, 2 >::n_raw_hexs ( const unsigned int  ) const

Definition at line 15266 of file tria.cc.

◆ n_active_quads() [5/8]

unsigned int Triangulation< 1, 2 >::n_active_quads ( const unsigned int  ) const

Definition at line 15274 of file tria.cc.

◆ n_active_quads() [6/8]

unsigned int Triangulation< 1, 2 >::n_active_quads ( ) const

Definition at line 15282 of file tria.cc.

◆ n_quads() [7/8]

unsigned int Triangulation< 1, 3 >::n_quads ( ) const

Definition at line 15290 of file tria.cc.

◆ n_quads() [8/8]

unsigned int Triangulation< 1, 3 >::n_quads ( const unsigned int  ) const

Definition at line 15298 of file tria.cc.

◆ n_raw_quads() [5/9]

unsigned int Triangulation< 1, 3 >::n_raw_quads ( const unsigned int  ) const

Definition at line 15306 of file tria.cc.

◆ n_raw_hexs() [4/5]

unsigned int Triangulation< 1, 3 >::n_raw_hexs ( const unsigned int  ) const

Definition at line 15314 of file tria.cc.

◆ n_active_quads() [7/8]

unsigned int Triangulation< 1, 3 >::n_active_quads ( const unsigned int  ) const

Definition at line 15322 of file tria.cc.

◆ n_active_quads() [8/8]

unsigned int Triangulation< 1, 3 >::n_active_quads ( ) const

Definition at line 15330 of file tria.cc.

◆ n_raw_quads() [6/9]

unsigned int Triangulation< 2, 2 >::n_raw_quads ( const unsigned int  level) const

Definition at line 15359 of file tria.cc.

◆ n_raw_quads() [7/9]

unsigned int Triangulation< 2, 3 >::n_raw_quads ( const unsigned int  level) const

Definition at line 15369 of file tria.cc.

◆ n_raw_quads() [8/9]

unsigned int Triangulation< 3, 3 >::n_raw_quads ( const unsigned int  ) const

Definition at line 15378 of file tria.cc.

◆ n_raw_quads() [9/9]

unsigned int Triangulation< 3, 3 >::n_raw_quads ( ) const

Definition at line 15398 of file tria.cc.

◆ n_hexs() [3/4]

unsigned int Triangulation< 3, 3 >::n_hexs ( ) const

Definition at line 15472 of file tria.cc.

◆ n_hexs() [4/4]

unsigned int Triangulation< 3, 3 >::n_hexs ( const unsigned int  level) const

Definition at line 15481 of file tria.cc.

◆ n_raw_hexs() [5/5]

unsigned int Triangulation< 3, 3 >::n_raw_hexs ( const unsigned int  level) const

Definition at line 15492 of file tria.cc.

◆ n_active_hexs() [3/4]

unsigned int Triangulation< 3, 3 >::n_active_hexs ( ) const

Definition at line 15501 of file tria.cc.

◆ n_active_hexs() [4/4]

unsigned int Triangulation< 3, 3 >::n_active_hexs ( const unsigned int  level) const

Definition at line 15510 of file tria.cc.

◆ max_adjacent_cells() [2/4]

unsigned int Triangulation< 1, 1 >::max_adjacent_cells ( ) const

Definition at line 15539 of file tria.cc.

◆ max_adjacent_cells() [3/4]

unsigned int Triangulation< 1, 2 >::max_adjacent_cells ( ) const

Definition at line 15548 of file tria.cc.

◆ max_adjacent_cells() [4/4]

unsigned int Triangulation< 1, 3 >::max_adjacent_cells ( ) const

Definition at line 15556 of file tria.cc.

◆ prepare_coarsening_and_refinement() [2/4]

bool Triangulation< 1, 1 >::prepare_coarsening_and_refinement ( )

Definition at line 16458 of file tria.cc.

◆ prepare_coarsening_and_refinement() [3/4]

bool Triangulation< 1, 2 >::prepare_coarsening_and_refinement ( )

Definition at line 16476 of file tria.cc.

◆ prepare_coarsening_and_refinement() [4/4]

bool Triangulation< 1, 3 >::prepare_coarsening_and_refinement ( )

Definition at line 16494 of file tria.cc.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Subscribes a user of the object by storing the pointer validity. The subscriber may be identified by text supplied as identifier.

Definition at line 135 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 155 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
inlineinherited

Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.

Definition at line 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 203 of file subscriptor.cc.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
privatenoexceptinherited

Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 52 of file subscriptor.cc.

Friends And Related Symbol Documentation

◆ TriaAccessorBase

template<int dim, int spacedim = dim>
template<int , int , int >
friend class TriaAccessorBase
friend

Definition at line 4558 of file tria.h.

◆ TriaAccessor

template<int dim, int spacedim = dim>
template<int , int , int >
friend class TriaAccessor
friend

Definition at line 4560 of file tria.h.

◆ TriaAccessor< 0, 1, spacedim >

template<int dim, int spacedim = dim>
friend class TriaAccessor< 0, 1, spacedim >
friend

Definition at line 4560 of file tria.h.

◆ CellAccessor< dim, spacedim >

template<int dim, int spacedim = dim>
friend class CellAccessor< dim, spacedim >
friend

Definition at line 4560 of file tria.h.

◆ ::internal::TriaAccessorImplementation::Implementation

template<int dim, int spacedim = dim>
friend struct ::internal::TriaAccessorImplementation::Implementation
friend

Definition at line 4565 of file tria.h.

◆ ::internal::TriangulationImplementation::Implementation

template<int dim, int spacedim = dim>
friend struct ::internal::TriangulationImplementation::Implementation
friend

Definition at line 4567 of file tria.h.

◆ ::internal::TriangulationImplementation::ImplementationMixedMesh

template<int dim, int spacedim = dim>
friend struct ::internal::TriangulationImplementation::ImplementationMixedMesh
friend

Definition at line 4568 of file tria.h.

◆ ::internal::TriangulationImplementation::TriaObjects

template<int dim, int spacedim = dim>
friend class ::internal::TriangulationImplementation::TriaObjects
friend

Definition at line 4571 of file tria.h.

Member Data Documentation

◆ dimension

template<int dim, int spacedim = dim>
constexpr unsigned int Triangulation< dim, spacedim >::dimension = dim
staticconstexpr

Make the dimension available in function templates.

Definition at line 1757 of file tria.h.

◆ space_dimension

template<int dim, int spacedim = dim>
constexpr unsigned int Triangulation< dim, spacedim >::space_dimension = spacedim
staticconstexpr

Make the space-dimension available in function templates.

Definition at line 1762 of file tria.h.

◆ CELL_PERSIST

template<int dim, int spacedim = dim>
constexpr auto Triangulation< dim, spacedim >::CELL_PERSIST
staticconstexpr
Initial value:
Deprecated:
This is an alias for backward compatibility. Use CellStatus directly.

Definition at line 2250 of file tria.h.

◆ CELL_REFINE

template<int dim, int spacedim = dim>
constexpr auto Triangulation< dim, spacedim >::CELL_REFINE
staticconstexpr
Initial value:
Deprecated:
This is an alias for backward compatibility. Use CellStatus directly.

Definition at line 2257 of file tria.h.

◆ CELL_COARSEN

template<int dim, int spacedim = dim>
constexpr auto Triangulation< dim, spacedim >::CELL_COARSEN
staticconstexpr
Initial value:
Deprecated:
This is an alias for backward compatibility. Use CellStatus directly.

Definition at line 2264 of file tria.h.

◆ CELL_INVALID

template<int dim, int spacedim = dim>
constexpr auto Triangulation< dim, spacedim >::CELL_INVALID
staticconstexpr
Initial value:
Deprecated:
This is an alias for backward compatibility. Use CellStatus directly.

Definition at line 2271 of file tria.h.

◆ signals

template<int dim, int spacedim = dim>
Signals Triangulation< dim, spacedim >::signals
mutable

Signals for the various actions that a triangulation can do to itself.

Definition at line 2529 of file tria.h.

◆ cell_attached_data

template<int dim, int spacedim = dim>
internal::CellAttachedData<dim, spacedim> Triangulation< dim, spacedim >::cell_attached_data

Definition at line 3876 of file tria.h.

◆ local_cell_relations

template<int dim, int spacedim = dim>
std::vector<typename internal::CellAttachedDataSerializer<dim, spacedim>:: cell_relation_t> Triangulation< dim, spacedim >::local_cell_relations
protected

Vector of pairs, each containing a deal.II cell iterator and its respective CellStatus. To update its contents, use the update_cell_relations() member function.

Definition at line 3926 of file tria.h.

◆ data_serializer

template<int dim, int spacedim = dim>
internal::CellAttachedDataSerializer<dim, spacedim> Triangulation< dim, spacedim >::data_serializer
protected

Definition at line 3928 of file tria.h.

◆ smooth_grid

template<int dim, int spacedim = dim>
MeshSmoothing Triangulation< dim, spacedim >::smooth_grid
protected

Do some smoothing in the process of refining the triangulation. See the general doc of this class for more information about this.

Definition at line 4018 of file tria.h.

◆ reference_cells

template<int dim, int spacedim = dim>
std::vector<ReferenceCell> Triangulation< dim, spacedim >::reference_cells
protected

Vector caching all reference-cell types of the given triangulation (also in the distributed case).

Definition at line 4024 of file tria.h.

◆ policy

template<int dim, int spacedim = dim>
std::unique_ptr< ::internal::TriangulationImplementation::Policy<dim, spacedim> > Triangulation< dim, spacedim >::policy
private

Policy with the Triangulation-specific tasks related to creation, refinement, and coarsening.

Definition at line 4076 of file tria.h.

◆ periodic_face_pairs_level_0

template<int dim, int spacedim = dim>
std::vector<GridTools::PeriodicFacePair<cell_iterator> > Triangulation< dim, spacedim >::periodic_face_pairs_level_0
private

If add_periodicity() is called, this variable stores the given periodic face pairs on level 0 for later access during the identification of ghost cells for the multigrid hierarchy and for setting up the periodic_face_map.

Definition at line 4085 of file tria.h.

◆ periodic_face_map

template<int dim, int spacedim = dim>
std::map<std::pair<cell_iterator, unsigned int>, std::pair<std::pair<cell_iterator, unsigned int>, std::bitset<3> > > Triangulation< dim, spacedim >::periodic_face_map
private

If add_periodicity() is called, this variable stores the active periodic face pairs.

Definition at line 4093 of file tria.h.

◆ levels

template<int dim, int spacedim = dim>
std::vector< std::unique_ptr<::internal::TriangulationImplementation::TriaLevel> > Triangulation< dim, spacedim >::levels
private

Array of pointers pointing to the objects storing the cell data on the different levels.

Definition at line 4465 of file tria.h.

◆ faces

template<int dim, int spacedim = dim>
std::unique_ptr<::internal::TriangulationImplementation::TriaFaces> Triangulation< dim, spacedim >::faces
private

Pointer to the faces of the triangulation. In 1d this contains nothing, in 2d it contains data concerning lines and in 3d quads and lines. All of these have no level and are therefore treated separately.

Definition at line 4473 of file tria.h.

◆ vertices

template<int dim, int spacedim = dim>
std::vector<Point<spacedim> > Triangulation< dim, spacedim >::vertices
private

Array of the vertices of this triangulation.

Definition at line 4479 of file tria.h.

◆ vertices_used

template<int dim, int spacedim = dim>
std::vector<bool> Triangulation< dim, spacedim >::vertices_used
private

Array storing a bit-pattern which vertices are used.

Definition at line 4484 of file tria.h.

◆ manifolds

template<int dim, int spacedim = dim>
std::map<types::manifold_id, std::unique_ptr<const Manifold<dim, spacedim> > > Triangulation< dim, spacedim >::manifolds
private

Collection of manifold objects. We store only objects, which are not of type FlatManifold.

Definition at line 4491 of file tria.h.

◆ anisotropic_refinement

template<int dim, int spacedim = dim>
bool Triangulation< dim, spacedim >::anisotropic_refinement
private

Flag indicating whether anisotropic refinement took place.

Definition at line 4496 of file tria.h.

◆ check_for_distorted_cells

template<int dim, int spacedim = dim>
const bool Triangulation< dim, spacedim >::check_for_distorted_cells
private

A flag that determines whether we are to check for distorted cells upon creation and refinement of a mesh.

Definition at line 4503 of file tria.h.

◆ number_cache

template<int dim, int spacedim = dim>
::internal::TriangulationImplementation::NumberCache<dim> Triangulation< dim, spacedim >::number_cache
private

Cache to hold the numbers of lines, quads, hexes, etc. These numbers are set at the end of the refinement and coarsening functions and enable faster access later on. In the old days, whenever one wanted to access one of these numbers, one had to perform a loop over all lines, e.g., and count the elements until we hit the end iterator. This is time consuming and since access to the number of lines etc is a rather frequent operation, this was not an optimal solution.

Definition at line 4514 of file tria.h.

◆ vertex_to_boundary_id_map_1d

template<int dim, int spacedim = dim>
std::unique_ptr<std::map<unsigned int, types::boundary_id> > Triangulation< dim, spacedim >::vertex_to_boundary_id_map_1d
private

A map that relates the number of a boundary vertex to the boundary indicator. This field is only used in 1d. We have this field because we store boundary indicator information with faces in 2d and higher where we have space in the structures that store data for faces, but in 1d there is no such space for faces.

The field is declared as a pointer for a rather mundane reason: all other fields of this class that can be modified by the TriaAccessor hierarchy are pointers, and so these accessor classes store a const pointer to the triangulation. We could no longer do so for TriaAccessor<0,1,spacedim> if this field (that can be modified by TriaAccessor::set_boundary_id) were not a pointer.

Definition at line 4531 of file tria.h.

◆ vertex_to_manifold_id_map_1d

template<int dim, int spacedim = dim>
std::unique_ptr<std::map<unsigned int, types::manifold_id> > Triangulation< dim, spacedim >::vertex_to_manifold_id_map_1d
private

A map that relates the number of a boundary vertex to the manifold indicator. This field is only used in 1d. We have this field because we store manifold indicator information with faces in 2d and higher where we have space in the structures that store data for faces, but in 1d there is no such space for faces.

Note
Manifold objects are pretty useless for points since they are neither refined nor are their interiors mapped. We nevertheless allow storing manifold ids for points to be consistent in dimension-independent programs.

The field is declared as a pointer for a rather mundane reason: all other fields of this class that can be modified by the TriaAccessor hierarchy are pointers, and so these accessor classes store a const pointer to the triangulation. We could no longer do so for TriaAccessor<0,1,spacedim> if this field (that can be modified by TriaAccessor::set_manifold_id) were not a pointer.

Definition at line 4554 of file tria.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
mutableprivateinherited

Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).

The creator (and owner) of an object is counted in the map below if HE manages to supply identification.

We use the mutable keyword in order to allow subscription to constant objects also.

This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic class template.

Definition at line 218 of file subscriptor.h.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 224 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.

Definition at line 240 of file subscriptor.h.

◆ object_info

const std::type_info* Subscriptor::object_info
mutableprivateinherited

Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.

Definition at line 248 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
staticprivateinherited

A mutex used to ensure data consistency when accessing the mutable members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers().

Definition at line 271 of file subscriptor.h.


The documentation for this class was generated from the following files: