Loading [MathJax]/extensions/tex2jax.js
3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
singleshapelocalrules.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_SINGLE_SHAPE_LOCAL_RULES_HH
25#define DUMUX_PNM_2P_SINGLE_SHAPE_LOCAL_RULES_HH
26
31
33
38template<class ScalarType,
39 class BaseLaw,
40 class Regularization = Dumux::FluidMatrix::NoRegularization>
41class SingleShapeTwoPLocalRules : public Dumux::FluidMatrix::Adapter<SingleShapeTwoPLocalRules<ScalarType, BaseLaw, Regularization>, Dumux::FluidMatrix::PcKrSw>
42{
43public:
44
45 using Scalar = ScalarType;
46
47 using BasicParams = typename BaseLaw::template Params<Scalar>;
48 using RegularizationParams = typename Regularization::template Params<Scalar>;
49
50 static constexpr bool supportsMultipleGeometries()
51 { return false; }
52
53 template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
54 void updateParams(const SpatialParams& spatialParams,
55 const Element& element,
56 const SubControlVolume& scv,
57 const ElemSol& elemSol)
58 {
59 basicParams_.update(spatialParams, element, scv, elemSol);
60 regularization_.updateParams(spatialParams, element, scv, elemSol);
61 }
62
66 static constexpr int numFluidPhases()
67 { return 2; }
68
72 static constexpr bool isRegularized()
73 { return !std::is_same<Regularization, Dumux::FluidMatrix::NoRegularization>::value; }
74
79 SingleShapeTwoPLocalRules(const BasicParams& baseParams = {},
80 const RegularizationParams& regParams = {},
81 const std::string& paramGroup = "")
82 : basicParams_(baseParams)
83
84 {
85 if constexpr (isRegularized())
86 regularization_.init(this, regParams, paramGroup);
87 }
88
92 template<bool enableRegularization = isRegularized()>
93 Scalar pc(const Scalar sw) const
94 {
95 if constexpr (enableRegularization)
96 {
97 const auto regularized = regularization_.pc(sw);
98 if (regularized)
99 return regularized.value();
100 }
101
102 return BaseLaw::pc(sw, basicParams_);
103 }
104
108 template<bool enableRegularization = isRegularized()>
109 Scalar dpc_dsw(const Scalar sw) const
110 {
111 if constexpr (enableRegularization)
112 {
113 const auto regularized = regularization_.dpc_dsw(sw);
114 if (regularized)
115 return regularized.value();
116 }
117
118 return BaseLaw::dpc_dsw(sw, basicParams_);
119 }
120
124 template<bool enableRegularization = isRegularized()>
125 Scalar sw(const Scalar pc) const
126 {
127 if constexpr (enableRegularization)
128 {
129 const auto regularized = regularization_.sw(pc);
130 if (regularized)
131 return regularized.value();
132 }
133
134 return BaseLaw::sw(pc, basicParams_);
135 }
136
140 template<bool enableRegularization = isRegularized()>
141 Scalar dsw_dpc(const Scalar pc) const
142 {
143 if constexpr (enableRegularization)
144 {
145 const auto regularized = regularization_.dsw_dpc(pc);
146 if (regularized)
147 return regularized.value();
148 }
149
150 return BaseLaw::dsw_dpc(pc, basicParams_);
151 }
152
157 template<bool enableRegularization = isRegularized()>
158 Scalar krw(const Scalar sw) const
159 { return 1.0; }
160
165 template<bool enableRegularization = isRegularized()>
167 { return 0; }
168
173 template<bool enableRegularization = isRegularized()>
174 Scalar krn(const Scalar sw) const
175 { return 1.0; }
176
181 template<bool enableRegularization = isRegularized()>
183 { return 0.0; }
184
189 {
190 return basicParams_ == o.basicParams_
191 && regularization_ == o.regularization_;
192 }
193
198 { return basicParams_; }
199
200
201private:
202 BasicParams basicParams_;
203 Regularization regularization_;
204};
205
206} // end namespace Dumux::PoreNetwork::FluidMatrix
207
208#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
A tag to turn off regularization and it's overhead.
This file contains functions related to calculate pore-body properties.
Definition: localrulesforplatonicbody.hh:38
A tag to turn off regularization and it's overhead.
Definition: noregularization.hh:34
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:67
Base class for all standard pore-local pc-Sw curves.
Definition: singleshapelocalrules.hh:42
void updateParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: singleshapelocalrules.hh:54
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: singleshapelocalrules.hh:109
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition: singleshapelocalrules.hh:158
static constexpr bool supportsMultipleGeometries()
Definition: singleshapelocalrules.hh:50
ScalarType Scalar
Definition: singleshapelocalrules.hh:45
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: singleshapelocalrules.hh:141
Scalar sw(const Scalar pc) const
The saturation-capillary-pressure curve.
Definition: singleshapelocalrules.hh:125
Scalar dkrn_dsw(const Scalar sw) const
The derivative of the relative permeability for the non-wetting phase w.r.t. saturation.
Definition: singleshapelocalrules.hh:182
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: singleshapelocalrules.hh:93
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition: singleshapelocalrules.hh:166
typename BaseLaw::template Params< Scalar > BasicParams
Definition: singleshapelocalrules.hh:47
typename Regularization::template Params< Scalar > RegularizationParams
Definition: singleshapelocalrules.hh:48
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: singleshapelocalrules.hh:72
const BasicParams & basicParams() const
Return the base law's parameters.
Definition: singleshapelocalrules.hh:197
bool operator==(const SingleShapeTwoPLocalRules &o) const
Equality comparison with another instance.
Definition: singleshapelocalrules.hh:188
SingleShapeTwoPLocalRules(const BasicParams &baseParams={}, const RegularizationParams &regParams={}, const std::string &paramGroup="")
Construct from parameter structs.
Definition: singleshapelocalrules.hh:79
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition: singleshapelocalrules.hh:174
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: singleshapelocalrules.hh:66