version 3.8
plotpckrsw.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_IO_PLOT_PC_KR_SW_HH
13#define DUMUX_IO_PLOT_PC_KR_SW_HH
14
15#include <string>
16#include <tuple>
17#include <algorithm>
18#include <dumux/common/math.hh>
19
20namespace Dumux {
21
22namespace Detail {
23template<class Function, class Range>
24Range evalFunctionForRange(const Function& f, const Range& range)
25{
26 Range result = range;
27 std::transform(range.begin(), range.end(), result.begin(), [&](auto x){ return f(x); });
28 return result;
29}
30} // end namespace Detail
31
36template<class PcKrSw, class V>
37auto samplePcSw(const PcKrSw& curve, const V& sw)
38{
39 return Detail::evalFunctionForRange([&](const auto s){ return curve.pc(s); }, sw);
40}
41
46template<class PcKrSw, class V>
47auto samplePcSwDerivative(const PcKrSw& curve, const V& sw)
48{
49 return Detail::evalFunctionForRange([&](const auto s){ return curve.dpc_dsw(s); }, sw);
50}
51
56template<class PcKrSw, class V>
57auto samplePcSwInverseDerivative(const PcKrSw& curve, const V& pc)
58{
59 return Detail::evalFunctionForRange([&](const auto p){ return curve.dsw_dpc(p); }, pc);
60}
61
66template<class PcKrSw, class V>
67auto sampleLog10PcSw(const PcKrSw& curve, const V& sw)
68{
69 return Detail::evalFunctionForRange([&](const auto s){ using std::log10; return log10(curve.pc(s)); }, sw);
70}
71
76template<class PcKrSw, class V>
77auto sampleRelPerms(const PcKrSw& curve, const V& sw)
78{
79 return std::make_tuple(
80 Detail::evalFunctionForRange([&](const auto s){ return curve.krw(s); }, sw),
81 Detail::evalFunctionForRange([&](const auto s){ return curve.krn(s); }, sw)
82 );
83}
84
89template<class PcKrSw, class V>
90auto sampleRelPermDerivatives(const PcKrSw& curve, const V& sw)
91{
92 return std::make_tuple(
93 Detail::evalFunctionForRange([&](const auto s){ return curve.dkrw_dsw(s); }, sw),
94 Detail::evalFunctionForRange([&](const auto s){ return curve.dkrn_dsw(s); }, sw)
95 );
96}
97
98// forward declaration
99template<class S> class GnuplotInterface;
100
101namespace Detail {
102template<class S, class V>
104 const V& x, const V& y,
105 const std::string& curveName,
106 const std::string& curveOptions,
107 const std::string& xLabel,
108 const std::string& yLabel)
109{
110 gnuplot.setXlabel(xLabel);
111 gnuplot.setYlabel(yLabel);
112 gnuplot.addDataSetToPlot(x, y, curveName, curveOptions);
113}
114} // end namespace Detail
115
120namespace Gnuplot {
121
126template<class S, class V>
127void addPcSw(GnuplotInterface<S>& gnuplot, const V& sw, const V& pc,
128 const std::string& curveName = "pc-sw",
129 const std::string& curveOptions = "w l",
130 const std::string& xLabel = "wetting phase saturation [-]",
131 const std::string& yLabel = "capillary pressure [Pa]")
132{
133 Detail::addDataSetToGnuplot(gnuplot, sw, pc, curveName, curveOptions, xLabel, yLabel);
134}
135
140template<class S, class V>
141void addPcSwDerivative(GnuplotInterface<S>& gnuplot, const V& sw, const V& dpc_dsw,
142 const std::string& curveName = "dpc-dsw",
143 const std::string& curveOptions = "w l",
144 const std::string& xLabel = "wetting phase saturation [-]",
145 const std::string& yLabel = "derivative of capillary pressure [Pa]")
146{
147 Detail::addDataSetToGnuplot(gnuplot, sw, dpc_dsw, curveName, curveOptions, xLabel, yLabel);
148}
149
154template<class S, class V>
155void addPcSwInverseDerivative(GnuplotInterface<S>& gnuplot, const V& sw, const V& dpc_dsw,
156 const std::string& curveName = "dsw-dpc",
157 const std::string& curveOptions = "w l",
158 const std::string& xLabel = "capillary pressure [Pa]",
159 const std::string& yLabel = "derivative of saturation [1/Pa]")
160{
161 Detail::addDataSetToGnuplot(gnuplot, sw, dpc_dsw, curveName, curveOptions, xLabel, yLabel);
162}
163
168template<class S, class V>
169void addLog10PcSw(GnuplotInterface<S>& gnuplot, const V& sw, const V& log10pc,
170 const std::string& curveName = "log10-pc-sw",
171 const std::string& curveOptions = "w l",
172 const std::string& xLabel = "wetting phase saturation [-]",
173 const std::string& yLabel = "log10 of capillary pressure [Pa]")
174{
175 Detail::addDataSetToGnuplot(gnuplot, sw, log10pc, curveName, curveOptions, xLabel, yLabel);
176}
177
182template<class S, class V>
183void addRelPerms(GnuplotInterface<S>& gnuplot, const V& sw, const V& krw, const V& krn,
184 const std::string& curveName = "relperm",
185 const std::string& curveOptions = "w l",
186 const std::string& xLabel = "wetting phase saturation [-]",
187 const std::string& yLabel = "relative permeability [-]")
188{
189 Detail::addDataSetToGnuplot(gnuplot, sw, krw, curveName + "_krw", curveOptions, xLabel, yLabel);
190 Detail::addDataSetToGnuplot(gnuplot, sw, krn, curveName + "_krn", curveOptions, xLabel, yLabel);
191}
192
197template<class S, class V>
198void addRelPermDerivatives(GnuplotInterface<S>& gnuplot, const V& sw, const V& krw, const V& krn,
199 const std::string& curveName = "relperm_dsw",
200 const std::string& curveOptions = "w l",
201 const std::string& xLabel = "wetting phase saturation [-]",
202 const std::string& yLabel = "derivative of the relative permeability [-]")
203{
204 Detail::addDataSetToGnuplot(gnuplot, sw, krw, curveName + "_krw", curveOptions, xLabel, yLabel);
205 Detail::addDataSetToGnuplot(gnuplot, sw, krn, curveName + "_krn", curveOptions, xLabel, yLabel);
206}
207
208} // end namespace Gnuplot
209} // end namespace Dumux
210
211#endif
Interface for passing data sets to a file and plotting them, if gnuplot is installed.
Definition: gnuplotinterface.hh:45
void setXlabel(const std::string &label)
Sets the label for the x-axis.
Definition: gnuplotinterface.hh:267
void addDataSetToPlot(const std::vector< Scalar > &x, const std::vector< Scalar > &y, const std::string &fileName, const std::string &options="with lines")
Adds a data set and writes a data file.
Definition: gnuplotinterface.hh:230
void setYlabel(const std::string &label)
Sets the label for the y-axis.
Definition: gnuplotinterface.hh:277
PcKrSw(T &&) -> PcKrSw< T >
Deduction guide for the PcKrSw class. Makes sure that PcKrSw stores a copy of T if the constructor is...
void addPcSw(GnuplotInterface< S > &gnuplot, const V &sw, const V &pc, const std::string &curveName="pc-sw", const std::string &curveOptions="w l", const std::string &xLabel="wetting phase saturation [-]", const std::string &yLabel="capillary pressure [Pa]")
Convenience function for adding material law quantities to gnuplot.
Definition: plotpckrsw.hh:127
auto samplePcSw(const PcKrSw &curve, const V &sw)
sample the pc-sw curve
Definition: plotpckrsw.hh:37
auto samplePcSwDerivative(const PcKrSw &curve, const V &sw)
sample the pc-sw curve derivative wrt sw
Definition: plotpckrsw.hh:47
void addPcSwInverseDerivative(GnuplotInterface< S > &gnuplot, const V &sw, const V &dpc_dsw, const std::string &curveName="dsw-dpc", const std::string &curveOptions="w l", const std::string &xLabel="capillary pressure [Pa]", const std::string &yLabel="derivative of saturation [1/Pa]")
Convenience function for adding material law quantities to gnuplot.
Definition: plotpckrsw.hh:155
auto sampleRelPerms(const PcKrSw &curve, const V &sw)
sample krw-sw and krn-sw curves
Definition: plotpckrsw.hh:77
void addRelPerms(GnuplotInterface< S > &gnuplot, const V &sw, const V &krw, const V &krn, const std::string &curveName="relperm", const std::string &curveOptions="w l", const std::string &xLabel="wetting phase saturation [-]", const std::string &yLabel="relative permeability [-]")
Convenience function for adding material law quantities to gnuplot.
Definition: plotpckrsw.hh:183
auto samplePcSwInverseDerivative(const PcKrSw &curve, const V &pc)
sample the sw-pc curve derivative wrt pc
Definition: plotpckrsw.hh:57
void addLog10PcSw(GnuplotInterface< S > &gnuplot, const V &sw, const V &log10pc, const std::string &curveName="log10-pc-sw", const std::string &curveOptions="w l", const std::string &xLabel="wetting phase saturation [-]", const std::string &yLabel="log10 of capillary pressure [Pa]")
Convenience function for adding material law quantities to gnuplot.
Definition: plotpckrsw.hh:169
auto sampleRelPermDerivatives(const PcKrSw &curve, const V &sw)
sample the derivatives of the krw-sw and krn-sw curves
Definition: plotpckrsw.hh:90
auto sampleLog10PcSw(const PcKrSw &curve, const V &sw)
sample sw-pc curve but return the log10 of the capillary pressure
Definition: plotpckrsw.hh:67
void addRelPermDerivatives(GnuplotInterface< S > &gnuplot, const V &sw, const V &krw, const V &krn, const std::string &curveName="relperm_dsw", const std::string &curveOptions="w l", const std::string &xLabel="wetting phase saturation [-]", const std::string &yLabel="derivative of the relative permeability [-]")
Convenience function for adding material law quantities to gnuplot.
Definition: plotpckrsw.hh:198
void addPcSwDerivative(GnuplotInterface< S > &gnuplot, const V &sw, const V &dpc_dsw, const std::string &curveName="dpc-dsw", const std::string &curveOptions="w l", const std::string &xLabel="wetting phase saturation [-]", const std::string &yLabel="derivative of capillary pressure [Pa]")
Convenience function for adding material law quantities to gnuplot.
Definition: plotpckrsw.hh:141
Define some often used mathematical functions.
Range evalFunctionForRange(const Function &f, const Range &range)
Definition: plotpckrsw.hh:24
void addDataSetToGnuplot(GnuplotInterface< S > &gnuplot, const V &x, const V &y, const std::string &curveName, const std::string &curveOptions, const std::string &xLabel, const std::string &yLabel)
Definition: plotpckrsw.hh:103
Definition: adapt.hh:17