Reference documentation for deal.II version Git 78150d9107 2021-02-24 10:44:07 -0700
\(\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 | Private Attributes | List of all members
PersistentTriangulation< dim, spacedim > Class Template Reference

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

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

Public Types

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

Public Member Functions

 PersistentTriangulation (const Triangulation< dim, spacedim > &coarse_grid)
 
 PersistentTriangulation (const PersistentTriangulation< dim, spacedim > &old_tria)
 
virtual ~PersistentTriangulation () override=default
 
virtual void execute_coarsening_and_refinement () override
 
void restore ()
 
void restore (const unsigned int step_no)
 
unsigned int n_refinement_steps () const
 
virtual void copy_triangulation (const Triangulation< dim, spacedim > &tria) 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 create_triangulation_compatibility (const std::vector< Point< spacedim >> &vertices, const std::vector< CellData< dim >> &cells, const SubCellData &subcelldata) override
 
virtual void write_flags (std::ostream &out) const
 
virtual void read_flags (std::istream &in)
 
virtual void clear_flags ()
 
virtual std::size_t memory_consumption () 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 std::vector< types::boundary_idget_boundary_ids () const
 
virtual std::vector< types::manifold_idget_manifold_ids () const
 
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 create_cell_iterator (const CellId &cell_id) const
 
Cell iterator functions returning ranges of iterators
IteratorRange< cell_iteratorcell_iterators () const
 
IteratorRange< active_cell_iteratoractive_cell_iterators () const
 
IteratorRange< cell_iteratorcell_iterators_on_level (const unsigned int level) const
 
IteratorRange< active_cell_iteratoractive_cell_iterators_on_level (const unsigned int level) const
 
Face iterator functions
face_iterator begin_face () const
 
active_face_iterator begin_active_face () const
 
face_iterator end_face () const
 
IteratorRange< active_face_iteratoractive_face_iterators () const
 
Vertex iterator functions
vertex_iterator begin_vertex () const
 
active_vertex_iterator begin_active_vertex () const
 
vertex_iterator end_vertex () const
 
Information about the triangulation
unsigned int n_lines () const
 
unsigned int n_lines (const unsigned int level) const
 
unsigned int n_active_lines () const
 
unsigned int n_active_lines (const unsigned int level) const
 
unsigned int n_quads () const
 
unsigned int n_quads (const unsigned int level) const
 
unsigned int n_active_quads () const
 
unsigned int n_active_quads (const unsigned int level) const
 
unsigned int n_hexs () const
 
unsigned int n_hexs (const unsigned int level) const
 
unsigned int n_active_hexs () const
 
unsigned int n_active_hexs (const unsigned int level) const
 
unsigned int n_cells () const
 
unsigned int n_cells (const unsigned int level) const
 
unsigned int n_active_cells () const
 
unsigned int n_active_cells (const unsigned int level) const
 
virtual types::global_cell_index n_global_active_cells () const
 
unsigned int n_faces () const
 
unsigned int n_active_faces () const
 
unsigned int n_levels () const
 
virtual unsigned int n_global_levels () const
 
virtual bool has_hanging_nodes () const
 
unsigned int n_vertices () const
 
const std::vector< Point< spacedim > > & get_vertices () const
 
unsigned int n_used_vertices () const
 
bool vertex_used (const unsigned int index) const
 
const std::vector< bool > & get_used_vertices () const
 
unsigned int max_adjacent_cells () const
 
virtual types::subdomain_id locally_owned_subdomain () const
 
Triangulation< dim, spacedim > & get_triangulation ()
 
const Triangulation< dim, spacedim > & get_triangulation () const
 
Internal information about the number of objects
unsigned int n_raw_lines () const
 
unsigned int n_raw_lines (const unsigned int level) const
 
unsigned int n_raw_quads () const
 
unsigned int n_raw_quads (const unsigned int level) const
 
unsigned int n_raw_hexs (const unsigned int level) const
 
unsigned int n_raw_cells (const unsigned int level) const
 
unsigned int n_raw_faces () const
 
template<class Archive >
void save (Archive &ar, const unsigned int version) const
 
template<class Archive >
void load (Archive &ar, const unsigned int version)
 
virtual void add_periodicity (const std::vector< GridTools::PeriodicFacePair< cell_iterator >> &)
 
const std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, std::bitset< 3 > > > & get_periodic_face_map () const
 
const std::vector< ReferenceCell > & get_reference_cells () const
 
bool all_reference_cells_are_hyper_cube () const
 
template<class Archive >
void serialize (Archive &archive, const unsigned int version)
 
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 ::ExceptionBaseExcTriaNotEmpty ()
 
static ::ExceptionBaseExcFlagsNotCleared ()
 
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 spacedimension = spacedim
 
static const unsigned int space_dimension = spacedim
 

Private Attributes

SmartPointer< const Triangulation< dim, spacedim >, PersistentTriangulation< dim, spacedim > > coarse_grid
 
std::vector< std::vector< bool > > refine_flags
 
std::vector< std::vector< bool > > coarsen_flags
 

Keeping up with what happens to a triangulation

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

Exceptions

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 ()
 
MeshSmoothing smooth_grid
 
std::vector< ReferenceCellreference_cells
 
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)
 
