24#ifndef DUMUX_PNM_THROAT_TRANSMISSIBILITY_2P_HH
25#define DUMUX_PNM_THROAT_TRANSMISSIBILITY_2P_HH
32namespace WettingLayerTransmissibility {
45 template<
class Scalar>
46 static Scalar
beta(
const Scalar alpha,
const Scalar theta,
const Scalar f = 0)
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);
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);
66template<
class Scalar,
class CreviceResistanceFactor = CreviceResistanceFactorZhou>
71 using NumCornerVector = Dune::ReservedVector<Scalar, 4>;
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,
82 const auto& spatialParams = problem.spatialParams();
84 const Scalar contactAngle = spatialParams.contactAngle(element, elemVolVars);
92 {
return beta_[cornerIdx]; }
95 NumCornerVector beta_;
101 template<
class Element,
class FVElementGeometry,
class FluxVariablesCache>
103 const FVElementGeometry& fvGeometry,
104 const typename FVElementGeometry::SubControlVolumeFace& scvf,
105 const FluxVariablesCache& fluxVarsCache)
107 const Scalar throatLength = fluxVarsCache.throatLength();
108 const Scalar rC = fluxVarsCache.curvatureRadius();
110 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
111 const auto shape = fvGeometry.gridGeometry().throatCrossSectionShape(eIdx);
117 result += fluxVarsCache.wettingLayerCrossSectionalArea(i) * rC*rC / (throatLength*fluxVarsCache.wettingLayerFlowVariables().creviceResistanceFactor(i));
119 assert(std::isnormal(result));
125namespace NonWettingPhaseTransmissibility {
128template<
class Scalar>
133 template<
class Element,
class FVElementGeometry,
class FluxVariablesCache>
135 const FVElementGeometry& fvGeometry,
136 const typename FVElementGeometry::SubControlVolumeFace& scvf,
137 const FluxVariablesCache& fluxVarsCache)
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));
152template<
class Scalar,
class SinglePhaseTransmissibilityLaw>
157 template<
class Element,
class FVElementGeometry,
class FluxVariablesCache>
159 const FVElementGeometry& fvGeometry,
160 const typename FVElementGeometry::SubControlVolumeFace& scvf,
161 const FluxVariablesCache& fluxVarsCache)
165 const auto nPhaseIdx = fluxVarsCache.nPhaseIdx();
166 const Scalar aNw = fluxVarsCache.throatCrossSectionalArea(nPhaseIdx);
167 const Scalar aTot = fluxVarsCache.throatCrossSectionalArea();
169 const Scalar result = SinglePhaseTransmissibilityLaw::singlePhaseTransmissibility(element, fvGeometry, scvf, fluxVarsCache)
172 assert(std::isnormal(result));
177template<
class Scalar>
187 template<
class Element,
class FVElementGeometry,
class FluxVariablesCache>
189 const FVElementGeometry& fvGeometry,
190 const typename FVElementGeometry::SubControlVolumeFace& scvf,
191 const FluxVariablesCache& fluxVarsCache)
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));
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
Definition: transmissibility2p.hh:35
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
Definition: transmissibility2p.hh:68
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
Definition: transmissibility2p.hh:70
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
Definition: transmissibility2p.hh:179
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