KASKADE 7 development version
Classes | Functions

Components for defining elastomechanical contact problems. More...

Classes

struct  Kaskade::MortarB
 Defines the sample-based contact constraint formulation. More...
 

Functions

template<class Space , class DisplacedFace , class Displacement >
auto Kaskade::getContactConstraint (Space const &space, BoundaryLocator< typename Space::GridView, Displacement > const &boundaryLocator, DisplacedFace const &face, Dune::FieldVector< typename Space::GridView::ctype, Space::GridView::Grid::dimension-1 > const &xi, double overlap)
 Computes a single multibody contact point constraint. More...
 
template<class Space1 , class DisplacedFace1 , class Displacement1 , int dimw1 = Space1::GridView::Grid::dimension, class Space2 , class Displacement2 , int dimw2 = Space2::GridView::Grid::dimension>
auto Kaskade::getContactConstraint (Space1 const &space1, BoundaryLocator< typename Space1::GridView, Displacement1, Space1::GridView::Grid::dimension > const &boundaryLocator1, Space2 const &space2, BoundaryLocator< typename Space2::GridView, Displacement2, Space2::GridView::Grid::dimension > const &boundaryLocator2, DisplacedFace1 const &face, Dune::FieldVector< typename Space1::GridView::ctype, Space1::GridView::Grid::dimension-1 > const &xi, double overlap, int which=0)
 Computes a single multibody contact point constraint. More...
 
template<int dim>
std::vector< Dune::QuadraturePoint< double, dim > > Kaskade::constraintPositions (Dune::GeometryType gt, int nodesPerEdge)
 Computes contact sample points on faces. More...
 
template<class Space , class Displacement >
auto Kaskade::getContactConstraints (Space const &space, BoundaryLocator< typename Space::GridView, Displacement > const &boundaryLocator, int pointsPerEdge, double overlap=std::numeric_limits< double >::lowest(), Mortar< typename Space::Grid::ctype, Space::Grid::dimension > const &mortar=ContactConstraintsDetail::defaultMortar< typename Space::Grid::ctype, Space::Grid::dimension >())
 Computes a complete set of pointwise multibody contact constraints. More...
 
template<class Space1 , class Displacement1 , class Space2 , class Displacement2 , bool flattened = false>
auto Kaskade::getContactConstraints (Space1 const &space1, BoundaryLocator< typename Space1::GridView, Displacement1, Space1::GridView::Grid::dimension > const &boundaryLocator1, Space2 const &space2, BoundaryLocator< typename Space2::GridView, Displacement2, Space1::GridView::Grid::dimension > const &boundaryLocator2, int order, double overlap=std::numeric_limits< double >::lowest(), bool symmetric=true)
 Computes a complete pointwise twobody contact constraint. More...
 
template<class GridView , class Displacement >
double Kaskade::getContactArea (BoundaryLocator< GridView, Displacement > const &boundaryLocator1, BoundaryLocator< GridView, Displacement > const &boundaryLocator2, int order, double overlap=std::numeric_limits< double >::lowest(), double maxGap=1e-4)
 Computes the contact area for the first body. More...
 
template<class Displacement , class BoundaryDisplacement >
double Kaskade::contactLinesearch (Displacement const &u, BoundaryLocator< typename Displacement::GridView, BoundaryDisplacement > &boundaryLocator, Displacement const &direction, int pointsPerEdge, double trialStepSize=0.2, double overlap=std::numeric_limits< double >::lowest(), double targetOverlap=0, Mortar< typename Displacement::ctype, Displacement::GridView::dimension > const &mortar=ContactConstraintsDetail::defaultMortar< typename Displacement::ctype, Displacement::GridView::dimension >())
 Computes a feasible step size for contact constraints. More...
 
template<class Displacement1 , class BoundaryDisplacement1 , class Displacement2 , class BoundaryDisplacement2 >
double Kaskade::contactLinesearch (Displacement1 const &u1, BoundaryLocator< typename Displacement1::GridView, BoundaryDisplacement1 > &boundaryLocator1, Displacement1 const &direction1, Displacement2 const &u2, BoundaryLocator< typename Displacement2::GridView, BoundaryDisplacement2 > &boundaryLocator2, Displacement2 const &direction2, int order, double trialStepSize=0.2, double overlap=std::numeric_limits< double >::lowest(), double targetOverlap=0)
 Computes a feasible step size for contact constraints. More...
 

