version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWOP_MATERIAL_LAW_HH
14#define DUMUX_MATERIAL_FLUIDMATRIX_TWOP_MATERIAL_LAW_HH
15
20
21namespace Dumux::FluidMatrix {
22
40template<class ScalarType,
41 class BaseLaw,
42 class Regularization = NoRegularization,
43 class EffToAbsPolicy = TwoPEffToAbsDefaultPolicy>
44class TwoPMaterialLaw : public Adapter<TwoPMaterialLaw<ScalarType, BaseLaw, Regularization, EffToAbsPolicy>, PcKrSw>
45{
46public:
47
48 using Scalar = ScalarType;
49
50 using BasicParams = typename BaseLaw::template Params<Scalar>;
51 using EffToAbsParams = typename EffToAbsPolicy::template Params<Scalar>;
52 using RegularizationParams = typename Regularization::template Params<Scalar>;
53
54 using EffToAbs = EffToAbsPolicy;
55
59 static constexpr int numFluidPhases()
60 { return 2; }
61
65 static constexpr bool isRegularized()
66 { return !std::is_same<Regularization, NoRegularization>::value; }
67
72 TwoPMaterialLaw() = delete;
73
78 explicit TwoPMaterialLaw(const std::string& paramGroup)
79 : basicParams_(makeBasicParams(paramGroup))
80 , effToAbsParams_(makeEffToAbsParams(paramGroup))
81 {
82 if constexpr (isRegularized())
83 regularization_.init(this, paramGroup);
84 }
85
90 TwoPMaterialLaw(const BasicParams& baseParams,
92 const RegularizationParams& regParams = {})
93 : basicParams_(baseParams)
94 , effToAbsParams_(effToAbsParams)
95 {
96 if constexpr (isRegularized())
97 regularization_.init(this, baseParams, effToAbsParams, regParams);
98 }
99
103 template<bool enableRegularization = isRegularized()>
104 Scalar pc(const Scalar sw) const
105 {
106 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
107 if constexpr (enableRegularization)
108 {
109 const auto regularized = regularization_.pc(swe);
110 if (regularized)
111 return regularized.value();
112 }
113
114 return BaseLaw::pc(swe, basicParams_);
115 }
116
120 template<bool enableRegularization = isRegularized()>
121 Scalar dpc_dsw(const Scalar sw) const
122 {
123 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
124 if constexpr (enableRegularization)
125 {
126 const auto regularized = regularization_.dpc_dswe(swe);
127 if (regularized)
128 return regularized.value()*EffToAbs::dswe_dsw(effToAbsParams_);
129 }
130
131 return BaseLaw::dpc_dswe(swe, basicParams_)*EffToAbs::dswe_dsw(effToAbsParams_);
132 }
133
138 {
139 return BaseLaw::endPointPc(basicParams_);
140 }
141
145 template<bool enableRegularization = isRegularized()>
146 Scalar sw(const Scalar pc) const
147 {
148 if constexpr (enableRegularization)
149 {
150 const auto regularized = regularization_.swe(pc);
151 if (regularized)
152 return EffToAbs::sweToSw(regularized.value(), effToAbsParams_);
153 }
154
155 return EffToAbs::sweToSw(BaseLaw::swe(pc, basicParams_), effToAbsParams_);
156 }
157
161 template<bool enableRegularization = isRegularized()>
162 Scalar dsw_dpc(const Scalar pc) const
163 {
164 if constexpr (enableRegularization)
165 {
166 const auto regularized = regularization_.dswe_dpc(pc);
167 if (regularized)
168 return regularized.value()*EffToAbs::dsw_dswe(effToAbsParams_);
169 }
170
171 return BaseLaw::dswe_dpc(pc, basicParams_)*EffToAbs::dsw_dswe(effToAbsParams_);
172 }
173
177 template<bool enableRegularization = isRegularized()>
178 Scalar krw(const Scalar sw) const
179 {
180 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
181 if constexpr (enableRegularization)
182 {
183 const auto regularized = regularization_.krw(swe);
184 if (regularized)
185 return regularized.value();
186 }
187
188 return BaseLaw::krw(swe, basicParams_);
189 }
190
194 template<bool enableRegularization = isRegularized()>
196 {
197 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
198 if constexpr (enableRegularization)
199 {
200 const auto regularized = regularization_.dkrw_dswe(swe);
201 if (regularized)
202 return regularized.value()*EffToAbs::dswe_dsw(effToAbsParams_);
203 }
204
205 return BaseLaw::dkrw_dswe(swe, basicParams_)*EffToAbs::dswe_dsw(effToAbsParams_);
206 }
207
211 template<bool enableRegularization = isRegularized()>
212 Scalar krn(const Scalar sw) const
213 {
214 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
215 if constexpr (enableRegularization)
216 {
217 const auto regularized = regularization_.krn(swe);
218 if (regularized)
219 return regularized.value();
220 }
221
222 return BaseLaw::krn(swe, basicParams_);
223 }
224
228 template<bool enableRegularization = isRegularized()>
230 {
231 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
232 if constexpr (enableRegularization)
233 {
234 const auto regularized = regularization_.dkrn_dswe(swe);
235 if (regularized)
236 return regularized.value()*EffToAbs::dswe_dsw(effToAbsParams_);
237 }
238
239 return BaseLaw::dkrn_dswe(swe, basicParams_)*EffToAbs::dswe_dsw(effToAbsParams_);
240 }
241
245 bool operator== (const TwoPMaterialLaw& o) const
246 {
247 return basicParams_ == o.basicParams_
248 && effToAbsParams_ == o.effToAbsParams_
249 && regularization_ == o.regularization_;
250 }
251
256 static BasicParams makeBasicParams(const std::string& paramGroup)
257 {
258 return BaseLaw::template makeParams<Scalar>(paramGroup);
259 }
260
265 { return basicParams_; }
266
271 static EffToAbsParams makeEffToAbsParams(const std::string& paramGroup)
272 {
273 return EffToAbs::template makeParams<Scalar>(paramGroup);
274 }
275
280 { return effToAbsParams_; }
281
282private:
283 BasicParams basicParams_;
284 EffToAbsParams effToAbsParams_;
285 Regularization regularization_;
286};
287
288} // end namespace Dumux::FluidMatrix
289
290#endif
Wrapper class to implement regularized material laws (pc-sw, kr-sw) with a conversion policy between ...
Definition: materiallaw.hh:45
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:59
Scalar endPointPc() const
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: materiallaw.hh:137
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: materiallaw.hh:162
TwoPMaterialLaw(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: materiallaw.hh:78
TwoPMaterialLaw(const BasicParams &baseParams, const EffToAbsParams &effToAbsParams={}, const RegularizationParams &regParams={})
Construct from parameter structs.
Definition: materiallaw.hh:90
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: materiallaw.hh:65
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:229
static EffToAbsParams makeEffToAbsParams(const std::string &paramGroup)
Create the parameters of the EffToAbs policy using input file parameters.
Definition: materiallaw.hh:271
ScalarType Scalar
Definition: materiallaw.hh:48
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition: materiallaw.hh:212
const BasicParams & basicParams() const
Return the base law's parameters.
Definition: materiallaw.hh:264
static BasicParams makeBasicParams(const std::string &paramGroup)
Create the base law's parameters using input file parameters.
Definition: materiallaw.hh:256
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition: materiallaw.hh:195
typename Regularization::template Params< Scalar > RegularizationParams
Definition: materiallaw.hh:52
const EffToAbsParams & effToAbsParams() const
Return the parameters of the EffToAbs policy.
Definition: materiallaw.hh:279
bool operator==(const TwoPMaterialLaw &o) const
Equality comparison with another instance.
Definition: materiallaw.hh:245
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: materiallaw.hh:121
EffToAbsPolicy EffToAbs
Definition: materiallaw.hh:54
Scalar sw(const Scalar pc) const
The saturation-capillary pressure curve.
Definition: materiallaw.hh:146
typename EffToAbsPolicy::template Params< Scalar > EffToAbsParams
Definition: materiallaw.hh:51
typename BaseLaw::template Params< Scalar > BasicParams
Definition: materiallaw.hh:50
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: materiallaw.hh:104
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition: materiallaw.hh:178
This is a policy for 2p material laws how to convert absolute to relative saturations and vice versa.
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Definition: brookscorey.hh:23
A tag to turn off regularization and it's overhead.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:55