Reference documentation for deal.II version Git 16b398077e 2020-12-04 17:16:12 -0500
\(\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 Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
parallel::shared::Triangulation< dim, spacedim > Class Template Reference

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

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

Public Types

enum  Settings {
  partition_auto = 0x0, partition_metis = 0x1, partition_zorder = 0x2, partition_zoltan = 0x3,
  partition_custom_signal = 0x4, construct_multigrid_hierarchy = 0x8
}
 
using active_cell_iterator = typename ::Triangulation< dim, spacedim >::active_cell_iterator
 
using cell_iterator = typename ::Triangulation< dim, spacedim >::cell_iterator
 
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, const typename ::Triangulation< dim, spacedim >::MeshSmoothing=(::Triangulation< dim, spacedim >::none), const bool allow_artificial_cells=false, const Settings settings=partition_auto)
 
virtual ~Triangulation () override=default
 
virtual bool is_multilevel_hierarchy_constructed () const override
 
virtual void execute_coarsening_and_refinement () override
 
virtual void create_triangulation (const std::vector< Point< spacedim >> &vertices, const std::vector< CellData< dim >> &cells, const SubCellData &subcelldata) override
 
virtual void create_triangulation (const TriangulationDescription::Description< dim, spacedim > &construction_data) override
 
virtual void copy_triangulation (const ::Triangulation< dim, spacedim > &other_tria) override
 
template<class Archive >
void load (Archive &ar, const unsigned int version)
 
const std::vector< types::subdomain_id > & get_true_subdomain_ids_of_cells () const
 
const std::vector< types::subdomain_id > & get_true_level_subdomain_ids_of_cells (const unsigned int level) const
 
bool with_artificial_cells () const
 
virtual const MPI_Commget_communicator () const
 
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 std::size_t memory_consumption () 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 Utilities::MPI::Partitionerglobal_active_cell_index_partitioner () const
 
const Utilities::MPI::Partitionerglobal_level_cell_index_partitioner (const unsigned int level) const
 
virtual std::map< unsigned int, std::set<::types::subdomain_id > > compute_vertices_with_ghost_neighbors () const
 
virtual std::vector< types::boundary_idget_boundary_ids () const override
 
virtual std::vector< types::manifold_idget_manifold_ids () const override
 
virtual void clear ()
 
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_compatibility (const std::vector< Point< spacedim >> &vertices, const std::vector< CellData< dim >> &cells, const SubCellData &subcelldata)
 
void flip_all_direction_flags ()
 
template<>
bool prepare_coarsening_and_refinement ()
 
template<>
bool prepare_coarsening_and_refinement ()
 
template<>
bool prepare_coarsening_and_refinement ()
 
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
 
Mesh refinement
void set_all_refine_flags ()
 
void refine_global (const unsigned int times=1)
 
void coarsen_global (const unsigned int times=1)
 
virtual bool prepare_coarsening_and_refinement ()
 
History of a triangulation
void save_refine_flags (std::ostream &out) const
 
void save_refine_flags (std::vector< bool > &v) const
 
void load_refine_flags (std::istream &in)
 
void load_refine_flags (const std::vector< bool > &v)
 
void save_coarsen_flags (std::ostream &out) const
 
void save_coarsen_flags (std::vector< bool > &v) const
 
void load_coarsen_flags (std::istream &out)
 
void load_coarsen_flags (const std::vector< bool > &v)
 
bool get_anisotropic_refinement_flag () const
 
User data
void clear_user_flags ()
 
void save_user_flags (std::ostream &out) const
 
void save_user_flags (std::vector< bool > &v) const
 
void load_user_flags (std::istream &in)
 
void load_user_flags (const std::vector< bool > &v)
 
void clear_user_flags_line ()
 
void save_user_flags_line (std::ostream &out) const
 
void save_user_flags_line (std::vector< bool > &v) const
 
void load_user_flags_line (std::istream &in)
 
void load_user_flags_line (const std::vector< bool > &v)
 
void clear_user_flags_quad ()
 
void save_user_flags_quad (std::ostream &out) const
 
