25 #ifndef OPM_ECL_STONE2_MATERIAL_HPP 
   26 #define OPM_ECL_STONE2_MATERIAL_HPP 
   33 #include <opm/common/Exceptions.hpp> 
   34 #include <opm/common/ErrorMacros.hpp> 
   53 template <
class TraitsT,
 
   54           class GasOilMaterialLawT,
 
   55           class OilWaterMaterialLawT,
 
   56           class ParamsT = EclStone2MaterialParams<TraitsT,
 
   57                                                    typename GasOilMaterialLawT::Params,
 
   58                                                    typename OilWaterMaterialLawT::Params> >
 
   66     static_assert(TraitsT::numPhases == 3,
 
   67                   "The number of phases considered by this capillary pressure " 
   68                   "law is always three!");
 
   69     static_assert(GasOilMaterialLaw::numPhases == 2,
 
   70                   "The number of phases considered by the gas-oil capillary " 
   71                   "pressure law must be two!");
 
   72     static_assert(OilWaterMaterialLaw::numPhases == 2,
 
   73                   "The number of phases considered by the oil-water capillary " 
   74                   "pressure law must be two!");
 
   75     static_assert(std::is_same<
typename GasOilMaterialLaw::Scalar,
 
   76                                typename OilWaterMaterialLaw::Scalar>::value,
 
   77                   "The two two-phase capillary pressure laws must use the same " 
   78                   "type of floating point values.");
 
   80     static_assert(GasOilMaterialLaw::implementsTwoPhaseSatApi,
 
   81                   "The gas-oil material law must implement the two-phase saturation " 
   82                   "only API to for the default Ecl capillary pressure law!");
 
   83     static_assert(OilWaterMaterialLaw::implementsTwoPhaseSatApi,
 
   84                   "The oil-water material law must implement the two-phase saturation " 
   85                   "only API to for the default Ecl capillary pressure law!");
 
   87     typedef TraitsT Traits;
 
   88     typedef ParamsT Params;
 
   89     typedef typename Traits::Scalar Scalar;
 
   91     static const int numPhases = 3;
 
   92     static const int waterPhaseIdx = Traits::wettingPhaseIdx;
 
   93     static const int oilPhaseIdx = Traits::nonWettingPhaseIdx;
 
   94     static const int gasPhaseIdx = Traits::gasPhaseIdx;
 
   98     static const bool implementsTwoPhaseApi = 
false;
 
  102     static const bool implementsTwoPhaseSatApi = 
false;
 
  106     static const bool isSaturationDependent = 
true;
 
  110     static const bool isPressureDependent = 
false;
 
  114     static const bool isTemperatureDependent = 
false;
 
  118     static const bool isCompositionDependent = 
false;
 
  134     template <
class ContainerT, 
class Flu
idState>
 
  135     static void capillaryPressures(ContainerT &values,
 
  136                                    const Params ¶ms,
 
  137                                    const FluidState &state)
 
  139         typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
 
  140         values[gasPhaseIdx] = pcgn<FluidState, Evaluation>(params, state);
 
  141         values[oilPhaseIdx] = 0;
 
  142         values[waterPhaseIdx] = - pcnw<FluidState, Evaluation>(params, state);
 
  157     template <
class Flu
idState, 
class Evaluation = 
typename Flu
idState::Scalar>
 
  158     static Evaluation pcgn(
const Params ¶ms,
 
  159                            const FluidState &fs)
 
  163         const auto& Sw = 1.0 - FsToolbox::template toLhs<Evaluation>(fs.saturation(gasPhaseIdx));
 
  164         return GasOilMaterialLaw::twoPhaseSatPcnw(params.gasOilParams(), Sw);
 
  176     template <
class Flu
idState, 
class Evaluation = 
typename Flu
idState::Scalar>
 
  177     static Evaluation pcnw(
const Params ¶ms,
 
  178                            const FluidState &fs)
 
  180         typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
 
  182         const auto& Sw = FsToolbox::template toLhs<Evaluation>(fs.saturation(waterPhaseIdx));
 
  184         const auto& result = OilWaterMaterialLaw::twoPhaseSatPcnw(params.oilWaterParams(), Sw);
 
  192     template <
