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,
bool isNetwork>
61template <
class TypeTag>
64 GetPropType<TypeTag, Properties::GridGeometry>,
65 (GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension < GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld)>
72template<class GridGeometry>
73class TpfaForchheimersLawCacheFiller
75 using FVElementGeometry = typename GridGeometry::LocalView;
76 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
77 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
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)
92 scvfFluxVarsCache.updateAdvection(problem, element, fvGeometry, elemVolVars, scvf);
100template<class AdvectionType, class GridGeometry>
101class TpfaForchheimersLawCache
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>;
111 using Filler = TpfaForchheimersLawCacheFiller<GridGeometry>;
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)
120 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
121 harmonicMeanSqrtK_ = AdvectionType::calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf);
124 const Scalar& advectionTij() const
127 const DimWorldMatrix& harmonicMeanSqrtPermeability() const
128 { return harmonicMeanSqrtK_; }
132 DimWorldMatrix harmonicMeanSqrtK_;
139template<class ScalarType, class GridGeometry>
140class CCTpfaForchheimersLaw<ScalarType, GridGeometry, false>
142 using ThisType = CCTpfaForchheimersLaw<ScalarType, GridGeometry, 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;
150 static constexpr int dim = GridView::dimension;
151 static constexpr int dimWorld = GridView::dimensionworld;
153 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
154 using DimWorldVector = Dune::FieldVector<ScalarType, dimWorld>;
155 using DimWorldMatrix = Dune::FieldMatrix<ScalarType, dimWorld, dimWorld>;
157 using DarcysLaw = CCTpfaDarcysLaw<ScalarType, GridGeometry, false>;
161 using Scalar = ScalarType;
164 static const DiscretizationMethod discMethod = DiscretizationMethod::cctpfa;
167 using Cache = TpfaForchheimersLawCache<ThisType, GridGeometry>;
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,
238 const ElementFluxVarsCache& elemFluxVarsCache)
240 const auto velocity = velocity_(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
241 Scalar flux = velocity * scvf.unitOuterNormal();
242 flux *= Extrusion::area(scvf);
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)
255 return DarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
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)
270 Scalar harmonicMeanSqrtK(0.0);
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);
277 if (!scvf.boundary())
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());
286 harmonicMeanSqrtK = sqrtKi;
289 DimWorldMatrix result(0.0);
290 for (int i = 0; i < dimWorld; ++i)
291 result[i][i] = harmonicMeanSqrtK;
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)
308 DimWorldMatrix sqrtKi(0.0);
309 DimWorldMatrix sqrtKj(0.0);
310 DimWorldMatrix harmonicMeanSqrtK(0.0);
312 const auto insideScvIdx = scvf.insideScvIdx();
313 const auto& insideVolVars = elemVolVars[insideScvIdx];
314 const auto& Ki = getPermeability_(problem, insideVolVars, scvf.ipGlobal());
317 if (!isDiagonal_(Ki))
318 DUNE_THROW(Dune::NotImplemented, "Only diagonal permeability tensors are supported.");
320 for (int i = 0; i < dim; ++i)
321 sqrtKi[i][i] = sqrt(Ki[i][i]);
323 if (!scvf.boundary())
325 const auto outsideScvIdx = scvf.outsideScvIdx();
326 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
327 const auto& Kj = getPermeability_(problem, outsideVolVars, scvf.ipGlobal());
329 if (!isDiagonal_(Kj))
330 DUNE_THROW(Dune::NotImplemented, "Only diagonal permeability tensors are supported.");
332 for (int i = 0; i < dim; ++i)
333 sqrtKj[i][i] = sqrt(Kj[i][i]);
335 harmonicMeanSqrtK = problem.spatialParams().harmonicMean(sqrtKi, sqrtKj, scvf.unitOuterNormal());
338 harmonicMeanSqrtK = sqrtKi;
340 return harmonicMeanSqrtK;
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,
353 const ElementFluxVarsCache& elemFluxVarsCache)
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) );
360 darcyFlux *= darcyUpwindMobility;
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 Scalar forchheimerUpwindMobility = doUpwind_(scvf, elemVolVars, upwindTerm,
416 !std::signbit(velocity * scvf.unitOuterNormal()) );
419 if (forchheimerUpwindMobility > 0.0)
420 velocity /= forchheimerUpwindMobility;
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,
436 residual = oldVelocity;
439 residual -= darcyVelocity;
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);
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,
463 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.density(phaseIdx) / volVars.viscosity(phaseIdx); };
464 bool insideIsUpstream = !std::signbit(velocity * scvf.unitOuterNormal());
466 const Scalar absV = velocity.two_norm() ;
476 const Scalar rhoOverMu = doUpwind_(scvf, elemVolVars, upwindTerm, insideIsUpstream);
479 for (int i = 0; i < dim; i++)
481 for (int k = 0; k <= i; k++)
483 derivative[i][k] = sqrtK[i][i] * velocity[i] * velocity[k] * forchCoeff
487 derivative[k][i] = derivative[i][k];
494 for (int i = 0; i < dim; i++)
495 derivative[i][i] += (1.0 + (sqrtK[i][i]*forchCoeff * absV * rhoOverMu));
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)
504 static const Scalar upwindWeight = getParam<Scalar>("Flux.UpwindWeight");
506 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
507 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
509 if (!insideIsUpstream)
510 return (upwindWeight*upwindTerm(outsideVolVars)
511 + (1.0 - upwindWeight)*upwindTerm(insideVolVars));
513 return (upwindWeight*upwindTerm(insideVolVars)
514 + (1.0 - upwindWeight)*upwindTerm(outsideVolVars));
525 static bool isDiagonal_(const DimWorldMatrix & K)
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)) {
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(); }
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); }
557template<class ScalarType, class GridGeometry>
558class CCTpfaForchheimersLaw<ScalarType, GridGeometry, true>
560 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.