3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
transmissibility2p.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_2P_HH
25#define DUMUX_PNM_THROAT_TRANSMISSIBILITY_2P_HH
26
28#include "emptycache.hh"
29
30namespace Dumux::PoreNetwork {
31
32namespace WettingLayerTransmissibility {
33
35{
45 template<class Scalar>
46 static Scalar beta(const Scalar alpha, const Scalar theta, const Scalar f = 0)
47 {
48 using std::sin;
49 using std::cos;
50 using std::tan;
51 const Scalar sinAlpha = sin(alpha);
52 const Scalar sinSum = sin(alpha + theta);
53 const Scalar cosSum = cos(alpha + theta);
54 const Scalar phi1 = cosSum*cosSum + cosSum*sinSum*tan(alpha);
55 const Scalar phi2 = 1 - theta/(M_PI/2 - alpha);
56 const Scalar phi3 = cosSum / cos(alpha);
57 const Scalar B = (M_PI/2 - alpha)*tan(alpha);
58
59 Scalar result = 12*sinAlpha*sinAlpha*(1-B)*(1-B)*(phi1 - B*phi2)*(phi3 + f*B*phi2)*(phi3 + f*B*phi2);
60 result /= (1-sinAlpha)*(1-sinAlpha)*B*B*(phi1 - B*phi2)*(phi1 - B*phi2)*(phi1 - B*phi2);
61 return result;
62 }
63};
64
65
66template<class Scalar, class CreviceResistanceFactor = CreviceResistanceFactorZhou>
68{
70 {
71 using NumCornerVector = Dune::ReservedVector<Scalar, 4>;
72 public:
73 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
74 void fill(const Problem& problem,
75 const Element& element,
76 const FVElementGeometry& fvGeometry,
77 const typename FVElementGeometry::SubControlVolumeFace& scvf,
78 const ElementVolumeVariables& elemVolVars,
79 const FluxVariablesCache& fluxVarsCache,
80 const int phaseIdx)
81 {
82 const auto& spatialParams = problem.spatialParams();
83 const auto& cornerHalfAngles = spatialParams.cornerHalfAngles(element);
84 const Scalar contactAngle = spatialParams.contactAngle(element, elemVolVars);
85
86 beta_.clear(); beta_.resize(cornerHalfAngles.size());
87 for (int i = 0; i< cornerHalfAngles.size(); ++i)
88 beta_[i] = CreviceResistanceFactor::beta(cornerHalfAngles[i], contactAngle);
89 }
90
91 Scalar creviceResistanceFactor(const int cornerIdx) const
92 { return beta_[cornerIdx]; }
93
94 private:
95 NumCornerVector beta_;
96 };
97
101 template<class Element, class FVElementGeometry, class FluxVariablesCache>
102 static Scalar wettingLayerTransmissibility(const Element& element,
103 const FVElementGeometry& fvGeometry,
104 const typename FVElementGeometry::SubControlVolumeFace& scvf,
105 const FluxVariablesCache& fluxVarsCache)
106 {
107 const Scalar throatLength = fluxVarsCache.throatLength();
108 const Scalar rC = fluxVarsCache.curvatureRadius();
109
110 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
111 const auto shape = fvGeometry.gridGeometry().throatCrossSectionShape(eIdx);
112 const auto numCorners = Throat::numCorners(shape);
113
114 // treat the wetting film layer in each corner of the throat individually (might have different corner half-angle and beta)
115 Scalar result = 0.0;
116 for (int i = 0; i < numCorners; ++i)
117 result += fluxVarsCache.wettingLayerCrossSectionalArea(i) * rC*rC / (throatLength*fluxVarsCache.wettingLayerFlowVariables().creviceResistanceFactor(i));
118
119 assert(std::isnormal(result));
120 return result;
121 }
122};
123} // end namespace WettingLayerTransmissibility
124
125namespace NonWettingPhaseTransmissibility {
126
128template<class Scalar>
130{
132
133 template<class Element, class FVElementGeometry, class FluxVariablesCache>
134 static Scalar nonWettingPhaseTransmissibility(const Element& element,
135 const FVElementGeometry& fvGeometry,
136 const typename FVElementGeometry::SubControlVolumeFace& scvf,
137 const FluxVariablesCache& fluxVarsCache)
138 {
139 // Mogensen et al. (1999), does not really revover the single phase value
140 using std::sqrt;
141 const Scalar throatLength = fluxVarsCache.throatLength();
142 const auto nPhaseIdx = fluxVarsCache.nPhaseIdx();
143 const Scalar aNw = fluxVarsCache.throatCrossSectionalArea(nPhaseIdx);
144 const Scalar rEff = 0.5*(sqrt(aNw / M_PI) + fluxVarsCache.throatInscribedRadius());
145 const Scalar result = M_PI/(8*throatLength) * rEff*rEff*rEff*rEff;
146 assert(std::isnormal(result));
147 return result;
148 }
149};
150
152template<class Scalar, class SinglePhaseTransmissibilityLaw>
154{
156
157 template<class Element, class FVElementGeometry, class FluxVariablesCache>
158 static Scalar nonWettingPhaseTransmissibility(const Element& element,
159 const FVElementGeometry& fvGeometry,
160 const typename FVElementGeometry::SubControlVolumeFace& scvf,
161 const FluxVariablesCache& fluxVarsCache)
162 {
163 // Tora et al. (2012), also does not fully recover single-phase value, but is closer
164 using std::sqrt;
165 const auto nPhaseIdx = fluxVarsCache.nPhaseIdx();
166 const Scalar aNw = fluxVarsCache.throatCrossSectionalArea(nPhaseIdx);
167 const Scalar aTot = fluxVarsCache.throatCrossSectionalArea();
168
169 const Scalar result = SinglePhaseTransmissibilityLaw::singlePhaseTransmissibility(element, fvGeometry, scvf, fluxVarsCache)
170 * aNw / aTot;
171
172 assert(std::isnormal(result));
173 return result;
174 }
175};
176
177template<class Scalar>
179{
181
187 template<class Element, class FVElementGeometry, class FluxVariablesCache>
188 static Scalar nonWettingPhaseTransmissibility(const Element& element,
189 const FVElementGeometry& fvGeometry,
190 const typename FVElementGeometry::SubControlVolumeFace& scvf,
191 const FluxVariablesCache& fluxVarsCache)
192 {
193 // Tora et al. (2012), quite close for single-phase value of square
194 using std::sqrt;
195 const Scalar throatLength = fluxVarsCache.throatLength();
196 const auto nPhaseIdx = fluxVarsCache.nPhaseIdx();
197 const Scalar aNw = fluxVarsCache.throatCrossSectionalArea(nPhaseIdx);
198 const Scalar rEff = 0.5*(sqrt(aNw / M_PI) + fluxVarsCache.throatInscribedRadius());
199 const Scalar result = rEff*rEff*aNw / (8.0*throatLength);
200 assert(std::isnormal(result));
201 return result;
202 }
203};
204} // end namespace NonWettingPhaseTransmissibility
205} // end namespace Dumux::PoreNetwork
206
207#endif
An empty cache for transmissibility laws using only standard quantities.
This file contains functions related to calculate pore-throat properties.
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
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:215
Dune::ReservedVector< Scalar, 4 > cornerHalfAngles(Shape shape)
Returns the corner half angle.
Definition: throatproperties.hh:97
Definition: emptycache.hh:30
static Scalar beta(const Scalar alpha, const Scalar theta, const Scalar f=0)
Returns the crevice resistance factor used for calculating the w-phase conductance in an invaded pore...
Definition: transmissibility2p.hh:46
static Scalar wettingLayerTransmissibility(const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const FluxVariablesCache &fluxVarsCache)
Returns the integral conductivity of all wetting layers occupying the corners of a throat.
Definition: transmissibility2p.hh:102
void fill(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const FluxVariablesCache &fluxVarsCache, const int phaseIdx)
Definition: transmissibility2p.hh:74
Scalar creviceResistanceFactor(const int cornerIdx) const
Definition: transmissibility2p.hh:91
TODO: evalute and maybe remove.
Definition: transmissibility2p.hh:130
static Scalar nonWettingPhaseTransmissibility(const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const FluxVariablesCache &fluxVarsCache)
Definition: transmissibility2p.hh:134
TODO: evalute and maybe remove.
Definition: transmissibility2p.hh:154
static Scalar nonWettingPhaseTransmissibility(const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const FluxVariablesCache &fluxVarsCache)
Definition: transmissibility2p.hh:158
static Scalar nonWettingPhaseTransmissibility(const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const FluxVariablesCache &fluxVarsCache)
Returns the conductivity of a throat.
Definition: transmissibility2p.hh:188