3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
box/forchheimerslaw.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 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_BOX_FORCHHEIMERS_LAW_HH
25#define DUMUX_DISCRETIZATION_BOX_FORCHHEIMERS_LAW_HH
26
27#include <dune/common/fvector.hh>
28#include <dune/common/fmatrix.hh>
29
30#include <dumux/common/math.hh>
34
39
40namespace Dumux {
41
42// forward declarations
43template<class TypeTag, class ForchheimerVelocity, class DiscretizationMethod>
44class ForchheimersLawImplementation;
45
55template<class Scalar, class GridGeometry, class ForchheimerVelocity>
56class BoxForchheimersLaw;
57
62template <class TypeTag, class ForchheimerVelocity>
63class ForchheimersLawImplementation<TypeTag, ForchheimerVelocity, DiscretizationMethods::Box>
64: public BoxForchheimersLaw<GetPropType<TypeTag, Properties::Scalar>,
65 GetPropType<TypeTag, Properties::GridGeometry>,
66 ForchheimerVelocity>
67{};
68
73template<class ScalarType, class GridGeometry, class ForchheimerVelocity>
75{
77 using FVElementGeometry = typename GridGeometry::LocalView;
78 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
79 using Extrusion = Extrusion_t<GridGeometry>;
80 using GridView = typename GridGeometry::GridView;
81 using Element = typename GridView::template Codim<0>::Entity;
82
83 using DimWorldVector = typename ForchheimerVelocity::DimWorldVector;
84 using DimWorldMatrix = typename ForchheimerVelocity::DimWorldMatrix;
85
87
88public:
90 using Scalar = ScalarType;
91
94 static constexpr DiscretizationMethod discMethod{};
95
103 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
104 static Scalar flux(const Problem& problem,
105 const Element& element,
106 const FVElementGeometry& fvGeometry,
107 const ElementVolumeVariables& elemVolVars,
108 const SubControlVolumeFace& scvf,
109 int phaseIdx,
110 const ElementFluxVarsCache& elemFluxVarsCache)
111 {
112 const auto& fluxVarCache = elemFluxVarsCache[scvf];
113 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
114 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
115 const auto& insideVolVars = elemVolVars[insideScv];
116 const auto& outsideVolVars = elemVolVars[outsideScv];
117
118 auto insideK = insideVolVars.permeability();
119 auto outsideK = outsideVolVars.permeability();
120
121 // scale with correct extrusion factor
122 insideK *= insideVolVars.extrusionFactor();
123 outsideK *= outsideVolVars.extrusionFactor();
124
125 const auto K = faceTensorAverage(insideK, outsideK, scvf.unitOuterNormal());
126 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
127
128 const auto& shapeValues = fluxVarCache.shapeValues();
129
130 // evaluate gradP - rho*g at integration point
131 DimWorldVector gradP(0.0);
132 Scalar rho(0.0);
133 for (auto&& scv : scvs(fvGeometry))
134 {
135 const auto& volVars = elemVolVars[scv];
136
137 if (enableGravity)
138 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
139
140 // the global shape function gradient
141 gradP.axpy(volVars.pressure(phaseIdx), fluxVarCache.gradN(scv.indexInElement()));
142 }
143
144 if (enableGravity)
145 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
146
147 DimWorldVector darcyVelocity = mv(K, gradP);
148 darcyVelocity *= -1;
149
150 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.mobility(phaseIdx); };
151 darcyVelocity *= ForchheimerVelocity::UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm, darcyVelocity*scvf.unitOuterNormal(), phaseIdx);
152
153 const auto velocity = ForchheimerVelocity::velocity(fvGeometry,
154 elemVolVars,
155 scvf,
156 phaseIdx,
157 calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf),
158 darcyVelocity);
159
160 Scalar flux = velocity * scvf.unitOuterNormal();
161 flux *= Extrusion::area(fvGeometry, scvf);
162 return flux;
163 }
164
167 template<class Problem, class ElementVolumeVariables>
168 static Scalar calculateTransmissibility(const Problem& problem,
169 const Element& element,
170 const FVElementGeometry& fvGeometry,
171 const ElementVolumeVariables& elemVolVars,
172 const SubControlVolumeFace& scvf)
173 {
174 return DarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
175 }
176
181 template<class Problem, class ElementVolumeVariables>
182 static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem& problem,
183 const ElementVolumeVariables& elemVolVars,
184 const SubControlVolumeFace& scvf)
185 {
186 return ForchheimerVelocity::calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf);
187 }
188};
189
190} // end namespace Dumux
191
192#endif
Type traits.
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::DenseVector< V >::derived_type mv(const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V > &v)
Returns the result of the projection of a vector v with a Matrix M.
Definition: math.hh:812
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
Definition: method.hh:58
forward declare
Definition: forchheimerslaw_fwd.hh:39
Forchheimer's law for box scheme.
Definition: box/forchheimerslaw.hh:75
static Scalar calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: box/forchheimerslaw.hh:168
static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem &problem, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Returns the harmonic mean of and .
Definition: box/forchheimerslaw.hh:182
static constexpr DiscretizationMethod discMethod
state the discretization method this implementation belongs to
Definition: box/forchheimerslaw.hh:94
ScalarType Scalar
state the scalar type of the law
Definition: box/forchheimerslaw.hh:90
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, int phaseIdx, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the advective flux of a phase across the given sub-control volume face using the Forchheimer ...
Definition: box/forchheimerslaw.hh:104
Darcy's law for control-volume finite element schemes.
Definition: flux/cvfe/darcyslaw.hh:48
Forchheimer's law For a detailed description see dumux/flow/forchheimerslaw.hh.
Definition: forchheimervelocity.hh:51
static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem &problem, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Returns the harmonic mean of and .
Definition: forchheimervelocity.hh:138
Dune::FieldVector< Scalar, dimWorld > DimWorldVector
Definition: forchheimervelocity.hh:62
static DimWorldVector velocity(const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const DimWorldMatrix sqrtK, const DimWorldVector darcyVelocity)
Definition: forchheimervelocity.hh:67
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimWorldMatrix
Definition: forchheimervelocity.hh:63
Declares all properties used in Dumux.
Specialization of Darcy's Law for control-volume finite element schemes.