GasLiftGroupInfo: template Scalar type

This commit is contained in:
Arne Morten Kvarving 2024-02-19 15:46:29 +01:00
parent 37fa8cc5b8
commit 772b00bc19
13 changed files with 360 additions and 303 deletions

View File

@ -114,7 +114,7 @@ namespace Opm {
using GLiftProdWells = typename BlackoilWellModelGeneric<Scalar>::GLiftProdWells;
using GLiftWellStateMap =
typename BlackoilWellModelGeneric<Scalar>::GLiftWellStateMap;
using GLiftEclWells = typename GasLiftGroupInfo::GLiftEclWells;
using GLiftEclWells = typename GasLiftGroupInfo<Scalar>::GLiftEclWells;
using GLiftSyncGroups = typename GasLiftSingleWellGeneric::GLiftSyncGroups;
constexpr static std::size_t pressureVarIndex = GetPropType<TypeTag, Properties::Indices>::pressureSwitchIdx;
typedef typename BaseAuxiliaryModule<TypeTag>::NeighborSet NeighborSet;
@ -518,15 +518,19 @@ namespace Opm {
bool maybeDoGasLiftOptimize(DeferredLogger& deferred_logger);
void gasLiftOptimizationStage1(DeferredLogger& deferred_logger,
GLiftProdWells &prod_wells, GLiftOptWells &glift_wells,
GasLiftGroupInfo &group_info, GLiftWellStateMap &state_map);
GLiftProdWells& prod_wells,
GLiftOptWells& glift_wells,
GasLiftGroupInfo<Scalar>& group_info,
GLiftWellStateMap& state_map);
// cannot be const since it accesses the non-const WellState
void gasLiftOptimizationStage1SingleWell(WellInterface<TypeTag> *well,
DeferredLogger& deferred_logger,
GLiftProdWells &prod_wells, GLiftOptWells &glift_wells,
GasLiftGroupInfo &group_info, GLiftWellStateMap &state_map,
GLiftSyncGroups& groups_to_sync);
void gasLiftOptimizationStage1SingleWell(WellInterface<TypeTag>* well,
DeferredLogger& deferred_logger,
GLiftProdWells& prod_wells,
GLiftOptWells& glift_wells,
GasLiftGroupInfo<Scalar>& group_info,
GLiftWellStateMap& state_map,
GLiftSyncGroups& groups_to_sync);
void extractLegacyCellPvtRegionIndex_();

View File

@ -1397,7 +1397,7 @@ void BlackoilWellModelGeneric<Scalar>::
gasLiftOptimizationStage2(DeferredLogger& deferred_logger,
GLiftProdWells& prod_wells,
GLiftOptWells& glift_wells,
GasLiftGroupInfo& group_info,
GasLiftGroupInfo<Scalar>& group_info,
GLiftWellStateMap& glift_well_state_map,
const int episodeIndex)
{

View File

@ -53,7 +53,7 @@
namespace Opm {
class DeferredLogger;
class EclipseState;
class GasLiftGroupInfo;
template<class Scalar> class GasLiftGroupInfo;
class GasLiftSingleWellGeneric;
class GasLiftWellState;
class Group;
@ -387,7 +387,7 @@ protected:
void gasLiftOptimizationStage2(DeferredLogger& deferred_logger,
GLiftProdWells& prod_wells,
GLiftOptWells& glift_wells,
GasLiftGroupInfo& group_info,
GasLiftGroupInfo<Scalar>& group_info,
GLiftWellStateMap& map,
const int episodeIndex);

View File

@ -1340,8 +1340,10 @@ namespace Opm {
void
BlackoilWellModel<TypeTag>::
gasLiftOptimizationStage1(DeferredLogger& deferred_logger,
GLiftProdWells &prod_wells, GLiftOptWells &glift_wells,
GasLiftGroupInfo &group_info, GLiftWellStateMap &state_map)
GLiftProdWells& prod_wells,
GLiftOptWells &glift_wells,
GasLiftGroupInfo<Scalar>& group_info,
GLiftWellStateMap& state_map)
{
auto comm = simulator_.vanguard().grid().comm();
int num_procs = comm.size();
@ -1445,11 +1447,13 @@ namespace Opm {
template<typename TypeTag>
void
BlackoilWellModel<TypeTag>::
gasLiftOptimizationStage1SingleWell(WellInterface<TypeTag> *well,
DeferredLogger& deferred_logger,
GLiftProdWells &prod_wells, GLiftOptWells &glift_wells,
GasLiftGroupInfo &group_info, GLiftWellStateMap &state_map,
GLiftSyncGroups& sync_groups)
gasLiftOptimizationStage1SingleWell(WellInterface<TypeTag>* well,
DeferredLogger& deferred_logger,
GLiftProdWells& prod_wells,
GLiftOptWells& glift_wells,
GasLiftGroupInfo<Scalar>& group_info,
GLiftWellStateMap& state_map,
GLiftSyncGroups& sync_groups)
{
const auto& summary_state = simulator_.vanguard().summaryState();
std::unique_ptr<GasLiftSingleWell> glift

View File

@ -32,21 +32,20 @@
namespace Opm {
GasLiftGroupInfo::
GasLiftGroupInfo(
GLiftEclWells &ecl_wells,
const Schedule &schedule,
const SummaryState &summary_state,
const int report_step_idx,
const int iteration_idx,
const PhaseUsage &phase_usage,
DeferredLogger &deferred_logger,
WellState<double>& well_state,
const GroupState<double>& group_state,
const Communication &comm,
bool glift_debug
) :
GasLiftCommon(well_state, group_state, deferred_logger, comm, glift_debug)
template<class Scalar>
GasLiftGroupInfo<Scalar>::
GasLiftGroupInfo(GLiftEclWells& ecl_wells,
const Schedule& schedule,
const SummaryState& summary_state,
const int report_step_idx,
const int iteration_idx,
const PhaseUsage& phase_usage,
DeferredLogger& deferred_logger,
WellState<Scalar>& well_state,
const GroupState<Scalar>& group_state,
const Communication& comm,
bool glift_debug)
: GasLiftCommon<Scalar>(well_state, group_state, deferred_logger, comm, glift_debug)
, ecl_wells_{ecl_wells}
, schedule_{schedule}
, summary_state_{summary_state}
@ -54,69 +53,70 @@ GasLiftGroupInfo(
, iteration_idx_{iteration_idx}
, phase_usage_{phase_usage}
, glo_{schedule_.glo(report_step_idx_)}
{
}
{}
/****************************************
* Public methods in alphabetical order
****************************************/
double
GasLiftGroupInfo::
template<class Scalar>
Scalar GasLiftGroupInfo<Scalar>::
alqRate(const std::string& group_name)
{
auto& group_rate = this->group_rate_map_.at(group_name);
return group_rate.alq();
}
int
GasLiftGroupInfo::
template<class Scalar>
int GasLiftGroupInfo<Scalar>::
getGroupIdx(const std::string& group_name)
{
return this->group_idx_.at(group_name);
}
double
GasLiftGroupInfo::
template<class Scalar>
Scalar GasLiftGroupInfo<Scalar>::
gasRate(const std::string& group_name) const
{
auto& group_rate = this->group_rate_map_.at(group_name);
return group_rate.gasRate();
}
double
GasLiftGroupInfo::
template<class Scalar>
Scalar GasLiftGroupInfo<Scalar>::
gasPotential(const std::string& group_name) const
{
auto& group_rate = this->group_rate_map_.at(group_name);
return group_rate.gasPotential();
}
double
GasLiftGroupInfo::
template<class Scalar>
Scalar GasLiftGroupInfo<Scalar>::
waterPotential(const std::string& group_name) const
{
auto& group_rate = this->group_rate_map_.at(group_name);
return group_rate.waterPotential();
}
double
GasLiftGroupInfo::
template<class Scalar>
Scalar GasLiftGroupInfo<Scalar>::
oilPotential(const std::string& group_name) const
{
auto& group_rate = this->group_rate_map_.at(group_name);
return group_rate.oilPotential();
}
std::optional<double>
GasLiftGroupInfo::
template<class Scalar>
std::optional<Scalar>
GasLiftGroupInfo<Scalar>::
gasTarget(const std::string& group_name) const
{
auto& group_rate = this->group_rate_map_.at(group_name);
return group_rate.gasTarget();
}
double
GasLiftGroupInfo::
template<class Scalar>
Scalar GasLiftGroupInfo<Scalar>::
getRate(Rate rate_type, const std::string& group_name) const
{
switch (rate_type) {
@ -133,8 +133,9 @@ getRate(Rate rate_type, const std::string& group_name) const
throw std::runtime_error("This should not happen");
}
}
double
GasLiftGroupInfo::
template<class Scalar>
Scalar GasLiftGroupInfo<Scalar>::
getPotential(Rate rate_type, const std::string& group_name) const
{
switch (rate_type) {
@ -152,8 +153,9 @@ getPotential(Rate rate_type, const std::string& group_name) const
}
}
std::tuple<double, double, double, double>
GasLiftGroupInfo::
template<class Scalar>
std::tuple<Scalar, Scalar, Scalar, Scalar>
GasLiftGroupInfo<Scalar>::
getRates(const int group_idx) const
{
const auto& group_name = groupIdxToName(group_idx);
@ -161,8 +163,9 @@ getRates(const int group_idx) const
return std::make_tuple(rates.oilRate(), rates.gasRate(), rates.waterRate(), rates.alq());
}
std::optional<double>
GasLiftGroupInfo::
template<class Scalar>
std::optional<Scalar>
GasLiftGroupInfo<Scalar>::
getTarget(Rate rate_type, const std::string& group_name) const
{
switch (rate_type) {
@ -180,19 +183,21 @@ getTarget(Rate rate_type, const std::string& group_name) const
}
}
std::vector<std::pair<std::string,double>>&
GasLiftGroupInfo::
template<class Scalar>
std::vector<std::pair<std::string,Scalar>>&
GasLiftGroupInfo<Scalar>::
getWellGroups(const std::string& well_name)
{
assert(this->well_group_map_.count(well_name) == 1);
return this->well_group_map_[well_name];
}
template<class Scalar>
const std::string&
GasLiftGroupInfo::
GasLiftGroupInfo<Scalar>::
groupIdxToName(int group_idx) const
{
const std::string *group_name = nullptr;
const std::string* group_name = nullptr;
// TODO: An alternative to the below loop is to set up a reverse map from idx ->
// string, then we could in theory do faster lookup here..
for (const auto& [key, value] : this->group_idx_) {
@ -209,75 +214,80 @@ groupIdxToName(int group_idx) const
return *group_name;
}
bool
GasLiftGroupInfo::
template<class Scalar>
bool GasLiftGroupInfo<Scalar>::
hasAnyTarget(const std::string& group_name) const
{
return oilTarget(group_name) || gasTarget(group_name)
|| waterTarget(group_name) || liquidTarget(group_name);
}
bool
GasLiftGroupInfo::
template<class Scalar>
bool GasLiftGroupInfo<Scalar>::
hasWell(const std::string& well_name)
{
return this->well_group_map_.count(well_name) == 1;
}
void
GasLiftGroupInfo::
template<class Scalar>
void GasLiftGroupInfo<Scalar>::
initialize()
{
const auto& group = this->schedule_.getGroup("FIELD", this->report_step_idx_);
initializeGroupRatesRecursive_(group);
std::vector<std::string> group_names;
std::vector<double> group_efficiency;
std::vector<Scalar> group_efficiency;
initializeWell2GroupMapRecursive_(
group, group_names, group_efficiency, /*current efficiency=*/1.0);
}
std::optional<double>
GasLiftGroupInfo::
liquidTarget(const std::string &group_name) const
template<class Scalar>
std::optional<Scalar>
GasLiftGroupInfo<Scalar>::
liquidTarget(const std::string& group_name) const
{
auto& group_rate = this->group_rate_map_.at(group_name);
return group_rate.liquidTarget();
}
std::optional<double>
GasLiftGroupInfo::
template<class Scalar>
std::optional<Scalar>
GasLiftGroupInfo<Scalar>::
maxAlq(const std::string& group_name)
{
auto& group_rate = this->group_rate_map_.at(group_name);
return group_rate.maxAlq();
}
std::optional<double>
GasLiftGroupInfo::
template<class Scalar>
std::optional<Scalar>
GasLiftGroupInfo<Scalar>::
maxTotalGasRate(const std::string& group_name)
{
auto& group_rate = this->group_rate_map_.at(group_name);
return group_rate.maxTotalGasRate();
}
double
GasLiftGroupInfo::
oilRate(const std::string &group_name) const
template<class Scalar>
Scalar GasLiftGroupInfo<Scalar>::
oilRate(const std::string& group_name) const
{
auto& group_rate = this->group_rate_map_.at(group_name);
return group_rate.oilRate();
}
std::optional<double>
GasLiftGroupInfo::
oilTarget(const std::string &group_name) const
template<class Scalar>
std::optional<Scalar>
GasLiftGroupInfo<Scalar>::
oilTarget(const std::string& group_name) const
{
auto& group_rate = this->group_rate_map_.at(group_name);
return group_rate.oilTarget();
}
template<class Scalar>
const std::string
GasLiftGroupInfo::
GasLiftGroupInfo<Scalar>::
rateToString(Rate rate) {
switch (rate) {
case Rate::oil:
@ -293,34 +303,42 @@ rateToString(Rate rate) {
}
}
double
GasLiftGroupInfo::
waterRate(const std::string &group_name) const
template<class Scalar>
Scalar GasLiftGroupInfo<Scalar>::
waterRate(const std::string& group_name) const
{
auto& group_rate = this->group_rate_map_.at(group_name);
return group_rate.waterRate();
}
std::optional<double>
GasLiftGroupInfo::
waterTarget(const std::string &group_name) const
template<class Scalar>
std::optional<Scalar>
GasLiftGroupInfo<Scalar>::
waterTarget(const std::string& group_name) const
{
auto& group_rate = this->group_rate_map_.at(group_name);
return group_rate.waterTarget();
}
void
GasLiftGroupInfo::
update(
const std::string &group_name, double delta_oil, double delta_gas, double delta_water, double delta_alq)
template<class Scalar>
void GasLiftGroupInfo<Scalar>::
update(const std::string& group_name,
Scalar delta_oil,
Scalar delta_gas,
Scalar delta_water,
Scalar delta_alq)
{
auto& group_rate = this->group_rate_map_.at(group_name);
group_rate.update(delta_oil, delta_gas, delta_water, delta_alq);
}
void
GasLiftGroupInfo::
updateRate(int idx, double oil_rate, double gas_rate, double water_rate, double alq)
template<class Scalar>
void GasLiftGroupInfo<Scalar>::
updateRate(int idx,
Scalar oil_rate,
Scalar gas_rate,
Scalar water_rate,
Scalar alq)
{
const auto& group_name = groupIdxToName(idx);
auto& rates = this->group_rate_map_.at(group_name);
@ -331,10 +349,9 @@ updateRate(int idx, double oil_rate, double gas_rate, double water_rate, double
* Protected methods in alphabetical order
****************************************/
bool
GasLiftGroupInfo::
checkDoGasLiftOptimization_(const std::string &well_name)
template<class Scalar>
bool GasLiftGroupInfo<Scalar>::
checkDoGasLiftOptimization_(const std::string& well_name)
{
if (this->well_state_.gliftCheckAlqOscillation(well_name)) {
displayDebugMessage_(
@ -390,9 +407,9 @@ checkDoGasLiftOptimization_(const std::string &well_name)
}
}
bool
GasLiftGroupInfo::
checkNewtonIterationIdxOk_(const std::string &well_name)
template<class Scalar>
bool GasLiftGroupInfo<Scalar>::
checkNewtonIterationIdxOk_(const std::string& well_name)
{
if (this->glo_.all_newton()) {
const int nupcol = this->schedule_[this->report_step_idx_].nupcol();
@ -419,17 +436,20 @@ checkNewtonIterationIdxOk_(const std::string &well_name)
}
// This is called by each rank, but the value of "well_name" should be unique
// across ranks
void
GasLiftGroupInfo::
debugDisplayWellContribution_(
const std::string& gr_name, const std::string& well_name,
double eff_factor,
double well_oil_rate, double well_gas_rate, double well_water_rate,
double well_alq,
double oil_rate, double gas_rate, double water_rate,
double alq
) const
// across ranks
template<class Scalar>
void GasLiftGroupInfo<Scalar>::
debugDisplayWellContribution_(const std::string& gr_name,
const std::string& well_name,
Scalar eff_factor,
Scalar well_oil_rate,
Scalar well_gas_rate,
Scalar well_water_rate,
Scalar well_alq,
Scalar oil_rate,
Scalar gas_rate,
Scalar water_rate,
Scalar alq) const
{
const std::string msg = fmt::format("Group rate for {} : Well {} : "
"eff_factor = {}, oil_rate = {}, gas_rate = {}, water_rate = {}, "
@ -439,47 +459,49 @@ debugDisplayWellContribution_(
displayDebugMessage_(msg);
}
void
GasLiftGroupInfo::
debugDisplayUpdatedGroupRates(
const std::string& name,
double oil_rate, double gas_rate, double water_rate, double alq) const
template<class Scalar>
void GasLiftGroupInfo<Scalar>::
debugDisplayUpdatedGroupRates(const std::string& name,
Scalar oil_rate,
Scalar gas_rate,
Scalar water_rate,
Scalar alq) const
{
const std::string msg = fmt::format("Updated group info for {} : "
"oil_rate = {}, gas_rate = {}, water_rate = {}, alq = {}",
name, oil_rate, gas_rate, water_rate, alq);
displayDebugMessageOnRank0_(msg);
this->displayDebugMessageOnRank0_(msg);
}
void
GasLiftGroupInfo::
template<class Scalar>
void GasLiftGroupInfo<Scalar>::
debugEndInitializeGroup(const std::string& name) const
{
const std::string msg = fmt::format("Finished with group {} ...", name);
displayDebugMessageOnRank0_(msg);
this->displayDebugMessageOnRank0_(msg);
}
void
GasLiftGroupInfo::
template<class Scalar>
void GasLiftGroupInfo<Scalar>::
debugStartInitializeGroup(const std::string& name) const
{
const std::string msg = fmt::format("Initializing group {} ...", name);
displayDebugMessageOnRank0_(msg);
this->displayDebugMessageOnRank0_(msg);
}
void
GasLiftGroupInfo::
displayDebugMessage_(const std::string &msg) const
template<class Scalar>
void GasLiftGroupInfo<Scalar>::
displayDebugMessage_(const std::string& msg) const
{
if (this->debug) {
const std::string message = fmt::format("Init group info : {}", msg);
logMessage_(/*prefix=*/"GLIFT", message);
this->logMessage_(/*prefix=*/"GLIFT", message);
}
}
void
GasLiftGroupInfo::
displayDebugMessage_(const std::string &msg, const std::string &well_name)
template<class Scalar>
void GasLiftGroupInfo<Scalar>::
displayDebugMessage_(const std::string& msg, const std::string& well_name)
{
if (this->debug) {
const std::string message = fmt::format("Well {} : {}", well_name, msg);
@ -487,9 +509,9 @@ displayDebugMessage_(const std::string &msg, const std::string &well_name)
}
}
std::tuple<double, double, double, double, double, double>
GasLiftGroupInfo::
template<class Scalar>
std::tuple<Scalar, Scalar, Scalar, Scalar, Scalar, Scalar>
GasLiftGroupInfo<Scalar>::
getProducerWellRates_(const Well* well, int well_index)
{
const auto& pu = this->phase_usage_;
@ -509,21 +531,21 @@ getProducerWellRates_(const Well* well, int well_index)
: 0.0;
const auto controls = well->productionControls(this->summary_state_);
double oil_rate = oil_pot;
Scalar oil_rate = oil_pot;
if (controls.hasControl(Well::ProducerCMode::ORAT)) {
oil_rate = std::min(controls.oil_rate, oil_rate);
oil_rate = std::min(static_cast<Scalar>(controls.oil_rate), oil_rate);
}
double gas_rate = gas_pot;
Scalar gas_rate = gas_pot;
if (controls.hasControl(Well::ProducerCMode::GRAT)) {
gas_rate = std::min(controls.gas_rate, gas_rate);
gas_rate = std::min(static_cast<Scalar>(controls.gas_rate), gas_rate);
}
double water_rate = water_pot;
Scalar water_rate = water_pot;
if (controls.hasControl(Well::ProducerCMode::WRAT)) {
water_rate = std::min(controls.water_rate, water_rate);
water_rate = std::min(static_cast<Scalar>(controls.water_rate), water_rate);
}
if (controls.hasControl(Well::ProducerCMode::LRAT)) {
double liquid_rate = oil_rate + water_rate;
double liquid_rate_lim = std::min(controls.liquid_rate, liquid_rate);
Scalar liquid_rate = oil_rate + water_rate;
Scalar liquid_rate_lim = std::min(static_cast<Scalar>(controls.liquid_rate), liquid_rate);
water_rate = water_rate / liquid_rate * liquid_rate_lim;
oil_rate = oil_rate / liquid_rate * liquid_rate_lim;
}
@ -531,11 +553,12 @@ getProducerWellRates_(const Well* well, int well_index)
return {oil_rate, gas_rate, water_rate, oil_pot, gas_pot, water_pot};
}
std::tuple<double, double, double, double, double, double, double>
GasLiftGroupInfo::
initializeGroupRatesRecursive_(const Group &group)
template<class Scalar>
std::tuple<Scalar, Scalar, Scalar, Scalar, Scalar, Scalar, Scalar>
GasLiftGroupInfo<Scalar>::
initializeGroupRatesRecursive_(const Group& group)
{
std::array<double,7> rates{};
std::array<Scalar,7> rates{};
if (this->debug) debugStartInitializeGroup(group.name());
auto& [oil_rate, water_rate, gas_rate, oil_potential, water_potential, gas_potential, alq] = rates;
if (group.wellgroup()) {
@ -554,7 +577,7 @@ initializeGroupRatesRecursive_(const Group &group)
if (well->isProducer()) {
auto [sw_oil_rate, sw_gas_rate, sw_water_rate, sw_oil_pot, sw_gas_pot, sw_water_pot] = getProducerWellRates_(well, index);
auto sw_alq = this->well_state_.getALQ(well_name);
double factor = well->getEfficiencyFactor();
Scalar factor = well->getEfficiencyFactor();
oil_rate += (factor * sw_oil_rate);
gas_rate += (factor * sw_gas_rate);
water_rate += (factor * sw_water_rate);
@ -595,7 +618,7 @@ initializeGroupRatesRecursive_(const Group &group)
}
}
if (this->debug) debugEndInitializeGroup(group.name());
std::optional<double> oil_target, gas_target, water_target, liquid_target, max_total_gas, max_alq;
std::optional<Scalar> oil_target, gas_target, water_target, liquid_target, max_total_gas, max_alq;
const auto controls = group.productionControls(this->summary_state_);
if (group.has_control(Group::ProductionCMode::LRAT)) {
liquid_target = controls.liquid_target;
@ -610,21 +633,21 @@ initializeGroupRatesRecursive_(const Group &group)
water_target = controls.water_target;
}
if (this->glo_.has_group(group.name())) {
const auto &gl_group = this->glo_.group(group.name());
const auto& gl_group = this->glo_.group(group.name());
max_alq = gl_group.max_lift_gas();
max_total_gas = gl_group.max_total_gas();
}
if (oil_target || liquid_target || water_target || gas_target || max_total_gas || max_alq) {
updateGroupIdxMap_(group.name());
if(oil_target)
if (oil_target)
oil_rate = std::min(oil_rate, *oil_target);
if(gas_target)
if (gas_target)
gas_rate = std::min(gas_rate, *gas_target);
if(water_target)
if (water_target)
water_rate = std::min(water_rate, *water_target);
if(liquid_target) {
double liquid_rate = oil_rate + water_rate;
double liquid_rate_limited = std::min(liquid_rate, *liquid_target);
if (liquid_target) {
Scalar liquid_rate = oil_rate + water_rate;
Scalar liquid_rate_limited = std::min(liquid_rate, *liquid_target);
oil_rate = oil_rate / liquid_rate * liquid_rate_limited;
water_rate = water_rate / liquid_rate * liquid_rate_limited;
}
@ -641,15 +664,14 @@ initializeGroupRatesRecursive_(const Group &group)
return std::make_tuple(oil_rate, gas_rate, water_rate, oil_potential, gas_potential, water_potential, alq);
}
void
GasLiftGroupInfo::
initializeWell2GroupMapRecursive_(
const Group &group,
std::vector<std::string> &group_names,
std::vector<double> &group_efficiency,
double cur_efficiency)
template<class Scalar>
void GasLiftGroupInfo<Scalar>::
initializeWell2GroupMapRecursive_(const Group& group,
std::vector<std::string>& group_names,
std::vector<Scalar>& group_efficiency,
Scalar cur_efficiency)
{
double gfac = group.getGroupEfficiencyFactor();
Scalar gfac = group.getGroupEfficiencyFactor();
cur_efficiency = gfac * cur_efficiency;
for (auto &item : group_efficiency) {
item *= gfac;
@ -668,17 +690,17 @@ initializeWell2GroupMapRecursive_(
if (checkDoGasLiftOptimization_(well_name)) {
const auto &well = this->schedule_.getWell(
well_name, this->report_step_idx_);
double wfac = well.getEfficiencyFactor();
Scalar wfac = well.getEfficiencyFactor();
auto [itr, success] = this->well_group_map_.insert(
{well_name, /*empty vector*/ {}});
assert(success);
auto &vec = itr->second;
auto& vec = itr->second;
assert(group_names.size() == group_efficiency.size());
auto iter2 = group_efficiency.begin();
for (auto iter1 = group_names.begin();
iter1 != group_names.end(); ++iter1)
{
double efficiency = (*iter2) * wfac;
Scalar efficiency = (*iter2) * wfac;
vec.emplace_back(/*group_name=*/*iter1, efficiency);
++iter2;
}
@ -701,15 +723,13 @@ initializeWell2GroupMapRecursive_(
}
}
// TODO: It would be more efficient if the group idx map was build once
// per time step (or better: once per report step) and saved e.g. in
// the well state object, instead of rebuilding here for each of
// NUPCOL well iteration for each time step.
void
GasLiftGroupInfo::
updateGroupIdxMap_(const std::string &group_name)
template<class Scalar>
void GasLiftGroupInfo<Scalar>::
updateGroupIdxMap_(const std::string& group_name)
{
if (this->group_idx_.count(group_name) == 0) {
//auto [itr, success] =
@ -718,4 +738,6 @@ updateGroupIdxMap_(const std::string &group_name)
}
}
template class GasLiftGroupInfo<double>;
} // namespace Opm

View File

@ -29,8 +29,7 @@
#include <tuple>
#include <vector>
namespace Opm
{
namespace Opm {
class DeferredLogger;
class GasLiftOpt;
@ -41,7 +40,8 @@ class SummaryState;
class Well;
template<class Scalar> class WellState;
class GasLiftGroupInfo : public GasLiftCommon<double>
template<class Scalar>
class GasLiftGroupInfo : public GasLiftCommon<Scalar>
{
protected:
class GroupRates;
@ -51,7 +51,7 @@ protected:
// factors of the child groups of the group all the way down
// to the well group.
using Well2GroupMap =
std::map<std::string, std::vector<std::pair<std::string,double>>>;
std::map<std::string, std::vector<std::pair<std::string,Scalar>>>;
using GroupRateMap =
std::map<std::string, GroupRates>;
using GroupIdxMap = std::map<std::string, int>;
@ -62,128 +62,155 @@ protected:
static const int Water = BlackoilPhases::Aqua;
static const int Oil = BlackoilPhases::Liquid;
static const int Gas = BlackoilPhases::Vapour;
public:
enum class Rate {oil, gas, water, liquid};
using GLiftEclWells = std::map<std::string,std::pair<const Well *,int>>;
GasLiftGroupInfo(
GLiftEclWells& ecl_wells,
const Schedule& schedule,
const SummaryState& summary_state,
const int report_step_idx,
const int iteration_idx,
const PhaseUsage& phase_usage,
DeferredLogger& deferred_logger,
WellState<double>& well_state,
const GroupState<double>& group_state,
const Parallel::Communication& comm,
bool glift_debug
);
std::vector<std::pair<std::string,double>>& getWellGroups(
const std::string& well_name);
GasLiftGroupInfo(GLiftEclWells& ecl_wells,
const Schedule& schedule,
const SummaryState& summary_state,
const int report_step_idx,
const int iteration_idx,
const PhaseUsage& phase_usage,
DeferredLogger& deferred_logger,
WellState<Scalar>& well_state,
const GroupState<Scalar>& group_state,
const Parallel::Communication& comm,
bool glift_debug);
double alqRate(const std::string& group_name);
double gasRate(const std::string& group_name) const;
double gasPotential(const std::string& group_name) const;
double waterPotential(const std::string& group_name) const;
double oilPotential(const std::string& group_name) const;
std::vector<std::pair<std::string,Scalar>>&
getWellGroups(const std::string& well_name);
Scalar alqRate(const std::string& group_name);
Scalar gasRate(const std::string& group_name) const;
Scalar gasPotential(const std::string& group_name) const;
Scalar waterPotential(const std::string& group_name) const;
Scalar oilPotential(const std::string& group_name) const;
int getGroupIdx(const std::string& group_name);
double getRate(Rate rate_type, const std::string& group_name) const;
double getPotential(Rate rate_type, const std::string& group_name) const;
std::tuple<double,double,double,double> getRates(const int group_idx) const;
std::optional<double> gasTarget(const std::string& group_name) const;
std::optional<double> getTarget(
Rate rate_type, const std::string& group_name) const;
Scalar getRate(Rate rate_type, const std::string& group_name) const;
Scalar getPotential(Rate rate_type, const std::string& group_name) const;
std::tuple<Scalar,Scalar,Scalar,Scalar> getRates(const int group_idx) const;
std::optional<Scalar> gasTarget(const std::string& group_name) const;
std::optional<Scalar> getTarget(Rate rate_type, const std::string& group_name) const;
const std::string& groupIdxToName(int group_idx) const;
bool hasAnyTarget(const std::string& group_name) const;
bool hasWell(const std::string& well_name);
void initialize();
std::optional<double> liquidTarget(const std::string& group_name) const;
std::optional<double> maxAlq(const std::string& group_name);
std::optional<double> maxTotalGasRate(const std::string& group_name);
double oilRate(const std::string& group_name) const;
std::optional<double> oilTarget(const std::string& group_name) const;
std::optional<Scalar> liquidTarget(const std::string& group_name) const;
std::optional<Scalar> maxAlq(const std::string& group_name);
std::optional<Scalar> maxTotalGasRate(const std::string& group_name);
Scalar oilRate(const std::string& group_name) const;
std::optional<Scalar> oilTarget(const std::string& group_name) const;
static const std::string rateToString(Rate rate);
double waterRate(const std::string& group_name) const;
std::optional<double> waterTarget(const std::string& group_name) const;
Scalar waterRate(const std::string& group_name) const;
std::optional<Scalar> waterTarget(const std::string& group_name) const;
void update(const std::string& well_name,
double delta_oil, double delta_gas, double delta_water, double delta_alq);
void updateRate(int idx, double oil_rate, double gas_rate, double water_rate, double alq);
Scalar delta_oil,
Scalar delta_gas,
Scalar delta_water,
Scalar delta_alq);
void updateRate(int idx,
Scalar oil_rate,
Scalar gas_rate,
Scalar water_rate,
Scalar alq);
const Well2GroupMap& wellGroupMap() { return well_group_map_; }
protected:
bool checkDoGasLiftOptimization_(const std::string& well_name);
bool checkNewtonIterationIdxOk_(const std::string& well_name);
void debugDisplayWellContribution_(
const std::string& gr_name, const std::string& well_name,
double eff_factor,
double well_oil_rate, double well_gas_rate, double well_water_rate,
double well_alq,
double oil_rate, double gas_rate, double water_rate,
double alq
) const;
void debugDisplayWellContribution_(const std::string& gr_name,
const std::string& well_name,
Scalar eff_factor,
Scalar well_oil_rate,
Scalar well_gas_rate,
Scalar well_water_rate,
Scalar well_alq,
Scalar oil_rate,
Scalar gas_rate,
Scalar water_rate,
Scalar alq) const;
void debugDisplayUpdatedGroupRates(const std::string& name,
double oil_rate, double gas_rate, double water_rate, double alq) const;
Scalar oil_rate,
Scalar gas_rate,
Scalar water_rate,
Scalar alq) const;
void debugEndInitializeGroup(const std::string& name) const;
void debugStartInitializeGroup(const std::string& name) const;
void displayDebugMessage_(const std::string& msg) const override;
void displayDebugMessage_(const std::string& msg, const std::string& well_name);
std::tuple<double, double, double, double, double, double>
getProducerWellRates_(const Well* well, const int index);
std::tuple<double, double, double, double, double, double, double>
initializeGroupRatesRecursive_(const Group &group);
void initializeWell2GroupMapRecursive_(
const Group& group, std::vector<std::string>& group_names,
std::vector<double>& group_efficiency, double cur_efficiency);
void updateGroupIdxMap_(const std::string& group_name);
std::tuple<Scalar, Scalar, Scalar, Scalar, Scalar, Scalar>
getProducerWellRates_(const Well* well, const int index);
std::tuple<Scalar, Scalar, Scalar, Scalar, Scalar, Scalar, Scalar>
initializeGroupRatesRecursive_(const Group &group);
void initializeWell2GroupMapRecursive_(const Group& group,
std::vector<std::string>& group_names,
std::vector<Scalar>& group_efficiency,
Scalar cur_efficiency);
void updateGroupIdxMap_(const std::string& group_name);
class GroupRates {
public:
GroupRates( double oil_rate, double gas_rate, double water_rate, double alq,
double oil_potential, double gas_potential, double water_potential,
std::optional<double> oil_target,
std::optional<double> gas_target,
std::optional<double> water_target,
std::optional<double> liquid_target,
std::optional<double> total_gas,
std::optional<double> max_alq
) :
oil_rate_{oil_rate},
gas_rate_{gas_rate},
water_rate_{water_rate},
alq_{alq},
oil_potential_{oil_potential},
gas_potential_{gas_potential},
water_potential_{water_potential},
oil_target_{oil_target},
gas_target_{gas_target},
water_target_{water_target},
liquid_target_{liquid_target},
total_gas_{total_gas},
max_alq_{max_alq}
GroupRates(Scalar oil_rate,
Scalar gas_rate,
Scalar water_rate,
Scalar alq,
Scalar oil_potential,
Scalar gas_potential,
Scalar water_potential,
std::optional<Scalar> oil_target,
std::optional<Scalar> gas_target,
std::optional<Scalar> water_target,
std::optional<Scalar> liquid_target,
std::optional<Scalar> total_gas,
std::optional<Scalar> max_alq)
: oil_rate_{oil_rate}
, gas_rate_{gas_rate}
, water_rate_{water_rate}
, alq_{alq}
, oil_potential_{oil_potential}
, gas_potential_{gas_potential}
, water_potential_{water_potential}
, oil_target_{oil_target}
, gas_target_{gas_target}
, water_target_{water_target}
, liquid_target_{liquid_target}
, total_gas_{total_gas}
, max_alq_{max_alq}
{}
double alq() const { return alq_; }
void assign(double oil_rate, double gas_rate, double water_rate, double alq)
Scalar alq() const { return alq_; }
void assign(Scalar oil_rate,
Scalar gas_rate,
Scalar water_rate,
Scalar alq)
{
oil_rate_ = oil_rate;
gas_rate_ = gas_rate;
water_rate_ = water_rate;
alq_ = alq;
}
double gasRate() const { return gas_rate_; }
double waterRate() const { return water_rate_; }
std::optional<double> gasTarget() const { return gas_target_; }
std::optional<double> waterTarget() const { return water_target_; }
std::optional<double> maxAlq() const { return max_alq_; }
std::optional<double> maxTotalGasRate() const { return total_gas_; }
double oilRate() const { return oil_rate_; }
std::optional<double> oilTarget() const { return oil_target_; }
std::optional<double> liquidTarget() const { return liquid_target_; }
double oilPotential() const { return oil_potential_; }
double gasPotential() const { return gas_potential_; }
double waterPotential() const { return water_potential_; }
Scalar gasRate() const { return gas_rate_; }
Scalar waterRate() const { return water_rate_; }
std::optional<Scalar> gasTarget() const { return gas_target_; }
std::optional<Scalar> waterTarget() const { return water_target_; }
std::optional<Scalar> maxAlq() const { return max_alq_; }
std::optional<Scalar > maxTotalGasRate() const { return total_gas_; }
Scalar oilRate() const { return oil_rate_; }
std::optional<Scalar> oilTarget() const { return oil_target_; }
std::optional<Scalar> liquidTarget() const { return liquid_target_; }
Scalar oilPotential() const { return oil_potential_; }
Scalar gasPotential() const { return gas_potential_; }
Scalar waterPotential() const { return water_potential_; }
void update(double delta_oil, double delta_gas, double delta_water, double delta_alq)
void update(Scalar delta_oil,
Scalar delta_gas,
Scalar delta_water,
Scalar delta_alq)
{
oil_rate_ += delta_oil;
gas_rate_ += delta_gas;
@ -192,28 +219,29 @@ protected:
// Note. We don't updata the potentials at this point. They
// are only needed initially.
}
private:
double oil_rate_;
double gas_rate_;
double water_rate_;
double alq_;
double oil_potential_;
double gas_potential_;
double water_potential_;
std::optional<double> oil_target_;
std::optional<double> gas_target_;
std::optional<double> water_target_;
std::optional<double> liquid_target_;
std::optional<double> total_gas_;
std::optional<double> max_alq_;
Scalar oil_rate_;
Scalar gas_rate_;
Scalar water_rate_;
Scalar alq_;
Scalar oil_potential_;
Scalar gas_potential_;
Scalar water_potential_;
std::optional<Scalar> oil_target_;
std::optional<Scalar> gas_target_;
std::optional<Scalar> water_target_;
std::optional<Scalar> liquid_target_;
std::optional<Scalar> total_gas_;
std::optional<Scalar> max_alq_;
};
GLiftEclWells &ecl_wells_;
const Schedule &schedule_;
const SummaryState &summary_state_;
GLiftEclWells& ecl_wells_;
const Schedule& schedule_;
const SummaryState& summary_state_;
const int report_step_idx_;
const int iteration_idx_;
const PhaseUsage &phase_usage_;
const PhaseUsage& phase_usage_;
const GasLiftOpt& glo_;
GroupRateMap group_rate_map_;
Well2GroupMap well_group_map_;
@ -221,7 +249,6 @@ protected:
int next_group_idx_ = 0;
// Optimize only wells under THP control
bool optimize_only_thp_wells_ = false;
};
} // namespace Opm

View File

@ -48,7 +48,7 @@ namespace Opm
DeferredLogger &deferred_logger,
WellState<Scalar>& well_state,
const GroupState<Scalar>& group_state,
GasLiftGroupInfo &group_info,
GasLiftGroupInfo<Scalar>& group_info,
GLiftSyncGroups &sync_groups,
const Parallel::Communication& comm,
bool glift_debug

View File

@ -42,7 +42,7 @@ GasLiftSingleWellGeneric::GasLiftSingleWellGeneric(DeferredLogger& deferred_logg
const GroupState<double>& group_state,
const Well& ecl_well,
const SummaryState& summary_state,
GasLiftGroupInfo& group_info,
GasLiftGroupInfo<double>& group_info,
const PhaseUsage& phase_usage,
const Schedule& schedule,
const int report_step_idx,
@ -254,7 +254,7 @@ GasLiftSingleWellGeneric::checkGroupTargetsViolated(const BasicRates& rates, con
const std::string msg
= fmt::format("Group {} : {} rate {} exceeds target {}. Stopping iteration",
group_name,
GasLiftGroupInfo::rateToString(rate_type),
GasLiftGroupInfo<double>::rateToString(rate_type),
new_group_rate,
*target_opt);
displayDebugMessage_(msg);
@ -434,7 +434,7 @@ GasLiftSingleWellGeneric::debugShowLimitingTargets_(const LimitedRates& rates) c
if (rates.oil_is_limited) {
const std::string msg = fmt::format("oil rate {} is limited by {} target",
rates.oil,
GasLiftGroupInfo::rateToString(*(rates.oil_limiting_target)));
GasLiftGroupInfo<double>::rateToString(*(rates.oil_limiting_target)));
displayDebugMessage_(msg);
}
if (rates.gas_is_limited) {
@ -444,7 +444,7 @@ GasLiftSingleWellGeneric::debugShowLimitingTargets_(const LimitedRates& rates) c
if (rates.water_is_limited) {
const std::string msg = fmt::format("water rate {} is limited by {} target",
rates.water,
GasLiftGroupInfo::rateToString(*(rates.water_limiting_target)));
GasLiftGroupInfo<double>::rateToString(*(rates.water_limiting_target)));
displayDebugMessage_(msg);
}
} else {
@ -628,7 +628,7 @@ GasLiftSingleWellGeneric::getRateWithLimit_(Rate rate_type, const BasicRates& ra
if (new_rate > target) {
const std::string msg = fmt::format("limiting {} rate to target: "
"computed rate: {}, target: {}",
GasLiftGroupInfo::rateToString(rate_type),
GasLiftGroupInfo<double>::rateToString(rate_type),
new_rate,
target);
displayDebugMessage_(msg);
@ -663,7 +663,7 @@ GasLiftSingleWellGeneric::getRateWithLimit_(Rate rate_type, const BasicRates& ra
target_type = Rate::liquid;
const std::string msg = fmt::format("limiting {} rate to {} due to LRAT target: "
"computed LRAT: {}, target LRAT: {}",
GasLiftGroupInfo::rateToString(rate_type),
GasLiftGroupInfo<double>::rateToString(rate_type),
new_rate,
liq_rate,
liq_target);
@ -790,7 +790,7 @@ GasLiftSingleWellGeneric::getRateWithGroupLimit_(Rate rate_type, const double ne
if (this->debug) {
const std::string msg = fmt::format("limiting {} rate from {} to {} to meet group target {} "
"for group {}. Computed group rate was: {}",
GasLiftGroupInfo::rateToString(rate_type),
GasLiftGroupInfo<double>::rateToString(rate_type),
new_rate,
limited_rate,
gr_target,
@ -1455,7 +1455,7 @@ GasLiftSingleWellGeneric::debugInfoGroupRatesExceedTarget(Rate rate_type,
{
const std::string msg = fmt::format("{} rate for group {} exceeds target: "
"rate = {}, target = {}, the old rate is kept.",
GasLiftGroupInfo::rateToString(rate_type),
GasLiftGroupInfo<double>::rateToString(rate_type),
gr_name,
rate,
target);
@ -1706,13 +1706,13 @@ GasLiftSingleWellGeneric::OptimizeState::checkRatesViolated(const LimitedRates&
std::string target_type;
std::string rate_type;
if (rates.oil_is_limited) {
target_type = GasLiftGroupInfo::rateToString(*(rates.oil_limiting_target));
target_type = GasLiftGroupInfo<double>::rateToString(*(rates.oil_limiting_target));
rate_type = "oil";
} else if (rates.gas_is_limited) {
target_type = "gas";
rate_type = "gas";
} else if (rates.water_is_limited) {
target_type = GasLiftGroupInfo::rateToString(*(rates.water_limiting_target));
target_type = GasLiftGroupInfo<double>::rateToString(*(rates.water_limiting_target));
rate_type = "water";
}
const std::string msg = fmt::format("iteration {} : {} rate was limited due to {} {} target. "

View File

@ -56,7 +56,7 @@ protected:
public:
using GLiftSyncGroups = std::set<int>;
using Rate = GasLiftGroupInfo::Rate;
using Rate = GasLiftGroupInfo<double>::Rate;
struct GradInfo
{
GradInfo() { }
@ -108,7 +108,7 @@ protected:
const GroupState<double>& group_state,
const Well& ecl_well,
const SummaryState& summary_state,
GasLiftGroupInfo& group_info,
GasLiftGroupInfo<double>& group_info,
const PhaseUsage& phase_usage,
const Schedule& schedule,
const int report_step_idx,
@ -332,7 +332,7 @@ protected:
const Well& ecl_well_;
const SummaryState& summary_state_;
GasLiftGroupInfo& group_info_;
GasLiftGroupInfo<double>& group_info_;
const PhaseUsage& phase_usage_;
GLiftSyncGroups& sync_groups_;
const WellProductionControls controls_;

View File

@ -37,7 +37,7 @@ GasLiftSingleWell(const WellInterface<TypeTag> &well,
DeferredLogger &deferred_logger,
WellState<Scalar>& well_state,
const GroupState<Scalar>& group_state,
GasLiftGroupInfo &group_info,
GasLiftGroupInfo<Scalar>& group_info,
GLiftSyncGroups &sync_groups,
const Parallel::Communication& comm,
bool glift_debug

View File

@ -49,7 +49,7 @@ GasLiftStage2::GasLiftStage2(
const GroupState<double>& group_state,
GLiftProdWells &prod_wells,
GLiftOptWells &glift_wells,
GasLiftGroupInfo& group_info,
GasLiftGroupInfo<double>& group_info,
GLiftWellStateMap &state_map,
bool glift_debug
) :

View File

@ -66,7 +66,7 @@ public:
const GroupState<double>& group_state,
GLiftProdWells& prod_wells,
GLiftOptWells& glift_wells,
GasLiftGroupInfo& group_info,
GasLiftGroupInfo<double>& group_info,
GLiftWellStateMap& state_map,
bool glift_debug
);
@ -119,7 +119,7 @@ protected:
GLiftProdWells& prod_wells_;
GLiftOptWells& stage1_wells_;
GasLiftGroupInfo& group_info_;
GasLiftGroupInfo<double>& group_info_;
GLiftWellStateMap& well_state_map_;
int report_step_idx_;

View File

@ -126,7 +126,7 @@ BOOST_AUTO_TEST_CASE(G1)
using WellState = Opm::WellState<double>;
using StdWell = Opm::StandardWell<TypeTag>;
using GasLiftSingleWell = Opm::GasLiftSingleWell<TypeTag>;
using GasLiftGroupInfo = Opm::GasLiftGroupInfo;
using GasLiftGroupInfo = Opm::GasLiftGroupInfo<double>;
using GasLiftSingleWellGeneric = Opm::GasLiftSingleWellGeneric;
using GLiftEclWells = typename GasLiftGroupInfo::GLiftEclWells;
const std::string filename = "GLIFT1.DATA";