3.5-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
47template<class Scalar, class GridGeometry, class FluxVariables>
49{
50 using FVElementGeometry = typename GridGeometry::LocalView;
51 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
52 using GridView = typename GridGeometry::GridView;
53 using Element = typename GridView::template Codim<0>::Entity;
54
55 static constexpr int dim = GridView::dimension;
56 static constexpr int dimWorld = GridView::dimensionworld;
57 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
58 public:
60 using DimWorldVector = Dune::FieldVector<Scalar, dimWorld>;
61 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
62
122 template<class ElementVolumeVariables>
123 static DimWorldVector velocity(const FVElementGeometry& fvGeometry,
124 const ElementVolumeVariables& elemVolVars,
125 const SubControlVolumeFace& scvf,
126 const int phaseIdx,
127 const DimWorldMatrix sqrtK,
128 const DimWorldVector darcyVelocity)
129 {
130 const auto& problem = elemVolVars.gridVolVars().problem();
131
132 // Obtain the Forchheimer coefficient from the spatial parameters
133 const Scalar forchCoeff = problem.spatialParams().forchCoeff(scvf);
134
135 // Initial guess of velocity: Darcy relation
136 DimWorldVector velocity = darcyVelocity;
137
138 DimWorldVector deltaV(0.0); // the change in velocity between Newton iterations
139 DimWorldVector residual(10e10); // the residual (function value that is to be minimized)
140 DimWorldMatrix gradF(0.0); // slope of equation that is to be solved
141
142 // Search by means of the Newton method for a root of Forchheimer equation
143 static const Scalar epsilon = getParamFromGroup<Scalar>(problem.paramGroup(), "Forchheimer.NewtonTolerance", 1e-12);
144 static const std::size_t maxNumIter = getParamFromGroup<std::size_t>(problem.paramGroup(), "Forchheimer.MaxIterations", 30);
145 for (int k = 0; residual.two_norm() > epsilon ; ++k)
146 {
147 if (k >= maxNumIter)
148 DUNE_THROW(NumericalProblem, "could not determine forchheimer velocity within "
149 << k << " iterations");
150
151 // calculate current value of Forchheimer equation
152 forchheimerResidual_(residual,
153 forchCoeff,
154 sqrtK,
155 darcyVelocity,
156 velocity,
157 elemVolVars,
158 scvf,
159 phaseIdx);
160
161 // newton's method: slope (gradF) and current value (residual) of function is known,
162 forchheimerDerivative_(gradF,
163 forchCoeff,
164 sqrtK,
165 velocity,
166 elemVolVars,
167 scvf,
168 phaseIdx);
169
170 // solve for change in velocity ("x-Axis intercept")
171 gradF.solve(deltaV, residual);
172 velocity -= deltaV;
173 }
174
175 // The fluxvariables expect a value on which an upwinding of the mobility is performed.
176 // We therefore divide by the mobility here.
177 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.mobility(phaseIdx); };
178 const Scalar forchheimerUpwindMobility = UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm, velocity * scvf.unitOuterNormal(), phaseIdx);
179
180 // Do not divide by zero. If the mobility is zero, the resulting flux with upwinding will be zero anyway.
181 if (forchheimerUpwindMobility > 0.0)
182 velocity /= forchheimerUpwindMobility;
183
184 return velocity;
185 }
186
191 template<class Problem, class ElementVolumeVariables,
192 bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value,
193 std::enable_if_t<scalarPerm, int> = 0>
195 const ElementVolumeVariables& elemVolVars,
196 const SubControlVolumeFace& scvf)
197 {
198 using std::sqrt;
199 Scalar harmonicMeanSqrtK(0.0);
200
201 const auto insideScvIdx = scvf.insideScvIdx();
202 const auto& insideVolVars = elemVolVars[insideScvIdx];
203 const Scalar Ki = getPermeability_(problem, insideVolVars, scvf.ipGlobal());
204 const Scalar sqrtKi = sqrt(Ki);
205
206 if (!scvf.boundary())
207 {
208 const auto outsideScvIdx = scvf.outsideScvIdx();
209 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
210 const Scalar Kj = getPermeability_(problem, outsideVolVars, scvf.ipGlobal());
211 const Scalar sqrtKj = sqrt(Kj);
212 harmonicMeanSqrtK = faceTensorAverage(sqrtKi, sqrtKj, scvf.unitOuterNormal());
213 }
214 else
215 harmonicMeanSqrtK = sqrtKi;
216
217 // Convert to tensor
218 DimWorldMatrix result(0.0);
219 for (int i = 0; i < dimWorld; ++i)
220 result[i][i] = harmonicMeanSqrtK;
221
222 return result;
223 }
224
229 template<class Problem, class ElementVolumeVariables,
230 bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value,
231 std::enable_if_t<!scalarPerm, int> = 0>
233 const ElementVolumeVariables& elemVolVars,
234 const SubControlVolumeFace& scvf)
235 {
236 using std::sqrt;
237 DimWorldMatrix sqrtKi(0.0);
238 DimWorldMatrix sqrtKj(0.0);
239 DimWorldMatrix harmonicMeanSqrtK(0.0);
240
241 const auto insideScvIdx = scvf.insideScvIdx();
242 const auto& insideVolVars = elemVolVars[insideScvIdx];
243 const auto& Ki = getPermeability_(problem, insideVolVars, scvf.ipGlobal());
244 // Make sure the permeability matrix does not have off-diagonal entries
245 // TODO implement method using eigenvalues and eigenvectors for general tensors
246 if (!isDiagonal_(Ki))
247 DUNE_THROW(Dune::NotImplemented, "Only diagonal permeability tensors are supported.");
248
249 for (int i = 0; i < dim; ++i)
250 sqrtKi[i][i] = sqrt(Ki[i][i]);
251
252 if (!scvf.boundary())
253 {
254 const auto outsideScvIdx = scvf.outsideScvIdx();
255 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
256 const auto& Kj = getPermeability_(problem, outsideVolVars, scvf.ipGlobal());
257 // Make sure the permeability matrix does not have off-diagonal entries
258 if (!isDiagonal_(Kj))
259 DUNE_THROW(Dune::NotImplemented, "Only diagonal permeability tensors are supported.");
260
261 for (int i = 0; i < dim; ++i)
262 sqrtKj[i][i] = sqrt(Kj[i][i]);
263
264 harmonicMeanSqrtK = faceTensorAverage(sqrtKi, sqrtKj, scvf.unitOuterNormal());
265 }
266 else
267 harmonicMeanSqrtK = sqrtKi;
268
269 return harmonicMeanSqrtK;
270 }
271
272private:
273
275 template<class ElementVolumeVariables>
276 static void forchheimerResidual_(DimWorldVector& residual,
277 const Scalar forchCoeff,
278 const DimWorldMatrix& sqrtK,
279 const DimWorldVector& darcyVelocity,
280 const DimWorldVector& oldVelocity,
281 const ElementVolumeVariables& elemVolVars,
282 const SubControlVolumeFace& scvf,
283 const int phaseIdx)
284 {
285 residual = oldVelocity;
286 // residual += k_r/mu K grad p
287 residual -= darcyVelocity;
288 // residual += c_F rho / mu abs(v) sqrt(K) v
289 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.density(phaseIdx) / volVars.viscosity(phaseIdx); };
290 const Scalar rhoOverMu = UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm, oldVelocity * scvf.unitOuterNormal(), phaseIdx);
291 sqrtK.usmv(forchCoeff * rhoOverMu * oldVelocity.two_norm(), oldVelocity, residual);
292 }
293
295 template<class ElementVolumeVariables>
296 static void forchheimerDerivative_(DimWorldMatrix& derivative,
297 const Scalar forchCoeff,
298 const DimWorldMatrix& sqrtK,
300 const ElementVolumeVariables& elemVolVars,
301 const SubControlVolumeFace& scvf,
302 const int phaseIdx)
303 {
304 // Initialize because for low velocities, we add and do not set the entries.
305 derivative = 0.0;
306
307 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.density(phaseIdx) / volVars.viscosity(phaseIdx); };
308
309 const Scalar absV = velocity.two_norm();
310 // This part of the derivative is only calculated if the velocity is sufficiently small:
311 // do not divide by zero!
312 // This should not be very bad because the derivative is only intended to give an
313 // approximation of the gradient of the
314 // function that goes into the Newton scheme.
315 // This is important in the case of a one-phase region in two-phase flow. The non-existing
316 // phase has a velocity of zero (kr=0).
317 // f = sqrtK * vTimesV (this is a matrix)
318 // f *= forchCoeff density / |velocity| / viscosity (multiply all entries with scalars)
319 const Scalar rhoOverMu = UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm, velocity * scvf.unitOuterNormal(), phaseIdx);
320 if (absV > 1e-20)
321 {
322 for (int i = 0; i < dim; i++)
323 {
324 for (int k = 0; k <= i; k++)
325 {
326 derivative[i][k] = sqrtK[i][i] * velocity[i] * velocity[k] * forchCoeff
327 / absV * rhoOverMu;
328
329 if (k < i)
330 derivative[k][i] = derivative[i][k];
331 }
332 }
333 }
334
335 // add on the main diagonal:
336 // 1 + sqrtK_i forchCoeff density |v| / viscosity
337 for (int i = 0; i < dim; i++)
338 derivative[i][i] += (1.0 + (sqrtK[i][i]*forchCoeff * absV * rhoOverMu));
339 }
340
349 static bool isDiagonal_(const DimWorldMatrix & K)
350 {
351 using std::abs;
352 for (int i = 0; i < dim; i++)
353 for (int k = 0; k < dim; k++)
354 if ((i != k) && (abs(K[i][k]) >= 1e-25))
355 return false;
356 return true;
357 }
358
359 template<class Problem, class VolumeVariables>
360 static decltype(auto) getPermeability_(const Problem& problem,
361 const VolumeVariables& volVars,
362 const GlobalPosition& scvfIpGlobal)
363 {
364 if constexpr (Problem::SpatialParams::evaluatePermeabilityAtScvfIP())
365 return problem.spatialParams().permeabilityAtPos(scvfIpGlobal);
366 else
367 return volVars.permeability();
368 }
369};
370
371} // end namespace Dumux
372
373#endif
Type traits.
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A free function to average a Tensor at an interface.
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
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 This file contains the calculation of the Forchheimer velocity for a given Darcy ve...
Definition: forchheimervelocity.hh:49
typename FluxVariables::UpwindScheme UpwindScheme
Definition: forchheimervelocity.hh:59
static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem &problem, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Returns the harmonic mean of and .
Definition: forchheimervelocity.hh:194
Dune::FieldVector< Scalar, dimWorld > DimWorldVector
Definition: forchheimervelocity.hh:60
static DimWorldVector velocity(const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const DimWorldMatrix sqrtK, const DimWorldVector darcyVelocity)
Compute the Forchheimer velocity of a phase at the given sub-control volume face using the Forchheime...
Definition: forchheimervelocity.hh:123
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimWorldMatrix
Definition: forchheimervelocity.hh:61
Declares all properties used in Dumux.