class ContainerT, 
class Flu
idState>
 
  193     static void saturations(ContainerT& ,
 
  197         OPM_THROW(std::logic_error, 
"Not implemented: saturations()");
 
  203     template <
class Flu
idState, 
class Evaluation = 
typename Flu
idState::Scalar>
 
  204     static Evaluation Sg(
const Params& ,
 
  207         OPM_THROW(std::logic_error, 
"Not implemented: Sg()");
 
  213     template <
class Flu
idState, 
class Evaluation = 
typename Flu
idState::Scalar>
 
  214     static Evaluation Sn(
const Params& ,
 
  217         OPM_THROW(std::logic_error, 
"Not implemented: Sn()");
 
  223     template <
class Flu
idState, 
class Evaluation = 
typename Flu
idState::Scalar>
 
  224     static Evaluation Sw(
const Params& ,
 
  227         OPM_THROW(std::logic_error, 
"Not implemented: Sw()");
 
  245     template <
class ContainerT, 
class Flu
idState>
 
  246     static void relativePermeabilities(ContainerT &values,
 
  247                                        const Params ¶ms,
 
  248                                        const FluidState &fluidState)
 
  250         typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
 
  252         values[waterPhaseIdx] = krw<FluidState, Evaluation>(params, fluidState);
 
  253         values[oilPhaseIdx] = krn<FluidState, Evaluation>(params, fluidState);
 
  254         values[gasPhaseIdx] = krg<FluidState, Evaluation>(params, fluidState);
 
  260     template <
class Flu
idState, 
class Evaluation = 
typename Flu
idState::Scalar>
 
  261     static Evaluation krg(
const Params ¶ms,
 
  262                           const FluidState &fluidState)
 
  264         typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
 
  266         const Evaluation& Sw = 1 - FsToolbox::template toLhs<Evaluation>(fluidState.saturation(gasPhaseIdx));
 
  267         return GasOilMaterialLaw::twoPhaseSatKrn(params.gasOilParams(), Sw);
 
  273     template <
class Flu
idState, 
class Evaluation = 
typename Flu
idState::Scalar>
 
  274     static Evaluation krw(
const Params ¶ms,
 
  275                           const FluidState &fluidState)
 
  277         typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
 
  279         const Evaluation& Sw = FsToolbox::template toLhs<Evaluation>(fluidState.saturation(waterPhaseIdx));
 
  280         return OilWaterMaterialLaw::twoPhaseSatKrw(params.oilWaterParams(), Sw);
 
  286     template <
class Flu
idState, 
class Evaluation = 
typename Flu
idState::Scalar>
 
  287     static Evaluation krn(
const Params ¶ms,
 
  288                           const FluidState &fluidState)
 
  290         typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
 
  292         Scalar Swco = params.Swl();
 
  293         const Evaluation& Sw = FsToolbox::template toLhs<Evaluation>(fluidState.saturation(waterPhaseIdx));
 
  294         const Evaluation& Sg = FsToolbox::template toLhs<Evaluation>(fluidState.saturation(gasPhaseIdx));
 
  296         Scalar krocw = OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Swco);
 
  297         Evaluation krow = OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Sw);
 
  298         Evaluation krw = OilWaterMaterialLaw::twoPhaseSatKrw(params.oilWaterParams(), Sw);
 
  299         Evaluation krg = GasOilMaterialLaw::twoPhaseSatKrn(params.gasOilParams(), 1 - Sg);
 
  300         Evaluation krog = GasOilMaterialLaw::twoPhaseSatKrw(params.gasOilParams(), 1 - Sg);
 
  302         return krocw*((krow/krocw + krw)*(krog/krocw + krg) - krw - krg);
 
  312     template <
class Flu
idState>
 
  313     static void updateHysteresis(Params ¶ms, 
const FluidState &fluidState)
 
  315         typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
 
  317         Scalar Sw = FsToolbox::value(fluidState.saturation(waterPhaseIdx));
 
  318         Scalar Sg = FsToolbox::value(fluidState.saturation(gasPhaseIdx));
 
  320         params.oilWaterParams().update(Sw, Sw, Sw);
 
  321         params.gasOilParams().update(1 - Sg, 1 - Sg, 1 - Sg);
 
bool CheckDefined(const T &value OPM_UNUSED)
Make valgrind complain if any of the memory occupied by an object is undefined. 
Definition: Valgrind.hpp:74
 
Definition: Air_Mesitylene.hpp:31
 
Some templates to wrap the valgrind client request macros. 
 
OilWaterMaterialLawT OilWaterMaterialLaw
Definition: EclStone2Material.hpp:63
 
Implements the second phase capillary pressure/relperm law suggested by Stone as used by the ECLipse ...
Definition: EclStone2Material.hpp:59
 
GasOilMaterialLawT GasOilMaterialLaw
Definition: EclStone2Material.hpp:62
 
Default implementation for the parameters required by the three-phase capillary pressure/relperm Ston...