3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porenetwork/2p/fluxvariablescache.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_PNM_2P_FLUXVARIABLESCACHE_HH
25#define DUMUX_PNM_2P_FLUXVARIABLESCACHE_HH
26
27#include <array>
28#include <algorithm>
29#include <dune/common/reservedvector.hh>
31
32namespace Dumux::PoreNetwork {
33
39template<class AdvectionType, int maxNumCorners = 4>
41{
42 using Scalar = typename AdvectionType::Scalar;
43 static constexpr auto numPhases = 2;
44 using NumCornerVector = Dune::ReservedVector<Scalar, maxNumCorners>;
45
46public:
47
48 template<class Problem, class Element, class FVElementGeometry,
49 class ElementVolumeVariables, class SubControlVolumeFace>
50 void update(const Problem& problem,
51 const Element& element,
52 const FVElementGeometry& fvGeometry,
53 const ElementVolumeVariables& elemVolVars,
54 const SubControlVolumeFace& scvf,
55 bool invaded)
56 {
57 const auto eIdx = fvGeometry.gridGeometry().elementMapper().index(element);
58 throatCrossSectionShape_ = fvGeometry.gridGeometry().throatCrossSectionShape(eIdx);
59 throatShapeFactor_ = fvGeometry.gridGeometry().throatShapeFactor(eIdx);
60 pc_ = std::max(elemVolVars[0].capillaryPressure(), elemVolVars[1].capillaryPressure());
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);
65 invaded_ = invaded;
66 poreToPoreDistance_ = element.geometry().volume();
67
68 // get the non-wetting phase index
69 using FluidSystem = typename ElementVolumeVariables::VolumeVariables::FluidSystem;
70 const auto& spatialParams = problem.spatialParams();
71 nPhaseIdx_ = 1 - spatialParams.template wettingPhase<FluidSystem>(element, elemVolVars);
72
73 // take the average surface tension of both adjacent pores TODO: is this correct?
74 surfaceTension_ = 0.5*(elemVolVars[0].surfaceTension() + elemVolVars[1].surfaceTension());
75
76 const auto& cornerHalfAngles = spatialParams.cornerHalfAngles(element);
77 wettingLayerArea_.clear(); wettingLayerArea_.resize(cornerHalfAngles.size());
78
79 if (invaded) // two-phase flow
80 {
81 const Scalar theta = spatialParams.contactAngle(element, elemVolVars);
82 for (int i = 0; i< cornerHalfAngles.size(); ++i)
84
85 throatCrossSectionalArea_[wPhaseIdx()] = std::accumulate(wettingLayerArea_.begin(), wettingLayerArea_.end(), 0.0);
86 throatCrossSectionalArea_[nPhaseIdx()] = spatialParams.throatCrossSectionalArea(element, elemVolVars) - throatCrossSectionalArea_[wPhaseIdx()];
87 }
88 else // single-phase flow
89 {
90 for (int i = 0; i< cornerHalfAngles.size(); ++i)
91 wettingLayerArea_[i] = 0.0;
92
93 throatCrossSectionalArea_[wPhaseIdx()] = spatialParams.throatCrossSectionalArea(element, elemVolVars);
94 throatCrossSectionalArea_[nPhaseIdx()] = 0.0;
95 }
96
97 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
98 {
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);
102 }
103
104 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
105 transmissibility_[phaseIdx] = AdvectionType::calculateTransmissibility(problem, element, fvGeometry, scvf, elemVolVars, *this, phaseIdx);
106 }
107
112 { return throatCrossSectionShape_; }
113
117 Scalar throatShapeFactor() const
118 { return throatShapeFactor_; }
119
123 Scalar transmissibility(const int phaseIdx) const
124 { return transmissibility_[phaseIdx]; }
125
129 Scalar throatCrossSectionalArea(const int phaseIdx) const
130 { return throatCrossSectionalArea_[phaseIdx]; }
131
136 { return throatCrossSectionalArea_[0] + throatCrossSectionalArea_[1]; }
137
141 Scalar throatLength() const
142 { return throatLength_; }
143
148 { return throatInscribedRadius_; }
149
153 Scalar pcEntry() const
154 { return pcEntry_; }
155
159 Scalar pcSnapoff() const
160 { return pcSnapoff_; }
161
165 Scalar pc() const
166 { return pc_; }
167
171 Scalar surfaceTension() const
172 { return surfaceTension_; }
173
177 bool invaded() const
178 { return invaded_; }
179
183 Scalar curvatureRadius() const
184 { return surfaceTension_ / pc_;}
185
190 Scalar wettingLayerCrossSectionalArea(const int cornerIdx) const
191 { return wettingLayerArea_[cornerIdx]; }
192
196 std::size_t wPhaseIdx() const
197 { return 1 - nPhaseIdx_; }
198
202 std::size_t nPhaseIdx() const
203 { return nPhaseIdx_; }
204
208 const auto& singlePhaseFlowVariables() const
209 { return singlePhaseCache_; }
210
215 { return nonWettingPhaseCache_; }
216
220 const auto& wettingLayerFlowVariables() const
221 { return wettingLayerCache_; }
222
226 Scalar poreToPoreDistance() const
227 { return poreToPoreDistance_; }
228
229private:
230 Throat::Shape throatCrossSectionShape_;
231 Scalar throatShapeFactor_;
232 std::array<Scalar, numPhases> transmissibility_;
233 std::array<Scalar, numPhases> throatCrossSectionalArea_;
234 Scalar throatLength_;
235 Scalar throatInscribedRadius_;
236 Scalar pcEntry_;
237 Scalar pcSnapoff_;
238 Scalar pc_;
239 Scalar surfaceTension_;
240 bool invaded_;
241 NumCornerVector wettingLayerArea_;
242 std::size_t nPhaseIdx_;
243 Scalar poreToPoreDistance_;
244
245 typename AdvectionType::Transmissibility::SinglePhaseCache singlePhaseCache_;
246 typename AdvectionType::Transmissibility::NonWettingPhaseCache nonWettingPhaseCache_;
247 typename AdvectionType::Transmissibility::WettingLayerCache wettingLayerCache_;
248};
249
250} // end Dumux::PoreNetwork
251
252#endif
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