Some improvements to runtime type conversion when evaluating UDQ

This commit is contained in:
Joakim Hove
2019-12-13 15:59:04 +01:00
parent fdd159cc75
commit af1b140723
7 changed files with 164 additions and 40 deletions

View File

@@ -24,8 +24,36 @@
namespace Opm {
/*
The UDQ variables can be of of many different types. In addition they can be
either scalars or vector sets. The arch example of a vector set is well
variables - in the expressions:
UDQ
DEFINE WUBHP WBHP * 1.15 /
DEFINE WUORAT 1000 /
/
we define two UDQ values 'WUBHP' and 'WUORAT'. Both of these UDQ values will
apply to all wells; the WUBHP vector will correspond to the normal BHP scaled
up 15%, the WUORAT has the scalar value 1000 for all the wells. The well sets
can be qualified with a wellname, if the wellname has a wildcard we will get a
well set, if the wellname is fully qualified we have a scalar:
UDQ
DEFINE WUWCT WWCT 'OP*' /
DEFINE FUORAT WOPR 'OPX' * 100 /
/
Here the UDQ WUCWT corresponds to the well WWCT for all wells matching the
name 'OP*', and it is undefined for the remaing wells. The UDQ FUORAT is a
scalar, given by the WOPR of well 'OPX' - multiplied by 100.
There are clearly rules for how the different variable types can be combined
in expressions, and what will be resulting type from an expression -
unfortunately that is not yet very well implemented in the opm codebase. In
UDQParser.cpp there is a function static_type_check and in UDQDefine there is
a function dynamic_type_check - these functions try to verfiy that the type
conversions are legitimate, but currently they are woefully inadequate.
*/
enum class UDQVarType {

View File

@@ -78,8 +78,29 @@ UDQASTNode::UDQASTNode(UDQTokenType type_arg,
}
namespace {
bool is_scalar(UDQVarType var_type, const std::vector<std::string>& selector)
{
if (var_type == UDQVarType::SCALAR)
return true;
if (var_type == UDQVarType::FIELD_VAR)
return true;
if (var_type == UDQVarType::WELL_VAR) {
if (selector.empty())
return false;
return (selector[0].find("*") == std::string::npos);
}
if (var_type == UDQVarType::GROUP_VAR) {
if (selector.empty())
return false;
return (selector[0].find("*") == std::string::npos);
}
}
}
UDQASTNode::UDQASTNode(UDQTokenType type_arg,
const std::string& string_value_arg,
@@ -94,6 +115,13 @@ UDQASTNode::UDQASTNode(UDQTokenType type_arg,
if (type_arg == UDQTokenType::ecl_expr)
this->var_type = UDQ::targetType(string_value);
if (this->var_type == UDQVarType::CONNECTION_VAR ||
this->var_type == UDQVarType::REGION_VAR ||
this->var_type == UDQVarType::SEGMENT_VAR ||
this->var_type == UDQVarType::AQUIFER_VAR ||
this->var_type == UDQVarType::BLOCK_VAR)
throw std::logic_error("UDQ variable of type: " + UDQ::typeName(this->var_type) + " not yet supported in flow");
}
@@ -105,12 +133,12 @@ UDQSet UDQASTNode::eval(UDQVarType target_type, const UDQContext& context) const
const auto& wells = context.wells();
if (this->selector.size() > 0) {
int fnmatch_flags = 0;
const std::string& well_pattern = this->selector[0];
if (well_pattern.find("*") == std::string::npos)
return UDQSet::wells(this->string_value, wells, context.get_well_var(well_pattern, this->string_value));
return UDQSet::scalar(this->string_value, context.get_well_var(well_pattern, this->string_value));
else {
auto res = UDQSet::wells(this->string_value, wells);
int fnmatch_flags = 0;
for (const auto& well : wells) {
if (fnmatch(well_pattern.c_str(), well.c_str(), fnmatch_flags) == 0) {
if (context.has_well_var(well, this->string_value))
@@ -130,24 +158,14 @@ UDQSet UDQASTNode::eval(UDQVarType target_type, const UDQContext& context) const
}
if (this->var_type == UDQVarType::GROUP_VAR) {
const auto& groups = context.groups();
if (this->selector.size() > 0) {
int fnmatch_flags = 0;
const std::string& group_pattern = this->selector[0];
if (group_pattern.find("*") == std::string::npos)
return UDQSet::groups(this->string_value, groups, context.get_group_var(group_pattern, this->string_value));
else {
auto res = UDQSet::groups(this->string_value, groups);
for (const auto& group : groups) {
if (fnmatch(group_pattern.c_str(), group.c_str(), fnmatch_flags) == 0) {
if (context.has_group_var(group, this->string_value))
res.assign(group, context.get_group_var(group, this->string_value));
}
}
return res;
}
return UDQSet::scalar(this->string_value, context.get_group_var(group_pattern, this->string_value));
else
throw std::logic_error("Group names with wildcards is not yet supported");
} else {
const auto& groups = context.groups();
auto res = UDQSet::groups(this->string_value, groups);
for (const auto& group : groups) {
if (context.has_group_var(group, this->string_value))
@@ -156,6 +174,11 @@ UDQSet UDQASTNode::eval(UDQVarType target_type, const UDQContext& context) const
return res;
}
}
if (this->var_type == UDQVarType::FIELD_VAR)
return UDQSet::scalar(this->string_value, context.get(this->string_value));
throw std::logic_error("Should not be here: var_type: " + UDQ::typeName(this->var_type));
}

View File

@@ -49,8 +49,8 @@ private:
void func_tokens(std::set<UDQTokenType>& tokens) const;
std::string string_value;
double scalar_value;
std::vector<std::string> selector;
double scalar_value;
std::vector<UDQASTNode> arglist;
};

View File

@@ -262,9 +262,8 @@ bool static_type_check(UDQVarType lhs, UDQVarType rhs) {
UDQASTNode UDQParser::parse(const UDQParams& udq_params, UDQVarType target_type, const std::string& target_var, const std::vector<std::string>& tokens, const ParseContext& parseContext, ErrorGuard& errors)
{
UDQParser parser(udq_params, tokens);
UDQParser parser(udq_params, target_type, tokens);
parser.next();
auto tree = parser.parse_cmp();
auto current = parser.current();

View File

@@ -66,9 +66,10 @@ public:
static UDQASTNode parse(const UDQParams& udq_params, UDQVarType target_type, const std::string& target_var, const std::vector<std::string>& tokens, const ParseContext& parseContext, ErrorGuard& errors);
private:
UDQParser(const UDQParams& udq_params1, const std::vector<std::string>& tokens1) :
UDQParser(const UDQParams& udq_params1, UDQVarType ttype, const std::vector<std::string>& tokens1) :
udq_params(udq_params1),
udqft(UDQFunctionTable(udq_params)),
target_type(ttype),
tokens(tokens1)
{}
@@ -85,6 +86,7 @@ private:
const UDQParams& udq_params;
UDQFunctionTable udqft;
UDQVarType target_type;
std::vector<std::string> tokens;
ssize_t current_pos = -1;
};

View File

@@ -381,9 +381,37 @@ UDQScalar operator/(double lhs, const UDQScalar& rhs) {
return result;
}
/*-----------------------------------------------------------------*/
namespace {
/*
If one result set is scalar and the other represents a set of wells/groups,
the scalar result is promoted to a set of the right type.
*/
UDQSet udq_cast(const UDQSet& lhs, const UDQSet& rhs)
{
if (lhs.size() != rhs.size()) {
if (lhs.var_type() != UDQVarType::SCALAR) {
printf("s1: %ld s2: %ld \n", lhs.size(), rhs.size());
throw std::logic_error("Type/size mismatch");
}
if (rhs.var_type() == UDQVarType::WELL_VAR)
return UDQSet::wells(lhs.name(), rhs.wgnames(), lhs[0].value());
if (rhs.var_type() == UDQVarType::GROUP_VAR)
return UDQSet::groups(lhs.name(), rhs.wgnames(), lhs[0].value());
throw std::logic_error("Don't have a clue");
} else
return lhs;
}
}
UDQSet operator+(const UDQSet&lhs, const UDQSet& rhs) {
UDQSet sum = lhs;
UDQSet sum = udq_cast(lhs, rhs);
sum += rhs;
return sum;
}
@@ -401,9 +429,9 @@ UDQSet operator+(double lhs, const UDQSet& rhs) {
}
UDQSet operator-(const UDQSet&lhs, const UDQSet& rhs) {
UDQSet sum = lhs;
sum -= rhs;
return sum;
UDQSet diff = udq_cast(lhs, rhs);
diff -= rhs;
return diff;
}
UDQSet operator-(const UDQSet&lhs, double rhs) {
@@ -419,15 +447,15 @@ UDQSet operator-(double lhs, const UDQSet& rhs) {
}
UDQSet operator*(const UDQSet&lhs, const UDQSet& rhs) {
UDQSet sum = lhs;
sum *= rhs;
return sum;
UDQSet prod = udq_cast(lhs, rhs);
prod *= rhs;
return prod;
}
UDQSet operator*(const UDQSet&lhs, double rhs) {
UDQSet sum = lhs;
sum *= rhs;
return sum;
UDQSet prod = lhs;
prod *= rhs;
return prod;
}
UDQSet operator*(double lhs, const UDQSet& rhs) {
@@ -437,15 +465,15 @@ UDQSet operator*(double lhs, const UDQSet& rhs) {
}
UDQSet operator/(const UDQSet&lhs, const UDQSet& rhs) {
UDQSet sum = lhs;
sum /= rhs;
return sum;
UDQSet frac = udq_cast(lhs, rhs);
frac /= rhs;
return frac;
}
UDQSet operator/(const UDQSet&lhs, double rhs) {
UDQSet sum = lhs;
sum /= rhs;
return sum;
UDQSet frac = lhs;
frac /= rhs;
return frac;
}
UDQSet operator/(double lhs, const UDQSet&rhs) {

View File

@@ -34,6 +34,52 @@ Copyright 2018 Statoil ASA.
using namespace Opm;
BOOST_AUTO_TEST_CASE(TEST)
{
UDQParams udqp;
UDQFunctionTable udqft;
UDQDefine def1(udqp, "WUWI3", {"GOPR" , "MAU", "*", "2.0", "*", "0.25", "*", "10"});
UDQDefine def2(udqp, "WUWI3", {"2.0", "*", "0.25", "*", "3"});
UDQDefine def3(udqp, "WUWI3", {"GOPR" , "FIELD", "-", "2.0", "*", "3"});
UDQDefine def4(udqp, "WUWI3", {"FOPR" , "/", "2"});
SummaryState st(std::chrono::system_clock::now());
UDQContext context(udqft, st);
st.update_group_var("MAU", "GOPR", 4);
st.update_group_var("XXX", "GOPR", 5);
st.update_group_var("FIELD", "GOPR", 6);
st.update_well_var("P1", "WBHP", 0.5);
st.update_well_var("P2", "WBHP", 0.5);
st.update("FOPR", 2);
auto res1 = def1.eval(context);
BOOST_CHECK_EQUAL( res1["P1"].value(), 20 );
BOOST_CHECK_EQUAL( res1["P2"].value(), 20 );
auto res2 = def2.eval(context);
BOOST_CHECK_EQUAL( res2["P1"].value(), 1.50 );
BOOST_CHECK_EQUAL( res2["P2"].value(), 1.50 );
auto res3 = def3.eval(context);
BOOST_CHECK_EQUAL( res3["P1"].value(), 0.00 );
BOOST_CHECK_EQUAL( res3["P2"].value(), 0.00 );
auto res4 = def4.eval(context);
BOOST_CHECK_EQUAL( res4["P1"].value(), 1.00 );
BOOST_CHECK_EQUAL( res4["P2"].value(), 1.00 );
/*
This expression has a well set as target type, and involves group with
wildcard that is not supported by flow.
*/
BOOST_CHECK_THROW( UDQDefine(udqp, "WUWI2", {"GOPR", "G*", "*", "2.0"}), std::logic_error);
/*
UDQVarType == BLOCK is not yet supported.
*/
BOOST_CHECK_THROW( UDQDefine(udqp, "WUWI2", {"BPR", "1","1", "1", "*", "2.0"}), std::logic_error);
}
Schedule make_schedule(const std::string& input) {
Parser parser;
@@ -453,7 +499,7 @@ BOOST_AUTO_TEST_CASE(UDQ_SET) {
s1.assign(0,0.0);
{
UDQSet s2("NAME", 6);
BOOST_REQUIRE_THROW(s1 + s2, std::invalid_argument);
BOOST_REQUIRE_THROW(s1 + s2, std::logic_error);
}
{
UDQSet s2("NAME", 5);
@@ -568,9 +614,7 @@ BOOST_AUTO_TEST_CASE(CMP_FUNCTIONS) {
arg2.assign(2, 2.5);
arg2.assign(4, 4.0);
BOOST_CHECK_THROW(UDQBinaryFunction::EQ(0.25, arg1, arg3), std::invalid_argument);
BOOST_CHECK_THROW(UDQBinaryFunction::EQ(0.25, arg1, arg3), std::logic_error);
{
auto result = UDQBinaryFunction::EQ(0, arg1, arg2);