27#ifndef OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP
28#define OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP
31#include <opm/common/TimingMacros.hpp>
38#include <opm/common/utility/gpuDecorators.hpp>
51template <
class TraitsT,
class ParamsT = PiecewiseLinearTwoPhaseMaterialParams<TraitsT> >
54 using ValueVector =
typename ParamsT::ValueVector;
64 using Scalar =
typename Traits::Scalar;
69 "The piecewise linear two-phase capillary pressure law only"
70 "applies to the case of two fluid phases");
99 template <
class Container,
class Flu
idState>
102 OPM_TIMEFUNCTION_LOCAL();
103 using Evaluation =
typename std::remove_reference<
decltype(values[0])>::type;
105 values[Traits::wettingPhaseIdx] = 0.0;
106 values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fs);
113 template <
class Container,
class Flu
idState>
115 {
throw std::logic_error(
"Not implemented: saturations()"); }
120 template <
class Container,
class Flu
idState>
123 OPM_TIMEFUNCTION_LOCAL();
124 using Evaluation =
typename std::remove_reference<
decltype(values[0])>::type;
126 values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fs);
127 values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fs);
133 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
134 OPM_HOST_DEVICE
static Evaluation
pcnw(
const Params& params,
const FluidState& fs)
136 OPM_TIMEFUNCTION_LOCAL();
138 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
146 template <
class Evaluation>
149 OPM_TIMEFUNCTION_LOCAL();
150 return eval_(params.SwPcwnSamples(), params.pcwnSamples(),
Sw);
153 template <
class Evaluation>
154 static Evaluation twoPhaseSatPcnwInv(
const Params& params,
const Evaluation&
pcnw)
155 {
return eval_(params.pcwnSamples(), params.SwPcwnSamples(),
pcnw); }
160 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
161 OPM_HOST_DEVICE
static Evaluation
Sw(
const Params& ,
const FluidState& )
162 {
throw std::logic_error(
"Not implemented: Sw()"); }
164 template <
class Evaluation>
165 OPM_HOST_DEVICE
static Evaluation twoPhaseSatSw(
const Params& ,
const Evaluation& )
166 {
throw std::logic_error(
"Not implemented: twoPhaseSatSw()"); }
172 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
173 OPM_HOST_DEVICE
static Evaluation
Sn(
const Params& params,
const FluidState& fs)
174 {
return 1 - Sw<FluidState, Scalar>(params, fs); }
176 template <
class Evaluation>
177 OPM_HOST_DEVICE
static Evaluation twoPhaseSatSn(
const Params& params,
const Evaluation& pC)
178 {
return 1 - twoPhaseSatSw(params, pC); }
184 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
185 OPM_HOST_DEVICE
static Evaluation
krw(
const Params& params,
const FluidState& fs)
187 OPM_TIMEFUNCTION_LOCAL();
189 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
191 return twoPhaseSatKrw(params, sw);
194 template <
class Evaluation>
195 OPM_HOST_DEVICE
static Evaluation twoPhaseSatKrw(
const Params& params,
const Evaluation&
Sw)
197 OPM_TIMEFUNCTION_LOCAL();
198 return eval_(params.SwKrwSamples(), params.krwSamples(),
Sw);
201 template <
class Evaluation>
202 OPM_HOST_DEVICE
static Evaluation twoPhaseSatKrwInv(
const Params& params,
const Evaluation&
krw)
204 OPM_TIMEFUNCTION_LOCAL();
205 return eval_(params.krwSamples(), params.SwKrwSamples(),
krw);
212 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
213 OPM_HOST_DEVICE
static Evaluation
krn(
const Params& params,
const FluidState& fs)
215 OPM_TIMEFUNCTION_LOCAL();
217 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
219 return twoPhaseSatKrn(params, sw);
222 template <
class Evaluation>
223 OPM_HOST_DEVICE
static Evaluation twoPhaseSatKrn(
const Params& params,
const Evaluation&
Sw)
225 OPM_TIMEFUNCTION_LOCAL();
226 return eval_(params.SwKrnSamples(), params.krnSamples(),
Sw);
229 template <
class Evaluation>
230 OPM_HOST_DEVICE
static Evaluation twoPhaseSatKrnInv(
const Params& params,
const Evaluation&
krn)
231 {
return eval_(params.krnSamples(), params.SwKrnSamples(),
krn); }
233 template <
class Evaluation>
234 OPM_HOST_DEVICE
static size_t findSegmentIndex(
const ValueVector& xValues,
const Evaluation& x){
235 return findSegmentIndex_(xValues, scalarValue(x));
238 template <
class Evaluation>
239 OPM_HOST_DEVICE
static size_t findSegmentIndexDescending(
const ValueVector& xValues,
const Evaluation& x){
240 return findSegmentIndexDescending_(xValues, scalarValue(x));
243 template <
class Evaluation>
244 OPM_HOST_DEVICE
static Evaluation eval(
const ValueVector& xValues,
const ValueVector& yValues,
const Evaluation& x,
unsigned segIdx){
245 Scalar x0 = xValues[segIdx];
246 Scalar x1 = xValues[segIdx + 1];
248 Scalar y0 = yValues[segIdx];
249 Scalar y1 = yValues[segIdx + 1];
251 Scalar m = (y1 - y0)/(x1 - x0);
253 return y0 + (x - x0)*m;
257 template <
class Evaluation>
258 OPM_HOST_DEVICE
static Evaluation eval_(
const ValueVector& xValues,
259 const ValueVector& yValues,
262 OPM_TIMEFUNCTION_LOCAL();
263 if (xValues.front() < xValues.back())
264 return evalAscending_(xValues, yValues, x);
265 return evalDescending_(xValues, yValues, x);
268 template <
class Evaluation>
269 OPM_HOST_DEVICE
static Evaluation evalAscending_(
const ValueVector& xValues,
270 const ValueVector& yValues,
273 OPM_TIMEFUNCTION_LOCAL();
274 if (x <= xValues.front())
275 return yValues.front();
276 if (x >= xValues.back())
277 return yValues.back();
279 size_t segIdx = findSegmentIndex_(xValues, scalarValue(x));
281 return eval(xValues, yValues, x, segIdx);
284 template <
class Evaluation>
285 OPM_HOST_DEVICE
static Evaluation evalDescending_(
const ValueVector& xValues,
286 const ValueVector& yValues,
289 OPM_TIMEFUNCTION_LOCAL();
290 if (x >= xValues.front())
291 return yValues.front();
292 if (x <= xValues.back())
293 return yValues.back();
295 size_t segIdx = findSegmentIndexDescending_(xValues, scalarValue(x));
297 return eval(xValues, yValues, x, segIdx);
300 template <
class Evaluation>
301 OPM_HOST_DEVICE
static Evaluation evalDeriv_(
const ValueVector& xValues,
302 const ValueVector& yValues,
305 OPM_TIMEFUNCTION_LOCAL();
306 if (x <= xValues.front())
308 if (x >= xValues.back())
311 size_t segIdx = findSegmentIndex_(xValues, scalarValue(x));
313 Scalar x0 = xValues[segIdx];
314 Scalar x1 = xValues[segIdx + 1];
316 Scalar y0 = yValues[segIdx];
317 Scalar y1 = yValues[segIdx + 1];
319 return (y1 - y0)/(x1 - x0);
321 template<
class ScalarT>
322 OPM_HOST_DEVICE
static size_t findSegmentIndex_(
const ValueVector& xValues,
const ScalarT& x)
324 OPM_TIMEFUNCTION_LOCAL();
325 assert(xValues.size() > 1);
326 size_t n = xValues.size() - 1;
327 if (xValues.back() <= x)
329 else if (x <= xValues.front())
333 size_t lowIdx = 0, highIdx = n;
334 while (lowIdx + 1 < highIdx) {
335 size_t curIdx = (lowIdx + highIdx)/2;
336 if (xValues[curIdx] < x)
345 OPM_HOST_DEVICE
static size_t findSegmentIndexDescending_(
const ValueVector& xValues,
Scalar x)
347 OPM_TIMEFUNCTION_LOCAL();
348 assert(xValues.size() > 1);
349 size_t n = xValues.size() - 1;
350 if (x <= xValues.back())
352 else if (xValues.front() <= x)
356 size_t lowIdx = 0, highIdx = n;
357 while (lowIdx + 1 < highIdx) {
358 size_t curIdx = (lowIdx + highIdx)/2;
359 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:53
TraitsT Traits
The traits class for this material law.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:58
ParamsT Params
The type of the parameter objects for this law.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:61
static OPM_HOST_DEVICE void capillaryPressures(Container &values, const Params ¶ms, const FluidState &fs)
The capillary pressure-saturation curve.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:100
static OPM_HOST_DEVICE Evaluation Sn(const Params ¶ms, const FluidState &fs)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:173
static constexpr int numPhases
The number of fluid phases.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:67
static OPM_HOST_DEVICE void relativePermeabilities(Container &values, const Params ¶ms, const FluidState &fs)
The relative permeabilities.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:121
static constexpr bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition PiecewiseLinearTwoPhaseMaterial.hpp:78
static OPM_HOST_DEVICE Evaluation twoPhaseSatPcnw(const Params ¶ms, const Evaluation &Sw)
The saturation-capillary pressure curve.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:147
typename Traits::Scalar Scalar
The type of the scalar values for this law.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:64
static OPM_HOST_DEVICE Evaluation krw(const Params ¶ms, const FluidState &fs)
The relative permeability for the wetting phase of the porous medium.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:185
static constexpr bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:94
static constexpr bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:82
static OPM_HOST_DEVICE void saturations(Container &, const Params &, const FluidState &)
The saturations of the fluid phases starting from their pressure differences.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:114
static OPM_HOST_DEVICE Evaluation krn(const Params ¶ms, const FluidState &fs)
The relative permeability for the non-wetting phase of the porous medium.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:213
static constexpr bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:90
static constexpr bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:86
static OPM_HOST_DEVICE Evaluation Sw(const Params &, const FluidState &)
The saturation-capillary pressure curve.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:161
static constexpr bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:74
static OPM_HOST_DEVICE Evaluation pcnw(const Params ¶ms, const FluidState &fs)
The capillary pressure-saturation curve.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:134
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30