KASKADE 7 development version
shellKinematics.hh
Go to the documentation of this file.
1/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2/* */
3/* This file is part of the library KASKADE 7 */
4/* https://www.zib.de/research/projects/kaskade7-finite-element-toolbox */
5/* */
6/* Copyright (C) 2019-2019 Zuse Institute Berlin */
7/* */
8/* KASKADE 7 is distributed under the terms of the ZIB Academic License. */
9/* see $KASKADE/academic.txt */
10/* */
11/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
12
13#ifndef SHELLKINEMATICS
14#define SHELLKINEMATICS
15
16#include "dune/common/fmatrix.hh"
17#include "dune/common/fvector.hh"
18
20#include "fem/fixdune.hh"
22
23
25{
44 template <int dim=2, class Real=double>
46 {
47 public:
52
54 using MatrixPair = std::pair<Matrix,Matrix>;
55
60
65 Tensor<Real,dim+1,dim,dim> const& ddphi_ddx_,
66 Real t, Dune::FieldVector<Real,dim> const& dt_dx_);
67
72 Tensor<Real,dim+1,dim,dim> const& ddphi_ddx_,
73 Real t, Dune::FieldVector<Real,dim> const& dt_dx_);
74
78 MatrixPair d0() const;
79
90
104 VariationalArg<Real,dim,1> const& w, int j) const;
105
106 private:
110 Real t;
112 };
113
126 template <int dim=2, class Real=double>
128 {
129 public:
134
135 ShellKinematics(Vector const& phi, Dune::FieldMatrix<Real,dim+1,dim> const& dphi_dx, Real t);
136
141
142 Vector d0() const
143 {
144 }
145
153 template <int row>
154 auto d1(VariationalArg<Real,dim,1> const& v) const;
155
156 private:
159 Real t;
160 };
161
162 // ----------------------------------------------------------------------------------------------
163
195 template <class Energy=MaterialLaws::StVenantKirchhoff<3,double>,
196 class Real=double>
198 {
199 public:
203 static int const dim = Energy::dim-1;
204
205 ShellEnergy(Energy const& W, Real tau);
206
215 ShellEnergy(Energy const& W, Real tau,
217 Tensor<Real,dim+1,dim,dim> const& ddphi_ddx,
218 Real t, Dune::FieldVector<Real,dim> const& dt_dx)
219 : ShellEnergy(W,tau)
220 {
221 setLinearizationPoint(dphi_dx,ddphi_ddx,t,dt_dx);
222 }
223
228 Tensor<Real,dim+1,dim,dim> const& ddphi_ddx,
229 Real t, Dune::FieldVector<Real,dim> const& dt_dx);
230
234 Real d0() const;
235
247 template <int row>
249
254 template <int row, int col>
256 VariationalArg<Real,dim,1> const& w) const;
257
258 private:
259 Energy W;
260 Real tau;
262
263 using Matrix = typename ShellKinematicsDerivative<dim,Real>::Matrix;
264
265 Matrix A, B, E, BtB, BtA, BtA_AtB;
266
267 // compute first directional derivative of energy, given the first directional derivatives of contributing matrices A and B
268 Real d1(std::pair<Matrix,Matrix> const& dAdB) const;
269
270 // compute second directional derivative of energy, given the first and second directional derivatives of contributing matrices A and B
271 // first derivatives need to be given in both directions, whereas the second derivative needs only one value
272 // (order of differentiation irrelevant)
273 Real d2(std::pair<Matrix,Matrix> const& dAdB1, std::pair<Matrix,Matrix> const& dAdB2,
274 std::pair<Matrix,Matrix> const& ddAddB) const;
275 };
276}
277
278#endif
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 .
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.
Definition: tensor.hh:86
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.