Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup GeomechanicsTests | ||
10 | * \brief Definition of a test problem for the linear elastic model. | ||
11 | */ | ||
12 | #ifndef DUMUX_ELASTICPROBLEM_HH | ||
13 | #define DUMUX_ELASTICPROBLEM_HH | ||
14 | |||
15 | #include <dune/common/fmatrix.hh> | ||
16 | |||
17 | #include <dumux/common/math.hh> | ||
18 | #include <dumux/common/properties.hh> | ||
19 | #include <dumux/common/parameters.hh> | ||
20 | #include <dumux/common/boundarytypes.hh> | ||
21 | #include <dumux/common/numeqvector.hh> | ||
22 | |||
23 | #include <dumux/geomechanics/fvproblem.hh> | ||
24 | |||
25 | namespace Dumux { | ||
26 | |||
27 | /*! | ||
28 | * \ingroup GeomechanicsTests | ||
29 | * \brief Problem definition for the deformation of an elastic body. | ||
30 | */ | ||
31 | template<class TypeTag> | ||
32 | 1 | class ElasticProblem : public GeomechanicsFVProblem<TypeTag> | |
33 | { | ||
34 | using ParentType = GeomechanicsFVProblem<TypeTag>; | ||
35 | |||
36 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
37 | using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices; | ||
38 | using BoundaryTypes = Dumux::BoundaryTypes<GetPropType<TypeTag, Properties::ModelTraits>::numEq()>; | ||
39 | using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>; | ||
40 | using NumEqVector = Dumux::NumEqVector<PrimaryVariables>; | ||
41 | |||
42 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
43 | using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView; | ||
44 | using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView; | ||
45 | |||
46 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
47 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
48 | using SubControlVolume = typename GridGeometry::SubControlVolume; | ||
49 | using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace; | ||
50 | |||
51 | using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView; | ||
52 | using Element = typename GridView::template Codim<0>::Entity; | ||
53 | using GlobalPosition = typename Element::Geometry::GlobalCoordinate; | ||
54 | |||
55 | static constexpr Scalar pi = M_PI; | ||
56 | static constexpr int dim = GridView::dimension; | ||
57 | static constexpr int dimWorld = GridView::dimensionworld; | ||
58 | using GradU = Dune::FieldMatrix<Scalar, dim, dimWorld>; | ||
59 | |||
60 | public: | ||
61 | 1 | ElasticProblem(std::shared_ptr<const GridGeometry> gridGeometry) | |
62 |
5/14✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 1 times.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
|
3 | : ParentType(gridGeometry) {} |
63 | |||
64 | //! Evaluates the initial value for a control volume. | ||
65 | ✗ | PrimaryVariables initialAtPos(const GlobalPosition& globalPos) const | |
66 | ✗ | { return PrimaryVariables(0.0); } | |
67 | |||
68 | //! Evaluates the boundary conditions for a Dirichlet boundary segment. | ||
69 | ✗ | PrimaryVariables dirichletAtPos(const GlobalPosition& globalPos) const | |
70 | 304 | { return PrimaryVariables(0.0); } | |
71 | |||
72 | /*! | ||
73 | * \brief Specifies which kind of boundary condition should be | ||
74 | * used for which equation on a given boundary segment. | ||
75 | * | ||
76 | * \param globalPos The global position | ||
77 | */ | ||
78 | 152 | BoundaryTypes boundaryTypesAtPos(const GlobalPosition& globalPos) const | |
79 | { | ||
80 | 152 | BoundaryTypes values; | |
81 | 152 | values.setAllDirichlet(); | |
82 | |||
83 |
12/12✓ Branch 0 taken 40 times.
✓ Branch 1 taken 112 times.
✓ Branch 2 taken 40 times.
✓ Branch 3 taken 112 times.
✓ Branch 4 taken 38 times.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 38 times.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 36 times.
✓ Branch 9 taken 2 times.
✓ Branch 10 taken 36 times.
✓ Branch 11 taken 2 times.
|
304 | if (globalPos[0] < eps_ && globalPos[1] > eps_ && globalPos[1] < 1.0 - eps_) |
84 | 36 | values.setNeumann(Indices::momentum(/*x-dir*/0)); | |
85 |
12/12✓ Branch 0 taken 40 times.
✓ Branch 1 taken 112 times.
✓ Branch 2 taken 40 times.
✓ Branch 3 taken 112 times.
✓ Branch 4 taken 38 times.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 38 times.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 36 times.
✓ Branch 9 taken 2 times.
✓ Branch 10 taken 36 times.
✓ Branch 11 taken 2 times.
|
304 | if (globalPos[1] > 1.0 - eps_ && globalPos[0] > eps_ && globalPos[0] < 1.0 - eps_) |
86 | 36 | values.setNeumann(Indices::momentum(/*y-dir*/1)); | |
87 | |||
88 | 152 | return values; | |
89 | } | ||
90 | |||
91 | //! Evaluates the boundary conditions for a Neumann boundary segment. | ||
92 | 648 | NumEqVector neumann(const Element& element, | |
93 | const FVElementGeometry& fvGeometry, | ||
94 | const ElementVolumeVariables& elemVolvars, | ||
95 | const ElementFluxVariablesCache& elemFluxVarsCache, | ||
96 | const SubControlVolumeFace& scvf) const | ||
97 | { | ||
98 | 1296 | GradU gradU = exactGradient(scvf.ipGlobal()); | |
99 | |||
100 | 648 | GradU epsilon; | |
101 |
2/2✓ Branch 0 taken 1296 times.
✓ Branch 1 taken 648 times.
|
1944 | for (int i = 0; i < dim; ++i) |
102 |
2/2✓ Branch 0 taken 2592 times.
✓ Branch 1 taken 1296 times.
|
3888 | for (int j = 0; j < dimWorld; ++j) |
103 | 18144 | epsilon[i][j] = 0.5*(gradU[i][j] + gradU[j][i]); | |
104 | |||
105 | 1296 | const auto& lameParams = this->spatialParams().lameParamsAtPos(scvf.ipGlobal()); | |
106 | 648 | GradU sigma(0.0); | |
107 | const auto traceEpsilon = trace(epsilon); | ||
108 |
2/2✓ Branch 0 taken 1296 times.
✓ Branch 1 taken 648 times.
|
1944 | for (int i = 0; i < dim; ++i) |
109 | { | ||
110 | 2592 | sigma[i][i] = lameParams.lambda()*traceEpsilon; | |
111 |
2/2✓ Branch 0 taken 2592 times.
✓ Branch 1 taken 1296 times.
|
3888 | for (int j = 0; j < dimWorld; ++j) |
112 | 12960 | sigma[i][j] += 2.0*lameParams.mu()*epsilon[i][j]; | |
113 | } | ||
114 | |||
115 | 648 | NumEqVector values; | |
116 | 1296 | sigma.mv(scvf.unitOuterNormal(), values); | |
117 | 648 | return values; | |
118 | } | ||
119 | |||
120 | /*! | ||
121 | * \brief Evaluates the source term for all phases within a given | ||
122 | * sub-control volume. | ||
123 | */ | ||
124 | 7200 | NumEqVector source(const Element& element, | |
125 | const FVElementGeometry& fvGeometry, | ||
126 | const ElementVolumeVariables& elemVolVars, | ||
127 | const SubControlVolume& scv) const | ||
128 | { | ||
129 | using std::sin; | ||
130 | using std::cos; | ||
131 | |||
132 | 7200 | const auto ipGlobal = scv.center(); | |
133 | 7200 | const auto x = ipGlobal[0]; | |
134 | 7200 | const auto y = ipGlobal[1]; | |
135 | |||
136 | // the lame parameters (we know they only depend on the position here) | ||
137 | 21600 | const auto& lameParams = this->spatialParams().lameParamsAtPos(scv.center()); | |
138 | 7200 | const auto lambda = lameParams.lambda(); | |
139 | 7200 | const auto mu = lameParams.mu(); | |
140 | |||
141 | // precalculated products | ||
142 | 7200 | const Scalar pi_2 = 2.0*pi; | |
143 | 7200 | const Scalar pi_2_square = pi_2*pi_2; | |
144 | 7200 | const Scalar cos_2pix = cos(pi_2*x); | |
145 | 7200 | const Scalar sin_2pix = sin(pi_2*x); | |
146 | 7200 | const Scalar cos_2piy = cos(pi_2*y); | |
147 | 7200 | const Scalar sin_2piy = sin(pi_2*y); | |
148 | |||
149 | 7200 | const Scalar dE11_dx = -2.0*sin_2piy; | |
150 | 7200 | const Scalar dE22_dx = pi_2_square*cos_2pix*cos_2piy; | |
151 | 7200 | const Scalar dE11_dy = pi_2*(1.0-2.0*x)*cos_2piy; | |
152 | 7200 | const Scalar dE22_dy = -1.0*pi_2_square*sin_2pix*sin_2piy; | |
153 | 7200 | const Scalar dE12_dy = 0.5*pi_2_square*(cos_2pix*cos_2piy - (x-x*x)*sin_2piy); | |
154 | 7200 | const Scalar dE21_dx = 0.5*((1.0-2*x)*pi_2*cos_2piy - pi_2_square*sin_2pix*sin_2piy); | |
155 | |||
156 | // compute exact divergence of sigma | ||
157 | 7200 | PrimaryVariables divSigma(0.0); | |
158 | 14400 | divSigma[Indices::momentum(/*x-dir*/0)] = lambda*(dE11_dx + dE22_dx) + 2*mu*(dE11_dx + dE12_dy); | |
159 | 14400 | divSigma[Indices::momentum(/*y-dir*/1)] = lambda*(dE11_dy + dE22_dy) + 2*mu*(dE21_dx + dE22_dy); | |
160 | 7200 | return divSigma; | |
161 | } | ||
162 | |||
163 | /*! | ||
164 | * \brief Evaluates the exact displacement to this problem at a given position. | ||
165 | */ | ||
166 | 121 | PrimaryVariables exactSolution(const GlobalPosition& globalPos) const | |
167 | { | ||
168 | using std::sin; | ||
169 | |||
170 | 121 | const auto x = globalPos[0]; | |
171 | 121 | const auto y = globalPos[1]; | |
172 | |||
173 | 121 | PrimaryVariables exact(0.0); | |
174 | 242 | exact[Indices::momentum(/*x-dir*/0)] = (x-x*x)*sin(2*pi*y); | |
175 | 242 | exact[Indices::momentum(/*y-dir*/1)] = sin(2*pi*x)*sin(2*pi*y); | |
176 | 121 | return exact; | |
177 | } | ||
178 | |||
179 | /*! | ||
180 | * \brief Evaluates the exact displacement gradient to this problem at a given position. | ||
181 | */ | ||
182 | 648 | GradU exactGradient(const GlobalPosition& globalPos) const | |
183 | { | ||
184 | using std::sin; | ||
185 | using std::cos; | ||
186 | |||
187 | 648 | const auto x = globalPos[0]; | |
188 | 648 | const auto y = globalPos[1]; | |
189 | |||
190 | static constexpr int xIdx = Indices::momentum(/*x-dir*/0); | ||
191 | static constexpr int yIdx = Indices::momentum(/*y-dir*/1); | ||
192 | |||
193 | 648 | GradU exactGrad(0.0); | |
194 | 1944 | exactGrad[xIdx][xIdx] = (1-2*x)*sin(2*pi*y); | |
195 | 1944 | exactGrad[xIdx][yIdx] = (x - x*x)*2*pi*cos(2*pi*y); | |
196 | 1944 | exactGrad[yIdx][xIdx] = 2*pi*cos(2*pi*x)*sin(2*pi*y); | |
197 | 1944 | exactGrad[yIdx][yIdx] = 2*pi*sin(2*pi*x)*cos(2*pi*y); | |
198 | 648 | return exactGrad; | |
199 | } | ||
200 | |||
201 | private: | ||
202 | static constexpr Scalar eps_ = 3e-6; | ||
203 | }; | ||
204 | |||
205 | } // end namespace Dumux | ||
206 | |||
207 | #endif | ||
208 |