void update_periodic_face_map ()
 
virtual void update_reference_cells ()
 

Detailed Description

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

This class handles the history of a triangulation and can rebuild it after it was deleted some time before. Its main purpose is support for time- dependent problems where one frequently deletes a triangulation due to memory pressure and later wants to rebuild it; this class has all the information to rebuild it exactly as it was before including the mapping of cell numbers to the geometrical cells.

Basically, this is a drop-in replacement for the triangulation. Since it is derived from the Triangulation class, it shares all the functionality, but it overrides some virtual functions and adds some functions, too. The main change to the base class is that it overrides the execute_coarsening_and_refinement function, where the new version first stores all refinement and coarsening flags and only then calls the respective function of the base class. The stored flags may later be used to restore the grid just as it was before. Some other functions have been extended slightly as well, see their documentation for more information.

We note that since the triangulation is created in exactly the same state as it was before, other objects working on it should result in the same state as well. This holds in particular for the DoFHandler object, which will assign the same degrees of freedom to the original cells and the ones after reconstruction of the triangulation. You can therefore safely use data vectors computed on the original grid on the reconstructed grid as well.

Usage

You can use objects of this class almost in the same way as objects of the Triangulation class. One of the few differences is that you can only construct such an object by giving a coarse grid to the constructor. The coarse grid will be used to base the triangulation on, and therefore the lifetime of the coarse grid has to be longer than the lifetime of the object of this class.

Basically, usage looks like this:

... // initialize coarse grid
PersistentTriangulation<dim> grid (coarse_grid);
for (...)
{
// restore grid from coarse grid
// and stored refinement flags
grid.restore ();
... // do something with the grid
... // flag some cells for refinement
// or coarsening
grid.execute_coarsening_and_refinement ();
// actually refine grid and store
// the flags
... // so something more with the grid
grid.clear (); // delete the grid, but keep the
// refinement flags for later use
// in grid.restore() above
... // do something where the grid
// is not needed anymore, e.g.
// working with another grid
};

Note that initially, the PersistentTriangulation object does not constitute a triangulation; it only becomes one after restore is first called. Note also that the execute_coarsening_and_refinement stores all necessary flags for later reconstruction using the restore function. Triangulation::clear() resets the underlying triangulation to a virgin state, but does not affect the stored refinement flags needed for later reconstruction and does also not touch the coarse grid which is used within restore().

Definition at line 108 of file persistent_tria.h.

Member Typedef Documentation

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

Member Enumeration Documentation

◆ MeshSmoothing

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

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

Enumerator
none 

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

limit_level_difference_at_vertices 

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

limit_level_difference_at_vertices.png

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

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

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

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

eliminate_unrefined_islands 

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

patch_level_1 

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

coarsest_level_1 

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

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

allow_anisotropic_smoothing 

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

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

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

eliminate_refined_inner_islands 

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

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

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

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

eliminate_refined_boundary_islands 

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

do_not_produce_unrefined_islands 

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

smoothing_on_refinement 

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

smoothing_on_coarsening 

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

maximum_smoothing 

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

Definition at line 1142 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 2023 of file tria.h.

Constructor & Destructor Documentation

◆ PersistentTriangulation() [1/2]

template<int dim, int spacedim>
PersistentTriangulation< dim, spacedim >::PersistentTriangulation ( const Triangulation< dim, spacedim > &  coarse_grid)

Build up the triangulation from the coarse grid in future. Copy smoothing flags, etc from that grid as well. Note that the initial state of the triangulation is empty, until restore_grid is called for the first time.

The coarse grid must persist until the end of this object, since it will be used upon reconstruction of the grid.

Definition at line 34 of file persistent_tria.cc.

◆ PersistentTriangulation() [2/2]

template<int dim, int spacedim>
PersistentTriangulation< dim, spacedim >::PersistentTriangulation ( const PersistentTriangulation< dim, spacedim > &  old_tria)

Copy constructor. This operation is only allowed, if the triangulation underlying the object to be copied is presently empty. Refinement flags as well as the pointer to the coarse grid are copied, however.

Definition at line 42 of file persistent_tria.cc.

◆ ~PersistentTriangulation()

template<int dim, int spacedim = dim>
virtual PersistentTriangulation< dim, spacedim >::~PersistentTriangulation ( )
overridevirtualdefault

