Add well_var methods to SummaryState

This commit is contained in:
Joakim Hove 2018-11-09 10:57:33 +01:00
parent 1ba51534ae
commit ca9da29c53
4 changed files with 79 additions and 12 deletions

View File

@ -23,6 +23,8 @@
#include <string>
#include <unordered_map>
#include <ert/ecl/smspec_node.h>
namespace Opm{
@ -31,19 +33,51 @@ namespace Opm{
computed, ready to use summary values. The values will typically be used by
the UDQ, WTEST and ACTIONX calculations. Observe that all value *have been
converted to the correct output units*.
The main key used to access the content of this container is the eclipse style
colon separated string - i.e. 'WWCT:OPX' to get the watercut in well 'OPX'.
The main usage of the SummaryState class is a temporary holding ground while
assembling data for the summary output, but it is also used as a context
object when evaulating the condition in ACTIONX keywords. For that reason some
of the data is duplicated both in the general structure and a specialized
structure:
SummaryState st;
st.add_well_var("OPX", "WWCT", 0.75);
st.add("WGOR:OPY", 120);
// The WWCT:OPX key has been added with the specialized add_well_var()
// method and this data is available both with the general
// st.has("WWCT:OPX") and the specialized st.has_well_var("OPX", "WWCT");
st.has("WWCT:OPX") => True
st.has_well_var("OPX", "WWCT") => True
// The WGOR:OPY key is added with the general add("WGOR:OPY") and is *not*
// accessible through the specialized st.has_well_var("OPY", "WGOR").
st.has("WGOR:OPY") => True
st.has_well_var("OPY", "WGOR") => False
*/
class SummaryState {
public:
typedef std::unordered_map<std::string, double>::const_iterator const_iterator;
void add(const std::string& key, double value);
double get(const std::string&) const;
bool has(const std::string& key) const;
void add(const std::string& key, double value);
void add(const smspec_node_type * node_ptr, double value);
void add_well_var(const std::string& well, const std::string& var, double value);
bool has_well_var(const std::string& well, const std::string& var) const;
double get_well_var(const std::string& well, const std::string& var) const;
const_iterator begin() const;
const_iterator end() const;
private:
std::unordered_map<std::string,double> values;
std::unordered_map<std::string, std::unordered_map<std::string, double>> well_values;
};
}

View File

@ -1310,7 +1310,7 @@ Summary::Summary( const EclipseState& st,
for (const auto& pair : this->handlers->handlers) {
const auto * nodeptr = pair.first;
if (smspec_node_is_total(nodeptr))
this->prev_state.add(smspec_node_get_gen_key1(nodeptr), 0);
this->prev_state.add(nodeptr, 0);
}
}
@ -1410,7 +1410,6 @@ void Summary::add_timestep( int report_step,
for( auto& f : this->handlers->handlers ) {
const int num = smspec_node_get_num( f.first );
const auto* genkey = smspec_node_get_gen_key1( f.first );
const auto schedule_wells = find_wells( schedule, f.first, sim_step, this->regionCache );
auto eff_factors = well_efficiency_factors( f.first, schedule, schedule_wells, sim_step );
@ -1425,22 +1424,22 @@ void Summary::add_timestep( int report_step,
eff_factors});
double unit_applied_val = es.getUnits().from_si( val.unit, val.value );
if (smspec_node_is_total(f.first))
if (smspec_node_is_total(f.first)) {
const auto* genkey = smspec_node_get_gen_key1( f.first );
unit_applied_val += this->prev_state.get(genkey);
}
st.add(genkey, unit_applied_val);
st.add(f.first, unit_applied_val);
}
for( const auto& value_pair : single_values ) {
const std::string key = value_pair.first;
const auto node_pair = this->handlers->single_value_nodes.find( key );
if (node_pair != this->handlers->single_value_nodes.end()) {
const auto * nodeptr = node_pair->second;
const auto * genkey = smspec_node_get_gen_key1( nodeptr );
const auto unit = single_values_units.at( key );
double si_value = value_pair.second;
double output_value = es.getUnits().from_si(unit , si_value );
st.add(genkey, output_value);
st.add(node_pair->second, output_value);
}
}
@ -1450,13 +1449,12 @@ void Summary::add_timestep( int report_step,
const auto node_pair = this->handlers->region_nodes.find( std::make_pair(key, reg+1) );
if (node_pair != this->handlers->region_nodes.end()) {
const auto * nodeptr = node_pair->second;
const auto* genkey = smspec_node_get_gen_key1( nodeptr );
const auto unit = region_units.at( key );
assert (smspec_node_get_num( nodeptr ) - 1 == static_cast<int>(reg));
double si_value = value_pair.second[reg];
double output_value = es.getUnits().from_si(unit , si_value );
st.add(genkey, output_value);
st.add(nodeptr, output_value);
}
}
}
@ -1466,11 +1464,10 @@ void Summary::add_timestep( int report_step,
const auto node_pair = this->handlers->block_nodes.find( key );
if (node_pair != this->handlers->block_nodes.end()) {
const auto * nodeptr = node_pair->second;
const auto * genkey = smspec_node_get_gen_key1( nodeptr );
const auto unit = block_units.at( key.first );
double si_value = value_pair.second;
double output_value = es.getUnits().from_si(unit , si_value );
st.add(genkey, output_value);
st.add(nodeptr, output_value);
}
}

View File

@ -21,6 +21,14 @@
#include <opm/parser/eclipse/EclipseState/Schedule/SummaryState.hpp>
namespace Opm{
void SummaryState::add(const smspec_node_type * node_ptr, double value) {
if (smspec_node_get_var_type(node_ptr) == ECL_SMSPEC_WELL_VAR)
this->add_well_var(smspec_node_get_wgname(node_ptr),
smspec_node_get_keyword(node_ptr),
value);
else
this->add(smspec_node_get_gen_key1(node_ptr), value);
}
void SummaryState::add(const std::string& key, double value) {
this->values[key] = value;
@ -40,6 +48,28 @@ namespace Opm{
return iter->second;
}
void SummaryState::add_well_var(const std::string& well, const std::string& var, double value) {
this->add(var + ":" + well, value);
this->well_values[well][var] = value;
}
bool SummaryState::has_well_var(const std::string& well, const std::string& var) const {
const auto& well_iter = this->well_values.find(well);
if (well_iter == this->well_values.end())
return false;
const auto& var_iter = well_iter->second.find(var);
if (var_iter == well_iter->second.end())
return false;
return true;
}
double SummaryState::get_well_var(const std::string& well, const std::string& var) const {
return this->well_values.at(well).at(var);
}
SummaryState::const_iterator SummaryState::begin() const {
return this->values.begin();
}

View File

@ -1306,6 +1306,12 @@ BOOST_AUTO_TEST_CASE(Test_SummaryState) {
BOOST_CHECK_THROW(st.get("NO_SUCH_KEY"), std::invalid_argument);
BOOST_CHECK(st.has("WWCT:OP_2"));
BOOST_CHECK(!st.has("NO_SUCH_KEY"));
st.add_well_var("OP1", "WWCT", 0.75);
BOOST_CHECK( st.has_well_var("OP1", "WWCT"));
BOOST_CHECK_EQUAL( st.get_well_var("OP1", "WWCT"), 0.75);
BOOST_CHECK_EQUAL( st.get_well_var("OP1", "WWCT"), st.get("WWCT:OP1"));
}
BOOST_AUTO_TEST_SUITE_END()