Reference documentation for deal.II version Git c92c73f660 2021-06-12 09:30:03 +0200
\(\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\}}\)
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Types | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
parallel::fullydistributed::Triangulation< dim, spacedim > Class Template Reference

#include <deal.II/distributed/fully_distributed_tria.h>

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

Public Types

using cell_iterator = typename ::Triangulation< dim, spacedim >::cell_iterator
 
using active_cell_iterator = typename ::Triangulation< dim, spacedim >::active_cell_iterator
 
using CellStatus = typename ::Triangulation< dim, spacedim >::CellStatus
 
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 level_cell_iterator = cell_iterator
 
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
 

Public Member Functions

 Triangulation (const MPI_Comm &mpi_communicator)
 
virtual ~Triangulation ()=default
 
void create_triangulation (const TriangulationDescription::Description< dim, spacedim > &construction_data) override
 
virtual void create_triangulation (const std::vector< Point< spacedim >> &vertices, const std::vector<::CellData< dim >> &cells, const SubCellData &subcelldata) override
 
void copy_triangulation (const ::Triangulation< dim, spacedim > &other_tria) override
 
void set_partitioner (const std::function< void(::Triangulation< dim, spacedim > &, const unsigned int)> &partitioner, const TriangulationDescription::Settings &settings)
 
virtual void execute_coarsening_and_refinement () override
 
virtual bool prepare_coarsening_and_refinement () override
 
virtual bool has_hanging_nodes () const override
 
virtual std::size_t memory_consumption () const override
 
virtual bool is_multilevel_hierarchy_constructed () const override
 
virtual void save (const std::string &filename) const override
 
virtual void load (const std::string &filename, const bool autopartition=false) override
 
virtual void clear () override
 
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)
 
virtual MPI_Comm get_communicator () const override
 
virtual void copy_triangulation (const ::Triangulation< dim, spacedim > &old_tria) override
 
virtual void copy_triangulation (const Triangulation< dim, spacedim > &other_tria)
 
unsigned int n_locally_owned_active_cells () const
 
virtual types::global_cell_index n_global_active_cells () const override
 
virtual unsigned int n_global_levels () const override
 
types::subdomain_id locally_owned_subdomain () const override
 
const std::set< types::subdomain_id > & ghost_owners () const
 
const std::set< types::subdomain_id > & level_ghost_owners () const
 
const std::weak_ptr< const Utilities::MPI::Partitionerglobal_active_cell_index_partitioner () const
 
const std::weak_ptr< const Utilities::MPI::Partitionerglobal_level_cell_index_partitioner (const unsigned int level) const
 
virtual std::vector< types::boundary_idget_boundary_ids () const override
 
virtual std::vector< types::manifold_idget_manifold_ids () const override
 
void communicate_locally_moved_vertices (const std::vector< bool > &vertex_locally_moved)
 
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 void create_triangulation (const std::vector< Point< spacedim >> &vertices, const std::vector< CellData< dim >> &cells, const SubCellData &subcelldata)
 
virtual void create_triangulation_compatibility (const std::vector< Point< spacedim >> &vertices, const std::vector< CellData< dim >> &cells, const SubCellData &subcelldata)
 
void flip_all_direction_flags ()
 
template<>
unsigned int n_quads () const
 
template<>
unsigned int n_quads (const unsigned int) const
 
template<>
unsigned int n_quads () const
 
template<>
unsigned int n_quads (const unsigned int) const
 
template<>
unsigned int n_quads () const
 
template<>
unsigned int n_quads (const unsigned int) const
 
template<>
unsigned int n_active_quads (const unsigned int) const
 
template<>
unsigned int n_active_quads () const
 
template<>
unsigned int n_active_quads (const unsigned int) const
 
template<>
unsigned int n_active_quads () const
 
template<>
unsigned int n_active_quads (const unsigned int) const
 
template<>
unsigned int n_active_quads () const
 
template<>
unsigned int n_hexs () const
 
template<>
unsigned int n_hexs (const unsigned int level) const
 
template<>
unsigned int n_active_hexs () const
 
template<>
unsigned int n_active_hexs (const unsigned int level) const
 
template<>
unsigned int max_adjacent_cells () const
 
template<>
unsigned int max_adjacent_cells () const
 
template<>
unsigned int max_adjacent_cells () const
 
template<>
unsigned int n_raw_quads (const unsigned int) const
 
template<>
unsigned int n_raw_quads (const unsigned int) const
 
template<>
unsigned int n_raw_quads (const unsigned int) const
 
template<>
unsigned int n_raw_quads (const unsigned int level) const
 
template<>
unsigned int n_raw_quads (const unsigned int level) const
 
template<>
unsigned int n_raw_quads (const unsigned int) const
 
template<>
unsigned int n_raw_quads () const
 
template<>
unsigned int n_raw_hexs (const unsigned int) const
 
template<>
unsigned int n_raw_hexs (const unsigned int) const
 
template<>
unsigned int n_raw_hexs (const unsigned int) const
 
template<>
unsigned int n_raw_hexs (const unsigned int level) const
 
Internal information about the number of objects
template<class Archive >
void save (Archive &ar, const unsigned int version) const
 
template<class Archive >
void load (Archive &ar, const unsigned int version)
 
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
 
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
 
template<class Archive >
void serialize (Archive &archive, const unsigned int version)
 
Mesh refinement
void set_all_refine_flags ()
 
void refine_global (const unsigned int times=1)
 
void coarsen_global (const unsigned int times=1)
 
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
 
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
 
unsigned int n_active_cells (const unsigned int level) const
 
unsigned int n_faces () const
 
unsigned int n_active_faces () const
 
unsigned int n_levels () 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
 
