24#ifndef DUMUX_PLOT_FLUID_MATRIX_LAW_HH
25#define DUMUX_PLOT_FLUID_MATRIX_LAW_HH
27#warning "This header is deprecated and will be removed after 3.3. Use new 2p material laws and plot tools from io/plotpckrsw.hh"
36template<
class Scalar>
class GnuplotInterface;
42template<
class Scalar,
class MaterialLaw>
45 using MaterialLawParams =
typename MaterialLaw::Params;
64 const MaterialLawParams ¶ms,
65 Scalar lowerSat = 0.0,
66 Scalar upperSat = 1.0,
67 std::string curveName =
"pc-Sw",
68 std::string curveOptions =
"w l")
70 std::vector<Scalar> sw(numIntervals_+1);
71 std::vector<Scalar> pc(numIntervals_+1);
72 Scalar satInterval = upperSat - lowerSat;
74 for (
int i = 0; i <= numIntervals_; i++)
76 sw[i] = lowerSat + satInterval * Scalar(i) / Scalar(numIntervals_);
77 pc[i] = MaterialLaw::pc(params, sw[i]);
80 gnuplot.
setXlabel(
"wetting phase saturation [-]");
81 gnuplot.
setYlabel(
"capillary pressure [Pa]");
96 const MaterialLawParams ¶ms,
97 Scalar lowerSat = 0.0,
98 Scalar upperSat = 1.0,
99 std::string curveName =
"log10_pc-Sw.dat",
100 std::string curveOptions =
"w l")
102 std::vector<Scalar> sw(numIntervals_+1);
103 std::vector<Scalar> log10pc(numIntervals_+1);
104 Scalar satInterval = upperSat - lowerSat;
106 for (
int i = 0; i <= numIntervals_; i++)
108 sw[i] = lowerSat + satInterval * Scalar(i) / Scalar(numIntervals_);
109 const Scalar pc = std::log10(MaterialLaw::pc(params, sw[i]));
110 log10pc[i] = std::isnan(pc) ? 0.0 : pc;
113 gnuplot.
setXlabel(
"wetting phase saturation [-]");
114 gnuplot.
setYlabel(
"log10 of capillary pressure [Pa]");
129 const MaterialLawParams ¶ms,
130 Scalar lowerpc = 0.0,
131 Scalar upperpc = 5000.0,
132 std::string curveName =
"Sw-pc",
133 std::string curveOptions =
"w l")
135 std::vector<Scalar> sw;
136 std::vector<Scalar> pc;
137 Scalar pcInterval = upperpc - lowerpc;
139 Scalar pcTemp, swTemp = 0.0;
140 for (
int i = 0; i <= numIntervals_; i++)
142 pcTemp = lowerpc + pcInterval * Scalar(i) / Scalar(numIntervals_);
143 swTemp = MaterialLaw::sw(params, pcTemp);
144 if (checkValues_(pcTemp, swTemp))
146 pc.push_back(pcTemp);
147 sw.push_back(swTemp);
151 gnuplot.
setXlabel(
"capillary pressure [Pa]");
152 gnuplot.
setYlabel(
"wetting phase saturation [-]");
167 const MaterialLawParams ¶ms,
169 Scalar upperpc = 1e9,
170 std::string curveName =
"Sw-pc",
171 std::string curveOptions =
"w l")
173 std::vector<Scalar> sw;
174 std::vector<Scalar> pc;
175 Scalar pcInterval = std::log10(upperpc) - std::log10(lowerpc);
177 Scalar pcTemp, swTemp = 0.0;
178 for (
int i = 0; i <= numIntervals_; i++)
180 pcTemp = std::log10(lowerpc) + pcInterval * Scalar(i) / Scalar(numIntervals_);
181 swTemp = MaterialLaw::sw(params, std::pow(10, pcTemp));
182 pc.push_back(pcTemp);
183 sw.push_back(swTemp);
186 gnuplot.
setXlabel(
"log10 capillary pressure [Pa]");
187 gnuplot.
setYlabel(
"wetting phase saturation [-]");
202 const MaterialLawParams ¶ms,
203 Scalar lowerSat = 0.0,
204 Scalar upperSat = 1.0,
205 std::string curveName =
"dpcdsw",
206 std::string curveOptions =
"w l")
208 std::vector<Scalar> sw;
209 std::vector<Scalar> dpcdsw;
210 Scalar satInterval = upperSat - lowerSat;
212 Scalar swTemp, dpcdswTemp = 0.0;
213 for (
int i = 0; i <= numIntervals_; i++)
215 swTemp = lowerSat + satInterval * Scalar(i) / Scalar(numIntervals_);
216 dpcdswTemp = MaterialLaw::dpc_dsw(params, swTemp);
217 if (checkValues_(swTemp, dpcdsw))
219 sw.push_back(swTemp);
220 dpcdsw.push_back(dpcdswTemp);
224 gnuplot.
setXlabel(
"wetting phase saturation [-]");
225 gnuplot.
setYlabel(
"gradient of the pc-Sw curve [Pa]");
240 const MaterialLawParams ¶ms,
241 Scalar lowerpc = 0.0,
242 Scalar upperpc = 5000.0,
243 std::string curveName =
"dswdpc",
244 std::string curveOptions =
"w l")
246 std::vector<Scalar> pc;
247 std::vector<Scalar> dswdpc;
248 Scalar pcInterval = upperpc - lowerpc;
250 Scalar dswdpcTemp, pcTemp = 0.0;
251 for (
int i = 0; i <= numIntervals_; i++)
253 pcTemp = lowerpc + pcInterval * Scalar(i) / Scalar(numIntervals_);
254 dswdpcTemp = MaterialLaw::dsw_dpc(params, pcTemp);
255 if (checkValues_(pcTemp, dswdpcTemp))
257 pc.push_back(pcTemp);
258 dswdpc.push_back(dswdpcTemp);
262 gnuplot.
setXlabel(
"capillary pressure [Pa]");
263 gnuplot.
setYlabel(
"gradient of the Sw-pc curve [1/Pa]");
278 const MaterialLawParams ¶ms,
279 Scalar lowerSat = 0.0,
280 Scalar upperSat = 1.0,
281 std::string curveName =
"kr",
282 std::string curveOptions =
"w l")
284 std::vector<Scalar> sw;
285 std::vector<Scalar> krw;
286 std::vector<Scalar> krn;
287 Scalar satInterval = upperSat - lowerSat;
289 Scalar swTemp, krwTemp, krnTemp = 0.0;
290 for (
int i = 0; i <= numIntervals_; i++)
292 swTemp = lowerSat + satInterval * Scalar(i) / Scalar(numIntervals_);
293 krwTemp = MaterialLaw::krw(params, swTemp);
294 krnTemp = MaterialLaw::krn(params, swTemp);
295 if (checkValues_(swTemp, krwTemp) && checkValues_(swTemp, krnTemp))
297 sw.push_back(swTemp);
298 krw.push_back(krwTemp);
299 krn.push_back(krnTemp);
303 gnuplot.
setXlabel(
"wetting phase saturation [-]");
304 gnuplot.
setYlabel(
"relative permeability [-]");
320 const MaterialLawParams ¶ms,
321 Scalar lowerSat = 0.0,
322 Scalar upperSat = 1.0,
323 std::string curveName =
"dkrndsw",
324 std::string curveOptions =
"w l")
326 std::vector<Scalar> sw;
327 std::vector<Scalar> dkrw_dsw;
328 std::vector<Scalar> dkrn_dsw;
329 Scalar satInterval = upperSat - lowerSat;
331 Scalar swTemp, dkrwdswTemp, dkrndswTemp = 0.0;
332 for (
int i = 0; i <= numIntervals_; i++)
334 swTemp = lowerSat + satInterval * Scalar(i) / Scalar(numIntervals_);
335 dkrwdswTemp = MaterialLaw::dkrw_dsw(params, swTemp);
336 dkrndswTemp = MaterialLaw::dkrn_dsw(params, swTemp);
337 if (checkValues_(swTemp, dkrwdswTemp) && checkValues_(swTemp, dkrndswTemp))
339 sw.push_back(swTemp);
340 dkrw_dsw.push_back(dkrwdswTemp);
341 dkrn_dsw.push_back(dkrndswTemp);
345 gnuplot.
setXlabel(
"wetting phase saturation [-]");
346 gnuplot.
setYlabel(
"gradient of the kr-Sw function [-]");
347 gnuplot.
addDataSetToPlot(sw, dkrw_dsw, curveName +
"_dkrw_dsw", curveOptions);
348 gnuplot.
addDataSetToPlot(sw, dkrn_dsw, curveName +
"_dkrn_dsw", curveOptions);
358 bool checkValues_(Scalar value1, Scalar value2)
362 return !isnan(value1) && !isinf(value1)
363 && !isnan(value2) && !isinf(value2);
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
Interface for plotting the two-phase fluid-matrix-interaction laws.
Definition: plotmateriallaw.hh:44
void addLog10PcSwCurve(GnuplotInterface< Scalar > &gnuplot, const MaterialLawParams ¶ms, Scalar lowerSat=0.0, Scalar upperSat=1.0, std::string curveName="log10_pc-Sw.dat", std::string curveOptions="w l")
Add a capillary pressure-saturation data set to the plot.
Definition: plotmateriallaw.hh:95
void addswpccurve(GnuplotInterface< Scalar > &gnuplot, const MaterialLawParams ¶ms, Scalar lowerpc=0.0, Scalar upperpc=5000.0, std::string curveName="Sw-pc", std::string curveOptions="w l")
Add a saturation-capillary pressure data set to the plot.
Definition: plotmateriallaw.hh:128
PlotMaterialLaw()
Constructor.
Definition: plotmateriallaw.hh:49
void addSwLog10PcCurve(GnuplotInterface< Scalar > &gnuplot, const MaterialLawParams ¶ms, Scalar lowerpc=1, Scalar upperpc=1e9, std::string curveName="Sw-pc", std::string curveOptions="w l")
Add a saturation-capillary pressure data set to the plot.
Definition: plotmateriallaw.hh:166
void adddswdpccurve(GnuplotInterface< Scalar > &gnuplot, const MaterialLawParams ¶ms, Scalar lowerpc=0.0, Scalar upperpc=5000.0, std::string curveName="dswdpc", std::string curveOptions="w l")
Add a saturation-capillary pressure gradient data set to the plot.
Definition: plotmateriallaw.hh:239
void adddkrdswcurves(GnuplotInterface< Scalar > &gnuplot, const MaterialLawParams ¶ms, Scalar lowerSat=0.0, Scalar upperSat=1.0, std::string curveName="dkrndsw", std::string curveOptions="w l")
Add relative permeabilities gradients data sets to the plot.
Definition: plotmateriallaw.hh:319
void adddpcdswcurve(GnuplotInterface< Scalar > &gnuplot, const MaterialLawParams ¶ms, Scalar lowerSat=0.0, Scalar upperSat=1.0, std::string curveName="dpcdsw", std::string curveOptions="w l")
Add a capillary pressure-saturation gradient data set to the plot.
Definition: plotmateriallaw.hh:201
void addpcswcurve(GnuplotInterface< Scalar > &gnuplot, const MaterialLawParams ¶ms, Scalar lowerSat=0.0, Scalar upperSat=1.0, std::string curveName="pc-Sw", std::string curveOptions="w l")
Add a capillary pressure-saturation data set to the plot.
Definition: plotmateriallaw.hh:63
void addkrcurves(GnuplotInterface< Scalar > &gnuplot, const MaterialLawParams ¶ms, Scalar lowerSat=0.0, Scalar upperSat=1.0, std::string curveName="kr", std::string curveOptions="w l")
Add relative permeabilities data sets to the plot.
Definition: plotmateriallaw.hh:277