13#ifndef SHELLKINEMATICS
14#define SHELLKINEMATICS
16#include "dune/common/fmatrix.hh"
17#include "dune/common/fvector.hh"
44 template <
int dim=2,
class Real=
double>
126 template <
int dim=2,
class Real=
double>
195 template <
class Energy=MaterialLaws::StVenantKirchhoff<3,
double>,
203 static int const dim = Energy::dim-1;
254 template <
int row,
int col>
265 Matrix A, B, E, BtB, BtA, BtA_AtB;
268 Real
d1(std::pair<Matrix,Matrix>
const& dAdB)
const;
273 Real
d2(std::pair<Matrix,Matrix>
const& dAdB1, std::pair<Matrix,Matrix>
const& dAdB2,
274 std::pair<Matrix,Matrix>
const& ddAddB)
const;
A stored energy density class for Taylor-based shell models.
Dune::FieldMatrix< Real, row==0?dim+1:1, col==0?dim+1:1 > d2(VariationalArg< Real, dim, 1 > const &v, VariationalArg< Real, dim, 1 > const &w) const
Evaulates the second parametric directional derivative of the marginal stored energy density wrt cha...
static int const dim
The dimension of the parameter domain.
Real d0() const
Evaluates the marginal stored energy density .
Dune::FieldVector< Real, row==0?dim+1:1 > d1(VariationalArg< Real, dim, 1 > const &v) const
Evaulates the parametric directional derivative of the marginal stored energy density wrt changes....
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.
ShellEnergy(Energy const &W, Real tau)
ShellEnergy(Energy const &W, Real tau, 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.
A lower-dimensional kinematics ansatz useful for shell models.
MatrixPair d1(VariationalArg< Real, dim, 1 > const &v, int i) const
Computes the directional parametrization derivative of .
ShellKinematicsDerivative()
Constructor defining a (pretty useless) default linearization point.
MatrixPair d0() const
Evaluates making up .
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.
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.
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 .
std::pair< Matrix, Matrix > MatrixPair
A lower-dimensional kinematics ansatz useful for shell models.
ShellKinematicsDerivative< dim, Real > derivative() const
Returns an object for evaluating and its parametrization derivatives.
auto d1(VariationalArg< Real, dim, 1 > const &v) const
Computes the directional parametrization derivative.
ShellKinematics(Vector const &phi, Dune::FieldMatrix< Real, dim+1, dim > const &dphi_dx, Real t)
A class for representing tensors of arbitrary static rank and extents.
This file contains various utility functions that augment the basic functionality of Dune.
A class that stores values, gradients, and Hessians of evaluated shape functions.