27#ifndef OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP
28#define OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP
48template <
class TraitsT,
class ParamsT = PiecewiseLinearTwoPhaseMaterialParams<TraitsT> >
51 using ValueVector =
typename ParamsT::ValueVector;
61 using Scalar =
typename Traits::Scalar;
66 "The piecewise linear two-phase capillary pressure law only"
67 "applies to the case of two fluid phases");
96 template <
class Container,
class Flu
idState>
99 using Evaluation =
typename std::remove_reference<
decltype(values[0])>::type;
101 values[Traits::wettingPhaseIdx] = 0.0;
102 values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fs);
109 template <
class Container,
class Flu
idState>
111 {
throw std::logic_error(
"Not implemented: saturations()"); }
116 template <
class Container,
class Flu
idState>
119 using Evaluation =
typename std::remove_reference<
decltype(values[0])>::type;
121 values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fs);
122 values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fs);
128 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
129 static Evaluation
pcnw(
const Params& params,
const FluidState& fs)
132 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
140 template <
class Evaluation>
142 {
return eval_(params.SwPcwnSamples(), params.pcnwSamples(),
Sw); }
144 template <
class Evaluation>
146 {
return eval_(params.pcnwSamples(), params.SwPcwnSamples(),
pcnw); }
151 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
152 static Evaluation
Sw(
const Params& ,
const FluidState& )
153 {
throw std::logic_error(
"Not implemented: Sw()"); }
155 template <
class Evaluation>
157 {
throw std::logic_error(
"Not implemented: twoPhaseSatSw()"); }
163 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
164 static Evaluation
Sn(
const Params& params,
const FluidState& fs)
165 {
return 1 - Sw<FluidState, Scalar>(params, fs); }
167 template <
class Evaluation>
175 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
176 static Evaluation
krw(
const Params& params,
const FluidState& fs)
179 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
184 template <
class Evaluation>
186 {
return eval_(params.SwKrwSamples(), params.krwSamples(),
Sw); }
188 template <
class Evaluation>
190 {
return eval_(params.krwSamples(), params.SwKrwSamples(),
krw); }
196 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
197 static Evaluation
krn(
const Params& params,
const FluidState& fs)
200 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
205 template <
class Evaluation>
207 {
return eval_(params.SwKrnSamples(), params.krnSamples(),
Sw); }
209 template <
class Evaluation>
211 {
return eval_(params.krnSamples(), params.SwKrnSamples(),
krn); }
214 template <
class Evaluation>
215 static Evaluation eval_(
const ValueVector& xValues,
216 const ValueVector& yValues,
219 if (xValues.front() < xValues.back())
220 return evalAscending_(xValues, yValues, x);
221 return evalDescending_(xValues, yValues, x);
224 template <
class Evaluation>
225 static Evaluation evalAscending_(
const ValueVector& xValues,
226 const ValueVector& yValues,
229 if (x <= xValues.front())
230 return yValues.front();
231 if (x >= xValues.back())
232 return yValues.back();
234 size_t segIdx = findSegmentIndex_(xValues,
scalarValue(x));
236 Scalar x0 = xValues[segIdx];
237 Scalar x1 = xValues[segIdx + 1];
239 Scalar y0 = yValues[segIdx];
240 Scalar y1 = yValues[segIdx + 1];
242 Scalar m = (y1 - y0)/(x1 - x0);
244 return y0 + (x - x0)*m;
247 template <
class Evaluation>
248 static Evaluation evalDescending_(
const ValueVector& xValues,
249 const ValueVector& yValues,
252 if (x >= xValues.front())
253 return yValues.front();
254 if (x <= xValues.back())
255 return yValues.back();
257 size_t segIdx = findSegmentIndexDescending_(xValues,
scalarValue(x));
259 Scalar x0 = xValues[segIdx];
260 Scalar x1 = xValues[segIdx + 1];
262 Scalar y0 = yValues[segIdx];
263 Scalar y1 = yValues[segIdx + 1];
265 Scalar m = (y1 - y0)/(x1 - x0);
267 return y0 + (x - x0)*m;
270 template <
class Evaluation>
271 static Evaluation evalDeriv_(
const ValueVector& xValues,
272 const ValueVector& yValues,
275 if (x <= xValues.front())
277 if (x >= xValues.back())
280 size_t segIdx = findSegmentIndex_(xValues,
scalarValue(x));
282 Scalar x0 = xValues[segIdx];
283 Scalar x1 = xValues[segIdx + 1];
285 Scalar y0 = yValues[segIdx];
286 Scalar y1 = yValues[segIdx + 1];
288 return (y1 - y0)/(x1 - x0);
291 static size_t findSegmentIndex_(
const ValueVector& xValues,
Scalar x)
293 assert(xValues.size() > 1);
294 size_t n = xValues.size() - 1;
295 if (xValues.back() <= x)
297 else if (x <= xValues.front())
301 size_t lowIdx = 0, highIdx = n;
302 while (lowIdx + 1 < highIdx) {
303 size_t curIdx = (lowIdx + highIdx)/2;
304 if (xValues[curIdx] < x)
313 static size_t findSegmentIndexDescending_(
const ValueVector& xValues,
Scalar x)
315 assert(xValues.size() > 1);
316 size_t n = xValues.size() - 1;
317 if (x <= xValues.back())
319 else if (xValues.front() <= x)
323 size_t lowIdx = 0, highIdx = n;
324 while (lowIdx + 1 < highIdx) {
325 size_t curIdx = (lowIdx + highIdx)/2;
326 if (xValues[curIdx] >= x)
Implementation of a tabulated, piecewise linear capillary pressure law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:50
TraitsT Traits
The traits class for this material law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:55
ParamsT Params
The type of the parameter objects for this law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:58
static Evaluation twoPhaseSatSw(const Params &, const Evaluation &)
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:156
static constexpr int numPhases
The number of fluid phases.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:64
static Evaluation krn(const Params ¶ms, const FluidState &fs)
The relative permeability for the non-wetting phase of the porous medium.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:197
static Evaluation twoPhaseSatPcnw(const Params ¶ms, const Evaluation &Sw)
The saturation-capillary pressure curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:141
static Evaluation krw(const Params ¶ms, const FluidState &fs)
The relative permeability for the wetting phase of the porous medium.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:176
static constexpr bool implementsTwoPhaseSatApi
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:75
static Evaluation twoPhaseSatKrn(const Params ¶ms, const Evaluation &Sw)
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:206
static void relativePermeabilities(Container &values, const Params ¶ms, const FluidState &fs)
The relative permeabilities.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:117
typename Traits::Scalar Scalar
The type of the scalar values for this law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:61
static Evaluation Sn(const Params ¶ms, const FluidState &fs)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:164
static Evaluation twoPhaseSatKrw(const Params ¶ms, const Evaluation &Sw)
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:185
static constexpr bool isCompositionDependent
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:91
static constexpr bool isSaturationDependent
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:79
static void saturations(Container &, const Params &, const FluidState &)
The saturations of the fluid phases starting from their pressure differences.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:110
static constexpr bool isTemperatureDependent
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:87
static constexpr bool isPressureDependent
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:83
static Evaluation twoPhaseSatPcnwInv(const Params ¶ms, const Evaluation &pcnw)
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:145
static constexpr bool implementsTwoPhaseApi
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:71
static Evaluation pcnw(const Params ¶ms, const FluidState &fs)
The capillary pressure-saturation curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:129
static Evaluation twoPhaseSatSn(const Params ¶ms, const Evaluation &pC)
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:168
static Evaluation twoPhaseSatKrnInv(const Params ¶ms, const Evaluation &krn)
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:210
static void capillaryPressures(Container &values, const Params ¶ms, const FluidState &fs)
The capillary pressure-saturation curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:97
static Evaluation twoPhaseSatKrwInv(const Params ¶ms, const Evaluation &krw)
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:189
static Evaluation Sw(const Params &, const FluidState &)
The saturation-capillary pressure curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:152
Definition: Air_Mesitylene.hpp:34
auto scalarValue(const Evaluation &val) -> decltype(MathToolbox< Evaluation >::scalarValue(val))
Definition: MathToolbox.hpp:335