26 #ifndef OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP
27 #define OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP
31 #include <opm/common/ErrorMacros.hpp>
32 #include <opm/common/Exceptions.hpp>
50 template <
class TraitsT,
class ParamsT = PiecewiseLinearTwoPhaseMaterialParams<TraitsT> >
53 typedef typename ParamsT::ValueVector ValueVector;
63 typedef typename Traits::Scalar
Scalar;
67 static_assert(numPhases == 2,
68 "The piecewise linear two-phase capillary pressure law only"
69 "applies to the case of two fluid phases");
73 static const bool implementsTwoPhaseApi =
true;
77 static const bool implementsTwoPhaseSatApi =
true;
81 static const bool isSaturationDependent =
true;
85 static const bool isPressureDependent =
false;
89 static const bool isTemperatureDependent =
false;
93 static const bool isCompositionDependent =
false;
98 template <
class Container,
class Flu
idState>
99 static void capillaryPressures(Container &values,
const Params ¶ms,
const FluidState &fs)
101 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
103 values[Traits::wettingPhaseIdx] = 0.0;
104 values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fs);
111 template <
class Container,
class Flu
idState>
112 static void saturations(Container& ,
const Params& ,
const FluidState& )
113 { OPM_THROW(std::logic_error,
"Not implemented: saturations()"); }
118 template <
class Container,
class Flu
idState>
119 static void relativePermeabilities(Container &values,
const Params ¶ms,
const FluidState &fs)
121 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
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 ¶ms,
const FluidState &fs)
133 typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
135 FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
137 return twoPhaseSatPcnw(params, Sw);
143 template <
class Evaluation>
144 static Evaluation twoPhaseSatPcnw(
const Params ¶ms,
const Evaluation& Sw)
145 {
return eval_(params.SwPcwnSamples(), params.pcnwSamples(), Sw); }
147 template <
class Evaluation>
148 static Evaluation twoPhaseSatPcnwInv(
const Params ¶ms,
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 ¶ms,
const FluidState &fs)
168 {
return 1 - Sw<FluidState, Scalar>(params, fs); }
170 template <
class Evaluation>
171 static Evaluation twoPhaseSatSn(
const Params ¶ms,
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 ¶ms,
const FluidState &fs)
181 typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
183 FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
185 return twoPhaseSatKrw(params, Sw);
188 template <
class Evaluation>
189 static Evaluation twoPhaseSatKrw(
const Params ¶ms,
const Evaluation& Sw)
190 {
return eval_(params.SwKrwSamples(), params.krwSamples(), Sw); }
192 template <
class Evaluation>
193 static Evaluation twoPhaseSatKrwInv(
const Params ¶ms,
const Evaluation& krw)
194 {
return eval_(params.krwSamples(), params.SwKrwSamples(), krw); }
200 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
201 static Evaluation krn(
const Params ¶ms,
const FluidState &fs)
203 typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
205 FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
207 return twoPhaseSatKrn(params, Sw);
210 template <
class Evaluation>
211 static Evaluation twoPhaseSatKrn(
const Params ¶ms,
const Evaluation& Sw)
212 {
return eval_(params.SwKrnSamples(), params.krnSamples(), Sw); }
214 template <
class Evaluation>
215 static Evaluation twoPhaseSatKrnInv(
const Params ¶ms,
const Evaluation& krn)
216 {
return eval_(params.krnSamples(), params.SwKrnSamples(), krn); }
219 template <
class Evaluation>
220 static Evaluation eval_(
const ValueVector &xValues,
221 const ValueVector &yValues,
224 if (xValues.front() < xValues.back())
225 return evalAscending_(xValues, yValues, x);
226 return evalDescending_(xValues, yValues, x);
229 template <
class Evaluation>
230 static Evaluation evalAscending_(
const ValueVector &xValues,
231 const ValueVector &yValues,
234 typedef MathToolbox<Evaluation> Toolbox;
236 if (x <= xValues.front())
237 return yValues.front();
238 if (x >= xValues.back())
239 return yValues.back();
241 size_t segIdx = findSegmentIndex_(xValues, Toolbox::value(x));
243 Scalar x0 = xValues[segIdx];
244 Scalar x1 = xValues[segIdx + 1];
246 Scalar y0 = yValues[segIdx];
247 Scalar y1 = yValues[segIdx + 1];
249 Scalar m = (y1 - y0)/(x1 - x0);
251 return y0 + (x - x0)*m;
254 template <
class Evaluation>
255 static Evaluation evalDescending_(
const ValueVector &xValues,
256 const ValueVector &yValues,
259 typedef MathToolbox<Evaluation> Toolbox;
261 if (x >= xValues.front())
262 return yValues.front();
263 if (x <= xValues.back())
264 return yValues.back();
266 size_t segIdx = findSegmentIndexDescending_(xValues, Toolbox::value(x));
268 Scalar x0 = xValues[segIdx];
269 Scalar x1 = xValues[segIdx + 1];
271 Scalar y0 = yValues[segIdx];
272 Scalar y1 = yValues[segIdx + 1];
274 Scalar m = (y1 - y0)/(x1 - x0);
276 return y0 + (x - x0)*m;
279 template <
class Evaluation>
280 static Evaluation evalDeriv_(
const ValueVector &xValues,
281 const ValueVector &yValues,
284 typedef MathToolbox<Evaluation> Toolbox;
286 if (Toolbox::value(x) <= xValues.front())
287 return Toolbox::createConstant(0.0);
288 if (Toolbox::value(x) >= xValues.back())
289 return Toolbox::createConstant(0.0);
291 size_t segIdx = findSegmentIndex_(xValues, Toolbox::value(x));
293 Scalar x0 = xValues[segIdx];
294 Scalar x1 = xValues[segIdx + 1];
296 Scalar y0 = yValues[segIdx];
297 Scalar y1 = yValues[segIdx + 1];
299 return Toolbox::createConstant((y1 - y0)/(x1 - x0));
302 static size_t findSegmentIndex_(
const ValueVector &xValues, Scalar x)
304 assert(xValues.size() > 1);
305 size_t n = xValues.size() - 1;
306 if (xValues.back() <= x)
308 else if (x <= xValues.front())
312 size_t lowIdx = 0, highIdx = n;
313 while (lowIdx + 1 < highIdx) {
314 size_t curIdx = (lowIdx + highIdx)/2;
315 if (xValues[curIdx] < x)
324 static size_t findSegmentIndexDescending_(
const ValueVector &xValues, Scalar x)
326 assert(xValues.size() > 1);
327 size_t n = xValues.size() - 1;
328 if (x <= xValues.back())
330 else if (xValues.front() <= x)
334 size_t lowIdx = 0, highIdx = n;
335 while (lowIdx + 1 < highIdx) {
336 size_t curIdx = (lowIdx + highIdx)/2;
337 if (xValues[curIdx] >= x)
TraitsT Traits
The traits class for this material law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:57
ParamsT Params
The type of the parameter objects for this law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:60
Definition: Air_Mesitylene.hpp:31
static const int numPhases
The number of fluid phases.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:66
Traits::Scalar Scalar
The type of the scalar values for this law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:63
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