KASKADE 7 development version
Public Types | Public Member Functions | List of all members
Kaskade::Mortar< Scalar, dim > Class Template Referenceabstract

Abstract base class for computation of constraints from constraint samples. More...

#include <contactConstraints.hh>

Detailed Description

template<class Scalar, int dim>
class Kaskade::Mortar< Scalar, dim >

Abstract base class for computation of constraints from constraint samples.

Let \( x \) denote the vector of displacements. On a boundary face, a set of \( n \) constraint samples is defined at local positions \( \xi_i \), leading to (pointwise) nonpenetration conditions of the form

\[ \sum_{k\in K_i} G_{ik} x_k \le g_i \]

where \( g_i \) is the gap between opposing boundary points and \( G_i \) is the impact of the displacements \( x_k \) on this gap, i.e. \( G_{ik} = \phi_k(\xi_i) \), where \( \phi_k \) is the shape function associated with \( x_k \). Note that the entries \( x_k \) may be (column) vector-valued, and the corresponding matrix entries \( G_{ik} \) are (row) vector-valued.

Now the actual constraints of the form \( Bx \le b \) can be defined by certain linear combinations of the pointwise constraints. Taking them one by one as they are results in pointwise constraints (with Dirac weights),

\[ Gx \le g, \]

see DiracMortar. If the pointwise constraints are multiplied by test functions \( \psi_j \) and integrated (summed) over the boundary face, different constraints result (see, e.g., BezierMortar):

\[ \forall j: \quad \sum_{i\in I} \psi_j(\xi_i) \sum_{k\in K_i} G_{ik} x_k \le \sum_{i\in I} \psi_j(\xi_i) g_i. \]

In any case, linear inequalities of the form \( Bx \le b \) are defined.

This class defines an interface for computing \( B \) and \( b \) from \( G \) and \( g \).

Definition at line 265 of file contactConstraints.hh.

Inheritance diagram for Kaskade::Mortar< Scalar, dim >:
Kaskade::BezierMortar< Scalar, dim > Kaskade::DiracMortar< Scalar, dim >

Public Types

using Entry = Dune::FieldMatrix< Scalar, 1, dim >
 The type of matrix entries in the constraint matrix \( B \) in \( Bx \le b \). More...
 
using Row = std::vector< std::pair< size_t, Entry > >
 

Public Member Functions

virtual ~Mortar ()
 
virtual int size (int n) const =0
 The number of resulting constraints if n samples are provided. More...
 
virtual void updateConstraints (int n, Dune::FieldVector< Scalar, dim-1 > const &xi, Row &Gi, Scalar gi, std::vector< Row > &B, std::vector< Scalar > &b) const =0
 Updates the given set of constraints (by \( B \) and \( b \)) when a new constraint sample becomes available. More...
 

Member Typedef Documentation

◆ Entry

template<class Scalar , int dim>
using Kaskade::Mortar< Scalar, dim >::Entry = Dune::FieldMatrix<Scalar,1,dim>

The type of matrix entries in the constraint matrix \( B \) in \( Bx \le b \).

The entries are row vectors with the same length as the primal variable \( x \) entries.

Definition at line 274 of file contactConstraints.hh.

◆ Row

template<class Scalar , int dim>
using Kaskade::Mortar< Scalar, dim >::Row = std::vector<std::pair<size_t,Entry> >

Definition at line 275 of file contactConstraints.hh.

Constructor & Destructor Documentation

◆ ~Mortar()

template<class Scalar , int dim>
virtual Kaskade::Mortar< Scalar, dim >::~Mortar ( )
virtual

Member Function Documentation

◆ size()

template<class Scalar , int dim>
virtual int Kaskade::Mortar< Scalar, dim >::size ( int  n) const
pure virtual

The number of resulting constraints if n samples are provided.

Implemented in Kaskade::DiracMortar< Scalar, dim >, and Kaskade::BezierMortar< Scalar, dim >.

◆ updateConstraints()

template<class Scalar , int dim>
virtual void Kaskade::Mortar< Scalar, dim >::updateConstraints ( int  n,
Dune::FieldVector< Scalar, dim-1 > const &  xi,
Row Gi,
Scalar  gi,
std::vector< Row > &  B,
std::vector< Scalar > &  b 
) const
pure virtual

Updates the given set of constraints (by \( B \) and \( b \)) when a new constraint sample becomes available.

This method needs to be implemented in derived classes.

Parameters
nThe total number of samples that may be provided (actually there may be fewer, because some points may lack a contact partner).
xiThe face-local coordinate of the current contact sample.
GiA sparse representation of the matrix row \( G_i \).
giThe gap \( g_i \)
BThe sparse matrix \( B \) to be updated. The update process can modify values as well as change the shape of the matrix.
bThe bounds vector \( b \)

If B and b are empty on entry, a new face is started. The containers are resized as needed.

Implemented in Kaskade::DiracMortar< Scalar, dim >, and Kaskade::BezierMortar< Scalar, dim >.


The documentation for this class was generated from the following file: