24#ifndef DUMUX_DISCRETIZATION_BOX_FORCHHEIMERS_LAW_HH
25#define DUMUX_DISCRETIZATION_BOX_FORCHHEIMERS_LAW_HH
27#include <dune/common/fvector.hh>
28#include <dune/common/fmatrix.hh>
43template<
class TypeTag,
class ForchheimerVelocity,
class DiscretizationMethod>
44class ForchheimersLawImplementation;
55template<
class Scalar,
class Gr
idGeometry,
class ForchheimerVelocity>
56class BoxForchheimersLaw;
62template <
class TypeTag,
class ForchheimerVelocity>
65 GetPropType<TypeTag, Properties::GridGeometry>,
73template<
class ScalarType,
class Gr
idGeometry,
class ForchheimerVelocity>
77 using FVElementGeometry =
typename GridGeometry::LocalView;
78 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
80 using GridView =
typename GridGeometry::GridView;
81 using Element =
typename GridView::template Codim<0>::Entity;
103 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache>
105 const Element& element,
106 const FVElementGeometry& fvGeometry,
107 const ElementVolumeVariables& elemVolVars,
108 const SubControlVolumeFace& scvf,
110 const ElementFluxVarsCache& elemFluxVarsCache)
112 const auto& fluxVarCache = elemFluxVarsCache[scvf];
113 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
114 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
115 const auto& insideVolVars = elemVolVars[insideScv];
116 const auto& outsideVolVars = elemVolVars[outsideScv];
118 auto insideK = insideVolVars.permeability();
119 auto outsideK = outsideVolVars.permeability();
122 insideK *= insideVolVars.extrusionFactor();
123 outsideK *= outsideVolVars.extrusionFactor();
126 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(),
"Problem.EnableGravity");
128 const auto& shapeValues = fluxVarCache.shapeValues();
131 DimWorldVector gradP(0.0);
133 for (
auto&& scv : scvs(fvGeometry))
135 const auto& volVars = elemVolVars[scv];
138 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
141 gradP.axpy(volVars.pressure(phaseIdx), fluxVarCache.gradN(scv.indexInElement()));
145 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
147 DimWorldVector darcyVelocity =
mv(K, gradP);
150 auto upwindTerm = [phaseIdx](
const auto& volVars){
return volVars.mobility(phaseIdx); };
151 darcyVelocity *= ForchheimerVelocity::UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm, darcyVelocity*scvf.unitOuterNormal(), phaseIdx);
160 Scalar flux = velocity * scvf.unitOuterNormal();
161 flux *= Extrusion::area(scvf);
167 template<
class Problem,
class ElementVolumeVariables>
169 const Element& element,
170 const FVElementGeometry& fvGeometry,
171 const ElementVolumeVariables& elemVolVars,
172 const SubControlVolumeFace& scvf)
174 return DarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
181 template<
class Problem,
class ElementVolumeVariables>
183 const ElementVolumeVariables& elemVolVars,
184 const SubControlVolumeFace& scvf)
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
A free function to average a Tensor at an interface.
Dune::DenseVector< V >::derived_type mv(const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V > &v)
Returns the result of the projection of a vector v with a Matrix M.
Definition: math.hh:812
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
Scalar faceTensorAverage(const Scalar T1, const Scalar T2, const Dune::FieldVector< Scalar, dim > &normal)
Average of a discontinuous scalar field at discontinuity interface (for compatibility reasons with th...
Definition: facetensoraverage.hh:41
Darcy's law for the box scheme.
Definition: flux/box/darcyslaw.hh:65
forward declare
Definition: forchheimerslaw_fwd.hh:39
Forchheimer's law for box scheme.
Definition: box/forchheimerslaw.hh:75
static Scalar calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: box/forchheimerslaw.hh:168
static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem &problem, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Returns the harmonic mean of and .
Definition: box/forchheimerslaw.hh:182
static constexpr DiscretizationMethod discMethod
state the discretization method this implementation belongs to
Definition: box/forchheimerslaw.hh:94
ScalarType Scalar
state the scalar type of the law
Definition: box/forchheimerslaw.hh:90
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, int phaseIdx, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the advective flux of a phase across the given sub-control volume face using the Forchheimer ...
Definition: box/forchheimerslaw.hh:104
Forchheimer's law This file contains the calculation of the Forchheimer velocity for a given Darcy ve...
Definition: forchheimervelocity.hh:49
static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem &problem, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Returns the harmonic mean of and .
Definition: forchheimervelocity.hh:194
Dune::FieldVector< Scalar, dimWorld > DimWorldVector
Definition: forchheimervelocity.hh:60
static DimWorldVector velocity(const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const DimWorldMatrix sqrtK, const DimWorldVector darcyVelocity)
Compute the Forchheimer velocity of a phase at the given sub-control volume face using the Forchheime...
Definition: forchheimervelocity.hh:123
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimWorldMatrix
Definition: forchheimervelocity.hh:61
Declares all properties used in Dumux.
Specialization of Darcy's Law for the box method.