25#ifndef DUMUX_PNM_THROAT_TRANSMISSIBILITY_2P_HH
26#define DUMUX_PNM_THROAT_TRANSMISSIBILITY_2P_HH
33namespace WettingLayerTransmissibility {
46 template<
class Scalar>
47 static Scalar
beta(
const Scalar alpha,
const Scalar theta,
const Scalar f = 0)
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);
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);
67template<
class Scalar,
class CreviceResistanceFactor = CreviceResistanceFactorZhou>
72 using NumCornerVector = Dune::ReservedVector<Scalar, 4>;
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,
83 const auto& spatialParams = problem.spatialParams();
85 const Scalar contactAngle = spatialParams.contactAngle(element, elemVolVars);
93 {
return beta_[cornerIdx]; }
96 NumCornerVector beta_;
102 template<
class Element,
class FVElementGeometry,
class FluxVariablesCache>
104 const FVElementGeometry& fvGeometry,
105 const typename FVElementGeometry::SubControlVolumeFace& scvf,
106 const FluxVariablesCache& fluxVarsCache)
108 const Scalar throatLength = fluxVarsCache.throatLength();
109 const Scalar rC = fluxVarsCache.curvatureRadius();
111 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
112 const auto shape = fvGeometry.gridGeometry().throatCrossSectionShape(eIdx);
118 result += fluxVarsCache.wettingLayerCrossSectionalArea(i) * rC*rC / (throatLength*fluxVarsCache.wettingLayerFlowVariables().creviceResistanceFactor(i));
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;
151template<
class Scalar,
class SinglePhaseTransmissibilityLaw>
156 template<
class Element,
class FVElementGeometry,
class FluxVariablesCache>
158 const FVElementGeometry& fvGeometry,
159 const typename FVElementGeometry::SubControlVolumeFace& scvf,
160 const FluxVariablesCache& fluxVarsCache)
164 const auto nPhaseIdx = fluxVarsCache.nPhaseIdx();
165 const Scalar aNw = fluxVarsCache.throatCrossSectionalArea(nPhaseIdx);
166 const Scalar aTot = fluxVarsCache.throatCrossSectionalArea();
168 const Scalar result = SinglePhaseTransmissibilityLaw::singlePhaseTransmissibility(element, fvGeometry, scvf, fluxVarsCache)
175template<
class Scalar>
185 template<
class Element,
class FVElementGeometry,
class FluxVariablesCache>
187 const FVElementGeometry& fvGeometry,
188 const typename FVElementGeometry::SubControlVolumeFace& scvf,
189 const FluxVariablesCache& fluxVarsCache)
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);
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:177
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:232
Dune::ReservedVector< Scalar, 4 > cornerHalfAngles(Shape shape)
Returns the corner half angle.
Definition: throatproperties.hh:97
Definition: emptycache.hh:31
Definition: transmissibility2p.hh:36
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
Definition: transmissibility2p.hh:69
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
Definition: transmissibility2p.hh:71
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: evaluate 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: evaluate 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
Definition: transmissibility2p.hh:177
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