3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
materiallaw.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_MATERIAL_FLUIDMATRIX_TWOP_MATERIAL_LAW_HH
26#define DUMUX_MATERIAL_FLUIDMATRIX_TWOP_MATERIAL_LAW_HH
27
32
33namespace Dumux::FluidMatrix {
34
52template<class ScalarType,
53 class BaseLaw,
54 class Regularization = NoRegularization,
55 class EffToAbsPolicy = TwoPEffToAbsDefaultPolicy>
56class TwoPMaterialLaw : public Adapter<TwoPMaterialLaw<ScalarType, BaseLaw, Regularization, EffToAbsPolicy>, PcKrSw>
57{
58public:
59
60 using Scalar = ScalarType;
61
62 using BasicParams = typename BaseLaw::template Params<Scalar>;
63 using EffToAbsParams = typename EffToAbsPolicy::template Params<Scalar>;
64 using RegularizationParams = typename Regularization::template Params<Scalar>;
65
66 using EffToAbs = EffToAbsPolicy;
67
71 static constexpr int numFluidPhases()
72 { return 2; }
73
77 static constexpr bool isRegularized()
78 { return !std::is_same<Regularization, NoRegularization>::value; }
79
84 TwoPMaterialLaw() = delete;
85
90 explicit TwoPMaterialLaw(const std::string& paramGroup)
91 : basicParams_(makeBasicParams(paramGroup))
92 , effToAbsParams_(makeEffToAbsParams(paramGroup))
93 {
94 if constexpr (isRegularized())
95 regularization_.init(this, paramGroup);
96 }
97
102 TwoPMaterialLaw(const BasicParams& baseParams,
103 const EffToAbsParams& effToAbsParams = {},
104 const RegularizationParams& regParams = {})
105 : basicParams_(baseParams)
106 , effToAbsParams_(effToAbsParams)
107 {
108 if constexpr (isRegularized())
109 regularization_.init(this, baseParams, effToAbsParams, regParams);
110 }
111
115 template<bool enableRegularization = isRegularized()>
116 Scalar pc(const Scalar sw) const
117 {
118 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
119 if constexpr (enableRegularization)
120 {
121 const auto regularized = regularization_.pc(swe);
122 if (regularized)
123 return regularized.value();
124 }
125
126 return BaseLaw::pc(swe, basicParams_);
127 }
128
132 template<bool enableRegularization = isRegularized()>
133 Scalar dpc_dsw(const Scalar sw) const
134 {
135 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
136 if constexpr (enableRegularization)
137 {
138 const auto regularized = regularization_.dpc_dswe(swe);
139 if (regularized)
140 return regularized.value()*EffToAbs::dswe_dsw(effToAbsParams_);
141 }
142
143 return BaseLaw::dpc_dswe(swe, basicParams_)*EffToAbs::dswe_dsw(effToAbsParams_);
144 }
145
150 {
151 return BaseLaw::endPointPc(basicParams_);
152 }
153
157 template<bool enableRegularization = isRegularized()>
158 Scalar sw(const Scalar pc) const
159 {
160 if constexpr (enableRegularization)
161 {
162 const auto regularized = regularization_.swe(pc);
163 if (regularized)
164 return EffToAbs::sweToSw(regularized.value(), effToAbsParams_);
165 }
166
167 return EffToAbs::sweToSw(BaseLaw::swe(pc, basicParams_), effToAbsParams_);
168 }
169
173 template<bool enableRegularization = isRegularized()>
174 Scalar dsw_dpc(const Scalar pc) const
175 {
176 if constexpr (enableRegularization)
177 {
178 const auto regularized = regularization_.dswe_dpc(pc);
179 if (regularized)
180 return regularized.value()*EffToAbs::dsw_dswe(effToAbsParams_);
181 }
182
183 return BaseLaw::dswe_dpc(pc, basicParams_)*EffToAbs::dsw_dswe(effToAbsParams_);
184 }
185
189 template<bool enableRegularization = isRegularized()>
190 Scalar krw(const Scalar sw) const
191 {
192 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
193 if constexpr (enableRegularization)
194 {
195 const auto regularized = regularization_.krw(swe);
196 if (regularized)
197 return regularized.value();
198 }
199
200 return BaseLaw::krw(swe, basicParams_);
201 }
202
206 template<bool enableRegularization = isRegularized()>
208 {
209 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
210 if constexpr (enableRegularization)
211 {
212 const auto regularized = regularization_.dkrw_dswe(swe);
213 if (regularized)
214 return regularized.value()*EffToAbs::dswe_dsw(effToAbsParams_);
215 }
216
217 return BaseLaw::dkrw_dswe(swe, basicParams_)*EffToAbs::dswe_dsw(effToAbsParams_);
218 }
219
223 template<bool enableRegularization = isRegularized()>
224 Scalar krn(const Scalar sw) const
225 {
226 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
227 if constexpr (enableRegularization)
228 {
229 const auto regularized = regularization_.krn(swe);
230 if (regularized)
231 return regularized.value();
232 }
233
234 return BaseLaw::krn(swe, basicParams_);
235 }
236
240 template<bool enableRegularization = isRegularized()>
242 {
243 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
244 if constexpr (enableRegularization)
245 {
246 const auto regularized = regularization_.dkrn_dswe(swe);
247 if (regularized)
248 return regularized.value()*EffToAbs::dswe_dsw(effToAbsParams_);
249 }
250
251 return BaseLaw::dkrn_dswe(swe, basicParams_)*EffToAbs::dswe_dsw(effToAbsParams_);
252 }
253
257 bool operator== (const TwoPMaterialLaw& o) const
258 {
259 return basicParams_ == o.basicParams_
260 && effToAbsParams_ == o.effToAbsParams_
261 && regularization_ == o.regularization_;
262 }
263
268 static BasicParams makeBasicParams(const std::string& paramGroup)
269 {
270 return BaseLaw::template makeParams<Scalar>(paramGroup);
271 }
272
277 { return basicParams_; }
278
283 static EffToAbsParams makeEffToAbsParams(const std::string& paramGroup)
284 {
285 return EffToAbs::template makeParams<Scalar>(paramGroup);
286 }
287
292 { return effToAbsParams_; }
293
294private:
295 BasicParams basicParams_;
296 EffToAbsParams effToAbsParams_;
297 Regularization regularization_;
298};
299
300} // end namespace Dumux::FluidMatrix
301
302#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
This is a policy for 2p material laws how to convert absolute to relative saturations and vice versa.
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....
Definition: brookscorey.hh:35
Wrapper class to implement regularized material laws (pc-sw, kr-sw) with a conversion policy between ...
Definition: materiallaw.hh:57
TwoPMaterialLaw()=delete
Deleted default constructor (so we are never in an undefined state)
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: materiallaw.hh:71
Scalar endPointPc() const
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: materiallaw.hh:149
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: materiallaw.hh:174
TwoPMaterialLaw(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: materiallaw.hh:90
TwoPMaterialLaw(const BasicParams &baseParams, const EffToAbsParams &effToAbsParams={}, const RegularizationParams &regParams={})
Construct from parameter structs.
Definition: materiallaw.hh:102
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: materiallaw.hh:77
Scalar dkrn_dsw(const Scalar sw) const
The derivative of the relative permeability for the non-wetting phase w.r.t. saturation.
Definition: materiallaw.hh:241
static EffToAbsParams makeEffToAbsParams(const std::string &paramGroup)
Create the parameters of the EffToAbs policy using input file parameters.
Definition: materiallaw.hh:283
ScalarType Scalar
Definition: materiallaw.hh:60
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition: materiallaw.hh:224
const BasicParams & basicParams() const
Return the base law's parameters.
Definition: materiallaw.hh:276
static BasicParams makeBasicParams(const std::string &paramGroup)
Create the base law's parameters using input file parameters.
Definition: materiallaw.hh:268
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition: materiallaw.hh:207
typename Regularization::template Params< Scalar > RegularizationParams
Definition: materiallaw.hh:64
const EffToAbsParams & effToAbsParams() const
Return the parameters of the EffToAbs policy.
Definition: materiallaw.hh:291
bool operator==(const TwoPMaterialLaw &o) const
Equality comparison with another instance.
Definition: materiallaw.hh:257
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: materiallaw.hh:133
EffToAbsPolicy EffToAbs
Definition: materiallaw.hh:66
Scalar sw(const Scalar pc) const
The saturation-capillary pressure curve.
Definition: materiallaw.hh:158
typename EffToAbsPolicy::template Params< Scalar > EffToAbsParams
Definition: materiallaw.hh:63
typename BaseLaw::template Params< Scalar > BasicParams
Definition: materiallaw.hh:62
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: materiallaw.hh:116
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition: materiallaw.hh:190
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:67