Reference documentation for deal.II version 9.1.1
|
#include <deal.II/grid/tria_accessor.h>
Public Types | |
using | AccessorData = typename TriaAccessorBase< structdim, dim, spacedim >::AccessorData |
Public Types inherited from TriaAccessorBase< structdim, dim, spacedim > | |
using | LocalData = void * |
Public Member Functions | |
TriaAccessor (const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr) | |
template<int structdim2, int dim2, int spacedim2> | |
TriaAccessor (const InvalidAccessor< structdim2, dim2, spacedim2 > &) | |
template<int structdim2, int dim2, int spacedim2> | |
TriaAccessor (const TriaAccessor< structdim2, dim2, spacedim2 > &) | |
void | operator= (const TriaAccessor &)=delete |
bool | used () const |
template<> | |
double | extent_in_direction (const unsigned int axis) const |
template<> | |
double | extent_in_direction (const unsigned int axis) const |
template<> | |
double | extent_in_direction (const unsigned int axis) const |
Accessing sub-objects | |
TriaIterator< TriaAccessor< 0, dim, spacedim > > | vertex_iterator (const unsigned int i) const |
unsigned int | vertex_index (const unsigned int i) const |
Point< spacedim > & | vertex (const unsigned int i) const |
typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::line_iterator | line (const unsigned int i) const |
unsigned int | line_index (const unsigned int i) const |
typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::quad_iterator | quad (const unsigned int i) const |
unsigned int | quad_index (const unsigned int i) const |
Orientation of sub-objects | |
bool | face_orientation (const unsigned int face) const |
bool | face_flip (const unsigned int face) const |
bool | face_rotation (const unsigned int face) const |
bool | line_orientation (const unsigned int line) const |
Accessing children | |
bool | has_children () const |
unsigned int | n_children () const |
unsigned int | number_of_children () const |
unsigned int | max_refinement_depth () const |
TriaIterator< TriaAccessor< structdim, dim, spacedim > > | child (const unsigned int i) const |
TriaIterator< TriaAccessor< structdim, dim, spacedim > > | isotropic_child (const unsigned int i) const |
RefinementCase< structdim > | refinement_case () const |
int | child_index (const unsigned int i) const |
int | isotropic_child_index (const unsigned int i) const |
Dealing with boundary indicators | |
types::boundary_id | boundary_id () const |
void | set_boundary_id (const types::boundary_id) const |
void | set_all_boundary_ids (const types::boundary_id) const |
bool | at_boundary () const |
const Manifold< dim, spacedim > & | get_manifold () const |
Dealing with manifold indicators | |
types::manifold_id | manifold_id () const |
void | set_manifold_id (const types::manifold_id) const |
void | set_all_manifold_ids (const types::manifold_id) 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 |
Geometric information about an object | |
double | diameter () const |
std::pair< Point< spacedim >, double > | enclosing_ball () const |
BoundingBox< spacedim > | bounding_box () const |
double | extent_in_direction (const unsigned int axis) const |
double | minimum_vertex_distance () const |
Point< spacedim > | intermediate_point (const Point< structdim > &coordinates) const |
Point< structdim > | real_to_unit_cell_affine_approximation (const Point< spacedim > &point) const |
Point< spacedim > | center (const bool respect_manifold=false, const bool interpolate_from_surrounding=false) const |
Point< spacedim > | barycenter () const |
double | measure () const |
bool | is_translation_of (const TriaIterator< TriaAccessor< structdim, dim, spacedim >> &o) const |
Public Member Functions inherited from TriaAccessorBase< structdim, dim, spacedim > | |
void | operator= (const TriaAccessorBase *)=delete |
int | level () const |
int | index () const |
IteratorState::IteratorStates | state () const |
const Triangulation< dim, spacedim > & | get_triangulation () const |
Private Member Functions | |
void | set_boundary_id_internal (const types::boundary_id id) const |
void | set (const ::internal::TriangulationImplementation::TriaObject< structdim > &o) const |
void | set_line_orientation (const unsigned int line, const bool orientation) const |
void | set_face_orientation (const unsigned int face, const bool orientation) const |
void | set_face_flip (const unsigned int face, const bool flip) const |
void | set_face_rotation (const unsigned int face, const bool rotation) const |
void | set_used_flag () const |
void | clear_used_flag () const |
void | set_refinement_case (const RefinementCase< structdim > &ref_case) const |
void | clear_refinement_case () const |
void | set_children (const unsigned int i, const int index) const |
void | clear_children () const |
Friends | |
template<int , int > | |
class | Triangulation |
A class that provides access to objects in a triangulation such as its vertices, sub-objects, children, geometric information, etc. This class represents objects of dimension structdim
(i.e. 1 for lines, 2 for quads, 3 for hexes) in a triangulation of dimensionality dim
(i.e. 1 for a triangulation of lines, 2 for a triangulation of quads, and 3 for a triangulation of hexes) that is embedded in a space of dimensionality spacedim
(for spacedim==dim
the triangulation represents a domain in \(R^{dim}\), for spacedim>dim
the triangulation is of a manifold embedded in a higher dimensional space).
There is a specialization of this class for the case where structdim
equals zero, i.e., for vertices of a triangulation.
using TriaAccessor< structdim, dim, spacedim >::AccessorData = typename TriaAccessorBase<structdim, dim, spacedim>::AccessorData |
Propagate alias from base class to this class.
Definition at line 704 of file tria_accessor.h.
TriaAccessor< structdim, dim, spacedim >::TriaAccessor | ( | const Triangulation< dim, spacedim > * | parent = nullptr , |
const int | level = -1 , |
||
const int | index = -1 , |
||
const AccessorData * | local_data = nullptr |
||
) |
Constructor.
TriaAccessor< structdim, dim, spacedim >::TriaAccessor | ( | const InvalidAccessor< structdim2, dim2, spacedim2 > & | ) |
Conversion constructor. This constructor exists to make certain constructs simpler to write in dimension independent code. For example, it allows assigning a face iterator to a line iterator, an operation that is useful in 2d but doesn't make any sense in 3d. The constructor here exists for the purpose of making the code conform to C++ but it will unconditionally abort; in other words, assigning a face iterator to a line iterator is better put into an if-statement that checks that the dimension is two, and assign to a quad iterator in 3d (an operator that, without this constructor would be illegal if we happen to compile for 2d).
Definition at line 3610 of file tria_accessor.h.
TriaAccessor< structdim, dim, spacedim >::TriaAccessor | ( | const TriaAccessor< structdim2, dim2, spacedim2 > & | ) |
Another conversion operator between objects that don't make sense, just like the previous one.
Definition at line 3640 of file tria_accessor.h.
|
delete |
Copy operator. These operators are usually used in a context like iterator a,b; *a=*b;
. Presumably, the intent here is to copy the object pointed to by b
to the object pointed to by a
. However, the result of dereferencing an iterator is not an object but an accessor; consequently, this operation is not useful for iterators on triangulations. Consequently, this operator is declared as deleted and can not be used.
bool TriaAccessor< structdim, dim, spacedim >::used | ( | ) | const |
Test for the element being used or not. The return value is true
for all iterators that are either normal iterators or active iterators, only raw iterators can return false
. Since raw iterators are only used in the interiors of the library, you will not usually need this function.
TriaIterator<TriaAccessor<0, dim, spacedim> > TriaAccessor< structdim, dim, spacedim >::vertex_iterator | ( | const unsigned int | i | ) | const |
Pointer to the ith
vertex bounding this object. Throw an exception if dim=1
.
unsigned int TriaAccessor< structdim, dim, spacedim >::vertex_index | ( | const unsigned int | i | ) | const |
Return the global index of i-th vertex of the current object. The convention regarding the numbering of vertices is laid down in the documentation of the GeometryInfo class.
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< structdim, dim, spacedim >::vertex | ( | const unsigned int | i | ) | const |
Return a reference to the ith
vertex. The reference is not const, i.e., it is possible to call this function on the left hand side of an assignment, thereby moving the vertex of a cell within the triangulation. Of course, doing so requires that you ensure that the new location of the vertex remains useful – for example, avoiding inverted or otherwise distorted (see also this glossary entry).
typename ::internal::TriangulationImplementation:: Iterators<dim, spacedim>::line_iterator TriaAccessor< structdim, dim, spacedim >::line | ( | const unsigned int | i | ) | const |
Pointer to the ith
line bounding this object.
unsigned int TriaAccessor< structdim, dim, spacedim >::line_index | ( | const unsigned int | i | ) | const |
Line index of the ith
line bounding this object.
Implemented only for structdim>1
, otherwise an exception generated.
typename ::internal::TriangulationImplementation:: Iterators<dim, spacedim>::quad_iterator TriaAccessor< structdim, dim, spacedim >::quad | ( | const unsigned int | i | ) | const |
Pointer to the ith
quad bounding this object.
unsigned int TriaAccessor< structdim, dim, spacedim >::quad_index | ( | const unsigned int | i | ) | const |
Quad index of the ith
quad bounding this object.
Implemented only for structdim>2
, otherwise an exception generated.
bool TriaAccessor< structdim, dim, spacedim >::face_orientation | ( | const unsigned int | face | ) | const |
Return whether the face with index face
has its normal pointing in the standard direction (true
) or whether it is the opposite (false
). Which is the standard direction is documented with the GeometryInfo class. In 1d and 2d, this is always true
, but in 3d it may be different, see the respective discussion in the documentation of the GeometryInfo class.
This function is really only for internal use in the library unless you absolutely know what this is all about.
bool TriaAccessor< structdim, dim, spacedim >::face_flip | ( | const unsigned int | face | ) | const |
Return whether the face with index face
is rotated by 180 degrees (true
) or not (false
). In 1d and 2d, this is always false
, but in 3d it may be different, see the respective discussion in the documentation of the GeometryInfo class.
This function is really only for internal use in the library unless you absolutely know what this is all about.
bool TriaAccessor< structdim, dim, spacedim >::face_rotation | ( | const unsigned int | face | ) | const |
Return whether the face with index face
is rotated by 90 degrees (true
) or not (false
). In 1d and 2d, this is always false
, but in 3d it may be different, see the respective discussion in the documentation of the GeometryInfo class.
This function is really only for internal use in the library unless you absolutely know what this is all about.
bool TriaAccessor< structdim, dim, spacedim >::line_orientation | ( | const unsigned int | line | ) | const |
Return whether the line with index line
is oriented in standard direction. true
indicates, that the line is oriented from vertex 0 to vertex 1, whereas it is the other way around otherwise. In 1d and 2d, this is always true
, but in 3d it may be different, see the respective discussion in the documentation of the GeometryInfo class.
This function is really only for internal use in the library unless you absolutely know what this is all about.
bool TriaAccessor< structdim, dim, spacedim >::has_children | ( | ) | const |
Test whether the object has children.
unsigned int TriaAccessor< structdim, dim, spacedim >::n_children | ( | ) | const |
Return the number of immediate children of this object. The number of children of an unrefined cell is zero.
unsigned int TriaAccessor< structdim, dim, spacedim >::number_of_children | ( | ) | const |
Compute and return the number of active descendants of this objects. For example, if all of the eight children of a hex are further refined isotropically exactly once, the returned number will be 64, not 80.
If the present cell is not refined, one is returned.
If one considers a triangulation as a forest where the root of each tree are the coarse mesh cells and nodes have descendants (the children of a cell), then this function returns the number of terminal nodes in the sub-tree originating from the current object; consequently, if the current object is not further refined, the answer is one.
unsigned int TriaAccessor< structdim, dim, spacedim >::max_refinement_depth | ( | ) | const |
Return the number of times that this object is refined. Note that not all its children are refined that often (which is why we prepend max_
), the returned number is rather the maximum number of refinement in any branch of children of this object.
For example, if this object is refined, and one of its children is refined exactly one more time, then max_refinement_depth
should return 2.
If this object is not refined (i.e. it is active), then the return value is zero.
TriaIterator<TriaAccessor<structdim, dim, spacedim> > TriaAccessor< structdim, dim, spacedim >::child | ( | const unsigned int | i | ) | const |
Return an iterator to the ith
child.
TriaIterator<TriaAccessor<structdim, dim, spacedim> > TriaAccessor< structdim, dim, spacedim >::isotropic_child | ( | const unsigned int | i | ) | const |
Return an iterator to that object that is identical to the ith child for isotropic refinement. If the current object is refined isotropically, then the returned object is the ith child. If the current object is refined anisotropically, the returned child may in fact be a grandchild of the object, or may not exist at all (in which case an exception is generated).
RefinementCase<structdim> TriaAccessor< structdim, dim, spacedim >::refinement_case | ( | ) | const |
Return the RefinementCase of this cell.
int TriaAccessor< structdim, dim, spacedim >::child_index | ( | const unsigned int | i | ) | const |
Index of the ith
child. The level of the child is one higher than that of the present cell, if the children of a cell are accessed. The children of faces have no level. If the child does not exist, -1 is returned.
int TriaAccessor< structdim, dim, spacedim >::isotropic_child_index | ( | const unsigned int | i | ) | const |
Index of the ith
isotropic child. See the isotropic_child() function for a definition of this concept. If the child does not exist, -1 is returned.
types::boundary_id TriaAccessor< structdim, dim, spacedim >::boundary_id | ( | ) | const |
Return the boundary indicator of this object.
If the return value is the special value numbers::internal_face_boundary_id, then this object is in the interior of the domain.
void TriaAccessor< structdim, dim, spacedim >::set_boundary_id | ( | const types::boundary_id | ) | const |
Set the boundary indicator of the current object. The same applies as for the boundary_id() function.
This function only sets the boundary object of the current object itself, not the indicators of the ones that bound it. For example, in 3d, if this function is called on a face, then the boundary indicator of the 4 edges that bound the face remain unchanged. If you want to set the boundary indicators of face and edges at the same time, use the set_all_boundary_ids() function. You can see the result of not using the correct function in the results section of step-49.
void TriaAccessor< structdim, dim, spacedim >::set_all_boundary_ids | ( | const types::boundary_id | ) | const |
Do as set_boundary_id() but also set the boundary indicators of the objects that bound the current object. For example, in 3d, if set_boundary_id() is called on a face, then the boundary indicator of the 4 edges that bound the face remain unchanged. In contrast, if you call the current function, the boundary indicators of face and edges are all set to the given value.
This function is useful if you set boundary indicators of faces in 3d (in 2d, the function does the same as set_boundary_id()) and you do so because you want a curved boundary object to represent the part of the boundary that corresponds to the current face. In that case, the Triangulation class needs to figure out where to put new vertices upon mesh refinement, and higher order Mapping objects also need to figure out where new interpolation points for a curved boundary approximation should be. In either case, the two classes first determine where interpolation points on the edges of a boundary face should be, asking the boundary object, before asking the boundary object for the interpolation points corresponding to the interior of the boundary face. For this to work properly, it is not sufficient to have set the boundary indicator for the face alone, but you also need to set the boundary indicators of the edges that bound the face. This function does all of this at once. You can see the result of not using the correct function in the results section of step-49.
bool TriaAccessor< structdim, dim, spacedim >::at_boundary | ( | ) | const |
Return whether this object is at the boundary. Obviously, the use of this function is only possible for dim>structdim
; however, for dim==structdim
, an object is a cell and the CellAccessor class offers another possibility to determine whether a cell is at the boundary or not.
const Manifold<dim, spacedim>& TriaAccessor< structdim, dim, spacedim >::get_manifold | ( | ) | const |
Return a constant reference to the manifold object used for this object.
As explained in the Manifold description for triangulations module, the process involved in finding the appropriate manifold description involves querying both the manifold or boundary indicators. See there for more information.
types::manifold_id TriaAccessor< structdim, dim, spacedim >::manifold_id | ( | ) | const |
Return the manifold indicator of this object.
If the return value is the special value numbers::flat_manifold_id, then this object is associated with a standard Cartesian Manifold Description.
bool TriaAccessor< structdim, dim, spacedim >::user_flag_set | ( | ) | const |
Read the user flag. See GlossUserFlags for more information.
void TriaAccessor< structdim, dim, spacedim >::set_user_flag | ( | ) | const |
Set the user flag. See GlossUserFlags for more information.
void TriaAccessor< structdim, dim, spacedim >::clear_user_flag | ( | ) | const |
Clear the user flag. See GlossUserFlags for more information.
void TriaAccessor< structdim, dim, spacedim >::recursively_set_user_flag | ( | ) | const |
Set the user flag for this and all descendants. See GlossUserFlags for more information.
void TriaAccessor< structdim, dim, spacedim >::recursively_clear_user_flag | ( | ) | const |
Clear the user flag for this and all descendants. See GlossUserFlags for more information.
void TriaAccessor< structdim, dim, spacedim >::clear_user_data | ( | ) | const |
Reset the user data to zero, independent if pointer or index. See GlossUserData for more information.
void TriaAccessor< structdim, dim, spacedim >::set_user_pointer | ( | void * | p | ) | const |
Set the user pointer to p
.
See GlossUserData for more information.
void TriaAccessor< structdim, dim, spacedim >::clear_user_pointer | ( | ) | const |
Reset the user pointer to a nullptr
pointer. See GlossUserData for more information.
void* TriaAccessor< structdim, dim, 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.
void TriaAccessor< structdim, dim, 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.
void TriaAccessor< structdim, dim, 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.
void TriaAccessor< structdim, dim, spacedim >::set_user_index | ( | const unsigned int | p | ) | const |
Set the user index to p
.
void TriaAccessor< structdim, dim, spacedim >::clear_user_index | ( | ) | const |
Reset the user index to 0. See GlossUserData for more information.
unsigned int TriaAccessor< structdim, dim, spacedim >::user_index | ( | ) | const |
Access the value of the user index.
See GlossUserData for more information.
void TriaAccessor< structdim, dim, 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.
void TriaAccessor< structdim, dim, 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.
double TriaAccessor< structdim, dim, spacedim >::diameter | ( | ) | const |
Diameter of the object.
The diameter of an object is computed to be the largest diagonal. This is not necessarily the true diameter for objects that may use higher order mappings, but completely sufficient for most computations.
std::pair<Point<spacedim>, double> TriaAccessor< structdim, dim, spacedim >::enclosing_ball | ( | ) | const |
Return a pair of Point and double corresponding to the center and the radius of a reasonably small enclosing ball of the object.
The function implements Ritter's O(n) algorithm to get a reasonably small enclosing ball around the vertices of the object. The initial guess for the enclosing ball is taken to be the ball which contains the largest diagonal of the object as its diameter. Starting from such an initial guess, the algorithm tests whether all the vertices of the object (except the vertices of the largest diagonal) are geometrically within the ball. If any vertex (v) is found to be geometrically outside the ball, a new iterate of the ball is constructed by shifting its center and increasing the radius so as to geometrically enclose both the previous ball and the vertex (v). The algorithm terminates when all the vertices are geometrically inside the ball.
If a vertex (v) is geometrically inside a particular iterate of the ball, then it will continue to be so in the subsequent iterates of the ball (this is true by construction).
see this and [Ritter 1990]
BoundingBox< spacedim > TriaAccessor< structdim, dim, spacedim >::bounding_box | ( | ) | const |
Return the smallest bounding box that encloses the object.
Notice that this method is not aware of any mapping you may be using to do your computations. If you are using a mapping object that modifies the position of the vertices, like MappingQEulerian, or MappingFEField, then you should call the function Mapping::get_bounding_box() instead.
Definition at line 1477 of file tria_accessor.cc.
double TriaAccessor< structdim, dim, spacedim >::extent_in_direction | ( | const unsigned int | axis | ) | const |
Length of an object in the direction of the given axis, specified in the local coordinate system. See the documentation of GeometryInfo for the meaning and enumeration of the local axes.
Note that the "length" of an object can be interpreted in a variety of ways. Here, we choose it as the maximal length of any of the edges of the object that are parallel to the chosen axis on the reference cell.
Definition at line 1499 of file tria_accessor.cc.
double TriaAccessor< structdim, dim, spacedim >::minimum_vertex_distance | ( | ) | const |
Return the minimal distance between any two vertices.
Point< spacedim > TriaAccessor< structdim, dim, spacedim >::intermediate_point | ( | const Point< structdim > & | coordinates | ) | const |
Return a point belonging to the Manifold<dim,spacedim> where this object lives, given its parametric coordinates on the reference structdim
cell. This function queries the underlying manifold object, and can be used to obtain the exact geometrical location of arbitrary points on this object.
Notice that the argument coordinates
are the coordinates on the reference cell, given in reference coordinates. In other words, the argument provides a weighting between the different vertices. For example, for lines, calling this function with argument Point<1>(.5), is equivalent to asking the line for its center.
Definition at line 1606 of file tria_accessor.cc.
Point< structdim > TriaAccessor< structdim, dim, spacedim >::real_to_unit_cell_affine_approximation | ( | const Point< spacedim > & | point | ) | const |
This function computes a fast approximate transformation from the real to the unit cell by inversion of an affine approximation of the \(d\)-linear function from the reference \(d\)-dimensional cell.
The affine approximation of the unit to real cell mapping is found by a least squares fit of an affine function to the \(2^d\) vertices of the present object. For any valid mesh cell whose geometry is not degenerate, this operation results in a unique affine mapping. Thus, this function will return a finite result for all given input points, even in cases where the actual transformation by an actual bi-/trilinear or higher order mapping might be singular. Besides only approximating the mapping from the vertex points, this function also ignores the attached manifold descriptions. The result is only exact in case the transformation from the unit to the real cell is indeed affine, such as in one dimension or for Cartesian and affine (parallelogram) meshes in 2D/3D.
For exact transformations to the unit cell, use Mapping::transform_real_to_unit_cell().
Definition at line 1724 of file tria_accessor.cc.
Point< spacedim > TriaAccessor< structdim, dim, spacedim >::center | ( | const bool | respect_manifold = false , |
const bool | interpolate_from_surrounding = false |
||
) | const |
Center of the object. The center of an object is defined to be the average of the locations of the vertices, which is also where a \(Q_1\) mapping would map the center of the reference cell. However, you can also ask this function to instead return the average of the vertices as computed by the underlying Manifold object associated with the current object, by setting to true the optional parameter respect_manifold
. Manifolds would then typically pull back the coordinates of the vertices to a reference domain (not necessarily the reference cell), compute the average there, and then push forward the coordinates of the averaged point to the physical space again; the resulting point is guaranteed to lie within the manifold, even if the manifold is curved.
When the object uses a different manifold description as its surrounding, like when part of the bounding objects of this TriaAccessor use a non-flat manifold description but the object itself is flat, the result given by the TriaAccessor::center() function may not be accurate enough, even when parameter respect_manifold
is set to true. If you find this to be case, than you can further refine the computation of the center by setting to true the second additional parameter interpolate_from_surrounding
. This computes the location of the center by a so-called transfinite interpolation from the center of all the bounding objects. For a 2D object, it puts a weight of 1/2
on each of the four surrounding lines and a weight -1/4
on the four vertices. This corresponds to a linear interpolation between the descriptions of the four faces, subtracting the contribution of the vertices that is added twice when coming through both lines adjacent to the vertex. In 3D, the weights for faces are 1/2
, the weights for lines are -1/4
, and the weights for vertices are 1/8
. For further information, also confer to the TransfiniteInterpolationManifold class that is able to not only apply this beneficial description to a single cell but all children of a coarse cell.
Definition at line 1753 of file tria_accessor.cc.
Point< spacedim > TriaAccessor< structdim, dim, spacedim >::barycenter | ( | ) | const |
Return the barycenter (also called centroid) of the object. The barycenter for an object \(K\) of dimension \(d\) in \(D\) space dimensions is given by the \(D\)-dimensional vector \(\mathbf x_K\) defined by
\[ \mathbf x_K = \frac{1}{|K|} \int_K \mathbf x \; \textrm{d}x \]
where the measure of the object is given by
\[ |K| = \int_K \mathbf 1 \; \textrm{d}x. \]
This function assumes that \(K\) is mapped by a \(d\)-linear function from the reference \(d\)-dimensional cell. Then the integrals above can be pulled back to the reference cell and evaluated exactly (if through lengthy and, compared to the center() function, expensive computations).
Definition at line 1455 of file tria_accessor.cc.
double TriaAccessor< structdim, dim, spacedim >::measure | ( | ) | const |
Compute the dim-dimensional measure of the object. For a dim-dimensional cell in dim-dimensional space, this equals its volume. On the other hand, for a 2d cell in 3d space, or if the current object pointed to is a 2d face of a 3d cell in 3d space, then the function computes the area the object occupies. For a one-dimensional object, return its length.
The function only computes the measure of cells, faces or edges assumed to be represented by (bi-/tri-)linear mappings. In other words, it only takes into account the locations of the vertices that bound the current object but not how the interior of the object may actually be mapped. In most simple cases, this is exactly what you want. However, for objects that are not "straight", e.g. 2d cells embedded in 3d space as part of a triangulation of a curved domain, two-dimensional faces of 3d cells that are not just parallelograms, or for faces that are at the boundary of a domain that is not just bounded by straight line segments or planes, this function only computes the dim-dimensional measure of a (bi-/tri-)linear interpolation of the "real" object as defined by the manifold or boundary object describing the real geometry of the object in question. If you want to consider the "real" geometry, you will need to compute the measure by integrating a function equal to one over the object, which after applying quadrature equals the summing the JxW values returned by the FEValues or FEFaceValues object you will want to use for the integral.
Definition at line 1466 of file tria_accessor.cc.
bool TriaAccessor< structdim, dim, spacedim >::is_translation_of | ( | const TriaIterator< TriaAccessor< structdim, dim, spacedim >> & | o | ) | const |
Return true if the current object is a translation of the given argument.
|
private |
Like set_boundary_id but without checking for internal faces or invalid ids.
|
private |
Copy the data of the given object into the internal data structures of a triangulation.
Definition at line 1444 of file tria_accessor.cc.
|
private |
Set the flag indicating, what line_orientation()
will return.
It is only possible to set the line_orientation of faces in 3d (i.e. structdim==2 && dim==3
).
|
private |
Set whether the quad with index face
has its normal pointing in the standard direction (true
) or whether it is the opposite (false
). Which is the standard direction is documented with the GeometryInfo class.
This function is only for internal use in the library. Setting this flag to any other value than the one that the triangulation has already set is bound to bring you disaster.
|
private |
Set the flag indicating, what face_flip()
will return.
It is only possible to set the face_orientation of cells in 3d (i.e. structdim==3 && dim==3
).
|
private |
Set the flag indicating, what face_rotation()
will return.
It is only possible to set the face_orientation of cells in 3d (i.e. structdim==3 && dim==3
).
|
private |
Set the used
flag. Only for internal use in the library.
|
private |
Clear the used
flag. Only for internal use in the library.
|
private |
Set the RefinementCase<dim>
this TriaObject is refined with. Not defined for structdim=1
as lines are always refined resulting in 2 children lines (isotropic refinement).
You should know quite exactly what you are doing if you touch this function. It is exclusively for internal use in the library.
|
private |
Clear the RefinementCase<dim> of this TriaObject, i.e. reset it to RefinementCase<dim>::no_refinement.
You should know quite exactly what you are doing if you touch this function. It is exclusively for internal use in the library.
|
private |
Set the index of the ith child. Since the children come at least in pairs, we need to store the index of only every second child, i.e. of the even numbered children. Make sure, that the index of child i=0 is set first. Calling this function for odd numbered children is not allowed.
|
private |
Clear the child field, i.e. set it to a value which indicates that this cell has no children.
double TriaAccessor< 2, 2, 2 >::extent_in_direction | ( | const unsigned int | axis | ) | const |
Lines along x-axis, see GeometryInfo
Lines along y-axis
Definition at line 1532 of file tria_accessor.cc.
double TriaAccessor< 2, 2, 3 >::extent_in_direction | ( | const unsigned int | axis | ) | const |
Lines along x-axis, see GeometryInfo
Lines along y-axis
Definition at line 1546 of file tria_accessor.cc.
double TriaAccessor< 3, 3, 3 >::extent_in_direction | ( | const unsigned int | axis | ) | const |
Lines along x-axis, see GeometryInfo
Lines along y-axis
Lines along z-axis
Definition at line 1561 of file tria_accessor.cc.