Reference documentation for deal.II version 9.2.0
\(\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\}}\)
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
parallel::TriangulationBase< dim, spacedim > Class Template Referenceabstract

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

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

Classes

struct  NumberCache
 

Public Member Functions

 TriangulationBase (MPI_Comm mpi_communicator, const typename ::Triangulation< dim, spacedim >::MeshSmoothing smooth_grid=(::Triangulation< dim, spacedim >::none), const bool check_for_distorted_cells=false)
 
virtual ~TriangulationBase () override
 
virtual const MPI_Commget_communicator () const
 
virtual bool is_multilevel_hierarchy_constructed () const =0
 
virtual void copy_triangulation (const ::Triangulation< dim, spacedim > &old_tria) override
 
unsigned int n_locally_owned_active_cells () const
 
virtual types::global_cell_index n_global_active_cells () const override
 
virtual std::size_t memory_consumption () const override
 
virtual unsigned int n_global_levels () const override
 
types::subdomain_id locally_owned_subdomain () const override
 
const std::set< types::subdomain_id > & ghost_owners () const
 
const std::set< types::subdomain_id > & level_ghost_owners () const
 
virtual std::map< unsigned int, std::set<::types::subdomain_id > > compute_vertices_with_ghost_neighbors () const
 
virtual std::vector< types::boundary_idget_boundary_ids () const override
 
virtual std::vector< types::manifold_idget_manifold_ids () const override
 
- Public Member Functions inherited from Triangulation< dim, dim >
 Triangulation (const MeshSmoothing smooth_grid=none, const bool check_for_distorted_cells=false)
 
 Triangulation (const Triangulation< dim, spacedim > &)=delete
 
 Triangulation (Triangulation< dim, spacedim > &&tria) noexcept
 
Triangulationoperator= (Triangulation< dim, spacedim > &&tria) noexcept
 
virtual ~Triangulation () override
 
virtual void clear ()
 
virtual void set_mesh_smoothing (const MeshSmoothing mesh_smoothing)
 
virtual const MeshSmoothingget_mesh_smoothing () const
 
void set_manifold (const types::manifold_id number, const Manifold< dim, spacedim > &manifold_object)
 
void reset_manifold (const types::manifold_id manifold_number)
 
void reset_all_manifolds ()
 
void set_all_manifold_ids (const types::manifold_id number)
 
void set_all_manifold_ids_on_boundary (const types::manifold_id number)
 
void set_all_manifold_ids_on_boundary (const types::boundary_id b_id, const types::manifold_id number)
 
const Manifold< dim, spacedim > & get_manifold (const types::manifold_id number) const
 
virtual std::vector< types::boundary_idget_boundary_ids () const
 
virtual std::vector< types::manifold_idget_manifold_ids () const
 
virtual void copy_triangulation (const Triangulation< dim, spacedim > &other_tria)
 
virtual void create_triangulation (const std::vector< Point< spacedim >> &vertices, const std::vector< CellData< dim >> &cells, const SubCellData &subcelldata)
 
virtual void create_triangulation (const TriangulationDescription::Description< dim, spacedim > &construction_data)
 
virtual void create_triangulation_compatibility (const std::vector< Point< spacedim >> &vertices, const std::vector< CellData< dim >> &cells, const SubCellData &subcelldata)
 
void flip_all_direction_flags ()
 
bool prepare_coarsening_and_refinement ()
 
bool prepare_coarsening_and_refinement ()
 
bool prepare_coarsening_and_refinement ()
 
unsigned int n_quads () const
 
unsigned int n_quads (const unsigned int) const
 
unsigned int n_quads () const
 
unsigned int n_quads (const unsigned int) const
 
unsigned int n_quads () const
 
unsigned int n_quads (const unsigned int) const
 
unsigned int n_active_quads (const unsigned int) const
 
unsigned int n_active_quads () const
 
unsigned int n_active_quads (const unsigned int) const
 
unsigned int n_active_quads () const
 
unsigned int n_active_quads (const unsigned int) const
 
unsigned int n_active_quads () 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 max_adjacent_cells () const
 
unsigned int max_adjacent_cells () const
 
unsigned int max_adjacent_cells () const
 
unsigned int n_raw_lines (const unsigned int level) const
 
unsigned int n_raw_lines () const
 
unsigned int n_raw_lines (const unsigned int level) const
 
unsigned int n_raw_lines () const
 
unsigned int n_raw_lines (const unsigned int level) const
 
unsigned int n_raw_lines () const
 