Triangulation< dim, spacedim > & get_triangulation ()
 
const Triangulation< dim, spacedim > & get_triangulation () 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)
 

Public Attributes

Keeping up with what happens to a triangulation
Signals signals
 

Static Public Attributes

static const unsigned int dimension = dim
 
static const unsigned int space_dimension = spacedim
 

Protected Types

using cell_relation_t = typename std::pair< cell_iterator, CellStatus >
 

Protected Member Functions

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_number_cache ()
 
void update_reference_cells () override
 
void reset_global_cell_indices ()
 

Protected Attributes

std::vector< cell_relation_tlocal_cell_relations
 
CellAttachedData cell_attached_data
 
DataTransfer data_transfer
 
const MPI_Comm mpi_communicator
 
types::subdomain_id my_subdomain
 
types::subdomain_id n_subdomains
 
NumberCache number_cache
 

Private Member Functions

virtual unsigned int coarse_cell_id_to_coarse_cell_index (const types::coarse_cell_id coarse_cell_id) const override
 
virtual types::coarse_cell_id coarse_cell_index_to_coarse_cell_id (const unsigned int coarse_cell_index) const override
 
virtual void update_cell_relations () override
 

Private Attributes

TriangulationDescription::Settings settings
 
std::function< void(::Triangulation< dim, spacedim > &, const unsigned int)> partitioner
 
std::vector< std::pair< types::coarse_cell_id, unsigned int > > coarse_cell_id_to_coarse_cell_index_vector
 
std::vector< types::coarse_cell_idcoarse_cell_index_to_coarse_cell_id_vector
 
bool currently_processing_create_triangulation_for_internal_usage
 
bool currently_processing_prepare_coarsening_and_refinement_for_internal_usage
 

Exceptions

MeshSmoothing smooth_grid
 
std::vector< ReferenceCellreference_cells
 
void update_periodic_face_map ()
 
static ::ExceptionBaseExcInvalidLevel (int arg1, int arg2)
 
static ::ExceptionBaseExcTriangulationNotEmpty (int arg1, int arg2)
 
static ::ExceptionBaseExcGridReadError ()
 
static ::ExceptionBaseExcFacesHaveNoLevel ()
 
static ::ExceptionBaseExcEmptyLevel (int arg1)
 
static ::ExceptionBaseExcNonOrientableTriangulation ()
 
static ::ExceptionBaseExcBoundaryIdNotFound (types::boundary_id arg1)
 
static ::ExceptionBaseExcInconsistentCoarseningFlags ()
 
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)
 

Detailed Description

template<int dim, int spacedim = dim>
class parallel::fullydistributed::Triangulation< dim, spacedim >

A distributed triangulation with a distributed coarse grid.

The motivation for parallel::fullydistributed::Triangulation has its origins in the following observations about complex geometries and/or about given meshes created by an external mesh generator. We regard complex geometries as geometries that can be meshed only with a non-negligible number of coarse cells (>10,000):

To be able to construct a fully partitioned triangulation that distributes the coarse grid and gives flexibility regarding partitioning, the following ingredients are required:

The ingredients listed above are bundled in the struct TriangulationDescription::Description. The user has to fill this data structure - in a pre-processing step - before actually creating the triangulation. Predefined functions to create TriangulationDescription::Description can be found in the namespace TriangulationDescription::Utilities.

Once the TriangulationDescription::Description construction_data has been constructed, the triangulation tria can be created by calling tria.create_triangulation(construction_data);.

Note
This triangulation supports: 1D/2D/3D, hanging nodes, geometric multigrid, and periodicity.
You can create a triangulation with hanging nodes and multigrid levels with create_triangulation(). However, once it has been created, it cannot be altered anymore, i.e. you cannot coarsen or refine afterwards.
Currently only simple periodicity conditions (i.e. without offsets and rotation matrices - see also the documentation of GridTools::collect_periodic_faces()) are supported.

Definition at line 113 of file fully_distributed_tria.h.

Member Typedef Documentation

◆ cell_iterator

template<int dim, int spacedim = dim>
using parallel::fullydistributed::Triangulation< dim, spacedim >::cell_iterator = typename ::Triangulation<dim, spacedim>::cell_iterator

Definition at line 118 of file fully_distributed_tria.h.

◆ active_cell_iterator

template<int dim, int spacedim = dim>
using parallel::fullydistributed::Triangulation< dim, spacedim >::active_cell_iterator = typename ::Triangulation<dim, spacedim>::active_cell_iterator

Definition at line 121 of file fully_distributed_tria.h.

◆ CellStatus

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

Definition at line 124 of file fully_distributed_tria.h.

◆ cell_relation_t

template<int dim, int spacedim = dim>
using parallel::DistributedTriangulationBase< dim, spacedim >::cell_relation_t = typename std::pair<cell_iterator, CellStatus>
protectedinherited

Auxiliary data structure for assigning a CellStatus to a deal.II cell iterator. For an extensive description of the former, see the documentation for the member function register_data_attach().

Definition at line 689 of file tria_base.h.

◆ level_cell_iterator

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

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

Definition at line 1361 of file tria.h.

Member Enumeration Documentation

◆ MeshSmoothing

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

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:

limit_level_difference_at_vertices.png

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 1136 of file tria.h.

Constructor & Destructor Documentation

◆ Triangulation()

template<int dim, int spacedim>
Triangulation< dim, spacedim >::Triangulation ( const MPI_Comm mpi_communicator)
explicit

Constructor.

Parameters
mpi_communicatorThe MPI communicator to be used for the triangulation.

Definition at line 45 of file fully_distributed_tria.cc.

◆ ~Triangulation()

