27#ifndef OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP
28#define OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP
31#include <opm/common/TimingMacros.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 OPM_TIMEFUNCTION_LOCAL();
100 using Evaluation =
typename std::remove_reference<
decltype(values[0])>::type;
102 values[Traits::wettingPhaseIdx] = 0.0;
103 values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fs);
110 template <
class Container,
class Flu
idState>
112 {
throw std::logic_error(
"Not implemented: saturations()"); }
117 template <
class Container,
class Flu
idState>
120 OPM_TIMEFUNCTION_LOCAL();
121 using Evaluation =
typename std::remove_reference<
decltype(values[0])>::type;
123 values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fs);
124 values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fs);
130 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
131 static Evaluation
pcnw(
const Params& params,
const FluidState& fs)
133 OPM_TIMEFUNCTION_LOCAL();
135 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
143 template <
class Evaluation>
146 OPM_TIMEFUNCTION_LOCAL();
147 return eval_(params.SwPcwnSamples(), params.pcnwSamples(),
Sw);
150 template <
class Evaluation>
151 static Evaluation twoPhaseSatPcnwInv(
const Params& params,
const Evaluation&
pcnw)
152 {
return eval_(params.pcnwSamples(), params.SwPcwnSamples(),
pcnw); }
157 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
158 static Evaluation
Sw(
const Params& ,
const FluidState& )
159 {
throw std::logic_error(
"Not implemented: Sw()"); }
161 template <
class Evaluation>
162 static Evaluation twoPhaseSatSw(
const Params& ,
const Evaluation& )
163 {
throw std::logic_error(
"Not implemented: twoPhaseSatSw()"); }
169 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
170 static Evaluation
Sn(
const Params& params,
const FluidState& fs)
171 {
return 1 - Sw<FluidState, Scalar>(params, fs); }
173 template <
class Evaluation>
174 static Evaluation twoPhaseSatSn(
const Params& params,
const Evaluation& pC)
175 {
return 1 - twoPhaseSatSw(params, pC); }
181 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
182 static Evaluation
krw(
const Params& params,
const FluidState& fs)
184 OPM_TIMEFUNCTION_LOCAL();
186 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
188 return twoPhaseSatKrw(params,
Sw);
191 template <
class Evaluation>
192 static Evaluation twoPhaseSatKrw(
const Params& params,
const Evaluation&
Sw)
194 OPM_TIMEFUNCTION_LOCAL();
195 return eval_(params.SwKrwSamples(), params.krwSamples(),
Sw);
198 template <
class Evaluation>
199 static Evaluation twoPhaseSatKrwInv(
const Params& params,
const Evaluation&
krw)
201 OPM_TIMEFUNCTION_LOCAL();
202 return eval_(params.krwSamples(), params.SwKrwSamples(),
krw);
209 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
210 static Evaluation
krn(
const Params& params,
const FluidState& fs)
212 OPM_TIMEFUNCTION_LOCAL();
214 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
216 return twoPhaseSatKrn(params,
Sw);
219 template <
class Evaluation>
220 static Evaluation twoPhaseSatKrn(
const Params& params,
const Evaluation&
Sw)
222 OPM_TIMEFUNCTION_LOCAL();
223 return eval_(params.SwKrnSamples(), params.krnSamples(),
Sw);
226 template <
class Evaluation>
227 static Evaluation twoPhaseSatKrnInv(
const Params& params,
const Evaluation&
krn)
228 {
return eval_(params.krnSamples(), params.SwKrnSamples(),
krn); }
230 template <
class Evaluation>
231 static size_t findSegmentIndex(
const ValueVector& xValues,
const Evaluation& x){
232 return findSegmentIndex_(xValues, scalarValue(x));
235 template <
class Evaluation>
236 static size_t findSegmentIndexDescending(
const ValueVector& xValues,
const Evaluation& x){
237 return findSegmentIndexDescending_(xValues, scalarValue(x));
240 template <
class Evaluation>
241 static Evaluation eval(
const ValueVector& xValues,
const ValueVector& yValues,
const Evaluation& x,
unsigned segIdx){
242 Scalar x0 = xValues[segIdx];
243 Scalar x1 = xValues[segIdx + 1];
245 Scalar y0 = yValues[segIdx];
246 Scalar y1 = yValues[segIdx + 1];
248 Scalar m = (y1 - y0)/(x1 - x0);
250 return y0 + (x - x0)*m;
254 template <
class Evaluation>
255 static Evaluation eval_(
const ValueVector& xValues,
256 const ValueVector& yValues,
259 OPM_TIMEFUNCTION_LOCAL();
260 if (xValues.front() < xValues.back())
261 return evalAscending_(xValues, yValues, x);
262 return evalDescending_(xValues, yValues, x);
265 template <
class Evaluation>
266 static Evaluation evalAscending_(
const ValueVector& xValues,
267 const ValueVector& yValues,
270 OPM_TIMEFUNCTION_LOCAL();
271 if (x <= xValues.front())
272 return yValues.front();
273 if (x >= xValues.back())
274 return yValues.back();
276 size_t segIdx = findSegmentIndex_(xValues, scalarValue(x));
278 return eval(xValues, yValues, x, segIdx);
281 template <
class Evaluation>
282 static Evaluation evalDescending_(
const ValueVector& xValues,
283 const ValueVector& yValues,
286 OPM_TIMEFUNCTION_LOCAL();
287 if (x >= xValues.front())
288 return yValues.front();
289 if (x <= xValues.back())
290 return yValues.back();
292 size_t segIdx = findSegmentIndexDescending_(xValues, scalarValue(x));
294 return eval(xValues, yValues, x, segIdx);
297 template <
class Evaluation>
298 static Evaluation evalDeriv_(
const ValueVector& xValues,
299 const ValueVector& yValues,
302 OPM_TIMEFUNCTION_LOCAL();
303 if (x <= xValues.front())
305 if (x >= xValues.back())
308 size_t segIdx = findSegmentIndex_(xValues, scalarValue(x));
310 Scalar x0 = xValues[segIdx];
311 Scalar x1 = xValues[segIdx + 1];
313 Scalar y0 = yValues[segIdx];
314 Scalar y1 = yValues[segIdx + 1];
316 return (y1 - y0)/(x1 - x0);
318 template<
class ScalarT>
319 static size_t findSegmentIndex_(
const ValueVector& xValues,
const ScalarT& x)
321 OPM_TIMEFUNCTION_LOCAL();
322 assert(xValues.size() > 1);
323 size_t n = xValues.size() - 1;
324 if (xValues.back() <= x)
326 else if (x <= xValues.front())
330 size_t lowIdx = 0, highIdx = n;
331 while (lowIdx + 1 < highIdx) {
332 size_t curIdx = (lowIdx + highIdx)/2;
333 if (xValues[curIdx] < x)
342 static size_t findSegmentIndexDescending_(
const ValueVector& xValues,
Scalar x)
344 OPM_TIMEFUNCTION_LOCAL();
345 assert(xValues.size() > 1);
346 size_t n = xValues.size() - 1;
347 if (x <= xValues.back())
349 else if (xValues.front() <= x)
353 size_t lowIdx = 0, highIdx = n;
354 while (lowIdx + 1 < highIdx) {
355 size_t curIdx = (lowIdx + highIdx)/2;
356 if (xValues[curIdx] >= x)
Specification of the material parameters for a two-phase material law which uses a table and piecewis...
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 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:210
static Evaluation twoPhaseSatPcnw(const Params ¶ms, const Evaluation &Sw)
The saturation-capillary pressure curve.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:144
static Evaluation krw(const Params ¶ms, const FluidState &fs)
The relative permeability for the wetting phase of the porous medium.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:182
static constexpr bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition PiecewiseLinearTwoPhaseMaterial.hpp:75
static void relativePermeabilities(Container &values, const Params ¶ms, const FluidState &fs)
The relative permeabilities.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:118
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:170
static constexpr bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:91
static constexpr bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
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:111
static constexpr bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:87
static constexpr bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:83
static constexpr bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:71
static Evaluation pcnw(const Params ¶ms, const FluidState &fs)
The capillary pressure-saturation curve.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:131
static void capillaryPressures(Container &values, const Params ¶ms, const FluidState &fs)
The capillary pressure-saturation curve.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:97
static Evaluation Sw(const Params &, const FluidState &)
The saturation-capillary pressure curve.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:158
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30