Reference documentation for deal.II version 9.1.1
|
#include <deal.II/grid/tria.h>
Classes | |
struct | CellWeightSum |
struct | DistortedCellList |
struct | Signals |
Public Types | |
enum | MeshSmoothing { none = 0x0, limit_level_difference_at_vertices = 0x1, eliminate_unrefined_islands = 0x2, patch_level_1 = 0x4, coarsest_level_1 = 0x8, allow_anisotropic_smoothing = 0x10, eliminate_refined_inner_islands = 0x100, eliminate_refined_boundary_islands = 0x200, do_not_produce_unrefined_islands = 0x400, smoothing_on_refinement, smoothing_on_coarsening, maximum_smoothing = 0xffff ^ allow_anisotropic_smoothing } |
using | cell_iterator = TriaIterator< CellAccessor< dim, spacedim > > |
using | active_cell_iterator = TriaActiveIterator< CellAccessor< dim, spacedim > > |
using | face_iterator = TriaIterator< TriaAccessor< dim - 1, dim, spacedim > > |
using | active_face_iterator = TriaActiveIterator< TriaAccessor< dim - 1, dim, spacedim > > |
using | vertex_iterator = TriaIterator<::TriaAccessor< 0, dim, spacedim > > |
using | active_vertex_iterator = TriaActiveIterator<::TriaAccessor< 0, dim, spacedim > > |
using | line_iterator = typename IteratorSelector::line_iterator |
using | active_line_iterator = typename IteratorSelector::active_line_iterator |
using | quad_iterator = typename IteratorSelector::quad_iterator |
using | active_quad_iterator = typename IteratorSelector::active_quad_iterator |
using | hex_iterator = typename IteratorSelector::hex_iterator |
using | active_hex_iterator = typename IteratorSelector::active_hex_iterator |
Public Member Functions | |
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 () override |
virtual void | clear () |
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 | 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 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 (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) |
void | flip_all_direction_flags () |
Mesh refinement | |
void | set_all_refine_flags () |
void | refine_global (const unsigned int times=1) |
virtual void | execute_coarsening_and_refinement () |
virtual bool | prepare_coarsening_and_refinement () |
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 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 |
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 |
Internal information about the number of objects | |
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 |
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 (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) |
Static Public Attributes | |
static const unsigned int | dimension = dim |
static const unsigned int | space_dimension = spacedim |
Private Types | |
using | IteratorSelector = ::internal::TriangulationImplementation::Iterators< dim, spacedim > |
Private Member Functions | |
Line iterator functions for internal use | |
raw_line_iterator | begin_raw_line (const unsigned int level=0) const |
line_iterator | begin_line (const unsigned int level=0) const |
active_line_iterator | begin_active_line (const unsigned int level=0) const |
line_iterator | end_line () const |
Quad iterator functions for internal use | |
raw_quad_iterator | begin_raw_quad (const unsigned int level=0) const |
quad_iterator | begin_quad (const unsigned int level=0) const |
active_quad_iterator | begin_active_quad (const unsigned int level=0) const |
quad_iterator | end_quad () const |
Keeping up with what happens to a triangulation | |
enum | CellStatus { CELL_PERSIST, CELL_REFINE, CELL_COARSEN, CELL_INVALID } |
Signals | signals |
Exceptions | |
MeshSmoothing | smooth_grid |
std::vector< GridTools::PeriodicFacePair< cell_iterator > > | periodic_face_pairs_level_0 |
std::map< std::pair< cell_iterator, unsigned int >, std::pair< std::pair< cell_iterator, unsigned int >, std::bitset< 3 > > > | 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 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) |
void | update_periodic_face_map () |
Cell iterator functions for internal use | |
using | raw_cell_iterator = TriaRawIterator< CellAccessor< dim, spacedim > > |
using | raw_face_iterator = TriaRawIterator< TriaAccessor< dim - 1, dim, spacedim > > |
using | raw_vertex_iterator = TriaRawIterator<::TriaAccessor< 0, dim, spacedim > > |
using | raw_line_iterator = typename IteratorSelector::raw_line_iterator |
using | raw_quad_iterator = typename IteratorSelector::raw_quad_iterator |
using | raw_hex_iterator = typename IteratorSelector::raw_hex_iterator |
raw_cell_iterator | begin_raw (const unsigned int level=0) const |
raw_cell_iterator | end_raw (const unsigned int level) const |
Hex iterator functions for internal use | |
std::vector< std::unique_ptr< ::internal::TriangulationImplementation::TriaLevel< dim > > > | levels |
std::unique_ptr<::internal::TriangulationImplementation::TriaFaces< dim > > | faces |
std::vector< Point< spacedim > > | vertices |
std::vector< bool > | vertices_used |
std::map< types::manifold_id, std::unique_ptr< const Manifold< dim, spacedim > > > | manifold |
bool | anisotropic_refinement |
const bool | check_for_distorted_cells |
::internal::TriangulationImplementation::NumberCache< dim > | number_cache |
std::unique_ptr< std::map< unsigned int, types::boundary_id > > | vertex_to_boundary_id_map_1d |
std::unique_ptr< std::map< unsigned int, types::manifold_id > > | vertex_to_manifold_id_map_1d |
template<int , int , int > | |
class | TriaAccessorBase |
template<int , int , int > | |
class | TriaAccessor |
class | TriaAccessor< 0, 1, spacedim > |
class | CellAccessor< dim, spacedim > |
struct | ::internal::TriaAccessorImplementation::Implementation |
class | hp::DoFHandler< dim, spacedim > |
struct | ::internal::TriangulationImplementation::Implementation |
template<typename > | |
class | ::internal::TriangulationImplementation::TriaObjects |
raw_hex_iterator | begin_raw_hex (const unsigned int level=0) const |
hex_iterator | begin_hex (const unsigned int level=0) const |
active_hex_iterator | begin_active_hex (const unsigned int level=0) const |
hex_iterator | end_hex () const |
void | clear_despite_subscriptions () |
void | reset_active_cell_indices () |
DistortedCellList | execute_refinement () |
void | execute_coarsening () |
void | fix_coarsen_flags () |
Additional Inherited Members | |
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) |
Triangulations denote a hierarchy of levels of elements which together form a dim
-dimensional manifold in spacedim
spatial dimensions (if spacedim is not specified it takes the default value spacedim=dim
).
Thus, for example, an object of type Triangulation<1,1>
(or simply Triangulation<1>
since spacedim==dim
by default) is used to represent and handle the usual one-dimensional triangulation used in the finite element method (so, segments on a straight line). On the other hand, objects such as Triangulation<1,2>
or Triangulation<2,3>
(that are associated with curves in 2D or surfaces in 3D) are the ones one wants to use in the boundary element method.
This class is written to be as independent of the dimension as possible (thus the complex construction of the internal::TriangulationImplementation::TriaLevel classes) to allow code-sharing, to allow reducing the need to mirror changes in the code for one dimension to the code for other dimensions. Nonetheless, some of the functions are dependent of the dimension and there only exist specialized versions for distinct dimensions.
This class satisfies the MeshType concept requirements.
The actual data structure of a Triangulation object is rather complex and quite inconvenient if one attempted to operate on it directly, since data is spread over quite a lot of arrays and other places. However, there are ways powerful enough to work on these data structures without knowing their exact relations. deal.II uses class local alias (see below) to make things as easy and dimension independent as possible.
The Triangulation class provides iterators which enable looping over all cells without knowing the exact representation used to describe them. For more information see the documentation of TriaIterator
. Their names are alias imported from the Iterators class (thus making them local types to this class) and are as follows:
cell_iterator
: loop over all cells used in the Triangulation active_cell_iterator
: loop over all active cells For dim==1
, these iterators are mapped as follows:
while for dim==2
we have the additional face iterator:
By using the cell iterators, you can write code independent of the spatial dimension. The same applies for substructure iterators, where a substructure is defined as a face of a cell. The face of a cell is a vertex in 1D and a line in 2D; however, vertices are handled in a different way and therefore lines have no faces.
The Triangulation class offers functions like begin_active() which gives you an iterator to the first active cell. There are quite a lot of functions returning iterators. Take a look at the class doc to get an overview.
Usage of these iterators is similar to usage of standard container iterators. Some examples taken from the Triangulation source code follow (notice that in the last two examples the template parameter spacedim
has been omitted, so it takes the default value dim
).
Counting the number of cells on a specific level
Another way, which uses std::distance
, would be to write
Usage of a Triangulation is mainly done through the use of iterators. An example probably shows best how to use it:
There are several possibilities to create a triangulation:
The most common domains, such as hypercubes (i.e. lines, squares, cubes, etc), hyper-balls (circles, balls, ...) and some other, more weird domains such as the L-shape region and higher dimensional generalizations and others, are provided by the GridGenerator class which takes a triangulation and fills it by a division of the required domain.
Reading in a triangulation: By using an object of the GridIn class, you can read in fairly general triangulations. See there for more information. The mentioned class uses the interface described directly below to transfer the data into the triangulation.
Explicitly creating a triangulation: you can create a triangulation by providing a list of vertices and a list of cells. Each such cell consists of a vector storing the indices of the vertices of this cell in the vertex list. To see how this works, you can take a look at the GridIn<dim>::read_* functions. The appropriate function to be called is create_triangulation().
Creating the hierarchical information needed for this library from cells storing only vertex information can be quite a complex task. For example in 2D, we have to create lines between vertices (but only once, though there are two cells which link these two vertices) and we have to create neighborhood information. Grids being read in should therefore not be too large, reading refined grids would be inefficient (although there is technically no problem in reading grids with several 10.000 or 100.000 cells; the library can handle this without much problems). Apart from the performance aspect, refined grids do not lend too well to multigrid algorithms, since solving on the coarsest level is expensive. It is wiser in any case to read in a grid as coarse as possible and then do the needed refinement steps.
It is your duty to guarantee that cells have the correct orientation. To guarantee this, in the input vector keeping the cell list, the vertex indices for each cell have to be in a defined order, see the documentation of GeometryInfo<dim>. In one dimension, the first vertex index must refer to that vertex with the lower coordinate value. In 2D and 3D, the corresponding conditions are not easy to verify and no full attempt to do so is made. If you violate this condition, you may end up with matrix entries having the wrong sign (clockwise vertex numbering, which results in a negative area element) of with wrong matrix elements (twisted quadrilaterals, i.e. two vertices interchanged; this results in a wrong area element).
There are more subtle conditions which must be imposed upon the vertex numbering within cells. They do not only hold for the data read from an UCD or any other input file, but also for the data passed to create_triangulation(). See the documentation for the GridIn class for more details on this, and above all to the GridReordering class that explains many of the problems and an algorithm to reorder cells such that they satisfy the conditions outlined above.
Copying a triangulation: when computing on time dependent meshes or when using adaptive refinement, you will often want to create a new triangulation to be the same as another one. This is facilitated by the copy_triangulation
function.
It is guaranteed that vertex, line or cell numbers in the two triangulations are the same and that two iterators walking on the two triangulations visit matching cells if they are incremented in parallel. It may be conceivable to implement a clean-up in the copy operation, which eliminates holes of unused memory, re-joins scattered data and so on. In principle this would be a useful operation but guaranteeing some parallelism in the two triangulations seems more important since usually data will have to be transferred between the grids.
Finally, there is a special function for folks who like bad grids: distort_random(). It moves all the vertices in the grid a bit around by a random value, leaving behind a distorted mesh. Note that you should apply this function to the final mesh, since refinement smoothes the mesh a bit.
The function will make sure that vertices on restricted faces (hanging nodes) will end up in the correct place, i.e. in the middle of the two other vertices of the mother line, and the analogue in higher space dimensions (vertices on the boundary are not corrected, so don't distort boundary vertices in more than two space dimension, i.e. in dimensions where boundary vertices can be hanging nodes). Applying the algorithm has another drawback related to the placement of cells, however: the children of a cell will not occupy the same region of the domain as the mother cell does. While this is the usual behavior with cells at the boundary, here you may get into trouble when using multigrid algorithms or when transferring solutions from coarse to fine grids and back. In general, the use of this function is only safe if you only use the most refined level of the triangulation for computations.
Refinement of a triangulation may be done through several ways. The most low-level way is directly through iterators: let i
be an iterator to an active cell (i.e. the cell pointed to has no children), then the function call i->set_refine_flag()
marks the respective cell for refinement. Marking non-active cells results in an error.
After all the cells you wanted to mark for refinement, call execute_coarsening_and_refinement() to actually perform the refinement. This function itself first calls the prepare_coarsening_and_refinement
function to regularize the resulting triangulation: since a face between two adjacent cells may only be subdivided once (i.e. the levels of two adjacent cells may differ by one at most; it is not possible to have a cell refined twice while the neighboring one is not refined), some additional cells are flagged for refinement to smooth the grid. This enlarges the number of resulting cells but makes the grid more regular, thus leading to better approximation properties and, above all, making the handling of data structures and algorithms much easier. To be honest, this is mostly an algorithmic step than one needed by the finite element method.
To coarsen a grid, the same way as above is possible by using i->set_coarsen_flag
and calling execute_coarsening_and_refinement().
The reason for first coarsening, then refining is that the refinement usually adds some additional cells to keep the triangulation regular and thus satisfies all refinement requests, while the coarsening does not delete cells not requested for; therefore the refinement will often revert some effects of coarsening while the opposite is not true. The stated order of coarsening before refinement will thus normally lead to a result closer to the intended one.
Marking cells for refinement 'by hand' through iterators is one way to produce a new grid, especially if you know what kind of grid you are looking for, e.g. if you want to have a grid successively refined towards the boundary or always at the center (see the example programs, they do exactly these things). There are more advanced functions, however, which are more suitable for automatic generation of hierarchical grids in the context of a posteriori error estimation and adaptive finite elements. These functions can be found in the GridRefinement class.
Some degradation of approximation properties has been observed for grids which are too unstructured. Therefore, prepare_coarsening_and_refinement() which is automatically called by execute_coarsening_and_refinement() can do some smoothing of the triangulation. Note that mesh smoothing is only done for two or more space dimensions, no smoothing is available at present for one spatial dimension. In the following, let execute_*
stand for execute_coarsening_and_refinement().
For the purpose of smoothing, the Triangulation constructor takes an argument specifying whether a smoothing step shall be performed on the grid each time execute_*
is called. The default is that such a step not be done, since this results in additional cells being produced, which may not be necessary in all cases. If switched on, calling execute_*
results in flagging additional cells for refinement to avoid vertices as the ones mentioned. The algorithms for both regularization and smoothing of triangulations are described below in the section on technical issues. The reason why this parameter must be given to the constructor rather than to execute_*
is that it would result in algorithmic problems if you called execute_*
once without and once with smoothing, since then in some refinement steps would need to be refined twice.
The parameter taken by the constructor is an integer which may be composed bitwise by the constants defined in the enum MeshSmoothing (see there for the possibilities).
Each cell, face or edge stores information denoting the material or the part of the boundary that an object belongs to. The material id of a cell is typically used to identify which cells belong to a particular part of the domain, e.g., when you have different materials (steel, concrete, wood) that are all part of the same domain. One would then usually query the material id associated with a cell during assembly of the bilinear form, and use it to determine (e.g., by table lookup, or a sequence of if-else statements) what the correct material coefficients would be for that cell. See also this glossary entry.
This material_id may be set upon construction of a triangulation (through the CellData data structure), or later through use of cell iterators. For a typical use of this functionality, see the step-28 tutorial program. The functions of the GridGenerator namespace typically set the material ID of all cells to zero. When reading a triangulation through the GridIn class, different input file formats have different conventions, but typically either explicitly specify the material id, or if they don't, then GridIn simply sets them to zero. Because the material of a cell is intended to pertain to a particular region of the domain, material ids are inherited by child cells from their parent upon mesh refinement.
Boundary indicators on lower dimensional objects (these have no material id) indicate the number of a boundary component. The weak formulation of the partial differential equation may have different boundary conditions on different parts of the boundary. The boundary indicator can be used in creating the matrix or the right hand side vector to indicate these different parts of the model (this use is like the material id of cells). Boundary indicators may be in the range from zero to numbers::internal_face_boundary_id-1. The value numbers::internal_face_boundary_id is reserved to denote interior lines (in 2D) and interior lines and quads (in 3D), which do not have a boundary indicator. This way, a program can easily determine, whether such an object is at the boundary or not. Material indicators may be in the range from zero to numbers::invalid_material_id-1.
Lines in two dimensions and quads in three dimensions inherit their boundary indicator to their children upon refinement. You should therefore make sure that if you have different boundary parts, the different parts are separated by a vertex (in 2D) or a line (in 3D) such that each boundary line or quad has a unique boundary indicator.
By default (unless otherwise specified during creation of a triangulation), all parts of the boundary have boundary indicator zero. As a historical wart, this isn't true for 1d meshes, however: For these, leftmost vertices have boundary indicator zero while rightmost vertices have boundary indicator one. In either case, the boundary indicator of a face can be changed using a call of the kind cell->face(1)->set_boundary_id(42);
.
It is possible to reconstruct a grid from its refinement history, which can be stored and loaded through the save_refine_flags
and load_refine_flags
functions. Normally, the code will look like this:
If you want to re-create the grid from the stored information, you write:
The same scheme is employed for coarsening and the coarsening flags.
You may write other information to the output file between different sets of refinement information, as long as you read it upon re-creation of the grid. You should make sure that the other information in the new triangulation which is to be created from the saved flags, matches that of the old triangulation, for example the smoothing level; if not, the cells actually created from the flags may be other ones, since smoothing adds additional cells, but their number may be depending on the smoothing level.
There actually are two sets of save_*_flags
and load_*_flags
functions. One takes a stream as argument and reads/writes the information from/to the stream, thus enabling storing flags to files. The other set takes an argument of type vector<bool>
. This enables the user to temporarily store some flags, e.g. if another function needs them, and restore them afterwards.
A triangulation offers one bit per line, quad, etc for user flags. This field can be accessed as all other data using iterators. Normally, this user flag is used if an algorithm walks over all cells and needs information whether another cell, e.g. a neighbor, has already been processed. See the glossary for more information.
There is another set of user data, which can be either an unsigned int
or a void *
, for each line, quad, etc. You can access these through the functions listed under User data
in the accessor classes. Again, see the glossary for more information.
The value of these user indices or pointers is nullptr
by default. Note that the pointers are not inherited to children upon refinement. Still, after a remeshing they are available on all cells, where they were set on the previous mesh.
The usual warning about the missing type safety of void
pointers are obviously in place here; responsibility for correctness of types etc lies entirely with the user of the pointer.
deal.II implements all geometries (curved and otherwise) with classes inheriting from Manifold; see the documentation of Manifold, step-49, or the Manifold description for triangulations module for examples and a complete description of the algorithms. By default, all cells in a Triangulation have a flat geometry, meaning that all lines in the Triangulation are assumed to be straight. If a cell has a manifold_id that is not equal to numbers::flat_manifold_id then the Triangulation uses the associated Manifold object for computations on that cell (e.g., cell refinement). Here is a quick example, taken from the implementation of GridGenerator::hyper_ball(), that sets up a polar grid:
This will set up a grid where the boundary lines will be refined by performing calculations in polar coordinates. When the mesh is refined the cells adjacent to the boundary will use this new line midpoint (as well as the other three midpoints and original cell vertices) to calculate the cell midpoint with a transfinite interpolation: this propagates the curved boundary into the interior in a smooth way. It is possible to generate a better grid (which interpolates across all cells between two different Manifold descriptions, instead of just going one cell at a time) by using TransfiniteInterpolationManifold; see the documentation of that class for more information.
You should take note of one caveat: if you have concave boundaries, you must make sure that a new boundary vertex does not lie too much inside the cell which is to be refined. The reason is that the center vertex is placed at the point which is a weighted average of the vertices of the original cell, new face midpoints, and (in 3D) new line midpoints. Therefore if your new boundary vertex is too near the center of the old quadrilateral or hexahedron, the distance to the midpoint vertex will become too small, thus generating distorted cells. This issue is discussed extensively in distorted cells.
There are cases where one object would like to know whenever a triangulation is being refined, copied, or modified in a number of other ways. This could of course be achieved if, in your user code, you tell every such object whenever you are about to refine the triangulation, but this will get tedious and is error prone. The Triangulation class implements a more elegant way to achieve this: signals.
In essence, a signal is an object (a member of the Triangulation class) that another object can connect to. A connection is in essence that the connecting object passes a function object taking a certain number and kind of arguments. Whenever the owner of the signal wants to indicate a certain kind of event, it 'triggers' the signal, which in turn means that all connections of the signal are triggered: in other word, the function objects are executed and can take the action that is necessary.
As a simple example, the following code will print something to the output every time the triangulation has just been refined:
This code will produce output twice, once for each refinement cycle.
A more interesting application would be the following, akin to what the FEValues class does. This class stores a pointer to a triangulation and also an iterator to the cell last handled (so that it can compare the current cell with the previous one and, for example, decide that there is no need to re-compute the Jacobian matrix if the new cell is a simple translation of the previous one). However, whenever the triangulation is modified, the iterator to the previously handled cell needs to be invalidated since it now no longer points to any useful cell (or, at the very least, points to something that may not necessarily resemble the cells previously handled). The code would look something like this (the real code has some more error checking and has to handle the case that subsequent cells might actually belong to different triangulation, but that is of no concern to us here):
Here, whenever the triangulation is refined, it triggers the post- refinement signal which calls the function object attached to it. This function object is the member function FEValues<dim>::invalidate_previous_cell
where we have bound the single argument (the this
pointer of a member function that otherwise takes no arguments) to the this
pointer of the FEValues object. Note how here there is no need for the code that owns the triangulation and the FEValues object to inform the latter if the former is refined. (In practice, the function would want to connect to some of the other signals that the triangulation offers as well, in particular to creation and deletion signals.)
The Triangulation class has a variety of signals that indicate different actions by which the triangulation can modify itself and potentially require follow-up action elsewhere. Please refer to Triangulation::Signals for details.
Like many other classes in deal.II, the Triangulation class can stream its contents to an archive using BOOST's serialization facilities. The data so stored can later be retrieved again from the archive to restore the contents of this object. This facility is frequently used to save the state of a program to disk for possible later resurrection, often in the context of checkpoint/restart strategies for long running computations or on computers that aren't very reliable (e.g. on very large clusters where individual nodes occasionally fail and then bring down an entire MPI job).
For technical reasons, writing and restoring a Triangulation object is not trivial. The primary reason is that unlike many other objects, triangulations rely on many other objects to which they store pointers or with which they interface; for example, triangulations store pointers to objects describing boundaries and manifolds, and they have signals that store pointers to other objects so they can be notified of changes in the triangulation (see the section on signals in this introduction). Since these objects are owned by the user space (for example the user can create a custom manifold object), they may not be serializable. So in cases like this, boost::serialize can store a reference to an object instead of the pointer, but the reference will never be satisfied at write time because the object pointed to is not serialized. Clearly, at load time, boost::serialize will not know where to let the pointer point to because it never gets to re-create the object originally pointed to.
For these reasons, saving a triangulation to an archive does not store all information, but only certain parts. More specifically, the information that is stored is everything that defines the mesh such as vertex locations, vertex indices, how vertices are connected to cells, boundary indicators, subdomain ids, material ids, etc. On the other hand, the following information is not stored:
In a sense, this approach to serialization means that re-loading a triangulation is more akin to calling the Triangulation::create_triangulation function and filling it with some additional content, as that function also does not touch the signals and boundary objects that belong to this triangulation. In keeping with this analogy, the Triangulation::load function also triggers the same kinds of signal as Triangulation::create_triangulation.
We chose an inductive point of view: since upon creation of the triangulation all cells are on the same level, all regularity assumptions regarding the maximum difference in level of cells sharing a common face, edge or vertex hold. Since we use the regularization and smoothing in each step of the mesh history, when coming to the point of refining it further the assumptions also hold.
The regularization and smoothing is done in the prepare_coarsening_and_refinement
function, which is called by execute_coarsening_and_refinement
at the very beginning. It decides which additional cells to flag for refinement by looking at the old grid and the refinement flags for each cell.
Regularization: The algorithm walks over all cells checking whether the present cell is flagged for refinement and a neighbor of the present cell is refined once less than the present one. If so, flag the neighbor for refinement. Because of the induction above, there may be no neighbor with level two less than the present one.
The neighbor thus flagged for refinement may induce more cells which need to be refined. However, such cells which need additional refinement always are on one level lower than the present one, so we can get away with only one sweep over all cells if we do the loop in the reverse way, starting with those on the highest level. This way, we may flag additional cells on lower levels, but if these induce more refinement needed, this is performed later on when we visit them in out backward running loop.
limit_level_difference_at_vertices:
First a list is set up which stores for each vertex the highest level one of the adjacent cells belongs to. Now, since we did smoothing in the previous refinement steps also, each cell may only have vertices with levels at most one greater than the level of the present cell.
However, if we store the level plus one for cells marked for refinement, we may end up with cells which have vertices of level two greater than the cells level. We need to refine this cell also, and need thus also update the levels of its vertices. This itself may lead to cells needing refinement, but these are on lower levels, as above, which is why we may do all kinds of additional flagging in one loop only.
eliminate_unrefined_islands:
For each cell we count the number of neighbors which are refined or flagged for refinement. If this exceeds the number of neighbors which are not refined and not flagged for refinement, then the current cell is flagged for refinement. Since this may lead to cells on the same level which also will need refinement, we will need additional loops of regularization and smoothing over all cells until nothing changes any more.
eliminate_refined_*_islands
: This one does much the same as the above one, but for coarsening. If a cell is flagged for refinement or if all of its children are active and if the number of neighbors which are either active and not flagged for refinement, or not active but all children flagged for coarsening equals the total number of neighbors, then this cell's children are flagged for coarsening or (if this cell was flagged for refinement) the refine flag is cleared.
For a description of the distinction between the two versions of the flag see above in the section about mesh smoothing in the general part of this classes description.
The same applies as above: several loops may be necessary.
Regularization and smoothing are a bit complementary in that we check whether we need to set additional refinement flags when being on a cell flagged for refinement (regularization) or on a cell not flagged for refinement. This makes readable programming easier.
All the described algorithms apply only for more than one space dimension, since for one dimension no restrictions apply. It may be necessary to apply some smoothing for multigrid algorithms, but this has to be decided upon later.
It seems impossible to preserve constness
of a triangulation through iterator usage. Thus, if you declare pointers to a const
triangulation object, you should be well aware that you might involuntarily alter the data stored in the triangulation.
Definition at line 50 of file dof_handler.h.
|
private |
|
private |
Declare a number of iterator types for raw iterators, i.e., iterators that also iterate over holes in the list of cells left by cells that have been coarsened away in previous mesh refinement cycles.
Since users should never have to access these internal properties of how we store data, these iterator types are made private.
enum Triangulation::MeshSmoothing |
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. There 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. |
enum Triangulation::CellStatus |
Used to inform functions in derived classes how the cell with the given cell_iterator is going to change. Note that this may me different than the refine_flag() and coarsen_flag() in the cell_iterator in parallel calculations because of refinement constraints that this machine does not see.
Triangulation< dim, spacedim >::Triangulation | ( | const MeshSmoothing | smooth_grid = none , |
const bool | check_for_distorted_cells = false |
||
) |
Create an empty triangulation. Do not create any cells.
smooth_grid | Determines the level of smoothness of the mesh size function that should be enforced upon mesh refinement. |
check_for_distorted_cells | Determines whether the triangulation should check whether any of the cells that are created by create_triangulation() or execute_coarsening_and_refinement() are distorted (see distorted cells). If set, these two functions may throw an exception if they encounter distorted cells. |
|
delete |
Copy constructor.
You should really use the copy_triangulation
function, so this constructor is deleted. The reason for this is that we may want to use triangulation objects in collections. However, C++ containers require that the objects stored in them are copyable, so we need to provide a copy constructor. On the other hand, copying triangulations is so expensive that we do not want such objects copied by accident, for example in compiler-generated temporary objects. By defining a copy constructor but throwing an error, we satisfy the formal requirements of containers, but at the same time disallow actual copies. Finally, through the exception, one easily finds the places where code has to be changed to avoid copies.
|
noexcept |
|
overridevirtual |
Delete the object and all levels of the hierarchy.
Reimplemented in parallel::distributed::Triangulation< 1, spacedim >, parallel::distributed::Triangulation< dim, spacedim >, parallel::distributed::Triangulation< dim >, parallel::shared::Triangulation< dim, spacedim >, parallel::Triangulation< dim, spacedim >, and parallel::Triangulation< 1, spacedim >.
|
noexcept |
|
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 in parallel::distributed::Triangulation< dim, spacedim >, and parallel::distributed::Triangulation< dim >.
|
virtual |
|
virtual |
std::vector< types::boundary_id > Triangulation< dim, spacedim >::get_boundary_ids | ( | ) | const |
Return a vector containing all boundary indicators assigned to boundary faces 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).
|
virtual |
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 >.
|
virtual |
Create a triangulation from a list of vertices and a list of cells, each of the latter being a list of 1<<dim
vertex indices. The triangulation must be empty upon calling this function and the cell list should be useful (connected domain, etc.). The result of calling this function is a coarse mesh.
Material data for the cells is given within the cells
array, while boundary information is given in the subcelldata
field.
The numbering of vertices within the cells
array is subject to some constraints; see the general class documentation for this.
For conditions when this function can generate a valid triangulation, see the documentation of this class, and the GridIn and GridReordering class.
If the check_for_distorted_cells
flag was specified upon creation of this object, at the very end of its operation, the current function walks over all cells and verifies that none of the cells is deformed (see the entry on distorted cells in the glossary), where we call a cell deformed if the determinant of the Jacobian of the mapping from reference cell to real cell is negative at least at one of the vertices (this computation is done using the GeometryInfo::jacobian_determinants_at_vertices function). If there are deformed cells, this function throws an exception of kind DistortedCellList. Since this happens after all data structures have been set up, you can catch and ignore this exception if you know what you do – for example, it may be that the determinant is zero (indicating that you have collapsed edges in a cell) but that this is ok because you didn't intend to integrate on this cell anyway. On the other hand, deformed cells are often a sign of a mesh that is too coarse to resolve the geometry of the domain, and in this case ignoring the exception is probably unwise.
Reimplemented in parallel::distributed::Triangulation< dim, spacedim >, parallel::shared::Triangulation< dim, spacedim >, and PersistentTriangulation< dim, spacedim >.
|
virtual |
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 >.
void Triangulation< dim, spacedim >::flip_all_direction_flags | ( | ) |
Revert or flip the direction_flags of a dim<spacedim triangulation, see GlossDirectionFlag.
This function throws an exception if dim equals spacedim.
void Triangulation< dim, spacedim >::set_all_refine_flags | ( | ) |
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.
void Triangulation< dim, spacedim >::refine_global | ( | const unsigned int | times = 1 | ) |
Refine all cells times
times, by alternatingly calling set_all_refine_flags and execute_coarsening_and_refinement.
The latter function 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.
|
virtual |
Execute both refinement and coarsening of the triangulation.
The function resets all refinement and coarsening flags to false. It uses the user flags for internal purposes. They will therefore be overwritten by undefined content.
To allow user programs to fix up these cells if that is desired, this function after completing all other work may throw an exception of type DistortedCellList that contains a list of those cells that have been refined and have at least one child that is distorted. The function does not create such an exception if no cells have created distorted children. Note that for the check for distorted cells to happen, the check_for_distorted_cells
flag has to be specified upon creation of a triangulation object.
See the general docs for more information.
virtual
to allow derived classes to insert hooks, such as saving refinement flags and the like (see e.g. the PersistentTriangulation class). Reimplemented in parallel::distributed::Triangulation< dim, spacedim >, parallel::distributed::Triangulation< dim >, parallel::shared::Triangulation< dim, spacedim >, and PersistentTriangulation< dim, spacedim >.
|
virtual |
Do both preparation for refinement and coarsening as well as mesh smoothing.
Regarding the refinement process it fixes the closure of the refinement in dim>=2
(make sure that no two cells are adjacent with a refinement level differing with more than one), etc. It performs some mesh smoothing if the according flag was given to the constructor of this class. The function returns whether additional cells have been flagged for refinement.
See the general doc of this class for more information on smoothing upon refinement.
Regarding the coarsening part, flagging and deflagging cells in preparation of the actual coarsening step are done. This includes deleting coarsen flags from cells which may not be deleted (e.g. because one neighbor is more refined than the cell), doing some smoothing, etc.
The effect is that only those cells are flagged for coarsening which will actually be coarsened. This includes the fact that all flagged cells belong to parent cells of which all children are flagged.
The function returns whether some cells' flagging has been changed in the process.
This function uses the user flags, so store them if you still need them afterwards.
Reimplemented in parallel::distributed::Triangulation< dim, spacedim >, and parallel::distributed::Triangulation< dim >.
void Triangulation< dim, spacedim >::save_refine_flags | ( | std::ostream & | out | ) | const |
void Triangulation< dim, spacedim >::save_refine_flags | ( | std::vector< bool > & | v | ) | const |
void Triangulation< dim, spacedim >::load_refine_flags | ( | std::istream & | in | ) |
void Triangulation< dim, spacedim >::load_refine_flags | ( | const std::vector< bool > & | v | ) |
void Triangulation< dim, spacedim >::save_coarsen_flags | ( | std::ostream & | out | ) | const |
void Triangulation< dim, spacedim >::save_coarsen_flags | ( | std::vector< bool > & | v | ) | const |
void Triangulation< dim, spacedim >::load_coarsen_flags | ( | std::istream & | out | ) |
void Triangulation< dim, spacedim >::load_coarsen_flags | ( | const std::vector< bool > & | v | ) |
bool Triangulation< dim, spacedim >::get_anisotropic_refinement_flag | ( | ) | const |
void Triangulation< dim, spacedim >::clear_user_flags | ( | ) |
Clear all user flags. See also GlossUserFlags.
void Triangulation< dim, spacedim >::save_user_flags | ( | std::ostream & | out | ) | const |
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.
void Triangulation< dim, spacedim >::save_user_flags | ( | std::vector< bool > & | v | ) | const |
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.
void Triangulation< dim, spacedim >::load_user_flags | ( | std::istream & | in | ) |
Read the information stored by save_user_flags
. See also GlossUserFlags.
void Triangulation< dim, spacedim >::load_user_flags | ( | const std::vector< bool > & | v | ) |
Read the information stored by save_user_flags
. See also GlossUserFlags.
void Triangulation< dim, spacedim >::clear_user_flags_line | ( | ) |
Clear all user flags on lines. See also GlossUserFlags.
void Triangulation< dim, spacedim >::save_user_flags_line | ( | std::ostream & | out | ) | const |
Save the user flags on lines. See also GlossUserFlags.
void Triangulation< dim, spacedim >::save_user_flags_line | ( | std::vector< bool > & | v | ) | const |
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.
void Triangulation< dim, spacedim >::load_user_flags_line | ( | std::istream & | in | ) |
Load the user flags located on lines. See also GlossUserFlags.
void Triangulation< dim, spacedim >::load_user_flags_line | ( | const std::vector< bool > & | v | ) |
Load the user flags located on lines. See also GlossUserFlags.
void Triangulation< dim, spacedim >::clear_user_flags_quad | ( | ) |
Clear all user flags on quads. See also GlossUserFlags.
void Triangulation< dim, spacedim >::save_user_flags_quad | ( | std::ostream & | out | ) | const |
Save the user flags on quads. See also GlossUserFlags.
void Triangulation< dim, spacedim >::save_user_flags_quad | ( | std::vector< bool > & | v | ) | const |
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.
void Triangulation< dim, spacedim >::load_user_flags_quad | ( | std::istream & | in | ) |
Load the user flags located on quads. See also GlossUserFlags.
void Triangulation< dim, spacedim >::load_user_flags_quad | ( | const std::vector< bool > & | v | ) |
Load the user flags located on quads. See also GlossUserFlags.
void Triangulation< dim, spacedim >::clear_user_flags_hex | ( | ) |
Clear all user flags on quads. See also GlossUserFlags.
void Triangulation< dim, spacedim >::save_user_flags_hex | ( | std::ostream & | out | ) | const |
Save the user flags on hexs. See also GlossUserFlags.
void Triangulation< dim, spacedim >::save_user_flags_hex | ( | std::vector< bool > & | v | ) | const |
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.
void Triangulation< dim, spacedim >::load_user_flags_hex | ( | std::istream & | in | ) |
Load the user flags located on hexs. See also GlossUserFlags.
void Triangulation< dim, spacedim >::load_user_flags_hex | ( | const std::vector< bool > & | v | ) |
Load the user flags located on hexs. See also GlossUserFlags.
void Triangulation< dim, spacedim >::clear_user_data | ( | ) |
Clear all user pointers and indices and allow the use of both for next access. See also GlossUserData.
void Triangulation< dim, spacedim >::save_user_indices | ( | std::vector< unsigned int > & | v | ) | const |
Save all user indices. The output vector is resized if necessary. See also GlossUserData.
void Triangulation< dim, spacedim >::load_user_indices | ( | const std::vector< unsigned int > & | v | ) |
Read the information stored by save_user_indices(). See also GlossUserData.
void Triangulation< dim, spacedim >::save_user_pointers | ( | std::vector< void *> & | v | ) | const |
Save all user pointers. The output vector is resized if necessary. See also GlossUserData.
void Triangulation< dim, spacedim >::load_user_pointers | ( | const std::vector< void *> & | v | ) |
Read the information stored by save_user_pointers(). See also GlossUserData.
void Triangulation< dim, spacedim >::save_user_indices_line | ( | std::vector< unsigned int > & | v | ) | const |
Save the user indices on lines. The output vector is resized if necessary. See also GlossUserData.
void Triangulation< dim, spacedim >::load_user_indices_line | ( | const std::vector< unsigned int > & | v | ) |
Load the user indices located on lines. See also GlossUserData.
void Triangulation< dim, spacedim >::save_user_indices_quad | ( | std::vector< unsigned int > & | v | ) | const |
Save the user indices on quads. The output vector is resized if necessary. See also GlossUserData.
void Triangulation< dim, spacedim >::load_user_indices_quad | ( | const std::vector< unsigned int > & | v | ) |
Load the user indices located on quads. See also GlossUserData.
void Triangulation< dim, spacedim >::save_user_indices_hex | ( | std::vector< unsigned int > & | v | ) | const |
Save the user indices on hexes. The output vector is resized if necessary. See also GlossUserData.
void Triangulation< dim, spacedim >::load_user_indices_hex | ( | const std::vector< unsigned int > & | v | ) |
Load the user indices located on hexs. See also GlossUserData.
void Triangulation< dim, spacedim >::save_user_pointers_line | ( | std::vector< void *> & | v | ) | const |
Save the user indices on lines. The output vector is resized if necessary. See also GlossUserData.
void Triangulation< dim, spacedim >::load_user_pointers_line | ( | const std::vector< void *> & | v | ) |
Load the user pointers located on lines. See also GlossUserData.
void Triangulation< dim, spacedim >::save_user_pointers_quad | ( | std::vector< void *> & | v | ) | const |
Save the user pointers on quads. The output vector is resized if necessary. See also GlossUserData.
void Triangulation< dim, spacedim >::load_user_pointers_quad | ( | const std::vector< void *> & | v | ) |
Load the user pointers located on quads. See also GlossUserData.
void Triangulation< dim, spacedim >::save_user_pointers_hex | ( | std::vector< void *> & | v | ) | const |
Save the user pointers on hexes. The output vector is resized if necessary. See also GlossUserData.
void Triangulation< dim, spacedim >::load_user_pointers_hex | ( | const std::vector< void *> & | v | ) |
Load the user pointers located on hexs. See also GlossUserData.
Triangulation< dim, spacedim >::cell_iterator Triangulation< dim, spacedim >::begin | ( | const unsigned int | level = 0 | ) | const |
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. Triangulation< dim, spacedim >::active_cell_iterator Triangulation< dim, spacedim >::begin_active | ( | const unsigned int | level = 0 | ) | const |
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. Triangulation< dim, spacedim >::cell_iterator Triangulation< dim, spacedim >::end | ( | ) | const |
Triangulation< dim, spacedim >::cell_iterator Triangulation< dim, spacedim >::end | ( | const unsigned int | level | ) | const |
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. Triangulation< dim, spacedim >::active_cell_iterator Triangulation< dim, spacedim >::end_active | ( | const unsigned int | level | ) | const |
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. Triangulation< dim, spacedim >::cell_iterator Triangulation< dim, spacedim >::last | ( | ) | const |
Triangulation< dim, spacedim >::active_cell_iterator Triangulation< dim, spacedim >::last_active | ( | ) | const |
Triangulation< dim, spacedim >::face_iterator Triangulation< dim, spacedim >::begin_face | ( | ) | const |
Triangulation< dim, spacedim >::active_face_iterator Triangulation< dim, spacedim >::begin_active_face | ( | ) | const |
Triangulation< dim, spacedim >::face_iterator Triangulation< dim, spacedim >::end_face | ( | ) | const |
Triangulation< dim, spacedim >::vertex_iterator Triangulation< dim, spacedim >::begin_vertex | ( | ) | const |
Triangulation< dim, spacedim >::active_vertex_iterator Triangulation< dim, spacedim >::begin_active_vertex | ( | ) | const |
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.
Triangulation< dim, spacedim >::vertex_iterator Triangulation< dim, spacedim >::end_vertex | ( | ) | const |
unsigned int Triangulation< dim, spacedim >::n_lines | ( | ) | const |
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.
unsigned int Triangulation< dim, spacedim >::n_lines | ( | const unsigned int | level | ) | const |
unsigned int Triangulation< dim, spacedim >::n_active_lines | ( | ) | const |
unsigned int Triangulation< dim, spacedim >::n_active_lines | ( | const unsigned int | level | ) | const |
unsigned int Triangulation< dim, spacedim >::n_quads | ( | ) | const |
unsigned int Triangulation< dim, spacedim >::n_quads | ( | const unsigned int | level | ) | const |
unsigned int Triangulation< dim, spacedim >::n_active_quads | ( | ) | const |
unsigned int Triangulation< dim, spacedim >::n_active_quads | ( | const unsigned int | level | ) | const |
unsigned int Triangulation< dim, spacedim >::n_hexs | ( | ) | const |
unsigned int Triangulation< dim, spacedim >::n_hexs | ( | const unsigned int | level | ) | const |
unsigned int Triangulation< dim, spacedim >::n_active_hexs | ( | ) | const |
unsigned int Triangulation< dim, spacedim >::n_active_hexs | ( | const unsigned int | level | ) | const |
unsigned int Triangulation< dim, spacedim >::n_cells | ( | ) | const |
unsigned int Triangulation< dim, spacedim >::n_cells | ( | const unsigned int | level | ) | const |
unsigned int Triangulation< dim, spacedim >::n_active_cells | ( | ) | const |
Return the total number of active cells. Maps to n_active_lines()
in one space dimension and so on.
|
virtual |
Return the total number of active cells. For the current class, this is the same as n_active_cells(). However, the function may be overloaded in derived classes (e.g., in parallel::distributed::Triangulation) where it may return a value greater than the number of active cells reported by the triangulation object on the current processor.
Reimplemented in parallel::Triangulation< dim, spacedim >, and parallel::Triangulation< 1, spacedim >.
unsigned int Triangulation< dim, spacedim >::n_active_cells | ( | const unsigned int | level | ) | const |
unsigned int Triangulation< dim, spacedim >::n_faces | ( | ) | const |
unsigned int Triangulation< dim, spacedim >::n_active_faces | ( | ) | const |
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.
unsigned int Triangulation< dim, spacedim >::n_levels | ( | ) | const |
Return the number of levels in this triangulation.
|
virtual |
Return the number of levels in use. This function is equivalent to n_levels() for a serial Triangulation, but gives the maximum of n_levels() over all processors for a parallel::distributed::Triangulation and therefore can be larger than n_levels().
Reimplemented in parallel::Triangulation< dim, spacedim >, and parallel::Triangulation< 1, spacedim >.
|
virtual |
Return true if the triangulation has hanging nodes.
The function is made virtual since the result can be interpreted in different ways, depending on whether the triangulation lives only on a single processor, or may be distributed as done in the parallel::distributed::Triangulation class (see there for a description of what the function is supposed to do in the parallel context).
Reimplemented in parallel::distributed::Triangulation< dim, spacedim >, and parallel::distributed::Triangulation< dim >.
unsigned int Triangulation< dim, spacedim >::n_vertices | ( | ) | const |
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.
const std::vector<Point<spacedim> >& Triangulation< dim, spacedim >::get_vertices | ( | ) | const |
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().
unsigned int Triangulation< dim, spacedim >::n_used_vertices | ( | ) | const |
bool Triangulation< dim, spacedim >::vertex_used | ( | const unsigned int | index | ) | const |
Return true
if the vertex with this index
is used.
const std::vector< bool > & Triangulation< dim, spacedim >::get_used_vertices | ( | ) | const |
unsigned int Triangulation< dim, spacedim >::max_adjacent_cells | ( | ) | const |
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.
|
virtual |
This function always returns invalid_subdomain_id
but is there for compatibility with the derived parallel::distributed::Triangulation
class. For distributed parallel triangulations this function returns the subdomain id of those cells that are owned by the current processor.
Reimplemented in parallel::Triangulation< dim, spacedim >, and parallel::Triangulation< 1, spacedim >.
Triangulation< dim, spacedim > & Triangulation< dim, spacedim >::get_triangulation | ( | ) |
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).
const Triangulation< dim, spacedim > & Triangulation< dim, spacedim >::get_triangulation | ( | ) | const |
unsigned int Triangulation< dim, spacedim >::n_raw_lines | ( | ) | const |
Total number of lines, used or unused.
unsigned int Triangulation< dim, spacedim >::n_raw_lines | ( | const unsigned int | level | ) | const |
Number of lines, used or unused, on the given level.
unsigned int Triangulation< dim, spacedim >::n_raw_quads | ( | ) | const |
Total number of quads, used or unused.
unsigned int Triangulation< dim, spacedim >::n_raw_quads | ( | const unsigned int | level | ) | const |
Number of quads, used or unused, on the given level.
unsigned int Triangulation< dim, spacedim >::n_raw_hexs | ( | const unsigned int | level | ) | const |
Number of hexs, used or unused, on the given level.
unsigned int Triangulation< dim, spacedim >::n_raw_cells | ( | const unsigned int | level | ) | const |
Number of cells, used or unused, on the given level.
unsigned int Triangulation< dim, spacedim >::n_raw_faces | ( | ) | const |
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.
|
virtual |
Determine an estimate for the memory consumption (in bytes) of this object.
This function is made virtual, since a triangulation object might be accessed through a pointer to this base class, even if the actual object is a derived class.
Reimplemented in parallel::distributed::Triangulation< dim, spacedim >, parallel::distributed::Triangulation< dim >, PersistentTriangulation< dim, spacedim >, parallel::Triangulation< dim, spacedim >, and parallel::Triangulation< 1, spacedim >.
void Triangulation< dim, spacedim >::save | ( | Archive & | ar, |
const unsigned int | version | ||
) | const |
Write the data of this object to a stream for the purpose of serialization.
void Triangulation< dim, spacedim >::load | ( | Archive & | ar, |
const unsigned int | version | ||
) |
Read the data of this object from a stream for the purpose of serialization. Throw away the previous content.
|
virtual |
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.
const std::map< std::pair< typename Triangulation< dim, spacedim >::cell_iterator, unsigned int >, std::pair< std::pair< typename Triangulation< dim, spacedim >::cell_iterator, unsigned int >, std::bitset< 3 > > > & Triangulation< dim, spacedim >::get_periodic_face_map | ( | ) | const |
|
staticprotected |
Write a bool vector to the given stream, writing a pre- and a postfix magic number. The vector is written in an almost binary format, i.e. the bool flags are packed but the data is written as ASCII text.
The flags are stored in a binary format: for each true
, a 1
bit is stored, a 0
bit otherwise. The bits are stored as unsigned char
, thus avoiding endianness. They are written to out
in plain text, thus amounting to 3.6 bits in the output per bits in the input on the average. Other information (magic numbers and number of elements of the input vector) is stored as plain text as well. The format should therefore be interplatform compatible.
|
staticprotected |
|
protected |
|
private |
|
private |
|
private |
|
private |
Iterator to the first used line 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.
|
private |
Iterator to the first active line 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.
|
private |
|
private |
Iterator to the first quad, used or not, on the given level. If a level has no quads, a past-the-end iterator is returned. If quads are no cells, i.e. for \(dim>2\) no level argument must be given.
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.
|
private |
Iterator to the first used quad 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.
|
private |
Iterator to the first active quad 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.
|
private |
|
private |
Iterator to the first hex, used or not, on level level
. If a level has no hexes, a past-the-end iterator is returned.
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.
|
private |
Iterator to the first used hex 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.
|
private |
Iterator to the first active hex 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.
|
private |
|
private |
The (public) function clear() will only work when the triangulation is not subscribed to by other users. The clear_despite_subscriptions() function now allows the triangulation being cleared even when there are subscriptions.
Make sure, you know what you do, when calling this function, as its use is reasonable in very rare cases, only. For example, when the subscriptions were for the initially empty Triangulation and the Triangulation object wants to release its memory before throwing an assertion due to input errors (e.g. in the create_triangulation() function).
|
private |
|
private |
Refine all cells on all levels which were previously flagged for refinement.
Note, that this function uses the line->user_flags
for dim=2,3
and the quad->user_flags
for dim=3
.
The function returns a list of cells that have produced children that satisfy the criteria of distorted cells if the check_for_distorted_cells
flag was specified upon creation of this object, at
|
private |
|
private |
|
static |
|
static |
|
mutable |
|
protected |
|
private |
If add_periodicity() is called, this variable stores the given periodic face pairs on level 0 for later access during the identification of ghost cells for the multigrid hierarchy and for setting up the periodic_face_map.
|
private |
If add_periodicity() is called, this variable stores the active periodic face pairs.
|
private |
|
private |
|
private |
|
private |
|
private |
Collection of manifold objects. We store only objects, which are not of type FlatManifold.
|
private |
|
private |
|
private |
Cache to hold the numbers of lines, quads, hexes, etc. These numbers are set at the end of the refinement and coarsening functions and enable faster access later on. In the old days, whenever one wanted to access one of these numbers, one had to perform a loop over all lines, e.g., and count the elements until we hit the end iterator. This is time consuming and since access to the number of lines etc is a rather frequent operation, this was not an optimal solution.
|
private |
A map that relates the number of a boundary vertex to the boundary indicator. This field is only used in 1d. We have this field because we store boundary indicator information with faces in 2d and higher where we have space in the structures that store data for faces, but in 1d there is no such space for faces.
The field is declared as a pointer for a rather mundane reason: all other fields of this class that can be modified by the TriaAccessor hierarchy are pointers, and so these accessor classes store a const pointer to the triangulation. We could no longer do so for TriaAccessor<0,1,spacedim> if this field (that can be modified by TriaAccessor::set_boundary_id) were not a pointer.
|
private |
A map that relates the number of a boundary vertex to the manifold indicator. This field is only used in 1d. We have this field because we store manifold indicator information with faces in 2d and higher where we have space in the structures that store data for faces, but in 1d there is no such space for faces.
The field is declared as a pointer for a rather mundane reason: all other fields of this class that can be modified by the TriaAccessor hierarchy are pointers, and so these accessor classes store a const pointer to the triangulation. We could no longer do so for TriaAccessor<0,1,spacedim> if this field (that can be modified by TriaAccessor::set_manifold_id) were not a pointer.