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