KASKADE 7 development version
Public Member Functions | Static Public Attributes | List of all members
Kaskade::Elastomechanics::ShellEnergy< Energy, Real > Class Template Reference

A stored energy density class for Taylor-based shell models. More...

#include <shellKinematics.hh>

Detailed Description

template<class Energy = MaterialLaws::StVenantKirchhoff<3,double>, class Real = double>
class Kaskade::Elastomechanics::ShellEnergy< Energy, Real >

A stored energy density class for Taylor-based shell models.

Assume on the reference configuration \(\Omega = \omega\times\mathopen]-\tau,\tau\mathclose[ \) with \( \omega\subset\R^d \) is a deformation \( \Phi:\Omega\to\R^{d+1}\) defined and there is a stored energy density \( W: \R^{(d+1)\times(d+1)}_{\rm sym} \to \R \) formulated in terms of the Euler-Lagrange strain tensor \( E = \frac{1}{2}(C-I) \) with \( C = \Phi'^T \Phi' \). Integrating over \( t\in \mathopen]-\tau,\tau\mathclose[ \), we define the marginal stored energy density

\[ w:\omega\to\R, \quad w(x) = \int_{-\tau}^\tau W(E(x,t)) \, dt. \]

For Kirchhoff shells ( \( d=2\)) with \( \Phi(x,t) = \phi(x) + t(x)n(x) \) and \( n(x) = \phi_{x_0}(x) \times \phi_{x_1}(x) \), the derivative is of the form \( \Phi'(x,t) = A(x)+tB(x) \). Then, the marginal stored energy is

\[ w(x) = 2\tau W(E(x,0)) + \frac{\tau^3}{3} W'(E(x,0))[B^TB] + \frac{\tau^3}{12} W''(E(x,0))[B^TA+A^TB,B^TA+A^TB] + \mathcal{O}(\tau^5), \]

where \( E(x,0) = \frac{1}{2} (A(x)^TA(x)-I) \). The first order term represents in-plane strain energy contributions, and the third order term bending energy. If \( \tau \) is small, the remainder term of order \( \tau^5 \) can be neglected.

This class implements this marginal stored energy density and its derivatives wrt \( \phi \) and \( t \) (and their spatial derivatives).

Template Parameters
Energythe type of stored energy density \( W \) (currently, only 3D material laws can be used)
Realthe scalar field type of real numbers

Explicit instantiations are provided for dim=2, Real=double, and Energy=StVenantKirchhoff<3,double>. If you want to instantiate the ShellEnergy for different material laws, include fem/diffops/shellKinematics.hpp, where the template definition is provided.

Definition at line 197 of file shellKinematics.hh.

Public Member Functions

 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. 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...
 
Real d0 () const
 Evaluates the marginal stored energy density \( w(x) \). More...
 
template<int row>
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 \( w \) wrt changes. in \( \phi\) and \( t \). More...
 
template<int row, int col>
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 \( w \) wrt changes. in \( \phi\) and \( t \). More...
 

Static Public Attributes

static int const dim = Energy::dim-1
 The dimension of the parameter domain. More...
 

Constructor & Destructor Documentation

◆ ShellEnergy() [1/2]

template<class Energy = MaterialLaws::StVenantKirchhoff<3,double>, class Real = double>
Kaskade::Elastomechanics::ShellEnergy< Energy, Real >::ShellEnergy ( Energy const &  W,
Real  tau 
)

◆ ShellEnergy() [2/2]

template<class Energy = MaterialLaws::StVenantKirchhoff<3,double>, class Real = double>
Kaskade::Elastomechanics::ShellEnergy< Energy, Real >::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 
)
inline

Constructor defining a linearization point.

Parameters
Wthe volume stored energy density
dphi_dxthe spatial derivative of deformation \( \phi \)
ddphi_ddxthe second spatial derivative of deformation \( \phi \)
tthe local shell thickness \( t \)
dt_dxthe spatial derivative of the thickness \( t \)

Definition at line 215 of file shellKinematics.hh.

Member Function Documentation

◆ d0()

template<class Energy = MaterialLaws::StVenantKirchhoff<3,double>, class Real = double>
Real Kaskade::Elastomechanics::ShellEnergy< Energy, Real >::d0 ( ) const

Evaluates the marginal stored energy density \( w(x) \).

◆ d1()

template<class Energy = MaterialLaws::StVenantKirchhoff<3,double>, class Real = double>
template<int row>
Dune::FieldVector< Real, row==0?dim+1:1 > Kaskade::Elastomechanics::ShellEnergy< Energy, Real >::d1 ( VariationalArg< Real, dim, 1 > const &  v) const

Evaulates the parametric directional derivative of the marginal stored energy density \( w \) wrt changes. in \( \phi\) and \( t \).

Template Parameters
roweither 0 or 1. For row=0, the derivative wrt \( \phi \) is computed, for row=1 the derivative wrt \( t \).
Parameters
vthe direction in which to take the derivative
Returns
for row=0, a dim+1 vector \( r_i \) with the derivative values wrt \( \phi_i \), for row=1, a scalar, the derivative wrt \( t \)

◆ d2()

template<class Energy = MaterialLaws::StVenantKirchhoff<3,double>, class Real = double>
template<int row, int col>
Dune::FieldMatrix< Real, row==0?dim+1:1, col==0?dim+1:1 > Kaskade::Elastomechanics::ShellEnergy< Energy, Real >::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 \( w \) wrt changes. in \( \phi\) and \( t \).

◆ setLinearizationPoint()

template<class Energy = MaterialLaws::StVenantKirchhoff<3,double>, class Real = double>
void Kaskade::Elastomechanics::ShellEnergy< Energy, 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.

Referenced by Kaskade::Elastomechanics::ShellEnergy< Energy, Real >::ShellEnergy().

Member Data Documentation

◆ dim

template<class Energy = MaterialLaws::StVenantKirchhoff<3,double>, class Real = double>
int const Kaskade::Elastomechanics::ShellEnergy< Energy, Real >::dim = Energy::dim-1
static

The dimension of the parameter domain.

Definition at line 203 of file shellKinematics.hh.


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