version 3.10-dev
box/hookeslaw.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_DISCRETIZATION_BOX_HOOKES_LAW_HH
14#define DUMUX_DISCRETIZATION_BOX_HOOKES_LAW_HH
15
16#include <dune/common/fmatrix.hh>
17
21
22namespace Dumux {
23
30template<class ScalarType, class GridGeometry>
31class HookesLaw<ScalarType, GridGeometry, typename GridGeometry::DiscretizationMethod>
32{
33 using FVElementGeometry = typename GridGeometry::LocalView;
34 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
35 using Extrusion = Extrusion_t<GridGeometry>;
36
37 using GridView = typename GridGeometry::GridView;
38 using Element = typename GridView::template Codim<0>::Entity;
39
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");
43
44public:
46 using Scalar = ScalarType;
48 using StressTensor = Dune::FieldMatrix<Scalar, dim, dimWorld>;
50 using ForceVector = typename StressTensor::row_type;
52
54 // state the discretization method this implementation belongs to
55 static constexpr DiscretizationMethod discMethod{};
56
60 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
61 static ForceVector force(const Problem& problem,
62 const Element& element,
63 const FVElementGeometry& fvGeometry,
64 const ElementVolumeVariables& elemVolVars,
65 const SubControlVolumeFace& scvf,
66 const ElementFluxVarsCache& elemFluxVarCache)
67 {
68 const auto sigma = stressTensor(problem, element, fvGeometry, elemVolVars, elemFluxVarCache[scvf]);
69
70 ForceVector scvfForce(0.0);
71 sigma.mv(scvf.unitOuterNormal(), scvfForce);
72 scvfForce *= Extrusion::area(fvGeometry, scvf);
73
74 return scvfForce;
75 }
76
78 template<class Problem, class ElementVolumeVariables, class FluxVarsCache>
79 static StressTensor stressTensor(const Problem& problem,
80 const Element& element,
81 const FVElementGeometry& fvGeometry,
82 const ElementVolumeVariables& elemVolVars,
83 const FluxVarsCache& fluxVarCache)
84 {
85 const auto& lameParams = problem.spatialParams().lameParams(element, fvGeometry, elemVolVars, fluxVarCache);
86
87 // evaluate displacement gradient
88 StressTensor gradU(0.0);
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()));
92
93 // evaluate strain tensor
94 StressTensor epsilon;
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]);
98
99 // calculate sigma
100 StressTensor sigma(0.0);
101 const auto traceEpsilon = trace(epsilon);
102 for (int i = 0; i < dim; ++i)
103 {
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];
107 }
108
109 return sigma;
110 }
111};
112
113} // end namespace Dumux
114
115#endif
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
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
Definition: method.hh:46