move actionOnBrokenConstraints into BlackoilWellModelConstraints

This commit is contained in:
Arne Morten Kvarving 2022-10-24 09:36:05 +02:00
parent 82cf0deef4
commit e04f20af8e
4 changed files with 105 additions and 88 deletions

View File

@ -29,6 +29,7 @@
#include <opm/simulators/wells/WellGroupHelpers.hpp> #include <opm/simulators/wells/WellGroupHelpers.hpp>
#include <opm/simulators/wells/WellInterfaceGeneric.hpp> #include <opm/simulators/wells/WellInterfaceGeneric.hpp>
#include <sstream>
#include <stdexcept> #include <stdexcept>
namespace Opm { namespace Opm {
@ -374,4 +375,77 @@ checkGroupConstraints(const Group& group,
return violated; return violated;
} }
void BlackoilWellModelConstraints::
actionOnBrokenConstraints(const Group& group,
const Group::InjectionCMode& newControl,
const Phase& controlPhase,
GroupState& group_state,
DeferredLogger& deferred_logger) const
{
auto oldControl = wellModel_.groupState().injection_control(group.name(), controlPhase);
std::ostringstream ss;
if (oldControl != newControl) {
const std::string from = Group::InjectionCMode2String(oldControl);
ss << "Switching injection control mode for group "<< group.name()
<< " from " << Group::InjectionCMode2String(oldControl)
<< " to " << Group::InjectionCMode2String(newControl);
group_state.injection_control(group.name(), controlPhase, newControl);
}
if (!ss.str().empty() && wellModel_.comm().rank() == 0)
deferred_logger.info(ss.str());
}
void BlackoilWellModelConstraints::
actionOnBrokenConstraints(const Group& group,
const Group::ExceedAction& exceed_action,
const Group::ProductionCMode& newControl,
GroupState& group_state,
DeferredLogger& deferred_logger) const
{
const Group::ProductionCMode oldControl = wellModel_.groupState().production_control(group.name());
std::ostringstream ss;
switch(exceed_action) {
case Group::ExceedAction::NONE: {
if (oldControl != newControl && oldControl != Group::ProductionCMode::NONE) {
ss << "Group production exceed action is NONE for group " + group.name() + ". Nothing happens.";
}
break;
}
case Group::ExceedAction::CON: {
OPM_DEFLOG_THROW(std::runtime_error, "Group " + group.name() + "GroupProductionExceedLimit CON not implemented", deferred_logger);
break;
}
case Group::ExceedAction::CON_PLUS: {
OPM_DEFLOG_THROW(std::runtime_error, "Group " + group.name() + "GroupProductionExceedLimit CON_PLUS not implemented", deferred_logger);
break;
}
case Group::ExceedAction::WELL: {
OPM_DEFLOG_THROW(std::runtime_error, "Group " + group.name() + "GroupProductionExceedLimit WELL not implemented", deferred_logger);
break;
}
case Group::ExceedAction::PLUG: {
OPM_DEFLOG_THROW(std::runtime_error, "Group " + group.name() + "GroupProductionExceedLimit PLUG not implemented", deferred_logger);
break;
}
case Group::ExceedAction::RATE: {
if (oldControl != newControl) {
group_state.production_control(group.name(), newControl);
ss << "Switching production control mode for group "<< group.name()
<< " from " << Group::ProductionCMode2String(oldControl)
<< " to " << Group::ProductionCMode2String(newControl);
}
break;
}
default:
throw("Invalid procedure for maximum rate limit selected for group" + group.name());
}
if (!ss.str().empty() && wellModel_.comm().rank() == 0)
deferred_logger.debug(ss.str());
}
} }

View File

@ -31,6 +31,7 @@ namespace Opm {
class BlackoilWellModelGeneric; class BlackoilWellModelGeneric;
class DeferredLogger; class DeferredLogger;
class GroupState;
class SummaryState; class SummaryState;
/// Class for handling constraints for the blackoil well model. /// Class for handling constraints for the blackoil well model.
@ -62,6 +63,20 @@ public:
const int reportStepIdx, const int reportStepIdx,
DeferredLogger& deferred_logger) const; DeferredLogger& deferred_logger) const;
//! \brief Execute action for broken constraint for an injection well group.
void actionOnBrokenConstraints(const Group& group,
const Group::InjectionCMode& newControl,
const Phase& controlPhase,
GroupState& group_state,
DeferredLogger& deferred_logger) const;
//! \brief Execute action on broken constraint for a production well group.
void actionOnBrokenConstraints(const Group& group,
const Group::ExceedAction& exceed_action,
const Group::ProductionCMode& newControl,
GroupState& group_state,
DeferredLogger& deferred_logger) const;
private: private:
const BlackoilWellModelGeneric& wellModel_; //!< Reference to well model const BlackoilWellModelGeneric& wellModel_; //!< Reference to well model
}; };