template<int dim, int spacedim = dim>
virtual parallel::fullydistributed::Triangulation< dim, spacedim >::~Triangulation ( )
virtualdefault

Destructor.

Reimplemented from Triangulation< dim, spacedim >.

Member Function Documentation

◆ create_triangulation() [1/3]

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

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 GridReordering class.

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.
This is the function to be used instead of Triangulation::create_triangulation() for some of the other triangulations of deal.II.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 62 of file fully_distributed_tria.cc.

◆ create_triangulation() [2/3]

template<int dim, int spacedim>
void Triangulation< dim, spacedim >::create_triangulation ( const std::vector< Point< spacedim >> &  vertices,
const std::vector<::CellData< dim >> &  cells,
const SubCellData subcelldata 
)
overridevirtual
Note
This function is not implemented for this class and throws an assertion. Instead, use the other create_triangulation() function to create the triangulation.

Definition at line 219 of file fully_distributed_tria.cc.

◆ copy_triangulation() [1/3]

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

Implementation of the same function as in the base class.

Parameters
other_triaThe triangulation to be copied. It can be a serial Triangulation or a parallel::distributed::Triangulation. Both can have been refined already.
Note
This function uses the partitioner registered with set_partitioner().

Definition at line 248 of file fully_distributed_tria.cc.

◆ set_partitioner()

template<int dim, int spacedim>
void Triangulation< dim, spacedim >::set_partitioner ( const std::function< void(::Triangulation< dim, spacedim > &, const unsigned int)> &  partitioner,
const TriangulationDescription::Settings settings 
)

Register a partitioner, which is used within the method copy_triangulation.

Parameters
partitionerA partitioning function, which takes as input argument a reference to the triangulation to be partitioned and the number of partitions to be created. The function needs to set subdomain ids for each active cell of the given triangulation, with values between zero (inclusive) and the second argument to the function (exclusive).
settingsSee the description of the Settings enumerator.
Note
As a default, GridTools::partition_triangulation_zorder() is used as partitioner and data structures on multigrid levels are not set up.

Definition at line 302 of file fully_distributed_tria.cc.

◆ execute_coarsening_and_refinement()

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

Coarsen and refine the mesh according to refinement and coarsening flags set.

Note
Not implemented yet.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 315 of file fully_distributed_tria.cc.

◆ prepare_coarsening_and_refinement()

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

Override the implementation of prepare_coarsening_and_refinement from the base class.

Note
Not implemented yet.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 324 of file fully_distributed_tria.cc.

◆ has_hanging_nodes()

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

Return true if the triangulation has hanging nodes.

Note
Not implemented yet.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 338 of file fully_distributed_tria.cc.

◆ memory_consumption()

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

Return the local memory consumption in bytes.

Reimplemented from parallel::TriangulationBase< dim, spacedim >.

Definition at line 348 of file fully_distributed_tria.cc.

◆ is_multilevel_hierarchy_constructed()

template<int dim, int spacedim>
bool Triangulation< dim, spacedim >::is_multilevel_hierarchy_constructed ( ) const
overridevirtual

Return if multilevel hierarchy is supported and has been constructed.

Implements parallel::TriangulationBase< dim, spacedim >.

Definition at line 364 of file fully_distributed_tria.cc.

◆ save() [1/2]

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

Save the triangulation into the given file. This file needs to be reachable from all nodes in the computation on a shared network file system. See the SolutionTransfer class on how to store solution vectors into this file. Additional cell-based data can be saved using register_data_attach().

Implements parallel::DistributedTriangulationBase< dim, spacedim >.

Definition at line 426 of file fully_distributed_tria.cc.

◆ load() [1/2]

template<int dim, int spacedim>
void Triangulation< dim, spacedim >::load ( const std::string &  filename,
const bool  autopartition = false 
)
overridevirtual

Load the triangulation saved with save() back in. The mesh must be empty before calling this function.

You need to load with the same number of MPI processes that you saved with, hence autopartition is disabled.

Cell-based data that was saved with register_data_attach() can be read in with notify_ready_to_unpack() after calling load().

Implements parallel::DistributedTriangulationBase< dim, spacedim >.

Definition at line 561 of file fully_distributed_tria.cc.

◆ coarse_cell_id_to_coarse_cell_index()

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

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.
Parameters
coarse_cell_idUnique id of the coarse cell.
Returns
Index of the coarse cell within the current triangulation.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 375 of file fully_distributed_tria.cc.

◆ coarse_cell_index_to_coarse_cell_id()

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

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 from Triangulation< dim, spacedim >.

Definition at line 394 of file fully_distributed_tria.cc.

◆ update_cell_relations()

template<int dim, int spacedim>
void Triangulation< dim, spacedim >::update_cell_relations ( )
overrideprivatevirtual

Go through all active cells that are locally owned and record how they will change in the private member vector local_cell_relations.

As no adaptive mesh refinement is supported at the moment for this class, all cells will be flagged with the CellStatus CELL_PERSIST. These relations will currently only be used for serialization.

The stored vector will have a size equal to the number of locally owned active cells and will be ordered by the occurrence of those cells.

Implements parallel::DistributedTriangulationBase< dim, spacedim >.

Definition at line 410 of file fully_distributed_tria.cc.

◆ clear()

template<int dim, int spacedim>
void parallel::DistributedTriangulationBase< dim, spacedim >::clear ( )
overridevirtualinherited

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 from Triangulation< dim, spacedim >.

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

Definition at line 670 of file tria_base.cc.

◆ save() [2/2]

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

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() [2/2]

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

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.

◆ register_data_attach()

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

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:

  • CELL_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.
  • CELL_REFINE: 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 CELL_REFINE, whereas the others will be marked with 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 CELL_REFINE, 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.
  • CELL_COARSEN: 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.
  • CELL_INVALID: See CELL_REFINE.
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 CELL_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.

