3.4
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, DiscretizationMethod::box>
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 static constexpr DiscretizationMethod discMethod = DiscretizationMethod::box;
65
69 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
70 static ForceVector force(const Problem& problem,
71 const Element& element,
72 const FVElementGeometry& fvGeometry,
73 const ElementVolumeVariables& elemVolVars,
74 const SubControlVolumeFace& scvf,
75 const ElementFluxVarsCache& elemFluxVarCache)
76 {
77 const auto sigma = stressTensor(problem, element, fvGeometry, elemVolVars, elemFluxVarCache[scvf]);
78
79 ForceVector scvfForce(0.0);
80 sigma.mv(scvf.unitOuterNormal(), scvfForce);
81 scvfForce *= Extrusion::area(scvf);
82
83 return scvfForce;
84 }
85
87 template<class Problem, class ElementVolumeVariables, class FluxVarsCache>
88 static StressTensor stressTensor(const Problem& problem,
89 const Element& element,
90 const FVElementGeometry& fvGeometry,
91 const ElementVolumeVariables& elemVolVars,
92 const FluxVarsCache& fluxVarCache)
93 {
94 const auto& lameParams = problem.spatialParams().lameParams(element, fvGeometry, elemVolVars, fluxVarCache);
95
96 // evaluate displacement gradient
97 StressTensor gradU(0.0);
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()));
101
102 // evaluate strain tensor
103 StressTensor epsilon;
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]);
107
108 // calculate sigma
109 StressTensor sigma(0.0);
110 const auto traceEpsilon = trace(epsilon);
111 for (int i = 0; i < dim; ++i)
112 {
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];
116 }
117
118 return sigma;
119 }
120};
121
122} // end namespace Dumux
123
124#endif
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
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
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