/*
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
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
namespace {
using P2PCommunicatorType = Dune::Point2PointCommunicator;
using MessageBufferType = P2PCommunicatorType::MessageBufferType;
class PackUnpackXConn : public P2PCommunicatorType::DataHandleInterface
{
public:
using XConn = std::vector;
explicit PackUnpackXConn(const bool isOwner,
const XConn& local,
XConn& global);
// Pack all data associated with link.
void pack(const int link, MessageBufferType& buffer);
// Unpack all data associated with link.
void unpack([[maybe_unused]] const int link,
MessageBufferType& buffer);
private:
const XConn& local_;
XConn& global_;
};
PackUnpackXConn::PackUnpackXConn(const bool isOwner,
const XConn& local,
XConn& global)
: local_ (local)
, global_(global)
{
if (! isOwner) {
return;
}
this->global_.insert(this->global_.end(),
this->local_.begin(),
this->local_.end());
}
void PackUnpackXConn::pack(const int link,
MessageBufferType& buffer)
{
// We should only get one link
if (link != 0) {
throw std::logic_error {
"link in pack() does not match expected value 0"
};
}
// Write all local connection results
{
const auto nconn = this->local_.size();
buffer.write(nconn);
}
for (const auto& conn : this->local_) {
conn.write(buffer);
}
}
void PackUnpackXConn::unpack([[maybe_unused]] const int link,
MessageBufferType& buffer)
{
const auto nconn = [this, &buffer]()
{
auto nc = 0 * this->local_.size();
buffer.read(nc);
return nc;
}();
this->global_.reserve(this->global_.size() + nconn);
for (auto conn = 0*nconn; conn < nconn; ++conn) {
this->global_.emplace_back().read(buffer);
}
}
} // Anonymous namespace
namespace Opm {
WellState::WellState(const ParallelWellInfo& pinfo)
: phase_usage_{}
{
wells_.add("test4",
SingleWellState{"dummy", pinfo, false, 0.0, {}, phase_usage_, 0.0});
}
WellState WellState::serializationTestObject(const ParallelWellInfo& pinfo)
{
WellState result(PhaseUsage{});
result.alq_state = ALQState::serializationTestObject();
result.well_rates = {{"test2", {true, {1.0}}}, {"test3", {false, {2.0}}}};
result.wells_.add("test4", SingleWellState::serializationTestObject(pinfo));
return result;
}
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->wells_.clear();
{
// const int nw = wells->number_of_wells;
const int nw = wells_ecl.size();
// const int np = wells->number_of_phases;
for (int w = 0; w < nw; ++w) {
const Well& well = wells_ecl[w];
// Initialize bhp(), thp(), wellRates(), temperature().
initSingleWell(cellPressures, well, well_perf_data[w], parallel_well_info[w], summary_state);
}
}
}
void WellState::initSingleProducer(const Well& well,
const ParallelWellInfo& well_info,
double pressure_first_connection,
const std::vector& well_perf_data,
const SummaryState& summary_state) {
const auto& pu = this->phase_usage_;
const double temp = 273.15 + 15.56;
auto& ws = this->wells_.add(well.name(),
SingleWellState{well.name(),
well_info,
true,
pressure_first_connection,
well_perf_data,
pu,
temp});
// the rest of the code needs to executed even if ws.perf_data is empty
// as this does not say anything for the whole well if it is distributed.
// Hence never ever return here!
if (well.getStatus() == Well::Status::OPEN) {
ws.status = Well::Status::OPEN;
}
ws.update_producer_targets(well, summary_state);
}
void WellState::initSingleInjector(const Well& well,
const ParallelWellInfo& well_info,
double pressure_first_connection,
const std::vector& well_perf_data,
const SummaryState& summary_state) {
const auto& pu = this->phase_usage_;
const double temp = well.temperature();
auto& ws = this->wells_.add(well.name(), SingleWellState{well.name(),
well_info,
false,
pressure_first_connection,
well_perf_data,
pu,
temp});
// the rest of the code needs to executed even if ws.perf_data is empty
// as this does not say anything for the whole well if it is distributed.
// Hence never ever return here!
if (well.getStatus() == Well::Status::OPEN) {
ws.status = Well::Status::OPEN;
}
ws.update_injector_targets(well, summary_state);
}
void WellState::initSingleWell(const std::vector& cellPressures,
const Well& well,
const std::vector& well_perf_data,
const ParallelWellInfo& well_info,
const SummaryState& summary_state)
{
double pressure_first_connection = -1;
if (!well_perf_data.empty())
pressure_first_connection = cellPressures[well_perf_data[0].cell_index];
pressure_first_connection = well_info.broadcastFirstPerforationValue(pressure_first_connection);
if (well.isInjector()) {
this->initSingleInjector(well, well_info, pressure_first_connection,
well_perf_data, summary_state);
} else {
this->initSingleProducer(well, well_info, pressure_first_connection,
well_perf_data, summary_state);
}
}
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& wname : schedule.wellNames(report_step))
{
well_rates.insert({wname, std::make_pair(false, std::vector(this->numPhases()))});
}
for (const auto& winfo: parallel_well_info)
{
well_rates[winfo.get().name()].first = winfo.get().isOwner();
}
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;
{
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->well(wname).events = wg_events.at(wname);
}
}
for (int w = 0; w < nw; ++w) {
// Initialize perfphaserates_ to well
// rates divided by the number of perforations.
const auto& ecl_well = wells_ecl[w];
auto& ws = this->well(w);
auto& perf_data = ws.perf_data;
const int num_perf_this_well = perf_data.size();
const int global_num_perf_this_well = ecl_well.getConnections().num_open();
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) {
perf_data.phase_rates[this->numPhases()*perf + p] = ws.surface_rates[p] / double(global_num_perf_this_well);
}
}
perf_data.pressure[perf] = cellPressures[well_perf_data[w][perf].cell_index];
}
}
for (int w = 0; w < nw; ++w) {
auto& ws = this->well(w);
if (wells_ecl[w].isProducer()) {
const auto controls = wells_ecl[w].productionControls(summary_state);
if (controls.cmode == Well::ProducerCMode::GRUP && !wells_ecl[w].isAvailableForGroupControl()) {
ws.production_cmode = Well::ProducerCMode::BHP; // wells always has a BHP control
} else {
ws.production_cmode = controls.cmode;
}
}
else {
const auto controls = wells_ecl[w].injectionControls(summary_state);
if (controls.cmode == Well::InjectorCMode::GRUP && !wells_ecl[w].isAvailableForGroupControl()) {
ws.injection_cmode = Well::InjectorCMode::BHP; // wells always has a BHP control
} else {
ws.injection_cmode = controls.cmode;
}
}
}
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->size() > 0) {
for (int w = 0; w < nw; ++w) {
const Well& well = wells_ecl[w];
if (well.getStatus() == Well::Status::SHUT) {
continue;
}
auto& new_well = this->well(w);
const auto& old_index = prevState->index(well.name());
if (old_index.has_value()) {
const auto& prev_well = prevState->well(old_index.value());
new_well.init_timestep(prev_well);
if (prev_well.status == Well::Status::SHUT) {
// Well was shut in previous state, do not use its values.
continue;
}
if (new_well.producer != prev_well.producer)
// Well changed to/from injector from/to producer, do not use its privious values.
continue;
// If new target is set using WCONPROD, WCONINJE etc. we use the new control
if (!new_well.events.hasEvent(WellState::event_mask)) {
new_well.injection_cmode = prev_well.injection_cmode;
new_well.production_cmode = prev_well.production_cmode;
}
new_well.surface_rates = prev_well.surface_rates;
new_well.reservoir_rates = prev_well.reservoir_rates;
new_well.well_potentials = prev_well.well_potentials;
// perfPhaseRates
const int num_perf_old_well = prev_well.perf_data.size();
const int num_perf_this_well = new_well.perf_data.size();
const bool global_num_perf_same = (num_perf_this_well == num_perf_old_well);
// copy perforation rates when the number of
// perforations is equal, otherwise initialize
// perfphaserates to well rates divided by the
// number of perforations.
// TODO: we might still need the values from the prev_well if the connection structure changes
if (global_num_perf_same)
{
auto& perf_data = new_well.perf_data;
const auto& prev_perf_data = prev_well.perf_data;
perf_data.try_assign( prev_perf_data );
} else {
const int global_num_perf_this_well = well.getConnections().num_open();
auto& perf_data = new_well.perf_data;
auto& target_rates = perf_data.phase_rates;
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] = new_well.surface_rates[p] / double(global_num_perf_this_well);
}
}
}
// Productivity index.
new_well.productivity_index = prev_well.productivity_index;
}
// 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) {
new_well.thp = 0;
}
}
}
updateWellsDefaultALQ(wells_ecl);
}
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
{
auto send = std::set{};
auto recv = std::set{};
const auto isOwner = comm.rank() == 0;
if (isOwner) {
for (auto other = 0*comm.size() + 1; other < comm.size(); ++other) {
recv.insert(other);
}
}
else {
send.insert(0);
}
auto toOwnerComm = P2PCommunicatorType{ comm };
toOwnerComm.insertRequest(send, recv);
PackUnpackXConn lineariser { isOwner, from_connections, to_connections };
toOwnerComm.exchange(lineariser);
}
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();
data::Wells res;
for (std::size_t well_index = 0; well_index < this->size(); ++well_index) {
const auto& ws = this->well(well_index);
if ((ws.status == Well::Status::SHUT) && !wasDynamicallyClosed(well_index))
{
continue;
}
const auto& reservoir_rates = ws.reservoir_rates;
const auto& well_potentials = ws.well_potentials;
const auto& wpi = ws.productivity_index;
const auto& wv = ws.surface_rates;
const auto& wname = this->name(well_index);
auto dummyWell = data::Well{};
auto& well = ws.parallel_info.get().isOwner() ? res[wname] : dummyWell;
well.bhp = ws.bhp;
well.thp = ws.thp;
well.temperature = ws.temperature;
if (pu.phase_used[BlackoilPhases::Aqua]) {
well.rates.set(rt::wat, wv[ pu.phase_pos[BlackoilPhases::Aqua] ] );
well.rates.set(rt::reservoir_water, reservoir_rates[pu.phase_pos[BlackoilPhases::Aqua]]);
well.rates.set(rt::productivity_index_water, wpi[pu.phase_pos[BlackoilPhases::Aqua]]);
well.rates.set(rt::well_potential_water, well_potentials[pu.phase_pos[BlackoilPhases::Aqua]]);
}
if (pu.phase_used[BlackoilPhases::Liquid]) {
well.rates.set(rt::oil, wv[ pu.phase_pos[BlackoilPhases::Liquid] ] );
well.rates.set(rt::reservoir_oil, reservoir_rates[pu.phase_pos[BlackoilPhases::Liquid]]);
well.rates.set(rt::productivity_index_oil, wpi[pu.phase_pos[BlackoilPhases::Liquid]]);
well.rates.set(rt::well_potential_oil, well_potentials[pu.phase_pos[BlackoilPhases::Liquid]]);
}
if( pu.phase_used[BlackoilPhases::Vapour] ) {
well.rates.set(rt::gas, wv[ pu.phase_pos[BlackoilPhases::Vapour] ] );
well.rates.set(rt::reservoir_gas, reservoir_rates[pu.phase_pos[BlackoilPhases::Vapour]]);
well.rates.set(rt::productivity_index_gas, wpi[pu.phase_pos[BlackoilPhases::Vapour]]);
well.rates.set(rt::well_potential_gas, well_potentials[pu.phase_pos[BlackoilPhases::Vapour]]);
}
if (pu.has_solvent || pu.has_zFraction) {
well.rates.set(rt::solvent, ws.sum_solvent_rates());
}
if (pu.has_polymer) {
well.rates.set(rt::polymer, ws.sum_polymer_rates());
}
if (pu.has_brine) {
well.rates.set(rt::brine, ws.sum_brine_rates());
}
if (ws.producer) {
well.rates.set(rt::alq, getALQ(wname));
}
else {
well.rates.set(rt::alq, 0.0);
}
well.rates.set(rt::dissolved_gas,
ws.phase_mixing_rates[ws.dissolved_gas] +
ws.phase_mixing_rates[ws.dissolved_gas_in_water]);
well.rates.set(rt::vaporized_oil, ws.phase_mixing_rates[ws.vaporized_oil]);
well.rates.set(rt::vaporized_water, ws.phase_mixing_rates[ws.vaporized_water]);
{
auto& curr = well.current_control;
curr.isProducer = ws.producer;
curr.prod = ws.production_cmode;
curr.inj = ws.injection_cmode;
}
const auto& pwinfo = ws.parallel_info.get();
if (pwinfo.communication().size() == 1) {
reportConnections(well.connections, pu, well_index, globalCellIdxMap);
}
else {
std::vector connections;
reportConnections(connections, pu, well_index, globalCellIdxMap);
gatherVectorsOnRoot(connections, well.connections, pwinfo.communication());
}
const auto nseg = ws.segments.size();
for (auto seg_ix = 0*nseg; seg_ix < nseg; ++seg_ix) {
const auto seg_no = ws.segments.segment_number()[seg_ix];
well.segments[seg_no] = this->reportSegmentResults(well_index, seg_ix, seg_no);
}
}
return res;
}
void WellState::reportConnections(std::vector& connections,
const PhaseUsage &pu,
std::size_t well_index,
const int* globalCellIdxMap) const
{
using rt = data::Rates::opt;
const auto& perf_data = this->well(well_index).perf_data;
const int num_perf_well = perf_data.size();
connections.resize(num_perf_well);
const auto& perf_rates = perf_data.rates;
const auto& perf_pressure = perf_data.pressure;
for (int i = 0; i < num_perf_well; ++i) {
const auto active_index = perf_data.cell_index[i];
auto& connection = connections[ i ];
connection.index = globalCellIdxMap[active_index];
connection.pressure = perf_pressure[i];
connection.reservoir_rate = perf_rates[i];
connection.trans_factor = perf_data.connection_transmissibility_factor[i];
}
const int np = pu.num_phases;
std::vector< rt > phs( np );
std::vector