unsigned int n_raw_quads (const unsigned int) const
 
unsigned int n_raw_quads (const unsigned int) const
 
unsigned int n_raw_quads (const unsigned int) const
 
unsigned int n_raw_quads (const unsigned int level) const
 
unsigned int n_raw_quads (const unsigned int level) const
 
unsigned int n_raw_quads (const unsigned int) const
 
unsigned int n_raw_quads () const
 
unsigned int n_raw_hexs (const unsigned int) const
 
unsigned int n_raw_hexs (const unsigned int) const
 
unsigned int n_raw_hexs (const unsigned int) const
 
unsigned int n_raw_hexs (const unsigned int level) const
 
void set_all_refine_flags ()
 
void refine_global (const unsigned int times=1)
 
virtual void execute_coarsening_and_refinement ()
 
virtual bool prepare_coarsening_and_refinement ()
 
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
 
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 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
 
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 begin_face () const
 
active_face_iterator begin_active_face () const
 
face_iterator end_face () const
 
IteratorRange< active_face_iteratoractive_face_iterators () const
 
vertex_iterator begin_vertex () const
 
active_vertex_iterator begin_active_vertex () const
 
vertex_iterator end_vertex () const
 
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
 
unsigned int n_raw_lines () const
 
unsigned int n_raw_lines (const unsigned int level) const
 
unsigned int n_raw_quads () const
 
unsigned int n_raw_quads (const unsigned int level) const
 
unsigned int n_raw_hexs (const unsigned int level) const
 
unsigned int n_raw_cells (const unsigned int level) const
 
unsigned int n_raw_faces () const
 
virtual std::size_t memory_consumption () const
 
void save (Archive &ar, const unsigned int version) const
 
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
 
void serialize (Archive &archive, const unsigned int version)
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
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
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Protected Member Functions

virtual void update_number_cache ()
 
- Protected Member Functions inherited from Triangulation< dim, dim >
void update_periodic_face_map ()
 

Protected Attributes

MPI_Comm mpi_communicator
 
types::subdomain_id my_subdomain
 
types::subdomain_id n_subdomains
 
NumberCache number_cache
 
- Protected Attributes inherited from Triangulation< dim, dim >
MeshSmoothing smooth_grid
 

Additional Inherited Members

- Public Types inherited from Triangulation< dim, dim >
enum  MeshSmoothing
 
using cell_iterator = TriaIterator< CellAccessor< dim, spacedim > >
 
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
 
enum  CellStatus
 
- Static Public Member Functions inherited from Triangulation< dim, dim >
static ::ExceptionBaseExcInvalidLevel (int arg1, int arg2)
 
static ::ExceptionBaseExcTriangulationNotEmpty (int arg1, int arg2)
 
static ::ExceptionBaseExcGridReadError ()
 
static ::ExceptionBaseExcFacesHaveNoLevel ()
 
static ::ExceptionBaseExcEmptyLevel (int arg1)
 
static ::ExceptionBaseExcNonOrientableTriangulation ()
 
static ::ExceptionBaseExcBoundaryIdNotFound (types::boundary_id arg1)
 
static ::ExceptionBaseExcInconsistentCoarseningFlags ()
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
- Public Attributes inherited from Triangulation< dim, dim >
Signals signals
 
- Static Public Attributes inherited from Triangulation< dim, dim >
static const unsigned int dimension
 
static const unsigned int space_dimension
 
- Static Protected Member Functions inherited from Triangulation< dim, dim >
static void write_bool_vector (const unsigned int magic_number1, const std::vector< bool > &v, const unsigned int magic_number2, std::ostream &out)
 
static void read_bool_vector (const unsigned int magic_number1, std::vector< bool > &v, const unsigned int magic_number2, std::istream &in)
 

Detailed Description

template<int dim, int spacedim = dim>
class parallel::TriangulationBase< dim, spacedim >

This class describes the interface for all triangulation classes that work in parallel, namely parallel::distributed::Triangulation, parallel::fullydistributed::Triangulation, and parallel::shared::Triangulation.

It is, consequently, a class that can be used to test whether a pointer of reference to a triangulation object refers to a sequential triangulation, or whether the triangulation is in fact parallel. In other words, one could write a function like this:

template <int dim, int spacedim>
bool is_parallel (const ::Triangulation<dim,spacedim> &tria)
{
(&tria)
!= nullptr)
return true;
else
return false;
}