Detailed Description

Components for defining elastomechanical contact problems.

Function Documentation

◆ constraintPositions()

template<int dim>
std::vector< Dune::QuadraturePoint< double, dim > > Kaskade::constraintPositions ( Dune::GeometryType  gt,
int  nodesPerEdge 
)

Computes contact sample points on faces.

Equidistant contact sample points are computed.

Parameters
gtthe type of face geometry, usually line or triangle (only those are supported right now)
nodesPerEdgethe number of sample points per one-dimensional face subentity

For lines (2D contact), \( n \) sample points are created, whereas for triangles (3D contact) \( n(n+1)/2 \) points are created.

Returns
a vector of sample points with associated (constant) weights

Definition at line 428 of file contactConstraints.hh.

Referenced by Kaskade::getContactArea(), and Kaskade::getContactConstraints().

◆ contactLinesearch() [1/2]

template<class Displacement , class BoundaryDisplacement >
double Kaskade::contactLinesearch ( Displacement const &  u,
BoundaryLocator< typename Displacement::GridView, BoundaryDisplacement > &  boundaryLocator,
Displacement const &  direction,
int  pointsPerEdge,
double  trialStepSize = 0.2,
double  overlap = std::numeric_limits<double>::lowest(),
double  targetOverlap = 0,
Mortar< typename Displacement::ctype, Displacement::GridView::dimension > const &  mortar = ContactConstraintsDetail::defaultMortar<typename Displacement::ctype,Displacement::GridView::dimension>() 
)

Computes a feasible step size for contact constraints.

The contact constraints as, e.g., obtained from getContactConstraints, are linearized, in that they restrict the displacement of boundary points in normal direction. A solution of this linearized contact problem may be infeasible: other body parts may not have been detected as obstacles in normal direction, but now become relevant. Taking a full step often leads to infeasible deformation states with considerable interpenetration. Reducing the stepsize appropriately, i.e. taking the largest possible step such that the configuration remains feasible not only in the linearized, but in the full nonlinear setting, is a promising way.

This function computes a feasible step size that is nearly optimal.