Destructor.

Member Function Documentation

◆ execute_coarsening_and_refinement()

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

Overloaded version of the same function in the base class which stores the refinement and coarsening flags for later reconstruction of the triangulation and after that calls the respective function of the base class.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 59 of file persistent_tria.cc.

◆ restore() [1/2]

template<int dim, int spacedim>
void PersistentTriangulation< dim, spacedim >::restore ( )

Restore the grid according to the saved data. For this, the coarse grid is copied and the grid is stepwise rebuilt using the saved flags.

Note that this function will result in an error if the underlying triangulation is not empty, i.e. it will only succeed if this object is newly created or the clear() function of the base class was called on it before.

Repeatedly calls the restore(unsigned int) function in a loop over all refinement steps.

Definition at line 78 of file persistent_tria.cc.

◆ restore() [2/2]

template<int dim, int spacedim>
void PersistentTriangulation< dim, spacedim >::restore ( const unsigned int  step_no)

Differential restore. Performs the step_noth local refinement and coarsening step. Step 0 stands for the copying of the coarse grid.

This function will only succeed if the triangulation is in just the state it were if restore would have been called from step=0...step_no-1 before.

Definition at line 90 of file persistent_tria.cc.

◆ n_refinement_steps()

template<int dim, int spacedim>
unsigned int PersistentTriangulation< dim, spacedim >::n_refinement_steps ( ) const

Return the number of refinement and coarsening steps. This is given by the size of the refine_flags vector.

Definition at line 116 of file persistent_tria.cc.

◆ copy_triangulation()

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

Overload this function to use tria as a new coarse grid. The present triangulation and all refinement and coarsening flags storing its history are deleted, and the state of the underlying triangulation is reset to be empty, until restore_grid is called the next time.

The coarse grid must persist until the end of this object, since it will be used upon reconstruction of the grid.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 125 of file persistent_tria.cc.

◆ create_triangulation() [1/2]

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

Throw an error, since this function is not useful in the context of this class.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 138 of file persistent_tria.cc.

◆ create_triangulation() [2/2]

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

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

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

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

For conditions when this function can generate a valid triangulation, see the documentation of this class, and the GridIn and GridReordering class.

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

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

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 150 of file persistent_tria.cc.

◆ create_triangulation_compatibility()

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

An overload of the respective function of the base class.

Throw an error, since this function is not useful in the context of this class.

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 162 of file persistent_tria.cc.

◆ write_flags()

template<int dim, int spacedim>
void PersistentTriangulation< dim, spacedim >::write_flags ( std::ostream &  out) const
virtual

Write all refine and coarsen flags to the ostream out.

Definition at line 174 of file persistent_tria.cc.

◆ read_flags()

template<int dim, int spacedim>
void PersistentTriangulation< dim, spacedim >::read_flags ( std::istream &  in)
virtual

Reads all refine and coarsen flags that previously were written by write_flags(...). This is especially useful for rebuilding the triangulation after the end or breakdown of a program and its restart.

Definition at line 203 of file persistent_tria.cc.

◆ clear_flags()

template<int dim, int spacedim>
void PersistentTriangulation< dim, spacedim >::clear_flags ( )
virtual

Clear all flags. Retains the same coarse grid.

Definition at line 241 of file persistent_tria.cc.

◆ memory_consumption()

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

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

Reimplemented from Triangulation< dim, spacedim >.

Definition at line 251 of file persistent_tria.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::DistributedTriangulationBase< dim, spacedim >, parallel::DistributedTriangulationBase< 1, spacedim >, parallel::DistributedTriangulationBase< dim, dim >, parallel::distributed::Triangulation< dim, spacedim >, and parallel::distributed::Triangulation< dim >.

Definition at line 10180 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 10196 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 10208 of file tria.cc.

◆ get_boundary_ids()

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

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

See also
Glossary entry on boundary indicators

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

Definition at line 10345 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 10806 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 10818 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 10834 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 10847 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 14124 of file tria.cc.

◆ prepare_coarsening_and_refinement() [2/4]

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

Definition at line 13881 of file tria.cc.

◆ prepare_coarsening_and_refinement() [3/4]

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

Definition at line 13900 of file tria.cc.

◆ prepare_coarsening_and_refinement() [4/4]

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

Definition at line 13919 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 10884 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 10867 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 10898 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 10909 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 10953 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 10936 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 10967 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 10981 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 11001 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 11171 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 11182 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 11200 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 11231 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 11249 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 11082 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 11301 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 11286 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 11315 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 11329 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 11122 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 11424 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 11405 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 11438 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 11452 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 11162 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 11493 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 11474 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 11507 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 11521 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 11044 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 11543 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 11575 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 11799 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 11830 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 11655 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 11669 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 11683 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 11701 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 11718 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 11736 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 11867 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 11881 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 11896 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 11914 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 11931 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 11949 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 11990 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 12010 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 12101 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 12141 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 12170 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 12030 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 12055 of file tria.cc.

