3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
transmissibility1p.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_PNM_THROAT_TRANSMISSIBILITY_1P_HH
25#define DUMUX_PNM_THROAT_TRANSMISSIBILITY_1P_HH
26
29#include "emptycache.hh"
30
31namespace Dumux::PoreNetwork {
32
39template<class Scalar>
41{
42public:
43
45
46 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
47 static Scalar singlePhaseTransmissibility(const Problem& problem,
48 const Element& element,
49 const FVElementGeometry& fvGeometry,
50 const typename FVElementGeometry::SubControlVolumeFace& scvf,
51 const ElementVolumeVariables& elemVolVars,
52 const FluxVariablesCache& fluxVarsCache,
53 const int phaseIdx)
54 {
55 const auto shape = fluxVarsCache.throatCrossSectionShape();
56 const Scalar throatInscribedRadius = fluxVarsCache.throatInscribedRadius();
57 const Scalar throatLength = fluxVarsCache.throatLength();
58 const Scalar area = fluxVarsCache.throatCrossSectionalArea();
59 return singlePhaseTransmissibility(shape, throatInscribedRadius, throatLength, area);
60 }
61
64 const Scalar inscribedRadius,
65 const Scalar throatLength,
66 const Scalar area)
67 { return 1.0 / rHydThroat_(shape, inscribedRadius, throatLength, area); }
68
69protected:
70
71 static Scalar rHydThroat_(const Throat::Shape shape,
72 const Scalar radius,
73 const Scalar length,
74 const Scalar area)
75 {
76 switch(shape)
77 {
79 {
80 const Scalar sideLength = 2.0*radius;
81 return 28.4*length * 1.0/(sideLength*sideLength*sideLength*sideLength);
82 }
84 {
85 return 8.0/M_PI * length * 1.0/(radius*radius*radius*radius);
86 }
88 {
89 using std::sqrt;
90 static Scalar sqrt3 = sqrt(3.0);
91 const Scalar sideLength = 6.0/sqrt3 * radius;
92 return 320.0/sqrt3 * length * 1.0/(sideLength*sideLength*sideLength*sideLength);
93 }
95 {
96 // the distance between the two parallel plates
97 const Scalar width = 2*radius;
98 return 12.0/(width*width*width) * length;
99 }
101 {
102 // requires width >> height for good accuracy
103 Scalar height = 2.0*radius;
104 Scalar width = area/height;
105
106 using std::swap;
107 if (width < height)
108 swap(width, height);
109
110 return 12.0*length / (1.0 - 0.63*(height/width)) * 1.0/(height*height*height*width);
111 }
112 default: DUNE_THROW(Dune::InvalidStateException, "Throat geometry not supported");
113 }
114 }
115};
116
121template<class Scalar, bool considerPoreResistance = false, bool interpolateK = false>
123{
124 static_assert(!interpolateK, "Interpolation of k not implemented");
125public:
126
128
129 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
130 static Scalar singlePhaseTransmissibility(const Problem& problem,
131 const Element& element,
132 const FVElementGeometry& fvGeometry,
133 const typename FVElementGeometry::SubControlVolumeFace& scvf,
134 const ElementVolumeVariables& elemVolVars,
135 const FluxVariablesCache& fluxVarsCache,
136 const int phaseIdx)
137 {
138 const auto shapeFactor = fluxVarsCache.throatShapeFactor();
139 const Scalar area = fluxVarsCache.throatCrossSectionalArea();
140 const Scalar throatLength = fluxVarsCache.throatLength();
141 const Scalar throatTransmissibility = singlePhaseTransmissibility(shapeFactor, throatLength, area);
142
143 if constexpr (!considerPoreResistance)
144 return throatTransmissibility;
145 else
146 {
147 static const bool considerPoreResistanceOnRuntime = getParamFromGroup<bool>(problem.paramGroup(), "Transmissibility.ConsiderPoreResistance", true);
148 if (!considerPoreResistanceOnRuntime)
149 return throatTransmissibility;
150
151 const auto& scv0 = fvGeometry.scv(scvf.insideScvIdx());
152 const auto& scv1 = fvGeometry.scv(scvf.outsideScvIdx());
153
154 const auto& spatialParams = problem.spatialParams();
155 const auto elemSol = elementSolution(element, elemVolVars, fvGeometry);
156
157 // TODO maybe include this in fluxVarsCache if this is general enough
158 const Scalar poreLength0 = spatialParams.poreLength(element, scv0, elemSol);
159 const Scalar poreLength1 = spatialParams.poreLength(element, scv1, elemSol);
160
161 const auto poreShapeFactor0 = spatialParams.poreShapeFactor(element, scv0, elemSol);
162 const auto poreShapeFactor1 = spatialParams.poreShapeFactor(element, scv1, elemSol);
163
164 const auto poreCrossSectionalArea0 = spatialParams.poreCrossSectionalArea(element, scv0, elemSol);
165 const auto poreCrossSectionalArea1 = spatialParams.poreCrossSectionalArea(element, scv1, elemSol);
166
167 const Scalar poreTransmissibility0 = singlePhaseTransmissibility(poreShapeFactor0, poreLength0, poreCrossSectionalArea0);
168 const Scalar poreTransmissibility1 = singlePhaseTransmissibility(poreShapeFactor1, poreLength1, poreCrossSectionalArea1);
169
170 return 1 / (1.0/throatTransmissibility + 1.0/poreTransmissibility0 + 1.0/poreTransmissibility1);
171 }
172 }
173
175 static Scalar singlePhaseTransmissibility(const Scalar shapeFactor,
176 const Scalar length,
177 const Scalar area)
178 {
179 const Scalar k = k_(shapeFactor);
180 return k * area*area * shapeFactor / length;
181 }
182
183private:
184 static Scalar k_(const Scalar shapeFactor)
185 {
186 if (shapeFactor <= Throat::shapeFactorEquilateralTriangle<Scalar>())
187 return 0.6; // == 3/5
188 else if (shapeFactor <= Throat::shapeFactorSquare<Scalar>())
189 return 0.5623;
190 else // circle
191 return 0.5;
192 }
193
194 // TODO interpolation
195};
196
197
199template<class Scalar>
201{
202public:
203
205
206 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
207 static Scalar singlePhaseTransmissibility(const Problem& problem,
208 const Element& element,
209 const FVElementGeometry& fvGeometry,
210 const typename FVElementGeometry::SubControlVolumeFace& scvf,
211 const ElementVolumeVariables& elemVolVars,
212 const FluxVariablesCache& fluxVarsCache,
213 const int phaseIdx)
214 {
215 const Scalar throatInscribedRadius = fluxVarsCache.throatInscribedRadius();
216 const Scalar throatLength = fluxVarsCache.throatLength();
217 return singlePhaseTransmissibility(throatInscribedRadius, throatLength);
218 }
219
221 static Scalar singlePhaseTransmissibility(const Scalar radius,
222 const Scalar length)
223 {
224 const Scalar rEff= std::sqrt(4.0/M_PI)*radius ;
225 return M_PI/(8.0*length) *rEff*rEff*rEff*rEff ;
226 }
227};
228
229} // end namespace Dumux::Porenetwork
230
231#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
An empty cache for transmissibility laws using only standard quantities.
This file contains functions related to calculate pore-throat properties.
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethod::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
Definition: discretization/porenetwork/fvelementgeometry.hh:33
constexpr Shape shape(const Scalar shapeFactor) noexcept
Returns the shape for a given shape factor.
Definition: throatproperties.hh:176
Shape
Collection of different pore-throat shapes.
Definition: throatproperties.hh:38
Scalar shapeFactor(Shape shape, const Scalar inscribedRadius)
Returns the value of the shape factor for a given shape.
Definition: throatproperties.hh:161
Definition: emptycache.hh:30
Collection of single-phase flow throat transmissibilities based on Bruus, H. (2011)....
Definition: transmissibility1p.hh:41
static Scalar singlePhaseTransmissibility(const Throat::Shape shape, const Scalar inscribedRadius, const Scalar throatLength, const Scalar area)
Returns the conductivity of a throat when only one phase is present.
Definition: transmissibility1p.hh:63
static Scalar singlePhaseTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const FluxVariablesCache &fluxVarsCache, const int phaseIdx)
Definition: transmissibility1p.hh:47
static Scalar rHydThroat_(const Throat::Shape shape, const Scalar radius, const Scalar length, const Scalar area)
Definition: transmissibility1p.hh:71
Single-phase flow throat transmissibility based on Patzek & Silin (2001) https://doi....
Definition: transmissibility1p.hh:123
static Scalar singlePhaseTransmissibility(const Scalar shapeFactor, const Scalar length, const Scalar area)
Returns the conductivity of a throat when only one phase is present. See Patzek & Silin (2001)
Definition: transmissibility1p.hh:175
static Scalar singlePhaseTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const FluxVariablesCache &fluxVarsCache, const int phaseIdx)
Definition: transmissibility1p.hh:130
Used by Joeakar-Niasar.
Definition: transmissibility1p.hh:201
static Scalar singlePhaseTransmissibility(const Scalar radius, const Scalar length)
Returns the conductivity of a throat when only one phase is present.
Definition: transmissibility1p.hh:221
static Scalar singlePhaseTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const FluxVariablesCache &fluxVarsCache, const int phaseIdx)
Definition: transmissibility1p.hh:207