3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
forchheimervelocity.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 *****************************************************************************/
26#ifndef DUMUX_FLUX_FORCHHEIMER_VELOCITY_HH
27#define DUMUX_FLUX_FORCHHEIMER_VELOCITY_HH
28
29#include <dune/common/fvector.hh>
30#include <dune/common/fmatrix.hh>
31
32#include <dumux/common/math.hh>
36
38
39namespace Dumux {
40
49template<class Scalar, class GridGeometry, class FluxVariables>
51{
52 using FVElementGeometry = typename GridGeometry::LocalView;
53 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
54 using GridView = typename GridGeometry::GridView;
55 using Element = typename GridView::template Codim<0>::Entity;
56
57 static constexpr int dim = GridView::dimension;
58 static constexpr int dimWorld = GridView::dimensionworld;
59 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
60 public:
62 using DimWorldVector = Dune::FieldVector<Scalar, dimWorld>;
63 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
64
65
66 template<class ElementVolumeVariables>
67 static DimWorldVector velocity(const FVElementGeometry& fvGeometry,
68 const ElementVolumeVariables& elemVolVars,
69 const SubControlVolumeFace& scvf,
70 const int phaseIdx,
71 const DimWorldMatrix sqrtK,
72 const DimWorldVector darcyVelocity)
73 {
74 const auto& problem = elemVolVars.gridVolVars().problem();
75
76 // Obtain the Forchheimer coefficient from the spatial parameters
77 const Scalar forchCoeff = problem.spatialParams().forchCoeff(scvf);
78
79 // Initial guess of velocity: Darcy relation
80 DimWorldVector velocity = darcyVelocity;
81
82 DimWorldVector deltaV(0.0); // the change in velocity between Newton iterations
83 DimWorldVector residual(10e10); // the residual (function value that is to be minimized)
84 DimWorldMatrix gradF(0.0); // slope of equation that is to be solved
85
86 // Search by means of the Newton method for a root of Forchheimer equation
87 static const Scalar epsilon = getParamFromGroup<Scalar>(problem.paramGroup(), "Forchheimer.NewtonTolerance", 1e-12);
88 static const std::size_t maxNumIter = getParamFromGroup<std::size_t>(problem.paramGroup(), "Forchheimer.MaxIterations", 30);
89 for (int k = 0; residual.two_norm() > epsilon ; ++k)
90 {
91 if (k >= maxNumIter)
92 DUNE_THROW(NumericalProblem, "could not determine forchheimer velocity within "
93 << k << " iterations");
94
95 // calculate current value of Forchheimer equation
96 forchheimerResidual_(residual,
97 forchCoeff,
98 sqrtK,
99 darcyVelocity,
100 velocity,
101 elemVolVars,
102 scvf,
103 phaseIdx);
104
105 // newton's method: slope (gradF) and current value (residual) of function is known,
106 forchheimerDerivative_(gradF,
107 forchCoeff,
108 sqrtK,
109 velocity,
110 elemVolVars,
111 scvf,
112 phaseIdx);
113
114 // solve for change in velocity ("x-Axis intercept")
115 gradF.solve(deltaV, residual);
116 velocity -= deltaV;
117 }
118
119 // The fluxvariables expect a value on which an upwinding of the mobility is performed.
120 // We therefore divide by the mobility here.
121 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.mobility(phaseIdx); };
122 const Scalar forchheimerUpwindMobility = UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm, velocity * scvf.unitOuterNormal(), phaseIdx);
123
124 // Do not divide by zero. If the mobility is zero, the resulting flux with upwinding will be zero anyway.
125 if (forchheimerUpwindMobility > 0.0)
126 velocity /= forchheimerUpwindMobility;
127
128 return velocity;
129 }
130
135 template<class Problem, class ElementVolumeVariables,
136 bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value,
137 std::enable_if_t<scalarPerm, int> = 0>
139 const ElementVolumeVariables& elemVolVars,
140 const SubControlVolumeFace& scvf)
141 {
142 using std::sqrt;
143 Scalar harmonicMeanSqrtK(0.0);
144
145 const auto insideScvIdx = scvf.insideScvIdx();
146 const auto& insideVolVars = elemVolVars[insideScvIdx];
147 const Scalar Ki = getPermeability_(problem, insideVolVars, scvf.ipGlobal());
148 const Scalar sqrtKi = sqrt(Ki);
149
150 if (!scvf.boundary())
151 {
152 const auto outsideScvIdx = scvf.outsideScvIdx();
153 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
154 const Scalar Kj = getPermeability_(problem, outsideVolVars, scvf.ipGlobal());
155 const Scalar sqrtKj = sqrt(Kj);
156 harmonicMeanSqrtK = faceTensorAverage(sqrtKi, sqrtKj, scvf.unitOuterNormal());
157 }
158 else
159 harmonicMeanSqrtK = sqrtKi;
160
161 // Convert to tensor
162 DimWorldMatrix result(0.0);
163 for (int i = 0; i < dimWorld; ++i)
164 result[i][i] = harmonicMeanSqrtK;
165
166 return result;
167 }
168
173 template<class Problem, class ElementVolumeVariables,
174 bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value,
175 std::enable_if_t<!scalarPerm, int> = 0>
177 const ElementVolumeVariables& elemVolVars,
178 const SubControlVolumeFace& scvf)
179 {
180 using std::sqrt;
181 DimWorldMatrix sqrtKi(0.0);
182 DimWorldMatrix sqrtKj(0.0);
183 DimWorldMatrix harmonicMeanSqrtK(0.0);
184
185 const auto insideScvIdx = scvf.insideScvIdx();
186 const auto& insideVolVars = elemVolVars[insideScvIdx];
187 const auto& Ki = getPermeability_(problem, insideVolVars, scvf.ipGlobal());
188 // Make sure the permeability matrix does not have off-diagonal entries
189 // TODO implement method using eigenvalues and eigenvectors for general tensors
190 if (!isDiagonal_(Ki))
191 DUNE_THROW(Dune::NotImplemented, "Only diagonal permeability tensors are supported.");
192
193 for (int i = 0; i < dim; ++i)
194 sqrtKi[i][i] = sqrt(Ki[i][i]);
195
196 if (!scvf.boundary())
197 {
198 const auto outsideScvIdx = scvf.outsideScvIdx();
199 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
200 const auto& Kj = getPermeability_(problem, outsideVolVars, scvf.ipGlobal());
201 // Make sure the permeability matrix does not have off-diagonal entries
202 if (!isDiagonal_(Kj))
203 DUNE_THROW(Dune::NotImplemented, "Only diagonal permeability tensors are supported.");
204
205 for (int i = 0; i < dim; ++i)
206 sqrtKj[i][i] = sqrt(Kj[i][i]);
207
208 harmonicMeanSqrtK = faceTensorAverage(sqrtKi, sqrtKj, scvf.unitOuterNormal());
209 }
210 else
211 harmonicMeanSqrtK = sqrtKi;
212
213 return harmonicMeanSqrtK;
214 }
215
216private:
217
219 template<class ElementVolumeVariables>
220 static void forchheimerResidual_(DimWorldVector& residual,
221 const Scalar forchCoeff,
222 const DimWorldMatrix& sqrtK,
223 const DimWorldVector& darcyVelocity,
224 const DimWorldVector& oldVelocity,
225 const ElementVolumeVariables& elemVolVars,
226 const SubControlVolumeFace& scvf,
227 const int phaseIdx)
228 {
229 residual = oldVelocity;
230 // residual += k_r/mu K grad p
231 residual -= darcyVelocity;
232 // residual += c_F rho / mu abs(v) sqrt(K) v
233 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.density(phaseIdx) / volVars.viscosity(phaseIdx); };
234 const Scalar rhoOverMu = UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm, oldVelocity * scvf.unitOuterNormal(), phaseIdx);
235 sqrtK.usmv(forchCoeff * rhoOverMu * oldVelocity.two_norm(), oldVelocity, residual);
236 }
237
239 template<class ElementVolumeVariables>
240 static void forchheimerDerivative_(DimWorldMatrix& derivative,
241 const Scalar forchCoeff,
242 const DimWorldMatrix& sqrtK,
244 const ElementVolumeVariables& elemVolVars,
245 const SubControlVolumeFace& scvf,
246 const int phaseIdx)
247 {
248 // Initialize because for low velocities, we add and do not set the entries.
249 derivative = 0.0;
250
251 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.density(phaseIdx) / volVars.viscosity(phaseIdx); };
252
253 const Scalar absV = velocity.two_norm();
254 // This part of the derivative is only calculated if the velocity is sufficiently small:
255 // do not divide by zero!
256 // This should not be very bad because the derivative is only intended to give an
257 // approximation of the gradient of the
258 // function that goes into the Newton scheme.
259 // This is important in the case of a one-phase region in two-phase flow. The non-existing
260 // phase has a velocity of zero (kr=0).
261 // f = sqrtK * vTimesV (this is a matrix)
262 // f *= forchCoeff density / |velocity| / viscosity (multiply all entries with scalars)
263 const Scalar rhoOverMu = UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm, velocity * scvf.unitOuterNormal(), phaseIdx);
264 if (absV > 1e-20)
265 {
266 for (int i = 0; i < dim; i++)
267 {
268 for (int k = 0; k <= i; k++)
269 {
270 derivative[i][k] = sqrtK[i][i] * velocity[i] * velocity[k] * forchCoeff
271 / absV * rhoOverMu;
272
273 if (k < i)
274 derivative[k][i] = derivative[i][k];
275 }
276 }
277 }
278
279 // add on the main diagonal:
280 // 1 + sqrtK_i forchCoeff density |v| / viscosity
281 for (int i = 0; i < dim; i++)
282 derivative[i][i] += (1.0 + (sqrtK[i][i]*forchCoeff * absV * rhoOverMu));
283 }
284
293 static bool isDiagonal_(const DimWorldMatrix & K)
294 {
295 using std::abs;
296 for (int i = 0; i < dim; i++)
297 for (int k = 0; k < dim; k++)
298 if ((i != k) && (abs(K[i][k]) >= 1e-25))
299 return false;
300 return true;
301 }
302
303 template<class Problem, class VolumeVariables>
304 static decltype(auto) getPermeability_(const Problem& problem,
305 const VolumeVariables& volVars,
306 const GlobalPosition& scvfIpGlobal)
307 {
308 if constexpr (Problem::SpatialParams::evaluatePermeabilityAtScvfIP())
309 return problem.spatialParams().permeabilityAtPos(scvfIpGlobal);
310 else
311 return volVars.permeability();
312 }
313};
314
315} // end namespace Dumux
316
317#endif
A free function to average a Tensor at an interface.
Type traits.
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
UpwindSchemeImpl< GridGeometry, typename GridGeometry::DiscretizationMethod > UpwindScheme
The upwind scheme used for the advective fluxes. This depends on the chosen discretization method.
Definition: flux/upwindscheme.hh:42
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Scalar faceTensorAverage(const Scalar T1, const Scalar T2, const Dune::FieldVector< Scalar, dim > &normal)
Average of a discontinuous scalar field at discontinuity interface (for compatibility reasons with th...
Definition: facetensoraverage.hh:41
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
Forchheimer's law For a detailed description see dumux/flow/forchheimerslaw.hh.
Definition: forchheimervelocity.hh:51
typename FluxVariables::UpwindScheme UpwindScheme
Definition: forchheimervelocity.hh:61
static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem &problem, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Returns the harmonic mean of and .
Definition: forchheimervelocity.hh:138
Dune::FieldVector< Scalar, dimWorld > DimWorldVector
Definition: forchheimervelocity.hh:62
static DimWorldVector velocity(const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const DimWorldMatrix sqrtK, const DimWorldVector darcyVelocity)
Definition: forchheimervelocity.hh:67
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimWorldMatrix
Definition: forchheimervelocity.hh:63
Declares all properties used in Dumux.