3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_DISCRETIZATION_BOX_HOOKES_LAW_HH
26#define DUMUX_DISCRETIZATION_BOX_HOOKES_LAW_HH
27
28#include <dune/common/fmatrix.hh>
29
33
34namespace Dumux {
35
42template<class ScalarType, class GridGeometry>
43class HookesLaw<ScalarType, GridGeometry, typename GridGeometry::DiscretizationMethod>
44{
45 using FVElementGeometry = typename GridGeometry::LocalView;
46 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
47 using Extrusion = Extrusion_t<GridGeometry>;
48
49 using GridView = typename GridGeometry::GridView;
50 using Element = typename GridView::template Codim<0>::Entity;
51
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");
55
56public:
58 using Scalar = ScalarType;
60 using StressTensor = Dune::FieldMatrix<Scalar, dim, dimWorld>;
62 using ForceVector = typename StressTensor::row_type;
64
66 // state the discretization method this implementation belongs to
67 static constexpr DiscretizationMethod discMethod{};
68
72 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
73 static ForceVector force(const Problem& problem,
74 const Element& element,
75 const FVElementGeometry& fvGeometry,
76 const ElementVolumeVariables& elemVolVars,
77 const SubControlVolumeFace& scvf,
78 const ElementFluxVarsCache& elemFluxVarCache)
79 {
80 const auto sigma = stressTensor(problem, element, fvGeometry, elemVolVars, elemFluxVarCache[scvf]);
81
82 ForceVector scvfForce(0.0);
83 sigma.mv(scvf.unitOuterNormal(), scvfForce);
84 scvfForce *= Extrusion::area(scvf);
85
86 return scvfForce;
87 }
88
90 template<class Problem, class ElementVolumeVariables, class FluxVarsCache>
91 static StressTensor stressTensor(const Problem& problem,
92 const Element& element,
93 const FVElementGeometry& fvGeometry,
94 const ElementVolumeVariables& elemVolVars,
95 const FluxVarsCache& fluxVarCache)
96 {
97 const auto& lameParams = problem.spatialParams().lameParams(element, fvGeometry, elemVolVars, fluxVarCache);
98
99 // evaluate displacement gradient
100 StressTensor gradU(0.0);
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()));
104
105 // evaluate strain tensor
106 StressTensor epsilon;
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]);
110
111 // calculate sigma
112 StressTensor sigma(0.0);
113 const auto traceEpsilon = trace(epsilon);
114 for (int i = 0; i < dim; ++i)
115 {
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];
119 }
120
121 return sigma;
122 }
123};
124
125} // end namespace Dumux
126
127#endif
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
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
Definition: method.hh:73
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