3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porenetwork2p.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 *****************************************************************************/
25#ifndef DUMUX_PNM2P_SPATIAL_PARAMS_HH
26#define DUMUX_PNM2P_SPATIAL_PARAMS_HH
27
30
33
34#include "porenetworkbase.hh"
35
36namespace Dumux::PoreNetwork {
37
46template<class GridGeometry, class Scalar, class LocalRules, class Implementation>
48: public BaseSpatialParams<GridGeometry, Scalar, TwoPBaseSpatialParams<GridGeometry, Scalar, LocalRules, Implementation>>
49{
51 using GridView = typename GridGeometry::GridView;
52 using SubControlVolume = typename GridGeometry::SubControlVolume;
53 using Element = typename GridView::template Codim<0>::Entity;
54
55public:
56
57 TwoPBaseSpatialParams(std::shared_ptr<const GridGeometry> gridGeometry)
59 {
60 if (!gridGeometry->useSameGeometryForAllPores() && LocalRules::supportsMultipleGeometries())
61 DUNE_THROW(Dune::InvalidStateException, "Your MaterialLaw does not support multiple pore body shapes.");
62
63 if (this->gridGeometry().useSameShapeForAllThroats())
64 {
65 cornerHalfAngles_.resize(1);
66 const auto& shape = this->gridGeometry().throatCrossSectionShape(/*eIdx*/0);
67 cornerHalfAngles_[0] = Throat::cornerHalfAngles<Scalar>(shape);
68 }
69 else
70 {
71 cornerHalfAngles_.resize(this->gridGeometry().gridView().size(0));
72 for (auto&& element : elements(this->gridGeometry().gridView()))
73 {
74 const auto eIdx = this->gridGeometry().elementMapper().index(element);
75 const auto& shape = this->gridGeometry().throatCrossSectionShape(eIdx);
76 cornerHalfAngles_[eIdx] = Throat::cornerHalfAngles<Scalar>(shape);
77 }
78 }
79 }
80
84 template<class FS, class ElementVolumeVariables>
85 int wettingPhase(const Element&, const ElementVolumeVariables& elemVolVars) const
86 { return 0; }
87
91 template<class FS, class ElementSolutionVector>
92 int wettingPhase(const Element&, const SubControlVolume& scv, const ElementSolutionVector& elemSol) const
93 { return 0; }
94
102 template<class ElementVolumeVariables>
103 Scalar contactAngle(const Element& element,
104 const ElementVolumeVariables& elemVolVars) const
105 {
106 static const Scalar theta = getParam<Scalar>("SpatialParams.ContactAngle", 0.0);
107 return theta;
108 }
109
118 template<class ElementSolutionVector>
119 Scalar contactAngle(const Element& element,
120 const SubControlVolume& scv,
121 const ElementSolutionVector& elemSol) const
122 {
123 static const Scalar theta = getParam<Scalar>("SpatialParams.ContactAngle", 0.0);
124 return theta;
125 }
126
132 template<class ElementVolumeVariables>
133 const Scalar pcEntry(const Element& element, const ElementVolumeVariables& elemVolVars) const
134 {
135 const auto eIdx = this->gridGeometry().elementMapper().index(element);
136 // take the average of both adjacent pores TODO: is this correct?
137 const Scalar surfaceTension = 0.5*(elemVolVars[0].surfaceTension() + elemVolVars[1].surfaceTension());
139 this->asImp_().contactAngle(element, elemVolVars),
140 this->asImp_().throatInscribedRadius(element, elemVolVars),
141 this->gridGeometry().throatShapeFactor(eIdx));
142 }
143
149 template<class ElementVolumeVariables>
150 const Scalar pcSnapoff(const Element& element, const ElementVolumeVariables& elemVolVars) const
151 {
152 // take the average of both adjacent pores TODO: is this correct?
153 const Scalar surfaceTension = 0.5*(elemVolVars[0].surfaceTension() + elemVolVars[1].surfaceTension());
155 this->asImp_().contactAngle(element, elemVolVars),
156 this->asImp_().throatInscribedRadius(element, elemVolVars));
157 }
158
166 template<class ElementSolution>
167 Scalar surfaceTension(const Element& element,
168 const SubControlVolume& scv,
169 const ElementSolution& elemSol) const
170 {
171 static const Scalar gamma = getParam<Scalar>("SpatialParams.SurfaceTension", 0.0725); // default to surface tension of water/air
172 return gamma;
173 }
174
183 template<class ElementSolution>
184 auto fluidMatrixInteraction(const Element& element,
185 const SubControlVolume& scv,
186 const ElementSolution& elemSol) const
187 {
188 const auto params = typename LocalRules::BasicParams(*this, element, scv, elemSol);
189 return makeFluidMatrixInteraction(LocalRules(params, "SpatialParams"));
190 }
191
192 const Dune::ReservedVector<Scalar, 4>& cornerHalfAngles(const Element& element) const
193 {
194 if (this->gridGeometry().useSameShapeForAllThroats())
195 return cornerHalfAngles_[0];
196 else
197 {
198 const auto eIdx = this->gridGeometry().gridView().indexSet().index(element);
199 return cornerHalfAngles_[eIdx];
200 }
201 }
202
203private:
204 std::vector<Dune::ReservedVector<Scalar, 4>> cornerHalfAngles_;
205};
206
207// TODO docme
208template<class GridGeometry, class Scalar, class MaterialLawT>
209class TwoPDefaultSpatialParams : public TwoPBaseSpatialParams<GridGeometry, Scalar, MaterialLawT,
210 TwoPDefaultSpatialParams<GridGeometry, Scalar, MaterialLawT>>
211{
212 using ParentType = TwoPBaseSpatialParams<GridGeometry, Scalar, MaterialLawT,
214public:
216};
217
218} // namespace Dumux::PoreNetwork
219
220#endif
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Specification of threshold capillary pressures for the PNM.
The base class for spatial parameters for pore-network models.
This file contains functions related to calculate pore-body properties.
This file contains functions related to calculate pore-throat properties.
auto makeFluidMatrixInteraction(Laws &&... laws)
Helper function to create an FluidMatrixInteraction object containing an arbitrary number of fluid ma...
Definition: fluidmatrixinteraction.hh:51
Definition: discretization/porenetwork/fvelementgeometry.hh:33
constexpr Shape shape(const Scalar shapeFactor) noexcept
Returns the shape for a given shape factor.
Definition: throatproperties.hh:176
Base class for the finite volume geometry for porenetwork models.
Definition: discretization/porenetwork/gridgeometry.hh:488
static constexpr Scalar pcSnapoff(const Scalar surfaceTension, const Scalar contactAngle, const Scalar inscribedRadius) noexcept
The snap-off capillary pressure of a pore throat.
Definition: thresholdcapillarypressures.hh:36
static constexpr Scalar pcEntry(const Scalar surfaceTension, const Scalar contactAngle, const Scalar inscribedRadius, const Scalar shapeFactor) noexcept
The entry capillary pressure of a pore throat.
Definition: thresholdcapillarypressures.hh:54
The base class for spatial parameters for pore-network models.
Definition: porenetwork2p.hh:49
const Scalar pcSnapoff(const Element &element, const ElementVolumeVariables &elemVolVars) const
Return the element (throat) specific snap-off capillary pressure .
Definition: porenetwork2p.hh:150
int wettingPhase(const Element &, const SubControlVolume &scv, const ElementSolutionVector &elemSol) const
The index of the wetting phase within a pore body.
Definition: porenetwork2p.hh:92
int wettingPhase(const Element &, const ElementVolumeVariables &elemVolVars) const
The index of the wetting phase within a pore throat.
Definition: porenetwork2p.hh:85
Scalar contactAngle(const Element &element, const SubControlVolume &scv, const ElementSolutionVector &elemSol) const
The contact angle within a pore body .
Definition: porenetwork2p.hh:119
TwoPBaseSpatialParams(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: porenetwork2p.hh:57
Scalar contactAngle(const Element &element, const ElementVolumeVariables &elemVolVars) const
The contact angle within a pore throat .
Definition: porenetwork2p.hh:103
auto fluidMatrixInteraction(const Element &element, const SubControlVolume &scv, const ElementSolution &elemSol) const
Returns the parameter object for the pore-local pc-Sw law.
Definition: porenetwork2p.hh:184
const Dune::ReservedVector< Scalar, 4 > & cornerHalfAngles(const Element &element) const
Definition: porenetwork2p.hh:192
Scalar surfaceTension(const Element &element, const SubControlVolume &scv, const ElementSolution &elemSol) const
Returns the surface tension .
Definition: porenetwork2p.hh:167
const Scalar pcEntry(const Element &element, const ElementVolumeVariables &elemVolVars) const
Return the element (throat) specific entry capillary pressure .
Definition: porenetwork2p.hh:133
Definition: porenetwork2p.hh:211
The base class for spatial parameters for pore-network models.
Definition: porenetworkbase.hh:92
TwoPBaseSpatialParams< GridGeometry, Scalar, LocalRules, Implementation > & asImp_()
Definition: porenetworkbase.hh:308
Scalar throatInscribedRadius(const Element &element, const ElementVolumeVariables &elemVolVars) const
Inscribed radius of the throat . Can be solution-dependent.
Definition: porenetworkbase.hh:135
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: porenetworkbase.hh:203
const GridView & gridView() const
Returns a reference to the gridview.
Definition: porenetworkbase.hh:176