3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
smoothedlinearlaw.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_SMOOTHED_LINEAR_LAW_HH
26#define DUMUX_MATERIAL_FLUIDMATRIX_TWOP_SMOOTHED_LINEAR_LAW_HH
27
28
29#include <cmath>
30#include <algorithm>
31
36
37namespace Dumux::FluidMatrix {
38
54template<class ScalarType, class EffToAbsPolicy = TwoPEffToAbsDefaultPolicy>
55class SmoothedLinearLaw : public Adapter<SmoothedLinearLaw<ScalarType, EffToAbsPolicy>, PcKrSw>
56{
57
58public:
59 using Scalar = ScalarType;
60 using EffToAbsParams = typename EffToAbsPolicy::template Params<Scalar>;
61 using EffToAbs = EffToAbsPolicy;
62
66 static constexpr bool isRegularized()
67 { return false; }
68
72 static constexpr int numFluidPhases()
73 { return 2; }
74
79 struct Params
80 {
82 : pe_(pe), pcMax_(pcMax), krLowS_(krLowS), krHighS_(krHighS)
83 {}
84
85 Scalar pe() const { return pe_; }
86 void setPe(Scalar pe) { pe_ = pe; }
87
88 Scalar pcMax() const { return pcMax_; }
89 void setPcMax(Scalar pcMax) { pcMax_ = pcMax; }
90
91 Scalar krLowS() const { return krLowS_; }
92 void setKrLowS(Scalar krLowS) { krLowS_ = krLowS; }
93
94 Scalar krHighS() const { return krHighS_; }
95 void setKrHighS(Scalar krHighS) { krHighS_ = krHighS; }
96
97 bool operator== (const Params& p) const
98 {
99 return Dune::FloatCmp::eq(pe(), p.pe(), 1e-6)
100 && Dune::FloatCmp::eq(pcMax(), p.pcMax(), 1e-6)
101 && Dune::FloatCmp::eq(krLowS(), p.krLowS(), 1e-6)
102 && Dune::FloatCmp::eq(krHighS(), p.krHighS(), 1e-6);
103 }
104
105 private:
106 Scalar pe_, pcMax_, krLowS_, krHighS_;
107 };
108
114
119 explicit SmoothedLinearLaw(const std::string& paramGroup)
120 : SmoothedLinearLaw(makeParams(paramGroup), EffToAbs::template makeParams<Scalar>(paramGroup))
121 {}
122
128 const EffToAbsParams& effToAbsParams = {})
129 : params_(params)
130 , effToAbsParams_(effToAbsParams)
131 , splineM_((1.0 - ((1.0 - params_.krHighS()) + params_.krLowS())/2.0 )
132 / (1.0 - (1.0 - params_.krHighS()) - params_.krLowS()))
133 , splineLowS_(0.0, params_.krLowS(), // x1, x2
134 0.0, params_.krLowS()/2.0, // y1, y2
135 0.0, splineM_) // m1, m2
136 , splineHighS_(params_.krHighS(), 1.0, // x1, x2
137 1.0 - (1.0 - params_.krHighS())/2.0, 1.0, // y1, y2
138 splineM_, 0.0) // m1, m2
139 {}
140
145 static Params makeParams(const std::string& paramGroup)
146 {
147 const auto pe = getParamFromGroup<Scalar>(paramGroup, "SmoothedLinearLawPe");
148 const auto pcMax = getParamFromGroup<Scalar>(paramGroup, "SmoothedLinearLawPcMax");
149 const auto krLowS = getParamFromGroup<Scalar>(paramGroup, "SmoothedLinearLawKrLowS");
150 const auto krHighS = getParamFromGroup<Scalar>(paramGroup, "SmoothedLinearLawKrHighS");
151 return {pe, pcMax, krLowS, krHighS};
152 }
153
158 {
159 return (1.0 - swe)*(params_.pcMax() - params_.pe()) + params_.pe();
160 }
161
166 {
167 return 1.0 - (pc - params_.pe())/(params_.pcMax() - params_.pe());
168 }
169
174 { return params_.pe(); }
175
180 {
181 return params_.pe() - params_.pcMax();
182 }
183
188 {
189 return 1.0/(params_.pe() - params_.pcMax());
190 }
191
196 {
197 return relperm_(swe);
198 }
199
204 {
205 Scalar sne = 1.0 - swe;
206 return relperm_(sne);
207 }
208
213 {
214 return params_ == o.params_
215 && effToAbsParams_ == o.effToAbsParams_;
216 }
217
219 { return effToAbsParams_; }
220
221private:
222
223 Scalar relperm_(Scalar S) const
224 {
225 const Scalar lowS = params_.krLowS();
226 const Scalar highS = params_.krHighS();
227
228 // check whether the saturation is unpyhsical
229 if (S >= 1.0)
230 return 1.0;
231 else if (S <= 0.0)
232 return 0;
233 // check wether the permeability needs to be regularized
234 else if (S < lowS)
235 return splineLowS_.eval(S);
236 else if (S > highS)
237 return splineHighS_.eval(S);
238
239 // straight line for S \in [lowS, highS]
240 return lowS/2.0 + splineM_*(S - lowS);
241 }
242
243 Params params_;
244 EffToAbsParams effToAbsParams_;
245 Scalar splineM_;
246 Dumux::Spline<Scalar> splineLowS_;
247 Dumux::Spline<Scalar> splineHighS_;
248};
249} // end namespace Dumux::FluidMatrix
250
251#endif
Provides 3rd order polynomial splines.
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.
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Definition: brookscorey.hh:35
A 3rd order polynomial spline.
Definition: spline.hh:55
Implements a linear saturation-capillary pressure relation.
Definition: smoothedlinearlaw.hh:56
Scalar krw(Scalar swe) const
The relative permeability for the wetting phase.
Definition: smoothedlinearlaw.hh:195
EffToAbsPolicy EffToAbs
Definition: smoothedlinearlaw.hh:61
Scalar endPointPc() const
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: smoothedlinearlaw.hh:173
SmoothedLinearLaw(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: smoothedlinearlaw.hh:119
SmoothedLinearLaw()=delete
Deleted default constructor (so we are never in an undefined state)
static Params makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: smoothedlinearlaw.hh:145
SmoothedLinearLaw(const Params &params, const EffToAbsParams &effToAbsParams={})
Construct from parameter structs.
Definition: smoothedlinearlaw.hh:127
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: smoothedlinearlaw.hh:72
Scalar swe(Scalar pc) const
The inverse saturation-capillary pressure curve.
Definition: smoothedlinearlaw.hh:165
typename EffToAbsPolicy::template Params< Scalar > EffToAbsParams
Definition: smoothedlinearlaw.hh:60
const EffToAbsParams & effToAbsParams() const
Definition: smoothedlinearlaw.hh:218
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: smoothedlinearlaw.hh:66
bool operator==(const SmoothedLinearLaw< Scalar, EffToAbs > &o) const
Equality comparison with another instance.
Definition: smoothedlinearlaw.hh:212
Scalar pc(Scalar swe) const
The capillary pressure-saturation curve.
Definition: smoothedlinearlaw.hh:157
ScalarType Scalar
Definition: smoothedlinearlaw.hh:59
Scalar dswe_dpc(Scalar pc) const
The partial derivative of the effective saturation w.r.t. the capillary pressure.
Definition: smoothedlinearlaw.hh:187
Scalar dpc_dswe(Scalar swe) const
The partial derivative of the capillary pressure w.r.t. the effective saturation.
Definition: smoothedlinearlaw.hh:179
Scalar krn(Scalar swe) const
The relative permeability for the non-wetting phase.
Definition: smoothedlinearlaw.hh:203
The parameter type.
Definition: smoothedlinearlaw.hh:80
void setKrLowS(Scalar krLowS)
Definition: smoothedlinearlaw.hh:92
void setPe(Scalar pe)
Definition: smoothedlinearlaw.hh:86
Scalar krHighS() const
Definition: smoothedlinearlaw.hh:94
Scalar pcMax() const
Definition: smoothedlinearlaw.hh:88
Scalar krLowS() const
Definition: smoothedlinearlaw.hh:91
Scalar pe() const
Definition: smoothedlinearlaw.hh:85
Params(Scalar pe, Scalar pcMax, Scalar krLowS, Scalar krHighS)
Definition: smoothedlinearlaw.hh:81
void setPcMax(Scalar pcMax)
Definition: smoothedlinearlaw.hh:89
void setKrHighS(Scalar krHighS)
Definition: smoothedlinearlaw.hh:95
bool operator==(const Params &p) const
Definition: smoothedlinearlaw.hh:97
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:67