27 #ifndef OPM_SPLINE_TWO_PHASE_MATERIAL_HPP
28 #define OPM_SPLINE_TWO_PHASE_MATERIAL_HPP
43 template <
class TraitsT,
class ParamsT = SplineTwoPhaseMaterialParams<TraitsT> >
46 typedef typename ParamsT::SamplePoints SamplePoints;
56 typedef typename Traits::Scalar
Scalar;
61 "The piecewise linear two-phase capillary pressure law only"
62 "applies to the case of two fluid phases");
91 template <
class Container,
class Flu
idState>
94 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
96 values[Traits::wettingPhaseIdx] = 0.0;
97 values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fluidState);
104 template <
class Container,
class Flu
idState>
106 {
throw std::logic_error(
"Not implemented: saturations()"); }
111 template <
class Container,
class Flu
idState>
114 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
116 values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fluidState);
117 values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fluidState);
123 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
124 static Evaluation
pcnw(
const Params& params,
const FluidState& fluidState)
126 const Evaluation&
Sw =
127 decay<Evaluation>(fluidState.saturation(Traits::wettingPhaseIdx));
135 template <
class Evaluation>
139 const auto& pcnwSpline = params.pcnwSpline();
140 if (
Sw <= pcnwSpline.xAt(0))
141 return Evaluation(pcnwSpline.valueAt(0));
142 if (
Sw >= pcnwSpline.xAt(pcnwSpline.numSamples() - 1))
143 return Evaluation(pcnwSpline.valueAt(pcnwSpline.numSamples() - 1));
145 return pcnwSpline.eval(
Sw);
148 template <
class Evaluation>
149 static Evaluation twoPhaseSatPcnwInv(
const Params& params,
const Evaluation&
pcnw)
151 static const Evaluation nil(0.0);
154 const auto& pcnwSpline = params.pcnwSpline();
155 if (
pcnw >= pcnwSpline.valueAt(0))
156 return Evaluation(pcnwSpline.xAt(0));
157 if (
pcnw <= pcnwSpline.y(pcnwSpline.numSamples() - 1))
158 return Evaluation(pcnwSpline.xAt(pcnwSpline.numSamples() - 1));
162 return pcnwSpline.intersect(nil, nil, nil,
pcnw);
168 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
169 static Evaluation
Sw(
const Params& ,
const FluidState& )
170 {
throw std::logic_error(
"Not implemented: Sw()"); }
172 template <
class Evaluation>
173 static Evaluation twoPhaseSatSw(
const Params& ,
const Evaluation& )
174 {
throw std::logic_error(
"Not implemented: twoPhaseSatSw()"); }
180 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
181 static Evaluation
Sn(
const Params& params,
const FluidState& fluidState)
182 {
return 1 - Sw<FluidState, Evaluation>(params, fluidState); }
184 template <
class Evaluation>
185 static Evaluation twoPhaseSatSn(
const Params& params,
const Evaluation& pC)
186 {
return 1 - twoPhaseSatSw(params, pC); }
192 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
193 static Evaluation
krw(
const Params& params,
const FluidState& fluidState)
195 const Evaluation&
Sw =
196 decay<Evaluation>(fluidState.saturation(Traits::wettingPhaseIdx));
198 return twoPhaseSatKrw(params,
Sw);
201 template <
class Evaluation>
202 static Evaluation twoPhaseSatKrw(
const Params& params,
const Evaluation&
Sw)
204 const auto& krwSpline = params.krwSpline();
205 if (
Sw <= krwSpline.xAt(0))
206 return Evaluation(krwSpline.valueAt(0));
207 if (
Sw >= krwSpline.xAt(krwSpline.numSamples() - 1))
208 return Evaluation(krwSpline.valueAt(krwSpline.numSamples() - 1));
210 return krwSpline.eval(
Sw);
213 template <
class Evaluation>
214 static Evaluation twoPhaseSatKrwInv(
const Params& params,
const Evaluation&
krw)
216 static const Evaluation nil(0.0);
218 const auto& krwSpline = params.krwSpline();
219 if (
krw <= krwSpline.valueAt(0))
220 return Evaluation(krwSpline.xAt(0));
221 if (
krw >= krwSpline.valueAt(krwSpline.numSamples() - 1))
222 return Evaluation(krwSpline.xAt(krwSpline.numSamples() - 1));
224 return krwSpline.intersect(nil, nil, nil,
krw);
231 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
232 static Evaluation
krn(
const Params& params,
const FluidState& fluidState)
234 const Evaluation&
Sn =
235 decay<Evaluation>(fluidState.saturation(Traits::nonWettingPhaseIdx));
237 return twoPhaseSatKrn(params, 1.0 -
Sn);
240 template <
class Evaluation>
241 static Evaluation twoPhaseSatKrn(
const Params& params,
const Evaluation&
Sw)
243 const auto& krnSpline = params.krnSpline();
244 if (
Sw <= krnSpline.xAt(0))
245 return Evaluation(krnSpline.valueAt(0));
246 if (
Sw >= krnSpline.xAt(krnSpline.numSamples() - 1))
247 return Evaluation(krnSpline.valueAt(krnSpline.numSamples() - 1));
249 return krnSpline.eval(
Sw);
252 template <
class Evaluation>
253 static Evaluation twoPhaseSatKrnInv(
const Params& params,
const Evaluation&
krn)
255 static const Evaluation nil(0.0);
257 const auto& krnSpline = params.krnSpline();
258 if (
krn >= krnSpline.valueAt(0))
259 return Evaluation(krnSpline.xAt(0));
260 if (
krn <= krnSpline.valueAt(krnSpline.numSamples() - 1))
261 return Evaluation(krnSpline.xAt(krnSpline.numSamples() - 1));
263 return krnSpline.intersect(nil, nil, nil,
krn);
Specification of the material parameters for a two-phase material law which uses a table and spline-b...
Implementation of a tabulated capillary pressure and relperm law which uses spline curves as interpol...
Definition: SplineTwoPhaseMaterial.hpp:45
static const bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition: SplineTwoPhaseMaterial.hpp:66
static Evaluation krw(const Params ¶ms, const FluidState &fluidState)
The relative permeability for the wetting phase of the porous medium.
Definition: SplineTwoPhaseMaterial.hpp:193
static const bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition: SplineTwoPhaseMaterial.hpp:82
Traits::Scalar Scalar
The type of the scalar values for this law.
Definition: SplineTwoPhaseMaterial.hpp:56
static void capillaryPressures(Container &values, const Params ¶ms, const FluidState &fluidState)
The capillary pressure-saturation curve.
Definition: SplineTwoPhaseMaterial.hpp:92
static Evaluation Sw(const Params &, const FluidState &)
The saturation-capillary pressure curve.
Definition: SplineTwoPhaseMaterial.hpp:169
static void relativePermeabilities(Container &values, const Params ¶ms, const FluidState &fluidState)
The relative permeabilities.
Definition: SplineTwoPhaseMaterial.hpp:112
static void saturations(Container &, const Params &, const FluidState &)
The saturations of the fluid phases starting from their pressure differences.
Definition: SplineTwoPhaseMaterial.hpp:105
static Evaluation krn(const Params ¶ms, const FluidState &fluidState)
The relative permeability for the non-wetting phase of the porous medium.
Definition: SplineTwoPhaseMaterial.hpp:232
static const bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition: SplineTwoPhaseMaterial.hpp:70
static const bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition: SplineTwoPhaseMaterial.hpp:78
static Evaluation pcnw(const Params ¶ms, const FluidState &fluidState)
The capillary pressure-saturation curve.
Definition: SplineTwoPhaseMaterial.hpp:124
ParamsT Params
The type of the parameter objects for this law.
Definition: SplineTwoPhaseMaterial.hpp:53
static Evaluation Sn(const Params ¶ms, const FluidState &fluidState)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition: SplineTwoPhaseMaterial.hpp:181
static const bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition: SplineTwoPhaseMaterial.hpp:86
static const bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition: SplineTwoPhaseMaterial.hpp:74
TraitsT Traits
The traits class for this material law.
Definition: SplineTwoPhaseMaterial.hpp:50
static const int numPhases
The number of fluid phases.
Definition: SplineTwoPhaseMaterial.hpp:59
static Evaluation twoPhaseSatPcnw(const Params ¶ms, const Evaluation &Sw)
The saturation-capillary pressure curve.
Definition: SplineTwoPhaseMaterial.hpp:136