3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
flux/box/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_BOX_DARCYS_LAW_HH
29#define DUMUX_DISCRETIZATION_BOX_DARCYS_LAW_HH
30
31#warning "This header is deprecated and will be removed after release 3.6. Use flux/cvfe/darcyslaw.hh"
32
33#include <dumux/common/math.hh>
39
40namespace Dumux {
41
42// forward declaration
43template<class TypeTag, class DiscretizationMethod>
45
46// forward declaration
47template<class Scalar, class GridGeometry>
48class BoxDarcysLaw;
49
54template<class TypeTag>
55class DarcysLawImplementation<TypeTag, DiscretizationMethods::Box>
56: public BoxDarcysLaw<GetPropType<TypeTag, Properties::Scalar>, GetPropType<TypeTag, Properties::GridGeometry>>
57{ };
58
65template<class Scalar, class GridGeometry>
67{
68 using FVElementGeometry = typename GridGeometry::LocalView;
69 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
70 using Extrusion = Extrusion_t<GridGeometry>;
71 using GridView = typename GridGeometry::GridView;
72 using Element = typename GridView::template Codim<0>::Entity;
73
74 enum { dim = GridView::dimension};
75 enum { dimWorld = GridView::dimensionworld};
76
77public:
78
89 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
90 static Scalar flux(const Problem& problem,
91 const Element& element,
92 const FVElementGeometry& fvGeometry,
93 const ElementVolumeVariables& elemVolVars,
94 const SubControlVolumeFace& scvf,
95 const int phaseIdx,
96 const ElementFluxVarsCache& elemFluxVarCache)
97 {
98 const auto& fluxVarCache = elemFluxVarCache[scvf];
99 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
100 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
101 const auto& insideVolVars = elemVolVars[insideScv];
102 const auto& outsideVolVars = elemVolVars[outsideScv];
103
104 auto insideK = insideVolVars.permeability();
105 auto outsideK = outsideVolVars.permeability();
106
107 // scale with correct extrusion factor
108 insideK *= insideVolVars.extrusionFactor();
109 outsideK *= outsideVolVars.extrusionFactor();
110
111 const auto K = faceTensorAverage(insideK, outsideK, scvf.unitOuterNormal());
112 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
113
114 const auto& shapeValues = fluxVarCache.shapeValues();
115
116 // evaluate gradP - rho*g at integration point
117 Dune::FieldVector<Scalar, dimWorld> gradP(0.0);
118 Scalar rho(0.0);
119 for (auto&& scv : scvs(fvGeometry))
120 {
121 const auto& volVars = elemVolVars[scv];
122
123 if (enableGravity)
124 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
125
126 // the global shape function gradient
127 gradP.axpy(volVars.pressure(phaseIdx), fluxVarCache.gradN(scv.indexInElement()));
128 }
129
130 if (enableGravity)
131 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
132
133 // apply the permeability and return the flux
134 return -1.0*vtmv(scvf.unitOuterNormal(), K, gradP)*Extrusion::area(fvGeometry, scvf);
135 }
136
137 // compute transmissibilities ti for analytical jacobians
138 template<class Problem, class ElementVolumeVariables, class FluxVarCache>
139 static std::vector<Scalar> calculateTransmissibilities(const Problem& problem,
140 const Element& element,
141 const FVElementGeometry& fvGeometry,
142 const ElementVolumeVariables& elemVolVars,
143 const SubControlVolumeFace& scvf,
144 const FluxVarCache& fluxVarCache)
145 {
146 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
147 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
148 const auto& insideVolVars = elemVolVars[insideScv];
149 const auto& outsideVolVars = elemVolVars[outsideScv];
150
151 auto insideK = insideVolVars.permeability();
152 auto outsideK = outsideVolVars.permeability();
153
154 // scale with correct extrusion factor
155 insideK *= insideVolVars.extrusionFactor();
156 outsideK *= outsideVolVars.extrusionFactor();
157
158 const auto K = faceTensorAverage(insideK, outsideK, scvf.unitOuterNormal());
159
160 std::vector<Scalar> ti(fvGeometry.numScv());
161 for (const auto& scv : scvs(fvGeometry))
162 ti[scv.indexInElement()] =
163 -1.0*Extrusion::area(fvGeometry, scvf)*vtmv(scvf.unitOuterNormal(), K, fluxVarCache.gradN(scv.indexInElement()));
164
165 return ti;
166 }
167};
168
169} // end namespace Dumux
170
171#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
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:83
forward declaration of the method-specific implementation
Definition: flux/box/darcyslaw.hh:44
Darcy's law for the box scheme.
Definition: flux/box/darcyslaw.hh:67
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:139
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:90
Declares all properties used in Dumux.