25#ifndef DUMUX_PNM_THROAT_TRANSMISSIBILITY_1P_HH
26#define DUMUX_PNM_THROAT_TRANSMISSIBILITY_1P_HH
48 template<
class Problem,
class Element,
class FVElementGeometry,
class ElementVolumeVariables,
class FluxVariablesCache>
50 const Element& element,
51 const FVElementGeometry& fvGeometry,
52 const typename FVElementGeometry::SubControlVolumeFace& scvf,
53 const ElementVolumeVariables& elemVolVars,
54 const FluxVariablesCache& fluxVarsCache,
57 const auto shape = fluxVarsCache.throatCrossSectionShape();
58 const Scalar throatInscribedRadius = fluxVarsCache.throatInscribedRadius();
59 const Scalar throatLength = fluxVarsCache.throatLength();
60 const Scalar area = fluxVarsCache.throatCrossSectionalArea();
66 const Scalar inscribedRadius,
67 const Scalar throatLength,
82 const Scalar sideLength = 2.0*radius;
83 return 28.4*length * 1.0/(sideLength*sideLength*sideLength*sideLength);
87 return 8.0/M_PI * length * 1.0/(radius*radius*radius*radius);
92 static Scalar sqrt3 = sqrt(3.0);
93 const Scalar sideLength = 6.0/sqrt3 * radius;
94 return 320.0/sqrt3 * length * 1.0/(sideLength*sideLength*sideLength*sideLength);
99 const Scalar width = 2*radius;
100 return 12.0/(width*width*width) * length;
105 Scalar height = 2.0*radius;
106 Scalar width = area/height;
112 return 12.0*length / (1.0 - 0.63*(height/width)) * 1.0/(height*height*height*width);
114 default: DUNE_THROW(Dune::InvalidStateException,
"Throat geometry not supported");
124template<
class Scalar,
bool cons
iderPoreResistance = false,
bool interpolateK = false>
127 static_assert(!interpolateK,
"Interpolation of k not implemented");
132 template<
class Problem,
class Element,
class FVElementGeometry,
class ElementVolumeVariables,
class FluxVariablesCache>
134 const Element& element,
135 const FVElementGeometry& fvGeometry,
136 const typename FVElementGeometry::SubControlVolumeFace& scvf,
137 const ElementVolumeVariables& elemVolVars,
138 const FluxVariablesCache& fluxVarsCache,
141 assert(fluxVarsCache.throatCrossSectionShape() !=
Throat::Shape::twoPlates &&
"TwoPlates not supported. Use TransmissibilityBruus instead!");
143 const Scalar
shapeFactor = fluxVarsCache.throatShapeFactor();
144 const Scalar area = fluxVarsCache.throatCrossSectionalArea();
145 const Scalar throatLength = fluxVarsCache.throatLength();
148 if constexpr (!considerPoreResistance)
149 return throatTransmissibility;
152 static const bool considerPoreResistanceOnRuntime = getParamFromGroup<bool>(problem.paramGroup(),
"Transmissibility.ConsiderPoreResistance",
true);
153 if (!considerPoreResistanceOnRuntime)
154 return throatTransmissibility;
156 const auto& scv0 = fvGeometry.scv(scvf.insideScvIdx());
157 const auto& scv1 = fvGeometry.scv(scvf.outsideScvIdx());
159 const auto& spatialParams = problem.spatialParams();
160 const auto elemSol =
elementSolution(element, elemVolVars, fvGeometry);
163 const Scalar poreLength0 = spatialParams.poreLength(element, scv0, elemSol);
164 const Scalar poreLength1 = spatialParams.poreLength(element, scv1, elemSol);
166 const Scalar poreShapeFactor0 = spatialParams.poreShapeFactor(element, scv0, elemSol);
167 const Scalar poreShapeFactor1 = spatialParams.poreShapeFactor(element, scv1, elemSol);
169 const Scalar poreCrossSectionalArea0 = spatialParams.poreCrossSectionalArea(element, scv0, elemSol);
170 const Scalar poreCrossSectionalArea1 = spatialParams.poreCrossSectionalArea(element, scv1, elemSol);
175 return 1 / (1.0/throatTransmissibility + 1.0/poreTransmissibility0 + 1.0/poreTransmissibility1);
191 if (
shapeFactor <= Throat::shapeFactorEquilateralTriangle<Scalar>())
193 else if (
shapeFactor <= Throat::shapeFactorSquare<Scalar>())
204template<
class Scalar>
211 template<
class Problem,
class Element,
class FVElementGeometry,
class ElementVolumeVariables,
class FluxVariablesCache>
213 const Element& element,
214 const FVElementGeometry& fvGeometry,
215 const typename FVElementGeometry::SubControlVolumeFace& scvf,
216 const ElementVolumeVariables& elemVolVars,
217 const FluxVariablesCache& fluxVarsCache,
220 const Scalar throatInscribedRadius = fluxVarsCache.throatInscribedRadius();
221 const Scalar throatLength = fluxVarsCache.throatLength();
229 const Scalar rEff= std::sqrt(4.0/M_PI)*radius ;
230 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==DiscretizationMethods::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:34
constexpr Shape shape(const Scalar shapeFactor) noexcept
Returns the shape for a given shape factor.
Definition: throatproperties.hh:177
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:162
Definition: emptycache.hh:31
Collection of single-phase flow throat transmissibilities based on Bruus, H. (2011)....
Definition: transmissibility1p.hh:43
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:65
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:49
static Scalar rHydThroat_(const Throat::Shape shape, const Scalar radius, const Scalar length, const Scalar area)
Definition: transmissibility1p.hh:73
Single-phase flow throat transmissibility based on Patzek & Silin (2001) https://doi....
Definition: transmissibility1p.hh:126
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:180
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:133
Used by Joeakar-Niasar.
Definition: transmissibility1p.hh:206
static Scalar singlePhaseTransmissibility(const Scalar radius, const Scalar length)
Returns the conductivity of a throat when only one phase is present.
Definition: transmissibility1p.hh:226
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:212