24#ifndef DUMUX_PNM_1P_FLUXVARIABLESCACHE_HH
25#define DUMUX_PNM_1P_FLUXVARIABLESCACHE_HH
36template<
class AdvectionType>
39 using Scalar =
typename AdvectionType::Scalar;
42 template<
class Problem,
class Element,
class FVElementGeometry,
class ElementVolumeVariables>
44 const Element& element,
45 const FVElementGeometry& fvGeometry,
46 const ElementVolumeVariables& elemVolVars,
47 const typename FVElementGeometry::SubControlVolumeFace& scvf)
49 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
50 throatCrossSectionShape_ = fvGeometry.gridGeometry().throatCrossSectionShape(eIdx);
51 throatShapeFactor_ = fvGeometry.gridGeometry().throatShapeFactor(eIdx);
52 throatCrossSectionalArea_ = problem.spatialParams().throatCrossSectionalArea(element, elemVolVars);
53 throatLength_ = problem.spatialParams().throatLength(element, elemVolVars);
54 throatInscribedRadius_ = problem.spatialParams().throatInscribedRadius(element, elemVolVars);
55 poreToPoreDistance_ = element.geometry().volume();
57 cache_.fill(problem, element, fvGeometry, scvf, elemVolVars, *
this, 0);
58 transmissibility_ = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, scvf, elemVolVars, *
this, 0);
65 {
return throatCrossSectionShape_; }
71 {
return throatShapeFactor_; }
77 {
return transmissibility_; }
83 {
return throatCrossSectionalArea_; }
89 {
return throatLength_; }
95 {
return throatInscribedRadius_; }
101 {
return poreToPoreDistance_; }
111 Scalar throatShapeFactor_;
112 Scalar transmissibility_;
113 Scalar throatCrossSectionalArea_;
114 Scalar throatLength_;
115 Scalar throatInscribedRadius_;
116 Scalar poreToPoreDistance_;
118 typename AdvectionType::Transmissibility::SinglePhaseCache cache_;
This file contains functions related to calculate pore-throat properties.
Definition: discretization/porenetwork/fvelementgeometry.hh:33
Shape
Collection of different pore-throat shapes.
Definition: throatproperties.hh:38
Flux variables cache for the single-phase-flow PNM Store data required for flux calculation.
Definition: porenetwork/1p/fluxvariablescache.hh:38
const auto & singlePhaseFlowVariables() const
Returns the throats's cached flow variables for single-phase flow.
Definition: porenetwork/1p/fluxvariablescache.hh:106
Scalar transmissibility(const int phaseIdx=0) const
Returns the throats's transmissibility.
Definition: porenetwork/1p/fluxvariablescache.hh:76
Scalar poreToPoreDistance() const
Returns the throats's pore-to-pore-center distance.
Definition: porenetwork/1p/fluxvariablescache.hh:100
Scalar throatShapeFactor() const
Returns the throats's shape factor.
Definition: porenetwork/1p/fluxvariablescache.hh:70
Scalar throatInscribedRadius() const
Returns the throats's inscribed radius.
Definition: porenetwork/1p/fluxvariablescache.hh:94
Throat::Shape throatCrossSectionShape() const
Returns the throats's cross-sectional shape.
Definition: porenetwork/1p/fluxvariablescache.hh:64
Scalar throatLength() const
Returns the throats's length.
Definition: porenetwork/1p/fluxvariablescache.hh:88
Scalar throatCrossSectionalArea(const int phaseIdx=0) const
Returns the throats's cross-sectional area.
Definition: porenetwork/1p/fluxvariablescache.hh:82
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf)
Definition: porenetwork/1p/fluxvariablescache.hh:43