27 #ifndef OPM_BROOKS_COREY_HPP
28 #define OPM_BROOKS_COREY_HPP
33 #include <opm/common/ErrorMacros.hpp>
34 #include <opm/common/Exceptions.hpp>
53 template <
class TraitsT,
class ParamsT = BrooksCoreyParams<TraitsT> >
59 typedef typename Traits::Scalar
Scalar;
63 static_assert(numPhases == 2,
64 "The Brooks-Corey capillary pressure law only applies "
65 "to the case of two fluid phases");
69 static const bool implementsTwoPhaseApi =
true;
73 static const bool implementsTwoPhaseSatApi =
true;
77 static const bool isSaturationDependent =
true;
81 static const bool isPressureDependent =
false;
85 static const bool isTemperatureDependent =
false;
89 static const bool isCompositionDependent =
false;
91 static_assert(Traits::numPhases == 2,
92 "The number of fluid phases must be two if you want to use "
93 "this material law!");
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 &values,
const Params ¶ms,
const FluidState &fs)
114 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
116 values[Traits::wettingPhaseIdx] = Sw<FluidState, Evaluation>(params, fs);
117 values[Traits::nonWettingPhaseIdx] = 1 - values[Traits::wettingPhaseIdx];
130 template <
class Container,
class Flu
idState>
131 static void relativePermeabilities(Container &values,
const Params ¶ms,
const FluidState &fs)
133 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
135 values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fs);
136 values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fs);
152 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
153 static Evaluation pcnw(
const Params ¶ms,
const FluidState &fs)
155 typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
157 const Evaluation& Sw =
158 FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
160 assert(0 <= Sw && Sw <= 1);
162 return twoPhaseSatPcnw(params, Sw);
165 template <
class Evaluation>
166 static Evaluation twoPhaseSatPcnw(
const Params ¶ms,
const Evaluation& Sw)
168 typedef MathToolbox<Evaluation> Toolbox;
170 assert(0 <= Sw && Sw <= 1);
172 return params.entryPressure()*
Toolbox::pow(Sw, -1/params.lambda());
175 template <
class Evaluation>
176 static Evaluation twoPhaseSatPcnwInv(
const Params ¶ms,
const Evaluation& pcnw)
178 typedef MathToolbox<Evaluation> Toolbox;
182 return Toolbox::pow(params.entryPressure()/pcnw, -params.lambda());
197 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
198 static Evaluation Sw(
const Params ¶ms,
const FluidState &fs)
200 typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
203 FsToolbox::template toLhs<Evaluation>(fs.pressure(Traits::nonWettingPhaseIdx))
204 - FsToolbox::template toLhs<Evaluation>(fs.pressure(Traits::wettingPhaseIdx));
205 return twoPhaseSatSw(params, pC);
208 template <
class Evaluation>
209 static Evaluation twoPhaseSatSw(
const Params ¶ms,
const Evaluation& pc)
211 typedef MathToolbox<Evaluation> Toolbox;
215 return Toolbox::pow(pc/params.entryPressure(), -params.lambda());
222 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
223 static Evaluation Sn(
const Params ¶ms,
const FluidState &fs)
224 {
return 1 - Sw<FluidState, Evaluation>(params, fs); }
226 template <
class Evaluation>
227 static Evaluation twoPhaseSatSn(
const Params ¶ms,
const Evaluation& pc)
228 {
return 1 - twoPhaseSatSw(params, pc); }
237 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
238 static Evaluation krw(
const Params ¶ms,
const FluidState &fs)
240 typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
243 FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
245 return twoPhaseSatKrw(params, Sw);
248 template <
class Evaluation>
249 static Evaluation twoPhaseSatKrw(
const Params ¶ms,
const Evaluation& Sw)
251 typedef MathToolbox<Evaluation> Toolbox;
253 assert(0 <= Sw && Sw <= 1);
258 template <
class Evaluation>
259 static Evaluation twoPhaseSatKrwInv(
const Params ¶ms,
const Evaluation& krw)
261 typedef MathToolbox<Evaluation> Toolbox;
263 return Toolbox::pow(krw, 1.0/(2.0/params.lambda() + 3));
274 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
275 static Evaluation krn(
const Params ¶ms,
const FluidState &fs)
277 typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
279 const Evaluation& Sw =
280 1.0 - FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
282 return twoPhaseSatKrn(params, Sw);
285 template <
class Evaluation>
286 static Evaluation twoPhaseSatKrn(
const Params ¶ms,
const Evaluation& Sw)
288 typedef MathToolbox<Evaluation> Toolbox;
290 assert(0 <= Sw && Sw <= 1);
292 Scalar exponent = 2.0/params.lambda() + 1;
293 const Evaluation Sn = 1. - Sw;
297 template <
class Evaluation>
298 static Evaluation twoPhaseSatKrnInv(
const Params ¶ms,
const Evaluation& krn)
300 typedef MathToolbox<Evaluation> Toolbox;
304 Evaluation Sw = Toolbox::createConstant(0.5);
306 for (
int i = 0; i < 20; ++i) {
307 Evaluation f = krn - twoPhaseSatKrn(params, Sw);
308 Evaluation fStar = krn - twoPhaseSatKrn(params, Sw + eps);
309 Evaluation fPrime = (fStar - f)/eps;
311 Evaluation delta = f/fPrime;
314 Sw = Toolbox::createConstant(0.0);
319 OPM_THROW(NumericalProblem,
320 "Couldn't invert the Brooks-Corey non-wetting phase"
321 " relperm within 20 iterations");
326 #endif // BROOKS_COREY_HPP
ParamsT Params
Definition: BrooksCorey.hpp:58
Definition: Air_Mesitylene.hpp:31
TraitsT Traits
Definition: BrooksCorey.hpp:57
Evaluation< Scalar, VarSetTag, numVars > abs(const Evaluation< Scalar, VarSetTag, numVars > &)
Definition: Math.hpp:41
Specification of the material parameters for the Brooks-Corey constitutive relations.
Implementation of the Brooks-Corey capillary pressure <-> saturation relation.
Definition: BrooksCorey.hpp:54
Evaluation< Scalar, VarSetTag, numVars > pow(const Evaluation< Scalar, VarSetTag, numVars > &base, Scalar exp)
Definition: Math.hpp:312
static const int numPhases
The number of fluid phases to which this material law applies.
Definition: BrooksCorey.hpp:62
Traits::Scalar Scalar
Definition: BrooksCorey.hpp:59