3.5-git
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#include <algorithm>
32
33namespace Dumux::FluidMatrix {
34
46{
47public:
52 template<class Scalar>
53 struct Params
54 {
55 Params(Scalar pcEntry, Scalar pcMax)
56 : pcEntry_(pcEntry), pcMax_(pcMax)
57 {}
58
59 Scalar pcEntry() const { return pcEntry_; }
60 void setPcEntry(Scalar pcEntry) { pcEntry_ = pcEntry; }
61
62 Scalar pcMax() const { return pcMax_; }
63 void setPcMax(Scalar pcMax) { pcMax_ = pcMax; }
64
65 bool operator== (const Params& p) const
66 {
67 return Dune::FloatCmp::eq(pcEntry(), p.pcEntry(), 1e-6)
68 && Dune::FloatCmp::eq(pcMax(), p.pcMax(), 1e-6);
69 }
70
71 private:
72 Scalar pcEntry_, pcMax_;
73 };
74
79 template<class Scalar = double>
80 static Params<Scalar> makeParams(const std::string& paramGroup)
81 {
82 const auto pcEntry = getParamFromGroup<Scalar>(paramGroup, "LinearPcEntry");
83 const auto pcMax = getParamFromGroup<Scalar>(paramGroup, "LinearPcMax");
84 return {pcEntry, pcMax};
85 }
86
90 template<class Scalar>
91 static Scalar pc(Scalar swe, const Params<Scalar>& params)
92 {
93 return (1.0 - swe)*(params.pcMax() - params.pcEntry()) + params.pcEntry();
94 }
95
99 template<class Scalar>
100 static Scalar swe(Scalar pc, const Params<Scalar>& params)
101 {
102 return 1.0 - (pc - params.pcEntry())/(params.pcMax() - params.pcEntry());
103 }
104
108 template<class Scalar>
109 static Scalar endPointPc(const Params<Scalar>& params)
110 { return params.pcEntry(); }
111
115 template<class Scalar>
116 static Scalar dpc_dswe(Scalar swe, const Params<Scalar>& params)
117 {
118 return params.pcEntry() - params.pcMax();
119 }
120
124 template<class Scalar>
125 static Scalar dswe_dpc(Scalar pc, const Params<Scalar>& params)
126 {
127 return 1.0/(params.pcEntry() - params.pcMax());
128 }
129
133 template<class Scalar>
134 static Scalar krw(Scalar swe, const Params<Scalar>& params)
135 {
136 using std::clamp;
137 return clamp(swe, 0.0, 1.0);
138 }
139
143 template<class Scalar>
144 static Scalar dkrw_dswe(Scalar swe, const Params<Scalar>& params)
145 {
146 return 1.0;
147 }
148
152 template<class Scalar>
153 static Scalar krn(Scalar swe, const Params<Scalar>& params)
154 {
155 using std::clamp;
156 return clamp(1.0-swe, 0.0, 1.0);
157 }
158
162 template<class Scalar>
163 static Scalar dkrn_dswe(Scalar swe, const Params<Scalar>& params)
164 {
165 return -1.0;
166 }
167};
168
169template<typename Scalar = double>
171
172} // end namespace Dumux::FluidMatrix
173
174#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Implementation helper for capillary pressure and relative permeability <-> saturation relations for t...
A tag to turn off regularization and it's overhead.
Definition: brookscorey.hh:35
Linear capillary pressure and relative permeability <-> saturation relations.
Definition: linearmaterial.hh:46
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:116
static Scalar endPointPc(const Params< Scalar > &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: linearmaterial.hh:109
static Scalar dkrw_dswe(Scalar swe, const Params< Scalar > &params)
The derivative of the relative permeability.
Definition: linearmaterial.hh:144
static Scalar pc(Scalar swe, const Params< Scalar > &params)
The capillary pressure-saturation curve.
Definition: linearmaterial.hh:91
static Scalar krn(Scalar swe, const Params< Scalar > &params)
The relative permeability for the non-wetting phase.
Definition: linearmaterial.hh:153
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:125
static Scalar dkrn_dswe(Scalar swe, const Params< Scalar > &params)
The derivative of the relative permeability.
Definition: linearmaterial.hh:163
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: linearmaterial.hh:80
static Scalar krw(Scalar swe, const Params< Scalar > &params)
The relative permeability for the wetting phase.
Definition: linearmaterial.hh:134
static Scalar swe(Scalar pc, const Params< Scalar > &params)
The inverse saturation-capillary pressure curve.
Definition: linearmaterial.hh:100
The parameter type.
Definition: linearmaterial.hh:54
Scalar pcMax() const
Definition: linearmaterial.hh:62
Params(Scalar pcEntry, Scalar pcMax)
Definition: linearmaterial.hh:55
Scalar pcEntry() const
Definition: linearmaterial.hh:59
void setPcMax(Scalar pcMax)
Definition: linearmaterial.hh:63
bool operator==(const Params &p) const
Definition: linearmaterial.hh:65
void setPcEntry(Scalar pcEntry)
Definition: linearmaterial.hh:60
Wrapper class to implement regularized material laws (pc-sw, kr-sw) with a conversion policy between ...
Definition: materiallaw.hh:57