/* Copyright 2017 Statoil ASA. 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 3 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 . */ #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include namespace { ::Opm::ECLPVT::ConvertUnits deadOilUnitConverter(const ::Opm::ECLUnits::UnitSystem& usys) { using ToSI = ::Opm::ECLPVT::CreateUnitConverter::ToSI; // [ Po, 1/B, 1/(B*mu), d(1/B)/dPo, d(1/(B*mu))/dPo ] return ::Opm::ECLPVT::ConvertUnits { ToSI::pressure(usys), { ToSI::recipFvf(usys), ToSI::recipFvfVisc(usys), ToSI::recipFvfDerivPress(usys), ToSI::recipFvfViscDerivPress(usys) } }; } ::Opm::ECLPVT::ConvertUnits deadOilUnitConverter(const int usys) { const auto u = ::Opm::ECLUnits::createUnitSystem(usys); return deadOilUnitConverter(*u); } std::pair< ::Opm::ECLPVT::ConvertUnits::Converter, ::Opm::ECLPVT::ConvertUnits> liveOilUnitConverter(const int usys) { using ToSI = ::Opm::ECLPVT::CreateUnitConverter::ToSI; const auto u = ::Opm::ECLUnits::createUnitSystem(usys); // Key = Rs // Table = [ Po, 1/B, 1/(B*mu), d(1/B)/dPo, d(1/(B*mu))/dPo ] // = dead oil table format. return std::make_pair(ToSI::disGas(*u), deadOilUnitConverter(*u)); } } // --------------------------------------------------------------------- // Enable runtime selection of dead or live oil functions. class PVxOBase { public: virtual std::vector formationVolumeFactor(const std::vector& rs, const std::vector& po) const = 0; virtual std::vector viscosity(const std::vector& rs, const std::vector& po) const = 0; virtual std::vector getPvtCurve(const Opm::ECLPVT::RawCurve curve, const Opm::ECLUnits::UnitSystem& usys) const = 0; virtual std::unique_ptr clone() const = 0; }; // ===================================================================== class DeadOil : public PVxOBase { public: using ElemIt = ::Opm::ECLPVT::PVDx::ElemIt; using ConvertUnits = ::Opm::ECLPVT::ConvertUnits; DeadOil(ElemIt xBegin, ElemIt xEnd, const ConvertUnits& convert, std::vector& colIt) : interpolant_(xBegin, xEnd, convert, colIt) {} virtual std::vector formationVolumeFactor(const std::vector& /* rs */, const std::vector& po) const override { return this->interpolant_.formationVolumeFactor(po); } virtual std::vector viscosity(const std::vector& /* rs */, const std::vector& po) const override { return this->interpolant_.viscosity(po); } virtual std::vector getPvtCurve(const Opm::ECLPVT::RawCurve curve, const Opm::ECLUnits::UnitSystem& usys) const override { if (curve == Opm::ECLPVT::RawCurve::SaturatedState) { // Not applicable to dead oil. Return empty. return { Opm::FlowDiagnostics::Graph{} }; } auto pvtcurve = this->interpolant_.getPvtCurve(curve); const auto x_unit = usys.pressure(); const auto y_unit = (curve == ::Opm::ECLPVT::RawCurve::FVF) ? (usys.reservoirVolume() / usys.surfaceVolumeLiquid()) : usys.viscosity(); auto& x = pvtcurve.first; auto& y = pvtcurve.second; for (auto n = x.size(), i = 0*n; i < n; ++i) { x[i] = ::Opm::unit::convert::to(x[i], x_unit); y[i] = ::Opm::unit::convert::to(y[i], y_unit); } return { std::move(pvtcurve) }; } virtual std::unique_ptr clone() const override { return std::unique_ptr(new DeadOil(*this)); } private: ::Opm::ECLPVT::PVDx interpolant_; }; // ===================================================================== class LiveOil : public PVxOBase { public: using Extrap = ::Opm::Interp1D::PiecewisePolynomial:: ExtrapolationPolicy::Linearly; using SubtableInterpolant = ::Opm::Interp1D::PiecewisePolynomial:: Linear; LiveOil(std::vector key, std::vector propInterp) : key_ (key) // Copy , interp_(std::move(key), std::move(propInterp)) {} virtual std::vector formationVolumeFactor(const std::vector& rs, const std::vector& po) const override { // PKey Inner C0 C1 C2 C3 // Rs Po 1/B 1/(B*mu) d(1/B)/dPo d(1/(B*mu))/dPo // : : : : : const auto key = TableInterpolant::PrimaryKey { rs }; const auto x = TableInterpolant::InnerVariate { po }; return this->interp_.formationVolumeFactor(key, x); } virtual std::vector viscosity(const std::vector& rs, const std::vector& po) const override { // PKey Inner C0 C1 C2 C3 // Rs Po 1/B 1/(B*mu) d(1/B)/dPo d(1/(B*mu))/dPo // : : : : : const auto key = TableInterpolant::PrimaryKey { rs }; const auto x = TableInterpolant::InnerVariate { po }; return this->interp_.viscosity(key, x); } virtual std::vector getPvtCurve(const Opm::ECLPVT::RawCurve curve, const Opm::ECLUnits::UnitSystem& usys) const override { if (curve == ::Opm::ECLPVT::RawCurve::SaturatedState) { return this->saturatedStateCurve(usys); } return this->mainPvtCurve(curve, usys); } virtual std::unique_ptr clone() const override { return std::unique_ptr(new LiveOil(*this)); } private: using TableInterpolant = ::Opm::ECLPVT::PVTx; std::vector key_; TableInterpolant interp_; std::vector mainPvtCurve(const Opm::ECLPVT::RawCurve curve, const Opm::ECLUnits::UnitSystem& usys) const; std::vector saturatedStateCurve(const Opm::ECLUnits::UnitSystem& usys) const; }; std::vector LiveOil::mainPvtCurve(const Opm::ECLPVT::RawCurve curve, const Opm::ECLUnits::UnitSystem& usys) const { auto curves = this->interp_.getPvtCurve(curve); const auto x_unit = usys.pressure(); const auto y_unit = (curve == ::Opm::ECLPVT::RawCurve::FVF) ? (usys.reservoirVolume() / usys.surfaceVolumeLiquid()) : usys.viscosity(); for (auto& crv : curves) { auto& x = crv.first; auto& y = crv.second; for (auto n = x.size(), i = 0*n; i < n; ++i) { x[i] = ::Opm::unit::convert::to(x[i], x_unit); y[i] = ::Opm::unit::convert::to(y[i], y_unit); } } return curves; } std::vector LiveOil::saturatedStateCurve(const Opm::ECLUnits::UnitSystem& usys) const { const auto press_unit = usys.pressure(); const auto rs_unit = usys.dissolvedGasOilRat(); auto rs = this->key_; auto p = this->interp_.getSaturatedPoints(); if (rs.size() != p.size()) { throw std::invalid_argument { "Inconsistent table sizes of saturated gas function (RsSat(p))" }; } for (auto n = rs.size(), i = 0*n; i < n; ++i) { p [i] = ::Opm::unit::convert::to(p [i], press_unit); rs[i] = ::Opm::unit::convert::to(rs[i], rs_unit); } auto graph = Opm::FlowDiagnostics::Graph { std::move(p), std::move(rs) }; return { std::move(graph) }; } // ##################################################################### namespace { std::vector> createDeadOil(const ::Opm::ECLPropTableRawData& raw, const int usys) { using PVTInterp = std::unique_ptr; using ElmIt = ::Opm::ECLPropTableRawData::ElementIterator; assert ((raw.numPrimary == 1) && "Can't Create Dead Oil Function From Live Oil Table"); const auto cvrt = deadOilUnitConverter(usys); return ::Opm::MakeInterpolants::fromRawData(raw, [&cvrt](ElmIt xBegin, ElmIt xEnd, std::vector& colIt) { return PVTInterp{ new DeadOil(xBegin, xEnd, cvrt, colIt) }; }); } std::vector extractPrimaryKey(const ::Opm::ECLPropTableRawData& raw, const ::Opm::ECLPropTableRawData::SizeType t, const ::Opm::ECLPVT::ConvertUnits::Converter& cvrtKey) { auto key = std::vector{}; key.reserve(raw.numPrimary); for (auto begin = std::begin(raw.primaryKey) + (t + 0)*raw.numPrimary, end = begin + raw.numPrimary; begin != end; ++begin) { if (std::abs(*begin) < 1.0e20) { key.push_back(cvrtKey(*begin)); } } return key; } std::vector> createLiveOil(const ::Opm::ECLPropTableRawData& raw, const int usys) { auto ret = std::vector>{}; ret.reserve(raw.numTables); using Extrap = LiveOil::Extrap; using StI = LiveOil::SubtableInterpolant; using ElemIt = ::Opm::ECLPropTableRawData::ElementIterator; const auto cvrt = liveOilUnitConverter(usys); auto sti = ::Opm::MakeInterpolants::fromRawData(raw, [&cvrt](ElemIt xBegin, ElemIt xEnd, std::vector& colIt) -> StI { try { return StI(Extrap{}, xBegin, xEnd, colIt, cvrt.second.indep, cvrt.second.column); } catch (const std::invalid_argument&) { // No valid nodes. Return invalid. return StI(Extrap{}); } }); for (auto t = 0*raw.numTables; t < raw.numTables; ++t) { auto key = extractPrimaryKey(raw, t, cvrt.first); const auto begin = (t + 0)*raw.numPrimary; const auto end = begin + key.size(); ret.emplace_back(new LiveOil(std::move(key), { std::make_move_iterator(std::begin(sti) + begin), std::make_move_iterator(std::begin(sti) + end) })); } return ret; } std::vector> createPVTFunction(const ::Opm::ECLPropTableRawData& raw, const int usys) { if (raw.numPrimary == 0) { // Malformed Gas PVT table. throw std::invalid_argument { "Oil PVT Table Without Primary Lookup Key" }; } if (raw.numCols != 5) { throw std::invalid_argument { "PVT Table for Oil Must Have Five Columns" }; } if (raw.primaryKey.size() != (raw.numPrimary * raw.numTables)) { throw std::invalid_argument { "Size Mismatch in RS Nodes of PVT Table for Oil" }; } if (raw.data.size() != (raw.numPrimary * raw.numRows * raw.numCols * raw.numTables)) { throw std::invalid_argument { "Size Mismatch in Condensed Table Data " "of PVT Table for Oil" }; } if (raw.numPrimary == 1) { return createDeadOil(raw, usys); } return createLiveOil(raw, usys); } std::vector> clone(const std::vector>& src) { auto dest = std::vector>{}; dest.reserve(src.size()); for (const auto& p : src) { dest.push_back(p->clone()); } return dest; } } // ##################################################################### // ##################################################################### // ===================================================================== // Class ECLPVT::Oil::Impl // --------------------------------------------------------------------- class Opm::ECLPVT::Oil::Impl { private: using EvalPtr = std::unique_ptr; public: Impl(const ECLPropTableRawData& raw, const int usys, std::vector rhoS); Impl(const Impl& rhs); Impl(Impl&& rhs); Impl& operator=(const Impl& rhs); Impl& operator=(Impl&& rhs); using RegIdx = std::vector::size_type; std::vector formationVolumeFactor(const RegIdx region, const std::vector& rs, const std::vector& po) const; std::vector viscosity(const RegIdx region, const std::vector& rs, const std::vector& po) const; double surfaceMassDensity(const RegIdx region) const { this->validateRegIdx(region); return this->rhoS_[region]; } std::vector getPvtCurve(const RegIdx region, const RawCurve curve, const ECLUnits::UnitSystem& usys) const; private: std::vector eval_; std::vector rhoS_; void validateRegIdx(const RegIdx region) const; }; Opm::ECLPVT::Oil::Impl:: Impl(const ECLPropTableRawData& raw, const int usys, std::vector rhoS) : eval_(createPVTFunction(raw, usys)) , rhoS_(std::move(rhoS)) {} Opm::ECLPVT::Oil::Impl::Impl(const Impl& rhs) : eval_(clone(rhs.eval_)) , rhoS_(rhs.rhoS_) {} Opm::ECLPVT::Oil::Impl::Impl(Impl&& rhs) : eval_(std::move(rhs.eval_)) , rhoS_(std::move(rhs.rhoS_)) {} Opm::ECLPVT::Oil::Impl& Opm::ECLPVT::Oil::Impl::operator=(const Impl& rhs) { this->eval_ = clone(rhs.eval_); this->rhoS_ = rhs.rhoS_; return *this; } Opm::ECLPVT::Oil::Impl& Opm::ECLPVT::Oil::Impl::operator=(Impl&& rhs) { this->eval_ = std::move(rhs.eval_); this->rhoS_ = std::move(rhs.rhoS_); return *this; } std::vector Opm::ECLPVT::Oil::Impl:: formationVolumeFactor(const RegIdx region, const std::vector& rs, const std::vector& po) const { this->validateRegIdx(region); return this->eval_[region]->formationVolumeFactor(rs, po); } std::vector Opm::ECLPVT::Oil::Impl:: viscosity(const RegIdx region, const std::vector& rs, const std::vector& po) const { this->validateRegIdx(region); return this->eval_[region]->viscosity(rs, po); } std::vector Opm::ECLPVT::Oil::Impl:: getPvtCurve(const RegIdx region, const RawCurve curve, const ECLUnits::UnitSystem& usys) const { this->validateRegIdx(region); return this->eval_[region]->getPvtCurve(curve, usys); } void Opm::ECLPVT::Oil::Impl::validateRegIdx(const RegIdx region) const { if (region >= this->eval_.size()) { throw std::invalid_argument { "Region Index " + std::to_string(region) + " Outside Valid Range (0 .. " + std::to_string(this->eval_.size() - 1) + ')' }; } } // ====================================================================== // Class ECLPVT::Oil // ---------------------------------------------------------------------- Opm::ECLPVT::Oil::Oil(const ECLPropTableRawData& raw, const int usys, std::vector rhoS) : pImpl_(new Impl(raw, usys, std::move(rhoS))) {} Opm::ECLPVT::Oil::~Oil() {} Opm::ECLPVT::Oil::Oil(const Oil& rhs) : pImpl_(new Impl(*rhs.pImpl_)) {} Opm::ECLPVT::Oil::Oil(Oil&& rhs) : pImpl_(std::move(rhs.pImpl_)) {} Opm::ECLPVT::Oil& Opm::ECLPVT::Oil::operator=(const Oil& rhs) { this->pImpl_.reset(new Impl(*rhs.pImpl_)); return *this; } Opm::ECLPVT::Oil& Opm::ECLPVT::Oil::operator=(Oil&& rhs) { this->pImpl_ = std::move(rhs.pImpl_); return *this; } std::vector Opm::ECLPVT::Oil:: formationVolumeFactor(const int region, const DissolvedGas& rs, const OilPressure& po) const { return this->pImpl_->formationVolumeFactor(region, rs.data, po.data); } std::vector Opm::ECLPVT::Oil:: viscosity(const int region, const DissolvedGas& rs, const OilPressure& po) const { return this->pImpl_->viscosity(region, rs.data, po.data); } double Opm::ECLPVT::Oil::surfaceMassDensity(const int region) const { return this->pImpl_->surfaceMassDensity(region); } std::vector Opm::ECLPVT::Oil:: getPvtCurve(const RawCurve curve, const int region, const ECLUnits::UnitSystem& usys) const { return this->pImpl_->getPvtCurve(region, curve, usys); } // ===================================================================== std::unique_ptr Opm::ECLPVT::CreateOilPVTInterpolant:: fromECLOutput(const ECLInitFileData& init) { using OPtr = std::unique_ptr; const auto& ih = init.keywordData(INTEHEAD_KW); const auto iphs = static_cast(ih[INTEHEAD_PHASE_INDEX]); if ((iphs & (1u << 0)) == 0) { // Oil is not an active phase (unexpected). // Return sentinel (null) pointer. return OPtr{}; } auto raw = ::Opm::ECLPropTableRawData{}; const auto& tabdims = init.keywordData("TABDIMS"); const auto& tab = init.keywordData("TAB"); raw.numPrimary = tabdims[ TABDIMS_NRPVTO_ITEM ]; // #Rs nodes/full table raw.numRows = tabdims[ TABDIMS_NPPVTO_ITEM ]; // #Po nodes/sub-table raw.numCols = 5; // [ Po, 1/B, 1/(B*mu), d(1/B)/dPo, d(1/(B*mu))/dPo ] raw.numTables = tabdims[ TABDIMS_NTPVTO_ITEM ]; // # PVTO tables // Extract Primary Key (Rs) { const auto nTabElem = raw.numPrimary * raw.numTables; // Subtract one to account for 1-based indices. const auto start = tabdims[ TABDIMS_JBPVTO_OFFSET_ITEM ] - 1; raw.primaryKey.assign(&tab[start], &tab[start] + nTabElem); } // Extract Full Table { const auto nTabElem = raw.numPrimary * raw.numRows * raw.numCols * raw.numTables; // Subtract one to account for 1-based indices. const auto start = tabdims[ TABDIMS_IBPVTO_OFFSET_ITEM ] - 1; raw.data.assign(&tab[start], &tab[start] + nTabElem); } auto rhoS = surfaceMassDensity(init, ECLPhaseIndex::Liquid); return OPtr{ new Oil(raw, ih[ INTEHEAD_UNIT_INDEX ], std::move(rhoS)) }; }