void save_user_flags_quad (std::vector< bool > &v) const
 
void load_user_flags_quad (std::istream &in)
 
void load_user_flags_quad (const std::vector< bool > &v)
 
void clear_user_flags_hex ()
 
void save_user_flags_hex (std::ostream &out) const
 
void save_user_flags_hex (std::vector< bool > &v) const
 
void load_user_flags_hex (std::istream &in)
 
void load_user_flags_hex (const std::vector< bool > &v)
 
void clear_user_data ()
 
void save_user_indices (std::vector< unsigned int > &v) const
 
void load_user_indices (const std::vector< unsigned int > &v)
 
void save_user_pointers (std::vector< void *> &v) const
 
void load_user_pointers (const std::vector< void *> &v)
 
void save_user_indices_line (std::vector< unsigned int > &v) const
 
void load_user_indices_line (const std::vector< unsigned int > &v)
 
void save_user_indices_quad (std::vector< unsigned int > &v) const
 
void load_user_indices_quad (const std::vector< unsigned int > &v)
 
void save_user_indices_hex (std::vector< unsigned int > &v) const
 
void load_user_indices_hex (const std::vector< unsigned int > &v)
 
void save_user_pointers_line (std::vector< void *> &v) const
 
void load_user_pointers_line (const std::vector< void *> &v)
 
void save_user_pointers_quad (std::vector< void *> &v) const
 
void load_user_pointers_quad (const std::vector< void *> &v)
 
void save_user_pointers_hex (std::vector< void *> &v) const
 
void load_user_pointers_hex (const std::vector< void *> &v)
 
Cell iterator functions
cell_iterator begin (const unsigned int level=0) const
 
active_cell_iterator begin_active (const unsigned int level=0) const
 
cell_iterator end () const
 
cell_iterator end (const unsigned int level) const
 
active_cell_iterator end_active (const unsigned int level) const
 
cell_iterator last () const
 
active_cell_iterator last_active () const
 
Cell iterator 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
 
virtual bool has_hanging_nodes () const
 
unsigned int n_vertices () const
 
const std::vector< Point< spacedim > > & get_vertices () const
 
unsigned int n_used_vertices () const
 
bool vertex_used (const unsigned int index) const
 
const std::vector< bool > & get_used_vertices () const
 
unsigned int max_adjacent_cells () const
 
Triangulation< dim, spacedim > & get_triangulation ()
 
const Triangulation< dim, spacedim > & get_triangulation () const
 
Internal information about the number of objects
unsigned int n_raw_lines () const
 
unsigned int n_raw_lines (const unsigned int level) const
 
unsigned int n_raw_quads () const
 
unsigned int n_raw_quads (const unsigned int level) const
 
unsigned int n_raw_hexs (const unsigned int level) const
 
unsigned int n_raw_cells (const unsigned int level) const
 
unsigned int n_raw_faces () const
 
template<class Archive >
void save (Archive &ar, const unsigned int version) 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::Type > & get_reference_cell_types () const
 
bool all_reference_cell_types_are_hyper_cube () const
 
template<class Archive >
void serialize (Archive &archive, const unsigned int version)
 
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)
 

Static Public Attributes

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

Protected Member Functions

virtual void update_number_cache ()
 
void update_reference_cell_types () override
 
void reset_global_cell_indices ()
 

Protected Attributes

MPI_Comm mpi_communicator
 
types::subdomain_id my_subdomain
 
types::subdomain_id n_subdomains
 
NumberCache number_cache
 

Private Member Functions

void partition ()
 

Private Attributes

const Settings settings
 
const bool allow_artificial_cells
 
std::vector< types::subdomain_idtrue_subdomain_ids_of_cells
 
std::vector< std::vector< types::subdomain_id > > true_level_subdomain_ids_of_cells
 

Keeping up with what happens to a triangulation

enum  CellStatus { CELL_PERSIST, CELL_REFINE, CELL_COARSEN, CELL_INVALID }
 
Signals signals
 

Exceptions

MeshSmoothing smooth_grid
 
