2021-05-11 16:03:33 -05:00
/*
Copyright 2017 SINTEF Digital , Mathematics and Cybernetics .
Copyright 2017 Statoil ASA .
Copyright 2018 IRIS
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 < http : //www.gnu.org/licenses/>.
*/
# include <config.h>
# include <opm/simulators/wells/WellInterfaceFluidSystem.hpp>
# include <opm/material/fluidsystems/BlackOilFluidSystem.hpp>
2021-12-14 01:30:15 -06:00
# include <opm/input/eclipse/Schedule/Schedule.hpp>
2021-05-11 16:03:33 -05:00
# include <opm/simulators/utils/DeferredLogger.hpp>
2022-10-19 02:55:14 -05:00
# include <opm/simulators/wells/GroupState.hpp>
2021-05-11 16:03:33 -05:00
# include <opm/simulators/wells/ParallelWellInfo.hpp>
2022-10-19 02:55:14 -05:00
# include <opm/simulators/wells/RateConverter.hpp>
2021-09-20 04:16:32 -05:00
# include <opm/simulators/wells/SingleWellState.hpp>
2021-06-10 08:09:05 -05:00
# include <opm/simulators/wells/TargetCalculator.hpp>
2022-10-19 02:55:14 -05:00
# include <opm/simulators/wells/WellGroupHelpers.hpp>
# include <opm/simulators/wells/WellState.hpp>
2021-05-11 16:03:33 -05:00
# include <cassert>
# include <cmath>
namespace Opm
{
template < class FluidSystem >
WellInterfaceFluidSystem < FluidSystem > : :
WellInterfaceFluidSystem ( const Well & well ,
const ParallelWellInfo & parallel_well_info ,
const int time_step ,
const RateConverterType & rate_converter ,
const int pvtRegionIdx ,
const int num_components ,
const int num_phases ,
const int index_of_well ,
const std : : vector < PerforationData > & perf_data )
: WellInterfaceGeneric ( well , parallel_well_info , time_step ,
pvtRegionIdx , num_components , num_phases ,
2021-06-06 01:58:41 -05:00
index_of_well , perf_data )
2021-05-11 16:03:33 -05:00
, rateConverter_ ( rate_converter )
{
}
template < typename FluidSystem >
void
WellInterfaceFluidSystem < FluidSystem > : :
2021-09-20 04:16:32 -05:00
calculateReservoirRates ( SingleWellState & ws ) const
2021-05-11 16:03:33 -05:00
{
const int fipreg = 0 ; // not considering the region for now
const int np = number_of_phases_ ;
std : : vector < double > surface_rates ( np , 0.0 ) ;
for ( int p = 0 ; p < np ; + + p ) {
2021-08-24 04:49:03 -05:00
surface_rates [ p ] = ws . surface_rates [ p ] ;
2021-05-11 16:03:33 -05:00
}
std : : vector < double > voidage_rates ( np , 0.0 ) ;
rateConverter_ . calcReservoirVoidageRates ( fipreg , pvtRegionIdx_ , surface_rates , voidage_rates ) ;
2021-08-06 03:45:13 -05:00
ws . reservoir_rates = voidage_rates ;
2021-05-11 16:03:33 -05:00
}
2021-07-26 10:42:26 -05:00
2021-05-11 16:03:33 -05:00
template < typename FluidSystem >
2021-07-26 10:42:26 -05:00
Well : : ProducerCMode
2021-05-11 16:03:33 -05:00
WellInterfaceFluidSystem < FluidSystem > : :
2021-09-20 04:16:32 -05:00
activeProductionConstraint ( const SingleWellState & ws ,
2022-02-23 04:11:11 -06:00
const SummaryState & summaryState ,
DeferredLogger & deferred_logger ) const
2021-05-11 16:03:33 -05:00
{
2021-07-26 10:42:26 -05:00
const PhaseUsage & pu = this - > phaseUsage ( ) ;
const auto controls = this - > well_ecl_ . productionControls ( summaryState ) ;
2021-08-04 05:03:36 -05:00
const auto currentControl = ws . production_cmode ;
2021-07-26 10:42:26 -05:00
if ( controls . hasControl ( Well : : ProducerCMode : : BHP ) & & currentControl ! = Well : : ProducerCMode : : BHP ) {
const double bhp_limit = controls . bhp_limit ;
2021-08-03 13:05:14 -05:00
double current_bhp = ws . bhp ;
2021-07-26 10:42:26 -05:00
if ( bhp_limit > current_bhp )
return Well : : ProducerCMode : : BHP ;
}
2021-05-11 16:03:33 -05:00
2021-07-26 10:42:26 -05:00
if ( controls . hasControl ( Well : : ProducerCMode : : ORAT ) & & currentControl ! = Well : : ProducerCMode : : ORAT ) {
2021-08-24 04:49:03 -05:00
double current_rate = - ws . surface_rates [ pu . phase_pos [ BlackoilPhases : : Liquid ] ] ;
2021-07-26 10:42:26 -05:00
if ( controls . oil_rate < current_rate )
return Well : : ProducerCMode : : ORAT ;
}
2021-05-11 16:03:33 -05:00
2021-07-26 10:42:26 -05:00
if ( controls . hasControl ( Well : : ProducerCMode : : WRAT ) & & currentControl ! = Well : : ProducerCMode : : WRAT ) {
2021-08-24 04:49:03 -05:00
double current_rate = - ws . surface_rates [ pu . phase_pos [ BlackoilPhases : : Aqua ] ] ;
2021-07-26 10:42:26 -05:00
if ( controls . water_rate < current_rate )
return Well : : ProducerCMode : : WRAT ;
}
2021-05-11 16:03:33 -05:00
2021-07-26 10:42:26 -05:00
if ( controls . hasControl ( Well : : ProducerCMode : : GRAT ) & & currentControl ! = Well : : ProducerCMode : : GRAT ) {
2021-08-24 04:49:03 -05:00
double current_rate = - ws . surface_rates [ pu . phase_pos [ BlackoilPhases : : Vapour ] ] ;
2021-07-26 10:42:26 -05:00
if ( controls . gas_rate < current_rate )
return Well : : ProducerCMode : : GRAT ;
}
2021-05-11 16:03:33 -05:00
2021-07-26 10:42:26 -05:00
if ( controls . hasControl ( Well : : ProducerCMode : : LRAT ) & & currentControl ! = Well : : ProducerCMode : : LRAT ) {
2021-08-24 04:49:03 -05:00
double current_rate = - ws . surface_rates [ pu . phase_pos [ BlackoilPhases : : Liquid ] ] ;
current_rate - = ws . surface_rates [ pu . phase_pos [ BlackoilPhases : : Aqua ] ] ;
2022-09-28 08:18:23 -05:00
bool skip = false ;
if ( controls . liquid_rate = = controls . oil_rate ) {
const double current_water_rate = ws . surface_rates [ pu . phase_pos [ BlackoilPhases : : Aqua ] ] ;
if ( std : : abs ( current_water_rate ) < 1e-12 ) {
skip = true ;
deferred_logger . debug ( " LRAT_ORAT_WELL " , " Well " + this - > name ( ) + " The LRAT target is equal the ORAT target and the water rate is zero, skip checking LRAT " ) ;
}
}
if ( ! skip & & controls . liquid_rate < current_rate )
2021-07-26 10:42:26 -05:00
return Well : : ProducerCMode : : LRAT ;
}
2021-05-11 16:03:33 -05:00
2021-07-26 10:42:26 -05:00
if ( controls . hasControl ( Well : : ProducerCMode : : RESV ) & & currentControl ! = Well : : ProducerCMode : : RESV ) {
double current_rate = 0.0 ;
if ( pu . phase_used [ BlackoilPhases : : Aqua ] )
2021-08-06 03:45:13 -05:00
current_rate - = ws . reservoir_rates [ pu . phase_pos [ BlackoilPhases : : Aqua ] ] ;
2021-05-11 16:03:33 -05:00
2021-07-26 10:42:26 -05:00
if ( pu . phase_used [ BlackoilPhases : : Liquid ] )
2021-08-06 03:45:13 -05:00
current_rate - = ws . reservoir_rates [ pu . phase_pos [ BlackoilPhases : : Liquid ] ] ;
2021-05-11 16:03:33 -05:00
2021-07-26 10:42:26 -05:00
if ( pu . phase_used [ BlackoilPhases : : Vapour ] )
2021-08-06 03:45:13 -05:00
current_rate - = ws . reservoir_rates [ pu . phase_pos [ BlackoilPhases : : Vapour ] ] ;
2021-05-11 16:03:33 -05:00
2021-07-26 10:42:26 -05:00
if ( controls . prediction_mode & & controls . resv_rate < current_rate )
return Well : : ProducerCMode : : RESV ;
2021-05-11 16:03:33 -05:00
2021-07-26 10:42:26 -05:00
if ( ! controls . prediction_mode ) {
const int fipreg = 0 ; // not considering the region for now
const int np = number_of_phases_ ;
std : : vector < double > surface_rates ( np , 0.0 ) ;
if ( pu . phase_used [ BlackoilPhases : : Aqua ] )
surface_rates [ pu . phase_pos [ BlackoilPhases : : Aqua ] ] = controls . water_rate ;
if ( pu . phase_used [ BlackoilPhases : : Liquid ] )
surface_rates [ pu . phase_pos [ BlackoilPhases : : Liquid ] ] = controls . oil_rate ;
if ( pu . phase_used [ BlackoilPhases : : Vapour ] )
surface_rates [ pu . phase_pos [ BlackoilPhases : : Vapour ] ] = controls . gas_rate ;
std : : vector < double > voidage_rates ( np , 0.0 ) ;
rateConverter_ . calcReservoirVoidageRates ( fipreg , pvtRegionIdx_ , surface_rates , voidage_rates ) ;
double resv_rate = 0.0 ;
for ( int p = 0 ; p < np ; + + p )
resv_rate + = voidage_rates [ p ] ;
if ( resv_rate < current_rate )
return Well : : ProducerCMode : : RESV ;
2021-05-11 16:03:33 -05:00
}
2021-07-26 10:42:26 -05:00
}
2021-05-11 16:03:33 -05:00
2021-07-26 10:42:26 -05:00
if ( controls . hasControl ( Well : : ProducerCMode : : THP ) & & currentControl ! = Well : : ProducerCMode : : THP ) {
const auto & thp = getTHPConstraint ( summaryState ) ;
2021-08-03 13:05:14 -05:00
double current_thp = ws . thp ;
2022-04-25 01:09:54 -05:00
if ( thp > current_thp & & ! ws . trivial_target ) {
2022-04-08 01:43:13 -05:00
// If WVFPEXP item 4 is set to YES1 or YES2
// switching to THP is prevented if the well will
// produce at a higher rate with THP control
const auto & wvfpexp = this - > well_ecl_ . getWVFPEXP ( ) ;
2022-02-23 04:11:11 -06:00
bool rate_less_than_potential = true ;
2022-04-08 01:43:13 -05:00
if ( wvfpexp . prevent ( ) ) {
for ( int p = 0 ; p < number_of_phases_ ; + + p ) {
// Currently we use the well potentials here computed before the iterations.
// We may need to recompute the well potentials to get a more
// accurate check here.
rate_less_than_potential = rate_less_than_potential & & ( - ws . surface_rates [ p ] ) < = ws . well_potentials [ p ] ;
}
2022-02-23 04:11:11 -06:00
}
2022-04-08 01:43:13 -05:00
if ( ! wvfpexp . prevent ( ) | | ! rate_less_than_potential ) {
2022-02-23 04:11:11 -06:00
this - > operability_status_ . thp_limit_violated_but_not_switched = false ;
return Well : : ProducerCMode : : THP ;
} else {
this - > operability_status_ . thp_limit_violated_but_not_switched = true ;
2022-09-13 01:58:51 -05:00
deferred_logger . info ( " NOT_SWITCHING_TO_THP " ,
2022-02-23 04:11:11 -06:00
" The THP limit is violated for producer " +
this - > name ( ) +
" . But the rate will increase if switched to THP. " +
" The well is therefore kept at " + Well : : ProducerCMode2String ( currentControl ) ) ;
}
}
2021-05-11 16:03:33 -05:00
}
2021-08-04 05:03:36 -05:00
return currentControl ;
2021-07-26 10:42:26 -05:00
}
template < typename FluidSystem >
Well : : InjectorCMode
WellInterfaceFluidSystem < FluidSystem > : :
2021-09-20 04:16:32 -05:00
activeInjectionConstraint ( const SingleWellState & ws ,
2022-02-23 04:11:11 -06:00
const SummaryState & summaryState ,
DeferredLogger & deferred_logger ) const
2021-07-26 10:42:26 -05:00
{
const PhaseUsage & pu = this - > phaseUsage ( ) ;
const auto controls = this - > well_ecl_ . injectionControls ( summaryState ) ;
2021-08-04 05:03:36 -05:00
const auto currentControl = ws . injection_cmode ;
2021-07-26 10:42:26 -05:00
if ( controls . hasControl ( Well : : InjectorCMode : : BHP ) & & currentControl ! = Well : : InjectorCMode : : BHP )
{
const auto & bhp = controls . bhp_limit ;
2021-08-03 13:05:14 -05:00
double current_bhp = ws . bhp ;
2021-07-26 10:42:26 -05:00
if ( bhp < current_bhp )
return Well : : InjectorCMode : : BHP ;
}
2021-05-11 16:03:33 -05:00
2021-07-26 10:42:26 -05:00
if ( controls . hasControl ( Well : : InjectorCMode : : RATE ) & & currentControl ! = Well : : InjectorCMode : : RATE )
{
InjectorType injectorType = controls . injector_type ;
double current_rate = 0.0 ;
switch ( injectorType ) {
case InjectorType : : WATER :
2021-05-11 16:03:33 -05:00
{
2021-08-24 04:49:03 -05:00
current_rate = ws . surface_rates [ pu . phase_pos [ BlackoilPhases : : Aqua ] ] ;
2021-07-26 10:42:26 -05:00
break ;
2021-05-11 16:03:33 -05:00
}
2021-07-26 10:42:26 -05:00
case InjectorType : : OIL :
{
2021-08-24 04:49:03 -05:00
current_rate = ws . surface_rates [ pu . phase_pos [ BlackoilPhases : : Liquid ] ] ;
2021-07-26 10:42:26 -05:00
break ;
2021-05-11 16:03:33 -05:00
}
2021-07-26 10:42:26 -05:00
case InjectorType : : GAS :
{
2021-08-24 04:49:03 -05:00
current_rate = ws . surface_rates [ pu . phase_pos [ BlackoilPhases : : Vapour ] ] ;
2021-07-26 10:42:26 -05:00
break ;
2021-05-11 16:03:33 -05:00
}
2021-07-26 10:42:26 -05:00
default :
throw ( " Expected WATER, OIL or GAS as type for injectors " + this - > well_ecl_ . name ( ) ) ;
2021-05-11 16:03:33 -05:00
}
2021-07-26 10:42:26 -05:00
if ( controls . surface_rate < current_rate )
return Well : : InjectorCMode : : RATE ;
}
2021-05-11 16:03:33 -05:00
2021-07-26 10:42:26 -05:00
if ( controls . hasControl ( Well : : InjectorCMode : : RESV ) & & currentControl ! = Well : : InjectorCMode : : RESV )
{
double current_rate = 0.0 ;
if ( pu . phase_used [ BlackoilPhases : : Aqua ] )
2021-08-06 03:45:13 -05:00
current_rate + = ws . reservoir_rates [ pu . phase_pos [ BlackoilPhases : : Aqua ] ] ;
2021-05-11 16:03:33 -05:00
2021-07-26 10:42:26 -05:00
if ( pu . phase_used [ BlackoilPhases : : Liquid ] )
2021-08-06 03:45:13 -05:00
current_rate + = ws . reservoir_rates [ pu . phase_pos [ BlackoilPhases : : Liquid ] ] ;
2021-07-26 10:42:26 -05:00
if ( pu . phase_used [ BlackoilPhases : : Vapour ] )
2021-08-06 03:45:13 -05:00
current_rate + = ws . reservoir_rates [ pu . phase_pos [ BlackoilPhases : : Vapour ] ] ;
2021-07-26 10:42:26 -05:00
if ( controls . reservoir_rate < current_rate )
return Well : : InjectorCMode : : RESV ;
}
if ( controls . hasControl ( Well : : InjectorCMode : : THP ) & & currentControl ! = Well : : InjectorCMode : : THP )
{
const auto & thp = getTHPConstraint ( summaryState ) ;
2021-08-03 13:05:14 -05:00
double current_thp = ws . thp ;
2022-02-23 04:11:11 -06:00
if ( thp < current_thp ) {
bool rate_less_than_potential = true ;
for ( int p = 0 ; p < number_of_phases_ ; + + p ) {
// Currently we use the well potentials here computed before the iterations.
// We may need to recompute the well potentials to get a more
// accurate check here.
rate_less_than_potential = rate_less_than_potential & & ( ws . surface_rates [ p ] ) < = ws . well_potentials [ p ] ;
}
if ( ! rate_less_than_potential ) {
this - > operability_status_ . thp_limit_violated_but_not_switched = false ;
return Well : : InjectorCMode : : THP ;
} else {
this - > operability_status_ . thp_limit_violated_but_not_switched = true ;
deferred_logger . debug ( " NOT_SWITCHING_TO_THP " ,
" The THP limit is violated for injector " +
this - > name ( ) +
" . But the rate will increase if switched to THP. " +
" The well is therefore kept at " + Well : : InjectorCMode2String ( currentControl ) ) ;
}
}
2021-07-26 10:42:26 -05:00
}
2021-08-04 05:03:36 -05:00
return currentControl ;
2021-07-26 10:42:26 -05:00
}
template < typename FluidSystem >
bool
WellInterfaceFluidSystem < FluidSystem > : :
2021-09-20 04:16:32 -05:00
checkIndividualConstraints ( SingleWellState & ws ,
2022-02-23 04:11:11 -06:00
const SummaryState & summaryState ,
DeferredLogger & deferred_logger ) const
2021-07-26 10:42:26 -05:00
{
if ( this - > well_ecl_ . isProducer ( ) ) {
2022-02-23 04:11:11 -06:00
auto new_cmode = this - > activeProductionConstraint ( ws , summaryState , deferred_logger ) ;
2021-08-04 05:03:36 -05:00
if ( new_cmode ! = ws . production_cmode ) {
ws . production_cmode = new_cmode ;
2021-07-26 10:42:26 -05:00
return true ;
2021-05-11 16:03:33 -05:00
}
2021-07-26 10:42:26 -05:00
}
2021-05-11 16:03:33 -05:00
2021-07-26 10:42:26 -05:00
if ( this - > well_ecl_ . isInjector ( ) ) {
2022-02-23 04:11:11 -06:00
auto new_cmode = this - > activeInjectionConstraint ( ws , summaryState , deferred_logger ) ;
2021-08-04 05:03:36 -05:00
if ( new_cmode ! = ws . injection_cmode ) {
ws . injection_cmode = new_cmode ;
2021-07-26 10:42:26 -05:00
return true ;
}
2021-05-11 16:03:33 -05:00
}
return false ;
}
template < typename FluidSystem >
std : : pair < bool , double >
WellInterfaceFluidSystem < FluidSystem > : :
checkGroupConstraintsInj ( const Group & group ,
const WellState & well_state ,
const GroupState & group_state ,
const double efficiencyFactor ,
const Schedule & schedule ,
const SummaryState & summaryState ,
DeferredLogger & deferred_logger ) const
{
// Translate injector type from control to Phase.
const auto & well_controls = this - > well_ecl_ . injectionControls ( summaryState ) ;
auto injectorType = well_controls . injector_type ;
Phase injectionPhase ;
switch ( injectorType ) {
case InjectorType : : WATER :
{
injectionPhase = Phase : : WATER ;
break ;
}
case InjectorType : : OIL :
{
injectionPhase = Phase : : OIL ;
break ;
}
case InjectorType : : GAS :
{
injectionPhase = Phase : : GAS ;
break ;
}
default :
throw ( " Expected WATER, OIL or GAS as type for injector " + name ( ) ) ;
}
// Make conversion factors for RESV <-> surface rates.
std : : vector < double > resv_coeff ( phaseUsage ( ) . num_phases , 1.0 ) ;
2021-06-08 08:29:15 -05:00
rateConverter_ . calcInjCoeff ( 0 , pvtRegionIdx_ , resv_coeff ) ; // FIPNUM region 0 here, should use FIPNUM from WELSPECS.
2021-05-11 16:03:33 -05:00
2021-08-24 04:49:03 -05:00
const auto & ws = well_state . well ( this - > index_of_well_ ) ;
2021-05-11 16:03:33 -05:00
// Call check for the well's injection phase.
return WellGroupHelpers : : checkGroupConstraintsInj ( name ( ) ,
well_ecl_ . groupName ( ) ,
group ,
well_state ,
group_state ,
current_step_ ,
guide_rate_ ,
2021-08-24 04:49:03 -05:00
ws . surface_rates . data ( ) ,
2021-05-11 16:03:33 -05:00
injectionPhase ,
phaseUsage ( ) ,
efficiencyFactor ,
schedule ,
summaryState ,
resv_coeff ,
deferred_logger ) ;
}
template < typename FluidSystem >
std : : pair < bool , double >
WellInterfaceFluidSystem < FluidSystem > : :
checkGroupConstraintsProd ( const Group & group ,
const WellState & well_state ,
const GroupState & group_state ,
const double efficiencyFactor ,
const Schedule & schedule ,
const SummaryState & summaryState ,
DeferredLogger & deferred_logger ) const
{
// Make conversion factors for RESV <-> surface rates.
std : : vector < double > resv_coeff ( this - > phaseUsage ( ) . num_phases , 1.0 ) ;
rateConverter_ . calcCoeff ( 0 , pvtRegionIdx_ , resv_coeff ) ; // FIPNUM region 0 here, should use FIPNUM from WELSPECS.
2021-08-24 04:49:03 -05:00
const auto & ws = well_state . well ( this - > index_of_well_ ) ;
2021-05-11 16:03:33 -05:00
return WellGroupHelpers : : checkGroupConstraintsProd ( name ( ) ,
well_ecl_ . groupName ( ) ,
group ,
well_state ,
group_state ,
current_step_ ,
guide_rate_ ,
2021-08-24 04:49:03 -05:00
ws . surface_rates . data ( ) ,
2021-05-11 16:03:33 -05:00
phaseUsage ( ) ,
efficiencyFactor ,
schedule ,
summaryState ,
resv_coeff ,
deferred_logger ) ;
}
template < typename FluidSystem >
bool
WellInterfaceFluidSystem < FluidSystem > : :
checkGroupConstraints ( WellState & well_state ,
const GroupState & group_state ,
const Schedule & schedule ,
const SummaryState & summaryState ,
DeferredLogger & deferred_logger ) const
{
const auto & well = well_ecl_ ;
const int well_index = index_of_well_ ;
2021-08-04 05:03:36 -05:00
auto & ws = well_state . well ( well_index ) ;
2021-05-11 16:03:33 -05:00
if ( well . isInjector ( ) ) {
2021-08-04 05:03:36 -05:00
const auto currentControl = ws . injection_cmode ;
2021-05-11 16:03:33 -05:00
if ( currentControl ! = Well : : InjectorCMode : : GRUP ) {
// This checks only the first encountered group limit,
// in theory there could be several, and then we should
// test all but the one currently applied. At that point,
// this if-statement should be removed and we should always
// check, skipping over only the single group parent whose
// control is the active one for the well (if any).
const auto & group = schedule . getGroup ( well . groupName ( ) , current_step_ ) ;
const double efficiencyFactor = well . getEfficiencyFactor ( ) ;
const std : : pair < bool , double > group_constraint =
checkGroupConstraintsInj ( group , well_state , group_state , efficiencyFactor ,
schedule , summaryState , deferred_logger ) ;
// If a group constraint was broken, we set the current well control to
// be GRUP.
if ( group_constraint . first ) {
2021-08-04 05:03:36 -05:00
ws . injection_cmode = Well : : InjectorCMode : : GRUP ;
2021-05-11 16:03:33 -05:00
const int np = well_state . numPhases ( ) ;
for ( int p = 0 ; p < np ; + + p ) {
2021-08-24 04:49:03 -05:00
ws . surface_rates [ p ] * = group_constraint . second ;
2021-05-11 16:03:33 -05:00
}
}
return group_constraint . first ;
}
}
if ( well . isProducer ( ) ) {
2021-08-04 05:03:36 -05:00
const auto currentControl = ws . production_cmode ;
2021-05-11 16:03:33 -05:00
if ( currentControl ! = Well : : ProducerCMode : : GRUP ) {
// This checks only the first encountered group limit,
// in theory there could be several, and then we should
// test all but the one currently applied. At that point,
// this if-statement should be removed and we should always
// check, skipping over only the single group parent whose
// control is the active one for the well (if any).
const auto & group = schedule . getGroup ( well . groupName ( ) , current_step_ ) ;
const double efficiencyFactor = well . getEfficiencyFactor ( ) ;
const std : : pair < bool , double > group_constraint =
checkGroupConstraintsProd ( group , well_state , group_state , efficiencyFactor ,
schedule , summaryState , deferred_logger ) ;
// If a group constraint was broken, we set the current well control to
// be GRUP.
if ( group_constraint . first ) {
2021-08-04 05:03:36 -05:00
ws . production_cmode = Well : : ProducerCMode : : GRUP ;
2021-05-11 16:03:33 -05:00
const int np = well_state . numPhases ( ) ;
for ( int p = 0 ; p < np ; + + p ) {
2021-08-24 04:49:03 -05:00
ws . surface_rates [ p ] * = group_constraint . second ;
2021-05-11 16:03:33 -05:00
}
}
return group_constraint . first ;
}
}
return false ;
}
template < typename FluidSystem >
bool
WellInterfaceFluidSystem < FluidSystem > : :
checkConstraints ( WellState & well_state ,
const GroupState & group_state ,
const Schedule & schedule ,
const SummaryState & summaryState ,
DeferredLogger & deferred_logger ) const
{
2022-02-23 04:11:11 -06:00
const bool ind_broken = checkIndividualConstraints ( well_state . well ( this - > index_of_well_ ) , summaryState , deferred_logger ) ;
2021-05-11 16:03:33 -05:00
if ( ind_broken ) {
return true ;
} else {
return checkGroupConstraints ( well_state , group_state , schedule , summaryState , deferred_logger ) ;
}
}
2021-05-31 06:21:31 -05:00
template < typename FluidSystem >
int
WellInterfaceFluidSystem < FluidSystem > : :
flowPhaseToEbosPhaseIdx ( const int phaseIdx ) const
{
const auto & pu = this - > phaseUsage ( ) ;
if ( FluidSystem : : phaseIsActive ( FluidSystem : : waterPhaseIdx ) & & pu . phase_pos [ Water ] = = phaseIdx )
return FluidSystem : : waterPhaseIdx ;
if ( FluidSystem : : phaseIsActive ( FluidSystem : : oilPhaseIdx ) & & pu . phase_pos [ Oil ] = = phaseIdx )
return FluidSystem : : oilPhaseIdx ;
if ( FluidSystem : : phaseIsActive ( FluidSystem : : gasPhaseIdx ) & & pu . phase_pos [ Gas ] = = phaseIdx )
return FluidSystem : : gasPhaseIdx ;
// for other phases return the index
return phaseIdx ;
}
2021-06-10 08:09:05 -05:00
template < typename FluidSystem >
std : : optional < double >
WellInterfaceFluidSystem < FluidSystem > : :
getGroupInjectionTargetRate ( const Group & group ,
const WellState & well_state ,
const GroupState & group_state ,
const Schedule & schedule ,
const SummaryState & summaryState ,
const InjectorType & injectorType ,
double efficiencyFactor ,
DeferredLogger & deferred_logger ) const
{
// Setting some defaults to silence warnings below.
// Will be overwritten in the switch statement.
Phase injectionPhase = Phase : : WATER ;
switch ( injectorType ) {
case InjectorType : : WATER :
{
injectionPhase = Phase : : WATER ;
break ;
}
case InjectorType : : OIL :
{
injectionPhase = Phase : : OIL ;
break ;
}
case InjectorType : : GAS :
{
injectionPhase = Phase : : GAS ;
break ;
}
default :
// Should not be here.
assert ( false ) ;
}
auto currentGroupControl = group_state . injection_control ( group . name ( ) , injectionPhase ) ;
if ( currentGroupControl = = Group : : InjectionCMode : : FLD | |
currentGroupControl = = Group : : InjectionCMode : : NONE ) {
if ( ! group . injectionGroupControlAvailable ( injectionPhase ) ) {
// We cannot go any further up the hierarchy. This could
// be the FIELD group, or any group for which this has
// been set in GCONINJE or GCONPROD. If we are here
// anyway, it is likely that the deck set inconsistent
// requirements, such as GRUP control mode on a well with
// no appropriate controls defined on any of its
// containing groups. We will therefore use the wells' bhp
// limit equation as a fallback.
return std : : nullopt ;
} else {
// Inject share of parents control
const auto & parent = schedule . getGroup ( group . parent ( ) , currentStep ( ) ) ;
efficiencyFactor * = group . getGroupEfficiencyFactor ( ) ;
return getGroupInjectionTargetRate ( parent , well_state , group_state , schedule , summaryState , injectorType , efficiencyFactor , deferred_logger ) ;
}
}
const auto pu = phaseUsage ( ) ;
if ( ! group . isInjectionGroup ( ) ) {
return std : : nullopt ;
}
// If we are here, we are at the topmost group to be visited in the recursion.
// This is the group containing the control we will check against.
// Make conversion factors for RESV <-> surface rates.
std : : vector < double > resv_coeff ( pu . num_phases , 1.0 ) ;
rateConverter_ . calcCoeff ( 0 , pvtRegionIdx ( ) , resv_coeff ) ; // FIPNUM region 0 here, should use FIPNUM from WELSPECS.
double sales_target = 0 ;
if ( schedule [ currentStep ( ) ] . gconsale ( ) . has ( group . name ( ) ) ) {
const auto & gconsale = schedule [ currentStep ( ) ] . gconsale ( ) . get ( group . name ( ) , summaryState ) ;
sales_target = gconsale . sales_target ;
}
2021-09-01 08:04:54 -05:00
WellGroupHelpers : : InjectionTargetCalculator tcalc ( currentGroupControl , pu , resv_coeff , group . name ( ) , sales_target , group_state , injectionPhase , group . has_gpmaint_control ( injectionPhase , currentGroupControl ) , deferred_logger ) ;
2021-06-10 08:09:05 -05:00
WellGroupHelpers : : FractionCalculator fcalc ( schedule , well_state , group_state , currentStep ( ) , guideRate ( ) , tcalc . guideTargetMode ( ) , pu , false , injectionPhase ) ;
auto localFraction = [ & ] ( const std : : string & child ) {
2021-06-28 01:59:30 -05:00
return fcalc . localFraction ( child , child ) ; //Note child needs to be passed to always include since the global isGrup map is not updated yet.
2021-06-10 08:09:05 -05:00
} ;
auto localReduction = [ & ] ( const std : : string & group_name ) {
const std : : vector < double > & groupTargetReductions = group_state . injection_reduction_rates ( group_name ) ;
return tcalc . calcModeRateFromRates ( groupTargetReductions ) ;
} ;
const double orig_target = tcalc . groupTarget ( group . injectionControls ( injectionPhase , summaryState ) , deferred_logger ) ;
const auto chain = WellGroupHelpers : : groupChainTopBot ( name ( ) , group . name ( ) , schedule , currentStep ( ) ) ;
// Because 'name' is the last of the elements, and not an ancestor, we subtract one below.
const size_t num_ancestors = chain . size ( ) - 1 ;
double target = orig_target ;
for ( size_t ii = 0 ; ii < num_ancestors ; + + ii ) {
if ( ( ii = = 0 ) | | guideRate ( ) - > has ( chain [ ii ] , injectionPhase ) ) {
// Apply local reductions only at the control level
// (top) and for levels where we have a specified
// group guide rate.
target - = localReduction ( chain [ ii ] ) ;
}
target * = localFraction ( chain [ ii + 1 ] ) ;
}
return std : : max ( 0.0 , target / efficiencyFactor ) ;
}
template < typename FluidSystem >
double
WellInterfaceFluidSystem < FluidSystem > : :
getGroupProductionTargetRate ( const Group & group ,
const WellState & well_state ,
const GroupState & group_state ,
const Schedule & schedule ,
const SummaryState & summaryState ,
double efficiencyFactor ) const
{
const Group : : ProductionCMode & currentGroupControl = group_state . production_control ( group . name ( ) ) ;
if ( currentGroupControl = = Group : : ProductionCMode : : FLD | |
currentGroupControl = = Group : : ProductionCMode : : NONE ) {
if ( ! group . productionGroupControlAvailable ( ) ) {
return 1.0 ;
} else {
// Produce share of parents control
const auto & parent = schedule . getGroup ( group . parent ( ) , currentStep ( ) ) ;
efficiencyFactor * = group . getGroupEfficiencyFactor ( ) ;
return getGroupProductionTargetRate ( parent , well_state , group_state , schedule , summaryState , efficiencyFactor ) ;
}
}
const auto pu = phaseUsage ( ) ;
if ( ! group . isProductionGroup ( ) ) {
return 1.0 ;
}
// If we are here, we are at the topmost group to be visited in the recursion.
// This is the group containing the control we will check against.
// Make conversion factors for RESV <-> surface rates.
std : : vector < double > resv_coeff ( phaseUsage ( ) . num_phases , 1.0 ) ;
rateConverter_ . calcCoeff ( 0 , pvtRegionIdx ( ) , resv_coeff ) ; // FIPNUM region 0 here, should use FIPNUM from WELSPECS.
// gconsale may adjust the grat target.
// the adjusted rates is send to the targetCalculator
double gratTargetFromSales = 0.0 ;
if ( group_state . has_grat_sales_target ( group . name ( ) ) )
gratTargetFromSales = group_state . grat_sales_target ( group . name ( ) ) ;
2021-09-01 08:04:54 -05:00
WellGroupHelpers : : TargetCalculator tcalc ( currentGroupControl , pu , resv_coeff , gratTargetFromSales , group . name ( ) , group_state , group . has_gpmaint_control ( currentGroupControl ) ) ;
2021-06-10 08:09:05 -05:00
WellGroupHelpers : : FractionCalculator fcalc ( schedule , well_state , group_state , currentStep ( ) , guideRate ( ) , tcalc . guideTargetMode ( ) , pu , true , Phase : : OIL ) ;
auto localFraction = [ & ] ( const std : : string & child ) {
2021-06-28 01:59:30 -05:00
return fcalc . localFraction ( child , child ) ; //Note child needs to be passed to always include since the global isGrup map is not updated yet.
2021-06-10 08:09:05 -05:00
} ;
auto localReduction = [ & ] ( const std : : string & group_name ) {
const std : : vector < double > & groupTargetReductions = group_state . production_reduction_rates ( group_name ) ;
return tcalc . calcModeRateFromRates ( groupTargetReductions ) ;
} ;
const double orig_target = tcalc . groupTarget ( group . productionControls ( summaryState ) ) ;
const auto chain = WellGroupHelpers : : groupChainTopBot ( name ( ) , group . name ( ) , schedule , currentStep ( ) ) ;
// Because 'name' is the last of the elements, and not an ancestor, we subtract one below.
const size_t num_ancestors = chain . size ( ) - 1 ;
double target = orig_target ;
for ( size_t ii = 0 ; ii < num_ancestors ; + + ii ) {
if ( ( ii = = 0 ) | | guideRate ( ) - > has ( chain [ ii ] ) ) {
// Apply local reductions only at the control level
// (top) and for levels where we have a specified
// group guide rate.
target - = localReduction ( chain [ ii ] ) ;
}
target * = localFraction ( chain [ ii + 1 ] ) ;
}
// Avoid negative target rates coming from too large local reductions.
const double target_rate = std : : max ( 0.0 , target / efficiencyFactor ) ;
2021-08-24 04:49:03 -05:00
const auto & ws = well_state . well ( this - > index_of_well_ ) ;
const auto & rates = ws . surface_rates ;
2021-06-10 08:09:05 -05:00
const auto current_rate = - tcalc . calcModeRateFromRates ( rates ) ; // Switch sign since 'rates' are negative for producers.
double scale = 1.0 ;
2022-04-25 01:09:54 -05:00
if ( target_rate = = 0.0 ) {
return 0.0 ;
}
2021-06-10 08:09:05 -05:00
if ( current_rate > 1e-14 )
scale = target_rate / current_rate ;
return scale ;
}
2021-05-11 16:03:33 -05:00
template class WellInterfaceFluidSystem < BlackOilFluidSystem < double , BlackOilDefaultIndexTraits > > ;
} // namespace Opm