KASKADE 7 development version
|
#include <boundaryLocator.hh>
Supports finding boundary intersections by geometric queries.
GridView | |
Function | a vector valued function defined on the grid view. It has to support a subset of the FunctionSpaceElement interface, namely the methods order(cell) , value(cell,xi) , and derivative(cell,pos) . |
A displacement function defined on the given grid view can be provided. In that case, instead of the original grid the deformed grid, defined by the displacement, is considered.
Design rationale: Dune::GeometryGrid provides similar functionality, but incurs the (small) overhead of defining a complete grid structe and is restricted to piecewise linear deformations.
Explicit instantiations for Dune::UGGrid<2> and Dune::UGGrid<3> leaf grid views are provided in the file boundaryLocator.cpp, and are compiled into the library. For other grids views, include boundaryLocator.hpp, which contains the implementations.
Definition at line 75 of file boundaryLocator.hh.
Public Types | |
using | ctype = typename GridView::ctype |
using | Position = Dune::FieldVector< ctype, dimension_world > |
using | LocalPosition = Dune::FieldVector< ctype, dimension_world-1 > |
using | Face = std::conditional_t< dimension==dimension_world, typename GridView::Intersection, typename GridView::template Codim< 0 >::Entity > |
using | DisplacedFace = BoundaryFace< Grid, Face, Function, dimw > |
Public Member Functions | |
BoundaryLocator (GridView const &gridView) | |
Constructs a locator, optionally providing a displacement. More... | |
BoundaryLocator (GridView const &gridView, std::vector< bool > const &boundarySegments) | |
Constructor. More... | |
template<class... Args> | |
BoundaryLocator (GridView const &gridView, Args &&... args) | |
Constructs a locator, constructing the displacement function with the given argument. More... | |
~BoundaryLocator () | |
template<class Update > | |
void | updateDisplacement (Update const &update) |
Creates or updates the displacement function. More... | |
std::optional< std::tuple< typename BoundaryLocator< GridView, Function, dimw >::DisplacedFace, typename BoundaryLocator< GridView, Function, dimw >::LocalPosition, typename GridView::ctype > > | byLineSegment (Position const &from, Position const &to, ctype minAngle=-0.5) const |
Finds the first boundary face that is intersected by the given line segment from outside towards inside. More... | |
std::optional< std::tuple< typename BoundaryLocator< GridView, Function, dimw >::DisplacedFace, typename BoundaryLocator< GridView, Function, dimw >::LocalPosition, typename GridView::ctype > > | byRay (Position const &from, Position const &direction, ctype minAngle=-0.5) const |
Finds the first boundary face that is intersected by the given ray from outside towards inside. More... | |
ctype | diameter () const |
The diagonal length of the grid's cartesian bounding box. More... | |
Position | global (Face const &face, LocalPosition const &xi) const |
Returns the global position of a boundary point. More... | |
Position | unitOuterNormal (Face const &face, LocalPosition const &xi) const |
DEPRECATED use BoundaryFace::unitOuterNormal() instead. More... | |
Function const * | displacement () const |
Returns the current displacement of the mesh. More... | |
std::vector< Face > | faces () const |
Returns a sequence of all faces. More... | |
std::vector< DisplacedFace > | displacedFaces () const |
Returns a sequence of all deformed faces. More... | |
Static Public Attributes | |
static int const | dimension = GridView::Grid::dimension |
static int const | dimension_world = dimw |
using Kaskade::BoundaryLocator< GridView, Function, dimw >::ctype = typename GridView::ctype |
Definition at line 83 of file boundaryLocator.hh.
using Kaskade::BoundaryLocator< GridView, Function, dimw >::DisplacedFace = BoundaryFace<Grid,Face,Function,dimw> |
Definition at line 89 of file boundaryLocator.hh.
using Kaskade::BoundaryLocator< GridView, Function, dimw >::Face = std::conditional_t<dimension==dimension_world, typename GridView::Intersection, typename GridView::template Codim<0>::Entity> |
Definition at line 86 of file boundaryLocator.hh.
using Kaskade::BoundaryLocator< GridView, Function, dimw >::LocalPosition = Dune::FieldVector<ctype,dimension_world-1> |
Definition at line 85 of file boundaryLocator.hh.
using Kaskade::BoundaryLocator< GridView, Function, dimw >::Position = Dune::FieldVector<ctype,dimension_world> |
Definition at line 84 of file boundaryLocator.hh.
Kaskade::BoundaryLocator< GridView, Function, dimw >::BoundaryLocator | ( | GridView const & | gridView | ) |
Constructs a locator, optionally providing a displacement.
Kaskade::BoundaryLocator< GridView, Function, dimw >::BoundaryLocator | ( | GridView const & | gridView, |
std::vector< bool > const & | boundarySegments | ||
) |
Constructor.
boundarySegments | states which boundary faces (by their boundary segment index) are to be considered; if empty all are considered |
|
inline |
Constructs a locator, constructing the displacement function with the given argument.
Args | the argument type for the displacement function. A constructor Function(Args) must exist. |
Definition at line 111 of file boundaryLocator.hh.
Kaskade::BoundaryLocator< GridView, Function, dimw >::~BoundaryLocator | ( | ) |
std::optional< std::tuple< typename BoundaryLocator< GridView, Function, dimw >::DisplacedFace, typename BoundaryLocator< GridView, Function, dimw >::LocalPosition, typename GridView::ctype > > Kaskade::BoundaryLocator< GridView, Function, dimw >::byLineSegment | ( | Position const & | from, |
Position const & | to, | ||
ctype | minAngle = -0.5 |
||
) | const |
Finds the first boundary face that is intersected by the given line segment from outside towards inside.
"First" means that the intersection is closest to the start point of the line.
from | the origin of the line segment |
to | the end point of the line segment (shall not coincide with from) |
minAngle | the minimum angle for which an intersection is accepted. Only faces with unit outer normal \( n \) and \( n^T (\mathrm{to}-\mathrm{from}) / \| \mathrm{to}-\mathrm{from} \| \le \mathrm{minAngle} \) are considered. |
An (optional) triple (i,c,q) is returned. If valid, i contains the boundary face and c is the face-local coordinate of the intersection point. q is the scalar product between normalized ray direction and intersecting face unit outer normal.
Referenced by Kaskade::BoundaryLocator< GridView, Function, dimw >::byRay().
|
inline |
Finds the first boundary face that is intersected by the given ray from outside towards inside.
"First" means that the intersection is closest to the start point of the ray.
from | the origin of the line segment |
direction | the direction of the ray (shall not be zero) |
minAngle | the minimum angle for which an intersection is accepted. Only faces with unit outer normal \( n \) and \( n^T \mathrm{dir} / \| \mathrm{dir} \| \le \mathrm{minAngle} \) are considered. |
An (optional) triple (i,c,q) is returned. If valid, i contains the boundary face and c is the face-local coordinate of the intersection point. q is the scalar product between normalized ray direction and intersecting face unit outer normal.
Definition at line 187 of file boundaryLocator.hh.
Referenced by Kaskade::getContactConstraint().
|
inline |
The diagonal length of the grid's cartesian bounding box.
The diameter of the geometry comes handy if intersections with a ray instead of a line segment is desired. Taking an end point on the ray in diameter distance to the origin defines a line segment with all possible intersections.
Definition at line 199 of file boundaryLocator.hh.
std::vector< DisplacedFace > Kaskade::BoundaryLocator< GridView, Function, dimw >::displacedFaces | ( | ) | const |
Returns a sequence of all deformed faces.
The sequence of faces is extracted and returned by value, which is a considerable computational effort. The order of faces is unspecified, and subject to change on calling updateDisplacement().
Referenced by Kaskade::getContactArea(), and Kaskade::getContactConstraints().
|
inline |
Returns the current displacement of the mesh.
The returned pointer can be null, which signals a trivial zero displacement.
Definition at line 221 of file boundaryLocator.hh.
Referenced by Kaskade::contactLinesearch(), and Kaskade::getContactArea().
std::vector< Face > Kaskade::BoundaryLocator< GridView, Function, dimw >::faces | ( | ) | const |
Returns a sequence of all faces.
The sequence of faces is extracted and returned by value, which is a considerable computational effort. The order of faces is unspecified, and subject to change on calling updateDisplacement().
Referenced by Kaskade::BoundaryLocator< GridView, Function, dimw >::updateDisplacement().
Position Kaskade::BoundaryLocator< GridView, Function, dimw >::global | ( | Face const & | face, |
LocalPosition const & | xi | ||
) | const |
Returns the global position of a boundary point.
This respects a potentially given nontrivial displacement of the boundary.
Position Kaskade::BoundaryLocator< GridView, Function, dimw >::unitOuterNormal | ( | Face const & | face, |
LocalPosition const & | xi | ||
) | const |
DEPRECATED use BoundaryFace::unitOuterNormal() instead.
|
inline |
Creates or updates the displacement function.
Update | the argument type for constructing or updating the displacement. Function must have both a constructor Function(Update) and an assignment operator operator=(Update) . |
Calling this method recreates the spatial search tree, which can lead to a new tree structure. This affects the order of faces in the sequence returned by faces().
Definition at line 130 of file boundaryLocator.hh.
Referenced by Kaskade::contactLinesearch().
|
static |
Definition at line 81 of file boundaryLocator.hh.
|
static |
Definition at line 82 of file boundaryLocator.hh.