25#ifndef DUMUX_DISCRETIZATION_BOX_HOOKES_LAW_HH
26#define DUMUX_DISCRETIZATION_BOX_HOOKES_LAW_HH
28#include <dune/common/fmatrix.hh>
42template<
class ScalarType,
class Gr
idGeometry>
45 using FVElementGeometry =
typename GridGeometry::LocalView;
46 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
47 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
49 using GridView =
typename GridGeometry::GridView;
50 using Element =
typename GridView::template Codim<0>::Entity;
52 static constexpr int dim = GridView::dimension;
53 static constexpr int dimWorld = GridView::dimensionworld;
54 static_assert(dim == dimWorld,
"Hookes Law not implemented for network/surface grids");
67 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache>
69 const Element& element,
70 const FVElementGeometry& fvGeometry,
71 const ElementVolumeVariables& elemVolVars,
72 const SubControlVolumeFace& scvf,
73 const ElementFluxVarsCache& elemFluxVarCache)
75 const auto sigma = stressTensor(problem, element, fvGeometry, elemVolVars, elemFluxVarCache[scvf]);
78 sigma.mv(scvf.unitOuterNormal(), scvfForce);
79 scvfForce *= scvf.area();
85 template<
class Problem,
class ElementVolumeVariables,
class FluxVarsCache>
87 const Element& element,
88 const FVElementGeometry& fvGeometry,
89 const ElementVolumeVariables& elemVolVars,
90 const FluxVarsCache& fluxVarCache)
92 const auto& lameParams = problem.spatialParams().lameParams(element, fvGeometry, elemVolVars, fluxVarCache);
96 for (
int dir = 0; dir < dim; ++dir)
97 for (
const auto& scv : scvs(fvGeometry))
98 gradU[dir].axpy(elemVolVars[scv.indexInElement()].displacement(dir), fluxVarCache.gradN(scv.indexInElement()));
102 for (
int i = 0; i < dim; ++i)
103 for (
int j = 0; j < dimWorld; ++j)
104 epsilon[i][j] = 0.5*(gradU[i][j] + gradU[j][i]);
108 const auto traceEpsilon =
trace(epsilon);
109 for (
int i = 0; i < dim; ++i)
111 sigma[i][i] = lameParams.lambda()*traceEpsilon;
112 for (
int j = 0; j < dimWorld; ++j)
113 sigma[i][j] += 2.0*lameParams.mu()*epsilon[i][j];
Define some often used mathematical functions.
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
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
ScalarType Scalar
export the type used for scalar values
Definition: box/hookeslaw.hh:58
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:86
Dune::FieldMatrix< Scalar, dim, dimWorld > StressTensor
export the type used for the stress tensor
Definition: box/hookeslaw.hh:60
typename StressTensor::row_type ForceVector
export the type used for force vectors
Definition: box/hookeslaw.hh:62
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:68
This computes the stress tensor and surface forces resulting from mechanical deformation.
Definition: hookeslaw.hh:39