diff --git a/opm/output/eclipse/VectorItems/well.hpp b/opm/output/eclipse/VectorItems/well.hpp
index 110f0d057..7b0887068 100644
--- a/opm/output/eclipse/VectorItems/well.hpp
+++ b/opm/output/eclipse/VectorItems/well.hpp
@@ -213,6 +213,7 @@ namespace Opm { namespace RestartIO { namespace Helpers { namespace VectorItems
GasFVF = 34, // Well's producing gas formation volume factor.
+ item36 = 35, // Unknown
item37 = 36, // Unknown
item38 = 37, // Unknown
diff --git a/src/opm/input/eclipse/Schedule/UDQ/UDQToken.cpp b/src/opm/input/eclipse/Schedule/UDQ/UDQToken.cpp
index 618553a9f..e9979c99c 100644
--- a/src/opm/input/eclipse/Schedule/UDQ/UDQToken.cpp
+++ b/src/opm/input/eclipse/Schedule/UDQ/UDQToken.cpp
@@ -16,20 +16,25 @@
You should have received a copy of the GNU General Public License
along with OPM. If not, see .
*/
-#include
-#include
#include
-namespace Opm {
+#include
+#include
+#include
+#include
+
+#include
namespace {
-std::string format_double(double d) {
- return fmt::format("{:g}", d);
-}
-}
-
+ std::string format_double(const double d)
+ {
+ // Use uppercase exponents for restart file compatibility.
+ return fmt::format("{:G}", d);
+ }
+} // namespace anonymous
+namespace Opm {
UDQToken::UDQToken(const std::string& string_token, UDQTokenType token_type_) :
token_type(token_type_)
@@ -74,5 +79,4 @@ std::string UDQToken::str() const {
return format_double(std::get(this->m_value));
}
-
-}
+} // namespace Opm
diff --git a/src/opm/output/eclipse/AggregateWellData.cpp b/src/opm/output/eclipse/AggregateWellData.cpp
index 5dc25cf51..d1a8de6ee 100644
--- a/src/opm/output/eclipse/AggregateWellData.cpp
+++ b/src/opm/output/eclipse/AggregateWellData.cpp
@@ -45,10 +45,9 @@
#include
#include
-#include
-
#include
#include
+#include
#include
#include
#include
@@ -58,6 +57,8 @@
#include
#include
+#include
+
namespace VI = Opm::RestartIO::Helpers::VectorItems;
// #####################################################################
@@ -806,12 +807,13 @@ namespace {
xWell[Ix::HistWatInjTotal] = get("WWITH");
xWell[Ix::HistGasInjTotal] = get("WGITH");
}
+
template
void assignProducer(const std::string& well,
const ::Opm::SummaryState& smry,
XWellArray& xWell)
{
- using Ix = ::Opm::RestartIO::Helpers::VectorItems::XWell::index;
+ using Ix = VI::XWell::index;
auto get = [&smry, &well](const std::string& vector)
{
@@ -832,6 +834,7 @@ namespace {
xWell[Ix::GORatio] = get("WGOR");
// Not fully characterised.
+ xWell[Ix::item36] = xWell[Ix::OilPrRate];
xWell[Ix::item37] = xWell[Ix::WatPrRate];
xWell[Ix::item38] = xWell[Ix::GasPrRate];
@@ -848,6 +851,16 @@ namespace {
{
using Ix = ::Opm::RestartIO::Helpers::VectorItems::XWell::index;
+ // Injection rates reported as negative.
+ xWell[Ix::OilPrRate] = -get("WOIR");
+ xWell[Ix::WatPrRate] = -get("WWIR");
+ xWell[Ix::GasPrRate] = -get("WGIR");
+
+ // Not fully characterised.
+ xWell[Ix::item36] = xWell[Ix::OilPrRate];
+ xWell[Ix::item37] = xWell[Ix::WatPrRate];
+ xWell[Ix::item38] = xWell[Ix::GasPrRate];
+
xWell[Ix::TubHeadPr] = get("WTHP");
xWell[Ix::FlowBHP] = get("WBHP");
}
@@ -861,20 +874,14 @@ namespace {
auto get = [&smry, &well](const std::string& vector)
{
- return smry.get_well_var(well, vector, 0);
+ return smry.get_well_var(well, vector, 0.0);
};
assignCommonInjector(get, xWell);
- // Injection rates reported as negative.
- xWell[Ix::WatPrRate] = -get("WWIR");
xWell[Ix::LiqPrRate] = xWell[Ix::WatPrRate];
- // Not fully characterised.
- xWell[Ix::item37] = xWell[Ix::WatPrRate];
-
xWell[Ix::PrimGuideRate] = xWell[Ix::PrimGuideRate_2] = -get("WWIGR");
-
xWell[Ix::WatVoidPrRate] = -get("WWVIR");
}
@@ -887,13 +894,12 @@ namespace {
auto get = [&smry, &well](const std::string& vector)
{
- return smry.get_well_var(well, vector, 0);
+ return smry.get_well_var(well, vector, 0.0);
};
assignCommonInjector(get, xWell);
// Injection rates reported as negative production rates.
- xWell[Ix::GasPrRate] = -get("WGIR");
xWell[Ix::VoidPrRate] = -get("WGVIR");
xWell[Ix::GasFVF] = (std::abs(xWell[Ix::GasPrRate]) > 0.0)
@@ -904,11 +910,7 @@ namespace {
xWell[Ix::GasFVF] = 0.0;
}
- // Not fully characterised.
- xWell[Ix::item38] = xWell[Ix::GasPrRate];
-
xWell[Ix::PrimGuideRate] = xWell[Ix::PrimGuideRate_2] = -get("WGIGR");
-
xWell[Ix::GasVoidPrRate] = xWell[Ix::VoidPrRate];
}
@@ -917,16 +919,17 @@ namespace {
const ::Opm::SummaryState& smry,
XWellArray& xWell)
{
- using Ix = ::Opm::RestartIO::Helpers::VectorItems::XWell::index;
+ using Ix = VI::XWell::index;
auto get = [&smry, &well](const std::string& vector)
{
- return smry.get_well_var(well, vector, 0);
+ return smry.get_well_var(well, vector, 0.0);
};
- xWell[Ix::TubHeadPr] = get("WTHP");
- xWell[Ix::FlowBHP] = get("WBHP");
+ assignCommonInjector(get, xWell);
+ // Injection rates reported as negative production rates.
+ xWell[Ix::VoidPrRate] = -get("WOVIR");
xWell[Ix::PrimGuideRate] = xWell[Ix::PrimGuideRate_2] = -get("WOIGR");
}
@@ -1091,13 +1094,13 @@ AggregateWellData(const std::vector& inteHead)
void
Opm::RestartIO::Helpers::AggregateWellData::
-captureDeclaredWellData(const Schedule& sched,
- const TracerConfig& tracers,
- const std::size_t sim_step,
+captureDeclaredWellData(const Schedule& sched,
+ const TracerConfig& tracers,
+ const std::size_t sim_step,
const ::Opm::Action::State& action_state,
- const Opm::WellTestState& wtest_state,
+ const Opm::WellTestState& wtest_state,
const ::Opm::SummaryState& smry,
- const std::vector& inteHead)
+ const std::vector& inteHead)
{
const auto& wells = sched.wellNames(sim_step);
const auto& step_glo = sched.glo(sim_step);
@@ -1106,7 +1109,7 @@ captureDeclaredWellData(const Schedule& sched,
{
//const auto grpNames = groupNames(sched.getGroups());
const auto groupMapNameIndex = IWell::currentGroupMapNameIndex(sched, sim_step, inteHead);
- auto msWellID = std::size_t{0};
+ auto msWellID = std::size_t{0};
wellLoop(wells, sched, sim_step, [&groupMapNameIndex, &msWellID, &step_glo, &wtest_state, &smry, &sched, &sim_step, this]
(const Well& well, const std::size_t wellID) -> void
@@ -1155,7 +1158,7 @@ Opm::RestartIO::Helpers::AggregateWellData::
captureDynamicWellData(const Opm::Schedule& sched,
const TracerConfig& tracers,
const std::size_t sim_step,
- const Opm::data::Wells& xw,
+ const Opm::data::Wells& xw,
const ::Opm::SummaryState& smry)
{
const auto& wells = sched.wellNames(sim_step);
diff --git a/tests/parser/UDQTests.cpp b/tests/parser/UDQTests.cpp
index 390726330..a81bdbeec 100644
--- a/tests/parser/UDQTests.cpp
+++ b/tests/parser/UDQTests.cpp
@@ -17,53 +17,65 @@ Copyright 2018 Statoil ASA.
#include
-#include
-#include
-
#include
-#include
-#include
-#include
-#include
-#include
-#include
-#include
-#include
-#include
+#include
+
#include
-#include
-#include
-#include
-#include
-#include
+#include
+
+#include
+
+#include
+#include
+#include
#include
+#include
+#include
+#include
#include
#include
-#include
+#include
#include
-#include
-#include
#include
-#include
+#include
+
+#include
+
+#include
+#include
+#include
+
+#include
+#include
+
+#include
+#include
+#include
+#include
+#include
+#include
+#include
using namespace Opm;
-Schedule make_schedule(const std::string& input) {
- Parser parser;
- auto python = std::make_shared();
+namespace {
+ Schedule make_schedule(const std::string& input) {
+ Parser parser;
+ auto python = std::make_shared();
- auto deck = parser.parseString(input);
- if (deck.hasKeyword("DIMENS")) {
- EclipseState es(deck);
- return Schedule(deck, es, python);
- } else {
- EclipseGrid grid(10,10,10);
- TableManager table ( deck );
- FieldPropsManager fp( deck, Phases{true, true, true}, grid, table);
- Runspec runspec (deck);
- return Schedule(deck, grid , fp, runspec, python);
+ auto deck = parser.parseString(input);
+ if (deck.hasKeyword("DIMENS")) {
+ EclipseState es(deck);
+ return Schedule(deck, es, python);
+ } else {
+ EclipseGrid grid(10,10,10);
+ TableManager table ( deck );
+ FieldPropsManager fp( deck, Phases{true, true, true}, grid, table);
+ Runspec runspec (deck);
+ return Schedule(deck, grid , fp, runspec, python);
+ }
}
-}
+} // namespace anonymous
BOOST_AUTO_TEST_CASE(TYPE_COERCION) {
BOOST_CHECK( UDQVarType::SCALAR == UDQ::coerce(UDQVarType::SCALAR, UDQVarType::SCALAR) );
@@ -1278,6 +1290,25 @@ BOOST_AUTO_TEST_CASE(UDQPARSE_TEST1) {
BOOST_CHECK_EQUAL( def2.input_string() , "2 * (1 + WBHP)");
}
+BOOST_AUTO_TEST_CASE(INPUT_STRING_SCIENTIFIC_NOTATION) {
+ const auto schedule = make_schedule(R"(
+SCHEDULE
+UDQ
+DEFINE FU_THREE (3000000 + FU_ONE*1500000 + 1000000*FU_TWO)/365 /
+/
+)");
+
+ const auto& udq = schedule.getUDQConfig(0);
+ const auto def = udq.definitions();
+
+ BOOST_CHECK_EQUAL(def.size(), 1ULL);
+
+ const auto expect_input_string = std::string {
+ "(3E+06 + FU_ONE * 1.5E+06 + 1E+06 * FU_TWO) / 365"
+ };
+
+ BOOST_CHECK_EQUAL(def[0].input_string(), expect_input_string);
+}
BOOST_AUTO_TEST_CASE(UDQ_PARSE_ERROR) {
UDQParams udqp;
@@ -1610,11 +1641,13 @@ BOOST_AUTO_TEST_CASE(IntegrationTest) {
}
}
-Schedule make_udq_schedule(const std::string& schedule_string) {
+namespace {
+ Schedule make_udq_schedule(const std::string& schedule_string) {
#include "data/integration_tests/udq2.data"
- deck_string += schedule_string;
- return make_schedule(deck_string);
-}
+ deck_string += schedule_string;
+ return make_schedule(deck_string);
+ }
+} // Namespace anonymous
BOOST_AUTO_TEST_CASE(IntegrationTest2) {
const std::string udq_string = R"(