Changes to Well status implementation
o The status of a well is maintened as a small object which is managed by a new std::shared_ptr member in the well objects. The consequence of this is that several well objects can share the same underlying status object. The advantage of this is that runtime (i.e. ACTIONX) updates of well status will affect the correct set of wells. o The general Schedule::updateWell() will use the DynamicState::upadte_equal()
This commit is contained in:
parent
85e2c641d3
commit
2e832f0fbe
@ -404,7 +404,7 @@ namespace Opm
|
||||
void updateGroup(std::shared_ptr<Group> group, std::size_t reportStep);
|
||||
bool checkGroups(const ParseContext& parseContext, ErrorGuard& errors);
|
||||
void updateUDQActive( std::size_t timeStep, std::shared_ptr<UDQActive> udq );
|
||||
bool updateWellStatus( const std::string& well, std::size_t reportStep , Well::Status status, bool update_connections, std::optional<KeywordLocation> = {});
|
||||
bool updateWellStatus( const std::string& well, std::size_t reportStep, bool runtime, Well::Status status, bool update_connections, std::optional<KeywordLocation> = {});
|
||||
void addWellToGroup( const std::string& group_name, const std::string& well_name , std::size_t timeStep);
|
||||
void iterateScheduleSection(std::shared_ptr<const Python> python, const std::string& input_path, const ParseContext& parseContext , ErrorGuard& errors, const SCHEDULESection& , const EclipseGrid& grid,
|
||||
const FieldPropsManager& fp);
|
||||
@ -460,7 +460,7 @@ namespace Opm
|
||||
|
||||
void applyEXIT(const DeckKeyword&, std::size_t currentStep);
|
||||
void applyMESSAGES(const DeckKeyword&, std::size_t currentStep);
|
||||
void applyWELOPEN(const DeckKeyword&, std::size_t currentStep, const ParseContext&, ErrorGuard&, const std::vector<std::string>& matching_wells = {});
|
||||
void applyWELOPEN(const DeckKeyword&, std::size_t currentStep, bool runtime, const ParseContext&, ErrorGuard&, const std::vector<std::string>& matching_wells = {});
|
||||
void applyWRFT(const DeckKeyword&, std::size_t currentStep);
|
||||
void applyWRFTPLT(const DeckKeyword&, std::size_t currentStep);
|
||||
|
||||
|
@ -79,6 +79,31 @@ public:
|
||||
static Status StatusFromString(const std::string& stringValue);
|
||||
|
||||
|
||||
struct WellStatus {
|
||||
Status status;
|
||||
|
||||
WellStatus() = default;
|
||||
|
||||
WellStatus(Status st) :
|
||||
status(st)
|
||||
{}
|
||||
|
||||
template<class Serializer>
|
||||
void serializeOp(Serializer& serializer)
|
||||
{
|
||||
serializer(status);
|
||||
}
|
||||
|
||||
bool operator==(const WellStatus& other) const {
|
||||
return this->status == other.status;
|
||||
}
|
||||
|
||||
static WellStatus serializeObject() {
|
||||
WellStatus ws(Well::Status::AUTO);
|
||||
return ws;
|
||||
}
|
||||
};
|
||||
|
||||
|
||||
|
||||
/*
|
||||
@ -545,7 +570,7 @@ public:
|
||||
void updateSegments(std::shared_ptr<WellSegments> segments_arg);
|
||||
bool updateConnections(std::shared_ptr<WellConnections> connections, bool force = false);
|
||||
bool updateConnections(std::shared_ptr<WellConnections> connections, const EclipseGrid& grid, const std::vector<int>& pvtnum);
|
||||
bool updateStatus(Status status, bool update_connections);
|
||||
bool updateStatus(Status status, bool runtime, bool update_connections);
|
||||
bool updateGroup(const std::string& group);
|
||||
bool updateWellGuideRate(bool available, double guide_rate, GuideRateTarget guide_phase, double scale_factor);
|
||||
bool updateWellGuideRate(double guide_rate);
|
||||
@ -650,7 +675,6 @@ private:
|
||||
GasInflowEquation gas_inflow = GasInflowEquation::STD; // Will NOT be loaded/assigned from restart file
|
||||
UnitSystem unit_system;
|
||||
double udq_undefined;
|
||||
Status status;
|
||||
WellType wtype;
|
||||
WellGuideRate guide_rate;
|
||||
double efficiency_factor;
|
||||
@ -669,6 +693,7 @@ private:
|
||||
std::shared_ptr<WellProductionProperties> production;
|
||||
std::shared_ptr<WellInjectionProperties> injection;
|
||||
std::shared_ptr<WellSegments> segments;
|
||||
std::shared_ptr<WellStatus> status;
|
||||
PAvg m_pavg;
|
||||
};
|
||||
|
||||
|
@ -870,7 +870,7 @@ namespace {
|
||||
const Well::Status status = Well::StatusFromString(record.getItem("STATUS").getTrimmedString(0));
|
||||
|
||||
for (const auto& well_name : well_names) {
|
||||
this->updateWellStatus( well_name , handlerContext.currentStep , status, false, handlerContext.keyword.location() );
|
||||
this->updateWellStatus( well_name , handlerContext.currentStep , false, status, false, handlerContext.keyword.location() );
|
||||
|
||||
std::optional<VFPProdTable::ALQ_TYPE> alq_type;
|
||||
auto& dynamic_state = this->wells_static.at(well_name);
|
||||
@ -922,7 +922,7 @@ namespace {
|
||||
"Well " + well2->name() + " is a history matched well with zero rate where crossflow is banned. " +
|
||||
"This well will be closed at " + std::to_string(m_timeMap.getTimePassedUntil(handlerContext.currentStep) / (60*60*24)) + " days";
|
||||
OpmLog::note(msg);
|
||||
this->updateWellStatus( well_name, handlerContext.currentStep, Well::Status::SHUT, false );
|
||||
this->updateWellStatus( well_name, handlerContext.currentStep, false, Well::Status::SHUT, false );
|
||||
}
|
||||
}
|
||||
}
|
||||
@ -939,13 +939,12 @@ namespace {
|
||||
const Well::Status status = Well::StatusFromString(record.getItem("STATUS").getTrimmedString(0));
|
||||
|
||||
for (const auto& well_name : well_names) {
|
||||
this->updateWellStatus(well_name, handlerContext.currentStep, status, false, handlerContext.keyword.location());
|
||||
bool update_well = this->updateWellStatus(well_name, handlerContext.currentStep, false, status, false, handlerContext.keyword.location());
|
||||
std::optional<VFPProdTable::ALQ_TYPE> alq_type;
|
||||
auto& dynamic_state = this->wells_static.at(well_name);
|
||||
auto well2 = std::make_shared<Well>(*dynamic_state[handlerContext.currentStep]);
|
||||
const bool switching_from_injector = !well2->isProducer();
|
||||
auto properties = std::make_shared<Well::WellProductionProperties>(well2->getProductionProperties());
|
||||
bool update_well = switching_from_injector;
|
||||
properties->clearControls();
|
||||
if (well2->isAvailableForGroupControl())
|
||||
properties->addProductionControl(Well::ProducerCMode::GRUP);
|
||||
@ -958,8 +957,10 @@ namespace {
|
||||
alq_type = this->getVFPProdTable(table_nr, handlerContext.currentStep).getALQType();
|
||||
properties->handleWCONPROD(alq_type, this->unit_system, well_name, record);
|
||||
|
||||
if (switching_from_injector)
|
||||
if (switching_from_injector) {
|
||||
properties->resetDefaultBHPLimit();
|
||||
update_well = true;
|
||||
}
|
||||
|
||||
if (well2->updateProduction(properties))
|
||||
update_well = true;
|
||||
@ -993,7 +994,7 @@ namespace {
|
||||
const Well::Status status = Well::StatusFromString(record.getItem("STATUS").getTrimmedString(0));
|
||||
|
||||
for (const auto& well_name : well_names) {
|
||||
this->updateWellStatus(well_name, handlerContext.currentStep, status, false, handlerContext.keyword.location());
|
||||
this->updateWellStatus(well_name, handlerContext.currentStep, false, status, false, handlerContext.keyword.location());
|
||||
|
||||
bool update_well = false;
|
||||
auto& dynamic_state = this->wells_static.at(well_name);
|
||||
@ -1026,14 +1027,14 @@ namespace {
|
||||
if (injection->surfaceInjectionRate.is<double>()) {
|
||||
if (injection->hasInjectionControl(Well::InjectorCMode::RATE) && injection->surfaceInjectionRate.zero()) {
|
||||
OpmLog::note(msg);
|
||||
this->updateWellStatus( well_name, handlerContext.currentStep, Well::Status::SHUT, false );
|
||||
this->updateWellStatus( well_name, handlerContext.currentStep, false, Well::Status::SHUT, false );
|
||||
}
|
||||
}
|
||||
|
||||
if (injection->reservoirInjectionRate.is<double>()) {
|
||||
if (injection->hasInjectionControl(Well::InjectorCMode::RESV) && injection->reservoirInjectionRate.zero()) {
|
||||
OpmLog::note(msg);
|
||||
this->updateWellStatus( well_name, handlerContext.currentStep, Well::Status::SHUT, false );
|
||||
this->updateWellStatus( well_name, handlerContext.currentStep, false, Well::Status::SHUT, false );
|
||||
}
|
||||
}
|
||||
}
|
||||
@ -1055,7 +1056,7 @@ namespace {
|
||||
const Well::Status status = Well::StatusFromString( record.getItem("STATUS").getTrimmedString(0));
|
||||
|
||||
for (const auto& well_name : well_names) {
|
||||
this->updateWellStatus(well_name, handlerContext.currentStep, status, false, handlerContext.keyword.location());
|
||||
this->updateWellStatus(well_name, handlerContext.currentStep, false, status, false, handlerContext.keyword.location());
|
||||
|
||||
bool update_well = false;
|
||||
auto& dynamic_state = this->wells_static.at(well_name);
|
||||
@ -1083,7 +1084,7 @@ namespace {
|
||||
"Well " + well_name + " is an injector with zero rate where crossflow is banned. " +
|
||||
"This well will be closed at " + std::to_string ( m_timeMap.getTimePassedUntil(handlerContext.currentStep) / (60*60*24) ) + " days";
|
||||
OpmLog::note(msg);
|
||||
this->updateWellStatus( well_name, handlerContext.currentStep, Well::Status::SHUT, false );
|
||||
this->updateWellStatus( well_name, handlerContext.currentStep, false, Well::Status::SHUT, false );
|
||||
}
|
||||
}
|
||||
}
|
||||
@ -1125,7 +1126,7 @@ namespace {
|
||||
}
|
||||
|
||||
void Schedule::handleWELOPEN (const HandlerContext& handlerContext, const ParseContext& parseContext, ErrorGuard& errors) {
|
||||
this->applyWELOPEN(handlerContext.keyword, handlerContext.currentStep, parseContext, errors);
|
||||
this->applyWELOPEN(handlerContext.keyword, handlerContext.currentStep, false, parseContext, errors);
|
||||
}
|
||||
|
||||
void Schedule::handleWELPI(const HandlerContext& handlerContext, const ParseContext& parseContext, ErrorGuard& errors) {
|
||||
|
@ -573,20 +573,20 @@ private:
|
||||
}
|
||||
|
||||
void Schedule::shut_well(const std::string& well_name, std::size_t report_step) {
|
||||
this->updateWellStatus(well_name, report_step, Well::Status::SHUT, true);
|
||||
this->updateWellStatus(well_name, report_step, true, Well::Status::SHUT, false);
|
||||
}
|
||||
|
||||
void Schedule::open_well(const std::string& well_name, std::size_t report_step) {
|
||||
this->updateWellStatus(well_name, report_step, Well::Status::OPEN, true);
|
||||
this->updateWellStatus(well_name, report_step, true, Well::Status::OPEN, false);
|
||||
}
|
||||
|
||||
void Schedule::stop_well(const std::string& well_name, std::size_t report_step) {
|
||||
this->updateWellStatus(well_name, report_step, Well::Status::STOP, true);
|
||||
this->updateWellStatus(well_name, report_step, true, Well::Status::STOP, false);
|
||||
}
|
||||
|
||||
void Schedule::updateWell(std::shared_ptr<Well> well, std::size_t reportStep) {
|
||||
auto& dynamic_state = this->wells_static.at(well->name());
|
||||
dynamic_state.update(reportStep, std::move(well));
|
||||
dynamic_state.update_equal(reportStep, std::move(well));
|
||||
}
|
||||
|
||||
|
||||
@ -594,7 +594,7 @@ private:
|
||||
Function is quite dangerous - because if this is called while holding a
|
||||
Well pointer that will go stale and needs to be refreshed.
|
||||
*/
|
||||
bool Schedule::updateWellStatus( const std::string& well_name, std::size_t reportStep , Well::Status status, bool update_connections, std::optional<KeywordLocation> location) {
|
||||
bool Schedule::updateWellStatus( const std::string& well_name, std::size_t reportStep , bool runtime, Well::Status status, bool update_connections, std::optional<KeywordLocation> location) {
|
||||
auto& dynamic_state = this->wells_static.at(well_name);
|
||||
auto well2 = std::make_shared<Well>(*dynamic_state[reportStep]);
|
||||
if (well2->getConnections().empty() && status == Well::Status::OPEN) {
|
||||
@ -608,14 +608,26 @@ private:
|
||||
return false;
|
||||
}
|
||||
|
||||
auto old_status = well2->getStatus();
|
||||
bool update = false;
|
||||
if (well2->updateStatus(status, update_connections)) {
|
||||
m_events.addEvent( ScheduleEvents::WELL_STATUS_CHANGE, reportStep );
|
||||
this->addWellGroupEvent( well2->name(), ScheduleEvents::WELL_STATUS_CHANGE, reportStep);
|
||||
if (well2->updateStatus(status, runtime, update_connections)) {
|
||||
this->updateWell(well2, reportStep);
|
||||
update = true;
|
||||
if (status == Well::Status::OPEN)
|
||||
this->rft_config.addWellOpen(well_name, reportStep);
|
||||
|
||||
/*
|
||||
The Well::updateStatus() will always return true because a new
|
||||
WellStatus object should be created. But the new object might have
|
||||
the same value as the previous object; therefor we need to check
|
||||
for an actual status change before we emit a WELL_STATUS_CHANGE
|
||||
event.
|
||||
*/
|
||||
if (old_status != status) {
|
||||
this->m_events.addEvent( ScheduleEvents::WELL_STATUS_CHANGE, reportStep );
|
||||
this->addWellGroupEvent( well2->name(), ScheduleEvents::WELL_STATUS_CHANGE, reportStep);
|
||||
}
|
||||
|
||||
update = true;
|
||||
}
|
||||
return update;
|
||||
}
|
||||
@ -654,7 +666,7 @@ private:
|
||||
}
|
||||
}
|
||||
|
||||
void Schedule::applyWELOPEN(const DeckKeyword& keyword, std::size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors, const std::vector<std::string>& matching_wells) {
|
||||
void Schedule::applyWELOPEN(const DeckKeyword& keyword, std::size_t currentStep, bool runtime, const ParseContext& parseContext, ErrorGuard& errors, const std::vector<std::string>& matching_wells) {
|
||||
|
||||
auto conn_defaulted = []( const DeckRecord& rec ) {
|
||||
auto defaulted = []( const DeckItem& item ) {
|
||||
@ -665,7 +677,6 @@ private:
|
||||
};
|
||||
|
||||
constexpr auto open = Well::Status::OPEN;
|
||||
bool action_mode = !matching_wells.empty();
|
||||
|
||||
for (const auto& record : keyword) {
|
||||
const auto& wellNamePattern = record.getItem( "WELL" ).getTrimmedString(0);
|
||||
@ -690,7 +701,7 @@ private:
|
||||
+ std::to_string( days ) + " days";
|
||||
OpmLog::note(msg);
|
||||
} else {
|
||||
this->updateWellStatus( wname, currentStep, well_status, false );
|
||||
this->updateWellStatus( wname, currentStep, runtime, well_status, false );
|
||||
if (well_status == open)
|
||||
this->rft_config.addWellOpen(wname, currentStep);
|
||||
}
|
||||
@ -705,7 +716,7 @@ private:
|
||||
{
|
||||
auto& dynamic_state = this->wells_static.at(wname);
|
||||
auto well_ptr = std::make_shared<Well>( *dynamic_state[currentStep] );
|
||||
if (well_ptr->handleWELOPEN(record, comp_status, action_mode))
|
||||
if (well_ptr->handleWELOPEN(record, comp_status, runtime))
|
||||
// The updateWell call breaks test at line 825 and 831 in ScheduleTests
|
||||
this->updateWell(std::move(well_ptr), currentStep);
|
||||
}
|
||||
@ -876,7 +887,7 @@ private:
|
||||
const std::string wname = well.name();
|
||||
|
||||
m_events.addEvent( ScheduleEvents::NEW_WELL , report_step );
|
||||
wellgroup_events.insert( std::make_pair(wname, Events(this->m_timeMap)));
|
||||
this->wellgroup_events.insert( std::make_pair(wname, Events(this->m_timeMap)));
|
||||
this->addWellGroupEvent(wname, ScheduleEvents::NEW_WELL, report_step);
|
||||
|
||||
well.setInsertIndex(this->wells_static.size());
|
||||
@ -1372,7 +1383,7 @@ private:
|
||||
"All completions in well " + well.name() + " is shut at " + std::to_string ( m_timeMap.getTimePassedUntil(timeStep) / (60*60*24) ) + " days. \n" +
|
||||
"The well is therefore also shut.";
|
||||
OpmLog::note(msg);
|
||||
this->updateWellStatus( well.name(), timeStep, Well::Status::SHUT, false);
|
||||
this->updateWellStatus( well.name(), timeStep, false, Well::Status::SHUT, false);
|
||||
}
|
||||
}
|
||||
}
|
||||
@ -1533,7 +1544,7 @@ private:
|
||||
this->updateUDQ(keyword, reportStep);
|
||||
|
||||
if (keyword.name() == "WELOPEN")
|
||||
this->applyWELOPEN(keyword, reportStep, parseContext, errors, result.wells());
|
||||
this->applyWELOPEN(keyword, reportStep, true, parseContext, errors, result.wells());
|
||||
|
||||
/*
|
||||
The WELPI functionality is implemented as a two-step process
|
||||
|
@ -140,7 +140,6 @@ Well::Well(const RestartIO::RstWell& rst_well,
|
||||
pvt_table(rst_well.pvt_table),
|
||||
unit_system(unit_system_arg),
|
||||
udq_undefined(udq_undefined_arg),
|
||||
status(status_from_int(rst_well.well_status)),
|
||||
wtype(rst_well.wtype),
|
||||
guide_rate(def_guide_rate),
|
||||
efficiency_factor(rst_well.efficiency_factor),
|
||||
@ -153,7 +152,8 @@ Well::Well(const RestartIO::RstWell& rst_well,
|
||||
tracer_properties(std::make_shared<WellTracerProperties>()),
|
||||
connections(std::make_shared<WellConnections>(order_from_int(rst_well.completion_ordering), headI, headJ)),
|
||||
production(std::make_shared<WellProductionProperties>(unit_system_arg, wname)),
|
||||
injection(std::make_shared<WellInjectionProperties>(unit_system_arg, wname))
|
||||
injection(std::make_shared<WellInjectionProperties>(unit_system_arg, wname)),
|
||||
status(std::make_shared<WellStatus>(status_from_int(rst_well.well_status)))
|
||||
{
|
||||
using CModeVal = ::Opm::RestartIO::Helpers::VectorItems::IWell::Value::WellCtrlMode;
|
||||
|
||||
@ -335,7 +335,6 @@ Well::Well(const std::string& wname_arg,
|
||||
gas_inflow(inflow_eq),
|
||||
unit_system(unit_system_arg),
|
||||
udq_undefined(udq_undefined_arg),
|
||||
status(Status::SHUT),
|
||||
wtype(wtype_arg),
|
||||
guide_rate({true, -1, Well::GuideRateTarget::UNDEFINED,ParserKeywords::WGRUPCON::SCALING_FACTOR::defaultValue}),
|
||||
efficiency_factor(1.0),
|
||||
@ -347,7 +346,8 @@ Well::Well(const std::string& wname_arg,
|
||||
tracer_properties(std::make_shared<WellTracerProperties>()),
|
||||
connections(std::make_shared<WellConnections>(ordering_arg, headI, headJ)),
|
||||
production(std::make_shared<WellProductionProperties>(unit_system, wname)),
|
||||
injection(std::make_shared<WellInjectionProperties>(unit_system, wname))
|
||||
injection(std::make_shared<WellInjectionProperties>(unit_system, wname)),
|
||||
status(std::make_shared<WellStatus>(Status::SHUT))
|
||||
{
|
||||
auto p = std::make_shared<WellProductionProperties>(this->unit_system, this->wname);
|
||||
p->whistctl_cmode = whistctl_cmode;
|
||||
@ -366,7 +366,7 @@ Well Well::serializeObject()
|
||||
result.ref_depth = 5;
|
||||
result.unit_system = UnitSystem::serializeObject();
|
||||
result.udq_undefined = 6.0;
|
||||
result.status = Status::SHUT;
|
||||
result.status = std::make_shared<WellStatus>(Status::SHUT);
|
||||
result.drainage_radius = 7.0;
|
||||
result.allow_cross_flow = true;
|
||||
result.automatic_shutin = false;
|
||||
@ -515,7 +515,7 @@ bool Well::updateWellProductivityIndex(const double prodIndex) {
|
||||
}
|
||||
|
||||
bool Well::updateHasProduced() {
|
||||
if (this->wtype.producer() && this->status == Status::OPEN) {
|
||||
if (this->wtype.producer() && this->getStatus() == Status::OPEN) {
|
||||
if (this->has_produced)
|
||||
return false;
|
||||
|
||||
@ -526,7 +526,7 @@ bool Well::updateHasProduced() {
|
||||
}
|
||||
|
||||
bool Well::updateHasInjected() {
|
||||
if (this->wtype.injector() && this->status == Status::OPEN) {
|
||||
if (this->wtype.injector() && this->getStatus() == Status::OPEN) {
|
||||
if (this->has_injected)
|
||||
return false;
|
||||
|
||||
@ -610,10 +610,108 @@ bool Well::updateHead(int I, int J) {
|
||||
}
|
||||
|
||||
|
||||
bool Well::updateStatus(Status well_state, bool update_connections) {
|
||||
bool update = false;
|
||||
/*
|
||||
The fileformat used by OPM/flow seems to work as an imperative programming
|
||||
language. The simulator can be percieved as a mutable imperative programming
|
||||
environment. The simulator program has an implicit DOM and the various
|
||||
keywords manipulate the elements in this DOM.
|
||||
|
||||
In opm flow the approach to the input file is not that of an imperative
|
||||
programming language, rather the entire input file is parsed and internalized
|
||||
into the mainly EclipseState and Schedule instances. For the most part this
|
||||
has worked out nicely, but some of the more advanced features of the simulator
|
||||
(notably ACTIONX) requires an interaction between the simulator and the
|
||||
Schedule datastructure which becomes awkward in the current implementation.
|
||||
E.g the complexity to shut/open a well is quite immense. To understand how
|
||||
this complexity arises it is important to understand:
|
||||
|
||||
1. How the DynamicState<T> class works.
|
||||
|
||||
2. How a new Well instance is created for each keyword which manipulates
|
||||
the well state.
|
||||
|
||||
3. How the well class uses pointer semantics to manage objects which should
|
||||
remane unchanged across several well keywords.
|
||||
|
||||
|
||||
START
|
||||
1 'JAN' 2000 /
|
||||
|
||||
SCHEDULE
|
||||
|
||||
WELSPECS
|
||||
W1 .... /
|
||||
/
|
||||
|
||||
WCONPROD
|
||||
W1 'OPEN' /
|
||||
/
|
||||
|
||||
DATES
|
||||
1 'FEB' 2000 /
|
||||
/
|
||||
|
||||
WELPI
|
||||
W1 1000 /
|
||||
/
|
||||
|
||||
DATES
|
||||
1 'MAR' 2000 /
|
||||
/
|
||||
|
||||
WCONPROD
|
||||
W1 'OPEN' /
|
||||
/
|
||||
|
||||
DATES
|
||||
1 'APR' 2000 /
|
||||
/
|
||||
|
||||
WELPI
|
||||
W1 1000 /
|
||||
/
|
||||
|
||||
0--------------------1--------------------2--------------------3-------------------->
|
||||
|
||||
[ W0 ---------------->
|
||||
|
|
||||
| [ W1 ---------------->
|
||||
| |
|
||||
| | [ W2 ---------------->
|
||||
| | |
|
||||
| | | [ W3 ---------------->
|
||||
| | | |
|
||||
| | | |
|
||||
\|/ | \|/ |
|
||||
| |
|
||||
[ WellStatus 0 ] <-----/ [ WellStatus 1] <------/
|
||||
|
||||
|
||||
This illustration shows "many things":
|
||||
|
||||
1. For each of the kewyords which manipulates wells a new well object are
|
||||
created. These are illustrated as W0, W1, W2 and W3. As illustrated the
|
||||
well objects have a validity in the time direction.
|
||||
|
||||
2. Each of the keywords which changes/sets the state of a well will create a
|
||||
new WellStatus object; these are illustrated as WellStatus 0 and WellStatus
|
||||
1. As we can see the WellStatus in general have different temporal ranges
|
||||
of validity than the well objects.
|
||||
|
||||
3. The main point of this complexity is to support runtime altering of the
|
||||
wells status - with ACTIONX or other means. If the runtime argument is true
|
||||
when calling Well::updateStatus() we update the wells status directly, and
|
||||
not go through creating a new WellStatus object. As a consequence the
|
||||
updated status will apply to all well/time points which share WellStatus
|
||||
object - this can even go backwards in time!
|
||||
*/
|
||||
|
||||
bool Well::updateStatus(Status well_state, bool runtime, bool update_connections) {
|
||||
if (update_connections) {
|
||||
Connection::State connection_state;
|
||||
if (runtime)
|
||||
throw std::logic_error("runtime and update_connections can not be combined");
|
||||
|
||||
|
||||
switch (well_state) {
|
||||
case Status::OPEN:
|
||||
@ -638,15 +736,15 @@ bool Well::updateStatus(Status well_state, bool update_connections) {
|
||||
new_connections->add(c);
|
||||
}
|
||||
|
||||
update = this->updateConnections(std::move(new_connections));
|
||||
this->updateConnections(std::move(new_connections));
|
||||
}
|
||||
|
||||
if (this->status != well_state) {
|
||||
this->status = well_state;
|
||||
update = true;
|
||||
}
|
||||
if (runtime)
|
||||
this->status->status = well_state;
|
||||
else
|
||||
this->status = std::make_shared<WellStatus>(well_state);
|
||||
|
||||
return update;
|
||||
return true;
|
||||
}
|
||||
|
||||
|
||||
@ -693,7 +791,7 @@ bool Well::updateConnections(std::shared_ptr<WellConnections> connections_arg, b
|
||||
if (force || *this->connections != *connections_arg) {
|
||||
this->connections = connections_arg;
|
||||
if (this->connections->empty())
|
||||
this->status = Status::SHUT;
|
||||
this->status = std::make_shared<WellStatus>(Well::Status::SHUT);
|
||||
|
||||
return true;
|
||||
}
|
||||
@ -957,7 +1055,7 @@ const Well::WellInjectionProperties& Well::getInjectionProperties() const {
|
||||
}
|
||||
|
||||
Well::Status Well::getStatus() const {
|
||||
return this->status;
|
||||
return this->status->status;
|
||||
}
|
||||
|
||||
const PAvg& Well::pavg() const {
|
||||
@ -1011,14 +1109,14 @@ int Well::fip_region_number() const {
|
||||
because there is some twisted logic aggregating connection changes over a
|
||||
complete report step.
|
||||
|
||||
However - when the WELOPEN is called as a ACTIONX action the full
|
||||
Schedule::iterateScheduleSection() is not run and the check if all connections
|
||||
is closed is not done. Therefor we have a action_mode flag here which makes
|
||||
sure to close the well in this case.
|
||||
However - when the WELOPEN is called runtime (typically as an ACTIONX action)
|
||||
the full Schedule::iterateScheduleSection() is not run and the check if all
|
||||
connections is closed is not done. Therefor we have a runtime flag here
|
||||
which makes sure to close the well in this case.
|
||||
*/
|
||||
|
||||
|
||||
bool Well::handleWELOPEN(const DeckRecord& record, Connection::State state_arg, bool action_mode) {
|
||||
bool Well::handleWELOPEN(const DeckRecord& record, Connection::State state_arg, bool runtime) {
|
||||
|
||||
auto match = [=]( const Connection &c) -> bool {
|
||||
if (!match_eq(c.getI(), record, "I" , -1)) return false;
|
||||
@ -1038,9 +1136,9 @@ bool Well::handleWELOPEN(const DeckRecord& record, Connection::State state_arg,
|
||||
|
||||
new_connections->add(c);
|
||||
}
|
||||
if (action_mode) {
|
||||
if (runtime) {
|
||||
if (new_connections->allConnectionsShut())
|
||||
this->status = Status::SHUT;
|
||||
this->updateStatus(Well::Status::SHUT, true, false);
|
||||
}
|
||||
|
||||
return this->updateConnections(std::move(new_connections));
|
||||
|
@ -425,6 +425,14 @@ DATES
|
||||
1 'FEB' 2015 /
|
||||
/
|
||||
|
||||
WCONPROD
|
||||
-- Item #:1 2 3 4 5 9
|
||||
'P1' 'OPEN' 'ORAT' 5000 4* 1000 /
|
||||
'P2' 'OPEN' 'ORAT' 5000 4* 1000 /
|
||||
'P3' 'OPEN' 'ORAT' 5000 4* 1000 /
|
||||
'P4' 'OPEN' 'ORAT' 5000 4* 1000 /
|
||||
/
|
||||
|
||||
DATES
|
||||
1 'MAR' 2015 /
|
||||
/
|
||||
|
@ -241,15 +241,11 @@ BOOST_AUTO_TEST_CASE(WELL_CLOSE_EXAMPLE) {
|
||||
BOOST_CHECK(w3.getStatus() == Well::Status::OPEN );
|
||||
}
|
||||
{
|
||||
const auto& w2_5 = td.schedule.getWell("P2", 5);
|
||||
const auto& w2_6 = td.schedule.getWell("P2", 6);
|
||||
BOOST_CHECK(w2_5.getStatus() == Well::Status::OPEN );
|
||||
BOOST_CHECK(w2_6.getStatus() == Well::Status::SHUT );
|
||||
}
|
||||
{
|
||||
const auto& w4_10 = td.schedule.getWell("P4", 10);
|
||||
const auto& w4_11 = td.schedule.getWell("P4", 11);
|
||||
BOOST_CHECK(w4_10.getStatus() == Well::Status::OPEN );
|
||||
BOOST_CHECK(w4_11.getStatus() == Well::Status::SHUT );
|
||||
}
|
||||
}
|
||||
@ -502,15 +498,11 @@ BOOST_AUTO_TEST_CASE(PYTHON_WELL_CLOSE_EXAMPLE) {
|
||||
BOOST_CHECK(w3.getStatus() == Well::Status::OPEN );
|
||||
}
|
||||
{
|
||||
const auto& w2_5 = td.schedule.getWell("P2", 5);
|
||||
const auto& w2_6 = td.schedule.getWell("P2", 6);
|
||||
BOOST_CHECK(w2_5.getStatus() == Well::Status::OPEN );
|
||||
BOOST_CHECK(w2_6.getStatus() == Well::Status::SHUT );
|
||||
}
|
||||
{
|
||||
const auto& w4_10 = td.schedule.getWell("P4", 10);
|
||||
const auto& w4_11 = td.schedule.getWell("P4", 11);
|
||||
BOOST_CHECK(w4_10.getStatus() == Well::Status::OPEN );
|
||||
BOOST_CHECK(w4_11.getStatus() == Well::Status::SHUT );
|
||||
}
|
||||
}
|
||||
|
@ -3233,9 +3233,8 @@ BOOST_AUTO_TEST_CASE(WELL_STATIC) {
|
||||
BOOST_CHECK(ws.updateRefDepth(1.0));
|
||||
BOOST_CHECK(!ws.updateRefDepth(1.0));
|
||||
|
||||
ws.updateStatus(Well::Status::OPEN, false);
|
||||
BOOST_CHECK(!ws.updateStatus(Well::Status::OPEN, false));
|
||||
BOOST_CHECK(ws.updateStatus(Well::Status::SHUT, false));
|
||||
ws.updateStatus(Well::Status::OPEN, false, false);
|
||||
ws.updateStatus(Well::Status::SHUT, false, false);
|
||||
|
||||
const auto& connections = ws.getConnections();
|
||||
BOOST_CHECK_EQUAL(connections.size(), 0U);
|
||||
@ -4353,3 +4352,111 @@ END
|
||||
}
|
||||
|
||||
|
||||
BOOST_AUTO_TEST_CASE(WELL_STATUS) {
|
||||
const std::string deck_string = R"(
|
||||
START
|
||||
7 OCT 2020 /
|
||||
|
||||
DIMENS
|
||||
10 10 3 /
|
||||
|
||||
GRID
|
||||
DXV
|
||||
10*100.0 /
|
||||
DYV
|
||||
10*100.0 /
|
||||
DZV
|
||||
3*10.0 /
|
||||
|
||||
DEPTHZ
|
||||
121*2000.0 /
|
||||
|
||||
PORO
|
||||
300*0.3 /
|
||||
|
||||
SCHEDULE
|
||||
WELSPECS -- 0
|
||||
'P1' 'G' 10 10 2005 'LIQ' /
|
||||
/
|
||||
|
||||
COMPDAT
|
||||
'P1' 9 9 1 1 'OPEN' 1* 32.948 0.311 3047.839 1* 1* 'X' 22.100 /
|
||||
/
|
||||
|
||||
WCONPROD
|
||||
'P1' 'OPEN' 'ORAT' 123.4 4* 50.0 /
|
||||
/
|
||||
|
||||
|
||||
TSTEP -- 1
|
||||
10 /
|
||||
|
||||
WELPI
|
||||
'P1' 200.0 /
|
||||
/
|
||||
|
||||
TSTEP -- 2
|
||||
10 /
|
||||
|
||||
WELOPEN
|
||||
'P1' SHUT /
|
||||
/
|
||||
|
||||
TSTEP -- 3,4,5
|
||||
10 10 10 /
|
||||
|
||||
WELOPEN
|
||||
'P1' OPEN /
|
||||
/
|
||||
|
||||
TSTEP -- 6,7,8
|
||||
10 10 10/
|
||||
|
||||
END
|
||||
|
||||
END
|
||||
)";
|
||||
|
||||
const auto deck = Parser{}.parseString(deck_string);
|
||||
const auto es = EclipseState{ deck };
|
||||
auto sched = Schedule{ deck, es };
|
||||
{
|
||||
const auto& well = sched.getWell("P1", 0);
|
||||
BOOST_CHECK( well.getStatus() == Well::Status::OPEN);
|
||||
}
|
||||
{
|
||||
const auto& well = sched.getWell("P1", 1);
|
||||
BOOST_CHECK( well.getStatus() == Well::Status::OPEN);
|
||||
}
|
||||
|
||||
{
|
||||
const auto& well = sched.getWell("P1", 2);
|
||||
BOOST_CHECK( well.getStatus() == Well::Status::SHUT);
|
||||
}
|
||||
{
|
||||
const auto& well = sched.getWell("P1", 5);
|
||||
BOOST_CHECK( well.getStatus() == Well::Status::OPEN);
|
||||
}
|
||||
|
||||
sched.shut_well("P1", 0);
|
||||
{
|
||||
const auto& well = sched.getWell("P1", 0);
|
||||
BOOST_CHECK( well.getStatus() == Well::Status::SHUT);
|
||||
}
|
||||
{
|
||||
const auto& well = sched.getWell("P1", 1);
|
||||
BOOST_CHECK( well.getStatus() == Well::Status::SHUT);
|
||||
}
|
||||
{
|
||||
const auto& well = sched.getWell("P1", 2);
|
||||
BOOST_CHECK( well.getStatus() == Well::Status::SHUT);
|
||||
}
|
||||
{
|
||||
const auto& well = sched.getWell("P1", 5);
|
||||
BOOST_CHECK( well.getStatus() == Well::Status::OPEN);
|
||||
}
|
||||
|
||||
//sched.open_well("P1", 2);
|
||||
}
|
||||
|
||||
|
||||
|
@ -79,13 +79,13 @@ BOOST_AUTO_TEST_CASE(WTEST_STATE2) {
|
||||
std::vector<Well> wells;
|
||||
wells.emplace_back("WELL_NAME", "A", 0, 0, 1, 1, 200., WellType(Phase::OIL), Well::ProducerCMode::NONE, Connection::Order::TRACK, us, 0., 1.0, true, true, 0, Well::GasInflowEquation::STD);
|
||||
{
|
||||
wells[0].updateStatus(Well::Status::SHUT, false);
|
||||
wells[0].updateStatus(Well::Status::SHUT, false, false);
|
||||
auto shut_wells = st.updateWells(wc, wells, 5000);
|
||||
BOOST_CHECK_EQUAL(shut_wells.size(), 0U);
|
||||
}
|
||||
|
||||
{
|
||||
wells[0].updateStatus(Well::Status::OPEN, false);
|
||||
wells[0].updateStatus(Well::Status::OPEN, false, false);
|
||||
auto shut_wells = st.updateWells(wc, wells, 5000);
|
||||
BOOST_CHECK_EQUAL( shut_wells.size(), 1U);
|
||||
}
|
||||
@ -116,12 +116,12 @@ BOOST_AUTO_TEST_CASE(WTEST_STATE) {
|
||||
|
||||
WellTestConfig wc;
|
||||
{
|
||||
wells[0].updateStatus(Well::Status::SHUT, false);
|
||||
wells[0].updateStatus(Well::Status::SHUT, false, false);
|
||||
auto shut_wells = st.updateWells(wc, wells, 110. * day);
|
||||
BOOST_CHECK_EQUAL(shut_wells.size(), 0U);
|
||||
}
|
||||
{
|
||||
wells[0].updateStatus(Well::Status::OPEN, false);
|
||||
wells[0].updateStatus(Well::Status::OPEN, false, false);
|
||||
auto shut_wells = st.updateWells(wc, wells, 110. * day);
|
||||
BOOST_CHECK_EQUAL(shut_wells.size(), 0U);
|
||||
}
|
||||
@ -152,10 +152,10 @@ BOOST_AUTO_TEST_CASE(WTEST_STATE) {
|
||||
wc.add_well("WELL_NAME", WellTestConfig::Reason::PHYSICAL, 1000. * day, 3, 0, 5);
|
||||
|
||||
|
||||
wells[0].updateStatus(Well::Status::SHUT, false);
|
||||
wells[0].updateStatus(Well::Status::SHUT, false, false);
|
||||
BOOST_CHECK_EQUAL( st.updateWells(wc, wells, 4100. * day).size(), 0U);
|
||||
|
||||
wells[0].updateStatus(Well::Status::OPEN, false);
|
||||
wells[0].updateStatus(Well::Status::OPEN, false, false);
|
||||
BOOST_CHECK_EQUAL( st.updateWells(wc, wells, 4100. * day).size(), 1U);
|
||||
|
||||
BOOST_CHECK_EQUAL( st.updateWells(wc, wells, 5200. * day).size(), 1U);
|
||||
@ -183,9 +183,9 @@ BOOST_AUTO_TEST_CASE(WTEST_STATE_COMPLETIONS) {
|
||||
const UnitSystem us{};
|
||||
std::vector<Well> wells;
|
||||
wells.emplace_back("WELL_NAME", "A", 0, 0, 1, 1, 200., WellType(Phase::OIL), Well::ProducerCMode::NONE, Connection::Order::TRACK, us, 0., 1.0, true, true, 0, Well::GasInflowEquation::STD);
|
||||
wells[0].updateStatus(Well::Status::OPEN, false);
|
||||
wells[0].updateStatus(Well::Status::OPEN, false, false);
|
||||
wells.emplace_back("WELLX", "A", 0, 0, 2, 2, 200., WellType(Phase::OIL), Well::ProducerCMode::NONE, Connection::Order::TRACK, us, 0., 1.0, true, true, 0, Well::GasInflowEquation::STD);
|
||||
wells[1].updateStatus(Well::Status::OPEN, false);
|
||||
wells[1].updateStatus(Well::Status::OPEN, false, false);
|
||||
|
||||
auto closed_completions = st.updateWells(wc, wells, 5000);
|
||||
BOOST_CHECK_EQUAL( closed_completions.size(), 0U);
|
||||
|
Loading…
Reference in New Issue
Block a user