28#ifndef DUMUX_DISCRETIZATION_BOX_DARCYS_LAW_HH
29#define DUMUX_DISCRETIZATION_BOX_DARCYS_LAW_HH
41template<
class TypeTag,
class DiscretizationMethod>
45template<
class Scalar,
class Gr
idGeometry>
52template<
class TypeTag>
54:
public BoxDarcysLaw<GetPropType<TypeTag, Properties::Scalar>, GetPropType<TypeTag, Properties::GridGeometry>>
63template<
class Scalar,
class Gr
idGeometry>
66 using FVElementGeometry =
typename GridGeometry::LocalView;
67 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
69 using GridView =
typename GridGeometry::GridView;
70 using Element =
typename GridView::template Codim<0>::Entity;
72 enum { dim = GridView::dimension};
73 enum { dimWorld = GridView::dimensionworld};
87 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache>
88 static Scalar
flux(
const Problem& problem,
89 const Element& element,
90 const FVElementGeometry& fvGeometry,
91 const ElementVolumeVariables& elemVolVars,
92 const SubControlVolumeFace& scvf,
94 const ElementFluxVarsCache& elemFluxVarCache)
96 const auto& fluxVarCache = elemFluxVarCache[scvf];
97 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
98 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
99 const auto& insideVolVars = elemVolVars[insideScv];
100 const auto& outsideVolVars = elemVolVars[outsideScv];
102 auto insideK = insideVolVars.permeability();
103 auto outsideK = outsideVolVars.permeability();
106 insideK *= insideVolVars.extrusionFactor();
107 outsideK *= outsideVolVars.extrusionFactor();
110 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(),
"Problem.EnableGravity");
112 const auto& shapeValues = fluxVarCache.shapeValues();
115 Dune::FieldVector<Scalar, dimWorld> gradP(0.0);
117 for (
auto&& scv : scvs(fvGeometry))
119 const auto& volVars = elemVolVars[scv];
122 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
125 gradP.axpy(volVars.pressure(phaseIdx), fluxVarCache.gradN(scv.indexInElement()));
129 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
132 return -1.0*
vtmv(scvf.unitOuterNormal(), K, gradP)*Extrusion::area(scvf);
136 template<
class Problem,
class ElementVolumeVariables,
class FluxVarCache>
138 const Element& element,
139 const FVElementGeometry& fvGeometry,
140 const ElementVolumeVariables& elemVolVars,
141 const SubControlVolumeFace& scvf,
142 const FluxVarCache& fluxVarCache)
144 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
145 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
146 const auto& insideVolVars = elemVolVars[insideScv];
147 const auto& outsideVolVars = elemVolVars[outsideScv];
149 auto insideK = insideVolVars.permeability();
150 auto outsideK = outsideVolVars.permeability();
153 insideK *= insideVolVars.extrusionFactor();
154 outsideK *= outsideVolVars.extrusionFactor();
158 std::vector<Scalar> ti(fvGeometry.numScv());
159 for (
const auto& scv : scvs(fvGeometry))
160 ti[scv.indexInElement()] =
161 -1.0*Extrusion::area(scvf)*
vtmv(scvf.unitOuterNormal(), K, fluxVarCache.gradN(scv.indexInElement()));
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::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:863
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
forward declaration of the method-specific implementation
Definition: flux/box/darcyslaw.hh:42
Darcy's law for the box scheme.
Definition: flux/box/darcyslaw.hh:65
static std::vector< Scalar > calculateTransmissibilities(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const FluxVarCache &fluxVarCache)
Definition: flux/box/darcyslaw.hh:137
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVarsCache &elemFluxVarCache)
Returns the advective flux of a fluid phase across the given sub-control volume face.
Definition: flux/box/darcyslaw.hh:88
Declares all properties used in Dumux.