3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
cctpfa/forchheimerslaw.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 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FORCHHEIMERS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FORCHHEIMERS_LAW_HH
26
27#include <dune/common/fvector.hh>
28#include <dune/common/fmatrix.hh>
29
30#include <dumux/common/math.hh>
34
38
39namespace Dumux {
40
41// forward declarations
42template<class TypeTag, DiscretizationMethod discMethod>
43class ForchheimersLawImplementation;
44
53template<class Scalar, class GridGeometry, bool isNetwork>
55
61template <class TypeTag>
63: public CCTpfaForchheimersLaw<GetPropType<TypeTag, Properties::Scalar>,
64 GetPropType<TypeTag, Properties::GridGeometry>,
65 (GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension < GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld)>
66{};
67
72template<class GridGeometry>
73class TpfaForchheimersLawCacheFiller
74{
75 using FVElementGeometry = typename GridGeometry::LocalView;
76 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
77 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
78
79public:
83 template<class FluxVariablesCache, class Problem, class ElementVolumeVariables, class FluxVariablesCacheFiller>
84 static void fill(FluxVariablesCache& scvfFluxVarsCache,
85 const Problem& problem,
86 const Element& element,
87 const FVElementGeometry& fvGeometry,
88 const ElementVolumeVariables& elemVolVars,
89 const SubControlVolumeFace& scvf,
90 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
91 {
92 scvfFluxVarsCache.updateAdvection(problem, element, fvGeometry, elemVolVars, scvf);
93 }
94};
95
100template<class AdvectionType, class GridGeometry>
101class TpfaForchheimersLawCache
102{
103 using Scalar = typename AdvectionType::Scalar;
104 using FVElementGeometry = typename GridGeometry::LocalView;
105 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
106 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
107 static constexpr int dimWorld = GridGeometry::GridView::dimensionworld;
108 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
109
110public:
111 using Filler = TpfaForchheimersLawCacheFiller<GridGeometry>;
112
113 template<class Problem, class ElementVolumeVariables>
114 void updateAdvection(const Problem& problem,
115 const Element& element,
116 const FVElementGeometry& fvGeometry,
117 const ElementVolumeVariables& elemVolVars,
118 const SubControlVolumeFace &scvf)
119 {
120 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
121 harmonicMeanSqrtK_ = AdvectionType::calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf);
122 }
123
124 const Scalar& advectionTij() const
125 { return tij_; }
126
127 const DimWorldMatrix& harmonicMeanSqrtPermeability() const
128 { return harmonicMeanSqrtK_; }
129
130private:
131 Scalar tij_;
132 DimWorldMatrix harmonicMeanSqrtK_;
133};
134
139template<class ScalarType, class GridGeometry>
140class CCTpfaForchheimersLaw<ScalarType, GridGeometry, /*isNetwork*/ false>
141{
142 using ThisType = CCTpfaForchheimersLaw<ScalarType, GridGeometry, /*isNetwork*/ false>;
143 using FVElementGeometry = typename GridGeometry::LocalView;
144 using SubControlVolume = typename GridGeometry::SubControlVolume;
145 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
146 using Extrusion = Extrusion_t<GridGeometry>;
147 using GridView = typename GridGeometry::GridView;
148 using Element = typename GridView::template Codim<0>::Entity;
149
150 static constexpr int dim = GridView::dimension;
151 static constexpr int dimWorld = GridView::dimensionworld;
152
153 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
154 using DimWorldVector = Dune::FieldVector<ScalarType, dimWorld>;
155 using DimWorldMatrix = Dune::FieldMatrix<ScalarType, dimWorld, dimWorld>;
156
157 using DarcysLaw = CCTpfaDarcysLaw<ScalarType, GridGeometry, /*isNetwork*/ false>;
158
159 public:
161 using Scalar = ScalarType;
162
164 static const DiscretizationMethod discMethod = DiscretizationMethod::cctpfa;
165
167 using Cache = TpfaForchheimersLawCache<ThisType, GridGeometry>;
168
231 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
232 static Scalar flux(const Problem& problem,
233 const Element& element,
234 const FVElementGeometry& fvGeometry,
235 const ElementVolumeVariables& elemVolVars,
236 const SubControlVolumeFace& scvf,
237 int phaseIdx,
238 const ElementFluxVarsCache& elemFluxVarsCache)
239 {
240 const auto velocity = velocity_(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
241 Scalar flux = velocity * scvf.unitOuterNormal();
242 flux *= Extrusion::area(scvf);
243 return flux;
244 }
245
248 template<class Problem, class ElementVolumeVariables>
249 static Scalar calculateTransmissibility(const Problem& problem,
250 const Element& element,
251 const FVElementGeometry& fvGeometry,
252 const ElementVolumeVariables& elemVolVars,
253 const SubControlVolumeFace& scvf)
254 {
255 return DarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
256 }
257
262 template<class Problem, class ElementVolumeVariables,
263 bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value,
264 std::enable_if_t<scalarPerm, int> = 0>
265 static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem& problem,
266 const ElementVolumeVariables& elemVolVars,
267 const SubControlVolumeFace& scvf)
268 {
269 using std::sqrt;
270 Scalar harmonicMeanSqrtK(0.0);
271
272 const auto insideScvIdx = scvf.insideScvIdx();
273 const auto& insideVolVars = elemVolVars[insideScvIdx];
274 const Scalar Ki = getPermeability_(problem, insideVolVars, scvf.ipGlobal());
275 const Scalar sqrtKi = sqrt(Ki);
276
277 if (!scvf.boundary())
278 {
279 const auto outsideScvIdx = scvf.outsideScvIdx();
280 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
281 const Scalar Kj = getPermeability_(problem, outsideVolVars, scvf.ipGlobal());
282 const Scalar sqrtKj = sqrt(Kj);
283 harmonicMeanSqrtK = problem.spatialParams().harmonicMean(sqrtKi, sqrtKj, scvf.unitOuterNormal());
284 }
285 else
286 harmonicMeanSqrtK = sqrtKi;
287
288 // Convert to tensor
289 DimWorldMatrix result(0.0);
290 for (int i = 0; i < dimWorld; ++i)
291 result[i][i] = harmonicMeanSqrtK;
292
293 return result;
294 }
295
300 template<class Problem, class ElementVolumeVariables,
301 bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value,
302 std::enable_if_t<!scalarPerm, int> = 0>
303 static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem& problem,
304 const ElementVolumeVariables& elemVolVars,
305 const SubControlVolumeFace& scvf)
306 {
307 using std::sqrt;
308 DimWorldMatrix sqrtKi(0.0);
309 DimWorldMatrix sqrtKj(0.0);
310 DimWorldMatrix harmonicMeanSqrtK(0.0);
311
312 const auto insideScvIdx = scvf.insideScvIdx();
313 const auto& insideVolVars = elemVolVars[insideScvIdx];
314 const auto& Ki = getPermeability_(problem, insideVolVars, scvf.ipGlobal());
315 // Make sure the permeability matrix does not have off-diagonal entries
316 // TODO implement method using eigenvalues and eigenvectors for general tensors
317 if (!isDiagonal_(Ki))
318 DUNE_THROW(Dune::NotImplemented, "Only diagonal permeability tensors are supported.");
319
320 for (int i = 0; i < dim; ++i)
321 sqrtKi[i][i] = sqrt(Ki[i][i]);
322
323 if (!scvf.boundary())
324 {
325 const auto outsideScvIdx = scvf.outsideScvIdx();
326 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
327 const auto& Kj = getPermeability_(problem, outsideVolVars, scvf.ipGlobal());
328 // Make sure the permeability matrix does not have off-diagonal entries
329 if (!isDiagonal_(Kj))
330 DUNE_THROW(Dune::NotImplemented, "Only diagonal permeability tensors are supported.");
331
332 for (int i = 0; i < dim; ++i)
333 sqrtKj[i][i] = sqrt(Kj[i][i]);
334
335 harmonicMeanSqrtK = problem.spatialParams().harmonicMean(sqrtKi, sqrtKj, scvf.unitOuterNormal());
336 }
337 else
338 harmonicMeanSqrtK = sqrtKi;
339
340 return harmonicMeanSqrtK;
341 }
342
343private:
344
346 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
347 static DimWorldVector velocity_(const Problem& problem,
348 const Element& element,
349 const FVElementGeometry& fvGeometry,
350 const ElementVolumeVariables& elemVolVars,
351 const SubControlVolumeFace& scvf,
352 int phaseIdx,
353 const ElementFluxVarsCache& elemFluxVarsCache)
354 {
355 // Get the volume flux based on Darcy's law. The value returned by this method needs to be multiplied with the
356 // mobility (upwinding).
357 Scalar darcyFlux = DarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
358 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.mobility(phaseIdx); };
359 const Scalar darcyUpwindMobility = doUpwind_(scvf, elemVolVars, upwindTerm, !std::signbit(darcyFlux) /*insideIsUpstream*/);
360 darcyFlux *= darcyUpwindMobility;
361
362 // Convert to Darcy velocity
363 DimWorldVector darcyVelocity = scvf.unitOuterNormal();
364 darcyVelocity *= darcyFlux / Extrusion::area(scvf);
365
366 // Get the harmonic mean of the square root of permeability
367 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
368 const auto& sqrtK = fluxVarsCache.harmonicMeanSqrtPermeability();
369
370 // Obtain the Forchheimer coefficient from the spatial parameters
371 const Scalar forchCoeff = problem.spatialParams().forchCoeff(scvf);
372
373 // Initial guess of velocity: Darcy relation
374 DimWorldVector velocity = darcyVelocity;
375
376 DimWorldVector deltaV(0.0); // the change in velocity between Newton iterations
377 DimWorldVector residual(10e10); // the residual (function value that is to be minimized)
378 DimWorldMatrix gradF(0.0); // slope of equation that is to be solved
379
380 // Search by means of the Newton method for a root of Forchheimer equation
381 static const Scalar epsilon = getParamFromGroup<Scalar>(problem.paramGroup(), "Forchheimer.NewtonTolerance", 1e-12);
382 static const std::size_t maxNumIter = getParamFromGroup<std::size_t>(problem.paramGroup(), "Forchheimer.MaxIterations", 30);
383 for (int k = 0; residual.two_norm() > epsilon ; ++k)
384 {
385 if (k >= maxNumIter)
386 DUNE_THROW(NumericalProblem, "could not determine forchheimer velocity within "
387 << k << " iterations");
388
389 // calculate current value of Forchheimer equation
390 forchheimerResidual_(residual,
391 forchCoeff,
392 sqrtK,
393 darcyVelocity,
394 velocity,
395 elemVolVars,
396 scvf,
397 phaseIdx);
398
399 // newton's method: slope (gradF) and current value (residual) of function is known,
400 forchheimerDerivative_(gradF,
401 forchCoeff,
402 sqrtK,
403 velocity,
404 elemVolVars,
405 scvf,
406 phaseIdx);
407
408 // solve for change in velocity ("x-Axis intercept")
409 gradF.solve(deltaV, residual);
410 velocity -= deltaV;
411 }
412
413 // The fluxvariables expect a value on which an upwinding of the mobility is performed.
414 // We therefore divide by the mobility here.
415 const Scalar forchheimerUpwindMobility = doUpwind_(scvf, elemVolVars, upwindTerm,
416 !std::signbit(velocity * scvf.unitOuterNormal()) /*insideIsUpstream*/);
417
418 // Do not divide by zero. If the mobility is zero, the resulting flux with upwinding will be zero anyway.
419 if (forchheimerUpwindMobility > 0.0)
420 velocity /= forchheimerUpwindMobility;
421
422 return velocity;
423 }
424
426 template<class ElementVolumeVariables>
427 static void forchheimerResidual_(DimWorldVector& residual,
428 const Scalar forchCoeff,
429 const DimWorldMatrix& sqrtK,
430 const DimWorldVector& darcyVelocity,
431 const DimWorldVector& oldVelocity,
432 const ElementVolumeVariables& elemVolVars,
433 const SubControlVolumeFace& scvf,
434 const int phaseIdx)
435 {
436 residual = oldVelocity;
437
438 // residual += k_r/mu K grad p
439 residual -= darcyVelocity;
440
441 // residual += c_F rho / mu abs(v) sqrt(K) v
442 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.density(phaseIdx) / volVars.viscosity(phaseIdx); };
443 bool insideIsUpstream = !std::signbit(oldVelocity * scvf.unitOuterNormal());
444 const Scalar rhoOverMu = doUpwind_(scvf, elemVolVars, upwindTerm, insideIsUpstream);
445 sqrtK.usmv(forchCoeff * rhoOverMu * oldVelocity.two_norm(), oldVelocity, residual);
446 }
447
449 template<class ElementVolumeVariables>
450 static void forchheimerDerivative_(DimWorldMatrix& derivative,
451 const Scalar forchCoeff,
452 const DimWorldMatrix& sqrtK,
453 const DimWorldVector& velocity,
454 const ElementVolumeVariables& elemVolVars,
455 const SubControlVolumeFace& scvf,
456 const int phaseIdx)
457 {
458
459
460 // Initialize because for low velocities, we add and do not set the entries.
461 derivative = 0.0;
462
463 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.density(phaseIdx) / volVars.viscosity(phaseIdx); };
464 bool insideIsUpstream = !std::signbit(velocity * scvf.unitOuterNormal());
465
466 const Scalar absV = velocity.two_norm() ;
467 // This part of the derivative is only calculated if the velocity is sufficiently small:
468 // do not divide by zero!
469 // This should not be very bad because the derivative is only intended to give an
470 // approximation of the gradient of the
471 // function that goes into the Newton scheme.
472 // This is important in the case of a one-phase region in two-phase flow. The non-existing
473 // phase has a velocity of zero (kr=0).
474 // f = sqrtK * vTimesV (this is a matrix)
475 // f *= forchCoeff density / |velocity| / viscosity (multiply all entries with scalars)
476 const Scalar rhoOverMu = doUpwind_(scvf, elemVolVars, upwindTerm, insideIsUpstream);
477 if (absV > 1e-20)
478 {
479 for (int i = 0; i < dim; i++)
480 {
481 for (int k = 0; k <= i; k++)
482 {
483 derivative[i][k] = sqrtK[i][i] * velocity[i] * velocity[k] * forchCoeff
484 / absV * rhoOverMu;
485
486 if (k < i)
487 derivative[k][i] = derivative[i][k];
488 }
489 }
490 }
491
492 // add on the main diagonal:
493 // 1 + sqrtK_i forchCoeff density |v| / viscosity
494 for (int i = 0; i < dim; i++)
495 derivative[i][i] += (1.0 + (sqrtK[i][i]*forchCoeff * absV * rhoOverMu));
496 }
497
498 template<class ElementVolumeVariables, class UpwindTermFunction>
499 static Scalar doUpwind_(const SubControlVolumeFace& scvf,
500 const ElementVolumeVariables& elemVolVars,
501 const UpwindTermFunction& upwindTerm,
502 const bool insideIsUpstream)
503 {
504 static const Scalar upwindWeight = getParam<Scalar>("Flux.UpwindWeight");
505
506 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
507 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
508
509 if (!insideIsUpstream) // if sign of flux is negative
510 return (upwindWeight*upwindTerm(outsideVolVars)
511 + (1.0 - upwindWeight)*upwindTerm(insideVolVars));
512 else
513 return (upwindWeight*upwindTerm(insideVolVars)
514 + (1.0 - upwindWeight)*upwindTerm(outsideVolVars));
515 }
516
525 static bool isDiagonal_(const DimWorldMatrix & K)
526 {
527 using std::abs;
528 for (int i = 0; i < dim; i++) {
529 for (int k = 0; k < dim; k++) {
530 if ((i != k) && (abs(K[i][k]) >= 1e-25)) {
531 return false;
532 }
533 }
534 }
535 return true;
536 }
537
538 template<class Problem, class VolumeVariables,
539 std::enable_if_t<!Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
540 static decltype(auto) getPermeability_(const Problem& problem,
541 const VolumeVariables& volVars,
542 const GlobalPosition& scvfIpGlobal)
543 { return volVars.permeability(); }
544
545 template<class Problem, class VolumeVariables,
546 std::enable_if_t<Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
547 static decltype(auto) getPermeability_(const Problem& problem,
548 const VolumeVariables& volVars,
549 const GlobalPosition& scvfIpGlobal)
550 { return problem.spatialParams().permeabilityAtPos(scvfIpGlobal); }
551};
552
557template<class ScalarType, class GridGeometry>
558class CCTpfaForchheimersLaw<ScalarType, GridGeometry, /*isNetwork*/ true>
559{
560 static_assert(AlwaysFalse<ScalarType>::value, "Forchheimer not implemented for network grids");
561};
562
563} // end namespace Dumux
564
565#endif
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Definition: adapt.hh:29
forward declare
Definition: forchheimerslaw.hh:38
Forchheimer's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: cctpfa/forchheimerslaw.hh:54
Declares all properties used in Dumux.
Darcy's law for cell-centered finite volume schemes with two-point flux approximation.