KASKADE 7 development version
|
A class providing information about domain boundary tangent vectors. More...
#include <boundaryInterpolation.hh>
A class providing information about domain boundary tangent vectors.
The surface normals are assumed to be continuous within each boundary face, but need not be continuous across non-feature boundary edges. Such a piecewise \( G^1 \)-continuous surface is determined by face-wise interpolation, where for each boundary edge the corresponding tangent vector of the smooth surface is provided.
Definition at line 523 of file boundaryInterpolation.hh.
Public Types | |
using | Vector = Dune::FieldVector< ctype, dim > |
Public Member Functions | |
template<class Grid > | |
BoundaryTangents (Grid const &grid, GeometricFeatureDetector< Index, ctype, dim > const &features) | |
Constructor. More... | |
Public Attributes | |
std::set< BoundaryEdge< Index, Vector > > | boundaryEdges2 |
Related Functions | |
(Note that these are not member functions.) | |
template<class GridView > | |
using | GridBoundaryTangents = BoundaryTangents< GridView::dimension, typename GridView::IndexSet::IndexType, typename GridView::ctype > |
A convenient type alias for BoundaryTangents if the grid view type is known. More... | |
using Kaskade::BoundaryTangents< dim, Index, ctype >::Vector = Dune::FieldVector<ctype,dim> |
Definition at line 526 of file boundaryInterpolation.hh.
|
inline |
Constructor.
Feature edges are determined according to two criteria:
featureEdges
list (the order of \( p_1, p_2 \) is not relevant).If any of these two criteria is satisfied, the edge is classified as feature edge.
If exactly two feature edges meet in one vertex, and their normalized emanating tangent vectors \( t_a, t_b \) are approximately parallel (i.e. \( \|t_a + t_b\| < \beta \), where \( \beta \) is the feature curve threshold), these two tangent vectors are modified such as to be collinear.
grid | |
features | an object telling whether an edge or a vertex are feature edges or vertices. The object need only exist during the constructor call. |
Definition at line 550 of file boundaryInterpolation.hh.
std::set<BoundaryEdge<Index,Vector> > Kaskade::BoundaryTangents< dim, Index, ctype >::boundaryEdges2 |
Definition at line 580 of file boundaryInterpolation.hh.
Referenced by Kaskade::BoundaryTangents< dim, Index, ctype >::BoundaryTangents().