Line |
Branch |
Exec |
Source |
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 |
|
|
// |
7 |
|
|
/*! |
8 |
|
|
* \file |
9 |
|
|
* \ingroup CCTpfaFlux |
10 |
|
|
* \brief Forchheimers's law for cell-centered finite volume schemes with two-point flux approximation |
11 |
|
|
*/ |
12 |
|
|
#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FORCHHEIMERS_LAW_HH |
13 |
|
|
#define DUMUX_DISCRETIZATION_CC_TPFA_FORCHHEIMERS_LAW_HH |
14 |
|
|
|
15 |
|
|
#include <dune/common/fvector.hh> |
16 |
|
|
#include <dune/common/fmatrix.hh> |
17 |
|
|
|
18 |
|
|
#include <dumux/common/math.hh> |
19 |
|
|
#include <dumux/common/parameters.hh> |
20 |
|
|
#include <dumux/common/properties.hh> |
21 |
|
|
#include <dumux/common/typetraits/typetraits.hh> |
22 |
|
|
|
23 |
|
|
#include <dumux/discretization/method.hh> |
24 |
|
|
#include <dumux/discretization/extrusion.hh> |
25 |
|
|
#include <dumux/flux/cctpfa/darcyslaw.hh> |
26 |
|
|
|
27 |
|
|
namespace Dumux { |
28 |
|
|
|
29 |
|
|
// forward declarations |
30 |
|
|
template<class TypeTag, class ForchheimerVelocity, class DiscretizationMethod> |
31 |
|
|
class ForchheimersLawImplementation; |
32 |
|
|
|
33 |
|
|
/*! |
34 |
|
|
* \ingroup CCTpfaFlux |
35 |
|
|
* \brief Forchheimer's law for cell-centered finite volume schemes with two-point flux approximation |
36 |
|
|
* \note Forchheimer's law is specialized for network and surface grids (i.e. if grid dim < dimWorld) |
37 |
|
|
* \tparam Scalar the scalar type for scalar physical quantities |
38 |
|
|
* \tparam GridGeometry the grid geometry |
39 |
|
|
* \tparam ForchheimerVelocity class for the calculation of the Forchheimer velocity |
40 |
|
|
* \tparam isNetwork whether we are computing on a network grid embedded in a higher world dimension |
41 |
|
|
*/ |
42 |
|
|
template<class Scalar, class GridGeometry, class ForchheimerVelocity, bool isNetwork> |
43 |
|
|
class CCTpfaForchheimersLaw; |
44 |
|
|
|
45 |
|
|
/*! |
46 |
|
|
* \ingroup CCTpfaFlux |
47 |
|
|
* \brief Forchheimer's law for cell-centered finite volume schemes with two-point flux approximation |
48 |
|
|
* \note Forchheimer's law is specialized for network and surface grids (i.e. if grid dim < dimWorld) |
49 |
|
|
*/ |
50 |
|
|
template <class TypeTag, class ForchheimerVelocity> |
51 |
|
|
class ForchheimersLawImplementation<TypeTag, ForchheimerVelocity, DiscretizationMethods::CCTpfa> |
52 |
|
|
: public CCTpfaForchheimersLaw<GetPropType<TypeTag, Properties::Scalar>, |
53 |
|
|
GetPropType<TypeTag, Properties::GridGeometry>, |
54 |
|
|
ForchheimerVelocity, |
55 |
|
|
(GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimension < GetPropType<TypeTag, Properties::GridGeometry>::GridView::dimensionworld)> |
56 |
|
|
{}; |
57 |
|
|
|
58 |
|
|
/*! |
59 |
|
|
* \ingroup CCTpfaFlux |
60 |
|
|
* \brief Class that fills the cache corresponding to tpfa Forchheimer's Law |
61 |
|
|
*/ |
62 |
|
|
template<class GridGeometry> |
63 |
|
|
class TpfaForchheimersLawCacheFiller |
64 |
|
|
{ |
65 |
|
|
using FVElementGeometry = typename GridGeometry::LocalView; |
66 |
|
|
using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace; |
67 |
|
|
using Element = typename GridGeometry::GridView::template Codim<0>::Entity; |
68 |
|
|
|
69 |
|
|
public: |
70 |
|
|
//! Function to fill a TpfaForchheimersLawCache of a given scvf |
71 |
|
|
//! This interface has to be met by any advection-related cache filler class |
72 |
|
|
//! TODO: Probably get cache type out of the filler |
73 |
|
|
template<class FluxVariablesCache, class Problem, class ElementVolumeVariables, class FluxVariablesCacheFiller> |
74 |
|
✗ |
static void fill(FluxVariablesCache& scvfFluxVarsCache, |
75 |
|
|
const Problem& problem, |
76 |
|
|
const Element& element, |
77 |
|
|
const FVElementGeometry& fvGeometry, |
78 |
|
|
const ElementVolumeVariables& elemVolVars, |
79 |
|
|
const SubControlVolumeFace& scvf, |
80 |
|
|
const FluxVariablesCacheFiller& fluxVarsCacheFiller) |
81 |
|
|
{ |
82 |
|
2720 |
scvfFluxVarsCache.updateAdvection(problem, element, fvGeometry, elemVolVars, scvf); |
83 |
|
✗ |
} |
84 |
|
|
}; |
85 |
|
|
|
86 |
|
|
/*! |
87 |
|
|
* \ingroup CCTpfaFlux |
88 |
|
|
* \brief The cache corresponding to tpfa Forchheimer's Law |
89 |
|
|
*/ |
90 |
|
|
template<class AdvectionType, class GridGeometry> |
91 |
|
2992 |
class TpfaForchheimersLawCache |
92 |
|
|
{ |
93 |
|
|
using Scalar = typename AdvectionType::Scalar; |
94 |
|
|
using FVElementGeometry = typename GridGeometry::LocalView; |
95 |
|
|
using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace; |
96 |
|
|
using Element = typename GridGeometry::GridView::template Codim<0>::Entity; |
97 |
|
|
static constexpr int dimWorld = GridGeometry::GridView::dimensionworld; |
98 |
|
|
using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>; |
99 |
|
|
|
100 |
|
|
public: |
101 |
|
|
using Filler = TpfaForchheimersLawCacheFiller<GridGeometry>; |
102 |
|
|
|
103 |
|
|
template<class Problem, class ElementVolumeVariables> |
104 |
|
2720 |
void updateAdvection(const Problem& problem, |
105 |
|
|
const Element& element, |
106 |
|
|
const FVElementGeometry& fvGeometry, |
107 |
|
|
const ElementVolumeVariables& elemVolVars, |
108 |
|
|
const SubControlVolumeFace &scvf) |
109 |
|
|
{ |
110 |
|
5440 |
tij_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf); |
111 |
|
5440 |
harmonicMeanSqrtK_ = AdvectionType::calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf); |
112 |
|
2720 |
} |
113 |
|
|
|
114 |
|
|
const Scalar& advectionTij() const |
115 |
|
✗ |
{ return tij_; } |
116 |
|
|
|
117 |
|
|
const DimWorldMatrix& harmonicMeanSqrtPermeability() const |
118 |
|
2208 |
{ return harmonicMeanSqrtK_; } |
119 |
|
|
|
120 |
|
|
private: |
121 |
|
|
Scalar tij_; |
122 |
|
|
DimWorldMatrix harmonicMeanSqrtK_; |
123 |
|
|
}; |
124 |
|
|
|
125 |
|
|
/*! |
126 |
|
|
* \ingroup CCTpfaFlux |
127 |
|
|
* \brief Specialization of the CCTpfaForchheimersLaw grids where dim=dimWorld |
128 |
|
|
*/ |
129 |
|
|
template<class ScalarType, class GridGeometry, class ForchheimerVelocity> |
130 |
|
|
class CCTpfaForchheimersLaw<ScalarType, GridGeometry, ForchheimerVelocity, /*isNetwork*/ false> |
131 |
|
|
{ |
132 |
|
|
using ThisType = CCTpfaForchheimersLaw<ScalarType, GridGeometry, ForchheimerVelocity, /*isNetwork*/ false>; |
133 |
|
|
using FVElementGeometry = typename GridGeometry::LocalView; |
134 |
|
|
using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace; |
135 |
|
|
using Extrusion = Extrusion_t<GridGeometry>; |
136 |
|
|
using GridView = typename GridGeometry::GridView; |
137 |
|
|
using Element = typename GridView::template Codim<0>::Entity; |
138 |
|
|
|
139 |
|
|
using DimWorldVector = typename ForchheimerVelocity::DimWorldVector; |
140 |
|
|
using DimWorldMatrix = typename ForchheimerVelocity::DimWorldMatrix; |
141 |
|
|
|
142 |
|
|
using DarcysLaw = CCTpfaDarcysLaw<ScalarType, GridGeometry, /*isNetwork*/ false>; |
143 |
|
|
|
144 |
|
|
public: |
145 |
|
|
//! state the scalar type of the law |
146 |
|
|
using Scalar = ScalarType; |
147 |
|
|
|
148 |
|
|
using DiscretizationMethod = DiscretizationMethods::CCTpfa; |
149 |
|
|
//! state the discretization method this implementation belongs to |
150 |
|
|
static constexpr DiscretizationMethod discMethod{}; |
151 |
|
|
|
152 |
|
|
//! state the type for the corresponding cache |
153 |
|
|
using Cache = TpfaForchheimersLawCache<ThisType, GridGeometry>; |
154 |
|
|
|
155 |
|
|
/*! \brief Compute the advective flux of a phase across |
156 |
|
|
* the given sub-control volume face using the Forchheimer equation. |
157 |
|
|
* |
158 |
|
|
* The flux is given in N*m, and can be converted |
159 |
|
|
* into a volume flux (m^3/s) or mass flux (kg/s) by applying an upwind scheme |
160 |
|
|
* for the mobility or the product of density and mobility, respectively. |
161 |
|
|
*/ |
162 |
|
|
template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache> |
163 |
|
2208 |
static Scalar flux(const Problem& problem, |
164 |
|
|
const Element& element, |
165 |
|
|
const FVElementGeometry& fvGeometry, |
166 |
|
|
const ElementVolumeVariables& elemVolVars, |
167 |
|
|
const SubControlVolumeFace& scvf, |
168 |
|
|
int phaseIdx, |
169 |
|
|
const ElementFluxVarsCache& elemFluxVarsCache) |
170 |
|
|
{ |
171 |
|
|
// Get the volume flux based on Darcy's law. The value returned by this method needs to be multiplied with the |
172 |
|
|
// mobility (upwinding). |
173 |
|
2208 |
Scalar darcyFlux = DarcysLaw::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache); |
174 |
|
8832 |
auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.mobility(phaseIdx); }; |
175 |
|
2208 |
DimWorldVector darcyVelocity = scvf.unitOuterNormal(); |
176 |
|
2208 |
darcyVelocity *= ForchheimerVelocity::UpwindScheme::apply(elemVolVars, scvf, upwindTerm, darcyFlux, phaseIdx); |
177 |
|
4416 |
darcyVelocity /= Extrusion::area(fvGeometry, scvf); |
178 |
|
|
|
179 |
|
4416 |
const auto velocity = ForchheimerVelocity::velocity(fvGeometry, |
180 |
|
|
elemVolVars, |
181 |
|
|
scvf, |
182 |
|
|
phaseIdx, |
183 |
|
2208 |
elemFluxVarsCache[scvf].harmonicMeanSqrtPermeability(), |
184 |
|
|
darcyVelocity); |
185 |
|
|
|
186 |
|
2208 |
Scalar flux = velocity * scvf.unitOuterNormal(); |
187 |
|
4416 |
flux *= Extrusion::area(fvGeometry, scvf); |
188 |
|
2208 |
return flux; |
189 |
|
|
} |
190 |
|
|
|
191 |
|
|
//! The flux variables cache has to be bound to an element prior to flux calculations |
192 |
|
|
//! During the binding, the transmissibility will be computed and stored using the method below. |
193 |
|
|
template<class Problem, class ElementVolumeVariables> |
194 |
|
|
static Scalar calculateTransmissibility(const Problem& problem, |
195 |
|
|
const Element& element, |
196 |
|
|
const FVElementGeometry& fvGeometry, |
197 |
|
|
const ElementVolumeVariables& elemVolVars, |
198 |
|
|
const SubControlVolumeFace& scvf) |
199 |
|
|
{ |
200 |
|
2720 |
return DarcysLaw::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf); |
201 |
|
|
} |
202 |
|
|
|
203 |
|
|
/*! \brief Returns the harmonic mean of \f$\sqrt{K_0}\f$ and \f$\sqrt{K_1}\f$. |
204 |
|
|
* |
205 |
|
|
* This is a specialization for scalar-valued permeabilities which returns a tensor with identical diagonal entries. |
206 |
|
|
*/ |
207 |
|
|
template<class Problem, class ElementVolumeVariables> |
208 |
|
|
static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem& problem, |
209 |
|
|
const ElementVolumeVariables& elemVolVars, |
210 |
|
|
const SubControlVolumeFace& scvf) |
211 |
|
|
{ |
212 |
|
2720 |
return ForchheimerVelocity::calculateHarmonicMeanSqrtPermeability(problem, elemVolVars, scvf); |
213 |
|
|
} |
214 |
|
|
}; |
215 |
|
|
|
216 |
|
|
/*! |
217 |
|
|
* \ingroup CCTpfaFlux |
218 |
|
|
* \brief Specialization of the CCTpfaForchheimersLaw grids where dim<dimWorld |
219 |
|
|
*/ |
220 |
|
|
template<class ScalarType, class GridGeometry, class ForchheimerVelocity> |
221 |
|
|
class CCTpfaForchheimersLaw<ScalarType, GridGeometry, ForchheimerVelocity, /*isNetwork*/ true> |
222 |
|
|
{ |
223 |
|
|
static_assert(AlwaysFalse<ScalarType>::value, "Forchheimer not implemented for network grids"); |
224 |
|
|
}; |
225 |
|
|
|
226 |
|
|
} // end namespace Dumux |
227 |
|
|
|
228 |
|
|
#endif |
229 |
|
|
|