3.3.0
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
73 struct Params
74 {
76 : pe_(pe), pcMax_(pcMax), krLowS_(krLowS), krHighS_(krHighS)
77 {}
78
79 Scalar pe() const { return pe_; }
80 void setPe(Scalar pe) { pe_ = pe; }
81
82 Scalar pcMax() const { return pcMax_; }
83 void setPcMax(Scalar pcMax) { pcMax_ = pcMax; }
84
85 Scalar krLowS() const { return krLowS_; }
86 void setKrLowS(Scalar krLowS) { krLowS_ = krLowS; }
87
88 Scalar krHighS() const { return krHighS_; }
89 void setKrHighS(Scalar krHighS) { krHighS_ = krHighS; }
90
91 bool operator== (const Params& p) const
92 {
93 return Dune::FloatCmp::eq(pe(), p.pe(), 1e-6)
94 && Dune::FloatCmp::eq(pcMax(), p.pcMax(), 1e-6)
95 && Dune::FloatCmp::eq(krLowS(), p.krLowS(), 1e-6)
96 && Dune::FloatCmp::eq(krHighS(), p.krHighS(), 1e-6);
97 }
98
99 private:
100 Scalar pe_, pcMax_, krLowS_, krHighS_;
101 };
102
108
113 explicit SmoothedLinearLaw(const std::string& paramGroup)
114 : SmoothedLinearLaw(makeParams(paramGroup), EffToAbs::template makeParams<Scalar>(paramGroup))
115 {}
116
122 const EffToAbsParams& effToAbsParams = {})
123 : params_(params)
124 , effToAbsParams_(effToAbsParams)
125 , splineM_((1.0 - ((1.0 - params_.krHighS()) + params_.krLowS())/2.0 )
126 / (1.0 - (1.0 - params_.krHighS()) - params_.krLowS()))
127 , splineLowS_(0.0, params_.krLowS(), // x1, x2
128 0.0, params_.krLowS()/2.0, // y1, y2
129 0.0, splineM_) // m1, m2
130 , splineHighS_(params_.krHighS(), 1.0, // x1, x2
131 1.0 - (1.0 - params_.krHighS())/2.0, 1.0, // y1, y2
132 splineM_, 0.0) // m1, m2
133 {}
134
139 static Params makeParams(const std::string& paramGroup)
140 {
141 const auto pe = getParamFromGroup<Scalar>(paramGroup, "SmoothedLinearLawPe");
142 const auto pcMax = getParamFromGroup<Scalar>(paramGroup, "SmoothedLinearLawPcMax");
143 const auto krLowS = getParamFromGroup<Scalar>(paramGroup, "SmoothedLinearLawKrLowS");
144 const auto krHighS = getParamFromGroup<Scalar>(paramGroup, "SmoothedLinearLawKrHighS");
145 return {pe, pcMax, krLowS, krHighS};
146 }
147
152 {
153 return (1.0 - swe)*(params_.pcMax() - params_.pe()) + params_.pe();
154 }
155
160 {
161 return 1.0 - (pc - params_.pe())/(params_.pcMax() - params_.pe());
162 }
163
168 { return params_.pe(); }
169
174 {
175 return params_.pe() - params_.pcMax();
176 }
177
182 {
183 return 1.0/(params_.pe() - params_.pcMax());
184 }
185
190 {
191 return relperm_(swe);
192 }
193
198 {
199 Scalar sne = 1.0 - swe;
200 return relperm_(sne);
201 }
202
207 {
208 return params_ == o.params_
209 && effToAbsParams_ == o.effToAbsParams_;
210 }
211
213 { return effToAbsParams_; }
214
215private:
216
217 Scalar relperm_(Scalar S) const
218 {
219 const Scalar lowS = params_.krLowS();
220 const Scalar highS = params_.krHighS();
221
222 // check whether the saturation is unpyhsical
223 if (S >= 1.0)
224 return 1.0;
225 else if (S <= 0.0)
226 return 0;
227 // check wether the permeability needs to be regularized
228 else if (S < lowS)
229 return splineLowS_.eval(S);
230 else if (S > highS)
231 return splineHighS_.eval(S);
232
233 // straight line for S \in [lowS, highS]
234 return lowS/2.0 + splineM_*(S - lowS);
235 }
236
237 Params params_;
238 EffToAbsParams effToAbsParams_;
239 Scalar splineM_;
240 Dumux::Spline<Scalar> splineLowS_;
241 Dumux::Spline<Scalar> splineHighS_;
242};
243} // end namespace Dumux::FluidMatrix
244
245#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:286
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:189
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:167
SmoothedLinearLaw(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: smoothedlinearlaw.hh:113
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:139
SmoothedLinearLaw(const Params &params, const EffToAbsParams &effToAbsParams={})
Construct from parameter structs.
Definition: smoothedlinearlaw.hh:121
Scalar swe(Scalar pc) const
The inverse saturation-capillary pressure curve.
Definition: smoothedlinearlaw.hh:159
typename EffToAbsPolicy::template Params< Scalar > EffToAbsParams
Definition: smoothedlinearlaw.hh:60
const EffToAbsParams & effToAbsParams() const
Definition: smoothedlinearlaw.hh:212
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:206
Scalar pc(Scalar swe) const
The capillary pressure-saturation curve.
Definition: smoothedlinearlaw.hh:151
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:181
Scalar dpc_dswe(Scalar swe) const
The partial derivative of the capillary pressure w.r.t. the effective saturation.
Definition: smoothedlinearlaw.hh:173
Scalar krn(Scalar swe) const
The relative permeability for the non-wetting phase.
Definition: smoothedlinearlaw.hh:197
The parameter type.
Definition: smoothedlinearlaw.hh:74
void setKrLowS(Scalar krLowS)
Definition: smoothedlinearlaw.hh:86
void setPe(Scalar pe)
Definition: smoothedlinearlaw.hh:80
Scalar krHighS() const
Definition: smoothedlinearlaw.hh:88
Scalar pcMax() const
Definition: smoothedlinearlaw.hh:82
Scalar krLowS() const
Definition: smoothedlinearlaw.hh:85
Scalar pe() const
Definition: smoothedlinearlaw.hh:79
Params(Scalar pe, Scalar pcMax, Scalar krLowS, Scalar krHighS)
Definition: smoothedlinearlaw.hh:75
void setPcMax(Scalar pcMax)
Definition: smoothedlinearlaw.hh:83
void setKrHighS(Scalar krHighS)
Definition: smoothedlinearlaw.hh:89
bool operator==(const Params &p) const
Definition: smoothedlinearlaw.hh:91
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:67