27#ifndef OPM_BROOKS_COREY_HPP
28#define OPM_BROOKS_COREY_HPP
53template <
class TraitsT,
class ParamsT = BrooksCoreyParams<TraitsT> >
57 typedef TraitsT Traits;
58 typedef ParamsT Params;
59 typedef typename Traits::Scalar Scalar;
64 "The Brooks-Corey capillary pressure law only applies "
65 "to the case of two fluid phases");
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>
101 typedef typename std::remove_reference<
decltype(values[0])>::type Evaluation;
103 values[Traits::wettingPhaseIdx] = 0.0;
111 template <
class Container,
class Flu
idState>
112 static void saturations(Container& values,
const Params& params,
const FluidState& fs)
114 typedef typename std::remove_reference<
decltype(values[0])>::type Evaluation;
117 values[Traits::nonWettingPhaseIdx] = 1.0 - values[Traits::wettingPhaseIdx];
130 template <
class Container,
class Flu
idState>
133 typedef typename std::remove_reference<
decltype(values[0])>::type Evaluation;
153 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
154 static Evaluation
pcnw(
const Params& params,
const FluidState& fs)
156 const Evaluation&
Sw =
157 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
159 assert(0.0 <=
Sw &&
Sw <= 1.0);
161 return twoPhaseSatPcnw(params,
Sw);
164 template <
class Evaluation>
165 static Evaluation twoPhaseSatPcnw(
const Params& params,
const Evaluation&
Sw)
167 assert(0.0 <=
Sw &&
Sw <= 1.0);
169 return params.entryPressure()*pow(
Sw, -1/params.lambda());
172 template <
class Evaluation>
173 static Evaluation twoPhaseSatPcnwInv(
const Params& params,
const Evaluation&
pcnw)
177 return pow(params.entryPressure()/
pcnw, -params.lambda());
193 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
194 static Evaluation
Sw(
const Params& params,
const FluidState& fs)
197 decay<Evaluation>(fs.pressure(Traits::nonWettingPhaseIdx))
198 - decay<Evaluation>(fs.pressure(Traits::wettingPhaseIdx));
199 return twoPhaseSatSw(params, pC);
202 template <
class Evaluation>
203 static Evaluation twoPhaseSatSw(
const Params& params,
const Evaluation& pc)
207 return pow(pc/params.entryPressure(), -params.lambda());
214 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
215 static Evaluation
Sn(
const Params& params,
const FluidState& fs)
218 template <
class Evaluation>
219 static Evaluation twoPhaseSatSn(
const Params& params,
const Evaluation& pc)
220 {
return 1.0 - twoPhaseSatSw(params, pc); }
230 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
231 static Evaluation
krw(
const Params& params,
const FluidState& fs)
234 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
236 return twoPhaseSatKrw(params, sw);
239 template <
class Evaluation>
240 static Evaluation twoPhaseSatKrw(
const Params& params,
const Evaluation&
Sw)
242 assert(0.0 <=
Sw &&
Sw <= 1.0);
244 return pow(
Sw, 2.0/params.lambda() + 3.0);
247 template <
class Evaluation>
248 static Evaluation twoPhaseSatKrwInv(
const Params& params,
const Evaluation&
krw)
250 return pow(
krw, 1.0/(2.0/params.lambda() + 3.0));
262 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
263 static Evaluation
krn(
const Params& params,
const FluidState& fs)
265 const Evaluation& sw =
266 1.0 - decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
268 return twoPhaseSatKrn(params, sw);
271 template <
class Evaluation>
272 static Evaluation twoPhaseSatKrn(
const Params& params,
const Evaluation&
Sw)
274 assert(0.0 <=
Sw &&
Sw <= 1.0);
276 Scalar exponent = 2.0/params.lambda() + 1.0;
277 const Evaluation sn = 1.0 -
Sw;
278 return sn*sn*(1. - pow(
Sw, exponent));
281 template <
class Evaluation>
282 static Evaluation twoPhaseSatKrnInv(
const Params& params,
const Evaluation&
krn)
288 for (
int i = 0; i < 20; ++i) {
289 Evaluation f =
krn - twoPhaseSatKrn(params, sw);
290 Evaluation fStar =
krn - twoPhaseSatKrn(params, sw + eps);
291 Evaluation fPrime = (fStar - f) / eps;
293 Evaluation delta = f / fPrime;
298 if (abs(delta) < 1e-10) {
303 throw NumericalProblem(
"Couldn't invert the Brooks-Corey non-wetting phase"
304 " relperm within 20 iterations");
Specification of the material parameters for the Brooks-Corey constitutive relations.
Provides the OPM specific exception classes.
Implementation of the Brooks-Corey capillary pressure <-> saturation relation.
Definition BrooksCorey.hpp:55
static void saturations(Container &values, const Params ¶ms, const FluidState &fs)
Calculate the saturations of the phases starting from their pressure differences.
Definition BrooksCorey.hpp:112
static const bool isPressureDependent
Definition BrooksCorey.hpp:81
static const bool isSaturationDependent
Definition BrooksCorey.hpp:77
static const int numPhases
Definition BrooksCorey.hpp:62
static void capillaryPressures(Container &values, const Params ¶ms, const FluidState &fs)
The capillary pressure-saturation curves.
Definition BrooksCorey.hpp:99
static void relativePermeabilities(Container &values, const Params ¶ms, const FluidState &fs)
The relative permeability-saturation curves.
Definition BrooksCorey.hpp:131
static const bool implementsTwoPhaseApi
Definition BrooksCorey.hpp:69
static Evaluation Sn(const Params ¶ms, const FluidState &fs)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition BrooksCorey.hpp:215
static const bool isTemperatureDependent
Definition BrooksCorey.hpp:85
static Evaluation Sw(const Params ¶ms, const FluidState &fs)
The saturation-capillary pressure curve according to Brooks & Corey.
Definition BrooksCorey.hpp:194
static Evaluation krw(const Params ¶ms, const FluidState &fs)
The relative permeability for the wetting phase of the medium implied by the Brooks-Corey parameteriz...
Definition BrooksCorey.hpp:231
static Evaluation krn(const Params ¶ms, const FluidState &fs)
The relative permeability for the non-wetting phase of the medium as implied by the Brooks-Corey para...
Definition BrooksCorey.hpp:263
static const bool isCompositionDependent
Definition BrooksCorey.hpp:89
static Evaluation pcnw(const Params ¶ms, const FluidState &fs)
The capillary pressure-saturation curve according to Brooks and Corey.
Definition BrooksCorey.hpp:154
static const bool implementsTwoPhaseSatApi
Definition BrooksCorey.hpp:73
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30