elasticity.hpp
Go to the documentation of this file.
1//==============================================================================
11//==============================================================================
12#ifndef ELASTICITY_HPP_
13#define ELASTICITY_HPP_
14
15#include <dune/common/fmatrix.hh>
16
17namespace Opm {
18namespace Elasticity {
19
21 template<class GridType>
23 public:
25 static const int dim = GridType::dimension;
26
28 typedef typename GridType::LeafGridView::ctype ctype;
29
32 explicit Elasticity(const GridType& gv_) : gv(gv_) {}
33
38 template<int funcdim>
39 void getBVector(Dune::FieldVector<ctype,funcdim>& BVector,
40 const Dune::FieldVector<ctype,dim>& point);
41
42 template<int components, int funcdim>
43 void getBmatrix(Dune::FieldMatrix<ctype,components,funcdim>& B,
44 const Dune::FieldVector<ctype,dim>& point,
45 const Dune::FieldMatrix<ctype,dim,dim>& Jinv);
46
52 template<int comp, int funcdim>
53 void getStiffnessMatrix(Dune::FieldMatrix<ctype,funcdim,funcdim>& A,
54 const Dune::FieldMatrix<ctype,comp,funcdim>& B,
55 const Dune::FieldMatrix<ctype,comp,comp>& C,
56 ctype detJW);
57
64 template<int comp, int funcdim>
65 void getStressVector(Dune::FieldVector<ctype,comp>& sigma,
66 const Dune::FieldVector<ctype,funcdim>& v,
67 const Dune::FieldVector<ctype,comp>& eps0,
68 const Dune::FieldMatrix<ctype,comp,funcdim>& B,
69 const Dune::FieldMatrix<ctype,comp,comp>& C);
70 protected:
72 const GridType& gv;
73};
74
80Dune::FieldVector<double,3> waveSpeeds(const Dune::FieldMatrix<double,6,6>& C, double phi,
81 double theta, double density);
82
83}
84}
85
86#include "elasticity_impl.hpp"
87
88#endif
Elasticity helper class.
Definition: elasticity.hpp:22
Elasticity(const GridType &gv_)
Default constructor.
Definition: elasticity.hpp:32
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 B matrix in a quadrature point.
GridType::LeafGridView::ctype ctype
A basic number.
Definition: elasticity.hpp:28
const GridType & gv
Const reference to our grid.
Definition: elasticity.hpp:72
static const int dim
The dimension of our grid.
Definition: elasticity.hpp:25
void getBmatrix(Dune::FieldMatrix< ctype, components, funcdim > &B, const Dune::FieldVector< ctype, dim > &point, const Dune::FieldMatrix< ctype, dim, dim > &Jinv)
Definition: elasticity_impl.hpp:22
Elasticity helper class - template implementations.
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_impl.hpp:116
Dune::FieldVector< ctype, comp > eps0
Definition: elasticity_upscale_impl.hpp:446
Definition: ImplicitAssembly.hpp:43