Template Parameters
Displacementa finite element function type for representing volume displacements
BoundaryDisplacementa boundary displacement
Parameters
uthe volume displacement, starting point for line search
boundaryLocatora spatial index for looking up ray-face intersections. Its internal boundary displacement will be modified during the line search
directionthe step to be taken (linearized contact solution)
pointsPerEdgethe number of contact sample points per edge
trialStepSizehow aggressively the method goes ahead. 1=fast but with risk of ending up infeasible, 0=very slow but safe. In the range ]0,1[. Use smaller values for highly nonlinear problems with large displacements.
overlapthe amount of penetration to be accepted when looking backwards for contact

Definition at line 855 of file contactConstraints.hh.

◆ contactLinesearch() [2/2]

template<class Displacement1 , class BoundaryDisplacement1 , class Displacement2 , class BoundaryDisplacement2 >
double Kaskade::contactLinesearch ( Displacement1 const &  u1,
BoundaryLocator< typename Displacement1::GridView, BoundaryDisplacement1 > &  boundaryLocator1,
Displacement1 const &  direction1,
Displacement2 const &  u2,
BoundaryLocator< typename Displacement2::GridView, BoundaryDisplacement2 > &  boundaryLocator2,
Displacement2 const &  direction2,
int  order,
double  trialStepSize = 0.2,
double  overlap = std::numeric_limits<double>::lowest(),
double  targetOverlap = 0 
)

Computes a feasible step size for contact constraints.

The contact constraints as, e.g., obtained from getContactConstraints, are linearized, in that they restrict the displacement of boundary points in normal direction. A solution of this linearized contact problem may be infeasible: other body parts may not have been detected as obstacles in normal direction, but now become relevant. Taking a full step often leads to infeasible deformation states with considerable interpenetration. Reducing the stepsize appropriately, i.e. taking the largest possible step such that the configuration remains feasible not only in the linearized, but in the full nonlinear setting, is a promising way.

This function computes a feasible step size that is nearly optimal.

This variant is using two BoundaryLocators to determine contact between objects discretized with two individual grids.

Template Parameters
Displacement1a FunctionSpaceElement (first object)
BoundaryDisplacement1a boundary displacement (first object)
Displacement2a FunctionSpaceElement (second object)
BoundaryDisplacement2a boundary displacement (second object)
Parameters
u1the current displacement of the first object, starting point for the step
boundaryLocator1the spatial index for looking up ray-face intersections, fitst object
direction1the step to be taken by object 1 (linearized contact solution)
u2the current displacement of the second object
boundaryLocator2the spatial index for looking up ray-face intersections, second object.
direction2the step to be taken by object 2 (linearized contact solution)
orderthe order of the quadrature rule that is used for selecting test points in the faces
trialStepSizehow aggressively the method goes ahead. 1=fast but with risk of ending up infeasible, 0=very slow but safe. In the range ]0,1[. Use smaller values for highly nonlinear problems with large displacements.
overlapthe amount of infeasibility that is to be accepted.

Definition at line 987 of file contactConstraints.hh.

◆ getContactArea()

template<class GridView , class Displacement >
double Kaskade::getContactArea ( BoundaryLocator< GridView, Displacement > const &  boundaryLocator1,
BoundaryLocator< GridView, Displacement > const &  boundaryLocator2,
int  order,
double  overlap = std::numeric_limits<double>::lowest(),
double  maxGap = 1e-4 
)

Computes the contact area for the first body.

The threshold value maxGap is not so easy to choose, maybe this is not the best way to compute the contact area.

Parameters
boundaryLocator1the spatial index of boundary faces of first body/grid
boundaryLocator2the spatial index of boundary faces of second body/grid
orderthe order of the quadrature rule that is used for selecting test points in the faces
overlapstates how far into the interior for an opposed face is searched. Negative value means search is only performed in the exterior with corresponding minimal distance. If not given a default value (small positive value) will be chosen.
maxGapthreshold value, distances below this value are treated as "in contact".
Returns
the area

Definition at line 786 of file contactConstraints.hh.

◆ getContactConstraint() [1/2]

template<class Space , class DisplacedFace , class Displacement >
auto Kaskade::getContactConstraint ( Space const &  space,
BoundaryLocator< typename Space::GridView, Displacement > const &  boundaryLocator,
DisplacedFace const &  face,
Dune::FieldVector< typename Space::GridView::ctype, Space::GridView::Grid::dimension-1 > const &  xi,
double  overlap 
)

Computes a single multibody contact point constraint.

For a given boundary point \( x \), the outer normal ray (in direction of the unit outer normal vector \( n\)) is tested for intersection with other boundary faces. More precisely, intersections are considered on the ray \( \{ x+tn \mid t\in[-\mathrm{overlap},\infty\mathclose[\} \). Positive values of overlap allow to detect contact partners in an infeasible, already overlapping configuration, which is common for penalty methods or may arise due to nonlinearity.

If an intersection is found at point \( y \), the (linearized) non-penetration condition for a displacement \( \delta u \) is

\[ n^T (\delta u(x)-\delta u(y)) \le |x-y|. \]

With a basis representation of \( \delta u = \sum_j \delta u_j \varphi_j\), this leads to a scalar algebraic constraint of the form

\[ \sum_j n^T(\varphi_j(x) -\varphi_j(y)) \delta u_j \le |x-y|. \]

Let \( J(z) := \{ j\in \mathbf{N} \mid z \in \mathrm{supp}\varphi_j \} \) denote the dofs associated to \( z \). The return value is the (optional) pair

\[ \left( \{ (j,n^T(\varphi_j(x) -\varphi_j(y))) \mid j \in J(x)\cup J(y) \}, |x-y| \right), \]

which can establish one row of a complete contact constraint of the form \( B\delta u \le g \).

Returns
an optional pair (vector[(i,c)],g).

In order to support nonlinear contact problems (with finite displacement), we consider an already deformed domain given in terms of a displacement \( u \), that is provided by the boundary locator. This displacement \( u \) can, but need not be a finite element function.

Parameters
spacethe FEFunctionSpace for which the constraints are to be computed
boundaryLocatora spatial index for looking up boundary intersections with rays
facethe current boundary face
xithe local coordinate in the face
overlapstates how far into the interior for an opposed face is searched. If std::numeric_limits<double>::lowest() is provided, a default value proportional to the cell diameter is chosen.

Definition at line 69 of file contactConstraints.hh.

◆ getContactConstraint() [2/2]

template<class Space1 , class DisplacedFace1 , class Displacement1 , int dimw1 = Space1::GridView::Grid::dimension, class Space2 , class Displacement2 , int dimw2 = Space2::GridView::Grid::dimension>
auto Kaskade::getContactConstraint ( Space1 const &  space1,
BoundaryLocator< typename Space1::GridView, Displacement1, Space1::GridView::Grid::dimension > const &  boundaryLocator1,
Space2 const &  space2,
BoundaryLocator< typename Space2::GridView, Displacement2, Space2::GridView::Grid::dimension > const &  boundaryLocator2,
DisplacedFace1 const &  face,
Dune::FieldVector< typename Space1::GridView::ctype, Space1::GridView::Grid::dimension-1 > const &  xi,
double  overlap,
int  which = 0 
)

Computes a single multibody contact point constraint.

This version is for problems where each body is represented by its own grid.

Parameters
space1the first grid's displacement FE-space with respect to which the contact constraint is to be computed.
boundaryLocator1boundary geometry for current point
space2the second grid's displacement FE-space with respect to which the contact constraint is to be computed
boundaryLocator2a spatial index for looking up boundary intersections with rays
facethe current boundary face
xithe local coordinate in the face
overlapstates how far into the interior for an opposed face is searched. Negative value means search is only performed in the exterior with corresponding minimal distance.

Definition at line 150 of file contactConstraints.hh.

Referenced by Kaskade::getContactArea(), and Kaskade::getContactConstraints().

◆ getContactConstraints() [1/2]

template<class Space , class Displacement >
auto Kaskade::getContactConstraints ( Space const &  space,
BoundaryLocator< typename Space::GridView, Displacement > const &  boundaryLocator,
int  pointsPerEdge,
double  overlap = std::numeric_limits<double>::lowest(),
Mortar< typename Space::Grid::ctype, Space::Grid::dimension > const &  mortar = ContactConstraintsDetail::defaultMortar<typename Space::Grid::ctype,Space::Grid::dimension>() 
)

Computes a complete set of pointwise multibody contact constraints.

This computes a linear inequality constraint of the form \( B \delta u \le b \). If satisfied by some small displacement \( \delta u \), this guarantees (up to linearization and discretization) mutual nonpenetration and self-nonpenetration.

Parameters
boundaryLocatorthe spatial index of boundary faces
pointsPerEdgethe number of points per edge of the boundary face
overlapstates how far into the interior for an opposed face is searched. If not given, a default value (small positive value) will be chosen in dependence of the face size.
mortarthe type of constraint agglomeration/averaging
Returns
a pair \( (B,b) \) of constraint matrix and right hand side
See also
getContactConstraint

Definition at line 504 of file contactConstraints.hh.

◆ getContactConstraints() [2/2]

template<class Space1 , class Displacement1 , class Space2 , class Displacement2 , bool flattened = false>
auto Kaskade::getContactConstraints ( Space1 const &  space1,
BoundaryLocator< typename Space1::GridView, Displacement1, Space1::GridView::Grid::dimension > const &  boundaryLocator1,
Space2 const &  space2,
BoundaryLocator< typename Space2::GridView, Displacement2, Space1::GridView::Grid::dimension > const &  boundaryLocator2,
int  order,
double  overlap = std::numeric_limits<double>::lowest(),
bool  symmetric = true 
)

Computes a complete pointwise twobody contact constraint.

This version is for problems where each body is represented by its own grid.

Parameters
boundaryLocator1the spatial index of boundary faces of first body/grid
boundaryLocator2the spatial index of boundary faces of second body/grid
orderthe order of the quadrature rule that is used for selecting test points in the faces
overlapstates how far into the interior for an opposed face is searched.
symmetricif true, this loops over both grids and collects pointwise constraints. If false, only boundary faces of first grid are considered Negative value means search is only performed in the exterior with corresponding minimal distance. If not given a default value (small positive value) will be chosen.
space1if given contact constraints are computed with respect to this space, otherwise the space from boundaryLocator1 is used
space2analogously to space1 but for the second grid
Returns
a three tuple \( (B1,B2,b) \) of constraint matrix blocks (belonging to grids) and right hand side

Definition at line 605 of file contactConstraints.hh.

Referenced by Kaskade::contactLinesearch().