View File

@ -723,7 +723,10 @@ checkGroupHigherConstraints(const Group& group,
deferred_logger); deferred_logger);
if (is_changed) { if (is_changed) {
switched_inj_groups_.insert({ {group.name(), phase}, Group::InjectionCMode2String(Group::InjectionCMode::FLD)}); switched_inj_groups_.insert({ {group.name(), phase}, Group::InjectionCMode2String(Group::InjectionCMode::FLD)});
actionOnBrokenConstraints(group, Group::InjectionCMode::FLD, phase, deferred_logger); BlackoilWellModelConstraints(*this).
actionOnBrokenConstraints(group, Group::InjectionCMode::FLD,
phase, this->groupState(),
deferred_logger);
WellGroupHelpers::updateWellRatesFromGroupTargetScale(scaling_factor, group, schedule(), reportStepIdx, /* isInjector */ true, this->groupState(), this->wellState()); WellGroupHelpers::updateWellRatesFromGroupTargetScale(scaling_factor, group, schedule(), reportStepIdx, /* isInjector */ true, this->groupState(), this->wellState());
changed = true; changed = true;
} }
@ -762,7 +765,11 @@ checkGroupHigherConstraints(const Group& group,
if (is_changed) { if (is_changed) {
switched_prod_groups_.insert({group.name(), Group::ProductionCMode2String(Group::ProductionCMode::FLD)}); switched_prod_groups_.insert({group.name(), Group::ProductionCMode2String(Group::ProductionCMode::FLD)});
const auto exceed_action = group.productionControls(summaryState_).exceed_action; const auto exceed_action = group.productionControls(summaryState_).exceed_action;
actionOnBrokenConstraints(group, exceed_action, Group::ProductionCMode::FLD, deferred_logger); BlackoilWellModelConstraints(*this).
actionOnBrokenConstraints(group, exceed_action,
Group::ProductionCMode::FLD,
this->groupState(),
deferred_logger);
WellGroupHelpers::updateWellRatesFromGroupTargetScale(scaling_factor, group, schedule(), reportStepIdx, /* isInjector */ false, this->groupState(), this->wellState()); WellGroupHelpers::updateWellRatesFromGroupTargetScale(scaling_factor, group, schedule(), reportStepIdx, /* isInjector */ false, this->groupState(), this->wellState());
changed = true; changed = true;
} }
@ -792,7 +799,9 @@ updateGroupIndividualControl(const Group& group,
if (changed_this.first != Group::InjectionCMode::NONE) if (changed_this.first != Group::InjectionCMode::NONE)
{ {
switched_inj_groups_.insert({{group.name(), phase}, Group::InjectionCMode2String(changed_this.first)}); switched_inj_groups_.insert({{group.name(), phase}, Group::InjectionCMode2String(changed_this.first)});
actionOnBrokenConstraints(group, changed_this.first, phase, deferred_logger); BlackoilWellModelConstraints(*this).
actionOnBrokenConstraints(group, changed_this.first, phase,
this->groupState(), deferred_logger);
WellGroupHelpers::updateWellRatesFromGroupTargetScale(changed_this.second, group, schedule(), reportStepIdx, /* isInjector */ false, this->groupState(), this->wellState()); WellGroupHelpers::updateWellRatesFromGroupTargetScale(changed_this.second, group, schedule(), reportStepIdx, /* isInjector */ false, this->groupState(), this->wellState());
changed = true; changed = true;
} }
@ -806,7 +815,10 @@ updateGroupIndividualControl(const Group& group,
if (changed_this.first != Group::ProductionCMode::NONE) if (changed_this.first != Group::ProductionCMode::NONE)
{ {
switched_prod_groups_.insert({group.name(), Group::ProductionCMode2String(changed_this.first)}); switched_prod_groups_.insert({group.name(), Group::ProductionCMode2String(changed_this.first)});
actionOnBrokenConstraints(group, controls.exceed_action, changed_this.first, deferred_logger); BlackoilWellModelConstraints(*this).
actionOnBrokenConstraints(group, controls.exceed_action,
changed_this.first,
this->groupState(), deferred_logger);
WellGroupHelpers::updateWellRatesFromGroupTargetScale(changed_this.second, group, schedule(), reportStepIdx, /* isInjector */ false, this->groupState(), this->wellState()); WellGroupHelpers::updateWellRatesFromGroupTargetScale(changed_this.second, group, schedule(), reportStepIdx, /* isInjector */ false, this->groupState(), this->wellState());
changed = true; changed = true;
} }
@ -815,81 +827,6 @@ updateGroupIndividualControl(const Group& group,
return changed; return changed;
} }
void
BlackoilWellModelGeneric::
actionOnBrokenConstraints(const Group& group,
const Group::ExceedAction& exceed_action,
const Group::ProductionCMode& newControl,
DeferredLogger& deferred_logger)
{
const Group::ProductionCMode oldControl = this->groupState().production_control(group.name());
std::ostringstream ss;
switch(exceed_action) {
case Group::ExceedAction::NONE: {
if (oldControl != newControl && oldControl != Group::ProductionCMode::NONE) {
ss << "Group production exceed action is NONE for group " + group.name() + ". Nothing happens.";
}
break;
}
case Group::ExceedAction::CON: {
OPM_DEFLOG_THROW(std::runtime_error, "Group " + group.name() + "GroupProductionExceedLimit CON not implemented", deferred_logger);
break;
}
case Group::ExceedAction::CON_PLUS: {
OPM_DEFLOG_THROW(std::runtime_error, "Group " + group.name() + "GroupProductionExceedLimit CON_PLUS not implemented", deferred_logger);
break;
}
case Group::ExceedAction::WELL: {
OPM_DEFLOG_THROW(std::runtime_error, "Group " + group.name() + "GroupProductionExceedLimit WELL not implemented", deferred_logger);
break;
}
case Group::ExceedAction::PLUG: {
OPM_DEFLOG_THROW(std::runtime_error, "Group " + group.name() + "GroupProductionExceedLimit PLUG not implemented", deferred_logger);
break;
}
case Group::ExceedAction::RATE: {
if (oldControl != newControl) {
this->groupState().production_control(group.name(), newControl);
ss << "Switching production control mode for group "<< group.name()
<< " from " << Group::ProductionCMode2String(oldControl)
<< " to " << Group::ProductionCMode2String(newControl);
}
break;
}
default:
throw("Invalid procedure for maximum rate limit selected for group" + group.name());
}
Parallel::Communication cc = comm_;
if (!ss.str().empty() && cc.rank() == 0)
deferred_logger.debug(ss.str());
}
void
BlackoilWellModelGeneric::
actionOnBrokenConstraints(const Group& group,
const Group::InjectionCMode& newControl,
const Phase& controlPhase,
DeferredLogger& deferred_logger)
{
auto oldControl = this->groupState().injection_control(group.name(), controlPhase);
std::ostringstream ss;
if (oldControl != newControl) {
const std::string from = Group::InjectionCMode2String(oldControl);
ss << "Switching injection control mode for group "<< group.name()
<< " from " << Group::InjectionCMode2String(oldControl)
<< " to " << Group::InjectionCMode2String(newControl);
this->groupState().injection_control(group.name(), controlPhase, newControl);
}
Parallel::Communication cc = comm_;
if (!ss.str().empty() && cc.rank() == 0)
deferred_logger.info(ss.str());
}
void void
BlackoilWellModelGeneric:: BlackoilWellModelGeneric::
updateEclWells(const int timeStepIdx, updateEclWells(const int timeStepIdx,

View File

@ -306,15 +306,6 @@ protected:
DeferredLogger& deferred_logger, DeferredLogger& deferred_logger,
const int reportStepIdx); const int reportStepIdx);
void actionOnBrokenConstraints(const Group& group,
const Group::ExceedAction& exceed_action,
const Group::ProductionCMode& newControl,
DeferredLogger& deferred_logger);
void actionOnBrokenConstraints(const Group& group,
const Group::InjectionCMode& newControl,
const Phase& controlPhase,
DeferredLogger& deferred_logger);
void updateAndCommunicateGroupData(const int reportStepIdx, void updateAndCommunicateGroupData(const int reportStepIdx,
const int iterationIdx); const int iterationIdx);