|
| Triangulation (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 |
|
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 () |
|
void | communicate_locally_moved_vertices (const std::vector< bool > &vertex_locally_moved) |
|
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 |
|
void | save (const std::string &filename) const |
|
void | load (const std::string &filename, const bool autopartition=true) |
|
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) |
|
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 typename ::internal::p4est::types< dim >::forest * | get_p4est () const |
|
virtual void | add_periodicity (const std::vector<::GridTools::PeriodicFacePair< cell_iterator >> &) override |
|
| DistributedTriangulationBase (MPI_Comm mpi_communicator, const typename ::Triangulation< dim, spacedim >::MeshSmoothing smooth_grid=(::Triangulation< dim, spacedim >::none), const bool check_for_distorted_cells=false) |
|
| TriangulationBase (MPI_Comm mpi_communicator, const typename ::Triangulation< dim, spacedim >::MeshSmoothing smooth_grid=(::Triangulation< dim, spacedim >::none), const bool check_for_distorted_cells=false) |
|
virtual | ~TriangulationBase () override |
|
virtual const MPI_Comm & | get_communicator () const |
|
virtual void | copy_triangulation (const ::Triangulation< dim, spacedim > &old_tria) override |
|
unsigned int | n_locally_owned_active_cells () const |
|
virtual types::global_cell_index | n_global_active_cells () const override |
|
virtual unsigned int | n_global_levels () const override |
|
types::subdomain_id | locally_owned_subdomain () const override |
|
const std::set< types::subdomain_id > & | ghost_owners () const |
|
const std::set< types::subdomain_id > & | level_ghost_owners () const |
|
virtual std::map< unsigned int, std::set<::types::subdomain_id > > | compute_vertices_with_ghost_neighbors () const |
|
virtual std::vector< types::boundary_id > | get_boundary_ids () const override |
|
virtual std::vector< types::manifold_id > | get_manifold_ids () const override |
|
| Triangulation (const MeshSmoothing smooth_grid=none, const bool check_for_distorted_cells=false) |
|
| Triangulation (const Triangulation< dim, spacedim > &)=delete |
|
| Triangulation (Triangulation< dim, spacedim > &&tria) noexcept |
|
Triangulation & | operator= (Triangulation< dim, spacedim > &&tria) noexcept |
|
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 std::vector< types::boundary_id > | get_boundary_ids () const |
|
virtual 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 () |
|
bool | prepare_coarsening_and_refinement () |
|
bool | prepare_coarsening_and_refinement () |
|
bool | prepare_coarsening_and_refinement () |
|
unsigned int | n_quads () const |
|
unsigned int | n_quads (const unsigned int) const |
|
unsigned int | n_quads () const |
|
unsigned int | n_quads (const unsigned int) const |
|
unsigned int | n_quads () const |
|
unsigned int | n_quads (const unsigned int) const |
|
unsigned int | n_active_quads (const unsigned int) const |
|
unsigned int | n_active_quads () const |
|
unsigned int | n_active_quads (const unsigned int) const |
|
unsigned int | n_active_quads () const |
|
unsigned int | n_active_quads (const unsigned int) const |
|
unsigned int | n_active_quads () const |
|
unsigned int | n_hexs () const |
|
unsigned int | n_hexs (const unsigned int level) const |
|
unsigned int | n_active_hexs () const |
|
unsigned int | n_active_hexs (const unsigned int level) const |
|
unsigned int | max_adjacent_cells () const |
|
unsigned int | max_adjacent_cells () const |
|
unsigned int | max_adjacent_cells () const |
|
unsigned int | n_raw_lines (const unsigned int level) const |
|
unsigned int | n_raw_lines () const |
|
unsigned int | n_raw_lines (const unsigned int level) const |
|
unsigned int | n_raw_lines () const |
|
unsigned int | n_raw_lines (const unsigned int level) const |
|
unsigned int | n_raw_lines () const |
|
unsigned int | n_raw_quads (const unsigned int) const |
|
unsigned int | n_raw_quads (const unsigned int) const |
|
unsigned int | n_raw_quads (const unsigned int) const |
|
unsigned int | n_raw_quads (const unsigned int level) const |
|
unsigned int | n_raw_quads (const unsigned int level) const |
|
unsigned int | n_raw_quads (const unsigned int) const |
|
unsigned int | n_raw_quads () const |
|
unsigned int | n_raw_hexs (const unsigned int) const |
|
unsigned int | n_raw_hexs (const unsigned int) const |
|
unsigned int | n_raw_hexs (const unsigned int) const |
|
unsigned int | n_raw_hexs (const unsigned int level) const |
|
void | set_all_refine_flags () |
|
void | refine_global (const unsigned int times=1) |
|
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 |
|
IteratorRange< active_face_iterator > | active_face_iterators () const |
|
vertex_iterator | begin_vertex () const |
|
active_vertex_iterator | begin_active_vertex () const |
|
vertex_iterator | end_vertex () const |
|
unsigned int | n_lines () const |
|
unsigned int | n_lines (const unsigned int level) const |
|
unsigned int | n_active_lines () const |
|
unsigned int | n_active_lines (const unsigned int level) const |
|
unsigned int | n_quads () const |
|
unsigned int | n_quads (const unsigned int level) const |
|
unsigned int | n_active_quads () const |
|
unsigned int | n_active_quads (const unsigned int level) const |
|
unsigned int | n_hexs () const |
|
unsigned int | n_hexs (const unsigned int level) const |
|
unsigned int | n_active_hexs () const |
|
unsigned int | n_active_hexs (const unsigned int level) const |
|
unsigned int | n_cells () const |
|
unsigned int | n_cells (const unsigned int level) const |
|
unsigned int | n_active_cells () const |
|
unsigned int | n_active_cells (const unsigned int level) const |
|
virtual types::global_cell_index | n_global_active_cells () const |
|
unsigned int | n_faces () const |
|
unsigned int | n_active_faces () const |
|
unsigned int | n_levels () const |
|
virtual unsigned int | n_global_levels () const |
|
virtual bool | has_hanging_nodes () const |
|
unsigned int | n_vertices () const |
|
const std::vector< Point< spacedim > > & | get_vertices () const |
|
unsigned int | n_used_vertices () const |
|
bool | vertex_used (const unsigned int index) const |
|
const std::vector< bool > & | get_used_vertices () const |
|
unsigned int | max_adjacent_cells () const |
|
virtual types::subdomain_id | locally_owned_subdomain () const |
|
Triangulation< dim, spacedim > & | get_triangulation () |
|
const Triangulation< dim, spacedim > & | get_triangulation () const |
|
unsigned int | n_raw_lines () const |
|
unsigned int | n_raw_lines (const unsigned int level) const |
|
unsigned int | n_raw_quads () const |
|
unsigned int | n_raw_quads (const unsigned int level) const |
|
unsigned int | n_raw_hexs (const unsigned int level) const |
|
unsigned int | n_raw_cells (const unsigned int level) const |
|
unsigned int | n_raw_faces () const |
|
virtual std::size_t | memory_consumption () const |
|
void | save (Archive &ar, const unsigned int version) const |
|
void | load (Archive &ar, const unsigned int version) |
|
virtual void | add_periodicity (const std::vector< GridTools::PeriodicFacePair< cell_iterator >> &) |
|
const std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, std::bitset< 3 > > > & | get_periodic_face_map () const |
|
void | serialize (Archive &archive, const unsigned int version) |
|
| Subscriptor () |
|
| Subscriptor (const Subscriptor &) |
|
| Subscriptor (Subscriptor &&) noexcept |
|
virtual | ~Subscriptor () |
|
Subscriptor & | operator= (const Subscriptor &) |
|
Subscriptor & | operator= (Subscriptor &&) noexcept |
|
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
|
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
|
unsigned int | n_subscriptions () const |
|
template<typename StreamType > |
void | list_subscribers (StreamType &stream) const |
|
void | list_subscribers () const |
|
template<class Archive > |
void | serialize (Archive &ar, const unsigned int version) |
|
template<int dim, int spacedim = dim>
class parallel::distributed::Triangulation< dim, spacedim >
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.
- Note
- This class does not support anisotropic refinement, because it relies on the p4est library that does not support this. Attempts to refine cells anisotropically will result in errors.
-
There is currently no support for distributing 1d triangulations.
Interaction with boundary description
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:
#include <functional>
template <int dim>
void set_boundary_ids (
{
}
template <int dim>
void
MyClass<dim>::create_coarse_mesh (
{
... create the coarse mesh ...
coarse_grid.
signals.post_refinement.connect(
[&coarse_grid](){
set_boundary_ids<dim>(coarse_grid);
});
}
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:
#include <functional>
template <int dim>
void
MyClass<dim>::set_boundary_ids (
{
}
template <int dim>
void
MyClass<dim>::create_coarse_mesh (
{
... create the coarse mesh ...
coarse_grid.
signals.post_refinement.connect(
[this, &coarse_grid]()
{
this->set_boundary_ids(coarse_grid);
});
}
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.
- Author
- Wolfgang Bangerth, Timo Heister 2008, 2009, 2010, 2011
Definition at line 242 of file tria.h.
template<int dim, int spacedim>
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 (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 (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
.
- Note
- If this function is used for serialization of data using save() and load(), then the cell status argument with which the callback is called will always be
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
).
- Note
- The purpose of this function is to register intent to attach data for a single, subsequent call to execute_coarsening_and_refinement() and notify_ready_to_unpack(), save(), load(). Consequently, notify_ready_to_unpack(), save(), and load() all forget the registered callbacks once these callbacks have been called, and you will have to re-register them with a triangulation if you want them to be active for another call to these functions.
Definition at line 4348 of file tria.cc.