27 #ifndef OPM_TWO_PHASE_LET_CURVES_HPP
28 #define OPM_TWO_PHASE_LET_CURVES_HPP
47 template <
class TraitsT,
class ParamsT = TwoPhaseLETCurvesParams<TraitsT> >
51 using Traits = TraitsT;
52 using Params = ParamsT;
53 using Scalar =
typename Traits::Scalar;
55 static_assert(Traits::numPhases == 2,
56 "The number of fluid phases must be two if you want to use "
57 "this material law!");
59 static constexpr Scalar eps = 1.0e-10;
91 template <
class Container,
class Flu
idState>
94 throw std::logic_error(
"The capillaryPressures(fs) method is not yet implemented");
101 template <
class Container,
class Flu
idState>
102 static void saturations(Container& ,
const Params& ,
const FluidState& )
104 throw std::logic_error(
"The saturations(fs) method is not yet implemented");
117 template <
class Container,
class Flu
idState>
120 throw std::logic_error(
"The relativePermeabilities(fs) method is not yet implemented");
128 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
129 static Evaluation
pcnw(
const Params& ,
const FluidState& )
131 throw std::logic_error(
"TwoPhaseLETCurves::pcnw"
132 " not implemented!");
136 template <
class Evaluation>
137 static Evaluation twoPhaseSatPcnw(
const Params& params,
const Evaluation& Sw)
139 Evaluation Ss = (Sw-params.Sminpc())/params.dSpc();
141 Ss -= (Opm::decay<Scalar>(Ss));
142 }
else if (Ss > 1.0) {
143 Ss -= (Opm::decay<Scalar>(Ss-1.0));
146 const Evaluation powS = Opm::pow(Ss,params.Tpc());
147 const Evaluation pow1mS = Opm::pow(1.0-Ss,params.Lpc());
149 const Evaluation F = pow1mS/(pow1mS+powS*params.Epc());
150 Evaluation tmp = params.Pct()+(params.Pcir()-params.Pct())*F;
155 template <
class Evaluation>
156 static Evaluation twoPhaseSatPcnwInv(
const Params& ,
const Evaluation&)
158 throw std::logic_error(
"TwoPhaseLETCurves::twoPhaseSatPcnwInv"
159 " not implemented!");
162 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
163 static Evaluation Sw(
const Params& ,
const FluidState& )
165 throw std::logic_error(
"The Sw(fs) method is not yet implemented");
168 template <
class Evaluation>
169 static Evaluation twoPhaseSatSw(
const Params& ,
const Evaluation& )
171 throw std::logic_error(
"The twoPhaseSatSw(fs) method is not yet implemented");
174 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
175 static Evaluation Sn(
const Params& ,
const FluidState& )
177 throw std::logic_error(
"The Sn(fs) method is not yet implemented");
180 template <
class Evaluation>
181 static Evaluation twoPhaseSatSn(
const Params& ,
const Evaluation& )
183 throw std::logic_error(
"The twoPhaseSatSn(fs) method is not yet implemented");
191 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
192 static Evaluation
krw(
const Params& ,
const FluidState& )
194 throw std::logic_error(
"TwoPhaseLETCurves::krw"
195 " not implemented!");
198 template <
class Evaluation>
199 static Evaluation twoPhaseSatKrw(
const Params& params,
const Evaluation& Sw)
201 return twoPhaseSatKrLET(Params::wIdx, params, Sw);
204 template <
class Evaluation>
205 static Evaluation twoPhaseSatKrLET(
const unsigned phaseIdx,
const Params& params,
const Evaluation& S)
207 Evaluation Ss = (S-params.Smin(phaseIdx))/params.dS(phaseIdx);
209 Ss -= (Opm::decay<Scalar>(Ss));
210 }
else if (Ss > 1.0) {
211 Ss -= (Opm::decay<Scalar>(Ss-1.0));
214 const Evaluation powS = Opm::pow(Ss,params.L(phaseIdx));
215 const Evaluation pow1mS = Opm::pow(1.0-Ss,params.T(phaseIdx));
217 const Evaluation tmp = params.Krt(phaseIdx)*powS/(powS+pow1mS*params.E(phaseIdx));
222 template <
class Evaluation>
223 static Evaluation twoPhaseSatKrwInv(
const Params& ,
const Evaluation& )
225 throw std::logic_error(
"TwoPhaseLETCurves::twoPhaseSatKrwInv"
226 " not implemented!");
235 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
236 static Evaluation
krn(
const Params& ,
const FluidState& )
238 throw std::logic_error(
"TwoPhaseLETCurves::krn"
239 " not implemented!");
242 template <
class Evaluation>
243 static Evaluation twoPhaseSatKrn(
const Params& params,
const Evaluation& Sw)
245 const Evaluation Sn = 1.0 - Sw;
247 return twoPhaseSatKrLET(Params::nwIdx, params, Sn);
250 template <
class Evaluation>
251 static Evaluation twoPhaseSatKrnInv(
const Params& params,
const Evaluation&
krn)
257 for (
int i = 0; i < 20; ++i) {
258 Evaluation f =
krn - twoPhaseSatKrn(params, Sw);
259 if (Opm::abs(f) < 1e-10)
261 Evaluation fStar =
krn - twoPhaseSatKrn(params, Sw + eps);
262 Evaluation fPrime = (fStar - f)/eps;
263 Evaluation delta = f/fPrime;
270 if (Opm::abs(delta) < 1e-10)
276 Evaluation fL =
krn - twoPhaseSatKrn(params, SL);
277 if (Opm::abs(fL) < eps)
280 Evaluation fR =
krn - twoPhaseSatKrn(params, SR);
281 if (Opm::abs(fR) < eps)
284 for (
int i = 0; i < 50; ++i) {
286 if (abs(SR-SL) < eps)
288 Evaluation fw =
krn - twoPhaseSatKrn(params, Sw);
289 if (Opm::abs(fw) < eps)
294 }
else if (fw * fL > 0) {
302 throw NumericalIssue(
"Couldn't invert the TwoPhaseLETCurves non-wetting phase"
303 " relperm within 20 newton iterations and 50 bisection iterations");
Provides the opm-material specific exception classes.
Implementation of the LET curve saturation functions.
Definition: TwoPhaseLETCurves.hpp:49
static Evaluation krn(const Params &, const FluidState &)
The relative permeability for the non-wetting phase of the medium as implied by the LET parameterizat...
Definition: TwoPhaseLETCurves.hpp:236
static Evaluation pcnw(const Params &, const FluidState &)
The capillary pressure-saturation curve.
Definition: TwoPhaseLETCurves.hpp:129
static constexpr bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition: TwoPhaseLETCurves.hpp:70
static Evaluation krw(const Params &, const FluidState &)
The relative permeability for the wetting phase of the medium implied by the LET parameterization.
Definition: TwoPhaseLETCurves.hpp:192
static constexpr bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition: TwoPhaseLETCurves.hpp:74
static constexpr int numPhases
The number of fluid phases to which this material law applies.
Definition: TwoPhaseLETCurves.hpp:62
static constexpr bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition: TwoPhaseLETCurves.hpp:86
static constexpr bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition: TwoPhaseLETCurves.hpp:82
static void saturations(Container &, const Params &, const FluidState &)
Calculate the saturations of the phases starting from their pressure differences.
Definition: TwoPhaseLETCurves.hpp:102
static constexpr bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition: TwoPhaseLETCurves.hpp:78
static void capillaryPressures(Container &, const Params &, const FluidState &)
The capillary pressure-saturation curves.
Definition: TwoPhaseLETCurves.hpp:92
static constexpr bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition: TwoPhaseLETCurves.hpp:66
static void relativePermeabilities(Container &, const Params &, const FluidState &)
The relative permeability-saturation curves.
Definition: TwoPhaseLETCurves.hpp:118