12#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FORCHHEIMERS_LAW_HH
13#define DUMUX_DISCRETIZATION_CC_TPFA_FORCHHEIMERS_LAW_HH
15#include <dune/common/fvector.hh>
16#include <dune/common/fmatrix.hh>
30template<
class TypeTag,
class ForchheimerVelocity,
class DiscretizationMethod>
31class ForchheimersLawImplementation;
42template<
class Scalar,
class Gr
idGeometry,
class ForchheimerVelocity,
bool isNetwork>
50template <
class TypeTag,
class ForchheimerVelocity>
53 GetPropType<TypeTag, Properties::GridGeometry>,
55 (GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension < GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld)>
62template<class GridGeometry>
63class TpfaForchheimersLawCacheFiller
65 using FVElementGeometry = typename GridGeometry::LocalView;
66 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
67 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
73 template<class FluxVariablesCache, class Problem, class ElementVolumeVariables, class FluxVariablesCacheFiller>
74 static void fill(FluxVariablesCache& scvfFluxVarsCache,
75 const Problem& problem,
76 const Element& element,
77 const FVElementGeometry& fvGeometry,
78 const ElementVolumeVariables& elemVolVars,
79 const SubControlVolumeFace& scvf,
80 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
82 scvfFluxVarsCache.updateAdvection(problem, element, fvGeometry, elemVolVars, scvf);
90template<class AdvectionType, class GridGeometry>
91class TpfaForchheimersLawCache
93 using Scalar = typename AdvectionType::Scalar;
94 using FVElementGeometry = typename GridGeometry::LocalView;
95 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
96 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
97 static constexpr int dimWorld = GridGeometry::GridView::dimensionworld;
98 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
101 using Filler = TpfaForchheimersLawCacheFiller<GridGeometry>;
103 template<class Problem, class ElementVolumeVariables>
104 void updateAdvection(const Problem& problem,
105 const Element& element,
106 const FVElementGeometry& fvGeometry,
107 const ElementVolumeVariables& elemVolVars,
108 const SubControlVolumeFace &scvf)
110 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
111 harmonicMeanSqrtK_ = AdvectionType::calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf);
114 const Scalar& advectionTij() const
117 const DimWorldMatrix& harmonicMeanSqrtPermeability() const
118 { return harmonicMeanSqrtK_; }
122 DimWorldMatrix harmonicMeanSqrtK_;
129template<class ScalarType, class GridGeometry, class ForchheimerVelocity>
130class CCTpfaForchheimersLaw<ScalarType, GridGeometry, ForchheimerVelocity, false>
132 using ThisType = CCTpfaForchheimersLaw<ScalarType, GridGeometry, ForchheimerVelocity, false>;
133 using FVElementGeometry = typename GridGeometry::LocalView;
134 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
135 using Extrusion = Extrusion_t<GridGeometry>;
136 using GridView = typename GridGeometry::GridView;
137 using Element = typename GridView::template Codim<0>::Entity;
139 using DimWorldVector = typename ForchheimerVelocity::DimWorldVector;
140 using DimWorldMatrix = typename ForchheimerVelocity::DimWorldMatrix;
142 using DarcysLaw = CCTpfaDarcysLaw<ScalarType, GridGeometry, false>;
146 using Scalar = ScalarType;
148 using DiscretizationMethod = DiscretizationMethods::CCTpfa;
150 static constexpr DiscretizationMethod discMethod{};
153 using Cache = TpfaForchheimersLawCache<ThisType, GridGeometry>;
162 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
163 static Scalar flux(const Problem& problem,
164 const Element& element,
165 const FVElementGeometry& fvGeometry,
166 const ElementVolumeVariables& elemVolVars,
167 const SubControlVolumeFace& scvf,
169 const ElementFluxVarsCache& elemFluxVarsCache)
173 Scalar darcyFlux = DarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
174 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.mobility(phaseIdx); };
175 DimWorldVector darcyVelocity = scvf.unitOuterNormal();
176 darcyVelocity *= ForchheimerVelocity::UpwindScheme::apply(elemVolVars, scvf, upwindTerm, darcyFlux, phaseIdx);
177 darcyVelocity /= Extrusion::area(fvGeometry, scvf);
179 const auto velocity = ForchheimerVelocity::velocity(fvGeometry,
183 elemFluxVarsCache[scvf].harmonicMeanSqrtPermeability(),
186 Scalar flux = velocity * scvf.unitOuterNormal();
187 flux *= Extrusion::area(fvGeometry, scvf);
193 template<class Problem, class ElementVolumeVariables>
194 static Scalar calculateTransmissibility(const Problem& problem,
195 const Element& element,
196 const FVElementGeometry& fvGeometry,
197 const ElementVolumeVariables& elemVolVars,
198 const SubControlVolumeFace& scvf)
200 return DarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
207 template<class Problem, class ElementVolumeVariables>
208 static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem& problem,
209 const ElementVolumeVariables& elemVolVars,
210 const SubControlVolumeFace& scvf)
212 return ForchheimerVelocity::calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf);
220template<class ScalarType, class GridGeometry, class ForchheimerVelocity>
221class CCTpfaForchheimersLaw<ScalarType, GridGeometry, ForchheimerVelocity, true>
223 static_assert(AlwaysFalse<ScalarType>::value, "Forchheimer not implemented for network grids");
Forchheimer's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: cctpfa/forchheimerslaw.hh:43
Forchheimer's law For a detailed description see dumux/flow/forchheimerslaw.hh.
Definition: forchheimervelocity.hh:39
forward declare
Definition: forchheimerslaw_fwd.hh:27
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
Darcy's law for cell-centered finite volume schemes with two-point flux approximation.
Define some often used mathematical functions.
The available discretization methods in Dumux.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.