Definition at line 753 of file tria_base.cc.

◆ notify_ready_to_unpack()

template<int dim, int spacedim>
void parallel::DistributedTriangulationBase< 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 
)
inherited

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 CELL_REFINE but no longer for callbacks with status CELL_COARSEN.

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().)

Definition at line 782 of file tria_base.cc.

◆ save_attached_data()

template<int dim, int spacedim>
void parallel::DistributedTriangulationBase< dim, spacedim >::save_attached_data ( const unsigned int  global_first_cell,
const unsigned int  global_num_cells,
const std::string &  filename 
) const
protectedinherited

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.

Definition at line 682 of file tria_base.cc.

◆ load_attached_data()

template<int dim, int spacedim>
void parallel::DistributedTriangulationBase< 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 
)
protectedinherited

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.

Definition at line 718 of file tria_base.cc.

◆ get_communicator()

template<int dim, int spacedim>
MPI_Comm parallel::TriangulationBase< dim, spacedim >::get_communicator ( ) const
overridevirtualinherited

Return MPI communicator used by this triangulation.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 140 of file tria_base.cc.

◆ copy_triangulation() [2/3]

template<int dim, int spacedim>
void parallel::TriangulationBase< dim, spacedim >::copy_triangulation ( const ::Triangulation< dim, spacedim > &  old_tria)
overridevirtualinherited

Implementation of the same function as in the base class.

Note
This function copies the cells, but not the communicator, of the source triangulation. In other words, the resulting triangulation will operate on the communicator it was constructed with.

Definition at line 65 of file tria_base.cc.

◆ copy_triangulation() [3/3]

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

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.

Keep in mind that this function also copies the pointer to the boundary descriptor previously set by the set_manifold function. You must therefore also guarantee that the Manifold objects describing the boundary have a lifetime at least as long as the copied triangulation.

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 >.

Definition at line 11511 of file tria.cc.

◆ n_locally_owned_active_cells()

template<int dim, int spacedim>
unsigned int parallel::TriangulationBase< dim, spacedim >::n_locally_owned_active_cells ( ) const
inherited

Return the number of active cells in the triangulation that are locally owned, i.e. that have a subdomain_id equal to locally_owned_subdomain(). Note that there may be more active cells in the triangulation stored on the present processor, such as for example ghost cells, or cells further away from the locally owned block of cells but that are needed to ensure that the triangulation that stores this processor's set of active cells still remains balanced with respect to the 2:1 size ratio of adjacent cells.

As a consequence of the remark above, the result of this function is always smaller or equal to the result of the function with the same name in the Triangulation base class, which includes the active ghost and artificial cells (see also GlossArtificialCell and GlossGhostCell).

Definition at line 119 of file tria_base.cc.

◆ n_global_active_cells()

template<int dim, int spacedim>
types::global_cell_index parallel::TriangulationBase< dim, spacedim >::n_global_active_cells ( ) const
overridevirtualinherited

Return the sum over all processors of the number of active cells owned by each processor. This equals the overall number of active cells in the triangulation.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 133 of file tria_base.cc.

◆ n_global_levels()

template<int dim, int spacedim>
unsigned int parallel::TriangulationBase< dim, spacedim >::n_global_levels ( ) const
overridevirtualinherited

