24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FORCHHEIMERS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FORCHHEIMERS_LAW_HH
27#include <dune/common/fvector.hh>
28#include <dune/common/fmatrix.hh>
42template<
class TypeTag, DiscretizationMethod discMethod>
43class ForchheimersLawImplementation;
53template<
class Scalar,
class Gr
idGeometry,
class UpwindScheme,
bool isNetwork>
61template <
class TypeTag>
64 GetPropType<TypeTag, Properties::GridGeometry>,
65 typename GetPropType<TypeTag, Properties::FluxVariables>::UpwindScheme,
66 (GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension < GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld)>
73template<class GridGeometry>
74class TpfaForchheimersLawCacheFiller
76 using FVElementGeometry = typename GridGeometry::LocalView;
77 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
78 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
84 template<class FluxVariablesCache, class Problem, class ElementVolumeVariables, class FluxVariablesCacheFiller>
85 static void fill(FluxVariablesCache& scvfFluxVarsCache,
86 const Problem& problem,
87 const Element& element,
88 const FVElementGeometry& fvGeometry,
89 const ElementVolumeVariables& elemVolVars,
90 const SubControlVolumeFace& scvf,
91 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
93 scvfFluxVarsCache.updateAdvection(problem, element, fvGeometry, elemVolVars, scvf);
101template<class AdvectionType, class GridGeometry>
102class TpfaForchheimersLawCache
104 using Scalar = typename AdvectionType::Scalar;
105 using FVElementGeometry = typename GridGeometry::LocalView;
106 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
107 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
108 static constexpr int dimWorld = GridGeometry::GridView::dimensionworld;
109 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
112 using Filler = TpfaForchheimersLawCacheFiller<GridGeometry>;
114 template<class Problem, class ElementVolumeVariables>
115 void updateAdvection(const Problem& problem,
116 const Element& element,
117 const FVElementGeometry& fvGeometry,
118 const ElementVolumeVariables& elemVolVars,
119 const SubControlVolumeFace &scvf)
121 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
122 harmonicMeanSqrtK_ = AdvectionType::calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf);
125 const Scalar& advectionTij() const
128 const DimWorldMatrix& harmonicMeanSqrtPermeability() const
129 { return harmonicMeanSqrtK_; }
133 DimWorldMatrix harmonicMeanSqrtK_;
140template<class ScalarType, class GridGeometry, class UpwindScheme>
141class CCTpfaForchheimersLaw<ScalarType, GridGeometry, UpwindScheme, false>
143 using ThisType = CCTpfaForchheimersLaw<ScalarType, GridGeometry, UpwindScheme, false>;
144 using FVElementGeometry = typename GridGeometry::LocalView;
145 using SubControlVolume = typename GridGeometry::SubControlVolume;
146 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
147 using Extrusion = Extrusion_t<GridGeometry>;
148 using GridView = typename GridGeometry::GridView;
149 using Element = typename GridView::template Codim<0>::Entity;
151 static constexpr int dim = GridView::dimension;
152 static constexpr int dimWorld = GridView::dimensionworld;
154 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
155 using DimWorldVector = Dune::FieldVector<ScalarType, dimWorld>;
156 using DimWorldMatrix = Dune::FieldMatrix<ScalarType, dimWorld, dimWorld>;
158 using DarcysLaw = CCTpfaDarcysLaw<ScalarType, GridGeometry, false>;
162 using Scalar = ScalarType;
165 static const DiscretizationMethod discMethod = DiscretizationMethod::cctpfa;
168 using Cache = TpfaForchheimersLawCache<ThisType, GridGeometry>;
232 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
233 static Scalar flux(const Problem& problem,
234 const Element& element,
235 const FVElementGeometry& fvGeometry,
236 const ElementVolumeVariables& elemVolVars,
237 const SubControlVolumeFace& scvf,
239 const ElementFluxVarsCache& elemFluxVarsCache)
241 const auto velocity = velocity_(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
242 Scalar flux = velocity * scvf.unitOuterNormal();
243 flux *= Extrusion::area(scvf);
249 template<class Problem, class ElementVolumeVariables>
250 static Scalar calculateTransmissibility(const Problem& problem,
251 const Element& element,
252 const FVElementGeometry& fvGeometry,
253 const ElementVolumeVariables& elemVolVars,
254 const SubControlVolumeFace& scvf)
256 return DarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
263 template<class Problem, class ElementVolumeVariables,
264 bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value,
265 std::enable_if_t<scalarPerm, int> = 0>
266 static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem& problem,
267 const ElementVolumeVariables& elemVolVars,
268 const SubControlVolumeFace& scvf)
271 Scalar harmonicMeanSqrtK(0.0);
273 const auto insideScvIdx = scvf.insideScvIdx();
274 const auto& insideVolVars = elemVolVars[insideScvIdx];
275 const Scalar Ki = getPermeability_(problem, insideVolVars, scvf.ipGlobal());
276 const Scalar sqrtKi = sqrt(Ki);
278 if (!scvf.boundary())
280 const auto outsideScvIdx = scvf.outsideScvIdx();
281 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
282 const Scalar Kj = getPermeability_(problem, outsideVolVars, scvf.ipGlobal());
283 const Scalar sqrtKj = sqrt(Kj);
284 harmonicMeanSqrtK = problem.spatialParams().harmonicMean(sqrtKi, sqrtKj, scvf.unitOuterNormal());
287 harmonicMeanSqrtK = sqrtKi;
290 DimWorldMatrix result(0.0);
291 for (int i = 0; i < dimWorld; ++i)
292 result[i][i] = harmonicMeanSqrtK;
301 template<class Problem, class ElementVolumeVariables,
302 bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value,
303 std::enable_if_t<!scalarPerm, int> = 0>
304 static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem& problem,
305 const ElementVolumeVariables& elemVolVars,
306 const SubControlVolumeFace& scvf)
309 DimWorldMatrix sqrtKi(0.0);
310 DimWorldMatrix sqrtKj(0.0);
311 DimWorldMatrix harmonicMeanSqrtK(0.0);
313 const auto insideScvIdx = scvf.insideScvIdx();
314 const auto& insideVolVars = elemVolVars[insideScvIdx];
315 const auto& Ki = getPermeability_(problem, insideVolVars, scvf.ipGlobal());
318 if (!isDiagonal_(Ki))
319 DUNE_THROW(Dune::NotImplemented, "Only diagonal permeability tensors are supported.");
321 for (int i = 0; i < dim; ++i)
322 sqrtKi[i][i] = sqrt(Ki[i][i]);
324 if (!scvf.boundary())
326 const auto outsideScvIdx = scvf.outsideScvIdx();
327 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
328 const auto& Kj = getPermeability_(problem, outsideVolVars, scvf.ipGlobal());
330 if (!isDiagonal_(Kj))
331 DUNE_THROW(Dune::NotImplemented, "Only diagonal permeability tensors are supported.");
333 for (int i = 0; i < dim; ++i)
334 sqrtKj[i][i] = sqrt(Kj[i][i]);
336 harmonicMeanSqrtK = problem.spatialParams().harmonicMean(sqrtKi, sqrtKj, scvf.unitOuterNormal());
339 harmonicMeanSqrtK = sqrtKi;
341 return harmonicMeanSqrtK;
347 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
348 static DimWorldVector velocity_(const Problem& problem,
349 const Element& element,
350 const FVElementGeometry& fvGeometry,
351 const ElementVolumeVariables& elemVolVars,
352 const SubControlVolumeFace& scvf,
354 const ElementFluxVarsCache& elemFluxVarsCache)
358 const auto flux = DarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
359 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.mobility(phaseIdx); };
360 const auto darcyFlux = UpwindScheme::apply(elemVolVars, scvf, upwindTerm, flux, phaseIdx);
363 DimWorldVector darcyVelocity = scvf.unitOuterNormal();
364 darcyVelocity *= darcyFlux / Extrusion::area(scvf);
367 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
368 const auto& sqrtK = fluxVarsCache.harmonicMeanSqrtPermeability();
371 const Scalar forchCoeff = problem.spatialParams().forchCoeff(scvf);
374 DimWorldVector velocity = darcyVelocity;
376 DimWorldVector deltaV(0.0);
377 DimWorldVector residual(10e10);
378 DimWorldMatrix gradF(0.0);
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)
386 DUNE_THROW(NumericalProblem, "could not determine forchheimer velocity within "
387 << k << " iterations");
390 forchheimerResidual_(residual,
400 forchheimerDerivative_(gradF,
409 gradF.solve(deltaV, residual);
415 const auto forchheimerUpwindMobility = UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm, velocity * scvf.unitOuterNormal(), phaseIdx);
418 if (forchheimerUpwindMobility > 0.0)
419 velocity /= forchheimerUpwindMobility;
425 template<class ElementVolumeVariables>
426 static void forchheimerResidual_(DimWorldVector& residual,
427 const Scalar forchCoeff,
428 const DimWorldMatrix& sqrtK,
429 const DimWorldVector& darcyVelocity,
430 const DimWorldVector& oldVelocity,
431 const ElementVolumeVariables& elemVolVars,
432 const SubControlVolumeFace& scvf,
435 residual = oldVelocity;
438 residual -= darcyVelocity;
441 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.density(phaseIdx) / volVars.viscosity(phaseIdx); };
442 const Scalar rhoOverMu = UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm, oldVelocity * scvf.unitOuterNormal(), phaseIdx);
443 sqrtK.usmv(forchCoeff * rhoOverMu * oldVelocity.two_norm(), oldVelocity, residual);
447 template<class ElementVolumeVariables>
448 static void forchheimerDerivative_(DimWorldMatrix& derivative,
449 const Scalar forchCoeff,
450 const DimWorldMatrix& sqrtK,
451 const DimWorldVector& velocity,
452 const ElementVolumeVariables& elemVolVars,
453 const SubControlVolumeFace& scvf,
461 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.density(phaseIdx) / volVars.viscosity(phaseIdx); };
463 const Scalar absV = velocity.two_norm() ;
473 const Scalar rhoOverMu = UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm, velocity * scvf.unitOuterNormal(), phaseIdx);
476 for (int i = 0; i < dim; i++)
478 for (int k = 0; k <= i; k++)
480 derivative[i][k] = sqrtK[i][i] * velocity[i] * velocity[k] * forchCoeff
484 derivative[k][i] = derivative[i][k];
491 for (int i = 0; i < dim; i++)
492 derivative[i][i] += (1.0 + (sqrtK[i][i]*forchCoeff * absV * rhoOverMu));
503 static bool isDiagonal_(const DimWorldMatrix & K)
506 for (int i = 0; i < dim; i++) {
507 for (int k = 0; k < dim; k++) {
508 if ((i != k) && (abs(K[i][k]) >= 1e-25)) {
516 template<class Problem, class VolumeVariables,
517 std::enable_if_t<!Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
518 static decltype(auto) getPermeability_(const Problem& problem,
519 const VolumeVariables& volVars,
520 const GlobalPosition& scvfIpGlobal)
521 { return volVars.permeability(); }
523 template<class Problem, class VolumeVariables,
524 std::enable_if_t<Problem::SpatialParams::evaluatePermeabilityAtScvfIP(), int> = 0>
525 static decltype(auto) getPermeability_(const Problem& problem,
526 const VolumeVariables& volVars,
527 const GlobalPosition& scvfIpGlobal)
528 { return problem.spatialParams().permeabilityAtPos(scvfIpGlobal); }
535template<class ScalarType, class GridGeometry, class UpwindScheme>
536class CCTpfaForchheimersLaw<ScalarType, GridGeometry, UpwindScheme, true>
538 static_assert(AlwaysFalse<ScalarType>::value, "Forchheimer not implemented for network grids");
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
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.