3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
28#ifndef DUMUX_DISCRETIZATION_CVFE_DARCYS_LAW_HH
29#define DUMUX_DISCRETIZATION_CVFE_DARCYS_LAW_HH
30
31#include <dumux/common/math.hh>
37
38namespace Dumux {
39
46template<class Scalar, class GridGeometry>
48{
49 using FVElementGeometry = typename GridGeometry::LocalView;
50 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
51 using Extrusion = Extrusion_t<GridGeometry>;
52 using GridView = typename GridGeometry::GridView;
53 using Element = typename GridView::template Codim<0>::Entity;
54
55 static constexpr int dim = GridView::dimension;
56 static constexpr int dimWorld = GridView::dimensionworld;
57
58public:
59
70 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
71 static Scalar flux(const Problem& problem,
72 const Element& element,
73 const FVElementGeometry& fvGeometry,
74 const ElementVolumeVariables& elemVolVars,
75 const SubControlVolumeFace& scvf,
76 const int phaseIdx,
77 const ElementFluxVarsCache& elemFluxVarCache)
78 {
79 const auto& fluxVarCache = elemFluxVarCache[scvf];
80 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
81 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
82 const auto& insideVolVars = elemVolVars[insideScv];
83 const auto& outsideVolVars = elemVolVars[outsideScv];
84
85 auto insideK = insideVolVars.permeability();
86 auto outsideK = outsideVolVars.permeability();
87
88 // scale with correct extrusion factor
89 insideK *= insideVolVars.extrusionFactor();
90 outsideK *= outsideVolVars.extrusionFactor();
91
92 const auto K = faceTensorAverage(insideK, outsideK, scvf.unitOuterNormal());
93 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
94
95 const auto& shapeValues = fluxVarCache.shapeValues();
96
97 // evaluate gradP - rho*g at integration point
98 Dune::FieldVector<Scalar, dimWorld> gradP(0.0);
99 Scalar rho(0.0);
100 for (auto&& scv : scvs(fvGeometry))
101 {
102 const auto& volVars = elemVolVars[scv];
103
104 if (enableGravity)
105 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
106
107 // the global shape function gradient
108 gradP.axpy(volVars.pressure(phaseIdx), fluxVarCache.gradN(scv.indexInElement()));
109 }
110
111 if (enableGravity)
112 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
113
114 // apply the permeability and return the flux
115 return -1.0*vtmv(scvf.unitOuterNormal(), K, gradP)*Extrusion::area(fvGeometry, scvf);
116 }
117
118 // compute transmissibilities ti for analytical Jacobians
119 template<class Problem, class ElementVolumeVariables, class FluxVarCache>
120 static std::vector<Scalar> calculateTransmissibilities(const Problem& problem,
121 const Element& element,
122 const FVElementGeometry& fvGeometry,
123 const ElementVolumeVariables& elemVolVars,
124 const SubControlVolumeFace& scvf,
125 const FluxVarCache& fluxVarCache)
126 {
127 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
128 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
129 const auto& insideVolVars = elemVolVars[insideScv];
130 const auto& outsideVolVars = elemVolVars[outsideScv];
131
132 auto insideK = insideVolVars.permeability();
133 auto outsideK = outsideVolVars.permeability();
134
135 // scale with correct extrusion factor
136 insideK *= insideVolVars.extrusionFactor();
137 outsideK *= outsideVolVars.extrusionFactor();
138
139 const auto K = faceTensorAverage(insideK, outsideK, scvf.unitOuterNormal());
140
141 std::vector<Scalar> ti(fvGeometry.numScv());
142 for (const auto& scv : scvs(fvGeometry))
143 ti[scv.indexInElement()] =
144 -1.0*Extrusion::area(fvGeometry, scvf)*vtmv(scvf.unitOuterNormal(), K, fluxVarCache.gradN(scv.indexInElement()));
145
146 return ti;
147 }
148};
149
150// forward declaration
151template<class TypeTag, class DiscretizationMethod>
152class DarcysLawImplementation;
153
158template<class TypeTag, class DM>
159class DarcysLawImplementation<TypeTag, DiscretizationMethods::CVFE<DM>>
160: public CVFEDarcysLaw<GetPropType<TypeTag, Properties::Scalar>, GetPropType<TypeTag, Properties::GridGeometry>>
161{};
162
163} // end namespace Dumux
164
165#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Define some often used mathematical functions.
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
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:44
Darcy's law for control-volume finite element schemes.
Definition: flux/cvfe/darcyslaw.hh:48
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:120
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:71
Declares all properties used in Dumux.