version 3.8
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// 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_TWOP_SMOOTHED_LINEAR_LAW_HH
14#define DUMUX_MATERIAL_FLUIDMATRIX_TWOP_SMOOTHED_LINEAR_LAW_HH
15
16
17#include <cmath>
18#include <algorithm>
19
24
25namespace Dumux::FluidMatrix {
26
42template<class ScalarType, class EffToAbsPolicy = TwoPEffToAbsDefaultPolicy>
43class SmoothedLinearLaw : public Adapter<SmoothedLinearLaw<ScalarType, EffToAbsPolicy>, PcKrSw>
44{
45
46public:
47 using Scalar = ScalarType;
48 using EffToAbsParams = typename EffToAbsPolicy::template Params<Scalar>;
49 using EffToAbs = EffToAbsPolicy;
50
54 static constexpr bool isRegularized()
55 { return false; }
56
60 static constexpr int numFluidPhases()
61 { return 2; }
62
67 struct Params
68 {
70 : pe_(pe), pcMax_(pcMax), krLowS_(krLowS), krHighS_(krHighS)
71 {}
72
73 Scalar pe() const { return pe_; }
74 void setPe(Scalar pe) { pe_ = pe; }
75
76 Scalar pcMax() const { return pcMax_; }
77 void setPcMax(Scalar pcMax) { pcMax_ = pcMax; }
78
79 Scalar krLowS() const { return krLowS_; }
80 void setKrLowS(Scalar krLowS) { krLowS_ = krLowS; }
81
82 Scalar krHighS() const { return krHighS_; }
83 void setKrHighS(Scalar krHighS) { krHighS_ = krHighS; }
84
85 bool operator== (const Params& p) const
86 {
87 return Dune::FloatCmp::eq(pe(), p.pe(), 1e-6)
88 && Dune::FloatCmp::eq(pcMax(), p.pcMax(), 1e-6)
89 && Dune::FloatCmp::eq(krLowS(), p.krLowS(), 1e-6)
90 && Dune::FloatCmp::eq(krHighS(), p.krHighS(), 1e-6);
91 }
92
93 private:
94 Scalar pe_, pcMax_, krLowS_, krHighS_;
95 };
96
102
107 explicit SmoothedLinearLaw(const std::string& paramGroup)
108 : SmoothedLinearLaw(makeParams(paramGroup), EffToAbs::template makeParams<Scalar>(paramGroup))
109 {}
110
116 const EffToAbsParams& effToAbsParams = {})
117 : params_(params)
118 , effToAbsParams_(effToAbsParams)
119 , splineM_((1.0 - ((1.0 - params_.krHighS()) + params_.krLowS())/2.0 )
120 / (1.0 - (1.0 - params_.krHighS()) - params_.krLowS()))
121 , splineLowS_(0.0, params_.krLowS(), // x1, x2
122 0.0, params_.krLowS()/2.0, // y1, y2
123 0.0, splineM_) // m1, m2
124 , splineHighS_(params_.krHighS(), 1.0, // x1, x2
125 1.0 - (1.0 - params_.krHighS())/2.0, 1.0, // y1, y2
126 splineM_, 0.0) // m1, m2
127 {}
128
133 static Params makeParams(const std::string& paramGroup)
134 {
135 const auto pe = getParamFromGroup<Scalar>(paramGroup, "SmoothedLinearLawPe");
136 const auto pcMax = getParamFromGroup<Scalar>(paramGroup, "SmoothedLinearLawPcMax");
137 const auto krLowS = getParamFromGroup<Scalar>(paramGroup, "SmoothedLinearLawKrLowS");
138 const auto krHighS = getParamFromGroup<Scalar>(paramGroup, "SmoothedLinearLawKrHighS");
139 return {pe, pcMax, krLowS, krHighS};
140 }
141
146 {
147 return (1.0 - swe)*(params_.pcMax() - params_.pe()) + params_.pe();
148 }
149
154 {
155 return 1.0 - (pc - params_.pe())/(params_.pcMax() - params_.pe());
156 }
157
162 { return params_.pe(); }
163
168 {
169 return params_.pe() - params_.pcMax();
170 }
171
176 {
177 return 1.0/(params_.pe() - params_.pcMax());
178 }
179
184 {
185 return relperm_(swe);
186 }
187
192 {
193 Scalar sne = 1.0 - swe;
194 return relperm_(sne);
195 }
196
201 {
202 return params_ == o.params_
203 && effToAbsParams_ == o.effToAbsParams_;
204 }
205
207 { return effToAbsParams_; }
208
209private:
210
211 Scalar relperm_(Scalar S) const
212 {
213 const Scalar lowS = params_.krLowS();
214 const Scalar highS = params_.krHighS();
215
216 // check whether the saturation is unpyhsical
217 if (S >= 1.0)
218 return 1.0;
219 else if (S <= 0.0)
220 return 0;
221 // check whether the permeability needs to be regularized
222 else if (S < lowS)
223 return splineLowS_.eval(S);
224 else if (S > highS)
225 return splineHighS_.eval(S);
226
227 // straight line for S \in [lowS, highS]
228 return lowS/2.0 + splineM_*(S - lowS);
229 }
230
231 Params params_;
232 EffToAbsParams effToAbsParams_;
233 Scalar splineM_;
234 Dumux::Spline<Scalar> splineLowS_;
235 Dumux::Spline<Scalar> splineHighS_;
236};
237} // end namespace Dumux::FluidMatrix
238
239#endif
Implements a linear saturation-capillary pressure relation.
Definition: smoothedlinearlaw.hh:44
Scalar krw(Scalar swe) const
The relative permeability for the wetting phase.
Definition: smoothedlinearlaw.hh:183
EffToAbsPolicy EffToAbs
Definition: smoothedlinearlaw.hh:49
Scalar endPointPc() const
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: smoothedlinearlaw.hh:161
SmoothedLinearLaw(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: smoothedlinearlaw.hh:107
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:133
SmoothedLinearLaw(const Params &params, const EffToAbsParams &effToAbsParams={})
Construct from parameter structs.
Definition: smoothedlinearlaw.hh:115
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: smoothedlinearlaw.hh:60
Scalar swe(Scalar pc) const
The inverse saturation-capillary pressure curve.
Definition: smoothedlinearlaw.hh:153
typename EffToAbsPolicy::template Params< Scalar > EffToAbsParams
Definition: smoothedlinearlaw.hh:48
const EffToAbsParams & effToAbsParams() const
Definition: smoothedlinearlaw.hh:206
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: smoothedlinearlaw.hh:54
bool operator==(const SmoothedLinearLaw< Scalar, EffToAbs > &o) const
Equality comparison with another instance.
Definition: smoothedlinearlaw.hh:200
Scalar pc(Scalar swe) const
The capillary pressure-saturation curve.
Definition: smoothedlinearlaw.hh:145
ScalarType Scalar
Definition: smoothedlinearlaw.hh:47
Scalar dswe_dpc(Scalar pc) const
The partial derivative of the effective saturation w.r.t. the capillary pressure.
Definition: smoothedlinearlaw.hh:175
Scalar dpc_dswe(Scalar swe) const
The partial derivative of the capillary pressure w.r.t. the effective saturation.
Definition: smoothedlinearlaw.hh:167
Scalar krn(Scalar swe) const
The relative permeability for the non-wetting phase.
Definition: smoothedlinearlaw.hh:191
A 3rd order polynomial spline.
Definition: spline.hh:43
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:23
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Provides 3rd order polynomial splines.
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:55
The parameter type.
Definition: smoothedlinearlaw.hh:68
void setKrLowS(Scalar krLowS)
Definition: smoothedlinearlaw.hh:80
void setPe(Scalar pe)
Definition: smoothedlinearlaw.hh:74
Scalar krHighS() const
Definition: smoothedlinearlaw.hh:82
Scalar pcMax() const
Definition: smoothedlinearlaw.hh:76
Scalar krLowS() const
Definition: smoothedlinearlaw.hh:79
Scalar pe() const
Definition: smoothedlinearlaw.hh:73
Params(Scalar pe, Scalar pcMax, Scalar krLowS, Scalar krHighS)
Definition: smoothedlinearlaw.hh:69
void setPcMax(Scalar pcMax)
Definition: smoothedlinearlaw.hh:77
void setKrHighS(Scalar krHighS)
Definition: smoothedlinearlaw.hh:83
bool operator==(const Params &p) const
Definition: smoothedlinearlaw.hh:85