KASKADE 7 development version
|
A lower-dimensional kinematics ansatz useful for shell models. More...
#include <shellKinematics.hh>
A lower-dimensional kinematics ansatz useful for shell models.
This represents the deformation \( \Phi: \Omega\subset \R^{d+1} \to \R^{d+1} \) on \(\Omega = \omega \times \R \) as
\[ \Phi(x,z) = \phi(x) + z t(x) n(x) \]
with \( \phi, n: \omega\to\R^{d+1} \), \( t:\omega\to\R \), and \( n^T\phi' = 0\).
dim | the dimension of the parametrization domain, usually 2 (or 1 for beams in 2D) |
Real | a real field type, usually double |
Definition at line 127 of file shellKinematics.hh.
Public Types | |
using | Vector = Dune::FieldVector< Real, dim+1 > |
The type of deformation vectors in \( \R^{d+1} \). More... | |
Public Member Functions | |
ShellKinematics (Vector const &phi, Dune::FieldMatrix< Real, dim+1, dim > const &dphi_dx, Real t) | |
ShellKinematicsDerivative< dim, Real > | derivative () const |
Returns an object for evaluating \( \Phi' \) and its parametrization derivatives. More... | |
Vector | d0 () const |
template<int row> | |
auto | d1 (VariationalArg< Real, dim, 1 > const &v) const |
Computes the directional parametrization derivative. More... | |
using Kaskade::Elastomechanics::ShellKinematics< dim, Real >::Vector = Dune::FieldVector<Real,dim+1> |
The type of deformation vectors in \( \R^{d+1} \).
Definition at line 133 of file shellKinematics.hh.
Kaskade::Elastomechanics::ShellKinematics< dim, Real >::ShellKinematics | ( | Vector const & | phi, |
Dune::FieldMatrix< Real, dim+1, dim > const & | dphi_dx, | ||
Real | t | ||
) |
|
inline |
Definition at line 142 of file shellKinematics.hh.
auto Kaskade::Elastomechanics::ShellKinematics< dim, Real >::d1 | ( | VariationalArg< Real, dim, 1 > const & | v | ) | const |
Computes the directional parametrization derivative.
\[ \frac{\partial\Phi}{\partial p} v, \quad p\in\{\phi,t\} \]
row | for row =0, the derivative wrt \( \phi \) is evaluated, for row =1 wrt \( t \) |
ShellKinematicsDerivative< dim, Real > Kaskade::Elastomechanics::ShellKinematics< dim, Real >::derivative | ( | ) | const |
Returns an object for evaluating \( \Phi' \) and its parametrization derivatives.