version 3.10-dev
flux/cvfe/darcyslaw.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
16#ifndef DUMUX_DISCRETIZATION_CVFE_DARCYS_LAW_HH
17#define DUMUX_DISCRETIZATION_CVFE_DARCYS_LAW_HH
18
19#include <dumux/common/math.hh>
25
26namespace Dumux {
27
34template<class Scalar, class GridGeometry>
36{
37 using FVElementGeometry = typename GridGeometry::LocalView;
38 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
39 using Extrusion = Extrusion_t<GridGeometry>;
40 using GridView = typename GridGeometry::GridView;
41 using Element = typename GridView::template Codim<0>::Entity;
42
43 static constexpr int dim = GridView::dimension;
44 static constexpr int dimWorld = GridView::dimensionworld;
45
46public:
47
58 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
59 static Scalar flux(const Problem& problem,
60 const Element& element,
61 const FVElementGeometry& fvGeometry,
62 const ElementVolumeVariables& elemVolVars,
63 const SubControlVolumeFace& scvf,
64 const int phaseIdx,
65 const ElementFluxVarsCache& elemFluxVarCache)
66 {
67 const auto& fluxVarCache = elemFluxVarCache[scvf];
68 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
69 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
70 const auto& insideVolVars = elemVolVars[insideScv];
71 const auto& outsideVolVars = elemVolVars[outsideScv];
72
73 auto insideK = insideVolVars.permeability();
74 auto outsideK = outsideVolVars.permeability();
75
76 // scale with correct extrusion factor
77 insideK *= insideVolVars.extrusionFactor();
78 outsideK *= outsideVolVars.extrusionFactor();
79
80 const auto K = faceTensorAverage(insideK, outsideK, scvf.unitOuterNormal());
81 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
82
83 const auto& shapeValues = fluxVarCache.shapeValues();
84
85 // evaluate gradP - rho*g at integration point
86 Dune::FieldVector<Scalar, dimWorld> gradP(0.0);
87 Scalar rho(0.0);
88 for (auto&& scv : scvs(fvGeometry))
89 {
90 const auto& volVars = elemVolVars[scv];
91
92 if (enableGravity)
93 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
94
95 // the global shape function gradient
96 gradP.axpy(volVars.pressure(phaseIdx), fluxVarCache.gradN(scv.indexInElement()));
97 }
98
99 if (enableGravity)
100 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
101
102 // apply the permeability and return the flux
103 return -1.0*vtmv(scvf.unitOuterNormal(), K, gradP)*Extrusion::area(fvGeometry, scvf);
104 }
105
106 // compute transmissibilities ti for analytical Jacobians
107 template<class Problem, class ElementVolumeVariables, class FluxVarCache>
108 static std::vector<Scalar> calculateTransmissibilities(const Problem& problem,
109 const Element& element,
110 const FVElementGeometry& fvGeometry,
111 const ElementVolumeVariables& elemVolVars,
112 const SubControlVolumeFace& scvf,
113 const FluxVarCache& fluxVarCache)
114 {
115 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
116 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
117 const auto& insideVolVars = elemVolVars[insideScv];
118 const auto& outsideVolVars = elemVolVars[outsideScv];
119
120 auto insideK = insideVolVars.permeability();
121 auto outsideK = outsideVolVars.permeability();
122
123 // scale with correct extrusion factor
124 insideK *= insideVolVars.extrusionFactor();
125 outsideK *= outsideVolVars.extrusionFactor();
126
127 const auto K = faceTensorAverage(insideK, outsideK, scvf.unitOuterNormal());
128
129 std::vector<Scalar> ti(fvGeometry.numScv());
130 for (const auto& scv : scvs(fvGeometry))
131 ti[scv.indexInElement()] =
132 -1.0*Extrusion::area(fvGeometry, scvf)*vtmv(scvf.unitOuterNormal(), K, fluxVarCache.gradN(scv.indexInElement()));
133
134 return ti;
135 }
136};
137
138// forward declaration
139template<class TypeTag, class DiscretizationMethod>
140class DarcysLawImplementation;
141
146template<class TypeTag, class DM>
147class DarcysLawImplementation<TypeTag, DiscretizationMethods::CVFE<DM>>
148: public CVFEDarcysLaw<GetPropType<TypeTag, Properties::Scalar>, GetPropType<TypeTag, Properties::GridGeometry>>
149{};
150
151} // end namespace Dumux
152
153#endif
Darcy's law for control-volume finite element schemes.
Definition: flux/cvfe/darcyslaw.hh:36
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/cvfe/darcyslaw.hh:108
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/cvfe/darcyslaw.hh:59
forward declaration of the method-specific implementation
Definition: flux/ccmpfa/darcyslaw.hh:27
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
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:880
Define some often used mathematical functions.
The available discretization methods in Dumux.
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
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:29
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.