12#ifndef ELASTICITY_HPP_
13#define ELASTICITY_HPP_
15#include <dune/common/fmatrix.hh>
16#include <dune/common/fvector.hh>
22 template<
class Gr
idType>
26 static const int dim = GridType::dimension;
29 typedef typename GridType::LeafGridView::ctype
ctype;
39 void getBVector(Dune::FieldVector<ctype,funcdim>& BVector,
40 const Dune::FieldVector<ctype,dim>& point);
46 template<
int components,
int funcdim>
47 void getBmatrix(Dune::FieldMatrix<ctype,components,funcdim>& B,
48 const Dune::FieldVector<ctype,dim>& point,
49 const Dune::FieldMatrix<ctype,dim,dim>& Jinv);
56 template<
int comp,
int funcdim>
58 const Dune::FieldMatrix<ctype,comp,funcdim>& B,
59 const Dune::FieldMatrix<ctype,comp,comp>& C,
68 template<
int comp,
int funcdim>
70 const Dune::FieldVector<ctype,funcdim>& v,
71 const Dune::FieldVector<ctype,comp>& eps0,
72 const Dune::FieldMatrix<ctype,comp,funcdim>& B,
73 const Dune::FieldMatrix<ctype,comp,comp>& C);
84Dune::FieldVector<double,3>
85waveSpeeds(
const Dune::FieldMatrix<double,6,6>& C,
86 double phi,
double theta,
double density);
Elasticity(const GridType &gv_)
Default constructor.
Definition elasticity.hpp:33
void getStressVector(Dune::FieldVector< ctype, comp > &sigma, const Dune::FieldVector< ctype, funcdim > &v, const Dune::FieldVector< ctype, comp > &eps0, const Dune::FieldMatrix< ctype, comp, funcdim > &B, const Dune::FieldMatrix< ctype, comp, comp > &C)
Return the stress vector in a quadrature point.
Definition elasticity_impl.hpp:107
void getStiffnessMatrix(Dune::FieldMatrix< ctype, funcdim, funcdim > &A, const Dune::FieldMatrix< ctype, comp, funcdim > &B, const Dune::FieldMatrix< ctype, comp, comp > &C, ctype detJW)
Return the stiffness matrix contributions in a quadrature point.
Definition elasticity_impl.hpp:91
void getBVector(Dune::FieldVector< ctype, funcdim > &BVector, const Dune::FieldVector< ctype, dim > &point)
Returns the vector of basis function values in a quadrature point.
GridType::LeafGridView::ctype ctype
A basic number.
Definition elasticity.hpp:29
const GridType & gv
Const reference to our grid.
Definition elasticity.hpp:76
static const int dim
The dimension of our grid.
Definition elasticity.hpp:26
void getBmatrix(Dune::FieldMatrix< ctype, components, funcdim > &B, const Dune::FieldVector< ctype, dim > &point, const Dune::FieldMatrix< ctype, dim, dim > &Jinv)
Returns the B matrix in a quadrature point.
Definition elasticity_impl.hpp:22
Dune::FieldVector< double, 3 > waveSpeeds(const Dune::FieldMatrix< double, 6, 6 > &C, double phi, double theta, double density)
Compute the elastic wave velocities.
Definition elasticity.cpp:24
Elasticity helper class - template implementations.
Inverting small matrices.
Definition ImplicitAssembly.hpp:43