std::vector< ReferenceCell::Typereference_cell_types
 
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::shared::Triangulation< dim, spacedim >

This class provides a parallel triangulation for which every processor knows about every cell of the global mesh (unlike for the parallel::distributed::Triangulation class) but in which cells are automatically partitioned when run with MPI so that each processor "owns" a subset of cells. The use of this class is demonstrated in step-18.

Different from the parallel::distributed::Triangulation and parallel::fullydistributed::Triangulation classes, this implies that the entire mesh is stored on each processor. While this is clearly a memory bottleneck that limits the use of this class to a few dozen or hundreds of MPI processes, the partitioning of the mesh can be used to partition work such as assembly or postprocessing between participating processors, and it can also be used to partition which processor stores which parts of matrices and vectors. As a consequence, using this class is often a gentler introduction to parallelizing a code than the more involved parallel::distributed::Triangulation class in which processors only know their own part of the mesh, but nothing about cells owned by other processors with the exception of a single layer of ghost cells around their own part of the domain.

The class is also useful in cases where compute time and memory considerations dictate that the program needs to be run in parallel, but where algorithmic concerns require that every processor knows about the entire mesh. An example could be where an application has to have both volume and surface meshes that can then both be partitioned independently, but where it is difficult to ensure that the locally owned set of surface mesh cells is adjacent to the locally owned set of volume mesh cells and the other way around. In such cases, knowing the entirety of both meshes ensures that assembly of coupling terms can be implemented without also implementing overly complicated schemes to transfer information about adjacent cells from processor to processor.

The partitioning of cells between processors is done internally based on a number of different possibilities. By passing appropriate flags to the constructor of this class (see the parallel::shared::Triangulation::Settings enum), it is possible to select different ways of partitioning the mesh, including ways that are dictated by the application and not by the desire to minimize the length of the interface between subdomains owned by processors (as is done by the METIS and Zoltan packages, both of which are options for partitioning). The DoFHandler class knows how to enumerate degrees of freedom in ways appropriate for the partitioned mesh.

Definition at line 101 of file shared_tria.h.

Member Typedef Documentation

◆ active_cell_iterator

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

Definition at line 106 of file shared_tria.h.

◆ cell_iterator

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

Definition at line 108 of file shared_tria.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 1363 of file tria.h.

Member Enumeration Documentation

◆ Settings

template<int dim, int spacedim = dim>
enum parallel::shared::Triangulation::Settings

Configuration flags for distributed Triangulations to be set in the constructor. Settings can be combined using bitwise OR.

The constructor requires that exactly one of partition_auto, partition_metis, partition_zorder, partition_zoltan and partition_custom_signal is set. If partition_auto is chosen, it will use partition_zoltan (if available), then partition_metis (if available) and finally partition_zorder.

Enumerator
partition_auto 

Choose the partitioner depending on the enabled dependencies that were found when configuring deal.II. In particular, if the Trilinos package Zoltan was found, then use the partition_zoltan strategy. If Zoltan was not found but the METIS package was found, then use the partition_metis strategy. If neither of these were found, then use the partition_zorder partitioning strategy.

partition_metis 

Use METIS partitioner to partition active cells.

partition_zorder 

Partition active cells with the same scheme used in the p4est library.

The term "Z-order" originates in the fact that cells are sorted using a space filling curve which in 2d connects the four children of a cell in the order bottom left, bottom right, top left, top right (i.e., with a curve that looks like a reverse "Z"), and does so recursively on all levels of a triangulation. This is also the order in which children are enumerated by the GeometryInfo class. The "Z-order" is also sometimes called "Morton ordering", see https://en.wikipedia.org/wiki/Z-order_curve .

See also
Z order glossary entry.
partition_zoltan 

Use Zoltan to partition active cells.

partition_custom_signal 

Partition cells using a custom, user defined function. This is accomplished by connecting the post_refinement signal to the triangulation whenever it is first created and passing the user defined function through the signal using std::bind. Here is an example:

