3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
linearmaterial.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_LINEAR_MATERIAL_HH
26#define DUMUX_MATERIAL_FLUIDMATRIX_LINEAR_MATERIAL_HH
27
28// remove from here after release 3.3 /////////////
30
31#include <algorithm>
32
33namespace Dumux {
34
47template <class ScalarT, class ParamsT = LinearMaterialParams<ScalarT> >
48class [[deprecated("Use new material laws and FluidMatrix::LinearMaterial instead!")]] LinearMaterial
49{
50public:
51 using Params = ParamsT;
52 using Scalar = typename Params::Scalar;
53
68 static Scalar pc(const Params &params, Scalar swe)
69 {
70 return (1 - swe)*(params.maxPc() - params.entryPc()) + params.entryPc();
71 }
72
87 static Scalar sw(const Params &params, Scalar pc)
88 {
89 return 1 - (pc - params.entryPc())/(params.maxPc() - params.entryPc());
90 }
91
99 static Scalar endPointPc(const Params &params)
100 { return params.entryPc(); }
101
117 static Scalar dpc_dswe(const Params &params, Scalar swe)
118 {
119 return - (params.maxPc() - params.entryPc());
120 }
121
132 static Scalar dsw_dpc(const Params &params, Scalar pc)
133 {
134 return - 1/(params.maxPc() - params.entryPc());
135 }
136
146 static Scalar krw(const Params &params, Scalar swe)
147 {
148 using std::clamp;
149 return clamp(swe, 0.0, 1.0);
150 }
151
161 static Scalar krn(const Params &params, Scalar swe)
162 {
163 using std::clamp;
164 return clamp(1.0-swe, 0.0, 1.0);
165 }
166};
167} // end namespace Dumux
168// remove until here after release 3.3 /////////////
169
170#include <algorithm>
174
175namespace Dumux::FluidMatrix {
176
188{
189public:
194 template<class Scalar>
195 struct Params
196 {
197 Params(Scalar pcEntry, Scalar pcMax)
198 : pcEntry_(pcEntry), pcMax_(pcMax)
199 {}
200
201 Scalar pcEntry() const { return pcEntry_; }
202 void setPcEntry(Scalar pcEntry) { pcEntry_ = pcEntry; }
203
204 Scalar pcMax() const { return pcMax_; }
205 void setPcMax(Scalar pcMax) { pcMax_ = pcMax; }
206
207 bool operator== (const Params& p) const
208 {
209 return Dune::FloatCmp::eq(pcEntry(), p.pcEntry(), 1e-6)
210 && Dune::FloatCmp::eq(pcMax(), p.pcMax(), 1e-6);
211 }
212
213 private:
214 Scalar pcEntry_, pcMax_;
215 };
216
221 template<class Scalar = double>
222 static Params<Scalar> makeParams(const std::string& paramGroup)
223 {
224 const auto pcEntry = getParamFromGroup<Scalar>(paramGroup, "LinearPcEntry");
225 const auto pcMax = getParamFromGroup<Scalar>(paramGroup, "LinearPcMax");
226 return {pcEntry, pcMax};
227 }
228
232 template<class Scalar>
233 static Scalar pc(Scalar swe, const Params<Scalar>& params)
234 {
235 return (1.0 - swe)*(params.pcMax() - params.pcEntry()) + params.pcEntry();
236 }
237
241 template<class Scalar>
242 static Scalar swe(Scalar pc, const Params<Scalar>& params)
243 {
244 return 1.0 - (pc - params.pcEntry())/(params.pcMax() - params.pcEntry());
245 }
246
250 template<class Scalar>
251 static Scalar endPointPc(const Params<Scalar>& params)
252 { return params.pcEntry(); }
253
257 template<class Scalar>
258 static Scalar dpc_dswe(Scalar swe, const Params<Scalar>& params)
259 {
260 return params.pcEntry() - params.pcMax();
261 }
262
266 template<class Scalar>
267 static Scalar dswe_dpc(Scalar pc, const Params<Scalar>& params)
268 {
269 return 1.0/(params.pcEntry() - params.pcMax());
270 }
271
275 template<class Scalar>
276 static Scalar krw(Scalar swe, const Params<Scalar>& params)
277 {
278 using std::clamp;
279 return clamp(swe, 0.0, 1.0);
280 }
281
285 template<class Scalar>
286 static Scalar dkrw_dswe(Scalar swe, const Params<Scalar>& params)
287 {
288 return 1.0;
289 }
290
294 template<class Scalar>
295 static Scalar krn(Scalar swe, const Params<Scalar>& params)
296 {
297 using std::clamp;
298 return clamp(1.0-swe, 0.0, 1.0);
299 }
300
304 template<class Scalar>
305 static Scalar dkrn_dswe(Scalar swe, const Params<Scalar>& params)
306 {
307 return -1.0;
308 }
309};
310
311template<typename Scalar = double>
313
314} // end namespace Dumux::FluidMatrix
315
316#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Parameters for the linear capillary pressure and relative permeability <-> saturation relations.
Implementation helper for capillary pressure and relative permeability <-> saturation relations for t...
A tag to turn off regularization and it's overhead.
Definition: adapt.hh:29
Definition: brookscorey.hh:286
Linear capillary pressure and relative permeability <-> saturation relations.
Definition: linearmaterial.hh:49
static Scalar krw(const Params &params, Scalar swe)
The relative permeability for the wetting phase.
Definition: linearmaterial.hh:146
static Scalar sw(const Params &params, Scalar pc)
The saturation-capillary pressure curve.
Definition: linearmaterial.hh:87
static Scalar dpc_dswe(const Params &params, Scalar swe)
Returns the partial derivative of the capillary pressure w.r.t. the effective saturation.
Definition: linearmaterial.hh:117
static Scalar endPointPc(const Params &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: linearmaterial.hh:99
static Scalar pc(const Params &params, Scalar swe)
The linear capillary pressure-saturation curve.
Definition: linearmaterial.hh:68
typename Params::Scalar Scalar
Definition: linearmaterial.hh:52
static Scalar dsw_dpc(const Params &params, Scalar pc)
Returns the partial derivative of the effective saturation to the capillary pressure.
Definition: linearmaterial.hh:132
static Scalar krn(const Params &params, Scalar swe)
The relative permeability for the nonwetting phase.
Definition: linearmaterial.hh:161
ParamsT Params
Definition: linearmaterial.hh:51
Linear capillary pressure and relative permeability <-> saturation relations.
Definition: linearmaterial.hh:188
static Scalar dpc_dswe(Scalar swe, const Params< Scalar > &params)
The partial derivative of the capillary pressure w.r.t. the effective saturation.
Definition: linearmaterial.hh:258
static Scalar endPointPc(const Params< Scalar > &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: linearmaterial.hh:251
static Scalar dkrw_dswe(Scalar swe, const Params< Scalar > &params)
The derivative of the relative permeability.
Definition: linearmaterial.hh:286
static Scalar pc(Scalar swe, const Params< Scalar > &params)
The capillary pressure-saturation curve.
Definition: linearmaterial.hh:233
static Scalar krn(Scalar swe, const Params< Scalar > &params)
The relative permeability for the non-wetting phase.
Definition: linearmaterial.hh:295
static Scalar dswe_dpc(Scalar pc, const Params< Scalar > &params)
The partial derivative of the effective saturation w.r.t. the capillary pressure.
Definition: linearmaterial.hh:267
static Scalar dkrn_dswe(Scalar swe, const Params< Scalar > &params)
The derivative of the relative permeability.
Definition: linearmaterial.hh:305
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: linearmaterial.hh:222
static Scalar krw(Scalar swe, const Params< Scalar > &params)
The relative permeability for the wetting phase.
Definition: linearmaterial.hh:276
static Scalar swe(Scalar pc, const Params< Scalar > &params)
The inverse saturation-capillary pressure curve.
Definition: linearmaterial.hh:242
The parameter type.
Definition: linearmaterial.hh:196
Scalar pcMax() const
Definition: linearmaterial.hh:204
Params(Scalar pcEntry, Scalar pcMax)
Definition: linearmaterial.hh:197
Scalar pcEntry() const
Definition: linearmaterial.hh:201
void setPcMax(Scalar pcMax)
Definition: linearmaterial.hh:205
bool operator==(const Params &p) const
Definition: linearmaterial.hh:207
void setPcEntry(Scalar pcEntry)
Definition: linearmaterial.hh:202
Wrapper class to implement regularized material laws (pc-sw, kr-sw) with a conversion policy between ...
Definition: materiallaw.hh:57