3.2-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
32
33namespace Dumux {
34
41template<class ScalarType, class GridGeometry>
42class HookesLaw<ScalarType, GridGeometry, DiscretizationMethod::box>
43{
44 using FVElementGeometry = typename GridGeometry::LocalView;
45 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
46
47 using GridView = typename GridGeometry::GridView;
48 using Element = typename GridView::template Codim<0>::Entity;
49
50 static constexpr int dim = GridView::dimension;
51 static constexpr int dimWorld = GridView::dimensionworld;
52 static_assert(dim == dimWorld, "Hookes Law not implemented for network/surface grids");
53
54public:
56 using Scalar = ScalarType;
58 using StressTensor = Dune::FieldMatrix<Scalar, dim, dimWorld>;
60 using ForceVector = typename StressTensor::row_type;
62 static constexpr DiscretizationMethod discMethod = DiscretizationMethod::box;
63
65 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
66 static ForceVector force(const Problem& problem,
67 const Element& element,
68 const FVElementGeometry& fvGeometry,
69 const ElementVolumeVariables& elemVolVars,
70 const SubControlVolumeFace& scvf,
71 const ElementFluxVarsCache& elemFluxVarCache)
72 {
73 const auto sigma = stressTensor(problem, element, fvGeometry, elemVolVars, elemFluxVarCache[scvf]);
74
75 ForceVector scvfForce(0.0);
76 sigma.mv(scvf.unitOuterNormal(), scvfForce);
77 scvfForce *= scvf.area();
78
79 return scvfForce;
80 }
81
83 template<class Problem, class ElementVolumeVariables, class FluxVarsCache>
84 static StressTensor stressTensor(const Problem& problem,
85 const Element& element,
86 const FVElementGeometry& fvGeometry,
87 const ElementVolumeVariables& elemVolVars,
88 const FluxVarsCache& fluxVarCache)
89 {
90 const auto& lameParams = problem.spatialParams().lameParams(element, fvGeometry, elemVolVars, fluxVarCache);
91
92 // evaluate displacement gradient
93 StressTensor gradU(0.0);
94 for (int dir = 0; dir < dim; ++dir)
95 for (const auto& scv : scvs(fvGeometry))
96 gradU[dir].axpy(elemVolVars[scv.indexInElement()].displacement(dir), fluxVarCache.gradN(scv.indexInElement()));
97
98 // evaluate strain tensor
99 StressTensor epsilon;
100 for (int i = 0; i < dim; ++i)
101 for (int j = 0; j < dimWorld; ++j)
102 epsilon[i][j] = 0.5*(gradU[i][j] + gradU[j][i]);
103
104 // calculate sigma
105 StressTensor sigma(0.0);
106 const auto traceEpsilon = trace(epsilon);
107 for (int i = 0; i < dim; ++i)
108 {
109 sigma[i][i] = lameParams.lambda()*traceEpsilon;
110 for (int j = 0; j < dimWorld; ++j)
111 sigma[i][j] += 2.0*lameParams.mu()*epsilon[i][j];
112 }
113
114 return sigma;
115 }
116};
117
118} // end namespace Dumux
119
120#endif
The available discretization methods in Dumux.
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:760
Definition: adapt.hh:29
ScalarType Scalar
export the type used for scalar values
Definition: box/hookeslaw.hh:56
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:84
Dune::FieldMatrix< Scalar, dim, dimWorld > StressTensor
export the type used for the stress tensor
Definition: box/hookeslaw.hh:58
typename StressTensor::row_type ForceVector
export the type used for force vectors
Definition: box/hookeslaw.hh:60
static ForceVector force(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVarsCache &elemFluxVarCache)
computes the force acting on a sub-control volume face
Definition: box/hookeslaw.hh:66
This computes the stress tensor and surface forces resulting from mechanical deformation.
Definition: hookeslaw.hh:39