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");
69 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache>
71 const Element& element,
72 const FVElementGeometry& fvGeometry,
73 const ElementVolumeVariables& elemVolVars,
74 const SubControlVolumeFace& scvf,
75 const ElementFluxVarsCache& elemFluxVarCache)
77 const auto sigma =
stressTensor(problem, element, fvGeometry, elemVolVars, elemFluxVarCache[scvf]);
80 sigma.mv(scvf.unitOuterNormal(), scvfForce);
81 scvfForce *= Extrusion::area(scvf);
87 template<
class Problem,
class ElementVolumeVariables,
class FluxVarsCache>
89 const Element& element,
90 const FVElementGeometry& fvGeometry,
91 const ElementVolumeVariables& elemVolVars,
92 const FluxVarsCache& fluxVarCache)
94 const auto& lameParams = problem.spatialParams().lameParams(element, fvGeometry, elemVolVars, fluxVarCache);
98 for (
int dir = 0; dir < dim; ++dir)
99 for (
const auto& scv : scvs(fvGeometry))
100 gradU[dir].axpy(elemVolVars[scv.indexInElement()].displacement(dir), fluxVarCache.gradN(scv.indexInElement()));
104 for (
int i = 0; i < dim; ++i)
105 for (
int j = 0; j < dimWorld; ++j)
106 epsilon[i][j] = 0.5*(gradU[i][j] + gradU[j][i]);
110 const auto traceEpsilon =
trace(epsilon);
111 for (
int i = 0; i < dim; ++i)
113 sigma[i][i] = lameParams.lambda()*traceEpsilon;
114 for (
int j = 0; j < dimWorld; ++j)
115 sigma[i][j] += 2.0*lameParams.mu()*epsilon[i][j];
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
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
@ box
Definition method.hh:38
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 constexpr DiscretizationMethod discMethod
state the discretization method this implementation belongs to
Definition box/hookeslaw.hh:64
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:88
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:70
This computes the stress tensor and surface forces resulting from mechanical deformation.
Definition hookeslaw.hh:39