13#ifndef DUMUX_PNM_THROAT_TRANSMISSIBILITY_1P_HH
14#define DUMUX_PNM_THROAT_TRANSMISSIBILITY_1P_HH
36 template<
class Problem,
class Element,
class FVElementGeometry,
class ElementVolumeVariables,
class FluxVariablesCache>
38 const Element& element,
39 const FVElementGeometry& fvGeometry,
40 const typename FVElementGeometry::SubControlVolumeFace& scvf,
41 const ElementVolumeVariables& elemVolVars,
42 const FluxVariablesCache& fluxVarsCache,
45 const auto shape = fluxVarsCache.throatCrossSectionShape();
46 const Scalar throatInscribedRadius = fluxVarsCache.throatInscribedRadius();
47 const Scalar throatLength = fluxVarsCache.throatLength();
48 const Scalar area = fluxVarsCache.throatCrossSectionalArea();
54 const Scalar inscribedRadius,
55 const Scalar throatLength,
70 const Scalar sideLength = 2.0*radius;
71 return 28.4*length * 1.0/(sideLength*sideLength*sideLength*sideLength);
75 return 8.0/M_PI * length * 1.0/(radius*radius*radius*radius);
80 static Scalar sqrt3 = sqrt(3.0);
81 const Scalar sideLength = 6.0/sqrt3 * radius;
82 return 320.0/sqrt3 * length * 1.0/(sideLength*sideLength*sideLength*sideLength);
87 const Scalar width = 2*radius;
88 return 12.0/(width*width*width) * length;
93 Scalar height = 2.0*radius;
94 Scalar width = area/height;
100 return 12.0*length / (1.0 - 0.63*(height/width)) * 1.0/(height*height*height*width);
102 default: DUNE_THROW(Dune::InvalidStateException,
"Throat geometry not supported");
112template<
class Scalar,
bool cons
iderPoreResistance = false,
bool interpolateK = false>
115 static_assert(!interpolateK,
"Interpolation of k not implemented");
120 template<
class Problem,
class Element,
class FVElementGeometry,
class ElementVolumeVariables,
class FluxVariablesCache>
122 const Element& element,
123 const FVElementGeometry& fvGeometry,
124 const typename FVElementGeometry::SubControlVolumeFace& scvf,
125 const ElementVolumeVariables& elemVolVars,
126 const FluxVariablesCache& fluxVarsCache,
129 assert(fluxVarsCache.throatCrossSectionShape() !=
Throat::Shape::twoPlates &&
"TwoPlates not supported. Use TransmissibilityBruus instead!");
131 const Scalar
shapeFactor = fluxVarsCache.throatShapeFactor();
132 const Scalar area = fluxVarsCache.throatCrossSectionalArea();
133 const Scalar throatLength = fluxVarsCache.throatLength();
136 if constexpr (!considerPoreResistance)
137 return throatTransmissibility;
140 static const bool considerPoreResistanceOnRuntime = getParamFromGroup<bool>(problem.paramGroup(),
"Transmissibility.ConsiderPoreResistance",
true);
141 if (!considerPoreResistanceOnRuntime)
142 return throatTransmissibility;
144 const auto& scv0 = fvGeometry.scv(scvf.insideScvIdx());
145 const auto& scv1 = fvGeometry.scv(scvf.outsideScvIdx());
147 const auto& spatialParams = problem.spatialParams();
148 const auto elemSol =
elementSolution(element, elemVolVars, fvGeometry);
151 const Scalar poreLength0 = spatialParams.poreLength(element, scv0, elemSol);
152 const Scalar poreLength1 = spatialParams.poreLength(element, scv1, elemSol);
154 const Scalar poreShapeFactor0 = spatialParams.poreShapeFactor(element, scv0, elemSol);
155 const Scalar poreShapeFactor1 = spatialParams.poreShapeFactor(element, scv1, elemSol);
157 const Scalar poreCrossSectionalArea0 = spatialParams.poreCrossSectionalArea(element, scv0, elemSol);
158 const Scalar poreCrossSectionalArea1 = spatialParams.poreCrossSectionalArea(element, scv1, elemSol);
163 return 1 / (1.0/throatTransmissibility + 1.0/poreTransmissibility0 + 1.0/poreTransmissibility1);
179 if (
shapeFactor <= Throat::shapeFactorEquilateralTriangle<Scalar>())
181 else if (
shapeFactor <= Throat::shapeFactorSquare<Scalar>())
192template<
class Scalar>
199 template<
class Problem,
class Element,
class FVElementGeometry,
class ElementVolumeVariables,
class FluxVariablesCache>
201 const Element& element,
202 const FVElementGeometry& fvGeometry,
203 const typename FVElementGeometry::SubControlVolumeFace& scvf,
204 const ElementVolumeVariables& elemVolVars,
205 const FluxVariablesCache& fluxVarsCache,
208 const Scalar throatInscribedRadius = fluxVarsCache.throatInscribedRadius();
209 const Scalar throatLength = fluxVarsCache.throatLength();
217 const Scalar rEff= std::sqrt(4.0/M_PI)*radius ;
218 return M_PI/(8.0*length) *rEff*rEff*rEff*rEff ;
Used by Joeakar-Niasar.
Definition: transmissibility1p.hh:194
static Scalar singlePhaseTransmissibility(const Scalar radius, const Scalar length)
Returns the conductivity of a throat when only one phase is present.
Definition: transmissibility1p.hh:214
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:200
Collection of single-phase flow throat transmissibilities based on Bruus, H. (2011)....
Definition: transmissibility1p.hh:31
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:53
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:37
static Scalar rHydThroat_(const Throat::Shape shape, const Scalar radius, const Scalar length, const Scalar area)
Definition: transmissibility1p.hh:61
Single-phase flow throat transmissibility based on Patzek & Silin (2001) https://doi....
Definition: transmissibility1p.hh:114
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:168
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:121
An empty cache for transmissibility laws using only standard quantities.
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
constexpr Shape shape(const Scalar shapeFactor) noexcept
Returns the shape for a given shape factor.
Definition: throatproperties.hh:165
Shape
Collection of different pore-throat shapes.
Definition: throatproperties.hh:26
Scalar shapeFactor(Shape shape, const Scalar inscribedRadius)
Returns the value of the shape factor for a given shape.
Definition: throatproperties.hh:150
Definition: discretization/porenetwork/fvelementgeometry.hh:24
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Definition: emptycache.hh:19
This file contains functions related to calculate pore-throat properties.