Reference documentation for deal.II version 9.0.0
|
#include <deal.II/grid/persistent_tria.h>
Public Member Functions | |
PersistentTriangulation (const Triangulation< dim, spacedim > &coarse_grid) | |
PersistentTriangulation (const PersistentTriangulation< dim, spacedim > &old_tria) | |
virtual | ~PersistentTriangulation ()=default |
virtual void | execute_coarsening_and_refinement () |
void | restore () |
void | restore (const unsigned int step_no) |
unsigned int | n_refinement_steps () const |
virtual void | copy_triangulation (const Triangulation< dim, spacedim > &tria) |
virtual void | create_triangulation (const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata) |
virtual void | create_triangulation_compatibility (const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata) |
virtual void | write_flags (std::ostream &out) const |
virtual void | read_flags (std::istream &in) |
virtual void | clear_flags () |
virtual std::size_t | memory_consumption () 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 > &)=delete | |
Triangulation (Triangulation< dim, spacedim > &&tria) noexcept | |
Triangulation & | operator= (Triangulation< dim, spacedim > &&tria) noexcept |
virtual | ~Triangulation () |
virtual void | clear () |
virtual void | set_mesh_smoothing (const MeshSmoothing mesh_smoothing) |
virtual const MeshSmoothing & | get_mesh_smoothing () const |
void | set_boundary (const types::manifold_id number, const Boundary< dim, spacedim > &boundary_object) |
void | set_boundary (const types::manifold_id number) |
void | set_manifold (const types::manifold_id number, const Manifold< dim, spacedim > &manifold_object) |
void | set_manifold (const types::manifold_id number) |
void | 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 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 |
void | flip_all_direction_flags () |
void | set_all_refine_flags () |
void | refine_global (const unsigned int times=1) |
virtual bool | prepare_coarsening_and_refinement () |
void | save_refine_flags (std::ostream &out) const |
void | save_refine_flags (std::vector< bool > &v) const |
void | load_refine_flags (std::istream &in) |
void | load_refine_flags (const std::vector< bool > &v) |
void | save_coarsen_flags (std::ostream &out) const |
void | save_coarsen_flags (std::vector< bool > &v) const |
void | load_coarsen_flags (std::istream &out) |
void | load_coarsen_flags (const std::vector< bool > &v) |
bool | get_anisotropic_refinement_flag () const |
void | clear_user_flags () |
void | save_user_flags (std::ostream &out) const |
void | save_user_flags (std::vector< bool > &v) const |
void | load_user_flags (std::istream &in) |
void | load_user_flags (const std::vector< bool > &v) |
void | clear_user_flags_line () |
void | save_user_flags_line (std::ostream &out) const |
void | save_user_flags_line (std::vector< bool > &v) const |
void | load_user_flags_line (std::istream &in) |
void | load_user_flags_line (const std::vector< bool > &v) |
void | clear_user_flags_quad () |
void | save_user_flags_quad (std::ostream &out) const |
void | save_user_flags_quad (std::vector< bool > &v) const |
void | load_user_flags_quad (std::istream &in) |
void | load_user_flags_quad (const std::vector< bool > &v) |
void | clear_user_flags_hex () |
void | save_user_flags_hex (std::ostream &out) const |
void | save_user_flags_hex (std::vector< bool > &v) const |
void | load_user_flags_hex (std::istream &in) |
void | load_user_flags_hex (const std::vector< bool > &v) |
void | clear_user_data () |
void | save_user_indices (std::vector< unsigned int > &v) const |
void | load_user_indices (const std::vector< unsigned int > &v) |
void | save_user_pointers (std::vector< void *> &v) const |
void | load_user_pointers (const std::vector< void *> &v) |
void | save_user_indices_line (std::vector< unsigned int > &v) const |
void | load_user_indices_line (const std::vector< unsigned int > &v) |
void | save_user_indices_quad (std::vector< unsigned int > &v) const |
void | load_user_indices_quad (const std::vector< unsigned int > &v) |
void | save_user_indices_hex (std::vector< unsigned int > &v) const |
void | load_user_indices_hex (const std::vector< unsigned int > &v) |
void | save_user_pointers_line (std::vector< void *> &v) const |
void | load_user_pointers_line (const std::vector< void *> &v) |
void | save_user_pointers_quad (std::vector< void *> &v) const |
void | load_user_pointers_quad (const std::vector< void *> &v) |
void | save_user_pointers_hex (std::vector< void *> &v) const |
void | load_user_pointers_hex (const std::vector< void *> &v) |
cell_iterator | begin (const unsigned int level=0) const |
active_cell_iterator | begin_active (const unsigned int level=0) const |
cell_iterator | end () const |
cell_iterator | end (const unsigned int level) const |
active_cell_iterator | end_active (const unsigned int level) const |
cell_iterator | last () const |
active_cell_iterator | last_active () const |
IteratorRange< cell_iterator > | cell_iterators () const |
IteratorRange< active_cell_iterator > | active_cell_iterators () const |
IteratorRange< cell_iterator > | cell_iterators_on_level (const unsigned int level) const |
IteratorRange< active_cell_iterator > | active_cell_iterators_on_level (const unsigned int level) const |
face_iterator | begin_face () const |
active_face_iterator | begin_active_face () const |
face_iterator | end_face () const |
vertex_iterator | begin_vertex () const |
active_vertex_iterator | begin_active_vertex () const |
vertex_iterator | end_vertex () const |
unsigned int | n_lines () const |
unsigned int | n_lines (const unsigned int level) const |
unsigned int | n_active_lines () const |
unsigned int | n_active_lines (const unsigned int level) const |
unsigned int | n_quads () const |
unsigned int | n_quads (const unsigned int level) const |
unsigned int | n_active_quads () const |
unsigned int | n_active_quads (const unsigned int level) const |
unsigned int | n_hexs () const |
unsigned int | n_hexs (const unsigned int level) const |
unsigned int | n_active_hexs () const |
unsigned int | n_active_hexs (const unsigned int level) const |
unsigned int | n_cells () const |
unsigned int | n_cells (const unsigned int level) const |
unsigned int | n_active_cells () const |
virtual types::global_dof_index | n_global_active_cells () const |
unsigned int | n_active_cells (const unsigned int level) const |
unsigned int | n_faces () const |
unsigned int | n_active_faces () const |
unsigned int | n_levels () const |
virtual 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 |
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 &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Member Functions | |
static ::ExceptionBase & | ExcTriaNotEmpty () |
static ::ExceptionBase & | ExcFlagsNotCleared () |
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, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Static Public Attributes | |
static const unsigned int | dimension = dim |
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 |
Private Attributes | |
SmartPointer< const Triangulation< dim, spacedim >, PersistentTriangulation< dim, spacedim > > | coarse_grid |
std::vector< std::vector< bool > > | refine_flags |
std::vector< std::vector< bool > > | coarsen_flags |
This class handles the history of a triangulation and can rebuild it after it was deleted some time before. Its main purpose is support for time- dependent problems where one frequently deletes a triangulation due to memory pressure and later wants to rebuild it; this class has all the information to rebuild it exactly as it was before including the mapping of cell numbers to the geometrical cells.
Basically, this is a drop-in replacement for the triangulation. Since it is derived from the Triangulation class, it shares all the functionality, but it overrides some virtual functions and adds some functions, too. The main change to the base class is that it overrides the execute_coarsening_and_refinement
function, where the new version first stores all refinement and coarsening flags and only then calls the respective function of the base class. The stored flags may later be used to restore the grid just as it was before. Some other functions have been extended slightly as well, see their documentation for more information.
We note that since the triangulation is created in exactly the same state as it was before, other objects working on it should result in the same state as well. This holds in particular for the DoFHandler object, which will assign the same degrees of freedom to the original cells and the ones after reconstruction of the triangulation. You can therefore safely use data vectors computed on the original grid on the reconstructed grid as well.
You can use objects of this class almost in the same way as objects of the Triangulation class. One of the few differences is that you can only construct such an object by giving a coarse grid to the constructor. The coarse grid will be used to base the triangulation on, and therefore the lifetime of the coarse grid has to be longer than the lifetime of the object of this class.
Basically, usage looks like this:
Note that initially, the PersistentTriangulation object does not constitute a triangulation; it only becomes one after restore
is first called. Note also that the execute_coarsening_and_refinement
stores all necessary flags for later reconstruction using the restore
function. Triangulation::clear() resets the underlying triangulation to a virgin state, but does not affect the stored refinement flags needed for later reconstruction and does also not touch the coarse grid which is used within restore().
Definition at line 106 of file persistent_tria.h.
PersistentTriangulation< dim, spacedim >::PersistentTriangulation | ( | const Triangulation< dim, spacedim > & | coarse_grid | ) |
Build up the triangulation from the coarse grid in future. Copy smoothing flags, etc from that grid as well. Note that the initial state of the triangulation is empty, until restore_grid
is called for the first time.
The coarse grid must persist until the end of this object, since it will be used upon reconstruction of the grid.
Definition at line 35 of file persistent_tria.cc.
PersistentTriangulation< dim, spacedim >::PersistentTriangulation | ( | const PersistentTriangulation< dim, spacedim > & | old_tria | ) |
Copy constructor. This operation is only allowed, if the triangulation underlying the object to be copied is presently empty. Refinement flags as well as the pointer to the coarse grid are copied, however.
Definition at line 44 of file persistent_tria.cc.
|
virtualdefault |
Destructor.
|
virtual |
Overloaded version of the same function in the base class which stores the refinement and coarsening flags for later reconstruction of the triangulation and after that calls the respective function of the base class.
Reimplemented from Triangulation< dim, spacedim >.
Definition at line 61 of file persistent_tria.cc.
void PersistentTriangulation< dim, spacedim >::restore | ( | ) |
Restore the grid according to the saved data. For this, the coarse grid is copied and the grid is stepwise rebuilt using the saved flags.
Note that this function will result in an error if the underlying triangulation is not empty, i.e. it will only succeed if this object is newly created or the clear()
function of the base class was called on it before.
Repeatedly calls the restore(unsigned int)
function in a loop over all refinement steps.
Definition at line 80 of file persistent_tria.cc.
void PersistentTriangulation< dim, spacedim >::restore | ( | const unsigned int | step_no | ) |
Differential restore. Performs the step_noth
local refinement and coarsening step. Step 0 stands for the copying of the coarse grid.
This function will only succeed if the triangulation is in just the state it were if restore would have been called from step=0...step_no-1
before.
Definition at line 92 of file persistent_tria.cc.
unsigned int PersistentTriangulation< dim, spacedim >::n_refinement_steps | ( | ) | const |
Return the number of refinement and coarsening steps. This is given by the size of the refine_flags
vector.
Definition at line 119 of file persistent_tria.cc.
|
virtual |
Overload this function to use tria
as a new coarse grid. The present triangulation and all refinement and coarsening flags storing its history are deleted, and the state of the underlying triangulation is reset to be empty, until restore_grid
is called the next time.
The coarse grid must persist until the end of this object, since it will be used upon reconstruction of the grid.
Reimplemented from Triangulation< dim, spacedim >.
Definition at line 128 of file persistent_tria.cc.
|
virtual |
Throw an error, since this function is not useful in the context of this class.
Reimplemented from Triangulation< dim, spacedim >.
Definition at line 140 of file persistent_tria.cc.
|
virtual |
An overload of the respective function of the base class.
Throw an error, since this function is not useful in the context of this class.
Reimplemented from Triangulation< dim, spacedim >.
Definition at line 151 of file persistent_tria.cc.
|
virtual |
Write all refine and coarsen flags to the ostream out
.
Definition at line 163 of file persistent_tria.cc.
|
virtual |
Reads all refine and coarsen flags that previously were written by write_flags(...)
. This is especially useful for rebuilding the triangulation after the end or breakdown of a program and its restart.
Definition at line 188 of file persistent_tria.cc.
|
virtual |
Clear all flags. Retains the same coarse grid.
Definition at line 222 of file persistent_tria.cc.
|
virtual |
Determine an estimate for the memory consumption (in bytes) of this object.
Reimplemented from Triangulation< dim, spacedim >.
Definition at line 232 of file persistent_tria.cc.
|
static |
Make the dimension available in function templates.
Definition at line 112 of file persistent_tria.h.
|
private |
This grid shall be used as coarse grid.
Definition at line 242 of file persistent_tria.h.
|
private |
Vectors holding the refinement and coarsening flags of the different sweeps on this time level. The vectors therefore hold the history of the grid.
Definition at line 249 of file persistent_tria.h.
|
private |
Definition at line 254 of file persistent_tria.h.