28#ifndef DUMUX_DISCRETIZATION_BOX_DARCYS_LAW_HH
29#define DUMUX_DISCRETIZATION_BOX_DARCYS_LAW_HH
40template<
class TypeTag, DiscretizationMethod discMethod>
41class DarcysLawImplementation;
44template<
class Scalar,
class Gr
idGeometry>
51template<
class TypeTag>
53:
public BoxDarcysLaw<GetPropType<TypeTag, Properties::Scalar>, GetPropType<TypeTag, Properties::GridGeometry>>
62template<
class Scalar,
class Gr
idGeometry>
65 using FVElementGeometry =
typename GridGeometry::LocalView;
66 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
68 using GridView =
typename GridGeometry::GridView;
69 using Element =
typename GridView::template Codim<0>::Entity;
71 enum { dim = GridView::dimension};
72 enum { dimWorld = GridView::dimensionworld};
86 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache>
87 static Scalar
flux(
const Problem& problem,
88 const Element& element,
89 const FVElementGeometry& fvGeometry,
90 const ElementVolumeVariables& elemVolVars,
91 const SubControlVolumeFace& scvf,
93 const ElementFluxVarsCache& elemFluxVarCache)
95 const auto& fluxVarCache = elemFluxVarCache[scvf];
96 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
97 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
98 const auto& insideVolVars = elemVolVars[insideScv];
99 const auto& outsideVolVars = elemVolVars[outsideScv];
101 auto insideK = insideVolVars.permeability();
102 auto outsideK = outsideVolVars.permeability();
105 insideK *= insideVolVars.extrusionFactor();
106 outsideK *= outsideVolVars.extrusionFactor();
108 const auto K = problem.spatialParams().harmonicMean(insideK, outsideK, scvf.unitOuterNormal());
109 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(),
"Problem.EnableGravity");
111 const auto& shapeValues = fluxVarCache.shapeValues();
114 Dune::FieldVector<Scalar, dimWorld> gradP(0.0);
116 for (
auto&& scv : scvs(fvGeometry))
118 const auto& volVars = elemVolVars[scv];
121 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
124 gradP.axpy(volVars.pressure(phaseIdx), fluxVarCache.gradN(scv.indexInElement()));
128 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
131 return -1.0*
vtmv(scvf.unitOuterNormal(), K, gradP)*Extrusion::area(scvf);
135 template<
class Problem,
class ElementVolumeVariables,
class FluxVarCache>
137 const Element& element,
138 const FVElementGeometry& fvGeometry,
139 const ElementVolumeVariables& elemVolVars,
140 const SubControlVolumeFace& scvf,
141 const FluxVarCache& fluxVarCache)
143 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
144 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
145 const auto& insideVolVars = elemVolVars[insideScv];
146 const auto& outsideVolVars = elemVolVars[outsideScv];
148 auto insideK = insideVolVars.permeability();
149 auto outsideK = outsideVolVars.permeability();
152 insideK *= insideVolVars.extrusionFactor();
153 outsideK *= outsideVolVars.extrusionFactor();
155 const auto K = problem.spatialParams().harmonicMean(insideK, outsideK, scvf.unitOuterNormal());
157 std::vector<Scalar> ti(fvGeometry.numScv());
158 for (
const auto& scv : scvs(fvGeometry))
159 ti[scv.indexInElement()] =
160 -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.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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
forward declaration of the method-specific implementation
Definition: flux/darcyslaw.hh:38
Darcy's law for the box scheme.
Definition: flux/box/darcyslaw.hh:64
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:136
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:87
Declares all properties used in Dumux.