28#ifndef DUMUX_DISCRETIZATION_BOX_DARCYS_LAW_HH
29#define DUMUX_DISCRETIZATION_BOX_DARCYS_LAW_HH
39template<
class TypeTag, DiscretizationMethod discMethod>
40class DarcysLawImplementation;
43template<
class Scalar,
class Gr
idGeometry>
50template<
class TypeTag>
52:
public BoxDarcysLaw<GetPropType<TypeTag, Properties::Scalar>, GetPropType<TypeTag, Properties::GridGeometry>>
61template<
class Scalar,
class Gr
idGeometry>
64 using FVElementGeometry =
typename GridGeometry::LocalView;
65 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
66 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
67 using GridView =
typename GridGeometry::GridView;
68 using Element =
typename GridView::template Codim<0>::Entity;
69 using CoordScalar =
typename GridView::ctype;
71 enum { dim = GridView::dimension};
72 enum { dimWorld = GridView::dimensionworld};
74 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
78 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache>
79 static Scalar
flux(
const Problem& problem,
80 const Element& element,
81 const FVElementGeometry& fvGeometry,
82 const ElementVolumeVariables& elemVolVars,
83 const SubControlVolumeFace& scvf,
85 const ElementFluxVarsCache& elemFluxVarCache)
87 const auto& fluxVarCache = elemFluxVarCache[scvf];
88 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
89 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
90 const auto& insideVolVars = elemVolVars[insideScv];
91 const auto& outsideVolVars = elemVolVars[outsideScv];
93 auto insideK = insideVolVars.permeability();
94 auto outsideK = outsideVolVars.permeability();
97 insideK *= insideVolVars.extrusionFactor();
98 outsideK *= outsideVolVars.extrusionFactor();
100 const auto K = problem.spatialParams().harmonicMean(insideK, outsideK, scvf.unitOuterNormal());
101 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(),
"Problem.EnableGravity");
103 const auto& shapeValues = fluxVarCache.shapeValues();
106 Dune::FieldVector<Scalar, dimWorld> gradP(0.0);
108 for (
auto&& scv : scvs(fvGeometry))
110 const auto& volVars = elemVolVars[scv];
113 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
116 gradP.axpy(volVars.pressure(phaseIdx), fluxVarCache.gradN(scv.indexInElement()));
120 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
123 return -1.0*
vtmv(scvf.unitOuterNormal(), K, gradP)*scvf.area();
127 template<
class Problem,
class ElementVolumeVariables,
class FluxVarCache>
129 const Element& element,
130 const FVElementGeometry& fvGeometry,
131 const ElementVolumeVariables& elemVolVars,
132 const SubControlVolumeFace& scvf,
133 const FluxVarCache& fluxVarCache)
135 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
136 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
137 const auto& insideVolVars = elemVolVars[insideScv];
138 const auto& outsideVolVars = elemVolVars[outsideScv];
140 auto insideK = insideVolVars.permeability();
141 auto outsideK = outsideVolVars.permeability();
144 insideK *= insideVolVars.extrusionFactor();
145 outsideK *= outsideVolVars.extrusionFactor();
147 const auto K = problem.spatialParams().harmonicMean(insideK, outsideK, scvf.unitOuterNormal());
149 std::vector<Scalar> ti(fvGeometry.numScv());
150 for (
const auto& scv : scvs(fvGeometry))
151 ti[scv.indexInElement()] =
152 -1.0*scvf.area()*
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.
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:840
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
forward declaration of the method-specific implementation
Definition: flux/darcyslaw.hh:38
Darcy's law for box schemes.
Definition: flux/box/darcyslaw.hh:63
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:128
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVarsCache &elemFluxVarCache)
Definition: flux/box/darcyslaw.hh:79
Declares all properties used in Dumux.