/*
Copyright 2014 SINTEF ICT, Applied Mathematics.
Copyright 2017 IRIS AS
This file is part of the Open Porous Media project (OPM).
OPM is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OPM is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with OPM. If not, see .
*/
#include
#include
#include
#include
#include
#include
#include
#include
namespace Opm
{
void WellState::base_init(const std::vector& cellPressures,
const std::vector& wells_ecl,
const std::vector& parallel_well_info,
const std::vector>& well_perf_data,
const SummaryState& summary_state)
{
// clear old name mapping
this->wellMap_.clear();
this->perfpress_.clear();
this->perfrates_.clear();
this->status_.clear();
this->well_perf_data_.clear();
this->parallel_well_info_.clear();
this->wellrates_.clear();
{
// const int nw = wells->number_of_wells;
const int nw = wells_ecl.size();
// const int np = wells->number_of_phases;
bhp_.resize(nw, 0.0);
thp_.resize(nw, 0.0);
temperature_.resize(nw, 273.15 + 15.56); // standard condition temperature
int connpos = 0;
for (int w = 0; w < nw; ++w) {
const Well& well = wells_ecl[w];
// Initialize bhp(), thp(), wellRates(), temperature().
initSingleWell(cellPressures, w, well, well_perf_data[w], parallel_well_info[w], summary_state);
// Setup wellname -> well index mapping.
const int num_perf_this_well = well_perf_data[w].size();
std::string name = well.name();
assert( name.size() > 0 );
mapentry_t& wellMapEntry = wellMap_[name];
wellMapEntry[ 0 ] = w;
wellMapEntry[ 1 ] = connpos;
wellMapEntry[ 2 ] = num_perf_this_well;
connpos += num_perf_this_well;
}
}
}
void WellState::initSingleWell(const std::vector& cellPressures,
const int w,
const Well& well,
const std::vector& well_perf_data,
const ParallelWellInfo* well_info,
const SummaryState& summary_state)
{
assert(well.isInjector() || well.isProducer());
// Set default zero initial well rates.
// May be overwritten below.
const auto& pu = this->phase_usage_;
const int np = pu.num_phases;
if ( well.isInjector() ) {
temperature_[w] = well.injectionControls(summary_state).temperature;
}
this->status_.add(well.name(), Well::Status::OPEN);
this->well_perf_data_.add(well.name(), well_perf_data);
this->parallel_well_info_.add(well.name(), well_info);
this->wellrates_.add(well.name(), std::vector(np, 0));
const int num_perf_this_well = well_info->communication().sum(well_perf_data_[w].size());
this->perfpress_.add(well.name(), std::vector(num_perf_this_well, -1e100));
this->perfrates_.add(well.name(), std::vector(num_perf_this_well, 0));
if ( num_perf_this_well == 0 ) {
// No perforations of the well. Initialize to zero.
bhp_[w] = 0.;
thp_[w] = 0.;
return;
}
const auto inj_controls = well.isInjector() ? well.injectionControls(summary_state) : Well::InjectionControls(0);
const auto prod_controls = well.isProducer() ? well.productionControls(summary_state) : Well::ProductionControls(0);
const bool is_bhp = well.isInjector() ? (inj_controls.cmode == Well::InjectorCMode::BHP)
: (prod_controls.cmode == Well::ProducerCMode::BHP);
const double bhp_limit = well.isInjector() ? inj_controls.bhp_limit : prod_controls.bhp_limit;
const bool is_grup = well.isInjector() ? (inj_controls.cmode == Well::InjectorCMode::GRUP)
: (prod_controls.cmode == Well::ProducerCMode::GRUP);
const double inj_surf_rate = well.isInjector() ? inj_controls.surface_rate : 0.0; // To avoid a "maybe-uninitialized" warning.
const double local_pressure = well_perf_data_[w].empty() ?
0 : cellPressures[well_perf_data_[w][0].cell_index];
const double global_pressure = well_info->broadcastFirstPerforationValue(local_pressure);
if (well.getStatus() == Well::Status::OPEN) {
this->status_[w] = Well::Status::OPEN;
}
if (well.getStatus() == Well::Status::STOP) {
// Stopped well:
// 1. Rates: zero well rates.
// 2. Bhp: assign bhp equal to bhp control, if
// applicable, otherwise assign equal to
// first perforation cell pressure.
if (is_bhp) {
bhp_[w] = bhp_limit;
} else {
bhp_[w] = global_pressure;
}
} else if (is_grup) {
// Well under group control.
// 1. Rates: zero well rates.
// 2. Bhp: initialize bhp to be a
// little above or below (depending on if
// the well is an injector or producer)
// pressure in first perforation cell.
const double safety_factor = well.isInjector() ? 1.01 : 0.99;
bhp_[w] = safety_factor * global_pressure;
} else {
// Open well, under own control:
// 1. Rates: initialize well rates to match
// controls if type is ORAT/GRAT/WRAT
// (producer) or RATE (injector).
// Otherwise, we cannot set the correct
// value here and initialize to zero rate.
auto& rates = this->wellrates_[w];
if (well.isInjector()) {
if (inj_controls.cmode == Well::InjectorCMode::RATE) {
switch (inj_controls.injector_type) {
case InjectorType::WATER:
assert(pu.phase_used[BlackoilPhases::Aqua]);
rates[pu.phase_pos[BlackoilPhases::Aqua]] = inj_surf_rate;
break;
case InjectorType::GAS:
assert(pu.phase_used[BlackoilPhases::Vapour]);
rates[pu.phase_pos[BlackoilPhases::Vapour]] = inj_surf_rate;
break;
case InjectorType::OIL:
assert(pu.phase_used[BlackoilPhases::Liquid]);
rates[pu.phase_pos[BlackoilPhases::Liquid]] = inj_surf_rate;
break;
case InjectorType::MULTI:
// Not currently handled, keep zero init.
break;
}
} else {
// Keep zero init.
}
} else {
assert(well.isProducer());
// Note negative rates for producing wells.
switch (prod_controls.cmode) {
case Well::ProducerCMode::ORAT:
assert(pu.phase_used[BlackoilPhases::Liquid]);
rates[pu.phase_pos[BlackoilPhases::Liquid]] = -prod_controls.oil_rate;
break;
case Well::ProducerCMode::WRAT:
assert(pu.phase_used[BlackoilPhases::Aqua]);
rates[pu.phase_pos[BlackoilPhases::Aqua]] = -prod_controls.water_rate;
break;
case Well::ProducerCMode::GRAT:
assert(pu.phase_used[BlackoilPhases::Vapour]);
rates[pu.phase_pos[BlackoilPhases::Vapour]] = -prod_controls.gas_rate;
break;
default:
// Keep zero init.
break;
}
}
// 2. Bhp: initialize bhp to be target pressure if
// bhp-controlled well, otherwise set to a
// little above or below (depending on if
// the well is an injector or producer)
// pressure in first perforation cell.
if (is_bhp) {
bhp_[w] = bhp_limit;
} else {
const double safety_factor = well.isInjector() ? 1.01 : 0.99;
bhp_[w] = safety_factor * global_pressure;
}
}
// 3. Thp: assign thp equal to thp target/limit, if such a limit exists,
// otherwise keep it zero.
const bool has_thp = well.isInjector() ? inj_controls.hasControl(Well::InjectorCMode::THP)
: prod_controls.hasControl(Well::ProducerCMode::THP);
const double thp_limit = well.isInjector() ? inj_controls.thp_limit : prod_controls.thp_limit;
if (has_thp) {
thp_[w] = thp_limit;
}
}
void WellState::init(const std::vector& cellPressures,
const Schedule& schedule,
const std::vector& wells_ecl,
const std::vector& parallel_well_info,
const int report_step,
const WellState* prevState,
const std::vector>& well_perf_data,
const SummaryState& summary_state)
{
// call init on base class
this->base_init(cellPressures, wells_ecl, parallel_well_info, well_perf_data, summary_state);
this->global_well_info = std::make_optional( schedule, report_step, wells_ecl );
for (const auto& winfo: parallel_well_info)
{
well_rates.insert({winfo->name(), std::make_pair(winfo->isOwner(), std::vector())});
}
const int nw = wells_ecl.size();
if( nw == 0 ) return ;
// Initialize perfphaserates_, which must be done here.
const auto& pu = this->phaseUsage();
const int np = pu.num_phases;
int nperf = 0;
for (const auto& wpd : well_perf_data) {
nperf += wpd.size();
}
well_reservoir_rates_.clear();
well_dissolved_gas_rates_.resize(nw, 0.0);
well_vaporized_oil_rates_.resize(nw, 0.0);
this->events_.clear();
{
const auto& wg_events = schedule[report_step].wellgroup_events();
for (const auto& ecl_well : wells_ecl) {
const auto& wname = ecl_well.name();
if (wg_events.has(wname))
this->events_.add( wname, wg_events.at(wname) );
else
this->events_.add( wname, Events() );
}
}
// Ensure that we start out with zero rates by default.
perfphaserates_.clear();
perfphaserates_.resize(nperf * this->numPhases(), 0.0);
// these are only used to monitor the injectivity
perf_water_throughput_.clear();
perf_water_throughput_.resize(nperf, 0.0);
perf_water_velocity_.clear();
perf_water_velocity_.resize(nperf, 0.0);
perf_skin_pressure_.clear();
perf_skin_pressure_.resize(nperf, 0.0);
num_perf_.resize(nw, 0);
first_perf_index_.resize(nw, 0);
first_perf_index_[0] = 0;
for (int w = 0; w < nw; ++w) {
// Initialize perfphaserates_ to well
// rates divided by the number of perforations.
const auto& wname = wells_ecl[w].name();
const auto& well_info = this->wellMap().at(wname);
const int connpos = well_info[1];
const int num_perf_this_well = well_info[2];
const int global_num_perf_this_well = parallel_well_info[w]->communication().sum(num_perf_this_well);
auto& perf_press = this->perfPress(w);
auto * phase_rates = &this->mutable_perfPhaseRates()[connpos * this->numPhases()];
for (int perf = 0; perf < num_perf_this_well; ++perf) {
if (wells_ecl[w].getStatus() == Well::Status::OPEN) {
for (int p = 0; p < this->numPhases(); ++p) {
phase_rates[this->numPhases()*perf + p] = wellRates(w)[p] / double(global_num_perf_this_well);
}
}
perf_press[perf] = cellPressures[well_perf_data[w][perf].cell_index];
}
num_perf_[w] = num_perf_this_well;
first_perf_index_[w] = connpos;
this->well_reservoir_rates_.add(wname, std::vector(np, 0));
}
is_producer_.clear();
for (int w = 0; w < nw; ++w) {
const auto& ecl_well = wells_ecl[w];
this->is_producer_.add( ecl_well.name(), ecl_well.isProducer());
}
current_injection_controls_.clear();
current_production_controls_.clear();
for (int w = 0; w < nw; ++w) {
const auto& wname = wells_ecl[w].name();
current_production_controls_.add(wname, Well::ProducerCMode::CMODE_UNDEFINED);
current_injection_controls_.add(wname, Well::InjectorCMode::CMODE_UNDEFINED);
if (wells_ecl[w].isProducer()) {
const auto controls = wells_ecl[w].productionControls(summary_state);
currentProductionControl(w, controls.cmode);
}
else {
const auto controls = wells_ecl[w].injectionControls(summary_state);
currentInjectionControl(w, controls.cmode);
}
}
perfRateSolvent_.clear();
perfRateSolvent_.resize(nperf, 0.0);
productivity_index_.resize(nw * this->numPhases(), 0.0);
conn_productivity_index_.resize(nperf * this->numPhases(), 0.0);
well_potentials_.resize(nw * this->numPhases(), 0.0);
perfRatePolymer_.clear();
perfRatePolymer_.resize(nperf, 0.0);
perfRateBrine_.clear();
perfRateBrine_.resize(nperf, 0.0);
for (int w = 0; w < nw; ++w) {
switch (wells_ecl[w].getStatus()) {
case Well::Status::SHUT:
this->shutWell(w);
break;
case Well::Status::STOP:
this->stopWell(w);
break;
default:
this->openWell(w);
break;
}
}
// intialize wells that have been there before
// order may change so the mapping is based on the well name
if (prevState && !prevState->wellMap().empty()) {
auto end = prevState->wellMap().end();
for (int w = 0; w < nw; ++w) {
const Well& well = wells_ecl[w];
if (well.getStatus() == Well::Status::SHUT) {
continue;
}
auto it = prevState->wellMap().find(well.name());
if (it != end)
{
const int newIndex = w;
const int oldIndex = it->second[ 0 ];
if (prevState->status_[oldIndex] == Well::Status::SHUT) {
// Well was shut in previous state, do not use its values.
continue;
}
if (is_producer_[newIndex] != prevState->is_producer_[oldIndex]) {
// Well changed to/from injector from/to producer, do not use its privious values.
continue;
}
// bhp
this->update_bhp( newIndex, prevState->bhp( oldIndex ));
// thp
this->update_thp( newIndex, prevState->thp( oldIndex ));
// If new target is set using WCONPROD, WCONINJE etc. we use the new control
if (!this->events_[w].hasEvent(WellState::event_mask)) {
current_injection_controls_[ newIndex ] = prevState->currentInjectionControl(oldIndex);
current_production_controls_[ newIndex ] = prevState->currentProductionControl(oldIndex);
}
wellRates(w) = prevState->wellRates(oldIndex);
wellReservoirRates(w) = prevState->wellReservoirRates(oldIndex);
// Well potentials
for( int i=0, idx=newIndex*np, oldidx=oldIndex*np; iwellPotentials()[ oldidx ];
}
// perfPhaseRates
const int oldPerf_idx_beg = (*it).second[ 1 ];
const int num_perf_old_well = (*it).second[ 2 ];
const auto new_iter = this->wellMap().find(well.name());
if (new_iter == this->wellMap().end()) {
throw std::logic_error {
well.name() + " is not in internal well map - "
"Bug in WellState"
};
}
const int connpos = new_iter->second[1];
const int num_perf_this_well = new_iter->second[2];
const int num_perf_changed = parallel_well_info[w]->communication()
.sum(static_cast(num_perf_old_well != num_perf_this_well));
const bool global_num_perf_same = num_perf_changed == 0;
// copy perforation rates when the number of
// perforations is equal, otherwise initialize
// perfphaserates to well rates divided by the
// number of perforations.
if (global_num_perf_same)
{
const auto * src_rates = &prevState->perfPhaseRates()[oldPerf_idx_beg* np];
auto * target_rates = &this->mutable_perfPhaseRates()[connpos*np];
for (int perf_index = 0; perf_index < num_perf_this_well; perf_index++) {
for (int p = 0; p < np; p++) {
target_rates[perf_index*np + p] = src_rates[perf_index*np + p];
}
}
} else {
const int global_num_perf_this_well = parallel_well_info[w]->communication().sum(num_perf_this_well);
auto * target_rates = &this->mutable_perfPhaseRates()[connpos*np];
for (int perf_index = 0; perf_index < num_perf_this_well; perf_index++) {
for (int p = 0; p < np; ++p) {
target_rates[perf_index*np + p] = wellRates(w)[p] / double(global_num_perf_this_well);
}
}
}
// perfPressures
if (global_num_perf_same)
{
auto& target_press = perfPress(w);
const auto& src_press = prevState->perfPress(well.name());
for (int perf = 0; perf < num_perf_this_well; ++perf)
{
target_press[perf] = src_press[perf];
}
}
// perfSolventRates
if (pu.has_solvent) {
if (global_num_perf_same)
{
int oldPerf_idx = oldPerf_idx_beg;
for (int perf = connpos; perf < connpos + num_perf_this_well; ++perf, ++oldPerf_idx )
{
perfRateSolvent()[ perf ] = prevState->perfRateSolvent()[ oldPerf_idx ];
}
}
}
// polymer injectivity related
//
// here we did not consider the case that we close
// some perforation during the running and also,
// wells can be shut and re-opened
if (pu.has_polymermw) {
if (global_num_perf_same)
{
int oldPerf_idx = oldPerf_idx_beg;
for (int perf = connpos; perf < connpos + num_perf_this_well; ++perf, ++oldPerf_idx )
{
perf_water_throughput_[ perf ] = prevState->perfThroughput()[ oldPerf_idx ];
perf_skin_pressure_[ perf ] = prevState->perfSkinPressure()[ oldPerf_idx ];
perf_water_velocity_[ perf ] = prevState->perfWaterVelocity()[ oldPerf_idx ];
}
}
}
// Productivity index.
{
auto* thisWellPI = &this ->productivityIndex()[newIndex*np + 0];
const auto* thatWellPI = &prevState->productivityIndex()[oldIndex*np + 0];
for (int p = 0; p < np; ++p) {
thisWellPI[p] = thatWellPI[p];
}
}
}
// If in the new step, there is no THP related
// target/limit anymore, its thp value should be set to
// zero.
const bool has_thp = well.isInjector()
? well.injectionControls (summary_state).hasControl(Well::InjectorCMode::THP)
: well.productionControls(summary_state).hasControl(Well::ProducerCMode::THP);
if (!has_thp) {
this->update_thp(w, 0.0);
}
}
}
{
// we need to create a trival segment related values to avoid there will be some
// multi-segment wells added later.
nseg_ = nw;
top_segment_index_.resize(nw);
seg_number_.resize(nw);
seg_press_.resize(nw);
for (int w = 0; w < nw; ++w) {
top_segment_index_[w] = w;
seg_number_[w] = 1; // Top segment is segment #1
this->seg_press_[w] = this->bhp(w);
}
//seg_rates_ = wellRates();
seg_rates_.assign(nw*np, 0);
seg_pressdrop_.assign(nw, 0.);
seg_pressdrop_hydorstatic_.assign(nw, 0.);
seg_pressdrop_friction_.assign(nw, 0.);
seg_pressdrop_acceleration_.assign(nw, 0.);
}
updateWellsDefaultALQ(wells_ecl);
do_glift_optimization_ = true;
}
void WellState::resize(const std::vector& wells_ecl,
const std::vector& parallel_well_info,
const Schedule& schedule,
const bool handle_ms_well,
const size_t numCells,
const std::vector>& well_perf_data,
const SummaryState& summary_state)
{
const std::vector tmp(numCells, 0.0); // <- UGLY HACK to pass the size
init(tmp, schedule, wells_ecl, parallel_well_info, 0, nullptr, well_perf_data, summary_state);
if (handle_ms_well) {
initWellStateMSWell(wells_ecl, nullptr);
}
}
const std::vector&
WellState::currentWellRates(const std::string& wellName) const
{
auto it = well_rates.find(wellName);
if (it == well_rates.end())
OPM_THROW(std::logic_error, "Could not find any rates for well " << wellName);
return it->second.second;
}
template
void WellState::gatherVectorsOnRoot(const std::vector& from_connections,
std::vector& to_connections,
const Communication& comm) const
{
int size = from_connections.size();
std::vector sizes;
std::vector displ;
if (comm.rank()==0){
sizes.resize(comm.size());
}
comm.gather(&size, sizes.data(), 1, 0);
if (comm.rank()==0){
displ.resize(comm.size()+1, 0);
std::partial_sum(sizes.begin(), sizes.end(), displ.begin()+1);
to_connections.resize(displ.back());
}
comm.gatherv(from_connections.data(), size, to_connections.data(),
sizes.data(), displ.data(), 0);
}
data::Wells
WellState::report(const int* globalCellIdxMap,
const std::function& wasDynamicallyClosed) const
{
if (this->numWells() == 0)
return {};
using rt = data::Rates::opt;
const auto& pu = this->phaseUsage();
const int np = pu.num_phases;
data::Wells res;
for( const auto& itr : this->wellMap() ) {
const auto well_index = itr.second[ 0 ];
if ((this->status_[well_index] == Well::Status::SHUT) &&
! wasDynamicallyClosed(well_index))
{
continue;
}
const auto& pwinfo = *this->parallel_well_info_[well_index];
using WellT = std::remove_reference_t;
WellT dummyWell; // dummy if we are not owner
auto& well = pwinfo.isOwner() ? res[ itr.first ] : dummyWell;
well.bhp = this->bhp(well_index);
well.thp = this->thp( well_index );
well.temperature = this->temperature( well_index );
const auto& wv = this->wellRates(well_index);
if( pu.phase_used[BlackoilPhases::Aqua] ) {
well.rates.set( rt::wat, wv[ pu.phase_pos[BlackoilPhases::Aqua] ] );
}
if( pu.phase_used[BlackoilPhases::Liquid] ) {
well.rates.set( rt::oil, wv[ pu.phase_pos[BlackoilPhases::Liquid] ] );
}
if( pu.phase_used[BlackoilPhases::Vapour] ) {
well.rates.set( rt::gas, wv[ pu.phase_pos[BlackoilPhases::Vapour] ] );
}
if (pwinfo.communication().size()==1)
{
reportConnections(well, pu, itr, globalCellIdxMap);
}
else
{
assert(pwinfo.communication().rank() != 0 || &dummyWell != &well);
// report the local connections
reportConnections(dummyWell, pu, itr, globalCellIdxMap);
// gather them to well on root.
gatherVectorsOnRoot(dummyWell.connections, well.connections,
pwinfo.communication());
}
}
std::vector