template <int dim>
void mypartition(parallel::shared::Triangulation<dim> &tria)
{
// user defined partitioning scheme: assign subdomain_ids
// round-robin in a mostly random way:
std::vector<unsigned int> assignment =
{0,0,1,2,0,0,2,1,0,2,2,1,2,2,0,0};
unsigned int index = 0;
for (const auto &cell : tria.active_cell_iterators())
cell->set_subdomain_id(assignment[(index++)%16]);
}
int main ()
{
...,
tria.signals.post_refinement.connect(std::bind(&mypartition<dim>,
std::ref(tria)));
}

An equivalent code using lambda functions would look like this:

int main ()
{
...,
tria.signals.post_refinement.connect (
[&tria]()
{
// user defined partitioning scheme as above
...
});
}
Note
If you plan to use a custom partition with geometric multigrid, you must manually partition the level cells in addition to the active cells.
construct_multigrid_hierarchy 

This flag needs to be set to use the geometric multigrid functionality. This option requires additional computation and communication.

Note: This flag should always be set alongside a flag for an active cell partitioning method.

Definition at line 123 of file shared_tria.h.

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

◆ CellStatus

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

Used to inform functions in derived classes how the cell with the given cell_iterator is going to change. Note that this may me different than the refine_flag() and coarsen_flag() in the cell_iterator in parallel calculations because of refinement constraints that this machine does not see.

Enumerator
CELL_PERSIST 

The cell will not be refined or coarsened and might or might not move to a different processor.

CELL_REFINE 

The cell will be or was refined.

CELL_COARSEN 

The children of this cell will be or were coarsened into this cell.

CELL_INVALID 

Invalid status. Will not occur for the user.

Definition at line 2017 of file tria.h.

Constructor & Destructor Documentation

◆ Triangulation()

template<int dim, int spacedim>
Triangulation< dim, spacedim >::Triangulation ( const MPI_Comm mpi_communicator,
const typename ::Triangulation< dim, spacedim >::MeshSmoothing  smooth_grid = (::Triangulation<dim, spacedim>::none),
const bool  allow_artificial_cells = false,
const Settings  settings = partition_auto 
)

Constructor.

If allow_aritifical_cells is true, this class will behave similar to parallel::distributed::Triangulation in that there will be locally owned, ghost and artificial cells.

Otherwise all non-locally owned cells are considered ghost.

Definition at line 39 of file shared_tria.cc.

◆ ~Triangulation()

template<int dim, int spacedim = dim>
virtual parallel::shared::Triangulation< dim, spacedim >::~Triangulation ( )
overridevirtualdefault

Destructor.

Reimplemented from Triangulation< dim, spacedim >.

Member Function Documentation

◆ 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 74 of file shared_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.

This step is equivalent to the Triangulation class with an addition of calling GridTools::partition_triangulation() at the end.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 357 of file shared_tria.cc.

◆ create_triangulation() [1/2]

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

Create a triangulation.

This function also partitions triangulation based on the MPI communicator provided to the constructor.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 368 of file shared_tria.cc.

◆ create_triangulation() [2/2]

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

Create a triangulation.

This function also partitions triangulation based on the MPI communicator provided to the constructor.

Note
Not implemented yet.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 394 of file shared_tria.cc.

◆ copy_triangulation() [1/3]

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

Copy other_tria to this triangulation.

This function also partitions triangulation based on the MPI communicator provided to the constructor.

Note
This function can not be used with parallel::distributed::Triangulation, since it only stores those cells that it owns, one layer of ghost cells around the ones it locally owns, and a number of artificial cells.

Definition at line 407 of file shared_tria.cc.

◆ load()

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

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

This function first does the same work as in Triangulation::load, then partitions the triangulation based on the MPI communicator provided to the constructor.

Definition at line 389 of file shared_tria.h.

◆ get_true_subdomain_ids_of_cells()

template<int dim, int spacedim>
const std::vector< types::subdomain_id > & Triangulation< dim, spacedim >::get_true_subdomain_ids_of_cells ( ) const

Return a vector of length Triangulation::n_active_cells() where each element stores the subdomain id of the owner of this cell. The elements of the vector are obviously the same as the subdomain ids for locally owned and ghost cells, but are also correct for artificial cells that do not store who the owner of the cell is in their subdomain_id field.

