3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
datasplinemateriallaw.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_MATERIAL_FLUIDMATRIX_DATA_SPLINE_MATERIAL_HH
25#define DUMUX_MATERIAL_FLUIDMATRIX_DATA_SPLINE_MATERIAL_HH
26
27#include <algorithm>
32
33namespace Dumux::FluidMatrix {
34
45template <class S, bool approximatePcSwInverse = false>
47: public Adapter<DataSplineTwoPMaterialLaw<S, approximatePcSwInverse>, PcKrSw>
48{
49public:
50
51 using Scalar = S;
52
56 static constexpr int numFluidPhases()
57 { return 2; }
58
59 static constexpr bool isRegularized()
60 { return false; }
61
67
71 explicit DataSplineTwoPMaterialLaw(const std::string& paramGroup)
72 {
73 using V = std::vector<Scalar>;
74 const std::string swPcGroup = paramGroup.empty() ? "Pc" : paramGroup + ".Pc";
75 const auto swPc = getParamFromGroup<V>(swPcGroup, "SwData");
76 const std::string krwPcGroup = paramGroup.empty() ? "Krw" : paramGroup + ".Krw";
77 const auto swKrw = getParamFromGroup<V>(krwPcGroup, "SwData");
78 const std::string krnPcGroup = paramGroup.empty() ? "Krn" : paramGroup + ".Krn";
79 const auto swKrn = getParamFromGroup<V>(krnPcGroup, "SwData");
80
81 const auto pc = getParamFromGroup<V>(paramGroup, "PcData");
82 const auto krw = getParamFromGroup<V>(paramGroup, "KrwData");
83 const auto krn = getParamFromGroup<V>(paramGroup, "KrnData");
84
85 updateData_(swPc, pc, swKrw, krw, swKrn, krn);
86 }
87
91 DataSplineTwoPMaterialLaw(const std::vector<Scalar>& swPc,
92 const std::vector<Scalar>& pc,
93 const std::vector<Scalar>& swKrw,
94 const std::vector<Scalar>& krw,
95 const std::vector<Scalar>& swKrn,
96 const std::vector<Scalar>& krn)
97 {
98 updateData_(swPc, pc, swKrw, krw, swKrn, krn);
99 }
100
104 Scalar pc(const Scalar sw) const
105 {
106 if constexpr (approximatePcSwInverse)
107 return pcSpline_.evalInverse(sw);
108 else
109 return pcSpline_.eval(sw);
110 }
111
115 Scalar dpc_dsw(const Scalar sw) const
116 {
117 if constexpr (approximatePcSwInverse)
118 return 1.0/pcSpline_.evalDerivative(pcSpline_.evalInverse(sw));
119 else
120 return pcSpline_.evalDerivative(sw);
121 }
122
126 Scalar sw(const Scalar pc) const
127 {
128 if constexpr (approximatePcSwInverse)
129 return pcSpline_.eval(pc);
130 else
131 return pcSpline_.evalInverse(pc);
132 }
133
137 Scalar dsw_dpc(const Scalar pc) const
138 {
139 if constexpr (approximatePcSwInverse)
140 return pcSpline_.evalDerivative(pc);
141 else
142 return 1.0/pcSpline_.evalDerivative(pcSpline_.evalInverse(pc));
143 }
144
148 Scalar krw(const Scalar sw) const
149 {
150 return krwSpline_.eval(sw);
151 }
152
157 {
158 return krwSpline_.evalDerivative(sw);
159 }
160
164 Scalar krn(const Scalar sw) const
165 {
166 return krnSpline_.eval(sw);
167 }
168
173 {
174 return krnSpline_.evalDerivative(sw);
175 }
176
177private:
178 void updateData_(const std::vector<Scalar>& swPc,
179 const std::vector<Scalar>& pc,
180 const std::vector<Scalar>& swKrw,
181 const std::vector<Scalar>& krw,
182 const std::vector<Scalar>& swKrn,
183 const std::vector<Scalar>& krn)
184 {
185 if constexpr (approximatePcSwInverse)
186 pcSpline_.updatePoints(pc, swPc);
187 else
188 pcSpline_.updatePoints(swPc, pc);
189
190 krwSpline_.updatePoints(swKrw, krw);
191 krnSpline_.updatePoints(swKrn, krn);
192 }
193
197};
198
199} // end namespace Dumux::FluidMatrix
200
201#endif
A monotone cubic spline.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Implementation helper for capillary pressure and relative permeability <-> saturation relations for t...
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Definition: brookscorey.hh:35
void updatePoints(const std::vector< Scalar > &x, const std::vector< Scalar > &y)
Create a monotone cubic spline from the control points (x[i], y[i])
Definition: monotonecubicspline.hh:75
Scalar evalInverse(const Scalar y) const
Evaluate the inverse function.
Definition: monotonecubicspline.hh:147
Scalar eval(const Scalar x) const
Evaluate the y value at a given x value.
Definition: monotonecubicspline.hh:116
Scalar evalDerivative(const Scalar x) const
Evaluate the first derivative dy/dx at a given x value.
Definition: monotonecubicspline.hh:131
Pc- and Kr-sw curves based on monotone splines through given data points.
Definition: datasplinemateriallaw.hh:48
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: datasplinemateriallaw.hh:104
DataSplineTwoPMaterialLaw()=delete
Deleted default constructor (so we are never in an undefined state)
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: datasplinemateriallaw.hh:56
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition: datasplinemateriallaw.hh:156
Scalar dkrn_dsw(const Scalar sw) const
The derivative of the relative permeability for the non-wetting phase w.r.t. saturation.
Definition: datasplinemateriallaw.hh:172
S Scalar
Definition: datasplinemateriallaw.hh:51
DataSplineTwoPMaterialLaw(const std::vector< Scalar > &swPc, const std::vector< Scalar > &pc, const std::vector< Scalar > &swKrw, const std::vector< Scalar > &krw, const std::vector< Scalar > &swKrn, const std::vector< Scalar > &krn)
Construct from user data.
Definition: datasplinemateriallaw.hh:91
static constexpr bool isRegularized()
Definition: datasplinemateriallaw.hh:59
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: datasplinemateriallaw.hh:115
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition: datasplinemateriallaw.hh:148
DataSplineTwoPMaterialLaw(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: datasplinemateriallaw.hh:71
Scalar sw(const Scalar pc) const
The saturation-capillary pressure curve.
Definition: datasplinemateriallaw.hh:126
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition: datasplinemateriallaw.hh:164
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: datasplinemateriallaw.hh:137
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:67