KASKADE 7 development version
Public Types | Public Member Functions | List of all members
Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real > Class Template Reference

A lower-dimensional kinematics ansatz useful for shell models. More...

#include <shellKinematics.hh>

Detailed Description

template<int dim = 2, class Real = double>
class Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real >

A lower-dimensional kinematics ansatz useful for shell models.

This represents the spatial derivative \( \Phi' = d\Phi / d[x,z] \) of the parametrized deformation \( \Phi: \Omega\subset \R^{d+1} \to \R^{d+1} \) on \(\Omega = \omega \times \R \) defined as

\[ \Phi(x,z;p) = \phi(x;p) + z t(x;p) n(x;p) \]

with \( \phi, n: \omega\to\R^{d+1} \), \( t:\omega\to\R \), and \( n^T\phi' = 0\). For the shell case dim=2, \( n = \phi_{x_1} \times \phi_{x_2} \) holds.

The spatial derivative is of the form \( \Phi'(x,z;p) = A(x;p) + z B(x;p) \) with \( A,B\in \R^{(d+1)\times(d+1)} \).

Template Parameters
dimthe dimension of the parametrization domain, usually 2 (or 1 for beams in 2D, currently not implemented)
Reala real field type, usually double

Definition at line 45 of file shellKinematics.hh.

Public Types

using Vector = Dune::FieldVector< Real, dim+1 >
 The type of deformation vectors in \( \R^{d+1} \). More...
 
using Matrix = Dune::FieldMatrix< Real, dim+1, dim+1 >
 
using MatrixPair = std::pair< Matrix, Matrix >
 

Public Member Functions

 ShellKinematicsDerivative ()
 Constructor defining a (pretty useless) default linearization point. More...
 
 ShellKinematicsDerivative (Dune::FieldMatrix< Real, dim+1, dim > const &dphi_dx_, Tensor< Real, dim+1, dim, dim > const &ddphi_ddx_, Real t, Dune::FieldVector< Real, dim > const &dt_dx_)
 Constructor defining a linearization point. More...
 
void setLinearizationPoint (Dune::FieldMatrix< Real, dim+1, dim > const &dphi_dx_, Tensor< Real, dim+1, dim, dim > const &ddphi_ddx_, Real t, Dune::FieldVector< Real, dim > const &dt_dx_)
 Defines a new evaluation/linearization point. More...
 
MatrixPair d0 () const
 Evaluates \( (A,B) \) making up \( \Phi' = A+zB \). More...
 
MatrixPair d1 (VariationalArg< Real, dim, 1 > const &v, int i) const
 Computes the directional parametrization derivative of \( \Phi'(x,z;p) \). More...
 
MatrixPair d2 (VariationalArg< Real, dim, 1 > const &v, int i, VariationalArg< Real, dim, 1 > const &w, int j) const
 Computes the second directional parametrization derivative of \( \Phi' \). More...
 

Member Typedef Documentation

◆ Matrix

template<int dim = 2, class Real = double>
using Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real >::Matrix = Dune::FieldMatrix<Real,dim+1,dim+1>

Definition at line 53 of file shellKinematics.hh.

◆ MatrixPair

template<int dim = 2, class Real = double>
using Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real >::MatrixPair = std::pair<Matrix,Matrix>

Definition at line 54 of file shellKinematics.hh.

◆ Vector

template<int dim = 2, class Real = double>
using Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real >::Vector = Dune::FieldVector<Real,dim+1>

The type of deformation vectors in \( \R^{d+1} \).

Definition at line 51 of file shellKinematics.hh.

Constructor & Destructor Documentation

◆ ShellKinematicsDerivative() [1/2]

template<int dim = 2, class Real = double>
Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real >::ShellKinematicsDerivative ( )

Constructor defining a (pretty useless) default linearization point.

◆ ShellKinematicsDerivative() [2/2]

template<int dim = 2, class Real = double>
Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real >::ShellKinematicsDerivative ( Dune::FieldMatrix< Real, dim+1, dim > const &  dphi_dx_,
Tensor< Real, dim+1, dim, dim > const &  ddphi_ddx_,
Real  t,
Dune::FieldVector< Real, dim > const &  dt_dx_ 
)

Constructor defining a linearization point.

Member Function Documentation

◆ d0()

template<int dim = 2, class Real = double>
MatrixPair Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real >::d0 ( ) const

Evaluates \( (A,B) \) making up \( \Phi' = A+zB \).

◆ d1()

template<int dim = 2, class Real = double>
MatrixPair Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real >::d1 ( VariationalArg< Real, dim, 1 > const &  v,
int  i 
) const

Computes the directional parametrization derivative of \( \Phi'(x,z;p) \).

\[ \frac{\partial (A,B)}{\partial p} v \]

Parameters
vthe parametrization variation direction in which to take the derivative
ifor i=0,...,dim-1, the derivative wrt to the i-th component of \( \phi \) is taken. For i=dim, the derivative wrt to \( t \) is computed.

◆ d2()

template<int dim = 2, class Real = double>
MatrixPair Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real >::d2 ( VariationalArg< Real, dim, 1 > const &  v,
int  i,
VariationalArg< Real, dim, 1 > const &  w,
int  j 
) const

Computes the second directional parametrization derivative of \( \Phi' \).

\[ \frac{\partial^2 (A,B)}{\partial p_i \partial p_j} [v,w], \quad p_i, p_j\in\{\phi,t\} \]

Parameters
vthe parametrization variation direction in which to take the first derivative
ifor i=0,...,dim-1, the derivative wrt to the i-th component of \( \phi \) is taken. For i=dim, the derivative wrt to \( t \) is computed.
wthe parametrization variation direction in which to take the second derivative
jfor j=0,...,dim-1, the derivative wrt to the j-th component of \( \phi \) is taken. For j=dim, the derivative wrt to \( t \) is computed.

◆ setLinearizationPoint()

template<int dim = 2, class Real = double>
void Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real >::setLinearizationPoint ( Dune::FieldMatrix< Real, dim+1, dim > const &  dphi_dx_,
Tensor< Real, dim+1, dim, dim > const &  ddphi_ddx_,
Real  t,
Dune::FieldVector< Real, dim > const &  dt_dx_ 
)

Defines a new evaluation/linearization point.


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