mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Support instantaneous flow rates in extended network (WEFAC and GEFAC item 3)
This commit is contained in:
@@ -32,6 +32,12 @@ const SupportedKeywordItems<std::string>&
|
||||
fullySupported()
|
||||
{
|
||||
static const SupportedKeywordItems<std::string> fully_supported_keywords_strings = {
|
||||
{
|
||||
"GEFAC",
|
||||
{
|
||||
{3,{true, is_bool_convertible {}, "GEFAC(GRPNETWK): String value must be convertible to bool."}}, // USE_GEFAC_IN_NETWORK
|
||||
},
|
||||
},
|
||||
{
|
||||
"NEXTSTEP",
|
||||
{
|
||||
@@ -44,6 +50,12 @@ fullySupported()
|
||||
{3,{true, allow_values<std::string> {"ORAT", "WRAT", "GRAT", "LRAT", "RESV", "BHP"}, "WCONHIST(TARGET): should be set to ORAT/WRAT/GRAT/LRAT/RESV or BHP"}}, // CMODE
|
||||
},
|
||||
},
|
||||
{
|
||||
"WEFAC",
|
||||
{
|
||||
{3,{true, is_bool_convertible {}, "WEFAC(WELNETWK): String value must be convertible to bool."}}, // USE_WEFAC_IN_NETWORK
|
||||
},
|
||||
},
|
||||
};
|
||||
|
||||
return fully_supported_keywords_strings;
|
||||
|
||||
@@ -111,12 +111,6 @@ partiallySupported()
|
||||
{8,{true, allow_values<std::string> {"NO"}, "GECON(ENDRUN): End run not implemented"}},
|
||||
},
|
||||
},
|
||||
{
|
||||
"GEFAC",
|
||||
{
|
||||
{3,{true, allow_values<std::string> {"YES"}, "GEFAC(GRPNETWK): Extended Network Model efficiency NO option not implemented"}}, // TRANSFER_EXT_NET
|
||||
},
|
||||
},
|
||||
{
|
||||
"GRIDOPTS",
|
||||
{
|
||||
@@ -273,12 +267,6 @@ partiallySupported()
|
||||
{5,{true, allow_values<std::string> {"NO"}, "WAGHYSTR(WATER_MODEL): only the NO option is supported – will STOP"}}, // WATER_MODEL
|
||||
},
|
||||
},
|
||||
{
|
||||
"WEFAC",
|
||||
{
|
||||
{3,{true, allow_values<std::string> {"YES"}, "WEFAC(WELNETWK): only the YES option is supported"}}, // EXTENDED_NETWORK_OPT
|
||||
},
|
||||
},
|
||||
{
|
||||
"WELSPECS",
|
||||
{
|
||||
|
||||
@@ -41,6 +41,7 @@ GroupState<Scalar> GroupState<Scalar>::serializationTestObject()
|
||||
{
|
||||
GroupState result(3);
|
||||
result.m_production_rates = {{"test1", {1.0, 2.0}}};
|
||||
result.m_network_production_rates={{"test1", {1.0, 20}}};
|
||||
result.production_controls = {{"test2", Group::ProductionCMode::LRAT}};
|
||||
result.prod_red_rates = {{"test3", {3.0, 4.0, 5.0}}};
|
||||
result.inj_red_rates = {{"test4", {6.0, 7.0}}};
|
||||
@@ -61,6 +62,7 @@ template<class Scalar>
|
||||
bool GroupState<Scalar>::operator==(const GroupState& other) const
|
||||
{
|
||||
return this->m_production_rates == other.m_production_rates &&
|
||||
this->m_network_production_rates == other.m_network_production_rates &&
|
||||
this->production_controls == other.production_controls &&
|
||||
this->prod_red_rates == other.prod_red_rates &&
|
||||
this->inj_red_rates == other.inj_red_rates &&
|
||||
@@ -83,6 +85,13 @@ bool GroupState<Scalar>::has_production_rates(const std::string& gname) const
|
||||
return (group_iter != this->m_production_rates.end());
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
bool GroupState<Scalar>::has_network_production_rates(const std::string& gname) const
|
||||
{
|
||||
auto group_iter = this->m_network_production_rates.find(gname);
|
||||
return (group_iter != this->m_network_production_rates.end());
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
void GroupState<Scalar>::update_production_rates(const std::string& gname,
|
||||
const std::vector<Scalar>& rates)
|
||||
@@ -93,6 +102,16 @@ void GroupState<Scalar>::update_production_rates(const std::string& gname,
|
||||
this->m_production_rates[gname] = rates;
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
void GroupState<Scalar>::update_network_production_rates(const std::string& gname,
|
||||
const std::vector<Scalar>& rates)
|
||||
{
|
||||
if (rates.size() != this->num_phases)
|
||||
throw std::logic_error("Wrong number of phases");
|
||||
|
||||
this->m_network_production_rates[gname] = rates;
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
const std::vector<Scalar>&
|
||||
GroupState<Scalar>::production_rates(const std::string& gname) const
|
||||
@@ -104,6 +123,17 @@ GroupState<Scalar>::production_rates(const std::string& gname) const
|
||||
return group_iter->second;
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
const std::vector<Scalar>&
|
||||
GroupState<Scalar>::network_production_rates(const std::string& gname) const
|
||||
{
|
||||
auto group_iter = this->m_network_production_rates.find(gname);
|
||||
if (group_iter == this->m_network_production_rates.end())
|
||||
throw std::logic_error("No such group");
|
||||
|
||||
return group_iter->second;
|
||||
}
|
||||
|
||||
//-------------------------------------------------------------------------
|
||||
|
||||
template<class Scalar>
|
||||
|
||||
@@ -50,9 +50,13 @@ public:
|
||||
bool operator==(const GroupState& other) const;
|
||||
|
||||
bool has_production_rates(const std::string& gname) const;
|
||||
bool has_network_production_rates(const std::string& gname) const;
|
||||
void update_production_rates(const std::string& gname,
|
||||
const std::vector<Scalar>& rates);
|
||||
void update_network_production_rates(const std::string& gname,
|
||||
const std::vector<Scalar>& rates);
|
||||
const std::vector<Scalar>& production_rates(const std::string& gname) const;
|
||||
const std::vector<Scalar>& network_production_rates(const std::string& gname) const;
|
||||
|
||||
void update_well_group_thp(const std::string& gname, const double& thp);
|
||||
Scalar well_group_thp(const std::string& gname) const;
|
||||
@@ -130,6 +134,7 @@ public:
|
||||
|
||||
auto forAllGroupData = [&](auto& func) {
|
||||
iterateContainer(m_production_rates, func);
|
||||
iterateContainer(m_network_production_rates, func);
|
||||
iterateContainer(prod_red_rates, func);
|
||||
iterateContainer(inj_red_rates, func);
|
||||
iterateContainer(inj_resv_rates, func);
|
||||
@@ -187,6 +192,7 @@ public:
|
||||
{
|
||||
serializer(num_phases);
|
||||
serializer(m_production_rates);
|
||||
serializer(m_network_production_rates);
|
||||
serializer(production_controls);
|
||||
serializer(group_thp);
|
||||
serializer(prod_red_rates);
|
||||
@@ -205,6 +211,7 @@ public:
|
||||
private:
|
||||
std::size_t num_phases{};
|
||||
std::map<std::string, std::vector<Scalar>> m_production_rates;
|
||||
std::map<std::string, std::vector<Scalar>> m_network_production_rates;
|
||||
std::map<std::string, Group::ProductionCMode> production_controls;
|
||||
std::map<std::string, std::vector<Scalar>> prod_red_rates;
|
||||
std::map<std::string, std::vector<Scalar>> inj_red_rates;
|
||||
|
||||
@@ -90,14 +90,15 @@ namespace Opm {
|
||||
const Opm::WellState<Scalar>& wellState,
|
||||
const int reportStepIdx,
|
||||
const int phasePos,
|
||||
const bool injector)
|
||||
const bool injector,
|
||||
const bool network)
|
||||
{
|
||||
|
||||
Scalar rate = 0.0;
|
||||
for (const std::string& groupName : group.groups()) {
|
||||
const auto& groupTmp = schedule.getGroup(groupName, reportStepIdx);
|
||||
const auto& gefac = groupTmp.getGroupEfficiencyFactor();
|
||||
rate += gefac * sumWellPhaseRates(res_rates, groupTmp, schedule, wellState, reportStepIdx, phasePos, injector);
|
||||
const auto& gefac = groupTmp.getGroupEfficiencyFactor(network);
|
||||
rate += gefac * sumWellPhaseRates(res_rates, groupTmp, schedule, wellState, reportStepIdx, phasePos, injector, network);
|
||||
}
|
||||
|
||||
// only sum satelite production once
|
||||
@@ -126,8 +127,7 @@ namespace Opm {
|
||||
if (wellEcl.getStatus() == Opm::Well::Status::SHUT)
|
||||
continue;
|
||||
|
||||
const Scalar factor = wellEcl.getEfficiencyFactor() *
|
||||
wellState[wellEcl.name()].efficiency_scaling_factor;
|
||||
Scalar factor = wellEcl.getEfficiencyFactor(network) * wellState[wellEcl.name()].efficiency_scaling_factor;
|
||||
const auto& ws = wellState.well(well_index.value());
|
||||
if (res_rates) {
|
||||
const auto& well_rates = ws.reservoir_rates;
|
||||
@@ -742,6 +742,17 @@ updateGroupProductionRates(const Group& group,
|
||||
rates[phase] = sumWellPhaseRates(false, group, schedule, wellState, reportStepIdx, phase, /*isInjector*/ false);
|
||||
}
|
||||
group_state.update_production_rates(group.name(), rates);
|
||||
|
||||
const auto& network = schedule[reportStepIdx].network();
|
||||
if (network.active()) {
|
||||
std::vector<Scalar> network_rates = rates;
|
||||
if (network.needs_instantaneous_rates(schedule, reportStepIdx)) {
|
||||
for (int phase = 0; phase < np; ++phase) {
|
||||
network_rates[phase] = sumWellPhaseRates(false, group, schedule, wellState, reportStepIdx, phase, /*isInjector*/ false, /*network*/ true);
|
||||
}
|
||||
}
|
||||
group_state.update_network_production_rates(group.name(), network_rates);
|
||||
}
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
@@ -927,24 +938,16 @@ computeNetworkPressures(const Network::ExtNetwork& network,
|
||||
node_inflows[node] = zero_rates;
|
||||
continue;
|
||||
}
|
||||
node_inflows[node] = group_state.network_production_rates(node);
|
||||
|
||||
node_inflows[node] = group_state.production_rates(node);
|
||||
// Add the ALQ amounts to the gas rates if requested.
|
||||
if (network.node(node).add_gas_lift_gas()) {
|
||||
const auto& group = schedule.getGroup(node, report_time_step);
|
||||
for (const std::string& wellname : group.wells()) {
|
||||
const Well& well = schedule.getWell(wellname, report_time_step);
|
||||
|
||||
if (well.isInjector() || !well_state.isOpen(wellname)) continue;
|
||||
// Here we use the efficiency unconditionally, but if WEFAC item 3
|
||||
// for the well is false (it defaults to true) then we should NOT use
|
||||
// the efficiency factor. Fixing this requires not only changing the
|
||||
// code here, but also:
|
||||
// - Adding a member to the well for this flag, and setting it in Schedule::handleWEFAC().
|
||||
// - Making the wells' maximum flows (i.e. not time-averaged by using a efficiency factor)
|
||||
// available and using those (for wells with WEFAC(3) true only) when accumulating group
|
||||
// rates, but ONLY for network calculations.
|
||||
const Scalar efficiency = well.getEfficiencyFactor() *
|
||||
well_state.getGlobalEfficiencyScalingFactor(wellname);
|
||||
const Scalar efficiency = well.getEfficiencyFactor(/*network*/ true) * well_state.getGlobalEfficiencyScalingFactor(wellname);
|
||||
node_inflows[node][BlackoilPhases::Vapour] += well_state.getALQ(wellname) * efficiency;
|
||||
}
|
||||
}
|
||||
@@ -962,13 +965,17 @@ computeNetworkPressures(const Network::ExtNetwork& network,
|
||||
// Add downbranch rates to upbranch.
|
||||
std::vector<Scalar>& up = node_inflows[(*upbranch).uptree_node()];
|
||||
const std::vector<Scalar>& down = node_inflows[node];
|
||||
// @TODO@ Also support NEFAC (for nodes that do not correspond to groups)
|
||||
double efficiency = 1.0;
|
||||
if (schedule.hasGroup(node, report_time_step)) {
|
||||
efficiency = schedule.getGroup(node, report_time_step).getGroupEfficiencyFactor(/*network*/ true);
|
||||
}
|
||||
if (up.empty()) {
|
||||
up = down;
|
||||
} else {
|
||||
assert (up.size() == down.size());
|
||||
for (std::size_t ii = 0; ii < up.size(); ++ii) {
|
||||
up[ii] += down[ii];
|
||||
}
|
||||
up = std::vector<Scalar>(down.size(), 0.0);
|
||||
}
|
||||
assert (up.size() == down.size());
|
||||
for (std::size_t ii = 0; ii < up.size(); ++ii) {
|
||||
up[ii] += (efficiency*down[ii]);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
@@ -56,7 +56,8 @@ public:
|
||||
const Opm::WellState<Scalar>& wellState,
|
||||
const int reportStepIdx,
|
||||
const int phasePos,
|
||||
const bool injector);
|
||||
const bool injector,
|
||||
const bool network = false);
|
||||
|
||||
static Scalar satelliteProduction(const ScheduleState& sched,
|
||||
const std::vector<std::string>& groups,
|
||||
|
||||
Reference in New Issue
Block a user