version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_PNM_THROAT_TRANSMISSIBILITY_2P_HH
14#define DUMUX_PNM_THROAT_TRANSMISSIBILITY_2P_HH
15
17#include "emptycache.hh"
18
19namespace Dumux::PoreNetwork {
20
21namespace WettingLayerTransmissibility {
22
24{
34 template<class Scalar>
35 static Scalar beta(const Scalar alpha, const Scalar theta, const Scalar f = 0)
36 {
37 using std::sin;
38 using std::cos;
39 using std::tan;
40 const Scalar sinAlpha = sin(alpha);
41 const Scalar sinSum = sin(alpha + theta);
42 const Scalar cosSum = cos(alpha + theta);
43 const Scalar phi1 = cosSum*cosSum + cosSum*sinSum*tan(alpha);
44 const Scalar phi2 = 1 - theta/(M_PI/2 - alpha);
45 const Scalar phi3 = cosSum / cos(alpha);
46 const Scalar B = (M_PI/2 - alpha)*tan(alpha);
47
48 Scalar result = 12*sinAlpha*sinAlpha*(1-B)*(1-B)*(phi1 - B*phi2)*(phi3 + f*B*phi2)*(phi3 + f*B*phi2);
49 result /= (1-sinAlpha)*(1-sinAlpha)*B*B*(phi1 - B*phi2)*(phi1 - B*phi2)*(phi1 - B*phi2);
50 return result;
51 }
52};
53
54
55template<class Scalar, class CreviceResistanceFactor = CreviceResistanceFactorZhou>
57{
59 {
60 using NumCornerVector = Dune::ReservedVector<Scalar, 4>;
61 public:
62 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables, class FluxVariablesCache>
63 void fill(const Problem& problem,
64 const Element& element,
65 const FVElementGeometry& fvGeometry,
66 const typename FVElementGeometry::SubControlVolumeFace& scvf,
67 const ElementVolumeVariables& elemVolVars,
68 const FluxVariablesCache& fluxVarsCache,
69 const int phaseIdx)
70 {
71 const auto& spatialParams = problem.spatialParams();
72 const auto& cornerHalfAngles = spatialParams.cornerHalfAngles(element);
73 const Scalar contactAngle = spatialParams.contactAngle(element, elemVolVars);
74
75 beta_.clear(); beta_.resize(cornerHalfAngles.size());
76 for (int i = 0; i< cornerHalfAngles.size(); ++i)
77 beta_[i] = CreviceResistanceFactor::beta(cornerHalfAngles[i], contactAngle);
78 }
79
80 Scalar creviceResistanceFactor(const int cornerIdx) const
81 { return beta_[cornerIdx]; }
82
83 private:
84 NumCornerVector beta_;
85 };
86
90 template<class Element, class FVElementGeometry, class FluxVariablesCache>
91 static Scalar wettingLayerTransmissibility(const Element& element,
92 const FVElementGeometry& fvGeometry,
93 const typename FVElementGeometry::SubControlVolumeFace& scvf,
94 const FluxVariablesCache& fluxVarsCache)
95 {
96 const Scalar throatLength = fluxVarsCache.throatLength();
97 const Scalar rC = fluxVarsCache.curvatureRadius();
98
99 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
100 const auto shape = fvGeometry.gridGeometry().throatCrossSectionShape(eIdx);
101 const auto numCorners = Throat::numCorners(shape);
102
103 // treat the wetting film layer in each corner of the throat individually (might have different corner half-angle and beta)
104 Scalar result = 0.0;
105 for (int i = 0; i < numCorners; ++i)
106 result += fluxVarsCache.wettingLayerCrossSectionalArea(i) * rC*rC / (throatLength*fluxVarsCache.wettingLayerFlowVariables().creviceResistanceFactor(i));
107
108 return result;
109 }
110};
111} // end namespace WettingLayerTransmissibility
112
113namespace NonWettingPhaseTransmissibility {
114
116template<class Scalar>
118{
120
121 template<class Element, class FVElementGeometry, class FluxVariablesCache>
122 static Scalar nonWettingPhaseTransmissibility(const Element& element,
123 const FVElementGeometry& fvGeometry,
124 const typename FVElementGeometry::SubControlVolumeFace& scvf,
125 const FluxVariablesCache& fluxVarsCache)
126 {
127 // Mogensen et al. (1999), does not really revover the single phase value
128 using std::sqrt;
129 const Scalar throatLength = fluxVarsCache.throatLength();
130 const auto nPhaseIdx = fluxVarsCache.nPhaseIdx();
131 const Scalar aNw = fluxVarsCache.throatCrossSectionalArea(nPhaseIdx);
132 const Scalar rEff = 0.5*(sqrt(aNw / M_PI) + fluxVarsCache.throatInscribedRadius());
133 const Scalar result = M_PI/(8*throatLength) * rEff*rEff*rEff*rEff;
134 return result;
135 }
136};
137
139template<class Scalar, class SinglePhaseTransmissibilityLaw>
141{
143
144 template<class Element, class FVElementGeometry, class FluxVariablesCache>
145 static Scalar nonWettingPhaseTransmissibility(const Element& element,
146 const FVElementGeometry& fvGeometry,
147 const typename FVElementGeometry::SubControlVolumeFace& scvf,
148 const FluxVariablesCache& fluxVarsCache)
149 {
150 // Tora et al. (2012), also does not fully recover single-phase value, but is closer
151 using std::sqrt;
152 const auto nPhaseIdx = fluxVarsCache.nPhaseIdx();
153 const Scalar aNw = fluxVarsCache.throatCrossSectionalArea(nPhaseIdx);
154 const Scalar aTot = fluxVarsCache.throatCrossSectionalArea();
155
156 const Scalar result = SinglePhaseTransmissibilityLaw::singlePhaseTransmissibility(element, fvGeometry, scvf, fluxVarsCache)
157 * aNw / aTot;
158
159 return result;
160 }
161};
162
163template<class Scalar>
165{
167
173 template<class Element, class FVElementGeometry, class FluxVariablesCache>
174 static Scalar nonWettingPhaseTransmissibility(const Element& element,
175 const FVElementGeometry& fvGeometry,
176 const typename FVElementGeometry::SubControlVolumeFace& scvf,
177 const FluxVariablesCache& fluxVarsCache)
178 {
179 // Tora et al. (2012), quite close for single-phase value of square
180 using std::sqrt;
181 const Scalar throatLength = fluxVarsCache.throatLength();
182 const auto nPhaseIdx = fluxVarsCache.nPhaseIdx();
183 const Scalar aNw = fluxVarsCache.throatCrossSectionalArea(nPhaseIdx);
184 const Scalar rEff = 0.5*(sqrt(aNw / M_PI) + fluxVarsCache.throatInscribedRadius());
185 const Scalar result = rEff*rEff*aNw / (8.0*throatLength);
186 return result;
187 }
188};
189} // end namespace NonWettingPhaseTransmissibility
190} // end namespace Dumux::PoreNetwork
191
192#endif
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:63
Scalar creviceResistanceFactor(const int cornerIdx) const
Definition: transmissibility2p.hh:80
An empty cache for transmissibility laws using only standard quantities.
constexpr Shape shape(const Scalar shapeFactor) noexcept
Returns the shape for a given shape factor.
Definition: throatproperties.hh:165
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:220
Dune::ReservedVector< Scalar, 4 > cornerHalfAngles(Shape shape)
Returns the corner half angle.
Definition: throatproperties.hh:85
Definition: discretization/porenetwork/fvelementgeometry.hh:24
Definition: emptycache.hh:19
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:174
TODO: evaluate and maybe remove.
Definition: transmissibility2p.hh:118
static Scalar nonWettingPhaseTransmissibility(const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const FluxVariablesCache &fluxVarsCache)
Definition: transmissibility2p.hh:122
TODO: evaluate and maybe remove.
Definition: transmissibility2p.hh:141
static Scalar nonWettingPhaseTransmissibility(const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const FluxVariablesCache &fluxVarsCache)
Definition: transmissibility2p.hh:145
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:35
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:91
This file contains functions related to calculate pore-throat properties.