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 SubControlVolumeFace =
typename GridGeometry::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");
72 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache>
74 const Element& element,
75 const FVElementGeometry& fvGeometry,
76 const ElementVolumeVariables& elemVolVars,
77 const SubControlVolumeFace& scvf,
78 const ElementFluxVarsCache& elemFluxVarCache)
80 const auto sigma = stressTensor(problem, element, fvGeometry, elemVolVars, elemFluxVarCache[scvf]);
83 sigma.mv(scvf.unitOuterNormal(), scvfForce);
84 scvfForce *= Extrusion::area(scvf);
90 template<
class Problem,
class ElementVolumeVariables,
class FluxVarsCache>
92 const Element& element,
93 const FVElementGeometry& fvGeometry,
94 const ElementVolumeVariables& elemVolVars,
95 const FluxVarsCache& fluxVarCache)
97 const auto& lameParams = problem.spatialParams().lameParams(element, fvGeometry, elemVolVars, fluxVarCache);
101 for (
int dir = 0; dir < dim; ++dir)
102 for (
const auto& scv : scvs(fvGeometry))
103 gradU[dir].axpy(elemVolVars[scv.indexInElement()].displacement(dir), fluxVarCache.gradN(scv.indexInElement()));
107 for (
int i = 0; i < dim; ++i)
108 for (
int j = 0; j < dimWorld; ++j)
109 epsilon[i][j] = 0.5*(gradU[i][j] + gradU[j][i]);
113 const auto traceEpsilon =
trace(epsilon);
114 for (
int i = 0; i < dim; ++i)
116 sigma[i][i] = lameParams.lambda()*traceEpsilon;
117 for (
int j = 0; j < dimWorld; ++j)
118 sigma[i][j] += 2.0*lameParams.mu()*epsilon[i][j];
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Hooke's law specialized for different discretization schemes. This computes the stress tensor and sur...
Dune::DenseMatrix< MatrixType >::field_type trace(const Dune::DenseMatrix< MatrixType > &M)
Trace of a dense matrix.
Definition: math.hh:783
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
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:91
ScalarType Scalar
export the type used for scalar values
Definition: box/hookeslaw.hh:58
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)
Returns the force (in Newton) acting on a sub-control volume face.
Definition: box/hookeslaw.hh:73
This computes the stress tensor and surface forces resulting from mechanical deformation.
Definition: hookeslaw_fwd.hh:39