3.5-git
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 *****************************************************************************/
25#ifndef DUMUX_PNM_THROAT_TRANSMISSIBILITY_2P_HH
26#define DUMUX_PNM_THROAT_TRANSMISSIBILITY_2P_HH
27
29#include "emptycache.hh"
30
31namespace Dumux::PoreNetwork {
32
33namespace WettingLayerTransmissibility {
34
36{
46 template<class Scalar>
47 static Scalar beta(const Scalar alpha, const Scalar theta, const Scalar f = 0)
48 {
49 using std::sin;
50 using std::cos;
51 using std::tan;
52 const Scalar sinAlpha = sin(alpha);
53 const Scalar sinSum = sin(alpha + theta);
54 const Scalar cosSum = cos(alpha + theta);
55 const Scalar phi1 = cosSum*cosSum + cosSum*sinSum*tan(alpha);
56 const Scalar phi2 = 1 - theta/(M_PI/2 - alpha);
57 const Scalar phi3 = cosSum / cos(alpha);
58 const Scalar B = (M_PI/2 - alpha)*tan(alpha);
59
60 Scalar result = 12*sinAlpha*sinAlpha*(1-B)*(1-B)*(phi1 - B*phi2)*(phi3 + f*B*phi2)*(phi3 + f*B*phi2);
61 result /= (1-sinAlpha)*(1-sinAlpha)*B*B*(phi1 - B*phi2)*(phi1 - B*phi2)*(phi1 - B*phi2);
62 return result;
63 }
64};
65
66
67template<class Scalar, class CreviceResistanceFactor = CreviceResistanceFactorZhou>
69{
71 {
72 using NumCornerVector = Dune::ReservedVector<Scalar, 4>;
73 public:
74 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
75 void fill(const Problem& problem,
76 const Element& element,
77 const FVElementGeometry& fvGeometry,
78 const typename FVElementGeometry::SubControlVolumeFace& scvf,
79 const ElementVolumeVariables& elemVolVars,
80 const FluxVariablesCache& fluxVarsCache,
81 const int phaseIdx)
82 {
83 const auto& spatialParams = problem.spatialParams();
84 const auto& cornerHalfAngles = spatialParams.cornerHalfAngles(element);
85 const Scalar contactAngle = spatialParams.contactAngle(element, elemVolVars);
86
87 beta_.clear(); beta_.resize(cornerHalfAngles.size());
88 for (int i = 0; i< cornerHalfAngles.size(); ++i)
89 beta_[i] = CreviceResistanceFactor::beta(cornerHalfAngles[i], contactAngle);
90 }
91
92 Scalar creviceResistanceFactor(const int cornerIdx) const
93 { return beta_[cornerIdx]; }
94
95 private:
96 NumCornerVector beta_;
97 };
98
102 template<class Element, class FVElementGeometry, class FluxVariablesCache>
103 static Scalar wettingLayerTransmissibility(const Element& element,
104 const FVElementGeometry& fvGeometry,
105 const typename FVElementGeometry::SubControlVolumeFace& scvf,
106 const FluxVariablesCache& fluxVarsCache)
107 {
108 const Scalar throatLength = fluxVarsCache.throatLength();
109 const Scalar rC = fluxVarsCache.curvatureRadius();
110
111 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
112 const auto shape = fvGeometry.gridGeometry().throatCrossSectionShape(eIdx);
113 const auto numCorners = Throat::numCorners(shape);
114
115 // treat the wetting film layer in each corner of the throat individually (might have different corner half-angle and beta)
116 Scalar result = 0.0;
117 for (int i = 0; i < numCorners; ++i)
118 result += fluxVarsCache.wettingLayerCrossSectionalArea(i) * rC*rC / (throatLength*fluxVarsCache.wettingLayerFlowVariables().creviceResistanceFactor(i));
119
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 return result;
147 }
148};
149
151template<class Scalar, class SinglePhaseTransmissibilityLaw>
153{
155
156 template<class Element, class FVElementGeometry, class FluxVariablesCache>
157 static Scalar nonWettingPhaseTransmissibility(const Element& element,
158 const FVElementGeometry& fvGeometry,
159 const typename FVElementGeometry::SubControlVolumeFace& scvf,
160 const FluxVariablesCache& fluxVarsCache)
161 {
162 // Tora et al. (2012), also does not fully recover single-phase value, but is closer
163 using std::sqrt;
164 const auto nPhaseIdx = fluxVarsCache.nPhaseIdx();
165 const Scalar aNw = fluxVarsCache.throatCrossSectionalArea(nPhaseIdx);
166 const Scalar aTot = fluxVarsCache.throatCrossSectionalArea();
167
168 const Scalar result = SinglePhaseTransmissibilityLaw::singlePhaseTransmissibility(element, fvGeometry, scvf, fluxVarsCache)
169 * aNw / aTot;
170
171 return result;
172 }
173};
174
175template<class Scalar>
177{
179
185 template<class Element, class FVElementGeometry, class FluxVariablesCache>
186 static Scalar nonWettingPhaseTransmissibility(const Element& element,
187 const FVElementGeometry& fvGeometry,
188 const typename FVElementGeometry::SubControlVolumeFace& scvf,
189 const FluxVariablesCache& fluxVarsCache)
190 {
191 // Tora et al. (2012), quite close for single-phase value of square
192 using std::sqrt;
193 const Scalar throatLength = fluxVarsCache.throatLength();
194 const auto nPhaseIdx = fluxVarsCache.nPhaseIdx();
195 const Scalar aNw = fluxVarsCache.throatCrossSectionalArea(nPhaseIdx);
196 const Scalar rEff = 0.5*(sqrt(aNw / M_PI) + fluxVarsCache.throatInscribedRadius());
197 const Scalar result = rEff*rEff*aNw / (8.0*throatLength);
198 return result;
199 }
200};
201} // end namespace NonWettingPhaseTransmissibility
202} // end namespace Dumux::PoreNetwork
203
204#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:34
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:31
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:47
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:103
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:75
Scalar creviceResistanceFactor(const int cornerIdx) const
Definition: transmissibility2p.hh:92
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:153
static Scalar nonWettingPhaseTransmissibility(const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const FluxVariablesCache &fluxVarsCache)
Definition: transmissibility2p.hh:157
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:186