/* Copyright (C) 2009-2013 by Andreas Lauser This file is part of the Open Porous Media project (OPM). OPM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 2 of the License, or (at your option) any later version. OPM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OPM. If not, see . */ /*! * \file * \copydoc Opm::PiecewiseLinearTwoPhaseMaterial */ #ifndef OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP #define OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_HPP #include "PiecewiseLinearTwoPhaseMaterialParams.hpp" #include #include #include #include #include #include namespace Opm { /*! * \ingroup FluidMatrixInteractions * * \brief Implementation of a tabulated, piecewise linear capillary * pressure law. * * It would be equally possible to use cubic splines, but since the * ECLIPSE reservoir simulator uses linear interpolation for capillary * pressure and relperm curves, we do the same. */ template > class PiecewiseLinearTwoPhaseMaterial : public TraitsT { typedef typename ParamsT::ValueVector ValueVector; public: //! The traits class for this material law typedef TraitsT Traits; //! The type of the parameter objects for this law typedef ParamsT Params; //! The type of the scalar values for this law typedef typename Traits::Scalar Scalar; //! The number of fluid phases static const int numPhases = Traits::numPhases; static_assert(numPhases == 2, "The piecewise linear two-phase capillary pressure law only" "applies to the case of two fluid phases"); //! Specify whether this material law implements the two-phase //! convenience API static const bool implementsTwoPhaseApi = true; //! Specify whether this material law implements the two-phase //! convenience API which only depends on the phase saturations static const bool implementsTwoPhaseSatApi = true; //! Specify whether the quantities defined by this material law //! are saturation dependent static const bool isSaturationDependent = true; //! Specify whether the quantities defined by this material law //! are dependent on the absolute pressure static const bool isPressureDependent = false; //! Specify whether the quantities defined by this material law //! are temperature dependent static const bool isTemperatureDependent = false; //! Specify whether the quantities defined by this material law //! are dependent on the phase composition static const bool isCompositionDependent = false; /*! * \brief The capillary pressure-saturation curve. */ template static void capillaryPressures(Container &values, const Params ¶ms, const FluidState &fs) { typedef typename std::remove_reference::type Evaluation; values[Traits::wettingPhaseIdx] = 0.0; // reference phase values[Traits::nonWettingPhaseIdx] = pcnw(params, fs); } /*! * \brief The saturations of the fluid phases starting from their * pressure differences. */ template static void saturations(Container &values, const Params ¶ms, const FluidState &fs) { OPM_THROW(std::logic_error, "Not implemented: saturations()"); } /*! * \brief The relative permeabilities */ template static void relativePermeabilities(Container &values, const Params ¶ms, const FluidState &fs) { typedef typename std::remove_reference::type Evaluation; values[Traits::wettingPhaseIdx] = krw(params, fs); values[Traits::nonWettingPhaseIdx] = krn(params, fs); } /*! * \brief The capillary pressure-saturation curve */ template static Evaluation pcnw(const Params ¶ms, const FluidState &fs) { typedef MathToolbox FsToolbox; const auto& Sw = FsToolbox::template toLhs(fs.saturation(Traits::wettingPhaseIdx)); return twoPhaseSatPcnw(params, Sw); } /*! * \brief The saturation-capillary pressure curve */ template static Evaluation twoPhaseSatPcnw(const Params ¶ms, const Evaluation& Sw) { return eval_(params.SwSamples(), params.pcnwSamples(), Sw); } /*! * \brief The saturation-capillary pressure curve */ template static Evaluation Sw(const Params ¶ms, const FluidState &fs) { OPM_THROW(std::logic_error, "Not implemented: Sw()"); } template static Evaluation twoPhaseSatSw(const Params ¶ms, const Evaluation& pC) { OPM_THROW(std::logic_error, "Not implemented: twoPhaseSatSw()"); } /*! * \brief Calculate the non-wetting phase saturations depending on * the phase pressures. */ template static Evaluation Sn(const Params ¶ms, const FluidState &fs) { return 1 - Sw(params, fs); } template static Evaluation twoPhaseSatSn(const Params ¶ms, const Evaluation& pC) { return 1 - twoPhaseSatSw(params, pC); } /*! * \brief The relative permeability for the wetting phase of the * porous medium */ template static Evaluation krw(const Params ¶ms, const FluidState &fs) { typedef MathToolbox FsToolbox; const auto& Sw = FsToolbox::template toLhs(fs.saturation(Traits::wettingPhaseIdx)); return twoPhaseSatKrw(params, Sw); } template static Evaluation twoPhaseSatKrw(const Params ¶ms, const Evaluation Sw) { typedef MathToolbox Toolbox; const auto& res = eval_(params.SwSamples(), params.krwSamples(), Sw); return Toolbox::max(0.0, Toolbox::min(1.0, res)); } /*! * \brief The relative permeability for the non-wetting phase * of the porous medium */ template static Evaluation krn(const Params ¶ms, const FluidState &fs) { typedef MathToolbox FsToolbox; const auto& Sw = FsToolbox::template toLhs(fs.saturation(Traits::wettingPhaseIdx)); return twoPhaseSatKrn(params, Sw); } template static Evaluation twoPhaseSatKrn(const Params ¶ms, const Evaluation& Sw) { typedef MathToolbox Toolbox; return Toolbox::max(0.0, Toolbox::min(1.0, eval_(params.SwSamples(), params.krnSamples(), Sw))); } private: template static Evaluation eval_(const ValueVector &xValues, const ValueVector &yValues, const Evaluation& x) { typedef MathToolbox Toolbox; if (Toolbox::value(x) < xValues.front()) return Toolbox::createConstant(yValues.front()); if (Toolbox::value(x) > xValues.back()) return Toolbox::createConstant(yValues.back()); int segIdx = findSegmentIndex_(xValues, Toolbox::value(x)); Scalar x0 = xValues[segIdx]; Scalar x1 = xValues[segIdx + 1]; Scalar y0 = yValues[segIdx]; Scalar y1 = yValues[segIdx + 1]; Scalar m = (y1 - y0)/(x1 - x0); return y0 + (x - x0)*m; } template static Evaluation evalDeriv_(const ValueVector &xValues, const ValueVector &yValues, const Evaluation& x) { typedef MathToolbox Toolbox; if (Toolbox::value(x) < xValues.front()) return Toolbox::createConstant(0.0); if (Toolbox::value(x) > xValues.back()) return Toolbox::createConstant(0.0); int segIdx = findSegmentIndex_(xValues, Toolbox::value(x)); Scalar x0 = xValues[segIdx]; Scalar x1 = xValues[segIdx + 1]; Scalar y0 = yValues[segIdx]; Scalar y1 = yValues[segIdx + 1]; return Toolbox::createConstant((y1 - y0)/(x1 - x0)); } static int findSegmentIndex_(const ValueVector &xValues, Scalar x) { int n = xValues.size() - 1; assert(n >= 1); // we need at least two sampling points! if (xValues[n] < x) return n - 1; else if (xValues[0] > x) return 0; // bisection int lowIdx = 0, highIdx = n; while (lowIdx + 1 < highIdx) { int curIdx = (lowIdx + highIdx)/2; if (xValues[curIdx] < x) lowIdx = curIdx; else highIdx = curIdx; } return lowIdx; } }; } // namespace Opm #endif