![]() |
Reference documentation for deal.II version Git 71165045b9 2021-01-23 20:21:06 -0500
|
#include <deal.II/distributed/tria.h>
Public Member Functions | |
Triangulation (const MPI_Comm &mpi_communicator, const typename ::Triangulation< dim, spacedim >::MeshSmoothing smooth_grid=(::Triangulation< dim, spacedim >::none), const Settings settings=default_setting) | |
virtual | ~Triangulation () override |
virtual void | clear () override |
bool | is_multilevel_hierarchy_constructed () const override |
void | execute_transfer (const typename ::internal::p4est::types< dim >::forest *parallel_forest, const typename ::internal::p4est::types< dim >::gloidx *previous_global_first_quadrant) |
virtual void | copy_triangulation (const ::Triangulation< dim, spacedim > &other_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 | execute_coarsening_and_refinement () override |
virtual bool | prepare_coarsening_and_refinement () override |
void | repartition () |
virtual bool | has_hanging_nodes () const override |
virtual std::size_t | memory_consumption () const override |
virtual std::size_t | memory_consumption_p4est () const |
void | write_mesh_vtk (const std::string &file_basename) const |
unsigned int | get_checksum () const |
virtual void | save (const std::string &filename) const override |
virtual void | load (const std::string &filename, const bool autopartition=true) override |
const std::vector< types::global_dof_index > & | get_p4est_tree_to_coarse_cell_permutation () const |
const std::vector< types::global_dof_index > & | get_coarse_cell_to_p4est_tree_permutation () const |
const ::internal::p4est::types< dim >::forest * | get_p4est () const |
virtual void | add_periodicity (const std::vector<::GridTools::PeriodicFacePair< cell_iterator >> &) override |
unsigned int | register_data_attach (const std::function< std::vector< char >(const cell_iterator &, const CellStatus)> &pack_callback, const bool returns_variable_size_data) |
void | notify_ready_to_unpack (const unsigned int handle, const std::function< void(const cell_iterator &, const CellStatus, const boost::iterator_range< std::vector< char >::const_iterator > &)> &unpack_callback) |
virtual const MPI_Comm & | get_communicator () const |
virtual void | copy_triangulation (const ::Triangulation< dim, spacedim > &old_tria) override |
virtual void | copy_triangulation (const Triangulation< dim, spacedim > &other_tria) |
unsigned int | n_locally_owned_active_cells () const |
virtual types::global_cell_index | n_global_active_cells () const override |
virtual unsigned int | n_global_levels () const override |
types::subdomain_id | locally_owned_subdomain () const override |
const std::set< types::subdomain_id > & | ghost_owners () const |
const std::set< types::subdomain_id > & | level_ghost_owners () const |
const Utilities::MPI::Partitioner & | global_active_cell_index_partitioner () const |
const Utilities::MPI::Partitioner & | global_level_cell_index_partitioner (const unsigned int level) const |
virtual std::map< unsigned int, std::set<::types::subdomain_id > > | compute_vertices_with_ghost_neighbors () const |
virtual std::vector< types::boundary_id > | get_boundary_ids () const override |
virtual std::vector< types::manifold_id > | get_manifold_ids () const override |
void | communicate_locally_moved_vertices (const std::vector< bool > &vertex_locally_moved) |
virtual void | set_mesh_smoothing (const MeshSmoothing mesh_smoothing) |
virtual const MeshSmoothing & | get_mesh_smoothing () const |
void | set_manifold (const types::manifold_id number, const Manifold< dim, spacedim > &manifold_object) |
void | reset_manifold (const types::manifold_id manifold_number) |
void | reset_all_manifolds () |
void | set_all_manifold_ids (const types::manifold_id number) |
void | set_all_manifold_ids_on_boundary (const types::manifold_id number) |
void | set_all_manifold_ids_on_boundary (const types::boundary_id b_id, const types::manifold_id number) |
const Manifold< dim, spacedim > & | get_manifold (const types::manifold_id number) const |
virtual void | create_triangulation_compatibility (const std::vector< Point< spacedim >> &vertices, const std::vector< CellData< dim >> &cells, const SubCellData &subcelldata) |
void | flip_all_direction_flags () |
template<> | |
unsigned int | n_quads () const |
template<> | |
unsigned int | n_quads (const unsigned int) const |
template<> | |
unsigned int | n_quads () const |
template<> | |
unsigned int | n_quads (const unsigned int) const |
template<> | |
unsigned int | n_quads () const |
template<> | |
unsigned int | n_quads (const unsigned int) const |
template<> | |
unsigned int | n_active_quads (const unsigned int) const |
template<> | |
unsigned int | n_active_quads () const |
template<> | |
unsigned int | n_active_quads (const unsigned int) const |
template<> | |
unsigned int | n_active_quads () const |
template<> | |
unsigned int | n_active_quads (const unsigned int) const |
template<> | |
unsigned int | n_active_quads () const |
template<> | |
unsigned int | n_hexs () const |
template<> | |
unsigned int | n_hexs (const unsigned int level) const |
template<> | |
unsigned int | n_active_hexs () const |
template<> | |
unsigned int | n_active_hexs (const unsigned int level) const |
template<> | |
unsigned int | max_adjacent_cells () const |
template<> | |
unsigned int | max_adjacent_cells () const |
template<> | |
unsigned int | max_adjacent_cells () const |
template<> | |
unsigned int | n_raw_quads (const unsigned int) const |
template<> | |
unsigned int | n_raw_quads (const unsigned int) const |
template<> | |
unsigned int | n_raw_quads (const unsigned int) const |
template<> | |
unsigned int | n_raw_quads (const unsigned int level) const |
template<> | |
unsigned int | n_raw_quads (const unsigned int level) const |
template<> | |
unsigned int | n_raw_quads (const unsigned int) const |
template<> | |
unsigned int | n_raw_quads () const |
template<> | |
unsigned int | n_raw_hexs (const unsigned int) const |
template<> | |
unsigned int | n_raw_hexs (const unsigned int) const |
template<> | |
unsigned int | n_raw_hexs (const unsigned int) const |
template<> | |
unsigned int | n_raw_hexs (const unsigned int level) const |
Internal information about the number of objects | |
template<class Archive > | |
void | save (Archive &ar, const unsigned int version) const |
template<class Archive > | |
void | load (Archive &ar, const unsigned int version) |
unsigned int | n_raw_lines () const |
unsigned int | n_raw_lines (const unsigned int level) const |
unsigned int | n_raw_quads () const |
unsigned int | n_raw_quads (const unsigned int level) const |
unsigned int | n_raw_hexs (const unsigned int level) const |
unsigned int | n_raw_cells (const unsigned int level) const |
unsigned int | n_raw_faces () const |
virtual void | add_periodicity (const std::vector< GridTools::PeriodicFacePair< cell_iterator >> &) |
const std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, std::bitset< 3 > > > & | get_periodic_face_map () const |
const std::vector< ReferenceCell::Type > & | get_reference_cell_types () const |
bool | all_reference_cell_types_are_hyper_cube () const |
template<class Archive > | |
void | serialize (Archive &archive, const unsigned int version) |
Mesh refinement | |
void | set_all_refine_flags () |
void | refine_global (const unsigned int times=1) |
void | coarsen_global (const unsigned int times=1) |
History of a triangulation | |
void | save_refine_flags (std::ostream &out) const |
void | save_refine_flags (std::vector< bool > &v) const |
void | load_refine_flags (std::istream &in) |
void | load_refine_flags (const std::vector< bool > &v) |
void | save_coarsen_flags (std::ostream &out) const |
void | save_coarsen_flags (std::vector< bool > &v) const |
void | load_coarsen_flags (std::istream &out) |
void | load_coarsen_flags (const std::vector< bool > &v) |
bool | get_anisotropic_refinement_flag () const |
User data | |
void | clear_user_flags () |
void | save_user_flags (std::ostream &out) const |
void | save_user_flags (std::vector< bool > &v) const |
void | load_user_flags (std::istream &in) |
void | load_user_flags (const std::vector< bool > &v) |
void | clear_user_flags_line () |
void | save_user_flags_line (std::ostream &out) const |
void | save_user_flags_line (std::vector< bool > &v) const |
void | load_user_flags_line (std::istream &in) |
void | load_user_flags_line (const std::vector< bool > &v) |
void | clear_user_flags_quad () |
void | save_user_flags_quad (std::ostream &out) const |
void | save_user_flags_quad (std::vector< bool > &v) const |
void | load_user_flags_quad (std::istream &in) |
void | load_user_flags_quad (const std::vector< bool > &v) |
void | clear_user_flags_hex () |
void | save_user_flags_hex (std::ostream &out) const |
void | save_user_flags_hex (std::vector< bool > &v) const |
void | load_user_flags_hex (std::istream &in) |
void | load_user_flags_hex (const std::vector< bool > &v) |
void | clear_user_data () |
void | save_user_indices (std::vector< unsigned int > &v) const |
void | load_user_indices (const std::vector< unsigned int > &v) |
void | save_user_pointers (std::vector< void *> &v) const |
void | load_user_pointers (const std::vector< void *> &v) |
void | save_user_indices_line (std::vector< unsigned int > &v) const |
void | load_user_indices_line (const std::vector< unsigned int > &v) |
void | save_user_indices_quad (std::vector< unsigned int > &v) const |
void | load_user_indices_quad (const std::vector< unsigned int > &v) |
void | save_user_indices_hex (std::vector< unsigned int > &v) const |
void | load_user_indices_hex (const std::vector< unsigned int > &v) |
void | save_user_pointers_line (std::vector< void *> &v) const |
void | load_user_pointers_line (const std::vector< void *> &v) |
void | save_user_pointers_quad (std::vector< void *> &v) const |
void | load_user_pointers_quad (const std::vector< void *> &v) |
void | save_user_pointers_hex (std::vector< void *> &v) const |
void | load_user_pointers_hex (const std::vector< void *> &v) |
Cell iterator functions | |
cell_iterator | begin (const unsigned int level=0) const |
active_cell_iterator | begin_active (const unsigned int level=0) const |
cell_iterator | end () const |
cell_iterator | end (const unsigned int level) const |
active_cell_iterator | end_active (const unsigned int level) const |
cell_iterator | last () const |
active_cell_iterator | last_active () const |
cell_iterator | create_cell_iterator (const CellId &cell_id) const |
Cell iterator functions returning ranges of iterators | |
IteratorRange< cell_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 functions | |
face_iterator | begin_face () const |
active_face_iterator | begin_active_face () const |
face_iterator | end_face () const |
IteratorRange< active_face_iterator > | active_face_iterators () const |
Vertex iterator functions | |
vertex_iterator | begin_vertex () const |
active_vertex_iterator | begin_active_vertex () const |
vertex_iterator | end_vertex () const |
Information about the triangulation | |
unsigned int | n_lines () const |
unsigned int | n_lines (const unsigned int level) const |
unsigned int | n_active_lines () const |
unsigned int | n_active_lines (const unsigned int level) const |
unsigned int | n_quads () const |
unsigned int | n_quads (const unsigned int level) const |
unsigned int | n_active_quads () const |
unsigned int | n_active_quads (const unsigned int level) const |
unsigned int | n_hexs () const |
unsigned int | n_hexs (const unsigned int level) const |
unsigned int | n_active_hexs () const |
unsigned int | n_active_hexs (const unsigned int level) const |
unsigned int | n_cells () const |
unsigned int | n_cells (const unsigned int level) const |
unsigned int | n_active_cells () const |
unsigned int | n_active_cells (const unsigned int level) const |
unsigned int | n_faces () const |
unsigned int | n_active_faces () const |
unsigned int | n_levels () const |
unsigned int | n_vertices () const |
const std::vector< Point< spacedim > > & | get_vertices () const |
unsigned int | n_used_vertices () const |
bool | vertex_used (const unsigned int index) const |
const std::vector< bool > & | get_used_vertices () const |
unsigned int | max_adjacent_cells () const |
Triangulation< dim, spacedim > & | get_triangulation () |
const Triangulation< dim, spacedim > & | get_triangulation () const |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
Static Public Member Functions | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Public Attributes | |
Keeping up with what happens to a triangulation | |
Signals | signals |
Static Public Attributes | |
static const unsigned int | dimension = dim |
static const unsigned int | space_dimension = spacedim |
Protected Types | |
using | cell_relation_t = typename std::tuple< CellStatus, cell_iterator > |
Protected Member Functions | |
void | save_attached_data (const unsigned int global_first_cell, const unsigned int global_num_cells, const std::string &filename) const |
void | load_attached_data (const unsigned int global_first_cell, const unsigned int global_num_cells, const unsigned int local_num_cells, const std::string &filename, const unsigned int n_attached_deserialize_fixed, const unsigned int n_attached_deserialize_variable) |
virtual void | update_number_cache () |
void | update_reference_cell_types () override |
void | reset_global_cell_indices () |
Private Member Functions | |
virtual void | update_cell_relations () override |
typename ::internal::p4est::types< dim >::tree * | init_tree (const int dealii_coarse_cell_index) const |
void | setup_coarse_cell_to_p4est_tree_permutation () |
void | copy_new_triangulation_to_p4est (std::integral_constant< int, 2 >) |
void | copy_new_triangulation_to_p4est (std::integral_constant< int, 3 >) |
void | copy_local_forest_to_triangulation () |
std::vector< unsigned int > | get_cell_weights () const |
std::vector< bool > | mark_locally_active_vertices_on_level (const int level) const |
virtual unsigned int | coarse_cell_id_to_coarse_cell_index (const types::coarse_cell_id coarse_cell_id) const override |
virtual types::coarse_cell_id | coarse_cell_index_to_coarse_cell_id (const unsigned int coarse_cell_index) const override |
Private Attributes | |
Settings | settings |
bool | triangulation_has_content |
typename ::internal::p4est::types< dim >::connectivity * | connectivity |
typename ::internal::p4est::types< dim >::forest * | parallel_forest |
typename ::internal::p4est::types< dim >::ghost * | parallel_ghost |
std::vector< types::global_dof_index > | coarse_cell_to_p4est_tree_permutation |
std::vector< types::global_dof_index > | p4est_tree_to_coarse_cell_permutation |
Friends | |
template<int , int , class > | |
class | ::FETools::internal::ExtrapolateImplementation |
Exceptions | |
MeshSmoothing | smooth_grid |
std::vector< ReferenceCell::Type > | reference_cell_types |
void | update_periodic_face_map () |
static ::ExceptionBase & | ExcInvalidLevel (int arg1, int arg2) |
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 ::ExceptionBase & | ExcInconsistentCoarseningFlags () |
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) |
This class acts like the Triangulation class, but it distributes the mesh across a number of different processors when using MPI. The class's interface does not add a lot to the Triangulation class but there are a number of difficult algorithms under the hood that ensure we always have a load-balanced, fully distributed mesh. Use of this class is explained in step-40, step-32, the Parallel computing with multiple processors using distributed memory documentation module, as well as the distributed_paper. See there for more information. This class satisfies the MeshType concept.
Refining and coarsening a distributed triangulation is a complicated process because cells may have to be migrated from one processor to another. On a single processor, materializing that part of the global mesh that we want to store here from what we have stored before therefore may involve several cycles of refining and coarsening the locally stored set of cells until we have finally gotten from the previous to the next triangulation. This process is described in more detail in the distributed_paper. Unfortunately, in this process, some information can get lost relating to flags that are set by user code and that are inherited from mother to child cell but that are not moved along with a cell if that cell is migrated from one processor to another.
An example are boundary indicators. Assume, for example, that you start with a single cell that is refined once globally, yielding four children. If you have four processors, each one owns one cell. Assume now that processor 1 sets the boundary indicators of the external boundaries of the cell it owns to 42. Since processor 0 does not own this cell, it doesn't set the boundary indicators of its ghost cell copy of this cell. Now, assume we do several mesh refinement cycles and end up with a configuration where this processor suddenly finds itself as the owner of this cell. If boundary indicator 42 means that we need to integrate Neumann boundary conditions along this boundary, then processor 0 will forget to do so because it has never set the boundary indicator along this cell's boundary to 42.
The way to avoid this dilemma is to make sure that things like setting boundary indicators or material ids is done immediately every time a parallel triangulation is refined. This is not necessary for sequential triangulations because, there, these flags are inherited from mother to child cell and remain with a cell even if it is refined and the children are later coarsened again, but this does not hold for distributed triangulations. It is made even more difficult by the fact that in the process of refining a parallel distributed triangulation, the triangulation may call Triangulation::execute_coarsening_and_refinement multiple times and this function needs to know about boundaries. In other words, it is not enough to just set boundary indicators on newly created faces only after calling distributed::parallel::TriangulationBase::execute_coarsening_and_refinement
: it actually has to happen while that function is still running.
The way to do this is by writing a function that sets boundary indicators and that will be called by the Triangulation class. The triangulation does not provide a pointer to itself to the function being called, nor any other information, so the trick is to get this information into the function. C++ provides a nice mechanism for this that is best explained using an example:
The object passed as argument to connect
is an object that can be called like a function with no arguments. It does so by wrapping a function that does, in fact, take an argument but this one argument is stored as a reference to the coarse grid triangulation when the lambda function is created. After each refinement step, the triangulation will then call the object so created which will in turn call set_boundary_ids<dim>
with the reference to the coarse grid as argument.
This approach can be generalized. In the example above, we have used a global function that will be called. However, sometimes it is necessary that this function is in fact a member function of the class that generates the mesh, for example because it needs to access run-time parameters. This can be achieved as follows: assuming the set_boundary_ids()
function has been declared as a (non- static, but possibly private) member function of the MyClass
class, then the following will work:
The lambda function above again is an object that can be called like a global function with no arguments, and this object in turn calls the current object's member function set_boundary_ids
with a reference to the triangulation to work on. Note that because the create_coarse_mesh
function is declared as const
, it is necessary that the set_boundary_ids
function is also declared const
.
Note:For reasons that have to do with the way the parallel::distributed::Triangulation is implemented, functions that have been attached to the post-refinement signal of the triangulation are called more than once, sometimes several times, every time the triangulation is actually refined.
using parallel::distributed::Triangulation< dim, spacedim >::CellStatus = typename ::Triangulation<dim, spacedim>::CellStatus |
|
protectedinherited |
This auxiliary data structure stores the relation between a deal.II cell and its current CellStatus. For an extensive description of the latter, see the documentation for the member function register_data_attach().
Definition at line 705 of file tria_base.h.
|
inherited |
enum parallel::distributed::Triangulation::Settings |
Configuration flags for distributed Triangulations to be set in the constructor. Settings can be combined using bitwise OR.
Enumerator | |
---|---|
default_setting | Default settings, other options are disabled. |
mesh_reconstruction_after_repartitioning | If set, the deal.II mesh will be reconstructed from the coarse mesh every time a repartitioning in p4est happens. This can be a bit more expensive, but guarantees the same memory layout and therefore cell ordering in the deal.II mesh. As assembly is done in the deal.II cell ordering, this flag is required to get reproducible behavior after snapshot/resume. |
construct_multigrid_hierarchy | This flags needs to be set to use the geometric multigrid functionality. This option requires additional computation and communication. |
no_automatic_repartitioning | Setting this flag will disable automatic repartitioning of the cells after a refinement cycle. It can be executed manually by calling repartition(). |
|
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: ![]() 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 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. |
|
explicit |
Constructor.
mpi_communicator | The MPI communicator to be used for the triangulation. |
smooth_grid | Degree and kind of mesh smoothing to be applied to the mesh. See the Triangulation class for a description of the kinds of smoothing operations that can be applied. |
settings | See the description of the Settings enumerator. Providing construct_multigrid_hierarchy enforces Triangulation::limit_level_difference_at_vertices for smooth_grid. |
check_for_distorted_cells
argument provided by the base class.
|
overridevirtual |
Destructor.
Reimplemented from Triangulation< dim, spacedim >.
|
overridevirtual |
Reset this triangulation into a virgin state by deleting all data.
Note that this operation is only allowed if no subscriptions to this object exist any more, such as DoFHandler objects using it.
Reimplemented from parallel::DistributedTriangulationBase< dim, spacedim >.
|
overridevirtual |
Return if multilevel hierarchy is supported and has been constructed.
Implements parallel::TriangulationBase< dim, spacedim >.
void Triangulation< dim, spacedim >::execute_transfer | ( | const typename ::internal::p4est::types< dim >::forest * | parallel_forest, |
const typename ::internal::p4est::types< dim >::gloidx * | previous_global_first_quadrant | ||
) |
Transfer data across forests.
Besides the actual parallel_forest
, which has been already refined and repartitioned, this function also needs information about its previous state, i.e. the locally owned intervals in p4est's sc_array of each processor. This information needs to be memcopyied out of the old p4est object and has to be provided via the parameter previous_global_first_quadrant
.
Data has to be previously packed with DistributedTriangulationBase::DataTransfer::pack_data().
|
overridevirtual |
Implementation of the same function as in the base class.
|
overridevirtual |
Create a triangulation as documented in the base class.
This function also sets up the various data structures necessary to distribute a mesh across a number of processors. This will be necessary once the mesh is being refined, though we will always keep the entire coarse mesh that is generated by this function on all processors.
Reimplemented from Triangulation< dim, spacedim >.
|
overridevirtual |
Create a triangulation as documented in the base class.
This function also sets up the various data structures necessary to distribute a mesh across a number of processors. This will be necessary once the mesh is being refined, though we will always keep the entire coarse mesh that is generated by this function on all processors.
Reimplemented from Triangulation< dim, spacedim >.
|
overridevirtual |
Coarsen and refine the mesh according to refinement and coarsening flags set.
Since the current processor only has control over those cells it owns (i.e. the ones for which cell->subdomain_id() == this->locally_owned_subdomain()
), refinement and coarsening flags are only respected for those locally owned cells. Flags may be set on other cells as well (and may often, in fact, if you call Triangulation::prepare_coarsening_and_refinement()) but will be largely ignored: the decision to refine the global mesh will only be affected by flags set on locally owned cells.
Reimplemented from Triangulation< dim, spacedim >.
|
overridevirtual |
Override the implementation of prepare_coarsening_and_refinement from the base class. This is necessary if periodic boundaries are enabled and the level difference over vertices over the periodic boundary must not be more than 2:1.
Reimplemented from Triangulation< dim, spacedim >.
void Triangulation< dim, spacedim >::repartition | ( | ) |
Manually repartition the active cells between processors. Normally this repartitioning will happen automatically when calling execute_coarsening_and_refinement() (or refine_global()) unless the no_automatic_repartitioning
is set in the constructor. Setting the flag and then calling repartition() gives the same result.
If you want to transfer data (using SolutionTransfer or manually with register_data_attach() and notify_ready_to_unpack()), you need to set it up twice: once when calling execute_coarsening_and_refinement(), which will handle coarsening and refinement but obviously won't ship any data between processors, and a second time when calling repartition(). Here, no coarsening and refinement will be done but information will be packed and shipped to different processors. In other words, you probably want to treat a call to repartition() in the same way as execute_coarsening_and_refinement() with respect to dealing with data movement (SolutionTransfer, etc.).
|
overridevirtual |
Return true if the triangulation has hanging nodes.
In the context of parallel distributed triangulations, every processor stores only that part of the triangulation it locally owns. However, it also stores the entire coarse mesh, and to guarantee the 2:1 relationship between cells, this may mean that there are hanging nodes between cells that are not locally owned or ghost cells (i.e., between ghost cells and artificial cells, or between artificial and artificial cells; see the glossary). One is not typically interested in this case, so the function returns whether there are hanging nodes between any two cells of the "global" mesh, i.e., the union of locally owned cells on all processors.
Reimplemented from Triangulation< dim, spacedim >.
|
overridevirtual |
Return the local memory consumption in bytes.
Reimplemented from parallel::TriangulationBase< dim, spacedim >.
|
virtual |
Return the local memory consumption contained in the p4est data structures alone. This is already contained in memory_consumption() but made available separately for debugging purposes.
void Triangulation< dim, spacedim >::write_mesh_vtk | ( | const std::string & | file_basename | ) | const |
unsigned int Triangulation< dim, spacedim >::get_checksum | ( | ) | const |
|
overridevirtual |
Save the refinement information from the coarse mesh into the given file. This file needs to be reachable from all nodes in the computation on a shared network file system. See the SolutionTransfer class on how to store solution vectors into this file. Additional cell-based data can be saved using DistributedTriangulationBase::DataTransfer::register_data_attach().
Implements parallel::DistributedTriangulationBase< dim, spacedim >.
|
overridevirtual |
Load the refinement information saved with save() back in. The mesh must contain the same coarse mesh that was used in save() before calling this function.
You do not need to load with the same number of MPI processes that you saved with. Rather, if a mesh is loaded with a different number of MPI processes than used at the time of saving, the mesh is repartitioned appropriately. Cell-based data that was saved with DistributedTriangulationBase::DataTransfer::register_data_attach() can be read in with DistributedTriangulationBase::DataTransfer::notify_ready_to_unpack() after calling load().
If you use p4est version > 0.3.4.2 the autopartition
flag tells p4est to ignore the partitioning that the triangulation had when it was saved and make it uniform upon loading. If autopartition
is set to false, the triangulation is only repartitioned if needed (i.e. if a different number of MPI processes is encountered).
Implements parallel::DistributedTriangulationBase< dim, spacedim >.
const std::vector< types::global_dof_index > & Triangulation< dim, spacedim >::get_p4est_tree_to_coarse_cell_permutation | ( | ) | const |
Return a permutation vector for the order the coarse cells are handed off to p4est. For example the value of the \(i\)th element in this vector is the index of the deal.II coarse cell (counting from begin(0)) that corresponds to the \(i\)th tree managed by p4est.
const std::vector< types::global_dof_index > & Triangulation< dim, spacedim >::get_coarse_cell_to_p4est_tree_permutation | ( | ) | const |
const ::internal::p4est::types< dim >::forest * Triangulation< dim, spacedim >::get_p4est | ( | ) | const |
|
overridevirtual |
In addition to the action in the base class Triangulation, this function joins faces in the p4est forest for periodic boundary conditions. As a result, each pair of faces will differ by at most one refinement level and ghost neighbors will be available across these faces.
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.
|
overrideprivatevirtual |
Go through all p4est trees and store the relations between the status of locally owned quadrants and cells in the private member local_cell_relations.
The stored vector will be ordered by the occurrence of quadrants in the corresponding local sc_array of the parallel_forest. p4est requires this specific ordering for its transfer functions.
Implements parallel::DistributedTriangulationBase< dim, spacedim >.
|
private |
|
private |
|
private |
Take the contents of a newly created triangulation we are attached to and copy it to p4est data structures.
This function exists in 2d and 3d variants.
|
private |
|
private |
|
private |
Internal function notifying all registered slots to provide their weights before repartitioning occurs. Called from execute_coarsening_and_refinement() and repartition().
|
private |
This method returns a bit vector of length tria.n_vertices() indicating the locally active vertices on a level, i.e., the vertices touched by the locally owned level cells for use in geometric multigrid (possibly including the vertices due to periodic boundary conditions) are marked by true.
Used by DoFHandler::Policy::ParallelDistributed.
ensure that if one of the two vertices on a periodic face is marked as active (i.e., belonging to an owned level cell), also the other one is active
|
overrideprivatevirtual |
Translate the unique id of a coarse cell to its index. See the glossary entry on coarse cell IDs for more information.
coarse_cell_id | Unique id of the coarse cell. |
Reimplemented from Triangulation< dim, spacedim >.
|
overrideprivatevirtual |
Translate the index of coarse cell to its unique id. See the glossary entry on coarse cell IDs for more information.
coarse_cell_index | Index of the coarse cell. |
Reimplemented from Triangulation< dim, spacedim >.
|
inherited |
Write the data of this object to a stream for the purpose of serialization using the BOOST serialization library.
|
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.
|
inherited |
Register a function that can be used to attach data of fixed size to cells. This is useful for two purposes: (i) Upon refinement and coarsening of a triangulation (e.g. in parallel::distributed::Triangulation::execute_coarsening_and_refinement()), one needs to be able to store one or more data vectors per cell that characterizes the solution values on the cell so that this data can then be transferred to the new owning processor of the cell (or its parent/children) when the mesh is re-partitioned; (ii) when serializing a computation to a file, it is necessary to attach data to cells so that it can be saved (e.g. in parallel::distributed::Triangulation::save()) along with the cell's other information and, if necessary, later be reloaded from disk with a different subdivision of cells among the processors.
The way this function works is that it allows any number of interest parties to register their intent to attach data to cells. One example of classes that do this is parallel::distributed::SolutionTransfer where each parallel::distributed::SolutionTransfer object that works on the current Triangulation object then needs to register its intent. Each of these parties registers a callback function (the first argument here, pack_callback
) that will be called whenever the triangulation's execute_coarsening_and_refinement() or save() functions are called.
The current function then returns an integer handle that corresponds to the number of data set that the callback provided here will attach. While this number could be given a precise meaning, this is not important: You will never actually have to do anything with this number except return it to the notify_ready_to_unpack() function. In other words, each interested party (i.e., the caller of the current function) needs to store their respective returned handle for later use when unpacking data in the callback provided to notify_ready_to_unpack().
Whenever pack_callback
is then called by execute_coarsening_and_refinement() or load() on a given cell, it receives a number of arguments. In particular, the first argument passed to the callback indicates the cell for which it is supposed to attach data. This is always an active cell.
The second, CellStatus, argument provided to the callback function will tell you if the given cell will be coarsened, refined, or will persist as is. (This status may be different than the refinement or coarsening flags set on that cell, to accommodate things such as the "one hanging node per edge" rule.). These flags need to be read in context with the p4est quadrant they belong to, as their relations are gathered in local_quadrant_cell_relations.
Specifically, the values for this argument mean the following:
CELL_PERSIST
: The cell won't be refined/coarsened, but might be moved to a different processor. If this is the case, the callback will want to pack up the data on this cell into an array and store it at the provided address for later unpacking wherever this cell may land.CELL_REFINE
: This cell will be refined into 4 or 8 cells (in 2d and 3d, respectively). However, because these children don't exist yet, you cannot access them at the time when the callback is called. Thus, in local_quadrant_cell_relations, the corresponding p4est quadrants of the children cells are linked to the deal.II cell which is going to be refined. To be specific, only the very first child is marked with CELL_REFINE
, whereas the others will be marked with CELL_INVALID
, which indicates that these cells will be ignored by default during the packing or unpacking process. This ensures that data is only transferred once onto or from the parent cell. If the callback is called with CELL_REFINE
, the callback will want to pack up the data on this cell into an array and store it at the provided address for later unpacking in a way so that it can then be transferred to the children of the cell that will then be available. In other words, if the data the callback will want to pack up corresponds to a finite element field, then the prolongation from parent to (new) children will have to happen during unpacking.CELL_COARSEN
: The children of this cell will be coarsened into the given cell. These children still exist, so if this is the value given to the callback as second argument, the callback will want to transfer data from the children to the current parent cell and pack it up so that it can later be unpacked again on a cell that then no longer has any children (and may also be located on a different processor). In other words, if the data the callback will want to pack up corresponds to a finite element field, then it will need to do the restriction from children to parent at this point.CELL_INVALID
: See CELL_REFINE
.CELL_PERSIST
.The callback function is expected to return a memory chunk of the format std::vector<char>
, representing the packed data on a certain cell.
The second parameter returns_variable_size_data
indicates whether the returned size of the memory region from the callback function varies by cell (=true
) or stays constant on each one throughout the whole domain (=false
).
Definition at line 702 of file tria_base.cc.
|
inherited |
This function is the opposite of register_data_attach(). It is called after the execute_coarsening_and_refinement() or save()/load() functions are done when classes and functions that have previously attached data to a triangulation for either transfer to other processors, across mesh refinement, or serialization of data to a file are ready to receive that data back. The important part about this process is that the triangulation cannot do this right away from the end of execute_coarsening_and_refinement() or load() via a previously attached callback function (as the register_data_attach() function does) because the classes that eventually want the data back may need to do some setup between the point in time where the mesh has been recreated and when the data can actually be received. An example is the parallel::distributed::SolutionTransfer class that can really only receive the data once not only the mesh is completely available again on the current processor, but only after a DoFHandler has been reinitialized and distributed degrees of freedom. In other words, there is typically a significant amount of set up that needs to happen in user space before the classes that can receive data attached to cell are ready to actually do so. When they are, they use the current function to tell the triangulation object that now is the time when they are ready by calling the current function.
The supplied callback function is then called for each newly locally owned cell. The first argument to the callback is an iterator that designates the cell; the second argument indicates the status of the cell in question; and the third argument localizes a memory area by two iterators that contains the data that was previously saved from the callback provided to register_data_attach().
The CellStatus will indicate if the cell was refined, coarsened, or persisted unchanged. The cell_iterator
argument to the callback will then either be an active, locally owned cell (if the cell was not refined), or the immediate parent if it was refined during execute_coarsening_and_refinement(). Therefore, contrary to during register_data_attach(), you can now access the children if the status is CELL_REFINE
but no longer for callbacks with status CELL_COARSEN
.
The first argument to this function, handle
, corresponds to the return value of register_data_attach(). (The precise meaning of what the numeric value of this handle is supposed to represent is neither important, nor should you try to use it for anything other than transmit information between a call to register_data_attach() to the corresponding call to notify_ready_to_unpack().)
Definition at line 731 of file tria_base.cc.
|
protectedinherited |
Save additional cell-attached data into the given file. The first arguments are used to determine the offsets where to write buffers to.
Called by save.
Definition at line 631 of file tria_base.cc.
|
protectedinherited |
Load additional cell-attached data from the given file, if any was saved. The first arguments are used to determine the offsets where to read buffers from.
Called by load.
Definition at line 667 of file tria_base.cc.
|
virtualinherited |
Return MPI communicator used by this triangulation.
Definition at line 139 of file tria_base.cc.
|
overridevirtualinherited |
Implementation of the same function as in the base class.
Definition at line 64 of file tria_base.cc.
|
virtualinherited |
Copy other_tria
to this triangulation. This operation is not cheap, so you should be careful with using this. We do not implement this function as a copy constructor, since it makes it easier to maintain collections of triangulations if you can assign them values later on.
Keep in mind that this function also copies the pointer to the boundary descriptor previously set by the set_manifold
function. You must therefore also guarantee that the Manifold objects describing the boundary have a lifetime at least as long as the copied triangulation.
This triangulation must be empty beforehand.
The function is made virtual
since some derived classes might want to disable or extend the functionality of this function.
Reimplemented in PersistentTriangulation< dim, spacedim >.
|
inherited |
Return the number of active cells in the triangulation that are locally owned, i.e. that have a subdomain_id equal to locally_owned_subdomain(). Note that there may be more active cells in the triangulation stored on the present processor, such as for example ghost cells, or cells further away from the locally owned block of cells but that are needed to ensure that the triangulation that stores this processor's set of active cells still remains balanced with respect to the 2:1 size ratio of adjacent cells.
As a consequence of the remark above, the result of this function is always smaller or equal to the result of the function with the same name in the Triangulation base class, which includes the active ghost and artificial cells (see also GlossArtificialCell and GlossGhostCell).
Definition at line 118 of file tria_base.cc.
|
overridevirtualinherited |
Return the sum over all processors of the number of active cells owned by each processor. This equals the overall number of active cells in the triangulation.
Reimplemented from Triangulation< dim, spacedim >.
Definition at line 132 of file tria_base.cc.
|
overridevirtualinherited |
Return the global maximum level. This may be bigger than the number Triangulation::n_levels() (a function in this class's base class) returns if the current processor only stores cells in parts of the domain that are not very refined, but if other processors store cells in more deeply refined parts of the domain.
Reimplemented from Triangulation< dim, spacedim >.
Definition at line 125 of file tria_base.cc.
|
overridevirtualinherited |
Return the subdomain id of those cells that are owned by the current processor. All cells in the triangulation that do not have this subdomain id are either owned by another processor or have children that only exist on other processors.
Reimplemented from Triangulation< dim, spacedim >.
Definition at line 319 of file tria_base.cc.
|
inherited |
Return a set of MPI ranks of the processors that have at least one ghost cell adjacent to the cells of the local processor. In other words, this is the set of subdomain_id() for all ghost cells.
The returned sets are symmetric, that is if i
is contained in the list of processor j
, then j
will also be contained in the list of processor i
.
Definition at line 328 of file tria_base.cc.
|
inherited |
Return a set of MPI ranks of the processors that have at least one level ghost cell adjacent to our cells used in geometric multigrid. In other words, this is the set of level_subdomain_id() for all level ghost cells.
The returned sets are symmetric, that is if i
is contained in the list of processor j
, then j
will also be contained in the list of processor i
.
Definition at line 337 of file tria_base.cc.
|
inherited |
Return partitioner for the global indices of the cells on the active level of the triangulation.
Definition at line 583 of file tria_base.cc.
|
inherited |
Return partitioner for the global indices of the cells on the given level
of the triangulation.
Definition at line 590 of file tria_base.cc.
|
virtualinherited |
Return a map that, for each vertex, lists all the processors whose subdomains are adjacent to that vertex.
Definition at line 346 of file tria_base.cc.
|
overridevirtualinherited |
Return a vector containing all boundary indicators assigned to boundary faces of active cells of this Triangulation object. Note, that each boundary indicator is reported only once. The size of the return vector will represent the number of different indicators (which is greater or equal one).
Reimplemented from Triangulation< dim, spacedim >.
Definition at line 356 of file tria_base.cc.
|
overridevirtualinherited |
Return a vector containing all manifold indicators assigned to the objects of the active cells of this Triangulation. Note, that each manifold indicator is reported only once. The size of the return vector will represent the number of different indicators (which is greater or equal one).
Reimplemented from Triangulation< dim, spacedim >.
Definition at line 367 of file tria_base.cc.
|
inherited |
When vertices have been moved locally, for example using code like
then this function can be used to update the location of vertices between MPI processes.
All the vertices that have been moved and might be in the ghost layer of a process have to be reported in the vertex_locally_moved
argument. This ensures that that part of the information that has to be send between processes is actually sent. Additionally, it is quite important that vertices on the boundary between processes are reported on exactly one process (e.g. the one with the highest id). Otherwise we could expect undesirable results if multiple processes move a vertex differently. A typical strategy is to let processor \(i\) move those vertices that are adjacent to cells whose owners include processor \(i\) but no other processor \(j\) with \(j<i\); in other words, for vertices at the boundary of a subdomain, the processor with the lowest subdomain id "owns" a vertex.
vertex_locally_moved
argument may not contain vertices that aren't at least on ghost cells.vertex_locally_moved | A bitmap indicating which vertices have been moved. The size of this array must be equal to Triangulation::n_vertices() and must be a subset of those vertices flagged by GridTools::get_locally_owned_vertices(). |
Definition at line 533 of file tria_base.cc.
|
protectedvirtualinherited |
Update the number_cache variable after mesh creation or refinement.
Definition at line 147 of file tria_base.cc.
|
overrideprotectedvirtualinherited |
Update the internal reference_cell_types vector.
Reimplemented from Triangulation< dim, spacedim >.
Definition at line 292 of file tria_base.cc.
|
protectedinherited |
Reset global active cell indices and global level cell indices.
Definition at line 378 of file tria_base.cc.
|
virtualinherited |
|
virtualinherited |
|
virtualinherited |
For backward compatibility, only. This function takes the cell data in the ordering as requested by deal.II versions up to 5.2, converts it to the new (lexicographic) ordering and calls create_triangulation().
Reimplemented in PersistentTriangulation< dim, spacedim >.
|
inherited |
Revert or flip the direction_flags of a dim<spacedim triangulation, see GlossDirectionFlag.
This function throws an exception if dim equals spacedim.
|
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.
|
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.
times > 1
) . See the section on signals in the general documentation of this class.
|
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.
times > 1
) . See the section on signals in the general documentation of this class.
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
Clear all user flags. See also GlossUserFlags.
|
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.
|
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.
|
inherited |
Read the information stored by save_user_flags
. See also GlossUserFlags.
|
inherited |
Read the information stored by save_user_flags
. See also GlossUserFlags.
|
inherited |
Clear all user flags on lines. See also GlossUserFlags.
|
inherited |
Save the user flags on lines. See also GlossUserFlags.
|
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.
|
inherited |
Load the user flags located on lines. See also GlossUserFlags.
|
inherited |
Load the user flags located on lines. See also GlossUserFlags.
|
inherited |
Clear all user flags on quads. See also GlossUserFlags.
|
inherited |
Save the user flags on quads. See also GlossUserFlags.
|
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.
|
inherited |
Load the user flags located on quads. See also GlossUserFlags.
|
inherited |
Load the user flags located on quads. See also GlossUserFlags.
|
inherited |
Clear all user flags on quads. See also GlossUserFlags.
|
inherited |
Save the user flags on hexs. See also GlossUserFlags.
|
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.
|
inherited |
Load the user flags located on hexs. See also GlossUserFlags.
|
inherited |
Load the user flags located on hexs. See also GlossUserFlags.
|
inherited |
Clear all user pointers and indices and allow the use of both for next access. See also GlossUserData.
|
inherited |
Save all user indices. The output vector is resized if necessary. See also GlossUserData.
|
inherited |
Read the information stored by save_user_indices(). See also GlossUserData.
|
inherited |
Save all user pointers. The output vector is resized if necessary. See also GlossUserData.
|
inherited |
Read the information stored by save_user_pointers(). See also GlossUserData.
|
inherited |
Save the user indices on lines. The output vector is resized if necessary. See also GlossUserData.
|
inherited |
Load the user indices located on lines. See also GlossUserData.
|
inherited |
Save the user indices on quads. The output vector is resized if necessary. See also GlossUserData.
|
inherited |
Load the user indices located on quads. See also GlossUserData.
|
inherited |
Save the user indices on hexes. The output vector is resized if necessary. See also GlossUserData.
|
inherited |
Load the user indices located on hexs. See also GlossUserData.
|
inherited |
Save the user indices on lines. The output vector is resized if necessary. See also GlossUserData.
|
inherited |
Load the user pointers located on lines. See also GlossUserData.
|
inherited |
Save the user pointers on quads. The output vector is resized if necessary. See also GlossUserData.
|
inherited |
Load the user pointers located on quads. See also GlossUserData.
|
inherited |
Save the user pointers on hexes. The output vector is resized if necessary. See also GlossUserData.
|
inherited |
Load the user pointers located on hexs. See also GlossUserData.
|
inherited |
Iterator to the first used cell on level level
.
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.
|
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
have zero iterations, as may be expected if there are no active cells on this level.
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.
|
inherited |
|
inherited |
Return an iterator which is the first iterator not on level. If level
is the last level, then this returns end()
.
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.
|
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()
.
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.
|
inherited |
|
inherited |
|
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.
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
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.
|
inherited |
|
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.
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
Return the total number of active cells. Maps to n_active_lines()
in one space dimension and so on.
|
inherited |
|
inherited |
|
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.
|
inherited |
Return the number of levels in this triangulation.
|
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.
|
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().
|
inherited |
|
inherited |
Return true
if the vertex with this index
is used.
|
inherited |
|
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.
|
inherited |
|
inherited |
|
inherited |
|
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).
|
inherited |
|
inherited |
Total number of lines, used or unused.
|
inherited |
Number of lines, used or unused, on the given level.
|
inherited |
Total number of quads, used or unused.
|
inherited |
Number of quads, used or unused, on the given level.
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
Number of hexs, used or unused, on the given level.
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
Number of cells, used or unused, on the given level.
|
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.
|
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.
|
inherited |
|
inherited |
|
inherited |
|
inherited |
Write and read the data of this object from a stream for the purpose of serialization. using the BOOST serialization library.
|
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.
|
inherited |
Unsubscribes a user from the object.
identifier
and the validity
pointer must be the same as the one supplied to subscribe(). Definition at line 156 of file subscriptor.cc.
|
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.
|
inlineinherited |
List the subscribers to the input stream
.
Definition at line 318 of file subscriptor.h.
|
inherited |
List the subscribers to deallog
.
Definition at line 204 of file subscriptor.cc.
|
friend |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
Two arrays that store which p4est tree corresponds to which coarse grid cell and vice versa. We need these arrays because p4est goes with the original order of coarse cells when it sets up its forest, and then applies the Morton ordering within each tree. But if coarse grid cells are badly ordered this may mean that individual parts of the forest stored on a local machine may be split across coarse grid cells that are not geometrically close. Consequently, we apply a hierarchical preordering according to SparsityTools::reorder_hierarchical() to ensure that the part of the forest stored by p4est is located on geometrically close coarse grid cells.
|
private |
|
protectedinherited |
Vector of tuples, which each contain a deal.II cell and their relation after refinement. To update its contents, use the compute_cell_relations member function.
The size of this vector is assumed to be equal to the number of locally owned quadrants in the parallel_forest object.
Definition at line 715 of file tria_base.h.
|
protectedinherited |
Definition at line 748 of file tria_base.h.
|
protectedinherited |
Definition at line 903 of file tria_base.h.
|
protectedinherited |
MPI communicator to be used for the triangulation. We create a unique communicator for this class, which is a duplicate of the one passed to the constructor.
Definition at line 311 of file tria_base.h.
|
protectedinherited |
The subdomain id to be used for the current processor. This is the MPI rank.
Definition at line 317 of file tria_base.h.
|
protectedinherited |
The total number of subdomains (or the size of the MPI communicator).
Definition at line 322 of file tria_base.h.
|
protectedinherited |
Definition at line 369 of file tria_base.h.
|
staticinherited |
|
staticinherited |
|
mutableinherited |