Definition at line 333 of file shared_tria.cc.

◆ get_true_level_subdomain_ids_of_cells()

template<int dim, int spacedim>
const std::vector< types::subdomain_id > & Triangulation< dim, spacedim >::get_true_level_subdomain_ids_of_cells ( const unsigned int  level) const

Return a vector of length Triangulation::n_cells(level) where each element stores the level subdomain id of the owner of this cell. The elements of the vector are obviously the same as the level subdomain ids for locally owned and ghost cells, but are also correct for artificial cells that do not store who the owner of the cell is in their level_subdomain_id field.

Definition at line 342 of file shared_tria.cc.

◆ with_artificial_cells()

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

Return allow_artificial_cells , namely true if artificial cells are allowed.

Definition at line 324 of file shared_tria.cc.

◆ partition()

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

This function calls GridTools::partition_triangulation () and if requested in the constructor of the class marks artificial cells.

Definition at line 83 of file shared_tria.cc.

◆ get_communicator()

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

Return MPI communicator used by this triangulation.

Definition at line 139 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 64 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 10312 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 118 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 132 of file tria_base.cc.

◆ memory_consumption()

template<int dim, int spacedim>
std::size_t parallel::TriangulationBase< dim, spacedim >::memory_consumption ( ) const
overridevirtualinherited

◆ 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 125 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 318 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 327 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 336 of file tria_base.cc.

◆ global_active_cell_index_partitioner()

template<int dim, int spacedim>
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 532 of file tria_base.cc.

◆ global_level_cell_index_partitioner()

template<int dim, int spacedim>
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 539 of file tria_base.cc.

◆ compute_vertices_with_ghost_neighbors()

template<int dim, int spacedim>
std::map< unsigned int, std::set<::types::subdomain_id > > parallel::TriangulationBase< dim, spacedim >::compute_vertices_with_ghost_neighbors ( ) const
virtualinherited

Return a map that, for each vertex, lists all the processors whose subdomains are adjacent to that vertex.

Deprecated:
Use GridTools::compute_vertices_with_ghost_neighbors() instead of parallel::TriangulationBase::compute_vertices_with_ghost_neighbors().

Definition at line 345 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 355 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 366 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 147 of file tria_base.cc.

◆ update_reference_cell_types()

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

Update the internal reference_cell_types vector.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 292 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 377 of file tria_base.cc.

◆ clear()

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

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

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

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

Definition at line 10093 of file tria.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 10109 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 10121 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 10379 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 10720 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 10732 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 10748 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 10761 of file tria.cc.

◆ prepare_coarsening_and_refinement() [1/4]

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

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

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

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

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

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

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

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

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

Definition at line 14010 of file tria.cc.

◆ prepare_coarsening_and_refinement() [2/4]

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

Definition at line 13767 of file tria.cc.

◆ prepare_coarsening_and_refinement() [3/4]

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

Definition at line 13786 of file tria.cc.

◆ prepare_coarsening_and_refinement() [4/4]

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

Definition at line 13805 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 10798 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 10781 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 10812 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 10823 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 10867 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 10850 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 10881 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 10895 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 10915 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 11085 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 11096 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 11114 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 11145 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 11163 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 10996 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 11215 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 11200 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 11229 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 11243 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 11036 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 11338 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 11319 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 11352 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 11366 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 11076 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 11407 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 11388 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 11421 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 11435 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 10958 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 11457 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 11489 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 11713 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 11744 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 11569 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 11583 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 11597 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 11615 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 11632 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 11650 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 11781 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 11795 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 11810 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 11828 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 11845 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 11863 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 11904 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 11924 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 11990 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 12030 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 12059 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 11944 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 11969 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 12132 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 12153 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 12174 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 12207 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 12224 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 12234 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 12730 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 12768 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 12778 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 12786 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 12943 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 12951 of file tria.cc.

◆ n_quads() [3/8]

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

Definition at line 12797 of file tria.cc.

◆ n_quads() [4/8]

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

Definition at line 12805 of file tria.cc.

