27 #ifndef OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP 28 #define OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP 32 #include <opm/common/ErrorMacros.hpp> 33 #include <opm/common/Exceptions.hpp> 51 template <
class TraitsT,
class ParamsT = PiecewiseLinearTwoPhaseMaterialParams<TraitsT> >
54 typedef typename ParamsT::ValueVector ValueVector;
64 typedef typename Traits::Scalar
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 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
104 values[Traits::wettingPhaseIdx] = 0.0;
105 values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fs);
112 template <
class Container,
class Flu
idState>
114 { OPM_THROW(std::logic_error,
"Not implemented: saturations()"); }
119 template <
class Container,
class Flu
idState>
122 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
124 values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fs);
125 values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fs);
131 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
132 static Evaluation
pcnw(
const Params& params,
const FluidState& fs)
135 Opm::decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
143 template <
class Evaluation>
145 {
return eval_(params.SwPcwnSamples(), params.pcnwSamples(),
Sw); }
147 template <
class Evaluation>
148 static Evaluation twoPhaseSatPcnwInv(
const Params& params,
const Evaluation&
pcnw)
149 {
return eval_(params.pcnwSamples(), params.SwPcwnSamples(),
pcnw); }
154 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
155 static Evaluation
Sw(
const Params& ,
const FluidState& )
156 { OPM_THROW(std::logic_error,
"Not implemented: Sw()"); }
158 template <
class Evaluation>
159 static Evaluation twoPhaseSatSw(
const Params& ,
const Evaluation& )
160 { OPM_THROW(std::logic_error,
"Not implemented: twoPhaseSatSw()"); }
166 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
167 static Evaluation
Sn(
const Params& params,
const FluidState& fs)
168 {
return 1 - Sw<FluidState, Scalar>(params, fs); }
170 template <
class Evaluation>
171 static Evaluation twoPhaseSatSn(
const Params& params,
const Evaluation& pC)
172 {
return 1 - twoPhaseSatSw(params, pC); }
178 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
179 static Evaluation
krw(
const Params& params,
const FluidState& fs)
182 Opm::decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
184 return twoPhaseSatKrw(params,
Sw);
187 template <
class Evaluation>
188 static Evaluation twoPhaseSatKrw(
const Params& params,
const Evaluation&
Sw)
189 {
return eval_(params.SwKrwSamples(), params.krwSamples(),
Sw); }
191 template <
class Evaluation>
192 static Evaluation twoPhaseSatKrwInv(
const Params& params,
const Evaluation&
krw)
193 {
return eval_(params.krwSamples(), params.SwKrwSamples(),
krw); }
199 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
200 static Evaluation
krn(
const Params& params,
const FluidState& fs)
203 Opm::decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
205 return twoPhaseSatKrn(params,
Sw);
208 template <
class Evaluation>
209 static Evaluation twoPhaseSatKrn(
const Params& params,
const Evaluation&
Sw)
210 {
return eval_(params.SwKrnSamples(), params.krnSamples(),
Sw); }
212 template <
class Evaluation>
213 static Evaluation twoPhaseSatKrnInv(
const Params& params,
const Evaluation&
krn)
214 {
return eval_(params.krnSamples(), params.SwKrnSamples(),
krn); }
217 template <
class Evaluation>
218 static Evaluation eval_(
const ValueVector& xValues,
219 const ValueVector& yValues,
222 if (xValues.front() < xValues.back())
223 return evalAscending_(xValues, yValues, x);
224 return evalDescending_(xValues, yValues, x);
227 template <
class Evaluation>
228 static Evaluation evalAscending_(
const ValueVector& xValues,
229 const ValueVector& yValues,
232 if (x <= xValues.front())
233 return yValues.front();
234 if (x >= xValues.back())
235 return yValues.back();
237 size_t segIdx = findSegmentIndex_(xValues, Opm::scalarValue(x));
239 Scalar x0 = xValues[segIdx];
240 Scalar x1 = xValues[segIdx + 1];
242 Scalar y0 = yValues[segIdx];
243 Scalar y1 = yValues[segIdx + 1];
245 Scalar m = (y1 - y0)/(x1 - x0);
247 return y0 + (x - x0)*m;
250 template <
class Evaluation>
251 static Evaluation evalDescending_(
const ValueVector& xValues,
252 const ValueVector& yValues,
255 if (x >= xValues.front())
256 return yValues.front();
257 if (x <= xValues.back())
258 return yValues.back();
260 size_t segIdx = findSegmentIndexDescending_(xValues, Opm::scalarValue(x));
262 Scalar x0 = xValues[segIdx];
263 Scalar x1 = xValues[segIdx + 1];
265 Scalar y0 = yValues[segIdx];
266 Scalar y1 = yValues[segIdx + 1];
268 Scalar m = (y1 - y0)/(x1 - x0);
270 return y0 + (x - x0)*m;
273 template <
class Evaluation>
274 static Evaluation evalDeriv_(
const ValueVector& xValues,
275 const ValueVector& yValues,
278 if (x <= xValues.front())
280 if (x >= xValues.back())
283 size_t segIdx = findSegmentIndex_(xValues, Opm::scalarValue(x));
285 Scalar x0 = xValues[segIdx];
286 Scalar x1 = xValues[segIdx + 1];
288 Scalar y0 = yValues[segIdx];
289 Scalar y1 = yValues[segIdx + 1];
291 return (y1 - y0)/(x1 - x0);
294 static size_t findSegmentIndex_(
const ValueVector& xValues,
Scalar x)
296 assert(xValues.size() > 1);
297 size_t n = xValues.size() - 1;
298 if (xValues.back() <= x)
300 else if (x <= xValues.front())
304 size_t lowIdx = 0, highIdx = n;
305 while (lowIdx + 1 < highIdx) {
306 size_t curIdx = (lowIdx + highIdx)/2;
307 if (xValues[curIdx] < x)
316 static size_t findSegmentIndexDescending_(
const ValueVector& xValues,
Scalar x)
318 assert(xValues.size() > 1);
319 size_t n = xValues.size() - 1;
320 if (x <= xValues.back())
322 else if (xValues.front() <= x)
326 size_t lowIdx = 0, highIdx = n;
327 while (lowIdx + 1 < highIdx) {
328 size_t curIdx = (lowIdx + highIdx)/2;
329 if (xValues[curIdx] >= x)
static Evaluation Sn(const Params ¶ms, const FluidState &fs)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:167
static const bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition...
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:94
static Evaluation krw(const Params ¶ms, const FluidState &fs)
The relative permeability for the wetting phase of the porous medium.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:179
static const bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:78
static const bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:82
static void saturations(Container &, const Params &, const FluidState &)
The saturations of the fluid phases starting from their pressure differences.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:113
ParamsT Params
The type of the parameter objects for this law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:61
Definition: Air_Mesitylene.hpp:33
static const bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:90
static Evaluation pcnw(const Params ¶ms, const FluidState &fs)
The capillary pressure-saturation curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:132
static Evaluation krn(const Params ¶ms, const FluidState &fs)
The relative permeability for the non-wetting phase of the porous medium.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:200
Implementation of a tabulated, piecewise linear capillary pressure law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:52
static Evaluation twoPhaseSatPcnw(const Params ¶ms, const Evaluation &Sw)
The saturation-capillary pressure curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:144
static Evaluation Sw(const Params &, const FluidState &)
The saturation-capillary pressure curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:155
static void capillaryPressures(Container &values, const Params ¶ms, const FluidState &fs)
The capillary pressure-saturation curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:100
static void relativePermeabilities(Container &values, const Params ¶ms, const FluidState &fs)
The relative permeabilities.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:120
static const int numPhases
The number of fluid phases.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:67
Specification of the material parameters for a two-phase material law which uses a table and piecewis...
Traits::Scalar Scalar
The type of the scalar values for this law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:64
static const bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure...
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:86
static const bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:74
TraitsT Traits
The traits class for this material law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:58