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 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)} \).
dim | the dimension of the parametrization domain, usually 2 (or 1 for beams in 2D, currently not implemented) |
Real | a 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... | |
using Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real >::Matrix = Dune::FieldMatrix<Real,dim+1,dim+1> |
Definition at line 53 of file shellKinematics.hh.
using Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real >::MatrixPair = std::pair<Matrix,Matrix> |
Definition at line 54 of file shellKinematics.hh.
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.
Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real >::ShellKinematicsDerivative | ( | ) |
Constructor defining a (pretty useless) default linearization point.
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.
MatrixPair Kaskade::Elastomechanics::ShellKinematicsDerivative< dim, Real >::d0 | ( | ) | const |
Evaluates \( (A,B) \) making up \( \Phi' = A+zB \).
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 \]
v | the parametrization variation direction in which to take the derivative |
i | for 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. |
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\} \]
v | the parametrization variation direction in which to take the first derivative |
i | for 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. |
w | the parametrization variation direction in which to take the second derivative |
j | for 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. |
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.