◆ n_quads() [5/8]

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

Definition at line 12846 of file tria.cc.

◆ n_quads() [6/8]

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

Definition at line 12854 of file tria.cc.

◆ n_quads() [7/8]

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

Definition at line 12894 of file tria.cc.

◆ n_quads() [8/8]

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

Definition at line 12902 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 13010 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 13018 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 12829 of file tria.cc.

◆ n_active_quads() [4/8]

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

Definition at line 12837 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 12878 of file tria.cc.

◆ n_active_quads() [6/8]

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

Definition at line 12886 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 12926 of file tria.cc.

◆ n_active_quads() [8/8]

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

Definition at line 12934 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 13029 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 13038 of file tria.cc.

◆ n_hexs() [3/4]

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

Definition at line 13072 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 13081 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 13055 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 13064 of file tria.cc.

◆ n_active_hexs() [3/4]

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

Definition at line 13101 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 13110 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 12578 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 12679 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 12586 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 12699 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 12602 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 12640 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.)

◆ has_hanging_nodes()

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

Return true if the triangulation has hanging nodes.

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

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

Definition at line 12718 of file tria.cc.

◆ 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 13121 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 13130 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 13164 of file tria.cc.

◆ max_adjacent_cells() [2/4]

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

Definition at line 13139 of file tria.cc.

◆ max_adjacent_cells() [3/4]

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

Definition at line 13148 of file tria.cc.

◆ max_adjacent_cells() [4/4]

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

Definition at line 13156 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 13206 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 13215 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 12754 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 12739 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 12991 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 12813 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 12862 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 12910 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 12962 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 12972 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 12981 of file tria.cc.

◆ n_raw_quads() [9/9]

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

Definition at line 13001 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 13047 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 12821 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 12870 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 12918 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 13092 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 12659 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 12621 of file tria.cc.

◆ save()

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.

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.

◆ 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 13224 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 13244 of file tria.cc.

◆ get_reference_cell_types()

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

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

Definition at line 13415 of file tria.cc.

◆ all_reference_cell_types_are_hyper_cube()

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

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

Definition at line 13424 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.

◆ 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 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

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

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 204 of file subscriptor.cc.

Member Data Documentation

◆ settings

template<int dim, int spacedim = dim>
const Settings parallel::shared::Triangulation< dim, spacedim >::settings
private

Settings

Definition at line 346 of file shared_tria.h.

◆ allow_artificial_cells

template<int dim, int spacedim = dim>
const bool parallel::shared::Triangulation< dim, spacedim >::allow_artificial_cells
private

A flag to decide whether or not artificial cells are allowed.

Definition at line 351 of file shared_tria.h.

◆ true_subdomain_ids_of_cells

template<int dim, int spacedim = dim>
std::vector<types::subdomain_id> parallel::shared::Triangulation< dim, spacedim >::true_subdomain_ids_of_cells
private

A vector containing subdomain IDs of cells obtained by partitioning using either zorder, METIS, or a user-defined partitioning scheme. In case allow_artificial_cells is false, this vector is consistent with IDs stored in cell->subdomain_id() of the triangulation class. When allow_artificial_cells is true, cells which are artificial will have cell->subdomain_id() == numbers::artificial;

The original partition information is stored to allow using sequential DoF distribution and partitioning functions with semi-artificial cells.

Definition at line 372 of file shared_tria.h.

◆ true_level_subdomain_ids_of_cells

template<int dim, int spacedim = dim>
std::vector<std::vector<types::subdomain_id> > parallel::shared::Triangulation< dim, spacedim >::true_level_subdomain_ids_of_cells
private

A vector containing level subdomain IDs of cells obtained by partitioning each level.

The original partition information is stored to allow using sequential DoF distribution and partitioning functions with semi-artificial cells.

Definition at line 383 of file shared_tria.h.

◆ mpi_communicator

template<int dim, int spacedim = dim>
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 255 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 261 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 266 of file tria_base.h.

◆ number_cache

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

Definition at line 313 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 1544 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 1549 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 2272 of file tria.h.


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