Reference documentation for deal.II version 9.6.0
|
#include <deal.II/grid/tria_accessor.h>
Public Types | |
enum | VertexKind { left_vertex , interior_vertex , right_vertex } |
using | AccessorData = void |
Public Member Functions | |
TriaAccessor (const Triangulation< 1, spacedim > *tria, const VertexKind vertex_kind, const unsigned int vertex_index) | |
TriaAccessor (const Triangulation< 1, spacedim > *tria=nullptr, const int=0, const int=0, const AccessorData *=nullptr) | |
template<int structdim2, int dim2, int spacedim2> | |
TriaAccessor (const TriaAccessor< structdim2, dim2, spacedim2 > &) | |
template<int structdim2, int dim2, int spacedim2> | |
TriaAccessor (const InvalidAccessor< structdim2, dim2, spacedim2 > &) | |
void | copy_from (const TriaAccessor &) |
void | copy_from (const TriaAccessorBase< 0, 1, spacedim > &) |
int | index () const |
const Triangulation< 1, spacedim > & | get_triangulation () const |
bool | at_boundary () const |
types::boundary_id | boundary_id () const |
const Manifold< 1, spacedim > & | get_manifold () const |
types::manifold_id | manifold_id () const |
bool | used () const |
ReferenceCell | reference_cell () const |
unsigned int | n_vertices () const |
unsigned int | n_lines () const |
std_cxx20::ranges::iota_view< unsigned int, unsigned int > | vertex_indices () const |
std_cxx20::ranges::iota_view< unsigned int, unsigned int > | line_indices () const |
Advancement of iterators | |
void | operator++ () const |
void | operator-- () const |
bool | operator== (const TriaAccessor &) const |
bool | operator!= (const TriaAccessor &) const |
bool | operator< (const TriaAccessor &other) const |
User data | |
bool | user_flag_set () const |
void | set_user_flag () const |
void | clear_user_flag () const |
void | recursively_set_user_flag () const |
void | recursively_clear_user_flag () const |
void | clear_user_data () const |
void | set_user_pointer (void *p) const |
void | clear_user_pointer () const |
void * | user_pointer () const |
void | recursively_set_user_pointer (void *p) const |
void | recursively_clear_user_pointer () const |
void | set_user_index (const unsigned int p) const |
void | clear_user_index () const |
unsigned int | user_index () const |
void | recursively_set_user_index (const unsigned int p) const |
void | recursively_clear_user_index () const |
Dealing with boundary indicators | |
void | set_boundary_id (const types::boundary_id) const |
void | set_manifold_id (const types::manifold_id) |
void | set_all_boundary_ids (const types::boundary_id) const |
void | set_all_manifold_ids (const types::manifold_id) |
Static Public Member Functions | |
static IteratorState::IteratorStates | state () |
static int | level () |
Orientation of sub-objects | |
static unsigned char | combined_face_orientation (const unsigned int face) |
Always return 0. | |
static bool | face_orientation (const unsigned int face) |
Always return false. | |
static bool | face_flip (const unsigned int face) |
Always return false. | |
static bool | face_rotation (const unsigned int face) |
Always return false. | |
static bool | line_orientation (const unsigned int line) |
Always return false. | |
Accessing children | |
static bool | has_children () |
static unsigned int | n_children () |
static unsigned int | n_active_descendants () |
static unsigned int | max_refinement_depth () |
static unsigned int | child_iterator_to_index (const TriaIterator< TriaAccessor< 0, 1, spacedim > > &) |
Return an invalid unsigned integer. | |
static TriaIterator< TriaAccessor< 0, 1, spacedim > > | child (const unsigned int) |
Return an invalid object. | |
static TriaIterator< TriaAccessor< 0, 1, spacedim > > | isotropic_child (const unsigned int) |
Return an invalid object. | |
static RefinementCase< 0 > | refinement_case () |
static int | child_index (const unsigned int i) |
Returns -1. | |
static int | isotropic_child_index (const unsigned int i) |
Returns -1. | |
Static Public Attributes | |
static constexpr unsigned int | space_dimension = spacedim |
static constexpr unsigned int | dimension = 1 |
static const unsigned int | structure_dimension = 0 |
Protected Attributes | |
const Triangulation< 1, spacedim > * | tria |
VertexKind | vertex_kind |
unsigned int | global_vertex_index |
Accessing sub-objects | |
unsigned int | vertex_index (const unsigned int i=0) const |
Point< spacedim > & | vertex (const unsigned int i=0) const |
Point< spacedim > | center () const |
static typename::internal::TriangulationImplementation::Iterators< 1, spacedim >::line_iterator | line (const unsigned int) |
static unsigned int | line_index (const unsigned int i) |
static typename::internal::TriangulationImplementation::Iterators< 1, spacedim >::quad_iterator | quad (const unsigned int i) |
static unsigned int | quad_index (const unsigned int i) |
This class is a specialization of TriaAccessor<structdim, dim, spacedim>
for the case that structdim
is zero and dim
is one. This class represents vertices in a one-dimensional triangulation that is embedded in a space of dimensionality spacedim
(for spacedim==dim==1
the triangulation represents a domain in \({\mathbb R}^\text{dim}\), for spacedim>dim==1
the triangulation is of a manifold embedded in a higher dimensional space).
The current specialization of the TriaAccessor<0,dim,spacedim> class for vertices of a one-dimensional triangulation exists since in the dim
== 1 case vertices are also faces.
Definition at line 2319 of file tria_accessor.h.
using TriaAccessor< 0, 1, spacedim >::AccessorData = void |
Pointer to internal data.
Definition at line 2346 of file tria_accessor.h.
enum TriaAccessor< 0, 1, spacedim >::VertexKind |
Whether the vertex represented here is at the left end of the domain, the right end, or in the interior.
Enumerator | |
---|---|
left_vertex | Left vertex. |
interior_vertex | Interior vertex. |
right_vertex | Right vertex. |
Definition at line 2352 of file tria_accessor.h.
TriaAccessor< 0, 1, spacedim >::TriaAccessor | ( | const Triangulation< 1, spacedim > * | tria, |
const VertexKind | vertex_kind, | ||
const unsigned int | vertex_index ) |
Constructor.
Since there is no mapping from vertices to cells, an accessor object for a point has no way to figure out whether it is at the boundary of the domain or not. Consequently, the second argument must be passed by the object that generates this accessor – e.g. a 1d cell that can figure out whether its left or right vertex are at the boundary.
The third argument is the global index of the vertex we point to.
TriaAccessor< 0, 1, spacedim >::TriaAccessor | ( | const Triangulation< 1, spacedim > * | tria = nullptr, |
const int | = 0, | ||
const int | = 0, | ||
const AccessorData * | = nullptr ) |
Constructor. This constructor exists in order to maintain interface compatibility with the other accessor classes. However, it doesn't do anything useful here and so may not actually be called.
TriaAccessor< 0, 1, spacedim >::TriaAccessor | ( | const TriaAccessor< structdim2, dim2, spacedim2 > & | ) |
Constructor. Should never be called and thus produces an error.
TriaAccessor< 0, 1, spacedim >::TriaAccessor | ( | const InvalidAccessor< structdim2, dim2, spacedim2 > & | ) |
Constructor. Should never be called and thus produces an error.
void TriaAccessor< 0, 1, spacedim >::copy_from | ( | const TriaAccessor< 0, 1, spacedim > & | ) |
Copy operator. Since this is only called from iterators, do not return anything, since the iterator will return itself.
void TriaAccessor< 0, 1, spacedim >::copy_from | ( | const TriaAccessorBase< 0, 1, spacedim > & | ) |
Copy operator. We need this function to support generic programming, but it just throws an exception because it cannot do the required operations.
|
static |
Return the state of the iterator. Since an iterator to points can not be incremented or decremented, its state remains constant, and in particular equal to IteratorState::valid.
|
static |
Level of this object. Vertices have no level, so this function always returns zero.
int TriaAccessor< 0, 1, spacedim >::index | ( | ) | const |
Index of this object. Returns the global index of the vertex this object points to.
const Triangulation< 1, spacedim > & TriaAccessor< 0, 1, spacedim >::get_triangulation | ( | ) | const |
Return a reference to the triangulation which the object pointed to by this class belongs to.
void TriaAccessor< 0, 1, spacedim >::operator++ | ( | ) | const |
This operator advances the iterator to the next element. For points, this operation is not defined, so you can't iterate over point iterators.
void TriaAccessor< 0, 1, spacedim >::operator-- | ( | ) | const |
This operator moves the iterator to the previous element. For points, this operation is not defined, so you can't iterate over point iterators.
bool TriaAccessor< 0, 1, spacedim >::operator== | ( | const TriaAccessor< 0, 1, spacedim > & | ) | const |
Compare for equality.
bool TriaAccessor< 0, 1, spacedim >::operator!= | ( | const TriaAccessor< 0, 1, spacedim > & | ) | const |
Compare for inequality.
bool TriaAccessor< 0, 1, spacedim >::operator< | ( | const TriaAccessor< 0, 1, spacedim > & | other | ) | const |
Comparison operator for accessors. This operator is used when comparing iterators into objects of a triangulation, for example when putting them into a std::map
.
This operator simply compares the global index of the vertex the current object points to.
unsigned int TriaAccessor< 0, 1, spacedim >::vertex_index | ( | const unsigned int | i = 0 | ) | const |
Return the global index of i-th vertex of the current object. If i is zero, this returns the index of the current point to which this object refers. Otherwise, it throws an exception.
Note that the returned value is only the index of the geometrical vertex. It has nothing to do with possible degrees of freedom associated with it. For this, see the DoFAccessor::vertex_dof_index
functions.
Point< spacedim > & TriaAccessor< 0, 1, spacedim >::vertex | ( | const unsigned int | i = 0 | ) | const |
Return a reference to the ith
vertex. If i is zero, this returns a reference to the current point to which this object refers. Otherwise, it throws an exception.
Point< spacedim > TriaAccessor< 0, 1, spacedim >::center | ( | ) | const |
Return the center of this object, which of course coincides with the location of the vertex this object refers to.
|
static |
Pointer to the ith
line bounding this object. Will point to an invalid object.
|
static |
Line index of the ith
line bounding this object.
Implemented only for structdim>1
, otherwise an exception generated.
|
static |
Pointer to the ith
quad bounding this object.
|
static |
Quad index of the ith
quad bounding this object.
Implemented only for structdim>2
, otherwise an exception generated.
bool TriaAccessor< 0, 1, spacedim >::at_boundary | ( | ) | const |
Return whether this point is at the boundary of the one-dimensional triangulation we deal with here.
types::boundary_id TriaAccessor< 0, 1, spacedim >::boundary_id | ( | ) | const |
Return the boundary indicator of this object. The convention for one dimensional triangulations is that left end vertices (of each line segment from which the triangulation may be constructed) have boundary indicator zero, and right end vertices have boundary indicator one, unless explicitly set differently.
If the return value is the special value numbers::internal_face_boundary_id, then this object is in the interior of the domain.
const Manifold< 1, spacedim > & TriaAccessor< 0, 1, spacedim >::get_manifold | ( | ) | const |
Return a constant reference to the manifold object used for this object.
types::manifold_id TriaAccessor< 0, 1, spacedim >::manifold_id | ( | ) | const |
Return the manifold indicator of this object.
bool TriaAccessor< 0, 1, spacedim >::user_flag_set | ( | ) | const |
Read the user flag. See GlossUserFlags for more information.
Definition at line 1736 of file tria_accessor.cc.
void TriaAccessor< 0, 1, spacedim >::set_user_flag | ( | ) | const |
Set the user flag. See GlossUserFlags for more information.
Definition at line 1747 of file tria_accessor.cc.
void TriaAccessor< 0, 1, spacedim >::clear_user_flag | ( | ) | const |
Clear the user flag. See GlossUserFlags for more information.
Definition at line 1757 of file tria_accessor.cc.
void TriaAccessor< 0, 1, spacedim >::recursively_set_user_flag | ( | ) | const |
Set the user flag for this and all descendants. See GlossUserFlags for more information.
Definition at line 1767 of file tria_accessor.cc.
void TriaAccessor< 0, 1, spacedim >::recursively_clear_user_flag | ( | ) | const |
Clear the user flag for this and all descendants. See GlossUserFlags for more information.
Definition at line 1780 of file tria_accessor.cc.
void TriaAccessor< 0, 1, spacedim >::clear_user_data | ( | ) | const |
Reset the user data to zero, independent if pointer or index. See GlossUserData for more information.
Definition at line 1793 of file tria_accessor.cc.
void TriaAccessor< 0, 1, spacedim >::set_user_pointer | ( | void * | p | ) | const |
Set the user pointer to p
.
See GlossUserData for more information.
Definition at line 1803 of file tria_accessor.cc.
void TriaAccessor< 0, 1, spacedim >::clear_user_pointer | ( | ) | const |
Reset the user pointer to nullptr
. See GlossUserData for more information.
Definition at line 1813 of file tria_accessor.cc.
void * TriaAccessor< 0, 1, spacedim >::user_pointer | ( | ) | const |
Access the value of the user pointer. It is in the responsibility of the user to make sure that the pointer points to something useful. You should use the new style cast operator to maintain a minimum of type safety, e.g.
A a=static_cast<A>(cell->user_pointer());
.See GlossUserData for more information.
Definition at line 1823 of file tria_accessor.cc.
void TriaAccessor< 0, 1, spacedim >::recursively_set_user_pointer | ( | void * | p | ) | const |
Set the user pointer of this object and all its children to the given value. This is useful for example if all cells of a certain subdomain, or all faces of a certain part of the boundary should have user pointers pointing to objects describing this part of the domain or boundary.
Note that the user pointer is not inherited under mesh refinement, so after mesh refinement there might be cells or faces that don't have user pointers pointing to the describing object. In this case, simply loop over all the elements of the coarsest level that has this information, and use this function to recursively set the user pointer of all finer levels of the triangulation.
See GlossUserData for more information.
Definition at line 1834 of file tria_accessor.cc.
void TriaAccessor< 0, 1, spacedim >::recursively_clear_user_pointer | ( | ) | const |
Clear the user pointer of this object and all of its descendants. The same holds as said for the recursively_set_user_pointer() function. See GlossUserData for more information.
Definition at line 1847 of file tria_accessor.cc.
void TriaAccessor< 0, 1, spacedim >::set_user_index | ( | const unsigned int | p | ) | const |
Set the user index to p
.
Definition at line 1860 of file tria_accessor.cc.
void TriaAccessor< 0, 1, spacedim >::clear_user_index | ( | ) | const |
Reset the user index to 0. See GlossUserData for more information.
Definition at line 1870 of file tria_accessor.cc.
unsigned int TriaAccessor< 0, 1, spacedim >::user_index | ( | ) | const |
Access the value of the user index.
See GlossUserData for more information.
Definition at line 1880 of file tria_accessor.cc.
void TriaAccessor< 0, 1, spacedim >::recursively_set_user_index | ( | const unsigned int | p | ) | const |
Set the user index of this object and all its children.
Note that the user index is not inherited under mesh refinement, so after mesh refinement there might be cells or faces that don't have the expected user indices. In this case, simply loop over all the elements of the coarsest level that has this information, and use this function to recursively set the user index of all finer levels of the triangulation.
See GlossUserData for more information.
Definition at line 1891 of file tria_accessor.cc.
void TriaAccessor< 0, 1, spacedim >::recursively_clear_user_index | ( | ) | const |
Clear the user index of this object and all of its descendants. The same holds as said for the recursively_set_user_index() function.
See GlossUserData for more information.
Definition at line 1904 of file tria_accessor.cc.
|
static |
Test whether the object has children. Always false.
|
static |
Return the number of immediate children of this object.This is always zero in dimension 0.
|
static |
Compute and return the number of active descendants of this objects. Always zero.
|
static |
Return the number of times that this object is refined. Always 0.
|
static |
Return an invalid unsigned integer.
|
static |
Return an invalid object.
|
static |
Return an invalid object.
|
static |
Always return no refinement.
|
static |
Returns -1.
|
static |
Returns -1.
void TriaAccessor< 0, 1, spacedim >::set_boundary_id | ( | const types::boundary_id | ) | const |
Set the boundary indicator. The same applies as for the boundary_id()
function.
void TriaAccessor< 0, 1, spacedim >::set_manifold_id | ( | const types::manifold_id | ) |
Set the manifold indicator of this vertex. This does nothing so far since manifolds are only used to refine and map objects, but vertices are not refined and the mapping is trivial. This function is here only to allow dimension independent programming.
void TriaAccessor< 0, 1, spacedim >::set_all_boundary_ids | ( | const types::boundary_id | ) | const |
Set the boundary indicator of this object and all of its lower-dimensional sub-objects. Since this object only represents a single vertex, there are no lower-dimensional object and this function is equivalent to calling set_boundary_id() with the same argument.
bool TriaAccessor< 0, 1, spacedim >::used | ( | ) | const |
Return whether the vertex pointed to here is used.
ReferenceCell TriaAccessor< 0, 1, spacedim >::reference_cell | ( | ) | const |
Reference cell type of the current object.
unsigned int TriaAccessor< 0, 1, spacedim >::n_vertices | ( | ) | const |
Number of vertices.
unsigned int TriaAccessor< 0, 1, spacedim >::n_lines | ( | ) | const |
Number of lines.
std_cxx20::ranges::iota_view< unsigned int, unsigned int > TriaAccessor< 0, 1, spacedim >::vertex_indices | ( | ) | const |
Return an object that can be thought of as an array containing all indices from zero to n_vertices().
std_cxx20::ranges::iota_view< unsigned int, unsigned int > TriaAccessor< 0, 1, spacedim >::line_indices | ( | ) | const |
Return an object that can be thought of as an array containing all indices from zero to n_lines().
|
staticconstexpr |
Dimension of the space the object represented by this accessor lives in. For example, if this accessor represents a quad that is part of a two-dimensional surface in four-dimensional space, then this value is four.
Definition at line 2327 of file tria_accessor.h.
|
staticconstexpr |
Dimensionality of the object that the thing represented by this accessor is part of. For example, if this accessor represents a line that is part of a hexahedron, then this value will be three.
Definition at line 2334 of file tria_accessor.h.
|
static |
Dimensionality of the current object represented by this accessor. For example, if it is line (irrespective of whether it is part of a 2d or 3d subobject), then this value equals 1.
Definition at line 2341 of file tria_accessor.h.
|
protected |
Pointer to the triangulation we operate on.
Definition at line 3031 of file tria_accessor.h.
|
protected |
Whether this is a left end, right end, or interior vertex. This information is provided by the cell at the time of creation.
Definition at line 3037 of file tria_accessor.h.
|
protected |
The global vertex index of the vertex this object corresponds to.
Definition at line 3042 of file tria_accessor.h.