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,
class ForchheimerVelocity,
class DiscretizationMethod>
43class ForchheimersLawImplementation;
54template<
class Scalar,
class Gr
idGeometry,
class ForchheimerVelocity,
bool isNetwork>
62template <
class TypeTag,
class ForchheimerVelocity>
65 GetPropType<TypeTag, Properties::GridGeometry>,
67 (GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension < GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld)>
74template<class GridGeometry>
75class TpfaForchheimersLawCacheFiller
77 using FVElementGeometry = typename GridGeometry::LocalView;
78 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
79 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
85 template<class FluxVariablesCache, class Problem, class ElementVolumeVariables, class FluxVariablesCacheFiller>
86 static void fill(FluxVariablesCache& scvfFluxVarsCache,
87 const Problem& problem,
88 const Element& element,
89 const FVElementGeometry& fvGeometry,
90 const ElementVolumeVariables& elemVolVars,
91 const SubControlVolumeFace& scvf,
92 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
94 scvfFluxVarsCache.updateAdvection(problem, element, fvGeometry, elemVolVars, scvf);
102template<class AdvectionType, class GridGeometry>
103class TpfaForchheimersLawCache
105 using Scalar = typename AdvectionType::Scalar;
106 using FVElementGeometry = typename GridGeometry::LocalView;
107 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
108 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
109 static constexpr int dimWorld = GridGeometry::GridView::dimensionworld;
110 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
113 using Filler = TpfaForchheimersLawCacheFiller<GridGeometry>;
115 template<class Problem, class ElementVolumeVariables>
116 void updateAdvection(const Problem& problem,
117 const Element& element,
118 const FVElementGeometry& fvGeometry,
119 const ElementVolumeVariables& elemVolVars,
120 const SubControlVolumeFace &scvf)
122 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
123 harmonicMeanSqrtK_ = AdvectionType::calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf);
126 const Scalar& advectionTij() const
129 const DimWorldMatrix& harmonicMeanSqrtPermeability() const
130 { return harmonicMeanSqrtK_; }
134 DimWorldMatrix harmonicMeanSqrtK_;
141template<class ScalarType, class GridGeometry, class ForchheimerVelocity>
142class CCTpfaForchheimersLaw<ScalarType, GridGeometry, ForchheimerVelocity, false>
144 using ThisType = CCTpfaForchheimersLaw<ScalarType, GridGeometry, ForchheimerVelocity, false>;
145 using FVElementGeometry = typename GridGeometry::LocalView;
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 using DimWorldVector = typename ForchheimerVelocity::DimWorldVector;
152 using DimWorldMatrix = typename ForchheimerVelocity::DimWorldMatrix;
154 using DarcysLaw = CCTpfaDarcysLaw<ScalarType, GridGeometry, false>;
158 using Scalar = ScalarType;
160 using DiscretizationMethod = DiscretizationMethods::CCTpfa;
162 static constexpr DiscretizationMethod discMethod{};
165 using Cache = TpfaForchheimersLawCache<ThisType, GridGeometry>;
174 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
175 static Scalar flux(const Problem& problem,
176 const Element& element,
177 const FVElementGeometry& fvGeometry,
178 const ElementVolumeVariables& elemVolVars,
179 const SubControlVolumeFace& scvf,
181 const ElementFluxVarsCache& elemFluxVarsCache)
185 Scalar darcyFlux = DarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
186 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.mobility(phaseIdx); };
187 DimWorldVector darcyVelocity = scvf.unitOuterNormal();
188 darcyVelocity *= ForchheimerVelocity::UpwindScheme::apply(elemVolVars, scvf, upwindTerm, darcyFlux, phaseIdx);
189 darcyVelocity /= Extrusion::area(fvGeometry, scvf);
191 const auto velocity = ForchheimerVelocity::velocity(fvGeometry,
195 elemFluxVarsCache[scvf].harmonicMeanSqrtPermeability(),
198 Scalar flux = velocity * scvf.unitOuterNormal();
199 flux *= Extrusion::area(fvGeometry, scvf);
205 template<class Problem, class ElementVolumeVariables>
206 static Scalar calculateTransmissibility(const Problem& problem,
207 const Element& element,
208 const FVElementGeometry& fvGeometry,
209 const ElementVolumeVariables& elemVolVars,
210 const SubControlVolumeFace& scvf)
212 return DarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
219 template<class Problem, class ElementVolumeVariables>
220 static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem& problem,
221 const ElementVolumeVariables& elemVolVars,
222 const SubControlVolumeFace& scvf)
224 return ForchheimerVelocity::calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf);
232template<class ScalarType, class GridGeometry, class ForchheimerVelocity>
233class CCTpfaForchheimersLaw<ScalarType, GridGeometry, ForchheimerVelocity, true>
235 static_assert(AlwaysFalse<ScalarType>::value, "Forchheimer not implemented for network grids");
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Define some often used mathematical functions.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
forward declare
Definition: forchheimerslaw_fwd.hh:39
Forchheimer's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: cctpfa/forchheimerslaw.hh:55
Forchheimer's law For a detailed description see dumux/flow/forchheimerslaw.hh.
Definition: forchheimervelocity.hh:51
Declares all properties used in Dumux.
Darcy's law for cell-centered finite volume schemes with two-point flux approximation.