27#ifndef OPM_TWO_PHASE_LET_CURVES_HPP
28#define OPM_TWO_PHASE_LET_CURVES_HPP
48template <
class TraitsT,
class ParamsT = TwoPhaseLETCurvesParams<TraitsT> >
52 using Traits = TraitsT;
53 using Params = ParamsT;
54 using Scalar =
typename Traits::Scalar;
56 static_assert(Traits::numPhases == 2,
57 "The number of fluid phases must be two if you want to use "
58 "this material law!");
60 static constexpr Scalar eps = 1.0e-10;
92 template <
class Container,
class Flu
idState>
95 throw std::logic_error(
"The capillaryPressures(fs) method is not yet implemented");
102 template <
class Container,
class Flu
idState>
103 static void saturations(Container& ,
const Params& ,
const FluidState& )
105 throw std::logic_error(
"The saturations(fs) method is not yet implemented");
118 template <
class Container,
class Flu
idState>
121 throw std::logic_error(
"The relativePermeabilities(fs) method is not yet implemented");
129 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
130 static Evaluation
pcnw(
const Params& ,
const FluidState& )
132 throw std::logic_error(
"TwoPhaseLETCurves::pcnw"
133 " not implemented!");
137 template <
class Evaluation>
138 static Evaluation twoPhaseSatPcnw(
const Params& params,
const Evaluation& Sw)
140 Evaluation Ss = (Sw-params.Sminpc())/params.dSpc();
142 Ss -= (Opm::decay<Scalar>(Ss));
143 }
else if (Ss > 1.0) {
144 Ss -= (Opm::decay<Scalar>(Ss-1.0));
147 const Evaluation powS = Opm::pow(Ss,params.Tpc());
148 const Evaluation pow1mS = Opm::pow(1.0-Ss,params.Lpc());
150 const Evaluation F = pow1mS/(pow1mS+powS*params.Epc());
151 Evaluation tmp = params.Pct()+(params.Pcir()-params.Pct())*F;
156 template <
class Evaluation>
157 static Evaluation twoPhaseSatPcnwInv(
const Params& ,
const Evaluation&)
159 throw std::logic_error(
"TwoPhaseLETCurves::twoPhaseSatPcnwInv"
160 " not implemented!");
163 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
164 static Evaluation Sw(
const Params& ,
const FluidState& )
166 throw std::logic_error(
"The Sw(fs) method is not yet implemented");
169 template <
class Evaluation>
170 static Evaluation twoPhaseSatSw(
const Params& ,
const Evaluation& )
172 throw std::logic_error(
"The twoPhaseSatSw(fs) method is not yet implemented");
175 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
176 static Evaluation Sn(
const Params& ,
const FluidState& )
178 throw std::logic_error(
"The Sn(fs) method is not yet implemented");
181 template <
class Evaluation>
182 static Evaluation twoPhaseSatSn(
const Params& ,
const Evaluation& )
184 throw std::logic_error(
"The twoPhaseSatSn(fs) method is not yet implemented");
192 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
193 static Evaluation
krw(
const Params& ,
const FluidState& )
195 throw std::logic_error(
"TwoPhaseLETCurves::krw"
196 " not implemented!");
199 template <
class Evaluation>
200 static Evaluation twoPhaseSatKrw(
const Params& params,
const Evaluation& Sw)
202 return twoPhaseSatKrLET(Params::wIdx, params, Sw);
205 template <
class Evaluation>
206 static Evaluation twoPhaseSatKrLET(
const unsigned phaseIdx,
const Params& params,
const Evaluation& S)
208 Evaluation Ss = (S-params.Smin(phaseIdx))/params.dS(phaseIdx);
210 Ss -= (Opm::decay<Scalar>(Ss));
211 }
else if (Ss > 1.0) {
212 Ss -= (Opm::decay<Scalar>(Ss-1.0));
215 const Evaluation powS = Opm::pow(Ss,params.L(phaseIdx));
216 const Evaluation pow1mS = Opm::pow(1.0-Ss,params.T(phaseIdx));
218 const Evaluation tmp = params.Krt(phaseIdx)*powS/(powS+pow1mS*params.E(phaseIdx));
223 template <
class Evaluation>
224 static Evaluation twoPhaseSatKrwInv(
const Params& ,
const Evaluation& )
226 throw std::logic_error(
"TwoPhaseLETCurves::twoPhaseSatKrwInv"
227 " not implemented!");
236 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
237 static Evaluation
krn(
const Params& ,
const FluidState& )
239 throw std::logic_error(
"TwoPhaseLETCurves::krn"
240 " not implemented!");
243 template <
class Evaluation>
244 static Evaluation twoPhaseSatKrn(
const Params& params,
const Evaluation& Sw)
246 const Evaluation Sn = 1.0 - Sw;
248 return twoPhaseSatKrLET(Params::nwIdx, params, Sn);
251 template <
class Evaluation>
252 static Evaluation twoPhaseSatKrnInv(
const Params& params,
const Evaluation&
krn)
258 for (
int i = 0; i < 20; ++i) {
259 Evaluation f =
krn - twoPhaseSatKrn(params, Sw);
260 if (Opm::abs(f) < 1e-10)
262 Evaluation fStar =
krn - twoPhaseSatKrn(params, Sw + eps);
263 Evaluation fPrime = (fStar - f)/eps;
264 Evaluation delta = f/fPrime;
271 if (Opm::abs(delta) < 1e-10)
277 Evaluation fL =
krn - twoPhaseSatKrn(params, SL);
278 if (Opm::abs(fL) < eps)
281 Evaluation fR =
krn - twoPhaseSatKrn(params, SR);
282 if (Opm::abs(fR) < eps)
285 for (
int i = 0; i < 50; ++i) {
287 if (abs(SR-SL) < eps)
289 Evaluation fw =
krn - twoPhaseSatKrn(params, Sw);
290 if (Opm::abs(fw) < eps)
295 }
else if (fw * fL > 0) {
303 throw NumericalProblem(
"Couldn't invert the TwoPhaseLETCurves non-wetting phase"
304 " relperm within 20 newton iterations and 50 bisection iterations");
Provides the OPM specific exception classes.
Implementation of the LET curve saturation functions.
Definition TwoPhaseLETCurves.hpp:50
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:237
static Evaluation pcnw(const Params &, const FluidState &)
The capillary pressure-saturation curve.
Definition TwoPhaseLETCurves.hpp:130
static constexpr bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition TwoPhaseLETCurves.hpp:71
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:193
static constexpr bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition TwoPhaseLETCurves.hpp:75
static constexpr int numPhases
The number of fluid phases to which this material law applies.
Definition TwoPhaseLETCurves.hpp:63
static constexpr bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition TwoPhaseLETCurves.hpp:87
static constexpr bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition TwoPhaseLETCurves.hpp:83
static void saturations(Container &, const Params &, const FluidState &)
Calculate the saturations of the phases starting from their pressure differences.
Definition TwoPhaseLETCurves.hpp:103
static constexpr bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition TwoPhaseLETCurves.hpp:79
static void capillaryPressures(Container &, const Params &, const FluidState &)
The capillary pressure-saturation curves.
Definition TwoPhaseLETCurves.hpp:93
static constexpr bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition TwoPhaseLETCurves.hpp:67
static void relativePermeabilities(Container &, const Params &, const FluidState &)
The relative permeability-saturation curves.
Definition TwoPhaseLETCurves.hpp:119
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30