Return the global maximum level. This may be bigger than the number Triangulation::n_levels() (a function in this class's base class) returns if the current processor only stores cells in parts of the domain that are not very refined, but if other processors store cells in more deeply refined parts of the domain.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 126 of file tria_base.cc.

◆ locally_owned_subdomain()

template<int dim, int spacedim>
types::subdomain_id parallel::TriangulationBase< dim, spacedim >::locally_owned_subdomain ( ) const
overridevirtualinherited

Return the subdomain id of those cells that are owned by the current processor. All cells in the triangulation that do not have this subdomain id are either owned by another processor or have children that only exist on other processors.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 320 of file tria_base.cc.

◆ ghost_owners()

template<int dim, int spacedim>
const std::set< types::subdomain_id > & parallel::TriangulationBase< dim, spacedim >::ghost_owners ( ) const
inherited

Return a set of MPI ranks of the processors that have at least one ghost cell adjacent to the cells of the local processor. In other words, this is the set of subdomain_id() for all ghost cells.

The returned sets are symmetric, that is if i is contained in the list of processor j, then j will also be contained in the list of processor i.

Definition at line 329 of file tria_base.cc.

◆ level_ghost_owners()

template<int dim, int spacedim>
const std::set< types::subdomain_id > & parallel::TriangulationBase< dim, spacedim >::level_ghost_owners ( ) const
inherited

Return a set of MPI ranks of the processors that have at least one level ghost cell adjacent to our cells used in geometric multigrid. In other words, this is the set of level_subdomain_id() for all level ghost cells.

The returned sets are symmetric, that is if i is contained in the list of processor j, then j will also be contained in the list of processor i.

Note
The level ghost owners can only be determined if the multigrid ownership has been assigned (by setting the construct_multigrid_hierarchy flag at construction time), otherwise the returned set will be empty.

Definition at line 338 of file tria_base.cc.

◆ global_active_cell_index_partitioner()

template<int dim, int spacedim>
const std::weak_ptr< const Utilities::MPI::Partitioner > parallel::TriangulationBase< dim, spacedim >::global_active_cell_index_partitioner ( ) const
inherited

Return partitioner for the global indices of the cells on the active level of the triangulation.

Definition at line 634 of file tria_base.cc.

◆ global_level_cell_index_partitioner()

template<int dim, int spacedim>
const std::weak_ptr< const Utilities::MPI::Partitioner > parallel::TriangulationBase< dim, spacedim >::global_level_cell_index_partitioner ( const unsigned int  level) const
inherited

Return partitioner for the global indices of the cells on the given level of the triangulation.

Definition at line 641 of file tria_base.cc.

◆ get_boundary_ids()

template<int dim, int spacedim>
std::vector< types::boundary_id > parallel::TriangulationBase< dim, spacedim >::get_boundary_ids ( ) const
overridevirtualinherited

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
Note
This function involves a global communication gathering all current IDs from all processes.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 347 of file tria_base.cc.

◆ get_manifold_ids()

template<int dim, int spacedim>
std::vector< types::manifold_id > parallel::TriangulationBase< dim, spacedim >::get_manifold_ids ( ) const
overridevirtualinherited

Return a vector containing all manifold indicators assigned to the objects of the active cells of this Triangulation. Note, that each manifold 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 manifold indicators
Note
This function involves a global communication gathering all current IDs from all processes.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 358 of file tria_base.cc.

◆ communicate_locally_moved_vertices()

template<int dim, int spacedim>
void parallel::TriangulationBase< dim, spacedim >::communicate_locally_moved_vertices ( const std::vector< bool > &  vertex_locally_moved)
inherited

When vertices have been moved locally, for example using code like

cell->vertex(0) = new_location;

then this function can be used to update the location of vertices between MPI processes.

All the vertices that have been moved and might be in the ghost layer of a process have to be reported in the vertex_locally_moved argument. This ensures that that part of the information that has to be send between processes is actually sent. Additionally, it is quite important that vertices on the boundary between processes are reported on exactly one process (e.g. the one with the highest id). Otherwise we could expect undesirable results if multiple processes move a vertex differently. A typical strategy is to let processor \(i\) move those vertices that are adjacent to cells whose owners include processor \(i\) but no other processor \(j\) with \(j<i\); in other words, for vertices at the boundary of a subdomain, the processor with the lowest subdomain id "owns" a vertex.

Note
It only makes sense to move vertices that are either located on locally owned cells or on cells in the ghost layer. This is because you can be sure that these vertices indeed exist on the finest mesh aggregated over all processors, whereas vertices on artificial cells but not at least in the ghost layer may or may not exist on the globally finest mesh. Consequently, the vertex_locally_moved argument may not contain vertices that aren't at least on ghost cells.
This function moves vertices in such a way that on every processor, the vertices of every locally owned and ghost cell is consistent with the corresponding location of these cells on other processors. On the other hand, the locations of artificial cells will in general be wrong since artificial cells may or may not exist on other processors and consequently it is not possible to determine their location in any way. This is not usually a problem since one never does anything on artificial cells. However, it may lead to problems if the mesh with moved vertices is refined in a later step. If that's what you want to do, the right way to do it is to save the offset applied to every vertex, call this function, and before refining or coarsening the mesh apply the opposite offset and call this function again.
Parameters
vertex_locally_movedA bitmap indicating which vertices have been moved. The size of this array must be equal to Triangulation::n_vertices() and must be a subset of those vertices flagged by GridTools::get_locally_owned_vertices().
See also
This function is used, for example, in GridTools::distort_random().

Definition at line 584 of file tria_base.cc.

◆ update_number_cache()

template<int dim, int spacedim>
void parallel::TriangulationBase< dim, spacedim >::update_number_cache ( )
protectedvirtualinherited

Update the number_cache variable after mesh creation or refinement.

Definition at line 148 of file tria_base.cc.

◆ update_reference_cells()

template<int dim, int spacedim>
void parallel::TriangulationBase< dim, spacedim >::update_reference_cells ( )
overrideprotectedvirtualinherited

Update the internal reference_cells vector.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 293 of file tria_base.cc.

◆ reset_global_cell_indices()

template<int dim, int spacedim>
void parallel::TriangulationBase< dim, spacedim >::reset_global_cell_indices ( )
protectedinherited

Reset global active cell indices and global level cell indices.

Definition at line 369 of file tria_base.cc.

◆ set_mesh_smoothing()

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

Set the mesh smoothing to mesh_smoothing. This overrides the MeshSmoothing given to the constructor. It is allowed to call this function only if the triangulation is empty.

Definition at line 11309 of file tria.cc.

◆ get_mesh_smoothing()

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

Return the mesh smoothing requirements that are obeyed.

Definition at line 11321 of file tria.cc.

◆ create_triangulation() [3/3]

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

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 GridReordering class.

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::distributed::Triangulation< dim, spacedim >, parallel::shared::Triangulation< dim, spacedim >, and PersistentTriangulation< dim, spacedim >.

Definition at line 11622 of file tria.cc.

◆ create_triangulation_compatibility()

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

For backward compatibility, only. This function takes the cell data in the ordering as requested by deal.II versions up to 5.2, converts it to the new (lexicographic) ordering and calls create_triangulation().

Note
This function internally calls create_triangulation and therefore can throw the same exception as the other function.

Reimplemented in PersistentTriangulation< dim, spacedim >.

Definition at line 11577 of file tria.cc.

◆ flip_all_direction_flags()

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

Revert or flip the direction_flags of a dim<spacedim triangulation, see GlossDirectionFlag.

This function throws an exception if dim equals spacedim.

Definition at line 11929 of file tria.cc.

◆ set_all_refine_flags()

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

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.

Definition at line 11941 of file tria.cc.

◆ refine_global()

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

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.

Definition at line 11957 of file tria.cc.

◆ coarsen_global()

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

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.

Definition at line 11970 of file tria.cc.

◆ save_refine_flags() [1/2]

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

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

Definition at line 12007 of file tria.cc.

◆ save_refine_flags() [2/2]

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

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

Definition at line 11990 of file tria.cc.

◆ load_refine_flags() [1/2]

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

Read the information stored by save_refine_flags.

Definition at line 12021 of file tria.cc.

◆ load_refine_flags() [2/2]

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

Read the information stored by save_refine_flags.

Definition at line 12032 of file tria.cc.

◆ save_coarsen_flags() [1/2]

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

Analogue to save_refine_flags.

Definition at line 12076 of file tria.cc.

◆ save_coarsen_flags() [2/2]

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

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

Definition at line 12059 of file tria.cc.

◆ load_coarsen_flags() [1/2]

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

Analogue to load_refine_flags.

Definition at line 12090 of file tria.cc.

◆ load_coarsen_flags() [2/2]

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

Analogue to load_refine_flags.

Definition at line 12104 of file tria.cc.

◆ get_anisotropic_refinement_flag()

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

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

Definition at line 12124 of file tria.cc.

◆ clear_user_flags()

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

Clear all user flags. See also GlossUserFlags.

Definition at line 12295 of file tria.cc.

◆ save_user_flags() [1/2]

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

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.

Definition at line 12306 of file tria.cc.

◆ save_user_flags() [2/2]

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

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.

Definition at line 12324 of file tria.cc.

◆ load_user_flags() [1/2]

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

Read the information stored by save_user_flags. See also GlossUserFlags.

Definition at line 12355 of file tria.cc.

◆ load_user_flags() [2/2]

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

Read the information stored by save_user_flags. See also GlossUserFlags.

Definition at line 12373 of file tria.cc.

◆ clear_user_flags_line()

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

Clear all user flags on lines. See also GlossUserFlags.

Definition at line 12206 of file tria.cc.

◆ save_user_flags_line() [1/2]

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

Save the user flags on lines. See also GlossUserFlags.

Definition at line 12425 of file tria.cc.

◆ save_user_flags_line() [2/2]

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

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.

Definition at line 12410 of file tria.cc.

◆ load_user_flags_line() [1/2]

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

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

Definition at line 12439 of file tria.cc.

◆ load_user_flags_line() [2/2]

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

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

Definition at line 12453 of file tria.cc.

◆ clear_user_flags_quad()

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

Clear all user flags on quads. See also GlossUserFlags.

Definition at line 12246 of file tria.cc.

◆ save_user_flags_quad() [1/2]

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

Save the user flags on quads. See also GlossUserFlags.

Definition at line 12548 of file tria.cc.

◆ save_user_flags_quad() [2/2]

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

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.

Definition at line 12529 of file tria.cc.

◆ load_user_flags_quad() [1/2]

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

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

Definition at line 12562 of file tria.cc.

◆ load_user_flags_quad() [2/2]

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

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

Definition at line 12576 of file tria.cc.

◆ clear_user_flags_hex()

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

Clear all user flags on quads. See also GlossUserFlags.

Definition at line 12286 of file tria.cc.

◆ save_user_flags_hex() [1/2]

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

Save the user flags on hexs. See also GlossUserFlags.

Definition at line 12617 of file tria.cc.

◆ save_user_flags_hex() [2/2]

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

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.

Definition at line 12598 of file tria.cc.

◆ load_user_flags_hex() [1/2]

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

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

Definition at line 12631 of file tria.cc.

◆ load_user_flags_hex() [2/2]

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

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

Definition at line 12645 of file tria.cc.

◆ clear_user_data()

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

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

Definition at line 12167 of file tria.cc.

◆ save_user_indices()

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

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

Definition at line 12667 of file tria.cc.

◆ load_user_indices()

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

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

Definition at line 12699 of file tria.cc.

◆ save_user_pointers()

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

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

Definition at line 12923 of file tria.cc.

◆ load_user_pointers()

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

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

Definition at line 12954 of file tria.cc.

◆ save_user_indices_line()

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

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

Definition at line 12779 of file tria.cc.

◆ load_user_indices_line()

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

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

Definition at line 12793 of file tria.cc.

◆ save_user_indices_quad()

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

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

Definition at line 12807 of file tria.cc.

◆ load_user_indices_quad()

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

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

Definition at line 12825 of file tria.cc.

◆ save_user_indices_hex()

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

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

Definition at line 12842 of file tria.cc.

◆ load_user_indices_hex()

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

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

Definition at line 12860 of file tria.cc.

◆ save_user_pointers_line()

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

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

Definition at line 12991 of file tria.cc.

◆ load_user_pointers_line()

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

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

Definition at line 13005 of file tria.cc.

◆ save_user_pointers_quad()

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

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

Definition at line 13020 of file tria.cc.

◆ load_user_pointers_quad()

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

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

Definition at line 13038 of file tria.cc.

◆ save_user_pointers_hex()

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

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

Definition at line 13055 of file tria.cc.

◆ load_user_pointers_hex()

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

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

Definition at line 13073 of file tria.cc.

◆ begin()

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

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.

Definition at line 13114 of file tria.cc.

◆ begin_active()

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

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)
{
...
}

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.