All parallel triangulations share certain traits, such as the fact that they communicate via MPI communicators or that they have locally owned, ghost, and possibly artificial cells. This class provides a number of member functions that allows querying some information about the triangulation that is independent of how exactly a parallel triangulation is implemented (i.e., which of the various classes derived from the current one it actually is).

Definition at line 78 of file tria_base.h.

Constructor & Destructor Documentation

◆ TriangulationBase()

template<int dim, int spacedim>
parallel::TriangulationBase< dim, spacedim >::TriangulationBase ( MPI_Comm  mpi_communicator,
const typename ::Triangulation< dim, spacedim >::MeshSmoothing  smooth_grid = (::Triangulation<dim, spacedim>::none),
const bool  check_for_distorted_cells = false 
)

Constructor.

Definition at line 43 of file tria_base.cc.

◆ ~TriangulationBase()

template<int dim, int spacedim>
parallel::TriangulationBase< dim, spacedim >::~TriangulationBase
overridevirtual

Destructor.

Definition at line 101 of file tria_base.cc.

Member Function Documentation

◆ get_communicator()

template<int dim, int spacedim>
const MPI_Comm & parallel::TriangulationBase< dim, spacedim >::get_communicator
virtual

Return MPI communicator used by this triangulation.

Definition at line 138 of file tria_base.cc.

◆ is_multilevel_hierarchy_constructed()

template<int dim, int spacedim = dim>
virtual bool parallel::TriangulationBase< dim, spacedim >::is_multilevel_hierarchy_constructed ( ) const
pure virtual

◆ copy_triangulation()

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

Implementation of the same function as in the base class.

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

Definition at line 63 of file tria_base.cc.

◆ n_locally_owned_active_cells()

template<int dim, int spacedim>
unsigned int parallel::TriangulationBase< dim, spacedim >::n_locally_owned_active_cells

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

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

Definition at line 117 of file tria_base.cc.

◆ n_global_active_cells()

template<int dim, int spacedim>
types::global_cell_index parallel::TriangulationBase< dim, spacedim >::n_global_active_cells
overridevirtual

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

Definition at line 131 of file tria_base.cc.

◆ memory_consumption()

template<int dim, int spacedim>
std::size_t parallel::TriangulationBase< dim, spacedim >::memory_consumption
overridevirtual

◆ n_global_levels()

template<int dim, int spacedim>
unsigned int parallel::TriangulationBase< dim, spacedim >::n_global_levels
overridevirtual

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

Definition at line 124 of file tria_base.cc.

◆ locally_owned_subdomain()

template<int dim, int spacedim>
types::subdomain_id parallel::TriangulationBase< dim, spacedim >::locally_owned_subdomain
override

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

Definition at line 288 of file tria_base.cc.

◆ ghost_owners()

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

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

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

Definition at line 297 of file tria_base.cc.

◆ level_ghost_owners()

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

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

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

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

Definition at line 306 of file tria_base.cc.

◆ compute_vertices_with_ghost_neighbors()

template<int dim, int spacedim>
std::map< unsigned int, std::set<::types::subdomain_id > > parallel::TriangulationBase< dim, spacedim >::compute_vertices_with_ghost_neighbors
virtual

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

Definition at line 315 of file tria_base.cc.

◆ get_boundary_ids()

template<int dim, int spacedim>
std::vector< types::boundary_id > parallel::TriangulationBase< dim, spacedim >::get_boundary_ids
overridevirtual

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

See also
Glossary entry on boundary indicators
Note
This function involves a global communication gathering all current IDs from all processes.

Definition at line 325 of file tria_base.cc.

◆ get_manifold_ids()

template<int dim, int spacedim>
std::vector< types::manifold_id > parallel::TriangulationBase< dim, spacedim >::get_manifold_ids
overridevirtual

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

See also
Glossary entry on manifold indicators
Note
This function involves a global communication gathering all current IDs from all processes.

Definition at line 336 of file tria_base.cc.

◆ update_number_cache()

template<int dim, int spacedim>
void parallel::TriangulationBase< dim, spacedim >::update_number_cache
protectedvirtual

Update the number_cache variable after mesh creation or refinement.

Definition at line 146 of file tria_base.cc.

Member Data Documentation

◆ mpi_communicator

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

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

Definition at line 240 of file tria_base.h.

◆ my_subdomain

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

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

Definition at line 246 of file tria_base.h.

◆ n_subdomains

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

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

Definition at line 251 of file tria_base.h.

◆ number_cache

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

Definition at line 288 of file tria_base.h.


The documentation for this class was generated from the following files:
parallel::TriangulationBase
Definition: tria_base.h:78