version 3.8
splinemateriallaw.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//
12#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWOP_SPLINE_MATERIAL_LAW_HH
13#define DUMUX_MATERIAL_FLUIDMATRIX_TWOP_SPLINE_MATERIAL_LAW_HH
14
15#include <memory> // unique_ptr
16
21
22namespace Dumux::FluidMatrix {
23
34template<class TwoPMaterialLaw, bool approximatePcSwInverse = false>
36: public TwoPMaterialLaw
37, public Adapter<SplineTwoPMaterialLaw<TwoPMaterialLaw, approximatePcSwInverse>, PcKrSw>
38{
39public:
41
45
49 static constexpr int numFluidPhases()
50 { return 2; }
51
56 static constexpr bool isRegularized()
57 { return true; }
58
64
69 explicit SplineTwoPMaterialLaw(const std::string& paramGroup)
70 : TwoPMaterialLaw(paramGroup)
71 {
72 const std::array<Scalar, 2> defaultInterval{{ 0.01, 1.0 }};
73 const auto sweInterval = getParamFromGroup<std::array<Scalar, 2>>(paramGroup, "SplineSweInterval", defaultInterval);
74 swInterval_ = {{ TwoPMaterialLaw::EffToAbs::sweToSw(sweInterval[0], this->effToAbsParams()),
75 TwoPMaterialLaw::EffToAbs::sweToSw(sweInterval[1], this->effToAbsParams()) }};
76 swIntervalPc_ = {{ TwoPMaterialLaw::pc(swInterval_[1]),
77 TwoPMaterialLaw::pc(swInterval_[0]) }};
78 numSwSamples_ = getParamFromGroup<std::size_t>(paramGroup, "SplineNumSwSamples", 30);
79
80 pcSpline_ = makeSweSpline_(
81 [&](const Scalar s){ return TwoPMaterialLaw::pc(s); },
82 approximatePcSwInverse
83 );
84
85 krwSpline_ = makeSweSpline_(
86 [&](const Scalar s){ return TwoPMaterialLaw::krw(s); }
87 );
88
89 krnSpline_ = makeSweSpline_(
90 [&](const Scalar s){ return TwoPMaterialLaw::krn(s); }
91 );
92 }
93
98 SplineTwoPMaterialLaw(const std::array<Scalar, 2>& sweInterval,
99 std::size_t numSwSamples,
100 TwoPMaterialLaw&& twoP)
101 : TwoPMaterialLaw(std::move(twoP))
102 , numSwSamples_(numSwSamples)
103 {
104 swInterval_ = {{ TwoPMaterialLaw::EffToAbs::sweToSw(sweInterval[0], this->effToAbsParams()),
105 TwoPMaterialLaw::EffToAbs::sweToSw(sweInterval[1], this->effToAbsParams()) }};
106 swIntervalPc_ = {{ TwoPMaterialLaw::pc(swInterval_[1]),
107 TwoPMaterialLaw::pc(swInterval_[0]) }};
108 }
109
113 Scalar pc(const Scalar sw) const
114 {
115 if (sw > swInterval_[0] && sw < swInterval_[1])
116 {
117 if constexpr (approximatePcSwInverse)
118 return pcSpline_->evalInverse(sw);
119 else
120 return pcSpline_->eval(sw);
121 }
122
123 return TwoPMaterialLaw::pc(sw);
124 }
125
129 Scalar dpc_dsw(const Scalar sw) const
130 {
131 if (sw > swInterval_[0] && sw < swInterval_[1])
132 {
133 if constexpr (approximatePcSwInverse)
134 return 1.0/pcSpline_->evalDerivative(pcSpline_->evalInverse(sw));
135 else
136 return pcSpline_->evalDerivative(sw);
137 }
138
140 }
141
145 Scalar sw(const Scalar pc) const
146 {
147 if (pc > swIntervalPc_[0] && pc < swIntervalPc_[1])
148 {
149 if constexpr (approximatePcSwInverse)
150 return pcSpline_->eval(pc);
151 else
152 return pcSpline_->evalInverse(pc);
153 }
154
155 return TwoPMaterialLaw::sw(pc);
156 }
157
161 Scalar dsw_dpc(const Scalar pc) const
162 {
163 if (pc > swIntervalPc_[0] && pc < swIntervalPc_[1])
164 {
165 if constexpr (approximatePcSwInverse)
166 return pcSpline_->evalDerivative(pc);
167 else
168 return 1.0/pcSpline_->evalDerivative(pcSpline_->evalInverse(pc));
169 }
170
172 }
173
177 Scalar krw(const Scalar sw) const
178 {
179 if (sw > swInterval_[0] && sw < swInterval_[1])
180 return krwSpline_->eval(sw);
181
182 return TwoPMaterialLaw::krw(sw);
183 }
184
189 {
190 if (sw > swInterval_[0] && sw < swInterval_[1])
191 return krwSpline_->evalDerivative(sw);
192
194 }
195
199 Scalar krn(const Scalar sw) const
200 {
201 if (sw > swInterval_[0] && sw < swInterval_[1])
202 return krnSpline_->eval(sw);
203
204 return TwoPMaterialLaw::krn(sw);
205 }
206
211 {
212 if (sw > swInterval_[0] && sw < swInterval_[1])
213 return krnSpline_->evalDerivative(sw);
214
216 }
217
218private:
219 template<class Function>
220 std::unique_ptr<MonotoneCubicSpline<Scalar>>
221 makeSweSpline_(const Function& f, bool invert = false) const
222 {
223 const auto sw = linspace(swInterval_[0], swInterval_[1], numSwSamples_);
224
225 auto values = sw;
226 std::transform(sw.begin(), sw.end(), values.begin(), f);
227
228 if (invert)
229 return std::make_unique<MonotoneCubicSpline<Scalar>>(values, sw);
230 else
231 return std::make_unique<MonotoneCubicSpline<Scalar>>(sw, values);
232 }
233
234 std::unique_ptr<MonotoneCubicSpline<Scalar>> pcSpline_;
235 std::unique_ptr<MonotoneCubicSpline<Scalar>> krwSpline_;
236 std::unique_ptr<MonotoneCubicSpline<Scalar>> krnSpline_;
237
238 std::array<Scalar, 2> swInterval_;
239 std::array<Scalar, 2> swIntervalPc_;
240 std::size_t numSwSamples_;
241};
242
243} // end namespace Dumux::FluidMatrix
244
245#endif
A spline approximation wrapper for 2p material laws.
Definition: splinemateriallaw.hh:38
static constexpr bool isRegularized()
We are always regularized in the sense that we replace the original curve by a cubic spline.
Definition: splinemateriallaw.hh:56
Scalar dkrn_dsw(const Scalar sw) const
The derivative of the relative permeability for the non-wetting phase w.r.t. saturation.
Definition: splinemateriallaw.hh:210
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: splinemateriallaw.hh:113
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: splinemateriallaw.hh:49
typename TwoPMaterialLaw::Scalar Scalar
Definition: splinemateriallaw.hh:40
SplineTwoPMaterialLaw(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: splinemateriallaw.hh:69
SplineTwoPMaterialLaw(const std::array< Scalar, 2 > &sweInterval, std::size_t numSwSamples, TwoPMaterialLaw &&twoP)
Construct from parameter structs.
Definition: splinemateriallaw.hh:98
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: splinemateriallaw.hh:161
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition: splinemateriallaw.hh:199
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition: splinemateriallaw.hh:177
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition: splinemateriallaw.hh:188
SplineTwoPMaterialLaw()=delete
Deleted default constructor (so we are never in an undefined state)
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: splinemateriallaw.hh:129
Scalar sw(const Scalar pc) const
The saturation-capillary pressure curve.
Definition: splinemateriallaw.hh:145
Wrapper class to implement regularized material laws (pc-sw, kr-sw) with a conversion policy between ...
Definition: materiallaw.hh:45
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: materiallaw.hh:162
Scalar dkrn_dsw(const Scalar sw) const
The derivative of the relative permeability for the non-wetting phase w.r.t. saturation.
Definition: materiallaw.hh:229
ScalarType Scalar
Definition: materiallaw.hh:48
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition: materiallaw.hh:212
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition: materiallaw.hh:195
typename Regularization::template Params< Scalar > RegularizationParams
Definition: materiallaw.hh:52
const EffToAbsParams & effToAbsParams() const
Return the parameters of the EffToAbs policy.
Definition: materiallaw.hh:279
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: materiallaw.hh:121
Scalar sw(const Scalar pc) const
The saturation-capillary pressure curve.
Definition: materiallaw.hh:146
typename EffToAbsPolicy::template Params< Scalar > EffToAbsParams
Definition: materiallaw.hh:51
typename BaseLaw::template Params< Scalar > BasicParams
Definition: materiallaw.hh:50
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: materiallaw.hh:104
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition: materiallaw.hh:178
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
std::vector< Scalar > linspace(const Scalar begin, const Scalar end, std::size_t samples, bool endPoint=true)
Generates linearly spaced vectors.
Definition: math.hh:582
Implementation helper for capillary pressure and relative permeability <-> saturation relations for t...
A monotone cubic spline.
Definition: brookscorey.hh:23
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:55