Definition at line 13134 of file tria.cc.

◆ end() [1/2]

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

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

Definition at line 13225 of file tria.cc.

◆ end() [2/2]

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

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.

Definition at line 13265 of file tria.cc.

◆ end_active()

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

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.

Definition at line 13294 of file tria.cc.

◆ last()

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

Return an iterator pointing to the last used cell.

Definition at line 13154 of file tria.cc.

◆ last_active()

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

Return an iterator pointing to the last active cell.

Definition at line 13179 of file tria.cc.

◆ create_cell_iterator()

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

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

If the given argument corresponds to a valid cell in this triangulation, this operation will always succeed 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 succeed for locally relevant cells, but may not for artificial cells that are less refined on the current processor.

Definition at line 13200 of file tria.cc.

◆ begin_face()

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

Iterator to the first used face.

Definition at line 13367 of file tria.cc.

◆ begin_active_face()

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

Iterator to the first active face.

Definition at line 13388 of file tria.cc.

◆ end_face()

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

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

Definition at line 13409 of file tria.cc.

◆ begin_vertex()

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

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

Definition at line 13442 of file tria.cc.

◆ begin_active_vertex()

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

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.

Definition at line 13459 of file tria.cc.

◆ end_vertex()

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

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.

Definition at line 13469 of file tria.cc.

