3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 *****************************************************************************/
25#ifndef DUMUX_PNM_2P_SINGLE_SHAPE_LOCAL_RULES_HH
26#define DUMUX_PNM_2P_SINGLE_SHAPE_LOCAL_RULES_HH
27
32
34
40template<class ScalarType,
41 class BaseLaw,
42 class Regularization = Dumux::FluidMatrix::NoRegularization>
43class SingleShapeTwoPLocalRules : public Dumux::FluidMatrix::Adapter<SingleShapeTwoPLocalRules<ScalarType, BaseLaw, Regularization>, Dumux::FluidMatrix::PcKrSw>
44{
45public:
46
47 using Scalar = ScalarType;
48
49 using BasicParams = typename BaseLaw::template Params<Scalar>;
50 using RegularizationParams = typename Regularization::template Params<Scalar>;
51
52 static constexpr bool supportsMultipleGeometries()
53 { return false; }
54
55 template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
56 void updateParams(const SpatialParams& spatialParams,
57 const Element& element,
58 const SubControlVolume& scv,
59 const ElemSol& elemSol)
60 {
61 basicParams_.update(spatialParams, element, scv, elemSol);
62 regularization_.updateParams(spatialParams, element, scv, elemSol);
63 }
64
68 static constexpr int numFluidPhases()
69 { return 2; }
70
74 static constexpr bool isRegularized()
75 { return !std::is_same<Regularization, Dumux::FluidMatrix::NoRegularization>::value; }
76
81 SingleShapeTwoPLocalRules(const BasicParams& baseParams = {},
82 const RegularizationParams& regParams = {},
83 const std::string& paramGroup = "")
84 : basicParams_(baseParams)
85
86 {
87 if constexpr (isRegularized())
88 regularization_.init(this, regParams, paramGroup);
89 }
90
94 template<bool enableRegularization = isRegularized()>
95 Scalar pc(const Scalar sw) const
96 {
97 if constexpr (enableRegularization)
98 {
99 const auto regularized = regularization_.pc(sw);
100 if (regularized)
101 return regularized.value();
102 }
103
104 return BaseLaw::pc(sw, basicParams_);
105 }
106
110 template<bool enableRegularization = isRegularized()>
111 Scalar dpc_dsw(const Scalar sw) const
112 {
113 if constexpr (enableRegularization)
114 {
115 const auto regularized = regularization_.dpc_dsw(sw);
116 if (regularized)
117 return regularized.value();
118 }
119
120 return BaseLaw::dpc_dsw(sw, basicParams_);
121 }
122
126 template<bool enableRegularization = isRegularized()>
127 Scalar sw(const Scalar pc) const
128 {
129 if constexpr (enableRegularization)
130 {
131 const auto regularized = regularization_.sw(pc);
132 if (regularized)
133 return regularized.value();
134 }
135
136 return BaseLaw::sw(pc, basicParams_);
137 }
138
142 template<bool enableRegularization = isRegularized()>
143 Scalar dsw_dpc(const Scalar pc) const
144 {
145 if constexpr (enableRegularization)
146 {
147 const auto regularized = regularization_.dsw_dpc(pc);
148 if (regularized)
149 return regularized.value();
150 }
151
152 return BaseLaw::dsw_dpc(pc, basicParams_);
153 }
154
159 template<bool enableRegularization = isRegularized()>
160 Scalar krw(const Scalar sw) const
161 { return 1.0; }
162
167 template<bool enableRegularization = isRegularized()>
169 { return 0; }
170
175 template<bool enableRegularization = isRegularized()>
176 Scalar krn(const Scalar sw) const
177 { return 1.0; }
178
183 template<bool enableRegularization = isRegularized()>
185 { return 0.0; }
186
191 {
192 return basicParams_ == o.basicParams_
193 && regularization_ == o.regularization_;
194 }
195
200 { return basicParams_; }
201
202
203private:
204 BasicParams basicParams_;
205 Regularization regularization_;
206};
207
208} // end namespace Dumux::PoreNetwork::FluidMatrix
209
210#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A tag to turn off regularization and it's overhead.
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
This file contains functions related to calculate pore-body properties.
Definition: localrulesforplatonicbody.hh:40
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:44
void updateParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: singleshapelocalrules.hh:56
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: singleshapelocalrules.hh:111
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition: singleshapelocalrules.hh:160
static constexpr bool supportsMultipleGeometries()
Definition: singleshapelocalrules.hh:52
ScalarType Scalar
Definition: singleshapelocalrules.hh:47
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: singleshapelocalrules.hh:143
Scalar sw(const Scalar pc) const
The saturation-capillary-pressure curve.
Definition: singleshapelocalrules.hh:127
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:184
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: singleshapelocalrules.hh:95
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition: singleshapelocalrules.hh:168
typename BaseLaw::template Params< Scalar > BasicParams
Definition: singleshapelocalrules.hh:49
typename Regularization::template Params< Scalar > RegularizationParams
Definition: singleshapelocalrules.hh:50
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: singleshapelocalrules.hh:74
const BasicParams & basicParams() const
Return the base law's parameters.
Definition: singleshapelocalrules.hh:199
bool operator==(const SingleShapeTwoPLocalRules &o) const
Equality comparison with another instance.
Definition: singleshapelocalrules.hh:190
SingleShapeTwoPLocalRules(const BasicParams &baseParams={}, const RegularizationParams &regParams={}, const std::string &paramGroup="")
Construct from parameter structs.
Definition: singleshapelocalrules.hh:81
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition: singleshapelocalrules.hh:176
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: singleshapelocalrules.hh:68
The base class for spatial parameters for pore-network models.
Definition: porenetwork/common/spatialparams.hh:45