3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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 *****************************************************************************/
24#ifndef DUMUX_IO_PLOT_PC_KR_SW_HH
25#define DUMUX_IO_PLOT_PC_KR_SW_HH
26
27#include <string>
28#include <tuple>
29#include <algorithm>
30#include <dumux/common/math.hh>
31
32namespace Dumux {
33
34namespace Detail {
35template<class Function, class Range>
36Range evalFunctionForRange(const Function& f, const Range& range)
37{
38 Range result = range;
39 std::transform(range.begin(), range.end(), result.begin(), [&](auto x){ return f(x); });
40 return result;
41}
42} // end namespace Detail
43
48template<class PcKrSw, class V>
49auto samplePcSw(const PcKrSw& curve, const V& sw)
50{
51 return Detail::evalFunctionForRange([&](const auto s){ return curve.pc(s); }, sw);
52}
53
58template<class PcKrSw, class V>
59auto samplePcSwDerivative(const PcKrSw& curve, const V& sw)
60{
61 return Detail::evalFunctionForRange([&](const auto s){ return curve.dpc_dsw(s); }, sw);
62}
63
68template<class PcKrSw, class V>
69auto samplePcSwInverseDerivative(const PcKrSw& curve, const V& pc)
70{
71 return Detail::evalFunctionForRange([&](const auto p){ return curve.dsw_dpc(p); }, pc);
72}
73
78template<class PcKrSw, class V>
79auto sampleLog10PcSw(const PcKrSw& curve, const V& sw)
80{
81 return Detail::evalFunctionForRange([&](const auto s){ using std::log10; return log10(curve.pc(s)); }, sw);
82}
83
88template<class PcKrSw, class V>
89auto sampleRelPerms(const PcKrSw& curve, const V& sw)
90{
91 return std::make_tuple(
92 Detail::evalFunctionForRange([&](const auto s){ return curve.krw(s); }, sw),
93 Detail::evalFunctionForRange([&](const auto s){ return curve.krn(s); }, sw)
94 );
95}
96
101template<class PcKrSw, class V>
102auto sampleRelPermDerivatives(const PcKrSw& curve, const V& sw)
103{
104 return std::make_tuple(
105 Detail::evalFunctionForRange([&](const auto s){ return curve.dkrw_dsw(s); }, sw),
106 Detail::evalFunctionForRange([&](const auto s){ return curve.dkrn_dsw(s); }, sw)
107 );
108}
109
110// forward declaration
111template<class S> class GnuplotInterface;
112
113namespace Detail {
114template<class S, class V>
116 const V& x, const V& y,
117 const std::string& curveName,
118 const std::string& curveOptions,
119 const std::string& xLabel,
120 const std::string& yLabel)
121{
122 gnuplot.setXlabel(xLabel);
123 gnuplot.setYlabel(yLabel);
124 gnuplot.addDataSetToPlot(x, y, curveName, curveOptions);
125}
126} // end namespace Detail
127
132namespace Gnuplot {
133
138template<class S, class V>
139void addPcSw(GnuplotInterface<S>& gnuplot, const V& sw, const V& pc,
140 const std::string& curveName = "pc-sw",
141 const std::string& curveOptions = "w l",
142 const std::string& xLabel = "wetting phase saturation [-]",
143 const std::string& yLabel = "capillary pressure [Pa]")
144{
145 Detail::addDataSetToGnuplot(gnuplot, sw, pc, curveName, curveOptions, xLabel, yLabel);
146}
147
152template<class S, class V>
153void addPcSwDerivative(GnuplotInterface<S>& gnuplot, const V& sw, const V& dpc_dsw,
154 const std::string& curveName = "dpc-dsw",
155 const std::string& curveOptions = "w l",
156 const std::string& xLabel = "wetting phase saturation [-]",
157 const std::string& yLabel = "derivative of capillary pressure [Pa]")
158{
159 Detail::addDataSetToGnuplot(gnuplot, sw, dpc_dsw, curveName, curveOptions, xLabel, yLabel);
160}
161
166template<class S, class V>
167void addPcSwInverseDerivative(GnuplotInterface<S>& gnuplot, const V& sw, const V& dpc_dsw,
168 const std::string& curveName = "dsw-dpc",
169 const std::string& curveOptions = "w l",
170 const std::string& xLabel = "capillary pressure [Pa]",
171 const std::string& yLabel = "derivative of saturation [1/Pa]")
172{
173 Detail::addDataSetToGnuplot(gnuplot, sw, dpc_dsw, curveName, curveOptions, xLabel, yLabel);
174}
175
180template<class S, class V>
181void addLog10PcSw(GnuplotInterface<S>& gnuplot, const V& sw, const V& log10pc,
182 const std::string& curveName = "log10-pc-sw",
183 const std::string& curveOptions = "w l",
184 const std::string& xLabel = "wetting phase saturation [-]",
185 const std::string& yLabel = "log10 of capillary pressure [Pa]")
186{
187 Detail::addDataSetToGnuplot(gnuplot, sw, log10pc, curveName, curveOptions, xLabel, yLabel);
188}
189
194template<class S, class V>
195void addRelPerms(GnuplotInterface<S>& gnuplot, const V& sw, const V& krw, const V& krn,
196 const std::string& curveName = "relperm",
197 const std::string& curveOptions = "w l",
198 const std::string& xLabel = "wetting phase saturation [-]",
199 const std::string& yLabel = "relative permeability [-]")
200{
201 Detail::addDataSetToGnuplot(gnuplot, sw, krw, curveName + "_krw", curveOptions, xLabel, yLabel);
202 Detail::addDataSetToGnuplot(gnuplot, sw, krn, curveName + "_krn", curveOptions, xLabel, yLabel);
203}
204
209template<class S, class V>
210void addRelPermDerivatives(GnuplotInterface<S>& gnuplot, const V& sw, const V& krw, const V& krn,
211 const std::string& curveName = "relperm_dsw",
212 const std::string& curveOptions = "w l",
213 const std::string& xLabel = "wetting phase saturation [-]",
214 const std::string& yLabel = "derivative of the relative permeability [-]")
215{
216 Detail::addDataSetToGnuplot(gnuplot, sw, krw, curveName + "_krw", curveOptions, xLabel, yLabel);
217 Detail::addDataSetToGnuplot(gnuplot, sw, krn, curveName + "_krn", curveOptions, xLabel, yLabel);
218}
219
220} // end namespace Gnuplot
221} // end namespace Dumux
222
223#endif
Define some often used mathematical functions.
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:139
auto samplePcSw(const PcKrSw &curve, const V &sw)
sample the pc-sw curve
Definition: plotpckrsw.hh:49
auto samplePcSwDerivative(const PcKrSw &curve, const V &sw)
sample the pc-sw curve derivative wrt sw
Definition: plotpckrsw.hh:59
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:167
auto sampleRelPerms(const PcKrSw &curve, const V &sw)
sample krw-sw and krn-sw curves
Definition: plotpckrsw.hh:89
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:195
auto samplePcSwInverseDerivative(const PcKrSw &curve, const V &pc)
sample the sw-pc curve derivative wrt pc
Definition: plotpckrsw.hh:69
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:181
auto sampleRelPermDerivatives(const PcKrSw &curve, const V &sw)
sample the derivatives of the krw-sw and krn-sw curves
Definition: plotpckrsw.hh:102
auto sampleLog10PcSw(const PcKrSw &curve, const V &sw)
sample sw-pc curve but return the log10 of the capillary pressure
Definition: plotpckrsw.hh:79
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:210
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:153
Definition: adapt.hh:29
Range evalFunctionForRange(const Function &f, const Range &range)
Definition: plotpckrsw.hh:36
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:115
Interface for passing data sets to a file and plotting them, if gnuplot is installed.
Definition: gnuplotinterface.hh:57
void setXlabel(const std::string &label)
Sets the label for the x-axis.
Definition: gnuplotinterface.hh:279
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:242
void setYlabel(const std::string &label)
Sets the label for the y-axis.
Definition: gnuplotinterface.hh:289