Reference documentation for deal.II version 8.5.1
|
#include <deal.II/distributed/shared_tria.h>
Public Member Functions | |
Triangulation (MPI_Comm mpi_communicator, const typename ::Triangulation< dim, spacedim >::MeshSmoothing=(::Triangulation< dim, spacedim >::none), const bool allow_artificial_cells=false) | |
virtual | ~Triangulation () |
virtual void | execute_coarsening_and_refinement () |
virtual void | create_triangulation (const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata) |
virtual void | copy_triangulation (const ::Triangulation< dim, spacedim > &other_tria) |
template<class Archive > | |
void | load (Archive &ar, const unsigned int version) |
const std::vector< types::subdomain_id > & | get_true_subdomain_ids_of_cells () const |
bool | with_artificial_cells () const |
Public Member Functions inherited from parallel::Triangulation< dim, spacedim > | |
Triangulation (MPI_Comm mpi_communicator, const typename ::Triangulation< dim, spacedim >::MeshSmoothing smooth_grid=(::Triangulation< dim, spacedim >::none), const bool check_for_distorted_cells=false) | |
virtual MPI_Comm | get_communicator () const |
virtual void | copy_triangulation (const ::Triangulation< dim, spacedim > &old_tria) |
const std::vector< unsigned int > & | n_locally_owned_active_cells_per_processor () const |
unsigned int | n_locally_owned_active_cells () const |
virtual types::global_dof_index | n_global_active_cells () const |
virtual std::size_t | memory_consumption () const |
virtual unsigned int | n_global_levels () const |
types::subdomain_id | locally_owned_subdomain () const |
const std::set< types::subdomain_id > & | ghost_owners () const |
const std::set< types::subdomain_id > & | level_ghost_owners () const |
Public Member Functions inherited from Triangulation< dim, spacedim > | |
Triangulation (const MeshSmoothing smooth_grid=none, const bool check_for_distorted_cells=false) | |
Triangulation (const Triangulation< dim, spacedim > &t) | |
Triangulation (Triangulation< dim, spacedim > &&tria) | |
Triangulation & | operator= (Triangulation< dim, spacedim > &&tria) |
virtual void | clear () |
virtual void | set_mesh_smoothing (const MeshSmoothing mesh_smoothing) |
virtual const MeshSmoothing & | get_mesh_smoothing () const |
void | set_boundary (const types::manifold_id number, const Boundary< dim, spacedim > &boundary_object) |
void | set_boundary (const types::manifold_id number) |
void | set_manifold (const types::manifold_id number, const Manifold< dim, spacedim > &manifold_object) |
void | set_manifold (const types::manifold_id number) |
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 Boundary< dim, spacedim > & | get_boundary (const types::manifold_id number) const |
const Manifold< dim, spacedim > & | get_manifold (const types::manifold_id number) const |
std::vector< types::boundary_id > | get_boundary_ids () const |
std::vector< types::manifold_id > | get_manifold_ids () const |
virtual void | copy_triangulation (const Triangulation< dim, spacedim > &other_tria) |
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 () |
void | set_all_refine_flags () |
void | refine_global (const unsigned int times=1) |
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_iterator > | cell_iterators () const |
IteratorRange< active_cell_iterator > | active_cell_iterators () const |
IteratorRange< cell_iterator > | cell_iterators_on_level (const unsigned int level) const |
IteratorRange< active_cell_iterator > | active_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 |
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 |
unsigned int | n_faces () const |
unsigned int | n_active_faces () const |
unsigned int | n_levels () const |
virtual bool | has_hanging_nodes () const |
unsigned int | n_vertices () const |
const std::vector< Point< spacedim > > & | get_vertices () const |
unsigned int | n_used_vertices () const |
bool | vertex_used (const unsigned int index) const |
const std::vector< bool > & | get_used_vertices () const |
unsigned int | max_adjacent_cells () const |
Triangulation< dim, spacedim > & | get_triangulation () |
const Triangulation< dim, spacedim > & | get_triangulation () const |
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 | 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 |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Member Functions | |
void | partition () |
Private Attributes | |
const bool | allow_artificial_cells |
std::vector< types::subdomain_id > | true_subdomain_ids_of_cells |
Additional Inherited Members | |
Public Types inherited from Triangulation< dim, spacedim > | |
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 } |
typedef TriaIterator< CellAccessor< dim, spacedim > > | cell_iterator |
typedef TriaActiveIterator< CellAccessor< dim, spacedim > > | active_cell_iterator |
enum | CellStatus { CELL_PERSIST, CELL_REFINE, CELL_COARSEN, CELL_INVALID } |
Static Public Member Functions inherited from Triangulation< dim, spacedim > | |
static ::ExceptionBase & | ExcInvalidLevel (int arg1) |
static ::ExceptionBase & | ExcTriangulationNotEmpty (int arg1, int arg2) |
static ::ExceptionBase & | ExcGridReadError () |
static ::ExceptionBase & | ExcFacesHaveNoLevel () |
static ::ExceptionBase & | ExcEmptyLevel (int arg1) |
static ::ExceptionBase & | ExcNonOrientableTriangulation () |
static ::ExceptionBase & | ExcBoundaryIdNotFound (types::boundary_id arg1) |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, char *arg2, std::string &arg3) |
static ::ExceptionBase & | ExcNoSubscriber (char *arg1, char *arg2) |
Public Attributes inherited from Triangulation< dim, spacedim > | |
Signals | signals |
Static Public Attributes inherited from Triangulation< dim, spacedim > | |
static const StraightBoundary< dim, spacedim > | straight_boundary = StraightBoundary<dim,spacedim>() |
static const unsigned int | dimension = dim |
static const unsigned int | space_dimension = spacedim |
Protected Member Functions inherited from parallel::Triangulation< dim, spacedim > | |
virtual void | update_number_cache () |
Protected Member Functions inherited from Triangulation< dim, spacedim > | |
void | update_periodic_face_map () |
Static Protected Member Functions inherited from Triangulation< dim, spacedim > | |
static void | write_bool_vector (const unsigned int magic_number1, const std::vector< bool > &v, const unsigned int magic_number2, std::ostream &out) |
static void | read_bool_vector (const unsigned int magic_number1, std::vector< bool > &v, const unsigned int magic_number2, std::istream &in) |
Protected Attributes inherited from parallel::Triangulation< dim, spacedim > | |
MPI_Comm | mpi_communicator |
types::subdomain_id | my_subdomain |
types::subdomain_id | n_subdomains |
Protected Attributes inherited from Triangulation< dim, spacedim > | |
MeshSmoothing | smooth_grid |
This is an extension of Triangulation class to automatically partition triangulation when run with MPI. Different from the parallel::distributed::Triangulation, the entire mesh is stored on each processor. However, cells are labeled according to the id of the processor which "owns" them. The partitioning is done automatically inside the DoFHandler by calling Metis. This enables distributing DoFs among processors and therefore splitting matrices and vectors across processors. The usage of this class is demonstrated in step-18.
Definition at line 68 of file shared_tria.h.
Triangulation< dim, spacedim >::Triangulation | ( | MPI_Comm | mpi_communicator, |
const typename ::Triangulation< dim, spacedim >::MeshSmoothing | smooth_grid = (::Triangulation<dim,spacedim>::none) , |
||
const bool | allow_artificial_cells = false |
||
) |
Constructor.
If allow_aritifical_cells
is true, this class will behave similar to parallel::distributed::Triangulation in that there will be locally owned, ghost and artificial cells.
Otherwise all non-locally owned cells are considered ghost.
Definition at line 37 of file shared_tria.cc.
|
virtual |
Destructor.
Reimplemented from parallel::Triangulation< dim, spacedim >.
Definition at line 110 of file shared_tria.cc.
|
virtual |
Coarsen and refine the mesh according to refinement and coarsening flags set.
This step is equivalent to the Triangulation class with an addition of calling GridTools::partition_triangulation() at the end.
Reimplemented from Triangulation< dim, spacedim >.
Definition at line 117 of file shared_tria.cc.
|
virtual |
Create a triangulation.
This function also partitions triangulation based on the MPI communicator provided to the constructor.
Reimplemented from Triangulation< dim, spacedim >.
Definition at line 128 of file shared_tria.cc.
|
virtual |
Copy other_tria
to this triangulation.
This function also partitions triangulation based on the MPI communicator provided to the constructor.
Definition at line 152 of file shared_tria.cc.
void Triangulation< dim, spacedim >::load | ( | Archive & | ar, |
const unsigned int | version | ||
) |
Read the data of this object from a stream for the purpose of serialization. Throw away the previous content.
This function first does the same work as in Triangulation::load, then partitions the triangulation based on the MPI communicator provided to the constructor.
Definition at line 181 of file shared_tria.h.
const std::vector< types::subdomain_id > & Triangulation< dim, spacedim >::get_true_subdomain_ids_of_cells | ( | ) | const |
Return a vector of length Triangulation::n_active_cells() where each element stores the subdomain id of the owner of this cell. The elements of the vector are obviously the same as the subdomain ids for locally owned and ghost cells, but are also correct for artificial cells that do not store who the owner of the cell is in their subdomain_id field.
Definition at line 102 of file shared_tria.cc.
bool Triangulation< dim, spacedim >::with_artificial_cells | ( | ) | const |
Return allow_artificial_cells , namely true if artificial cells are allowed.
Definition at line 93 of file shared_tria.cc.
|
private |
This function calls GridTools::partition_triangulation () and if requested in the constructor of the class marks artificial cells.
Definition at line 47 of file shared_tria.cc.
|
private |
A flag to decide whether or not artificial cells are allowed.
Definition at line 156 of file shared_tria.h.
|
private |
A vector containing subdomain IDs of cells obtained by partitioning using METIS. In case allow_artificial_cells is false, this vector is consistent with IDs stored in cell->subdomain_id() of the triangulation class. When allow_artificial_cells is true, cells which are artificial will have cell->subdomain_id() == numbers::artificial;
The original parition information is stored to allow using sequential DoF distribution and partitioning functions with semi-artificial cells.
Definition at line 175 of file shared_tria.h.