27 #ifndef OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP
28 #define OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP
49 template <
class TraitsT,
class ParamsT = PiecewiseLinearTwoPhaseMaterialParams<TraitsT> >
52 typedef typename ParamsT::ValueVector ValueVector;
62 typedef typename Traits::Scalar
Scalar;
67 "The piecewise linear two-phase capillary pressure law only"
68 "applies to the case of two fluid phases");
97 template <
class Container,
class Flu
idState>
100 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
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 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
122 values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fs);
123 values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fs);
129 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
130 static Evaluation
pcnw(
const Params& params,
const FluidState& fs)
133 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
141 template <
class Evaluation>
143 {
return eval_(params.SwPcwnSamples(), params.pcnwSamples(),
Sw); }
145 template <
class Evaluation>
146 static Evaluation twoPhaseSatPcnwInv(
const Params& params,
const Evaluation&
pcnw)
147 {
return eval_(params.pcnwSamples(), params.SwPcwnSamples(),
pcnw); }
152 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
153 static Evaluation
Sw(
const Params& ,
const FluidState& )
154 {
throw std::logic_error(
"Not implemented: Sw()"); }
156 template <
class Evaluation>
157 static Evaluation twoPhaseSatSw(
const Params& ,
const Evaluation& )
158 {
throw std::logic_error(
"Not implemented: twoPhaseSatSw()"); }
164 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
165 static Evaluation
Sn(
const Params& params,
const FluidState& fs)
166 {
return 1 - Sw<FluidState, Scalar>(params, fs); }
168 template <
class Evaluation>
169 static Evaluation twoPhaseSatSn(
const Params& params,
const Evaluation& pC)
170 {
return 1 - twoPhaseSatSw(params, pC); }
176 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
177 static Evaluation
krw(
const Params& params,
const FluidState& fs)
180 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
182 return twoPhaseSatKrw(params,
Sw);
185 template <
class Evaluation>
186 static Evaluation twoPhaseSatKrw(
const Params& params,
const Evaluation&
Sw)
187 {
return eval_(params.SwKrwSamples(), params.krwSamples(),
Sw); }
189 template <
class Evaluation>
190 static Evaluation twoPhaseSatKrwInv(
const Params& params,
const Evaluation&
krw)
191 {
return eval_(params.krwSamples(), params.SwKrwSamples(),
krw); }
197 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
198 static Evaluation
krn(
const Params& params,
const FluidState& fs)
201 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
203 return twoPhaseSatKrn(params,
Sw);
206 template <
class Evaluation>
207 static Evaluation twoPhaseSatKrn(
const Params& params,
const Evaluation&
Sw)
208 {
return eval_(params.SwKrnSamples(), params.krnSamples(),
Sw); }
210 template <
class Evaluation>
211 static Evaluation twoPhaseSatKrnInv(
const Params& params,
const Evaluation&
krn)
212 {
return eval_(params.krnSamples(), params.SwKrnSamples(),
krn); }
215 template <
class Evaluation>
216 static Evaluation eval_(
const ValueVector& xValues,
217 const ValueVector& yValues,
220 if (xValues.front() < xValues.back())
221 return evalAscending_(xValues, yValues, x);
222 return evalDescending_(xValues, yValues, x);
225 template <
class Evaluation>
226 static Evaluation evalAscending_(
const ValueVector& xValues,
227 const ValueVector& yValues,
230 if (x <= xValues.front())
231 return yValues.front();
232 if (x >= xValues.back())
233 return yValues.back();
235 size_t segIdx = findSegmentIndex_(xValues, scalarValue(x));
237 Scalar x0 = xValues[segIdx];
238 Scalar x1 = xValues[segIdx + 1];
240 Scalar y0 = yValues[segIdx];
241 Scalar y1 = yValues[segIdx + 1];
243 Scalar m = (y1 - y0)/(x1 - x0);
245 return y0 + (x - x0)*m;
248 template <
class Evaluation>
249 static Evaluation evalDescending_(
const ValueVector& xValues,
250 const ValueVector& yValues,
253 if (x >= xValues.front())
254 return yValues.front();
255 if (x <= xValues.back())
256 return yValues.back();
258 size_t segIdx = findSegmentIndexDescending_(xValues, scalarValue(x));
260 Scalar x0 = xValues[segIdx];
261 Scalar x1 = xValues[segIdx + 1];
263 Scalar y0 = yValues[segIdx];
264 Scalar y1 = yValues[segIdx + 1];
266 Scalar m = (y1 - y0)/(x1 - x0);
268 return y0 + (x - x0)*m;
271 template <
class Evaluation>
272 static Evaluation evalDeriv_(
const ValueVector& xValues,
273 const ValueVector& yValues,
276 if (x <= xValues.front())
278 if (x >= xValues.back())
281 size_t segIdx = findSegmentIndex_(xValues, scalarValue(x));
283 Scalar x0 = xValues[segIdx];
284 Scalar x1 = xValues[segIdx + 1];
286 Scalar y0 = yValues[segIdx];
287 Scalar y1 = yValues[segIdx + 1];
289 return (y1 - y0)/(x1 - x0);
292 static size_t findSegmentIndex_(
const ValueVector& xValues,
Scalar x)
294 assert(xValues.size() > 1);
295 size_t n = xValues.size() - 1;
296 if (xValues.back() <= x)
298 else if (x <= xValues.front())
302 size_t lowIdx = 0, highIdx = n;
303 while (lowIdx + 1 < highIdx) {
304 size_t curIdx = (lowIdx + highIdx)/2;
305 if (xValues[curIdx] < x)
314 static size_t findSegmentIndexDescending_(
const ValueVector& xValues,
Scalar x)
316 assert(xValues.size() > 1);
317 size_t n = xValues.size() - 1;
318 if (x <= xValues.back())
320 else if (xValues.front() <= x)
324 size_t lowIdx = 0, highIdx = n;
325 while (lowIdx + 1 < highIdx) {
326 size_t curIdx = (lowIdx + highIdx)/2;
327 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:51
static Evaluation krn(const Params ¶ms, const FluidState &fs)
The relative permeability for the non-wetting phase of the porous medium.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:198
static Evaluation twoPhaseSatPcnw(const Params ¶ms, const Evaluation &Sw)
The saturation-capillary pressure curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:142
TraitsT Traits
The traits class for this material law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:56
static Evaluation krw(const Params ¶ms, const FluidState &fs)
The relative permeability for the wetting phase of the porous medium.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:177
static const bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:72
static void relativePermeabilities(Container &values, const Params ¶ms, const FluidState &fs)
The relative permeabilities.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:118
static Evaluation Sn(const Params ¶ms, const FluidState &fs)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:165
static const bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:88
static void saturations(Container &, const Params &, const FluidState &)
The saturations of the fluid phases starting from their pressure differences.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:111
static const bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:92
static const bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:80
static const bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:76
static Evaluation pcnw(const Params ¶ms, const FluidState &fs)
The capillary pressure-saturation curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:130
static const bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:84
static const int numPhases
The number of fluid phases.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:65
Traits::Scalar Scalar
The type of the scalar values for this law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:62
ParamsT Params
The type of the parameter objects for this law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:59
static void capillaryPressures(Container &values, const Params ¶ms, const FluidState &fs)
The capillary pressure-saturation curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:98
static Evaluation Sw(const Params &, const FluidState &)
The saturation-capillary pressure curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:153