Elasticity helper class.
More...
#include <elasticity.hpp>
|
|
typedef GridType::LeafGridView::ctype | ctype |
| | A basic number.
|
|
| | Elasticity (const GridType &gv_) |
| | Default constructor.
|
| template<int funcdim> |
| void | getBVector (Dune::FieldVector< ctype, funcdim > &BVector, const Dune::FieldVector< ctype, dim > &point) |
| | Returns the vector of basis function values in a quadrature point.
|
| template<int components, int funcdim> |
| 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.
|
| template<int comp, int funcdim> |
| 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.
|
| template<int comp, int funcdim> |
| 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.
|
|
template<int funcdims> |
| void | getBVector (Dune::FieldVector< ctype, funcdims > &Bvector, const Dune::FieldVector< ctype, dim > &point) |
|
|
static const int | dim = GridType::dimension |
| | The dimension of our grid.
|
|
|
const GridType & | gv |
| | Const reference to our grid.
|
template<class
GridType>
class Opm::Elasticity::Elasticity< GridType >
Elasticity helper class.
◆ Elasticity()
Default constructor.
- Parameters
-
| [in] | gv_ | The grid we are doing the calculations on |
◆ getBmatrix()
template<int components, int funcdim>
Returns the B matrix in a quadrature point.
- Parameters
-
| [out] | B | The B matrix |
| [in] | point | (Reference) coordinates of quadrature point |
| [in] | Jinv | Jacobian matrix in quadrature point |
◆ getBVector()
Returns the vector of basis function values in a quadrature point.
- Parameters
-
| [out] | BVector | The vector of basis function values |
| [in] | point | (Reference) coordinates of quadrature point |
◆ getStiffnessMatrix()
template<int comp, int funcdim>
Return the stiffness matrix contributions in a quadrature point.
- Parameters
-
| [in] | B | The B matrix in the quadrature point |
| [in] | C | The constitutive matrix for the cell material |
| [in] | detJW | Det J times quadrature weight |
| [out] | A | The stiffness matrix contributions in the quadrature point |
◆ getStressVector()
template<int comp, int funcdim>
| void Opm::Elasticity::Elasticity< GridType >::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.
- Parameters
-
| [in] | v | The displacements in the quadrature point |
| [in] | eps0 | The load case vector |
| [in] | B | The B matrix in the quadrature point |
| [in] | C | The constitutive matrix for the cell material |
| [out] | sigma | The stress vector in the given quadrature point |
The documentation for this class was generated from the following files: