13#ifndef DUMUX_DISCRETIZATION_BOX_HOOKES_LAW_HH
14#define DUMUX_DISCRETIZATION_BOX_HOOKES_LAW_HH
16#include <dune/common/fmatrix.hh>
30template<
class ScalarType,
class Gr
idGeometry>
33 using FVElementGeometry =
typename GridGeometry::LocalView;
34 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
37 using GridView =
typename GridGeometry::GridView;
38 using Element =
typename GridView::template Codim<0>::Entity;
40 static constexpr int dim = GridView::dimension;
41 static constexpr int dimWorld = GridView::dimensionworld;
42 static_assert(dim == dimWorld,
"Hookes Law not implemented for network/surface grids");
60 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache>
62 const Element& element,
63 const FVElementGeometry& fvGeometry,
64 const ElementVolumeVariables& elemVolVars,
65 const SubControlVolumeFace& scvf,
66 const ElementFluxVarsCache& elemFluxVarCache)
68 const auto sigma = stressTensor(problem, element, fvGeometry, elemVolVars, elemFluxVarCache[scvf]);
71 sigma.mv(scvf.unitOuterNormal(), scvfForce);
72 scvfForce *= Extrusion::area(fvGeometry, scvf);
78 template<
class Problem,
class ElementVolumeVariables,
class FluxVarsCache>
80 const Element& element,
81 const FVElementGeometry& fvGeometry,
82 const ElementVolumeVariables& elemVolVars,
83 const FluxVarsCache& fluxVarCache)
85 const auto& lameParams = problem.spatialParams().lameParams(element, fvGeometry, elemVolVars, fluxVarCache);
89 for (
int dir = 0; dir < dim; ++dir)
90 for (
const auto& scv : scvs(fvGeometry))
91 gradU[dir].axpy(elemVolVars[scv.indexInElement()].displacement(dir), fluxVarCache.gradN(scv.indexInElement()));
95 for (
int i = 0; i < dim; ++i)
96 for (
int j = 0; j < dimWorld; ++j)
97 epsilon[i][j] = 0.5*(gradU[i][j] + gradU[j][i]);
101 const auto traceEpsilon =
trace(epsilon);
102 for (
int i = 0; i < dim; ++i)
104 sigma[i][i] = lameParams.lambda()*traceEpsilon;
105 for (
int j = 0; j < dimWorld; ++j)
106 sigma[i][j] += 2.0*lameParams.mu()*epsilon[i][j];
static StressTensor stressTensor(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const FluxVarsCache &fluxVarCache)
assembles the stress tensor at a given integration point
Definition: box/hookeslaw.hh:79
ScalarType Scalar
export the type used for scalar values
Definition: box/hookeslaw.hh:46
Dune::FieldMatrix< Scalar, dim, dimWorld > StressTensor
export the type used for the stress tensor
Definition: box/hookeslaw.hh:48
typename StressTensor::row_type ForceVector
export the type used for force vectors
Definition: box/hookeslaw.hh:50
static ForceVector force(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVarsCache &elemFluxVarCache)
Returns the force (in Newton) acting on a sub-control volume face.
Definition: box/hookeslaw.hh:61
This computes the stress tensor and surface forces resulting from mechanical deformation.
Definition: hookeslaw_fwd.hh:27
Helper classes to compute the integration elements.
Dune::DenseMatrix< MatrixType >::field_type trace(const Dune::DenseMatrix< MatrixType > &M)
Trace of a dense matrix.
Definition: math.hh:800
Stress-Strain relationship according to Hooke's law.
The available discretization methods in Dumux.
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166