/* Copyright 2015 SINTEF ICT, Applied Mathematics. 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 "config.h" #include #include #include #include #include #include namespace Opm { template Scalar VFPProdProperties:: thp(const int table_id, const Scalar aqua, const Scalar liquid, const Scalar vapour, const Scalar bhp_arg, const Scalar alq, const Scalar explicit_wfr, const Scalar explicit_gfr, const bool use_expvfp) const { const VFPProdTable& table = detail::getTable(m_tables, table_id); Scalar flo = detail::getFlo(table, aqua, liquid, vapour); Scalar wfr = detail::getWFR(table, aqua, liquid, vapour); Scalar gfr = detail::getGFR(table, aqua, liquid, vapour); if (use_expvfp || -flo < table.getFloAxis().front()) { wfr = explicit_wfr; gfr = explicit_gfr; } const std::vector thp_array = table.getTHPAxis(); int nthp = thp_array.size(); /** * Find the function bhp_array(thp) by creating a 1D view of the data * by interpolating for every value of thp. This might be somewhat * expensive, but let us assome that nthp is small. */ auto flo_i = VFPHelpers::findInterpData(-flo, table.getFloAxis()); auto wfr_i = VFPHelpers::findInterpData( wfr, table.getWFRAxis()); auto gfr_i = VFPHelpers::findInterpData( gfr, table.getGFRAxis()); auto alq_i = VFPHelpers::findInterpData( alq, table.getALQAxis()); std::vector bhp_array(nthp); for (int i = 0; i < nthp; ++i) { auto thp_i = VFPHelpers::findInterpData(thp_array[i], thp_array); bhp_array[i] = VFPHelpers::interpolate(table, flo_i, thp_i, wfr_i, gfr_i, alq_i).value; } return VFPHelpers::findTHP(bhp_array, thp_array, bhp_arg, /*find_largest*/ true); } template Scalar VFPProdProperties:: bhp(const int table_id, const Scalar aqua, const Scalar liquid, const Scalar vapour, const Scalar thp_arg, const Scalar alq, const Scalar explicit_wfr, const Scalar explicit_gfr, const bool use_expvfp) const { const VFPProdTable& table = detail::getTable(m_tables, table_id); detail::VFPEvaluation retval = VFPHelpers::bhp(table, aqua, liquid, vapour, thp_arg, alq, explicit_wfr, explicit_gfr, use_expvfp); return retval.value; } template const VFPProdTable& VFPProdProperties::getTable(const int table_id) const { return detail::getTable(m_tables, table_id); } template bool VFPProdProperties::hasTable(const int table_id) const { return detail::hasTable(m_tables, table_id); } template std::vector VFPProdProperties:: bhpwithflo(const std::vector& flos, const int table_id, const Scalar wfr, const Scalar gfr, const Scalar thp, const Scalar alq, const Scalar dp) const { // Get the table const VFPProdTable& table = detail::getTable(m_tables, table_id); const auto thp_i = VFPHelpers::findInterpData( thp, table.getTHPAxis()); // assume constant const auto wfr_i = VFPHelpers::findInterpData( wfr, table.getWFRAxis()); const auto gfr_i = VFPHelpers::findInterpData( gfr, table.getGFRAxis()); const auto alq_i = VFPHelpers::findInterpData( alq, table.getALQAxis()); //assume constant std::vector bhps(flos.size(), 0.); for (std::size_t i = 0; i < flos.size(); ++i) { // Value of FLO is negative in OPM for producers, but positive in VFP table const auto flo_i = VFPHelpers::findInterpData(-flos[i], table.getFloAxis()); const detail::VFPEvaluation bhp_val = VFPHelpers::interpolate(table, flo_i, thp_i, wfr_i, gfr_i, alq_i); // TODO: this kind of breaks the conventions for the functions here by putting dp within the function bhps[i] = bhp_val.value - dp; } return bhps; } template Scalar VFPProdProperties:: minimumBHP(const int table_id, const Scalar thp, const Scalar wfr, const Scalar gfr, const Scalar alq) const { // Get the table const VFPProdTable& table = detail::getTable(m_tables, table_id); const auto retval = VFPHelpers::getMinimumBHPCoordinate(table, thp, wfr, gfr, alq); // returned pair is (flo, bhp) return retval.second; } template void VFPProdProperties::addTable(const VFPProdTable& new_table) { this->m_tables.emplace( new_table.getTableNum(), new_table ); } template template EvalWell VFPProdProperties:: bhp(const int table_id, const EvalWell& aqua, const EvalWell& liquid, const EvalWell& vapour, const Scalar thp, const Scalar alq, const Scalar explicit_wfr, const Scalar explicit_gfr, const bool use_expvfp) const { //Get the table const VFPProdTable& table = detail::getTable(m_tables, table_id); EvalWell bhp = 0.0 * aqua; //Find interpolation variables EvalWell flo = detail::getFlo(table, aqua, liquid, vapour); EvalWell wfr = detail::getWFR(table, aqua, liquid, vapour); EvalWell gfr = detail::getGFR(table, aqua, liquid, vapour); if (use_expvfp || -flo.value() < table.getFloAxis().front()) { wfr = explicit_wfr; gfr = explicit_gfr; } //First, find the values to interpolate between //Value of FLO is negative in OPM for producers, but positive in VFP table auto flo_i = VFPHelpers::findInterpData(-flo.value(), table.getFloAxis()); auto thp_i = VFPHelpers::findInterpData( thp, table.getTHPAxis()); // assume constant auto wfr_i = VFPHelpers::findInterpData( wfr.value(), table.getWFRAxis()); auto gfr_i = VFPHelpers::findInterpData( gfr.value(), table.getGFRAxis()); auto alq_i = VFPHelpers::findInterpData( alq, table.getALQAxis()); //assume constant detail::VFPEvaluation bhp_val = VFPHelpers::interpolate(table, flo_i, thp_i, wfr_i, gfr_i, alq_i); bhp = (bhp_val.dwfr * wfr) + (bhp_val.dgfr * gfr) - (std::max(Scalar{0.0}, bhp_val.dflo) * flo); bhp.setValue(bhp_val.value); return bhp; } #define INSTANTIATE(T,...) \ template __VA_ARGS__ \ VFPProdProperties::bhp(const int, \ const __VA_ARGS__&, \ const __VA_ARGS__&, \ const __VA_ARGS__&, \ const T , \ const T , \ const T , \ const T , \ const bool) const; #define INSTANTIATE_TYPE(T) \ template class VFPProdProperties; \ INSTANTIATE(T,DenseAd::Evaluation) \ INSTANTIATE(T,DenseAd::Evaluation) \ INSTANTIATE(T,DenseAd::Evaluation) \ INSTANTIATE(T,DenseAd::Evaluation) \ INSTANTIATE(T,DenseAd::Evaluation) \ INSTANTIATE(T,DenseAd::Evaluation) \ INSTANTIATE(T,DenseAd::Evaluation) \ INSTANTIATE(T,DenseAd::Evaluation) \ INSTANTIATE(T,DenseAd::Evaluation) \ INSTANTIATE(T,DenseAd::Evaluation) \ INSTANTIATE(T,DenseAd::Evaluation) \ INSTANTIATE(T,DenseAd::Evaluation) \ INSTANTIATE(T,DenseAd::Evaluation) \ INSTANTIATE(T,DenseAd::Evaluation) \ INSTANTIATE(T,DenseAd::Evaluation) \ INSTANTIATE(T,DenseAd::Evaluation) INSTANTIATE_TYPE(double) #if FLOW_INSTANTIATE_FLOAT INSTANTIATE_TYPE(float) #endif }