13#ifndef DUMUX_PNM_THROAT_TRANSMISSIBILITY_2P_HH
14#define DUMUX_PNM_THROAT_TRANSMISSIBILITY_2P_HH
21namespace WettingLayerTransmissibility {
34 template<
class Scalar>
35 static Scalar
beta(
const Scalar alpha,
const Scalar theta,
const Scalar f = 0)
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);
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);
55template<
class Scalar,
class CreviceResistanceFactor = CreviceResistanceFactorZhou>
60 using NumCornerVector = Dune::ReservedVector<Scalar, 4>;
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,
71 const auto& spatialParams = problem.spatialParams();
73 const Scalar contactAngle = spatialParams.contactAngle(element, elemVolVars);
81 {
return beta_[cornerIdx]; }
84 NumCornerVector beta_;
90 template<
class Element,
class FVElementGeometry,
class FluxVariablesCache>
92 const FVElementGeometry& fvGeometry,
93 const typename FVElementGeometry::SubControlVolumeFace& scvf,
94 const FluxVariablesCache& fluxVarsCache)
96 const Scalar throatLength = fluxVarsCache.throatLength();
97 const Scalar rC = fluxVarsCache.curvatureRadius();
99 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
100 const auto shape = fvGeometry.gridGeometry().throatCrossSectionShape(eIdx);
106 result += fluxVarsCache.wettingLayerCrossSectionalArea(i) * rC*rC / (throatLength*fluxVarsCache.wettingLayerFlowVariables().creviceResistanceFactor(i));
113namespace NonWettingPhaseTransmissibility {
116template<
class Scalar>
121 template<
class Element,
class FVElementGeometry,
class FluxVariablesCache>
123 const FVElementGeometry& fvGeometry,
124 const typename FVElementGeometry::SubControlVolumeFace& scvf,
125 const FluxVariablesCache& fluxVarsCache)
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;
139template<
class Scalar,
class SinglePhaseTransmissibilityLaw>
144 template<
class Element,
class FVElementGeometry,
class FluxVariablesCache>
146 const FVElementGeometry& fvGeometry,
147 const typename FVElementGeometry::SubControlVolumeFace& scvf,
148 const FluxVariablesCache& fluxVarsCache)
152 const auto nPhaseIdx = fluxVarsCache.nPhaseIdx();
153 const Scalar aNw = fluxVarsCache.throatCrossSectionalArea(nPhaseIdx);
154 const Scalar aTot = fluxVarsCache.throatCrossSectionalArea();
156 const Scalar result = SinglePhaseTransmissibilityLaw::singlePhaseTransmissibility(element, fvGeometry, scvf, fluxVarsCache)
163template<
class Scalar>
173 template<
class Element,
class FVElementGeometry,
class FluxVariablesCache>
175 const FVElementGeometry& fvGeometry,
176 const typename FVElementGeometry::SubControlVolumeFace& scvf,
177 const FluxVariablesCache& fluxVarsCache)
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);
Definition: transmissibility2p.hh:59
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
Definition: transmissibility2p.hh:165
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
Definition: transmissibility2p.hh:24
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
Definition: transmissibility2p.hh:57
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.