Reference documentation for deal.II version 8.5.1
|
#include <deal.II/distributed/tria.h>
Public Member Functions | |
Triangulation (MPI_Comm mpi_communicator, const typename ::Triangulation< dim, spacedim >::MeshSmoothing smooth_grid=(::Triangulation< dim, spacedim >::none), const Settings settings=default_setting) | |
virtual | ~Triangulation () |
virtual void | clear () |
virtual void | copy_triangulation (const ::Triangulation< dim, spacedim > &other_tria) |
virtual void | create_triangulation (const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata) |
virtual void | execute_coarsening_and_refinement () |
virtual bool | prepare_coarsening_and_refinement () |
void | repartition () |
void | communicate_locally_moved_vertices (const std::vector< bool > &vertex_locally_moved) |
virtual bool | has_hanging_nodes () const |
virtual std::size_t | memory_consumption () const |
virtual std::size_t | memory_consumption_p4est () const |
void | write_mesh_vtk (const char *file_basename) const |
unsigned int | get_checksum () const |
void | save (const char *filename) const |
void | load (const char *filename, const bool autopartition=true) |
unsigned int | register_data_attach (const std::size_t size, const std_cxx11::function< void(const cell_iterator &, const CellStatus, void *)> &pack_callback) |
void | notify_ready_to_unpack (const unsigned int offset, const std_cxx11::function< void(const cell_iterator &, const CellStatus, const void *)> &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 |
virtual void | add_periodicity (const std::vector< GridTools::PeriodicFacePair< cell_iterator > > &) |
Public Member Functions inherited from parallel::Triangulation< dim, spacedim > | |
Triangulation (MPI_Comm mpi_communicator, const typename ::Triangulation< dim, spacedim >::MeshSmoothing smooth_grid=(::Triangulation< dim, spacedim >::none), const bool check_for_distorted_cells=false) | |
virtual MPI_Comm | get_communicator () const |
virtual void | copy_triangulation (const ::Triangulation< dim, spacedim > &old_tria) |
const std::vector< unsigned int > & | n_locally_owned_active_cells_per_processor () const |
unsigned int | n_locally_owned_active_cells () const |
virtual types::global_dof_index | n_global_active_cells () const |
virtual unsigned int | n_global_levels () const |
types::subdomain_id | locally_owned_subdomain () const |
const std::set< types::subdomain_id > & | ghost_owners () const |
const std::set< types::subdomain_id > & | level_ghost_owners () const |
Public Member Functions inherited from Triangulation< dim, spacedim > | |
Triangulation (const MeshSmoothing smooth_grid=none, const bool check_for_distorted_cells=false) | |
Triangulation (const Triangulation< dim, spacedim > &t) | |
Triangulation (Triangulation< dim, spacedim > &&tria) | |
Triangulation & | operator= (Triangulation< dim, spacedim > &&tria) |
virtual void | set_mesh_smoothing (const MeshSmoothing mesh_smoothing) |
virtual const MeshSmoothing & | get_mesh_smoothing () const |
void | set_boundary (const types::manifold_id number, const Boundary< dim, spacedim > &boundary_object) |
void | set_boundary (const types::manifold_id number) |
void | set_manifold (const types::manifold_id number, const Manifold< dim, spacedim > &manifold_object) |
void | set_manifold (const types::manifold_id number) |
void | set_all_manifold_ids (const types::manifold_id number) |
void | set_all_manifold_ids_on_boundary (const types::manifold_id number) |
void | set_all_manifold_ids_on_boundary (const types::boundary_id b_id, const types::manifold_id number) |
const Boundary< dim, spacedim > & | get_boundary (const types::manifold_id number) const |
const Manifold< dim, spacedim > & | get_manifold (const types::manifold_id number) const |
std::vector< types::boundary_id > | get_boundary_ids () const |
std::vector< types::manifold_id > | get_manifold_ids () const |
virtual void | copy_triangulation (const Triangulation< dim, spacedim > &other_tria) |
virtual void | create_triangulation_compatibility (const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata) |
void | flip_all_direction_flags () |
void | set_all_refine_flags () |
void | refine_global (const unsigned int times=1) |
void | save_refine_flags (std::ostream &out) const |
void | save_refine_flags (std::vector< bool > &v) const |
void | load_refine_flags (std::istream &in) |
void | load_refine_flags (const std::vector< bool > &v) |
void | save_coarsen_flags (std::ostream &out) const |
void | save_coarsen_flags (std::vector< bool > &v) const |
void | load_coarsen_flags (std::istream &out) |
void | load_coarsen_flags (const std::vector< bool > &v) |
bool | get_anisotropic_refinement_flag () const |
void | clear_user_flags () |
void | save_user_flags (std::ostream &out) const |
void | save_user_flags (std::vector< bool > &v) const |
void | load_user_flags (std::istream &in) |
void | load_user_flags (const std::vector< bool > &v) |
void | clear_user_flags_line () |
void | save_user_flags_line (std::ostream &out) const |
void | save_user_flags_line (std::vector< bool > &v) const |
void | load_user_flags_line (std::istream &in) |
void | load_user_flags_line (const std::vector< bool > &v) |
void | clear_user_flags_quad () |
void | save_user_flags_quad (std::ostream &out) const |
void | save_user_flags_quad (std::vector< bool > &v) const |
void | load_user_flags_quad (std::istream &in) |
void | load_user_flags_quad (const std::vector< bool > &v) |
void | clear_user_flags_hex () |
void | save_user_flags_hex (std::ostream &out) const |
void | save_user_flags_hex (std::vector< bool > &v) const |
void | load_user_flags_hex (std::istream &in) |
void | load_user_flags_hex (const std::vector< bool > &v) |
void | clear_user_data () |
void | save_user_indices (std::vector< unsigned int > &v) const |
void | load_user_indices (const std::vector< unsigned int > &v) |
void | save_user_pointers (std::vector< void *> &v) const |
void | load_user_pointers (const std::vector< void *> &v) |
void | save_user_indices_line (std::vector< unsigned int > &v) const |
void | load_user_indices_line (const std::vector< unsigned int > &v) |
void | save_user_indices_quad (std::vector< unsigned int > &v) const |
void | load_user_indices_quad (const std::vector< unsigned int > &v) |
void | save_user_indices_hex (std::vector< unsigned int > &v) const |
void | load_user_indices_hex (const std::vector< unsigned int > &v) |
void | save_user_pointers_line (std::vector< void *> &v) const |
void | load_user_pointers_line (const std::vector< void *> &v) |
void | save_user_pointers_quad (std::vector< void *> &v) const |
void | load_user_pointers_quad (const std::vector< void *> &v) |
void | save_user_pointers_hex (std::vector< void *> &v) const |
void | load_user_pointers_hex (const std::vector< void *> &v) |
cell_iterator | begin (const unsigned int level=0) const |
active_cell_iterator | begin_active (const unsigned int level=0) const |
cell_iterator | end () const |
cell_iterator | end (const unsigned int level) const |
active_cell_iterator | end_active (const unsigned int level) const |
cell_iterator | last () const |
active_cell_iterator | last_active () const |
IteratorRange< cell_iterator > | cell_iterators () const |
IteratorRange< active_cell_iterator > | active_cell_iterators () const |
IteratorRange< cell_iterator > | cell_iterators_on_level (const unsigned int level) const |
IteratorRange< active_cell_iterator > | active_cell_iterators_on_level (const unsigned int level) const |
face_iterator | begin_face () const |
active_face_iterator | begin_active_face () const |
face_iterator | end_face () const |
vertex_iterator | begin_vertex () const |
active_vertex_iterator | begin_active_vertex () const |
vertex_iterator | end_vertex () const |
unsigned int | n_lines () const |
unsigned int | n_lines (const unsigned int level) const |
unsigned int | n_active_lines () const |
unsigned int | n_active_lines (const unsigned int level) const |
unsigned int | n_quads () const |
unsigned int | n_quads (const unsigned int level) const |
unsigned int | n_active_quads () const |
unsigned int | n_active_quads (const unsigned int level) const |
unsigned int | n_hexs () const |
unsigned int | n_hexs (const unsigned int level) const |
unsigned int | n_active_hexs () const |
unsigned int | n_active_hexs (const unsigned int level) const |
unsigned int | n_cells () const |
unsigned int | n_cells (const unsigned int level) const |
unsigned int | n_active_cells () const |
unsigned int | n_active_cells (const unsigned int level) const |
unsigned int | n_faces () const |
unsigned int | n_active_faces () const |
unsigned int | n_levels () const |
unsigned int | n_vertices () const |
const std::vector< Point< spacedim > > & | get_vertices () const |
unsigned int | n_used_vertices () const |
bool | vertex_used (const unsigned int index) const |
const std::vector< bool > & | get_used_vertices () const |
unsigned int | max_adjacent_cells () const |
Triangulation< dim, spacedim > & | get_triangulation () |
const Triangulation< dim, spacedim > & | get_triangulation () const |
unsigned int | n_raw_lines () const |
unsigned int | n_raw_lines (const unsigned int level) const |
unsigned int | n_raw_quads () const |
unsigned int | n_raw_quads (const unsigned int level) const |
unsigned int | n_raw_hexs () const |
unsigned int | n_raw_hexs (const unsigned int level) const |
unsigned int | n_raw_cells (const unsigned int level) const |
unsigned int | n_raw_faces () const |
template<class Archive > | |
void | save (Archive &ar, const unsigned int version) const |
template<class Archive > | |
void | load (Archive &ar, const unsigned int version) |
virtual void | add_periodicity (const std::vector< GridTools::PeriodicFacePair< cell_iterator > > &) |
const std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, std::bitset< 3 > > > & | get_periodic_face_map () const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Member Functions | |
virtual void | update_number_cache () |
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 (::internal::int2type< 2 >) |
void | copy_local_forest_to_triangulation () |
void | attach_mesh_data () |
std::vector< unsigned int > | get_cell_weights () |
void | fill_vertices_with_ghost_neighbors (std::map< unsigned int, std::set<::types::subdomain_id > > &vertices_with_ghost_neighbors) |
void | fill_level_vertices_with_ghost_neighbors (const int level, std::map< unsigned int, std::set<::types::subdomain_id > > &vertices_with_ghost_neighbors) |
std::vector< bool > | mark_locally_active_vertices_on_level (const int level) const |
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 |
bool | refinement_in_progress |
unsigned int | attached_data_size |
unsigned int | n_attached_datas |
unsigned int | n_attached_deserialize |
callback_list_t | attached_data_pack_callbacks |
std::vector< types::global_dof_index > | coarse_cell_to_p4est_tree_permutation |
Additional Inherited Members | |
Static Public Member Functions inherited from Triangulation< dim, spacedim > | |
static ::ExceptionBase & | ExcInvalidLevel (int arg1) |
static ::ExceptionBase & | ExcTriangulationNotEmpty (int arg1, int arg2) |
static ::ExceptionBase & | ExcGridReadError () |
static ::ExceptionBase & | ExcFacesHaveNoLevel () |
static ::ExceptionBase & | ExcEmptyLevel (int arg1) |
static ::ExceptionBase & | ExcNonOrientableTriangulation () |
static ::ExceptionBase & | ExcBoundaryIdNotFound (types::boundary_id arg1) |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, char *arg2, std::string &arg3) |
static ::ExceptionBase & | ExcNoSubscriber (char *arg1, char *arg2) |
Public Attributes inherited from Triangulation< dim, spacedim > | |
Signals | signals |
Static Public Attributes inherited from Triangulation< dim, spacedim > | |
static const StraightBoundary< dim, spacedim > | straight_boundary = StraightBoundary<dim,spacedim>() |
static const unsigned int | dimension = dim |
static const unsigned int | space_dimension = spacedim |
Protected Member Functions inherited from Triangulation< dim, spacedim > | |
void | update_periodic_face_map () |
Static Protected Member Functions inherited from Triangulation< dim, spacedim > | |
static void | write_bool_vector (const unsigned int magic_number1, const std::vector< bool > &v, const unsigned int magic_number2, std::ostream &out) |
static void | read_bool_vector (const unsigned int magic_number1, std::vector< bool > &v, const unsigned int magic_number2, std::istream &in) |
Protected Attributes inherited from parallel::Triangulation< dim, spacedim > | |
MPI_Comm | mpi_communicator |
types::subdomain_id | my_subdomain |
types::subdomain_id | n_subdomains |
Protected Attributes inherited from Triangulation< dim, spacedim > | |
MeshSmoothing | smooth_grid |
This 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::Triangulation::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:
What the call to std_cxx11::bind
does is to produce an object that can be called like a function with no arguments. It does so by taking the address of a function that does, in fact, take an argument but permanently fix this one argument to a reference to the coarse grid triangulation. 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:
Here, like any other member function, set_boundary_ids
implicitly takes a pointer or reference to the object it belongs to as first argument. std::bind
again creates an object that can be called like a global function with no arguments, and this object in turn calls set_boundary_ids
with a pointer to the current object and 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.
Definition at line 46 of file p4est_wrappers.h.
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 behaviour 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. Note: geometric multigrid is still a work in progress. |
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(). |
Triangulation< dim, spacedim >::Triangulation | ( | MPI_Comm | mpi_communicator, |
const typename ::Triangulation< dim, spacedim >::MeshSmoothing | smooth_grid = (::Triangulation<dim,spacedim>::none) , |
||
const Settings | settings = default_setting |
||
) |
Constructor.
mpi_communicator | denotes 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.
|
virtual |
Destructor.
Reimplemented from parallel::Triangulation< dim, spacedim >.
|
virtual |
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 Triangulation< dim, spacedim >.
|
virtual |
Implementation of the same function as in the base class.
|
virtual |
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 >.
|
virtual |
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 >.
|
virtual |
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 be not 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.).
void Triangulation< dim, spacedim >::communicate_locally_moved_vertices | ( | const std::vector< bool > & | vertex_locally_moved | ) |
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(). |
|
virtual |
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 >.
|
virtual |
Return the local memory consumption in bytes.
Reimplemented from parallel::Triangulation< 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 char * | file_basename | ) | const |
unsigned int Triangulation< dim, spacedim >::get_checksum | ( | ) | const |
void Triangulation< dim, spacedim >::save | ( | const char * | filename | ) | const |
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 register_data_attach().
void Triangulation< dim, spacedim >::load | ( | const char * | filename, |
const bool | autopartition = true |
||
) |
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 register_data_attach() can be read in with 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).
unsigned int Triangulation< dim, spacedim >::register_data_attach | ( | const std::size_t | size, |
const std_cxx11::function< void(const cell_iterator &, const CellStatus, void *)> & | pack_callback | ||
) |
Register a function with the current Triangulation object that will be used to attach data to active cells before execute_coarsening_and_refinement(). In execute_coarsening_and_refinement() the Triangulation will call the given function pointer and provide size
bytes to store data. If necessary, this data will be transferred to the new owner of that cell during repartitioning the tree. See notify_ready_to_unpack() on how to retrieve the data.
Callers need to store the return value. It specifies an offset of the position at which data can later be retrieved during a call to notify_ready_to_unpack().
The CellStatus argument in the callback function will tell you if the given cell will be coarsened, refined, or will persist as is (this can be different than the coarsen and refine flags set by you). If it is
When unpacking the data with notify_ready_to_unpack() you can access the children of the cell if the status is CELL_REFINE but not for CELL_COARSEN. As a consequence you need to handle coarsening while packing and refinement during unpacking.
void Triangulation< dim, spacedim >::notify_ready_to_unpack | ( | const unsigned int | offset, |
const std_cxx11::function< void(const cell_iterator &, const CellStatus, const void *)> & | unpack_callback | ||
) |
The supplied callback function is called for each newly locally owned cell and corresponding data saved with register_data_attach(). This function needs to be called after execute_coarsening_and_refinement() with the offset returned by register_data_attach().
The CellStatus will indicate if the cell was refined, coarsened, or persisted unchanged. The cell_iterator will either by 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.
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 |
|
virtual |
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.
|
privatevirtual |
Override the function to update the number cache so we can fill data like level_ghost_owners
.
Reimplemented from parallel::Triangulation< 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 |
Internal function notifying all registered classes to attach their data before repartitioning occurs. Called from execute_coarsening_and_refinement().
|
private |
Internal function notifying all registered slots to provide their weights before repartitioning occurs. Called from execute_coarsening_and_refinement() and repartition().
|
private |
Fills a map that, for each vertex, lists all the processors whose subdomains are adjacent to that vertex. Used by DoFHandler::Policy::ParallelDistributed.
Determine the neighboring subdomains that are adjacent to each vertex. This is achieved via the p4est_iterate/p8est_iterate tool
|
private |
Fills a map that, for each vertex, lists all the processors whose subdomains are adjacent to that vertex on the given level for the multigrid hierarchy. Used by DoFHandler::Policy::ParallelDistributed.
Determine the neighboring subdomains that are adjacent to each vertex on the given multigrid level
|
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
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
A flag that indicates whether refinement of a triangulation is currently in progress. This flag is used to disambiguate whether a call to execute_coarsening_and_triangulation came from the outside or through a recursive call. While the first time we want to take over work to copy things from a refined p4est, the other times we don't want to get in the way as these latter calls to Triangulation::execute_coarsening_and_refinement() are simply there in order to re-create a triangulation that matches the p4est.
|
private |
number of bytes that get attached to the Triangulation through register_data_attach() for example SolutionTransfer.
|
private |
number of functions that get attached to the Triangulation through register_data_attach() for example SolutionTransfer.
|
private |
|
private |
List of callback functions registered by register_data_attach() that are going to be called for packing data.
|
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.