version 3.8
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// 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_DATA_SPLINE_MATERIAL_HH
13#define DUMUX_MATERIAL_FLUIDMATRIX_DATA_SPLINE_MATERIAL_HH
14
15#include <algorithm>
20
21namespace Dumux::FluidMatrix {
22
33template <class S, bool approximatePcSwInverse = false>
35: public Adapter<DataSplineTwoPMaterialLaw<S, approximatePcSwInverse>, PcKrSw>
36{
37public:
38
39 using Scalar = S;
40
44 static constexpr int numFluidPhases()
45 { return 2; }
46
47 static constexpr bool isRegularized()
48 { return false; }
49
55
59 explicit DataSplineTwoPMaterialLaw(const std::string& paramGroup)
60 {
61 using V = std::vector<Scalar>;
62 const std::string swPcGroup = paramGroup.empty() ? "Pc" : paramGroup + ".Pc";
63 const auto swPc = getParamFromGroup<V>(swPcGroup, "SwData");
64 const std::string krwPcGroup = paramGroup.empty() ? "Krw" : paramGroup + ".Krw";
65 const auto swKrw = getParamFromGroup<V>(krwPcGroup, "SwData");
66 const std::string krnPcGroup = paramGroup.empty() ? "Krn" : paramGroup + ".Krn";
67 const auto swKrn = getParamFromGroup<V>(krnPcGroup, "SwData");
68
69 const auto pc = getParamFromGroup<V>(paramGroup, "PcData");
70 const auto krw = getParamFromGroup<V>(paramGroup, "KrwData");
71 const auto krn = getParamFromGroup<V>(paramGroup, "KrnData");
72
73 updateData_(swPc, pc, swKrw, krw, swKrn, krn);
74 }
75
79 DataSplineTwoPMaterialLaw(const std::vector<Scalar>& swPc,
80 const std::vector<Scalar>& pc,
81 const std::vector<Scalar>& swKrw,
82 const std::vector<Scalar>& krw,
83 const std::vector<Scalar>& swKrn,
84 const std::vector<Scalar>& krn)
85 {
86 updateData_(swPc, pc, swKrw, krw, swKrn, krn);
87 }
88
92 Scalar pc(const Scalar sw) const
93 {
94 if constexpr (approximatePcSwInverse)
95 return pcSpline_.evalInverse(sw);
96 else
97 return pcSpline_.eval(sw);
98 }
99
103 Scalar dpc_dsw(const Scalar sw) const
104 {
105 if constexpr (approximatePcSwInverse)
106 return 1.0/pcSpline_.evalDerivative(pcSpline_.evalInverse(sw));
107 else
108 return pcSpline_.evalDerivative(sw);
109 }
110
114 Scalar sw(const Scalar pc) const
115 {
116 if constexpr (approximatePcSwInverse)
117 return pcSpline_.eval(pc);
118 else
119 return pcSpline_.evalInverse(pc);
120 }
121
125 Scalar dsw_dpc(const Scalar pc) const
126 {
127 if constexpr (approximatePcSwInverse)
128 return pcSpline_.evalDerivative(pc);
129 else
130 return 1.0/pcSpline_.evalDerivative(pcSpline_.evalInverse(pc));
131 }
132
136 Scalar krw(const Scalar sw) const
137 {
138 return krwSpline_.eval(sw);
139 }
140
145 {
146 return krwSpline_.evalDerivative(sw);
147 }
148
152 Scalar krn(const Scalar sw) const
153 {
154 return krnSpline_.eval(sw);
155 }
156
161 {
162 return krnSpline_.evalDerivative(sw);
163 }
164
165private:
166 void updateData_(const std::vector<Scalar>& swPc,
167 const std::vector<Scalar>& pc,
168 const std::vector<Scalar>& swKrw,
169 const std::vector<Scalar>& krw,
170 const std::vector<Scalar>& swKrn,
171 const std::vector<Scalar>& krn)
172 {
173 if constexpr (approximatePcSwInverse)
174 pcSpline_.updatePoints(pc, swPc);
175 else
176 pcSpline_.updatePoints(swPc, pc);
177
178 krwSpline_.updatePoints(swKrw, krw);
179 krnSpline_.updatePoints(swKrn, krn);
180 }
181
185};
186
187} // end namespace Dumux::FluidMatrix
188
189#endif
Pc- and Kr-sw curves based on monotone splines through given data points.
Definition: datasplinemateriallaw.hh:36
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: datasplinemateriallaw.hh:92
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:44
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition: datasplinemateriallaw.hh:144
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:160
S Scalar
Definition: datasplinemateriallaw.hh:39
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:79
static constexpr bool isRegularized()
Definition: datasplinemateriallaw.hh:47
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: datasplinemateriallaw.hh:103
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition: datasplinemateriallaw.hh:136
DataSplineTwoPMaterialLaw(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: datasplinemateriallaw.hh:59
Scalar sw(const Scalar pc) const
The saturation-capillary pressure curve.
Definition: datasplinemateriallaw.hh:114
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition: datasplinemateriallaw.hh:152
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: datasplinemateriallaw.hh:125
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:63
Scalar evalInverse(const Scalar y) const
Evaluate the inverse function.
Definition: monotonecubicspline.hh:135
Scalar eval(const Scalar x) const
Evaluate the y value at a given x value.
Definition: monotonecubicspline.hh:104
Scalar evalDerivative(const Scalar x) const
Evaluate the first derivative dy/dx at a given x value.
Definition: monotonecubicspline.hh:119
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
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