version 3.8
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// 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_LINEAR_MATERIAL_HH
14#define DUMUX_MATERIAL_FLUIDMATRIX_LINEAR_MATERIAL_HH
15
16#include <algorithm>
20
21namespace Dumux::FluidMatrix {
22
34{
35public:
40 template<class Scalar>
41 struct Params
42 {
43 Params(Scalar pcEntry, Scalar pcMax)
44 : pcEntry_(pcEntry), pcMax_(pcMax)
45 {}
46
47 Scalar pcEntry() const { return pcEntry_; }
48 void setPcEntry(Scalar pcEntry) { pcEntry_ = pcEntry; }
49
50 Scalar pcMax() const { return pcMax_; }
51 void setPcMax(Scalar pcMax) { pcMax_ = pcMax; }
52
53 bool operator== (const Params& p) const
54 {
55 return Dune::FloatCmp::eq(pcEntry(), p.pcEntry(), 1e-6)
56 && Dune::FloatCmp::eq(pcMax(), p.pcMax(), 1e-6);
57 }
58
59 private:
60 Scalar pcEntry_, pcMax_;
61 };
62
67 template<class Scalar = double>
68 static Params<Scalar> makeParams(const std::string& paramGroup)
69 {
70 const auto pcEntry = getParamFromGroup<Scalar>(paramGroup, "LinearPcEntry");
71 const auto pcMax = getParamFromGroup<Scalar>(paramGroup, "LinearPcMax");
72 return {pcEntry, pcMax};
73 }
74
78 template<class Scalar>
79 static Scalar pc(Scalar swe, const Params<Scalar>& params)
80 {
81 return (1.0 - swe)*(params.pcMax() - params.pcEntry()) + params.pcEntry();
82 }
83
87 template<class Scalar>
88 static Scalar swe(Scalar pc, const Params<Scalar>& params)
89 {
90 return 1.0 - (pc - params.pcEntry())/(params.pcMax() - params.pcEntry());
91 }
92
96 template<class Scalar>
97 static Scalar endPointPc(const Params<Scalar>& params)
98 { return params.pcEntry(); }
99
103 template<class Scalar>
104 static Scalar dpc_dswe(Scalar swe, const Params<Scalar>& params)
105 {
106 return params.pcEntry() - params.pcMax();
107 }
108
112 template<class Scalar>
113 static Scalar dswe_dpc(Scalar pc, const Params<Scalar>& params)
114 {
115 return 1.0/(params.pcEntry() - params.pcMax());
116 }
117
121 template<class Scalar>
122 static Scalar krw(Scalar swe, const Params<Scalar>& params)
123 {
124 using std::clamp;
125 return clamp(swe, 0.0, 1.0);
126 }
127
131 template<class Scalar>
132 static Scalar dkrw_dswe(Scalar swe, const Params<Scalar>& params)
133 {
134 return 1.0;
135 }
136
140 template<class Scalar>
141 static Scalar krn(Scalar swe, const Params<Scalar>& params)
142 {
143 using std::clamp;
144 return clamp(1.0-swe, 0.0, 1.0);
145 }
146
150 template<class Scalar>
151 static Scalar dkrn_dswe(Scalar swe, const Params<Scalar>& params)
152 {
153 return -1.0;
154 }
155};
156
157template<typename Scalar = double>
159
160} // end namespace Dumux::FluidMatrix
161
162#endif
Linear capillary pressure and relative permeability <-> saturation relations.
Definition: linearmaterial.hh:34
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:104
static Scalar endPointPc(const Params< Scalar > &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: linearmaterial.hh:97
static Scalar dkrw_dswe(Scalar swe, const Params< Scalar > &params)
The derivative of the relative permeability.
Definition: linearmaterial.hh:132
static Scalar pc(Scalar swe, const Params< Scalar > &params)
The capillary pressure-saturation curve.
Definition: linearmaterial.hh:79
static Scalar krn(Scalar swe, const Params< Scalar > &params)
The relative permeability for the non-wetting phase.
Definition: linearmaterial.hh:141
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:113
static Scalar dkrn_dswe(Scalar swe, const Params< Scalar > &params)
The derivative of the relative permeability.
Definition: linearmaterial.hh:151
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: linearmaterial.hh:68
static Scalar krw(Scalar swe, const Params< Scalar > &params)
The relative permeability for the wetting phase.
Definition: linearmaterial.hh:122
static Scalar swe(Scalar pc, const Params< Scalar > &params)
The inverse saturation-capillary pressure curve.
Definition: linearmaterial.hh:88
Wrapper class to implement regularized material laws (pc-sw, kr-sw) with a conversion policy between ...
Definition: materiallaw.hh:45
Implementation helper for capillary pressure and relative permeability <-> saturation relations for t...
Definition: brookscorey.hh:23
Scalar pcEntry(const Scalar surfaceTension, const Scalar contactAngle, const Scalar inscribedRadius, const Scalar shapeFactor) noexcept
The entry capillary pressure of a pore throat.
Definition: thresholdcapillarypressures.hh:82
A tag to turn off regularization and it's overhead.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The parameter type.
Definition: linearmaterial.hh:42
Scalar pcMax() const
Definition: linearmaterial.hh:50
Params(Scalar pcEntry, Scalar pcMax)
Definition: linearmaterial.hh:43
Scalar pcEntry() const
Definition: linearmaterial.hh:47
void setPcMax(Scalar pcMax)
Definition: linearmaterial.hh:51
bool operator==(const Params &p) const
Definition: linearmaterial.hh:53
void setPcEntry(Scalar pcEntry)
Definition: linearmaterial.hh:48