◆ n_lines() [1/2]

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

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.

Definition at line 13965 of file tria.cc.

◆ n_lines() [2/2]

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

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

Definition at line 14003 of file tria.cc.

◆ n_active_lines() [1/2]

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

Return the total number of active lines.

Definition at line 14013 of file tria.cc.

◆ n_active_lines() [2/2]

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

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

Definition at line 14021 of file tria.cc.

◆ n_quads() [1/8]

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

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

Definition at line 14178 of file tria.cc.

◆ n_quads() [2/8]

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

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

Definition at line 14186 of file tria.cc.

◆ n_quads() [3/8]

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

Definition at line 14032 of file tria.cc.

◆ n_quads() [4/8]

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

Definition at line 14040 of file tria.cc.

◆ n_quads() [5/8]

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

Definition at line 14081 of file tria.cc.

◆ n_quads() [6/8]

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

Definition at line 14089 of file tria.cc.

◆ n_quads() [7/8]

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

Definition at line 14129 of file tria.cc.

◆ n_quads() [8/8]

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

Definition at line 14137 of file tria.cc.

◆ n_active_quads() [1/8]

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

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

Definition at line 14245 of file tria.cc.

◆ n_active_quads() [2/8]

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

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

Definition at line 14253 of file tria.cc.

◆ n_active_quads() [3/8]

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

Definition at line 14064 of file tria.cc.

◆ n_active_quads() [4/8]

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

Definition at line 14072 of file tria.cc.

◆ n_active_quads() [5/8]

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

Definition at line 14113 of file tria.cc.

◆ n_active_quads() [6/8]

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

Definition at line 14121 of file tria.cc.

◆ n_active_quads() [7/8]

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

Definition at line 14161 of file tria.cc.

◆ n_active_quads() [8/8]

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

Definition at line 14169 of file tria.cc.

◆ n_hexs() [1/4]

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

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

Definition at line 14264 of file tria.cc.

◆ n_hexs() [2/4]

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

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

Definition at line 14273 of file tria.cc.

◆ n_hexs() [3/4]

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

Definition at line 14307 of file tria.cc.

◆ n_hexs() [4/4]

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

Definition at line 14316 of file tria.cc.

◆ n_active_hexs() [1/4]

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

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

Definition at line 14290 of file tria.cc.

◆ n_active_hexs() [2/4]

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

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

Definition at line 14299 of file tria.cc.

◆ n_active_hexs() [3/4]

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

Definition at line 14336 of file tria.cc.

◆ n_active_hexs() [4/4]

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

Definition at line 14345 of file tria.cc.

◆ n_cells() [1/2]

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

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

Definition at line 13813 of file tria.cc.

◆ n_cells() [2/2]

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

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.

Definition at line 13914 of file tria.cc.

◆ n_active_cells() [1/2]

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

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

Definition at line 13821 of file tria.cc.

◆ n_active_cells() [2/2]

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

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

Definition at line 13934 of file tria.cc.

◆ n_faces()

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

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.

Definition at line 13837 of file tria.cc.

◆ n_active_faces()

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

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.

Definition at line 13875 of file tria.cc.

◆ n_levels()

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

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_vertices()

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

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
inherited

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>
unsigned int Triangulation< dim, spacedim >::n_used_vertices ( ) const
inherited

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

Definition at line 14356 of file tria.cc.

◆ vertex_used()

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

Return true if the vertex with this index is used.

◆ get_used_vertices()

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

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

Definition at line 14365 of file tria.cc.

◆ max_adjacent_cells() [1/4]

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

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.

Definition at line 14399 of file tria.cc.

◆ max_adjacent_cells() [2/4]

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

Definition at line 14374 of file tria.cc.

◆ max_adjacent_cells() [3/4]

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

Definition at line 14383 of file tria.cc.

◆ max_adjacent_cells() [4/4]

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

Definition at line 14391 of file tria.cc.

◆ get_triangulation() [1/2]

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

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).

Definition at line 14441 of file tria.cc.

◆ get_triangulation() [2/2]

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

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

Definition at line 14450 of file tria.cc.

◆ n_raw_lines() [1/2]

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

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).

Definition at line 13989 of file tria.cc.

◆ n_raw_lines() [2/2]

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

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).

Definition at line 13974 of file tria.cc.

◆ n_raw_quads() [1/9]

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

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).

Definition at line 14226 of file tria.cc.

◆ n_raw_quads() [2/9]

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

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_quads() [3/9]

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

Definition at line 14048 of file tria.cc.

◆ n_raw_quads() [4/9]

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

Definition at line 14097 of file tria.cc.

◆ n_raw_quads() [5/9]

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

Definition at line 14145 of file tria.cc.

◆ n_raw_quads() [6/9]

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

Definition at line 14197 of file tria.cc.

◆ n_raw_quads() [7/9]

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

Definition at line 14207 of file tria.cc.

◆ n_raw_quads() [8/9]

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

Definition at line 14216 of file tria.cc.

◆ n_raw_quads() [9/9]

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

