24#ifndef DUMUX_MATERIAL_FLUIDMATRIX_DATA_SPLINE_MATERIAL_HH
25#define DUMUX_MATERIAL_FLUIDMATRIX_DATA_SPLINE_MATERIAL_HH
45template <
class S,
bool approximatePcSwInverse = false>
47:
public Adapter<DataSplineTwoPMaterialLaw<S, approximatePcSwInverse>, PcKrSw>
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");
81 const auto pc = getParamFromGroup<V>(paramGroup,
"PcData");
82 const auto krw = getParamFromGroup<V>(paramGroup,
"KrwData");
83 const auto krn = getParamFromGroup<V>(paramGroup,
"KrnData");
85 updateData_(swPc,
pc, swKrw,
krw, swKrn,
krn);
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)
98 updateData_(swPc,
pc, swKrw,
krw, swKrn,
krn);
106 if constexpr (approximatePcSwInverse)
109 return pcSpline_.
eval(
sw);
117 if constexpr (approximatePcSwInverse)
128 if constexpr (approximatePcSwInverse)
129 return pcSpline_.
eval(
pc);
139 if constexpr (approximatePcSwInverse)
150 return krwSpline_.
eval(
sw);
166 return krnSpline_.
eval(
sw);
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)
185 if constexpr (approximatePcSwInverse)
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 ¶mGroup)
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