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#include <dune/common/fvector.hh>
17
18namespace Opm {
19namespace Elasticity {
20
22 template<class GridType>
24 public:
26 static const int dim = GridType::dimension;
27
29 typedef typename GridType::LeafGridView::ctype ctype;
30
33 explicit Elasticity(const GridType& gv_) : gv(gv_) {}
34
38 template<int funcdim>
39 void getBVector(Dune::FieldVector<ctype,funcdim>& BVector,
40 const Dune::FieldVector<ctype,dim>& point);
41
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);
50
56 template<int comp, int funcdim>
57 void getStiffnessMatrix(Dune::FieldMatrix<ctype,funcdim,funcdim>& A,
58 const Dune::FieldMatrix<ctype,comp,funcdim>& B,
59 const Dune::FieldMatrix<ctype,comp,comp>& C,
60 ctype detJW);
61
68 template<int comp, int funcdim>
69 void getStressVector(Dune::FieldVector<ctype,comp>& sigma,
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);
74 protected:
76 const GridType& gv;
77};
78
84Dune::FieldVector<double,3>
85waveSpeeds(const Dune::FieldMatrix<double,6,6>& C,
86 double phi, double theta, double density);
87
88}
89}
90
91#include "elasticity_impl.hpp"
92
93#endif
Elasticity helper class.
Definition: elasticity.hpp:23
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
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.
Dune::FieldVector< ctype, comp > eps0
Definition: elasticity_upscale_impl.hpp:446
Definition: ImplicitAssembly.hpp:43