Definition at line 14236 of file tria.cc.

◆ n_raw_hexs() [1/5]

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

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).

Definition at line 14282 of file tria.cc.

◆ n_raw_hexs() [2/5]

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

Definition at line 14056 of file tria.cc.

◆ n_raw_hexs() [3/5]

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

Definition at line 14105 of file tria.cc.

◆ n_raw_hexs() [4/5]

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

Definition at line 14153 of file tria.cc.

◆ n_raw_hexs() [5/5]

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

Definition at line 14327 of file tria.cc.

◆ n_raw_cells()

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

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).

Definition at line 13894 of file tria.cc.

◆ n_raw_faces()

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

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).

Definition at line 13856 of file tria.cc.

◆ add_periodicity()

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

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.

Definition at line 14459 of file tria.cc.

◆ get_periodic_face_map()

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

Return the periodic_face_map.

Definition at line 14479 of file tria.cc.

◆ get_reference_cells()

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

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

Definition at line 14680 of file tria.cc.

◆ all_reference_cells_are_hyper_cube()

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

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

Definition at line 14689 of file tria.cc.

◆ serialize()

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

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

◆ 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 136 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 156 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 301 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 318 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 204 of file subscriptor.cc.

Member Data Documentation

◆ settings

template<int dim, int spacedim = dim>
TriangulationDescription::Settings parallel::fullydistributed::Triangulation< dim, spacedim >::settings
private

store the Settings.

Definition at line 284 of file fully_distributed_tria.h.

◆ partitioner

template<int dim, int spacedim = dim>
std::function<void(::Triangulation<dim, spacedim> &, const unsigned int)> parallel::fullydistributed::Triangulation< dim, spacedim >::partitioner
private

Partitioner used in copy_triangulation().

Definition at line 291 of file fully_distributed_tria.h.

◆ coarse_cell_id_to_coarse_cell_index_vector

template<int dim, int spacedim = dim>
std::vector<std::pair<types::coarse_cell_id, unsigned int> > parallel::fullydistributed::Triangulation< dim, spacedim >::coarse_cell_id_to_coarse_cell_index_vector
private

Sorted list of pairs of coarse-cell ids and their indices.

Definition at line 297 of file fully_distributed_tria.h.

◆ coarse_cell_index_to_coarse_cell_id_vector

template<int dim, int spacedim = dim>
std::vector<types::coarse_cell_id> parallel::fullydistributed::Triangulation< dim, spacedim >::coarse_cell_index_to_coarse_cell_id_vector
private

List of the coarse-cell id for each coarse cell (stored at cell->index()).

Definition at line 304 of file fully_distributed_tria.h.

◆ currently_processing_create_triangulation_for_internal_usage

template<int dim, int spacedim = dim>
bool parallel::fullydistributed::Triangulation< dim, spacedim >::currently_processing_create_triangulation_for_internal_usage
private

Boolean indicating that the function create_triangulation() was called for internal usage.

Definition at line 310 of file fully_distributed_tria.h.

◆ currently_processing_prepare_coarsening_and_refinement_for_internal_usage

template<int dim, int spacedim = dim>
bool parallel::fullydistributed::Triangulation< dim, spacedim >::currently_processing_prepare_coarsening_and_refinement_for_internal_usage
private

Boolean indicating that the function prepare_coarsening_and_refinement() was called for internal usage.

Definition at line 317 of file fully_distributed_tria.h.

◆ local_cell_relations

template<int dim, int spacedim = dim>
std::vector<cell_relation_t> parallel::DistributedTriangulationBase< dim, spacedim >::local_cell_relations
protectedinherited

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 696 of file tria_base.h.

◆ cell_attached_data

template<int dim, int spacedim = dim>
CellAttachedData parallel::DistributedTriangulationBase< dim, spacedim >::cell_attached_data
protectedinherited

Definition at line 729 of file tria_base.h.

◆ data_transfer

template<int dim, int spacedim = dim>
DataTransfer parallel::DistributedTriangulationBase< dim, spacedim >::data_transfer
protectedinherited

Definition at line 884 of file tria_base.h.

◆ mpi_communicator

template<int dim, int spacedim = dim>
const MPI_Comm parallel::TriangulationBase< dim, spacedim >::mpi_communicator
protectedinherited

MPI communicator to be used for the triangulation. We create a unique communicator for this class, which is a duplicate of the one passed to the constructor.

Definition at line 298 of file tria_base.h.

◆ my_subdomain

template<int dim, int spacedim = dim>
types::subdomain_id parallel::TriangulationBase< dim, spacedim >::my_subdomain
protectedinherited

The subdomain id to be used for the current processor. This is the MPI rank.

Definition at line 304 of file tria_base.h.

◆ n_subdomains

template<int dim, int spacedim = dim>
types::subdomain_id parallel::TriangulationBase< dim, spacedim >::n_subdomains
protectedinherited

The total number of subdomains (or the size of the MPI communicator).

Definition at line 309 of file tria_base.h.

◆ number_cache

template<int dim, int spacedim = dim>
NumberCache parallel::TriangulationBase< dim, spacedim >::number_cache
protectedinherited

Definition at line 358 of file tria_base.h.

◆ dimension

template<int dim, int spacedim = dim>
const unsigned int Triangulation< dim, spacedim >::dimension = dim
staticinherited

Make the dimension available in function templates.

Definition at line 1542 of file tria.h.

◆ space_dimension

template<int dim, int spacedim = dim>
const unsigned int Triangulation< dim, spacedim >::space_dimension = spacedim
staticinherited

Make the space-dimension available in function templates.

Definition at line 1547 of file tria.h.

◆ signals

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

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

Definition at line 2289 of file tria.h.


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