24#ifndef DUMUX_PNM_2P_FLUXVARIABLESCACHE_HH
25#define DUMUX_PNM_2P_FLUXVARIABLESCACHE_HH
29#include <dune/common/reservedvector.hh>
39template<
class AdvectionType,
int maxNumCorners = 4>
42 using Scalar =
typename AdvectionType::Scalar;
43 static constexpr auto numPhases = 2;
44 using NumCornerVector = Dune::ReservedVector<Scalar, maxNumCorners>;
48 template<
class Problem,
class Element,
class FVElementGeometry,
49 class ElementVolumeVariables,
class SubControlVolumeFace>
51 const Element& element,
52 const FVElementGeometry& fvGeometry,
53 const ElementVolumeVariables& elemVolVars,
54 const SubControlVolumeFace& scvf,
57 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
58 throatCrossSectionShape_ = fvGeometry.gridGeometry().throatCrossSectionShape(eIdx);
59 throatShapeFactor_ = fvGeometry.gridGeometry().throatShapeFactor(eIdx);
61 pcEntry_ = problem.spatialParams().pcEntry(element, elemVolVars);
62 pcSnapoff_ = problem.spatialParams().pcSnapoff(element, elemVolVars);
63 throatInscribedRadius_ = problem.spatialParams().throatInscribedRadius(element, elemVolVars);
64 throatLength_ = problem.spatialParams().throatLength(element, elemVolVars);
66 poreToPoreDistance_ = element.geometry().volume();
69 using FluidSystem =
typename ElementVolumeVariables::VolumeVariables::FluidSystem;
70 const auto& spatialParams = problem.spatialParams();
71 nPhaseIdx_ = 1 - spatialParams.template wettingPhase<FluidSystem>(element, elemVolVars);
74 surfaceTension_ = 0.5*(elemVolVars[0].surfaceTension() + elemVolVars[1].surfaceTension());
77 wettingLayerArea_.clear(); wettingLayerArea_.resize(
cornerHalfAngles.size());
81 const Scalar theta = spatialParams.contactAngle(element, elemVolVars);
85 throatCrossSectionalArea_[
wPhaseIdx()] = std::accumulate(wettingLayerArea_.begin(), wettingLayerArea_.end(), 0.0);
86 throatCrossSectionalArea_[
nPhaseIdx()] = spatialParams.throatCrossSectionalArea(element, elemVolVars) - throatCrossSectionalArea_[
wPhaseIdx()];
91 wettingLayerArea_[i] = 0.0;
93 throatCrossSectionalArea_[
wPhaseIdx()] = spatialParams.throatCrossSectionalArea(element, elemVolVars);
94 throatCrossSectionalArea_[
nPhaseIdx()] = 0.0;
97 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
99 singlePhaseCache_.fill(problem, element, fvGeometry, scvf, elemVolVars, *
this, phaseIdx);
100 nonWettingPhaseCache_.fill(problem, element, fvGeometry, scvf, elemVolVars, *
this, phaseIdx);
101 wettingLayerCache_.fill(problem, element, fvGeometry, scvf, elemVolVars, *
this, phaseIdx);
104 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
105 transmissibility_[phaseIdx] = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, scvf, elemVolVars, *
this, phaseIdx);
112 {
return throatCrossSectionShape_; }
118 {
return throatShapeFactor_; }
124 {
return transmissibility_[phaseIdx]; }
130 {
return throatCrossSectionalArea_[phaseIdx]; }
136 {
return throatCrossSectionalArea_[0] + throatCrossSectionalArea_[1]; }
142 {
return throatLength_; }
148 {
return throatInscribedRadius_; }
160 {
return pcSnapoff_; }
172 {
return surfaceTension_; }
184 {
return surfaceTension_ / pc_;}
191 {
return wettingLayerArea_[cornerIdx]; }
197 {
return 1 - nPhaseIdx_; }
203 {
return nPhaseIdx_; }
209 {
return singlePhaseCache_; }
215 {
return nonWettingPhaseCache_; }
221 {
return wettingLayerCache_; }
227 {
return poreToPoreDistance_; }
231 Scalar throatShapeFactor_;
232 std::array<Scalar, numPhases> transmissibility_;
233 std::array<Scalar, numPhases> throatCrossSectionalArea_;
234 Scalar throatLength_;
235 Scalar throatInscribedRadius_;
239 Scalar surfaceTension_;
241 NumCornerVector wettingLayerArea_;
242 std::size_t nPhaseIdx_;
243 Scalar poreToPoreDistance_;
245 typename AdvectionType::Transmissibility::SinglePhaseCache singlePhaseCache_;
246 typename AdvectionType::Transmissibility::NonWettingPhaseCache nonWettingPhaseCache_;
247 typename AdvectionType::Transmissibility::WettingLayerCache wettingLayerCache_;
This file contains functions related to calculate pore-throat properties.
Definition: discretization/porenetwork/fvelementgeometry.hh:33
std::string capillaryPressure() noexcept
I/O name of capillary pressure.
Definition: name.hh:135
constexpr Scalar wettingLayerCrossSectionalArea(const Scalar curvatureRadius, const Scalar contactAngle, const Scalar cornerHalfAngle) noexcept
Return the cross-sectional area of a wetting layer residing in a corner of a throat.
Definition: throatproperties.hh:238
Dune::ReservedVector< Scalar, 4 > cornerHalfAngles(Shape shape)
Returns the corner half angle.
Definition: throatproperties.hh:97
Shape
Collection of different pore-throat shapes.
Definition: throatproperties.hh:38
Flux variables cache for the two-phase-flow PNM Store data required for flux calculation.
Definition: porenetwork/2p/fluxvariablescache.hh:41
Scalar throatCrossSectionalArea() const
Returns the throats's total cross-sectional area.
Definition: porenetwork/2p/fluxvariablescache.hh:135
Scalar pc() const
Returns the capillary pressure within the throat.
Definition: porenetwork/2p/fluxvariablescache.hh:165
bool invaded() const
Returns true if the throat is invaded by the nonwetting phase.
Definition: porenetwork/2p/fluxvariablescache.hh:177
Scalar poreToPoreDistance() const
Returns the throats's pore-to-pore-center distance.
Definition: porenetwork/2p/fluxvariablescache.hh:226
Scalar wettingLayerCrossSectionalArea(const int cornerIdx) const
Returns the cross-sectional area of a wetting layer within one of the throat's corners.
Definition: porenetwork/2p/fluxvariablescache.hh:190
const auto & nonWettingPhaseFlowVariables() const
Returns the throats's cached flow variables for the nonwetting phase.
Definition: porenetwork/2p/fluxvariablescache.hh:214
Scalar curvatureRadius() const
Returns the curvature radius within the throat.
Definition: porenetwork/2p/fluxvariablescache.hh:183
Scalar throatCrossSectionalArea(const int phaseIdx) const
Returns the throats's cross-sectional area for a given phaseIdx.
Definition: porenetwork/2p/fluxvariablescache.hh:129
Scalar transmissibility(const int phaseIdx) const
Returns the throats's transmissibility.
Definition: porenetwork/2p/fluxvariablescache.hh:123
Scalar pcSnapoff() const
Returns the throats's snap-off capillary pressure.
Definition: porenetwork/2p/fluxvariablescache.hh:159
Scalar throatLength() const
Returns the throats's length.
Definition: porenetwork/2p/fluxvariablescache.hh:141
Throat::Shape throatCrossSectionShape() const
Returns the throats's cross-sectional shape.
Definition: porenetwork/2p/fluxvariablescache.hh:111
Scalar throatShapeFactor() const
Returns the throats's shape factor.
Definition: porenetwork/2p/fluxvariablescache.hh:117
const auto & singlePhaseFlowVariables() const
Returns the throats's cached flow variables for single-phase flow.
Definition: porenetwork/2p/fluxvariablescache.hh:208
Scalar pcEntry() const
Returns the throats's entry capillary pressure.
Definition: porenetwork/2p/fluxvariablescache.hh:153
Scalar throatInscribedRadius() const
Returns the throats's inscribed radius.
Definition: porenetwork/2p/fluxvariablescache.hh:147
std::size_t nPhaseIdx() const
Returns the index of the nonwetting phase.
Definition: porenetwork/2p/fluxvariablescache.hh:202
Scalar surfaceTension() const
Returns the surface tension within the throat.
Definition: porenetwork/2p/fluxvariablescache.hh:171
std::size_t wPhaseIdx() const
Returns the index of the wetting phase.
Definition: porenetwork/2p/fluxvariablescache.hh:196
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, bool invaded)
Definition: porenetwork/2p/fluxvariablescache.hh:50
const auto & wettingLayerFlowVariables() const
Returns the throats's cached flow variables for the wetting phase.
Definition: porenetwork/2p/fluxvariablescache.hh:220