3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
cctpfa/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_CC_TPFA_FORCHHEIMERS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FORCHHEIMERS_LAW_HH
26
27#include <dune/common/fvector.hh>
28#include <dune/common/fmatrix.hh>
29
30#include <dumux/common/math.hh>
34
38
39namespace Dumux {
40
41// forward declarations
42template<class TypeTag, class ForchheimerVelocity, class DiscretizationMethod>
43class ForchheimersLawImplementation;
44
54template<class Scalar, class GridGeometry, class ForchheimerVelocity, bool isNetwork>
56
62template <class TypeTag, class ForchheimerVelocity>
63class ForchheimersLawImplementation<TypeTag, ForchheimerVelocity, DiscretizationMethods::CCTpfa>
64: public CCTpfaForchheimersLaw<GetPropType<TypeTag, Properties::Scalar>,
65 GetPropType<TypeTag, Properties::GridGeometry>,
66 ForchheimerVelocity,
67 (GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension < GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld)>
68{};
69
74template<class GridGeometry>
75class TpfaForchheimersLawCacheFiller
76{
77 using FVElementGeometry = typename GridGeometry::LocalView;
78 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
79 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
80
81public:
85 template<class FluxVariablesCache, class Problem, class ElementVolumeVariables, class FluxVariablesCacheFiller>
86 static void fill(FluxVariablesCache& scvfFluxVarsCache,
87 const Problem& problem,
88 const Element& element,
89 const FVElementGeometry& fvGeometry,
90 const ElementVolumeVariables& elemVolVars,
91 const SubControlVolumeFace& scvf,
92 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
93 {
94 scvfFluxVarsCache.updateAdvection(problem, element, fvGeometry, elemVolVars, scvf);
95 }
96};
97
102template<class AdvectionType, class GridGeometry>
103class TpfaForchheimersLawCache
104{
105 using Scalar = typename AdvectionType::Scalar;
106 using FVElementGeometry = typename GridGeometry::LocalView;
107 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
108 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
109 static constexpr int dimWorld = GridGeometry::GridView::dimensionworld;
110 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
111
112public:
113 using Filler = TpfaForchheimersLawCacheFiller<GridGeometry>;
114
115 template<class Problem, class ElementVolumeVariables>
116 void updateAdvection(const Problem& problem,
117 const Element& element,
118 const FVElementGeometry& fvGeometry,
119 const ElementVolumeVariables& elemVolVars,
120 const SubControlVolumeFace &scvf)
121 {
122 tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
123 harmonicMeanSqrtK_ = AdvectionType::calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf);
124 }
125
126 const Scalar& advectionTij() const
127 { return tij_; }
128
129 const DimWorldMatrix& harmonicMeanSqrtPermeability() const
130 { return harmonicMeanSqrtK_; }
131
132private:
133 Scalar tij_;
134 DimWorldMatrix harmonicMeanSqrtK_;
135};
136
141template<class ScalarType, class GridGeometry, class ForchheimerVelocity>
142class CCTpfaForchheimersLaw<ScalarType, GridGeometry, ForchheimerVelocity, /*isNetwork*/ false>
143{
144 using ThisType = CCTpfaForchheimersLaw<ScalarType, GridGeometry, ForchheimerVelocity, /*isNetwork*/ false>;
145 using FVElementGeometry = typename GridGeometry::LocalView;
146 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
147 using Extrusion = Extrusion_t<GridGeometry>;
148 using GridView = typename GridGeometry::GridView;
149 using Element = typename GridView::template Codim<0>::Entity;
150
151 using DimWorldVector = typename ForchheimerVelocity::DimWorldVector;
152 using DimWorldMatrix = typename ForchheimerVelocity::DimWorldMatrix;
153
154 using DarcysLaw = CCTpfaDarcysLaw<ScalarType, GridGeometry, /*isNetwork*/ false>;
155
156 public:
158 using Scalar = ScalarType;
159
160 using DiscretizationMethod = DiscretizationMethods::CCTpfa;
162 static constexpr DiscretizationMethod discMethod{};
163
165 using Cache = TpfaForchheimersLawCache<ThisType, GridGeometry>;
166
174 template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache>
175 static Scalar flux(const Problem& problem,
176 const Element& element,
177 const FVElementGeometry& fvGeometry,
178 const ElementVolumeVariables& elemVolVars,
179 const SubControlVolumeFace& scvf,
180 int phaseIdx,
181 const ElementFluxVarsCache& elemFluxVarsCache)
182 {
183 // Get the volume flux based on Darcy's law. The value returned by this method needs to be multiplied with the
184 // mobility (upwinding).
185 Scalar darcyFlux = DarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
186 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.mobility(phaseIdx); };
187 DimWorldVector darcyVelocity = scvf.unitOuterNormal();
188 darcyVelocity *= ForchheimerVelocity::UpwindScheme::apply(elemVolVars, scvf, upwindTerm, darcyFlux, phaseIdx);
189 darcyVelocity /= Extrusion::area(scvf);
190
191 const auto velocity = ForchheimerVelocity::velocity(fvGeometry,
192 elemVolVars,
193 scvf,
194 phaseIdx,
195 elemFluxVarsCache[scvf].harmonicMeanSqrtPermeability(),
196 darcyVelocity);
197
198 Scalar flux = velocity * scvf.unitOuterNormal();
199 flux *= Extrusion::area(scvf);
200 return flux;
201 }
202
205 template<class Problem, class ElementVolumeVariables>
206 static Scalar calculateTransmissibility(const Problem& problem,
207 const Element& element,
208 const FVElementGeometry& fvGeometry,
209 const ElementVolumeVariables& elemVolVars,
210 const SubControlVolumeFace& scvf)
211 {
212 return DarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
213 }
214
219 template<class Problem, class ElementVolumeVariables>
220 static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem& problem,
221 const ElementVolumeVariables& elemVolVars,
222 const SubControlVolumeFace& scvf)
223 {
224 return ForchheimerVelocity::calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf);
225 }
226};
227
232template<class ScalarType, class GridGeometry, class ForchheimerVelocity>
233class CCTpfaForchheimersLaw<ScalarType, GridGeometry, ForchheimerVelocity, /*isNetwork*/ true>
234{
235 static_assert(AlwaysFalse<ScalarType>::value, "Forchheimer not implemented for network grids");
236};
237
238} // end namespace Dumux
239
240#endif
Type traits.
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Definition: adapt.hh:29
forward declare
Definition: forchheimerslaw_fwd.hh:39
Forchheimer's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: cctpfa/forchheimerslaw.hh:55
Forchheimer's law This file contains the calculation of the Forchheimer velocity for a given Darcy ve...
Definition: forchheimervelocity.hh:49
Declares all properties used in Dumux.
Darcy's law for cell-centered finite volume schemes with two-point flux approximation.