25#ifndef DUMUX_DISCRETIZATION_BOX_HOOKES_LAW_HH
26#define DUMUX_DISCRETIZATION_BOX_HOOKES_LAW_HH
28#include <dune/common/fmatrix.hh>
41template<
class ScalarType,
class Gr
idGeometry>
44 using FVElementGeometry =
typename GridGeometry::LocalView;
45 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
47 using GridView =
typename GridGeometry::GridView;
48 using Element =
typename GridView::template Codim<0>::Entity;
50 static constexpr int dim = GridView::dimension;
51 static constexpr int dimWorld = GridView::dimensionworld;
52 static_assert(dim == dimWorld,
"Hookes Law not implemented for network/surface grids");
65 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache>
67 const Element& element,
68 const FVElementGeometry& fvGeometry,
69 const ElementVolumeVariables& elemVolVars,
70 const SubControlVolumeFace& scvf,
71 const ElementFluxVarsCache& elemFluxVarCache)
73 const auto sigma = stressTensor(problem, element, fvGeometry, elemVolVars, elemFluxVarCache[scvf]);
76 sigma.mv(scvf.unitOuterNormal(), scvfForce);
77 scvfForce *= scvf.area();
83 template<
class Problem,
class ElementVolumeVariables,
class FluxVarsCache>
85 const Element& element,
86 const FVElementGeometry& fvGeometry,
87 const ElementVolumeVariables& elemVolVars,
88 const FluxVarsCache& fluxVarCache)
90 const auto& lameParams = problem.spatialParams().lameParams(element, fvGeometry, elemVolVars, fluxVarCache);
94 for (
int dir = 0; dir < dim; ++dir)
95 for (
const auto& scv : scvs(fvGeometry))
96 gradU[dir].axpy(elemVolVars[scv.indexInElement()].displacement(dir), fluxVarCache.gradN(scv.indexInElement()));
100 for (
int i = 0; i < dim; ++i)
101 for (
int j = 0; j < dimWorld; ++j)
102 epsilon[i][j] = 0.5*(gradU[i][j] + gradU[j][i]);
106 const auto traceEpsilon =
trace(epsilon);
107 for (
int i = 0; i < dim; ++i)
109 sigma[i][i] = lameParams.lambda()*traceEpsilon;
110 for (
int j = 0; j < dimWorld; ++j)
111 sigma[i][j] += 2.0*lameParams.mu()*epsilon[i][j];
The available discretization methods in Dumux.
Hooke's law specialized for different discretization schemes. This computes the stress tensor and sur...
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Dune::DenseMatrix< MatrixType >::field_type trace(const Dune::DenseMatrix< MatrixType > &M)
Trace of a dense matrix.
Definition: math.hh:760
ScalarType Scalar
export the type used for scalar values
Definition: box/hookeslaw.hh:56
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:84
Dune::FieldMatrix< Scalar, dim, dimWorld > StressTensor
export the type used for the stress tensor
Definition: box/hookeslaw.hh:58
typename StressTensor::row_type ForceVector
export the type used for force vectors
Definition: box/hookeslaw.hh:60
static ForceVector force(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVarsCache &elemFluxVarCache)
computes the force acting on a sub-control volume face
Definition: box/hookeslaw.hh:66
This computes the stress tensor and surface forces resulting from mechanical deformation.
Definition: hookeslaw.hh:39