◆ create_cell_iterator()

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

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

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

Definition at line 12076 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 12243 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 12264 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 12285 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 12318 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 12335 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 12345 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 12841 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 12879 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 12889 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 12897 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 13054 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 13062 of file tria.cc.

◆ n_quads() [3/8]

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

Definition at line 12908 of file tria.cc.

◆ n_quads() [4/8]

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

Definition at line 12916 of file tria.cc.

◆ n_quads() [5/8]

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

Definition at line 12957 of file tria.cc.

◆ n_quads() [6/8]

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

Definition at line 12965 of file tria.cc.

◆ n_quads() [7/8]

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

Definition at line 13005 of file tria.cc.

◆ n_quads() [8/8]

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

Definition at line 13013 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 13121 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 13129 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 12940 of file tria.cc.

◆ n_active_quads() [4/8]

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

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

◆ n_active_quads() [6/8]

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

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

◆ n_active_quads() [8/8]

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

Definition at line 13045 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 13140 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 13149 of file tria.cc.

◆ n_hexs() [3/4]

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

Definition at line 13183 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 13192 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 13166 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 13175 of file tria.cc.

◆ n_active_hexs() [3/4]

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

Definition at line 13212 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 13221 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 12689 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 12790 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 12697 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 12810 of file tria.cc.

◆ n_global_active_cells()

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

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

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

Definition at line 12704 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 12713 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 12751 of file tria.cc.

◆ n_levels()

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

Return the number of levels in this triangulation.

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

◆ n_global_levels()

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

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

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

◆ 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 12829 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 13232 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 13241 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 13275 of file tria.cc.

◆ max_adjacent_cells() [2/4]

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

Definition at line 13250 of file tria.cc.

◆ max_adjacent_cells() [3/4]

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

Definition at line 13259 of file tria.cc.

◆ max_adjacent_cells() [4/4]

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

Definition at line 13267 of file tria.cc.

◆ locally_owned_subdomain()

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

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

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

Definition at line 13308 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 13317 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 13326 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 12865 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 12850 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 13102 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 12924 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 12973 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 13021 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 13073 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 13083 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 13092 of file tria.cc.

◆ n_raw_quads() [9/9]

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

Definition at line 13112 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 13158 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 12932 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 12981 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 13029 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 13203 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 12770 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 12732 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 using the BOOST serialization library.

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

◆ load()

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

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

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

◆ 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 13335 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 13355 of file tria.cc.

◆ get_reference_cells()

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

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

Definition at line 13526 of file tria.cc.

◆ all_reference_cells_are_hyper_cube()

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

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

Definition at line 13535 of file tria.cc.

◆ serialize()

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

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

◆ subscribe()

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

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

Definition at line 136 of file subscriptor.cc.

◆ unsubscribe()

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

Unsubscribes a user from the object.

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

Definition at line 156 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
inlineinherited

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

Definition at line 301 of file subscriptor.h.

◆ list_subscribers() [1/2]

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

List the subscribers to the input stream.

Definition at line 318 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 204 of file subscriptor.cc.

Member Data Documentation

◆ dimension

template<int dim, int spacedim = dim>
const unsigned int PersistentTriangulation< dim, spacedim >::dimension = dim
static

Make the dimension available in function templates.

Definition at line 114 of file persistent_tria.h.

◆ spacedimension

template<int dim, int spacedim = dim>
const unsigned int PersistentTriangulation< dim, spacedim >::spacedimension = spacedim
static

Definition at line 115 of file persistent_tria.h.

◆ coarse_grid

template<int dim, int spacedim = dim>
SmartPointer<const Triangulation<dim, spacedim>, PersistentTriangulation<dim, spacedim> > PersistentTriangulation< dim, spacedim >::coarse_grid
private

This grid shall be used as coarse grid.

Definition at line 268 of file persistent_tria.h.

◆ refine_flags

template<int dim, int spacedim = dim>
std::vector<std::vector<bool> > PersistentTriangulation< dim, spacedim >::refine_flags
private

Vectors holding the refinement and coarsening flags of the different sweeps on this time level. The vectors therefore hold the history of the grid.

Definition at line 275 of file persistent_tria.h.

◆ coarsen_flags

template<int dim, int spacedim = dim>
std::vector<std::vector<bool> > PersistentTriangulation< dim, spacedim >::coarsen_flags
private

refine_flags

Definition at line 280 of file persistent_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 1553 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 2278 of file tria.h.


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