24#ifndef DUMUX_PNM_THROAT_TRANSMISSIBILITY_1P_HH
25#define DUMUX_PNM_THROAT_TRANSMISSIBILITY_1P_HH
46 template<
class Problem,
class Element,
class FVElementGeometry,
class ElementVolumeVariables,
class FluxVariablesCache>
48 const Element& element,
49 const FVElementGeometry& fvGeometry,
50 const typename FVElementGeometry::SubControlVolumeFace& scvf,
51 const ElementVolumeVariables& elemVolVars,
52 const FluxVariablesCache& fluxVarsCache,
55 const auto shape = fluxVarsCache.throatCrossSectionShape();
56 const Scalar throatInscribedRadius = fluxVarsCache.throatInscribedRadius();
57 const Scalar throatLength = fluxVarsCache.throatLength();
58 const Scalar area = fluxVarsCache.throatCrossSectionalArea();
64 const Scalar inscribedRadius,
65 const Scalar throatLength,
80 const Scalar sideLength = 2.0*radius;
81 return 28.4*length * 1.0/(sideLength*sideLength*sideLength*sideLength);
85 return 8.0/M_PI * length * 1.0/(radius*radius*radius*radius);
90 static Scalar sqrt3 = sqrt(3.0);
91 const Scalar sideLength = 6.0/sqrt3 * radius;
92 return 320.0/sqrt3 * length * 1.0/(sideLength*sideLength*sideLength*sideLength);
97 const Scalar width = 2*radius;
98 return 12.0/(width*width*width) * length;
103 Scalar height = 2.0*radius;
104 Scalar width = area/height;
110 return 12.0*length / (1.0 - 0.63*(height/width)) * 1.0/(height*height*height*width);
112 default: DUNE_THROW(Dune::InvalidStateException,
"Throat geometry not supported");
121template<
class Scalar,
bool cons
iderPoreResistance = false,
bool interpolateK = false>
124 static_assert(!interpolateK,
"Interpolation of k not implemented");
129 template<
class Problem,
class Element,
class FVElementGeometry,
class ElementVolumeVariables,
class FluxVariablesCache>
131 const Element& element,
132 const FVElementGeometry& fvGeometry,
133 const typename FVElementGeometry::SubControlVolumeFace& scvf,
134 const ElementVolumeVariables& elemVolVars,
135 const FluxVariablesCache& fluxVarsCache,
138 const auto shapeFactor = fluxVarsCache.throatShapeFactor();
139 const Scalar area = fluxVarsCache.throatCrossSectionalArea();
140 const Scalar throatLength = fluxVarsCache.throatLength();
143 if constexpr (!considerPoreResistance)
144 return throatTransmissibility;
147 static const bool considerPoreResistanceOnRuntime = getParamFromGroup<bool>(problem.paramGroup(),
"Transmissibility.ConsiderPoreResistance",
true);
148 if (!considerPoreResistanceOnRuntime)
149 return throatTransmissibility;
151 const auto& scv0 = fvGeometry.scv(scvf.insideScvIdx());
152 const auto& scv1 = fvGeometry.scv(scvf.outsideScvIdx());
154 const auto& spatialParams = problem.spatialParams();
155 const auto elemSol =
elementSolution(element, elemVolVars, fvGeometry);
158 const Scalar poreLength0 = spatialParams.poreLength(element, scv0, elemSol);
159 const Scalar poreLength1 = spatialParams.poreLength(element, scv1, elemSol);
161 const auto poreShapeFactor0 = spatialParams.poreShapeFactor(element, scv0, elemSol);
162 const auto poreShapeFactor1 = spatialParams.poreShapeFactor(element, scv1, elemSol);
164 const auto poreCrossSectionalArea0 = spatialParams.poreCrossSectionalArea(element, scv0, elemSol);
165 const auto poreCrossSectionalArea1 = spatialParams.poreCrossSectionalArea(element, scv1, elemSol);
170 return 1 / (1.0/throatTransmissibility + 1.0/poreTransmissibility0 + 1.0/poreTransmissibility1);
186 if (
shapeFactor <= Throat::shapeFactorEquilateralTriangle<Scalar>())
188 else if (
shapeFactor <= Throat::shapeFactorSquare<Scalar>())
199template<
class Scalar>
206 template<
class Problem,
class Element,
class FVElementGeometry,
class ElementVolumeVariables,
class FluxVariablesCache>
208 const Element& element,
209 const FVElementGeometry& fvGeometry,
210 const typename FVElementGeometry::SubControlVolumeFace& scvf,
211 const ElementVolumeVariables& elemVolVars,
212 const FluxVariablesCache& fluxVarsCache,
215 const Scalar throatInscribedRadius = fluxVarsCache.throatInscribedRadius();
216 const Scalar throatLength = fluxVarsCache.throatLength();
224 const Scalar rEff= std::sqrt(4.0/M_PI)*radius ;
225 return M_PI/(8.0*length) *rEff*rEff*rEff*rEff ;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
An empty cache for transmissibility laws using only standard quantities.
This file contains functions related to calculate pore-throat properties.
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethod::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
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
Shape
Collection of different pore-throat shapes.
Definition: throatproperties.hh:38
Scalar shapeFactor(Shape shape, const Scalar inscribedRadius)
Returns the value of the shape factor for a given shape.
Definition: throatproperties.hh:161
Definition: emptycache.hh:30
Collection of single-phase flow throat transmissibilities based on Bruus, H. (2011)....
Definition: transmissibility1p.hh:41
static Scalar singlePhaseTransmissibility(const Throat::Shape shape, const Scalar inscribedRadius, const Scalar throatLength, const Scalar area)
Returns the conductivity of a throat when only one phase is present.
Definition: transmissibility1p.hh:63
static Scalar singlePhaseTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const FluxVariablesCache &fluxVarsCache, const int phaseIdx)
Definition: transmissibility1p.hh:47
static Scalar rHydThroat_(const Throat::Shape shape, const Scalar radius, const Scalar length, const Scalar area)
Definition: transmissibility1p.hh:71
Single-phase flow throat transmissibility based on Patzek & Silin (2001) https://doi....
Definition: transmissibility1p.hh:123
static Scalar singlePhaseTransmissibility(const Scalar shapeFactor, const Scalar length, const Scalar area)
Returns the conductivity of a throat when only one phase is present. See Patzek & Silin (2001)
Definition: transmissibility1p.hh:175
static Scalar singlePhaseTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const FluxVariablesCache &fluxVarsCache, const int phaseIdx)
Definition: transmissibility1p.hh:130
Used by Joeakar-Niasar.
Definition: transmissibility1p.hh:201
static Scalar singlePhaseTransmissibility(const Scalar radius, const Scalar length)
Returns the conductivity of a throat when only one phase is present.
Definition: transmissibility1p.hh:221
static Scalar singlePhaseTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const FluxVariablesCache &fluxVarsCache, const int phaseIdx)
Definition: transmissibility1p.hh:207