2017-08-29 10:26:36 -05:00
/*
Copyright 2017 SINTEF Digital , Mathematics and Cybernetics .
Copyright 2017 Statoil ASA .
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/>.
*/
2024-02-23 01:45:25 -06:00
// Improve IDE experience
# ifndef OPM_MULTISEGMENTWELL_HEADER_INCLUDED
# define OPM_MULTISEGMENTWELL_IMPL_HEADER_INCLUDED
# include <config.h>
# include <opm/simulators/wells/MultisegmentWell.hpp>
# endif
2022-12-13 05:54:27 -06:00
# include <opm/common/Exceptions.hpp>
# include <opm/common/OpmLog/OpmLog.hpp>
2023-06-06 14:31:17 -05:00
# include <opm/input/eclipse/Schedule/MSW/Segment.hpp>
2022-12-13 05:54:27 -06:00
# include <opm/input/eclipse/Schedule/MSW/Valve.hpp>
2023-06-06 14:31:17 -05:00
# include <opm/input/eclipse/Schedule/MSW/WellSegments.hpp>
# include <opm/input/eclipse/Schedule/Well/Connection.hpp>
2023-06-21 04:48:57 -05:00
# include <opm/input/eclipse/Schedule/Well/WellConnections.hpp>
2023-06-06 14:31:17 -05:00
2023-01-16 05:21:29 -06:00
# include <opm/input/eclipse/Units/Units.hpp>
2017-08-29 10:26:36 -05:00
2022-12-20 05:42:09 -06:00
# include <opm/material/densead/EvaluationFormat.hpp>
2022-11-18 04:57:37 -06:00
# include <opm/simulators/wells/MultisegmentWellAssemble.hpp>
2022-10-19 02:55:14 -05:00
# include <opm/simulators/wells/WellBhpThpCalculator.hpp>
2019-05-07 06:06:02 -05:00
# include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>
2017-08-29 10:26:36 -05:00
2020-11-06 14:13:10 -06:00
# include <algorithm>
2023-08-15 02:32:10 -05:00
# include <cstddef>
# include <string>
2020-11-06 14:13:10 -06:00
2024-11-04 05:31:34 -06:00
# if COMPILE_GPU_BRIDGE && (HAVE_CUDA || HAVE_OPENCL)
# include <opm/simulators/linalg/gpubridge/WellContributions.hpp>
2021-06-01 08:49:24 -05:00
# endif
2017-08-29 10:26:36 -05:00
namespace Opm
{
2017-09-12 07:26:32 -05:00
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2017-08-29 10:26:36 -05:00
MultisegmentWell < TypeTag > : :
2020-10-09 08:09:28 -05:00
MultisegmentWell ( const Well & well ,
2024-02-20 08:35:13 -06:00
const ParallelWellInfo < Scalar > & pw_info ,
2020-10-09 08:09:28 -05:00
const int time_step ,
2017-11-23 01:37:30 -06:00
const ModelParameters & param ,
const RateConverterType & rate_converter ,
2017-11-30 09:31:48 -06:00
const int pvtRegionIdx ,
2019-10-23 02:09:45 -05:00
const int num_components ,
const int num_phases ,
const int index_of_well ,
2024-02-20 07:37:48 -06:00
const std : : vector < PerforationData < Scalar > > & perf_data )
2021-06-06 01:58:41 -05:00
: Base ( well , pw_info , time_step , param , rate_converter , pvtRegionIdx , num_components , num_phases , index_of_well , perf_data )
2024-09-30 01:55:00 -05:00
, MSWEval ( static_cast < WellInterfaceIndices < FluidSystem , Indices > & > ( * this ) , pw_info )
2022-04-06 03:09:21 -05:00
, regularize_ ( false )
2024-02-20 15:32:18 -06:00
, segment_fluid_initial_ ( this - > numberOfSegments ( ) , std : : vector < Scalar > ( this - > num_components_ , 0.0 ) )
2017-08-29 10:26:36 -05:00
{
2017-10-13 03:16:44 -05:00
// not handling solvent or polymer for now with multisegment well
2021-05-14 14:38:47 -05:00
if constexpr ( has_solvent ) {
2017-10-13 03:16:44 -05:00
OPM_THROW ( std : : runtime_error , " solvent is not supported by multisegment well yet " ) ;
}
2021-05-14 14:38:47 -05:00
if constexpr ( has_polymer ) {
2017-10-13 03:16:44 -05:00
OPM_THROW ( std : : runtime_error , " polymer is not supported by multisegment well yet " ) ;
}
2018-08-16 04:51:36 -05:00
2021-05-14 14:49:47 -05:00
if constexpr ( Base : : has_energy ) {
2018-08-16 04:51:36 -05:00
OPM_THROW ( std : : runtime_error , " energy is not supported by multisegment well yet " ) ;
2019-04-12 03:43:30 -05:00
}
2020-06-30 06:16:03 -05:00
2021-05-14 14:38:47 -05:00
if constexpr ( Base : : has_foam ) {
2020-06-30 06:16:03 -05:00
OPM_THROW ( std : : runtime_error , " foam is not supported by multisegment well yet " ) ;
}
2021-05-14 14:38:47 -05:00
if constexpr ( Base : : has_brine ) {
2020-06-30 06:16:03 -05:00
OPM_THROW ( std : : runtime_error , " brine is not supported by multisegment well yet " ) ;
}
2022-05-11 04:50:54 -05:00
if constexpr ( Base : : has_watVapor ) {
OPM_THROW ( std : : runtime_error , " water evaporation is not supported by multisegment well yet " ) ;
}
2022-04-21 04:03:11 -05:00
if ( this - > rsRvInj ( ) > 0 ) {
2023-01-02 08:24:35 -06:00
OPM_THROW ( std : : runtime_error ,
" dissolved gas/ vapporized oil in injected oil/gas not supported by multisegment well yet. "
" \n See (WCONINJE item 10 / WCONHIST item 8) " ) ;
2022-04-21 04:03:11 -05:00
}
2022-11-24 07:32:41 -06:00
if constexpr ( ! Indices : : oilEnabled & & Indices : : numPhases > 1 ) {
OPM_THROW ( std : : runtime_error , " water + gas case not supported by multisegment well yet " ) ;
}
2023-05-12 04:57:48 -05:00
this - > thp_update_iterations = true ;
2017-08-29 10:26:36 -05:00
}
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2017-08-29 10:26:36 -05:00
void
MultisegmentWell < TypeTag > : :
init ( const PhaseUsage * phase_usage_arg ,
2024-02-20 15:32:18 -06:00
const std : : vector < Scalar > & depth_arg ,
const Scalar gravity_arg ,
2022-04-12 01:44:52 -05:00
const std : : vector < Scalar > & B_avg ,
const bool changed_to_open_this_step )
2017-08-29 10:26:36 -05:00
{
2024-10-10 11:53:47 -05:00
Base : : init ( phase_usage_arg , depth_arg , gravity_arg , B_avg , changed_to_open_this_step ) ;
2017-08-29 10:26:36 -05:00
// TODO: for StandardWell, we need to update the perf depth here using depth_arg.
// for MultisegmentWell, it is much more complicated.
// It can be specified directly, it can be calculated from the segment depth,
// it can also use the cell center, which is the same for StandardWell.
// For the last case, should we update the depth with the depth_arg? For the
// future, it can be a source of wrong result with Multisegment well.
2017-08-31 11:22:48 -05:00
// An indicator from the opm-parser should indicate what kind of depth we should use here.
2017-09-01 10:37:51 -05:00
// \Note: we do not update the depth here. And it looks like for now, we only have the option to use
// specified perforation depth
2024-10-10 11:53:47 -05:00
this - > initMatrixAndVectors ( ) ;
2017-09-12 07:26:32 -05:00
2022-11-05 16:51:59 -05:00
// calculate the depth difference between the perforations and the perforated grid block
2024-10-15 07:11:34 -05:00
for ( int local_perf_index = 0 ; local_perf_index < this - > number_of_perforations_ ; + + local_perf_index ) {
// This variable loops over the number_of_perforations_ of *this* process, hence it is *local*.
const int cell_idx = this - > well_cells_ [ local_perf_index ] ;
this - > cell_perforation_depth_diffs_ [ local_perf_index ] = depth_arg [ cell_idx ] - this - > perf_depth_ [ local_perf_index ] ;
2017-09-12 07:26:32 -05:00
}
2017-08-29 10:26:36 -05:00
}
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2017-08-29 10:26:36 -05:00
void
MultisegmentWell < TypeTag > : :
2022-11-09 06:01:30 -06:00
initPrimaryVariablesEvaluation ( )
2017-08-29 10:26:36 -05:00
{
2022-11-07 23:38:12 -06:00
this - > primary_variables_ . init ( ) ;
2017-08-29 10:26:36 -05:00
}
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2017-08-29 10:26:36 -05:00
void
MultisegmentWell < TypeTag > : :
2024-03-19 06:50:34 -05:00
updatePrimaryVariables ( const Simulator & simulator ,
2024-02-17 11:13:46 -06:00
const WellState < Scalar > & well_state ,
2024-03-19 06:50:34 -05:00
DeferredLogger & deferred_logger )
2017-08-29 10:26:36 -05:00
{
2024-03-19 06:50:34 -05:00
const bool stop_or_zero_rate_target = this - > stoppedOrZeroRateTarget ( simulator , well_state , deferred_logger ) ;
2023-03-29 07:33:05 -05:00
this - > primary_variables_ . update ( well_state , stop_or_zero_rate_target ) ;
2017-08-29 10:26:36 -05:00
}
2021-06-03 08:34:14 -05:00
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2017-08-29 10:26:36 -05:00
void
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
updateWellStateWithTarget ( const Simulator & simulator ,
2024-02-17 11:13:46 -06:00
const GroupState < Scalar > & group_state ,
2024-02-17 11:13:46 -06:00
WellState < Scalar > & well_state ,
2021-05-05 04:22:44 -05:00
DeferredLogger & deferred_logger ) const
2017-08-29 10:26:36 -05:00
{
2024-02-06 04:55:07 -06:00
Base : : updateWellStateWithTarget ( simulator , group_state , well_state , deferred_logger ) ;
2020-12-14 08:55:14 -06:00
// scale segment rates based on the wellRates
// and segment pressure based on bhp
2022-12-19 08:52:36 -06:00
this - > scaleSegmentRatesWithWellRates ( this - > segments_ . inlets ( ) ,
this - > segments_ . perforations ( ) ,
2022-12-19 06:52:21 -06:00
well_state ) ;
2021-06-03 08:34:14 -05:00
this - > scaleSegmentPressuresWithBhp ( well_state ) ;
2017-08-29 10:26:36 -05:00
}
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2018-10-25 04:57:47 -05:00
ConvergenceReport
2017-08-29 10:26:36 -05:00
MultisegmentWell < TypeTag > : :
2024-03-19 06:50:34 -05:00
getWellConvergence ( const Simulator & /* simulator */ ,
2024-02-17 11:13:46 -06:00
const WellState < Scalar > & well_state ,
2024-02-20 15:32:18 -06:00
const std : : vector < Scalar > & B_avg ,
2021-06-03 08:34:14 -05:00
DeferredLogger & deferred_logger ,
const bool relax_tolerance ) const
2017-08-29 10:26:36 -05:00
{
2021-06-03 08:34:14 -05:00
return this - > MSWEval : : getWellConvergence ( well_state ,
B_avg ,
deferred_logger ,
2021-09-06 05:58:16 -05:00
this - > param_ . max_residual_allowed_ ,
this - > param_ . tolerance_wells_ ,
2021-09-06 03:08:28 -05:00
this - > param_ . relaxed_tolerance_flow_well_ ,
2021-09-06 05:58:16 -05:00
this - > param_ . tolerance_pressure_ms_wells_ ,
2021-09-06 03:08:28 -05:00
this - > param_ . relaxed_tolerance_pressure_ms_well_ ,
2023-11-08 14:07:47 -06:00
relax_tolerance ,
this - > wellIsStopped ( ) ) ;
2023-03-30 08:49:27 -05:00
2017-08-29 10:26:36 -05:00
}
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2017-08-29 10:26:36 -05:00
void
MultisegmentWell < TypeTag > : :
2017-09-01 10:37:51 -05:00
apply ( const BVector & x , BVector & Ax ) const
2017-08-29 10:26:36 -05:00
{
2022-11-11 14:41:24 -06:00
if ( ! this - > isOperableAndSolvable ( ) & & ! this - > wellIsStopped ( ) ) {
return ;
}
2020-11-27 00:57:55 -06:00
2022-11-11 14:41:24 -06:00
if ( this - > param_ . matrix_add_well_contributions_ ) {
2020-10-29 06:57:29 -05:00
// Contributions are already in the matrix itself
return ;
}
2017-09-13 10:48:08 -05:00
2022-11-11 14:41:24 -06:00
this - > linSys_ . apply ( x , Ax ) ;
2017-08-29 10:26:36 -05:00
}
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2017-08-29 10:26:36 -05:00
void
MultisegmentWell < TypeTag > : :
apply ( BVector & r ) const
{
2022-11-11 14:41:24 -06:00
if ( ! this - > isOperableAndSolvable ( ) & & ! this - > wellIsStopped ( ) ) {
return ;
}
2020-11-27 00:57:55 -06:00
2022-11-11 14:41:24 -06:00
this - > linSys_ . apply ( r ) ;
2017-08-29 10:26:36 -05:00
}
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2017-08-29 10:26:36 -05:00
void
MultisegmentWell < TypeTag > : :
2024-03-19 06:50:34 -05:00
recoverWellSolutionAndUpdateWellState ( const Simulator & simulator ,
2023-03-22 09:10:00 -05:00
const BVector & x ,
2024-02-17 11:13:46 -06:00
WellState < Scalar > & well_state ,
2022-11-09 06:01:30 -06:00
DeferredLogger & deferred_logger )
2017-08-29 10:26:36 -05:00
{
2022-11-11 14:41:24 -06:00
if ( ! this - > isOperableAndSolvable ( ) & & ! this - > wellIsStopped ( ) ) {
return ;
}
2020-11-27 00:57:55 -06:00
2017-09-18 09:31:11 -05:00
BVectorWell xw ( 1 ) ;
2022-11-11 14:41:24 -06:00
this - > linSys_ . recoverSolutionWell ( x , xw ) ;
2024-03-19 06:50:34 -05:00
updateWellState ( simulator , xw , well_state , deferred_logger ) ;
2017-08-29 10:26:36 -05:00
}
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2017-08-29 10:26:36 -05:00
void
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
computeWellPotentials ( const Simulator & simulator ,
2024-02-17 11:13:46 -06:00
const WellState < Scalar > & well_state ,
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > & well_potentials ,
2021-05-05 04:22:44 -05:00
DeferredLogger & deferred_logger )
2017-08-29 10:26:36 -05:00
{
2023-05-12 05:53:59 -05:00
const auto [ compute_potential , bhp_controlled_well ] =
2024-02-19 07:34:38 -06:00
this - > WellInterfaceGeneric < Scalar > : : computeWellPotentials ( well_potentials , well_state ) ;
2019-11-18 09:33:10 -06:00
2023-05-12 05:53:59 -05:00
if ( ! compute_potential ) {
2019-11-18 09:33:10 -06:00
return ;
}
2021-05-25 08:24:00 -05:00
2020-12-22 02:25:56 -06:00
debug_cost_counter_ = 0 ;
2023-10-26 10:28:05 -05:00
bool converged_implicit = false ;
if ( this - > param_ . local_well_solver_control_switching_ ) {
2024-02-06 04:55:07 -06:00
converged_implicit = computeWellPotentialsImplicit ( simulator , well_potentials , deferred_logger ) ;
2023-10-26 10:28:05 -05:00
if ( ! converged_implicit ) {
deferred_logger . debug ( " Implicit potential calculations failed for well "
+ this - > name ( ) + " , reverting to original aproach. " ) ;
}
}
if ( ! converged_implicit ) {
// does the well have a THP related constraint?
2024-02-06 04:55:07 -06:00
const auto & summaryState = simulator . vanguard ( ) . summaryState ( ) ;
2023-10-26 10:28:05 -05:00
if ( ! Base : : wellHasTHPConstraints ( summaryState ) | | bhp_controlled_well ) {
2024-02-06 04:55:07 -06:00
computeWellRatesAtBhpLimit ( simulator , well_potentials , deferred_logger ) ;
2023-10-26 10:28:05 -05:00
} else {
well_potentials = computeWellPotentialWithTHP (
2024-02-06 04:55:07 -06:00
well_state , simulator , deferred_logger ) ;
2023-10-26 10:28:05 -05:00
}
2019-11-18 09:33:10 -06:00
}
deferred_logger . debug ( " Cost in iterations of finding well potential for well "
2021-09-06 05:58:16 -05:00
+ this - > name ( ) + " : " + std : : to_string ( debug_cost_counter_ ) ) ;
2021-11-18 05:57:16 -06:00
2023-05-12 05:53:59 -05:00
this - > checkNegativeWellPotentials ( well_potentials ,
this - > param_ . check_well_operability_ ,
deferred_logger ) ;
2019-11-18 09:33:10 -06:00
}
2019-04-23 06:30:12 -05:00
2019-11-18 09:33:10 -06:00
template < typename TypeTag >
void
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
computeWellRatesAtBhpLimit ( const Simulator & simulator ,
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > & well_flux ,
2021-05-05 04:22:44 -05:00
DeferredLogger & deferred_logger ) const
2019-11-18 09:33:10 -06:00
{
2021-09-06 05:58:16 -05:00
if ( this - > well_ecl_ . isInjector ( ) ) {
2024-02-06 04:55:07 -06:00
const auto controls = this - > well_ecl_ . injectionControls ( simulator . vanguard ( ) . summaryState ( ) ) ;
computeWellRatesWithBhpIterations ( simulator , controls . bhp_limit , well_flux , deferred_logger ) ;
2019-11-18 09:33:10 -06:00
} else {
2024-02-06 04:55:07 -06:00
const auto controls = this - > well_ecl_ . productionControls ( simulator . vanguard ( ) . summaryState ( ) ) ;
computeWellRatesWithBhpIterations ( simulator , controls . bhp_limit , well_flux , deferred_logger ) ;
2019-04-23 06:30:12 -05:00
}
2017-08-29 10:26:36 -05:00
}
2019-04-23 06:30:12 -05:00
template < typename TypeTag >
void
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
computeWellRatesWithBhp ( const Simulator & simulator ,
2024-02-20 15:32:18 -06:00
const Scalar & bhp ,
std : : vector < Scalar > & well_flux ,
2021-05-05 04:22:44 -05:00
DeferredLogger & deferred_logger ) const
2019-04-23 06:30:12 -05:00
{
2021-09-13 03:58:50 -05:00
const int np = this - > number_of_phases_ ;
well_flux . resize ( np , 0.0 ) ;
const bool allow_cf = this - > getAllowCrossFlow ( ) ;
const int nseg = this - > numberOfSegments ( ) ;
2024-02-17 11:13:46 -06:00
const WellState < Scalar > & well_state = simulator . problem ( ) . wellModel ( ) . wellState ( ) ;
2021-09-13 03:58:50 -05:00
const auto & ws = well_state . well ( this - > indexOfWell ( ) ) ;
auto segments_copy = ws . segments ;
segments_copy . scale_pressure ( bhp ) ;
const auto & segment_pressure = segments_copy . pressure ;
for ( int seg = 0 ; seg < nseg ; + + seg ) {
2022-12-19 08:52:36 -06:00
for ( const int perf : this - > segments_ . perforations ( ) [ seg ] ) {
2024-10-15 07:11:34 -05:00
const int local_perf_index = this - > pw_info_ . globalToLocal ( perf ) ;
if ( local_perf_index < 0 ) // then the perforation is not on this process
continue ;
const int cell_idx = this - > well_cells_ [ local_perf_index ] ;
2024-02-06 04:55:07 -06:00
const auto & intQuants = simulator . model ( ) . intensiveQuantities ( cell_idx , /*timeIdx=*/ 0 ) ;
2021-09-13 03:58:50 -05:00
// flux for each perforation
std : : vector < Scalar > mob ( this - > num_components_ , 0. ) ;
2024-10-15 07:11:34 -05:00
getMobility ( simulator , local_perf_index , mob , deferred_logger ) ;
2024-02-20 15:32:18 -06:00
const Scalar trans_mult = simulator . problem ( ) . template wellTransMultiplier < Scalar > ( intQuants , cell_idx ) ;
2024-02-06 04:55:07 -06:00
const auto & wellstate_nupcol = simulator . problem ( ) . wellModel ( ) . nupcolWellState ( ) . well ( this - > index_of_well_ ) ;
2024-10-15 07:11:34 -05:00
const std : : vector < Scalar > Tw = this - > wellIndex ( local_perf_index , intQuants , trans_mult , wellstate_nupcol ) ;
2021-09-13 03:58:50 -05:00
const Scalar seg_pressure = segment_pressure [ seg ] ;
std : : vector < Scalar > cq_s ( this - > num_components_ , 0. ) ;
2023-05-07 14:56:09 -05:00
Scalar perf_press = 0.0 ;
2024-02-20 07:18:30 -06:00
PerforationRates < Scalar > perf_rates ;
2023-05-07 14:56:09 -05:00
computePerfRate ( intQuants , mob , Tw , seg , perf , seg_pressure ,
allow_cf , cq_s , perf_press , perf_rates , deferred_logger ) ;
2021-09-13 03:58:50 -05:00
for ( int p = 0 ; p < np ; + + p ) {
2024-02-06 05:12:13 -06:00
well_flux [ this - > modelCompIdxToFlowCompIdx ( p ) ] + = cq_s [ p ] ;
2021-09-13 03:58:50 -05:00
}
}
}
this - > parallel_well_info_ . communication ( ) . sum ( well_flux . data ( ) , well_flux . size ( ) ) ;
}
template < typename TypeTag >
void
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
computeWellRatesWithBhpIterations ( const Simulator & simulator ,
2021-09-13 03:58:50 -05:00
const Scalar & bhp ,
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > & well_flux ,
2021-09-13 03:58:50 -05:00
DeferredLogger & deferred_logger ) const
{
2022-11-05 16:51:59 -05:00
// creating a copy of the well itself, to avoid messing up the explicit information
2019-11-18 09:33:10 -06:00
// during this copy, the only information not copied properly is the well controls
MultisegmentWell < TypeTag > well_copy ( * this ) ;
2024-10-10 02:25:28 -05:00
well_copy . resetDampening ( ) ;
2019-11-18 09:33:10 -06:00
well_copy . debug_cost_counter_ = 0 ;
2019-05-24 09:45:27 -05:00
2019-04-23 06:30:12 -05:00
// store a copy of the well state, we don't want to update the real well state
2024-02-17 11:13:46 -06:00
WellState < Scalar > well_state_copy = simulator . problem ( ) . wellModel ( ) . wellState ( ) ;
2024-02-06 04:55:07 -06:00
const auto & group_state = simulator . problem ( ) . wellModel ( ) . groupState ( ) ;
2021-08-03 13:05:14 -05:00
auto & ws = well_state_copy . well ( this - > index_of_well_ ) ;
2019-11-18 09:33:10 -06:00
2019-12-05 05:40:35 -06:00
// Get the current controls.
2024-02-06 04:55:07 -06:00
const auto & summary_state = simulator . vanguard ( ) . summaryState ( ) ;
2019-12-05 05:40:35 -06:00
auto inj_controls = well_copy . well_ecl_ . isInjector ( )
? well_copy . well_ecl_ . injectionControls ( summary_state )
: Well : : InjectionControls ( 0 ) ;
auto prod_controls = well_copy . well_ecl_ . isProducer ( )
? well_copy . well_ecl_ . productionControls ( summary_state ) :
Well : : ProductionControls ( 0 ) ;
// Set current control to bhp, and bhp value in state, modify bhp limit in control object.
2019-11-18 09:33:10 -06:00
if ( well_copy . well_ecl_ . isInjector ( ) ) {
2019-12-05 05:40:35 -06:00
inj_controls . bhp_limit = bhp ;
2021-08-04 05:03:36 -05:00
ws . injection_cmode = Well : : InjectorCMode : : BHP ;
2019-11-18 09:33:10 -06:00
} else {
2019-12-05 05:40:35 -06:00
prod_controls . bhp_limit = bhp ;
2021-08-04 05:03:36 -05:00
ws . production_cmode = Well : : ProducerCMode : : BHP ;
2019-11-18 09:33:10 -06:00
}
2021-08-03 13:05:14 -05:00
ws . bhp = bhp ;
2020-12-16 05:58:17 -06:00
well_copy . scaleSegmentPressuresWithBhp ( well_state_copy ) ;
2020-11-27 00:57:55 -06:00
2020-12-14 08:55:14 -06:00
// initialized the well rates with the potentials i.e. the well rates based on bhp
2021-09-06 05:58:16 -05:00
const int np = this - > number_of_phases_ ;
2022-03-23 06:29:40 -05:00
bool trivial = true ;
2020-12-14 08:55:14 -06:00
for ( int phase = 0 ; phase < np ; + + phase ) {
2022-03-23 06:29:40 -05:00
trivial = trivial & & ( ws . well_potentials [ phase ] = = 0.0 ) ;
}
if ( ! trivial ) {
2024-02-20 15:32:18 -06:00
const Scalar sign = well_copy . well_ecl_ . isInjector ( ) ? 1.0 : - 1.0 ;
2022-03-23 06:29:40 -05:00
for ( int phase = 0 ; phase < np ; + + phase ) {
ws . surface_rates [ phase ] = sign * ws . well_potentials [ phase ] ;
}
2020-12-14 08:55:14 -06:00
}
2022-12-19 08:52:36 -06:00
well_copy . scaleSegmentRatesWithWellRates ( this - > segments_ . inlets ( ) ,
this - > segments_ . perforations ( ) ,
2022-12-19 06:52:21 -06:00
well_state_copy ) ;
2020-12-14 08:55:14 -06:00
2024-02-06 04:55:07 -06:00
well_copy . calculateExplicitQuantities ( simulator , well_state_copy , deferred_logger ) ;
const double dt = simulator . timeStepSize ( ) ;
2019-11-18 09:33:10 -06:00
// iterate to get a solution at the given bhp.
2024-02-06 04:55:07 -06:00
well_copy . iterateWellEqWithControl ( simulator , dt , inj_controls , prod_controls , well_state_copy , group_state ,
2020-06-02 07:54:45 -05:00
deferred_logger ) ;
2019-04-23 06:30:12 -05:00
// compute the potential and store in the flux vector.
2019-11-18 09:33:10 -06:00
well_flux . clear ( ) ;
2019-04-23 06:30:12 -05:00
well_flux . resize ( np , 0.0 ) ;
2021-09-06 05:58:16 -05:00
for ( int compIdx = 0 ; compIdx < this - > num_components_ ; + + compIdx ) {
2022-12-19 02:52:48 -06:00
const EvalWell rate = well_copy . primary_variables_ . getQs ( compIdx ) ;
2024-02-06 05:12:13 -06:00
well_flux [ this - > modelCompIdxToFlowCompIdx ( compIdx ) ] = rate . value ( ) ;
2019-04-23 06:30:12 -05:00
}
2019-11-18 09:33:10 -06:00
debug_cost_counter_ + = well_copy . debug_cost_counter_ ;
}
template < typename TypeTag >
2024-02-20 15:32:18 -06:00
std : : vector < typename MultisegmentWell < TypeTag > : : Scalar >
2019-11-18 09:33:10 -06:00
MultisegmentWell < TypeTag > : :
2024-02-17 11:13:46 -06:00
computeWellPotentialWithTHP ( const WellState < Scalar > & well_state ,
const Simulator & simulator ,
DeferredLogger & deferred_logger ) const
2019-11-18 09:33:10 -06:00
{
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > potentials ( this - > number_of_phases_ , 0.0 ) ;
2024-02-06 04:55:07 -06:00
const auto & summary_state = simulator . vanguard ( ) . summaryState ( ) ;
2019-04-23 06:30:12 -05:00
2021-09-06 05:58:16 -05:00
const auto & well = this - > well_ecl_ ;
2024-02-20 15:32:18 -06:00
if ( well . isInjector ( ) ) {
2024-02-06 04:55:07 -06:00
auto bhp_at_thp_limit = computeBhpAtThpLimitInj ( simulator , summary_state , deferred_logger ) ;
2019-11-18 09:33:10 -06:00
if ( bhp_at_thp_limit ) {
2021-09-06 05:58:16 -05:00
const auto & controls = well . injectionControls ( summary_state ) ;
2024-02-20 15:32:18 -06:00
const Scalar bhp = std : : min ( * bhp_at_thp_limit ,
static_cast < Scalar > ( controls . bhp_limit ) ) ;
2024-02-06 04:55:07 -06:00
computeWellRatesWithBhpIterations ( simulator , bhp , potentials , deferred_logger ) ;
2019-11-18 09:33:10 -06:00
deferred_logger . debug ( " Converged thp based potential calculation for well "
2021-09-06 05:58:16 -05:00
+ this - > name ( ) + " , at bhp = " + std : : to_string ( bhp ) ) ;
2019-11-18 09:33:10 -06:00
} else {
deferred_logger . warning ( " FAILURE_GETTING_CONVERGED_POTENTIAL " ,
" Failed in getting converged thp based potential calculation for well "
2021-09-06 05:58:16 -05:00
+ this - > name ( ) + " . Instead the bhp based value is used " ) ;
2021-09-06 05:58:16 -05:00
const auto & controls = well . injectionControls ( summary_state ) ;
2024-02-20 15:32:18 -06:00
const Scalar bhp = controls . bhp_limit ;
2024-02-06 04:55:07 -06:00
computeWellRatesWithBhpIterations ( simulator , bhp , potentials , deferred_logger ) ;
2019-11-18 09:33:10 -06:00
}
} else {
2022-02-07 04:28:35 -06:00
auto bhp_at_thp_limit = computeBhpAtThpLimitProd (
2024-02-06 04:55:07 -06:00
well_state , simulator , summary_state , deferred_logger ) ;
2019-11-18 09:33:10 -06:00
if ( bhp_at_thp_limit ) {
2021-09-06 05:58:16 -05:00
const auto & controls = well . productionControls ( summary_state ) ;
2024-02-20 15:32:18 -06:00
const Scalar bhp = std : : max ( * bhp_at_thp_limit ,
static_cast < Scalar > ( controls . bhp_limit ) ) ;
2024-02-06 04:55:07 -06:00
computeWellRatesWithBhpIterations ( simulator , bhp , potentials , deferred_logger ) ;
2019-11-18 09:33:10 -06:00
deferred_logger . debug ( " Converged thp based potential calculation for well "
2021-09-06 05:58:16 -05:00
+ this - > name ( ) + " , at bhp = " + std : : to_string ( bhp ) ) ;
2019-11-18 09:33:10 -06:00
} else {
deferred_logger . warning ( " FAILURE_GETTING_CONVERGED_POTENTIAL " ,
" Failed in getting converged thp based potential calculation for well "
2021-09-06 05:58:16 -05:00
+ this - > name ( ) + " . Instead the bhp based value is used " ) ;
2021-09-06 05:58:16 -05:00
const auto & controls = well . productionControls ( summary_state ) ;
2024-02-20 15:32:18 -06:00
const Scalar bhp = controls . bhp_limit ;
2024-02-06 04:55:07 -06:00
computeWellRatesWithBhpIterations ( simulator , bhp , potentials , deferred_logger ) ;
2019-11-18 09:33:10 -06:00
}
}
return potentials ;
2019-04-23 06:30:12 -05:00
}
2017-08-29 10:26:36 -05:00
2023-10-26 10:28:05 -05:00
template < typename TypeTag >
bool
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
computeWellPotentialsImplicit ( const Simulator & simulator ,
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > & well_potentials ,
2023-10-26 10:28:05 -05:00
DeferredLogger & deferred_logger ) const
{
// Create a copy of the well.
// TODO: check if we can avoid taking multiple copies. Call from updateWellPotentials
// is allready a copy, but not from other calls.
MultisegmentWell < TypeTag > well_copy ( * this ) ;
well_copy . debug_cost_counter_ = 0 ;
// store a copy of the well state, we don't want to update the real well state
2024-02-17 11:13:46 -06:00
WellState < Scalar > well_state_copy = simulator . problem ( ) . wellModel ( ) . wellState ( ) ;
2024-02-06 04:55:07 -06:00
const auto & group_state = simulator . problem ( ) . wellModel ( ) . groupState ( ) ;
2023-10-26 10:28:05 -05:00
auto & ws = well_state_copy . well ( this - > index_of_well_ ) ;
// get current controls
2024-02-06 04:55:07 -06:00
const auto & summary_state = simulator . vanguard ( ) . summaryState ( ) ;
2023-10-26 10:28:05 -05:00
auto inj_controls = well_copy . well_ecl_ . isInjector ( )
? well_copy . well_ecl_ . injectionControls ( summary_state )
: Well : : InjectionControls ( 0 ) ;
auto prod_controls = well_copy . well_ecl_ . isProducer ( )
? well_copy . well_ecl_ . productionControls ( summary_state )
: Well : : ProductionControls ( 0 ) ;
// prepare/modify well state and control
well_copy . prepareForPotentialCalculations ( summary_state , well_state_copy , inj_controls , prod_controls ) ;
well_copy . scaleSegmentPressuresWithBhp ( well_state_copy ) ;
// initialize rates from previous potentials
const int np = this - > number_of_phases_ ;
bool trivial = true ;
for ( int phase = 0 ; phase < np ; + + phase ) {
trivial = trivial & & ( ws . well_potentials [ phase ] = = 0.0 ) ;
}
if ( ! trivial ) {
2024-02-20 15:32:18 -06:00
const Scalar sign = well_copy . well_ecl_ . isInjector ( ) ? 1.0 : - 1.0 ;
2023-10-26 10:28:05 -05:00
for ( int phase = 0 ; phase < np ; + + phase ) {
ws . surface_rates [ phase ] = sign * ws . well_potentials [ phase ] ;
}
}
well_copy . scaleSegmentRatesWithWellRates ( this - > segments_ . inlets ( ) ,
this - > segments_ . perforations ( ) ,
well_state_copy ) ;
2024-02-06 04:55:07 -06:00
well_copy . calculateExplicitQuantities ( simulator , well_state_copy , deferred_logger ) ;
const double dt = simulator . timeStepSize ( ) ;
2023-10-26 10:28:05 -05:00
// solve equations
2023-11-07 10:29:50 -06:00
bool converged = false ;
if ( this - > well_ecl_ . isProducer ( ) & & this - > wellHasTHPConstraints ( summary_state ) ) {
2024-02-06 04:55:07 -06:00
converged = well_copy . solveWellWithTHPConstraint ( simulator , dt , inj_controls , prod_controls , well_state_copy , group_state , deferred_logger ) ;
2023-11-07 10:29:50 -06:00
} else {
2024-02-06 04:55:07 -06:00
converged = well_copy . iterateWellEqWithSwitching ( simulator , dt , inj_controls , prod_controls , well_state_copy , group_state , deferred_logger ) ;
2023-11-07 10:29:50 -06:00
}
2017-08-29 10:26:36 -05:00
2023-10-26 10:28:05 -05:00
// fetch potentials (sign is updated on the outside).
well_potentials . clear ( ) ;
well_potentials . resize ( np , 0.0 ) ;
for ( int compIdx = 0 ; compIdx < this - > num_components_ ; + + compIdx ) {
const EvalWell rate = well_copy . primary_variables_ . getQs ( compIdx ) ;
2024-02-06 05:12:13 -06:00
well_potentials [ this - > modelCompIdxToFlowCompIdx ( compIdx ) ] = rate . value ( ) ;
2023-10-26 10:28:05 -05:00
}
debug_cost_counter_ + = well_copy . debug_cost_counter_ ;
return converged ;
}
2017-08-29 10:26:36 -05:00
2023-04-21 02:10:10 -05:00
template < typename TypeTag >
void
MultisegmentWell < TypeTag > : :
2024-03-19 06:50:34 -05:00
solveEqAndUpdateWellState ( const Simulator & simulator ,
2024-02-17 11:13:46 -06:00
WellState < Scalar > & well_state ,
2023-04-21 02:10:10 -05:00
DeferredLogger & deferred_logger )
{
if ( ! this - > isOperableAndSolvable ( ) & & ! this - > wellIsStopped ( ) ) return ;
// We assemble the well equations, then we check the convergence,
// which is why we do not put the assembleWellEq here.
2023-12-04 16:15:08 -06:00
try {
const BVectorWell dx_well = this - > linSys_ . solve ( ) ;
2023-04-21 02:10:10 -05:00
2024-03-19 06:50:34 -05:00
updateWellState ( simulator , dx_well , well_state , deferred_logger ) ;
2023-12-04 16:15:08 -06:00
}
catch ( const NumericalProblem & exp ) {
// Add information about the well and log to deferred logger
// (Logging done inside of solve() method will only be seen if
// this is the process with rank zero)
deferred_logger . problem ( " In MultisegmentWell::solveEqAndUpdateWellState for well "
+ this - > name ( ) + " : " + exp . what ( ) ) ;
throw ;
}
2023-04-21 02:10:10 -05:00
}
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2017-09-05 08:16:04 -05:00
void
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
computePerfCellPressDiffs ( const Simulator & simulator )
2017-09-05 08:16:04 -05:00
{
2021-09-06 05:58:16 -05:00
for ( int perf = 0 ; perf < this - > number_of_perforations_ ; + + perf ) {
2017-09-29 02:39:32 -05:00
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > kr ( this - > number_of_phases_ , 0.0 ) ;
std : : vector < Scalar > density ( this - > number_of_phases_ , 0.0 ) ;
2017-09-29 02:39:32 -05:00
2021-09-06 05:58:16 -05:00
const int cell_idx = this - > well_cells_ [ perf ] ;
2024-02-06 04:55:07 -06:00
const auto & intQuants = simulator . model ( ) . intensiveQuantities ( cell_idx , /*timeIdx=*/ 0 ) ;
2017-09-29 02:39:32 -05:00
const auto & fs = intQuants . fluidState ( ) ;
2024-02-20 15:32:18 -06:00
Scalar sum_kr = 0. ;
2017-09-29 02:39:32 -05:00
2021-09-06 05:58:16 -05:00
const PhaseUsage & pu = this - > phaseUsage ( ) ;
2017-11-16 04:42:34 -06:00
if ( pu . phase_used [ Water ] ) {
const int water_pos = pu . phase_pos [ Water ] ;
2017-09-29 02:39:32 -05:00
kr [ water_pos ] = intQuants . relativePermeability ( FluidSystem : : waterPhaseIdx ) . value ( ) ;
sum_kr + = kr [ water_pos ] ;
density [ water_pos ] = fs . density ( FluidSystem : : waterPhaseIdx ) . value ( ) ;
}
2017-11-16 04:42:34 -06:00
if ( pu . phase_used [ Oil ] ) {
const int oil_pos = pu . phase_pos [ Oil ] ;
2017-09-29 02:39:32 -05:00
kr [ oil_pos ] = intQuants . relativePermeability ( FluidSystem : : oilPhaseIdx ) . value ( ) ;
sum_kr + = kr [ oil_pos ] ;
density [ oil_pos ] = fs . density ( FluidSystem : : oilPhaseIdx ) . value ( ) ;
}
2017-11-16 04:42:34 -06:00
if ( pu . phase_used [ Gas ] ) {
const int gas_pos = pu . phase_pos [ Gas ] ;
2017-09-29 02:39:32 -05:00
kr [ gas_pos ] = intQuants . relativePermeability ( FluidSystem : : gasPhaseIdx ) . value ( ) ;
sum_kr + = kr [ gas_pos ] ;
density [ gas_pos ] = fs . density ( FluidSystem : : gasPhaseIdx ) . value ( ) ;
}
assert ( sum_kr ! = 0. ) ;
// calculate the average density
2024-02-20 15:32:18 -06:00
Scalar average_density = 0. ;
2021-09-06 05:58:16 -05:00
for ( int p = 0 ; p < this - > number_of_phases_ ; + + p ) {
2017-09-29 02:39:32 -05:00
average_density + = kr [ p ] * density [ p ] ;
}
average_density / = sum_kr ;
2021-09-06 05:58:16 -05:00
this - > cell_perforation_pressure_diffs_ [ perf ] = this - > gravity_ * average_density * this - > cell_perforation_depth_diffs_ [ perf ] ;
2017-09-29 02:39:32 -05:00
}
2017-09-05 08:16:04 -05:00
}
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2017-09-05 08:16:04 -05:00
void
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
computeInitialSegmentFluids ( const Simulator & simulator )
2017-09-05 08:16:04 -05:00
{
2021-06-03 08:34:14 -05:00
for ( int seg = 0 ; seg < this - > numberOfSegments ( ) ; + + seg ) {
2019-02-12 05:48:51 -06:00
// TODO: trying to reduce the times for the surfaceVolumeFraction calculation
2024-02-20 15:32:18 -06:00
const Scalar surface_volume = getSegmentSurfaceVolume ( simulator , seg ) . value ( ) ;
2021-09-06 05:58:16 -05:00
for ( int comp_idx = 0 ; comp_idx < this - > num_components_ ; + + comp_idx ) {
2022-12-19 02:52:48 -06:00
segment_fluid_initial_ [ seg ] [ comp_idx ] = surface_volume * this - > primary_variables_ . surfaceVolumeFraction ( seg , comp_idx ) . value ( ) ;
2017-09-05 08:34:31 -05:00
}
}
2017-09-05 08:16:04 -05:00
}
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2017-09-12 10:32:27 -05:00
void
MultisegmentWell < TypeTag > : :
2024-03-19 06:50:34 -05:00
updateWellState ( const Simulator & simulator ,
2023-03-22 09:10:00 -05:00
const BVectorWell & dwells ,
2024-02-17 11:13:46 -06:00
WellState < Scalar > & well_state ,
2021-05-05 04:22:44 -05:00
DeferredLogger & deferred_logger ,
2024-02-20 15:32:18 -06:00
const Scalar relaxation_factor )
2017-09-12 10:32:27 -05:00
{
2021-09-29 09:01:16 -05:00
if ( ! this - > isOperableAndSolvable ( ) & & ! this - > wellIsStopped ( ) ) return ;
2020-11-27 00:57:55 -06:00
2024-02-20 15:32:18 -06:00
const Scalar dFLimit = this - > param_ . dwell_fraction_max_ ;
const Scalar max_pressure_change = this - > param_ . max_pressure_change_ms_wells_ ;
2024-03-19 06:50:34 -05:00
const bool stop_or_zero_rate_target =
this - > stoppedOrZeroRateTarget ( simulator , well_state , deferred_logger ) ;
2022-12-19 02:52:48 -06:00
this - > primary_variables_ . updateNewton ( dwells ,
relaxation_factor ,
dFLimit ,
2023-03-29 07:33:05 -05:00
stop_or_zero_rate_target ,
2024-03-02 09:21:53 -06:00
max_pressure_change ) ;
2017-09-15 09:28:10 -05:00
2024-03-19 06:50:34 -05:00
const auto & summary_state = simulator . vanguard ( ) . summaryState ( ) ;
this - > primary_variables_ . copyToWellState ( * this , getRefDensity ( ) ,
well_state ,
summary_state ,
deferred_logger ) ;
2023-07-04 11:40:13 -05:00
{
auto & ws = well_state . well ( this - > index_of_well_ ) ;
this - > segments_ . copyPhaseDensities ( ws . pu , ws . segments ) ;
}
2024-07-29 04:27:39 -05:00
Base : : calculateReservoirRates ( simulator . vanguard ( ) . eclState ( ) . runspec ( ) . co2Storage ( ) , well_state . well ( this - > index_of_well_ ) ) ;
2017-09-12 10:32:27 -05:00
}
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2017-09-05 04:09:58 -05:00
void
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
calculateExplicitQuantities ( const Simulator & simulator ,
2024-02-17 11:13:46 -06:00
const WellState < Scalar > & well_state ,
2021-05-05 04:22:44 -05:00
DeferredLogger & deferred_logger )
2017-09-05 04:09:58 -05:00
{
2024-03-19 06:50:34 -05:00
updatePrimaryVariables ( simulator , well_state , deferred_logger ) ;
2019-02-27 03:46:10 -06:00
initPrimaryVariablesEvaluation ( ) ;
2024-02-06 04:55:07 -06:00
computePerfCellPressDiffs ( simulator ) ;
computeInitialSegmentFluids ( simulator ) ;
2017-09-05 04:09:58 -05:00
}
2020-10-09 06:38:33 -05:00
template < typename TypeTag >
void
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
updateProductivityIndex ( const Simulator & simulator ,
2024-02-20 16:39:23 -06:00
const WellProdIndexCalculator < Scalar > & wellPICalc ,
2024-02-17 11:13:46 -06:00
WellState < Scalar > & well_state ,
2020-10-26 11:33:07 -05:00
DeferredLogger & deferred_logger ) const
2020-10-09 06:38:33 -05:00
{
2024-02-06 04:55:07 -06:00
auto fluidState = [ & simulator , this ] ( const int perf )
2020-10-09 06:38:33 -05:00
{
const auto cell_idx = this - > well_cells_ [ perf ] ;
2024-02-06 04:55:07 -06:00
return simulator . model ( )
2023-03-24 13:57:31 -05:00
. intensiveQuantities ( cell_idx , /*timeIdx=*/ 0 ) . fluidState ( ) ;
2020-10-09 06:38:33 -05:00
} ;
2020-10-26 11:33:07 -05:00
const int np = this - > number_of_phases_ ;
2024-02-20 15:32:18 -06:00
auto setToZero = [ np ] ( Scalar * x ) - > void
2020-10-09 06:38:33 -05:00
{
2020-10-26 11:33:07 -05:00
std : : fill_n ( x , np , 0.0 ) ;
} ;
2020-10-09 06:38:33 -05:00
2024-02-20 15:32:18 -06:00
auto addVector = [ np ] ( const Scalar * src , Scalar * dest ) - > void
2020-10-26 11:33:07 -05:00
{
std : : transform ( src , src + np , dest , dest , std : : plus < > { } ) ;
2020-10-09 06:38:33 -05:00
} ;
2021-08-05 05:17:57 -05:00
auto & ws = well_state . well ( this - > index_of_well_ ) ;
2021-08-06 06:14:00 -05:00
auto & perf_data = ws . perf_data ;
2021-06-07 02:10:53 -05:00
auto * connPI = perf_data . prod_index . data ( ) ;
2021-08-05 05:17:57 -05:00
auto * wellPI = ws . productivity_index . data ( ) ;
2020-10-09 06:38:33 -05:00
2020-10-26 11:33:07 -05:00
setToZero ( wellPI ) ;
2020-10-09 06:38:33 -05:00
2020-11-06 14:13:10 -06:00
const auto preferred_phase = this - > well_ecl_ . getPreferredPhase ( ) ;
2020-12-01 06:54:50 -06:00
auto subsetPerfID = 0 ;
2020-11-06 14:13:10 -06:00
2020-12-01 06:54:50 -06:00
for ( const auto & perf : * this - > perf_data_ ) {
auto allPerfID = perf . ecl_index ;
2020-10-09 06:38:33 -05:00
2024-02-20 15:32:18 -06:00
auto connPICalc = [ & wellPICalc , allPerfID ] ( const Scalar mobility ) - > Scalar
2020-11-06 14:13:10 -06:00
{
return wellPICalc . connectionProdIndStandard ( allPerfID , mobility ) ;
} ;
2021-09-13 02:36:16 -05:00
std : : vector < Scalar > mob ( this - > num_components_ , 0.0 ) ;
2024-10-15 07:11:34 -05:00
// The subsetPerfID loops over 0 .. this->perf_data_->size().
// *(this->perf_data_) contains info about the local processes only,
// hence subsetPerfID is a local perf id and we can call getMobility
// directly with that.
2024-02-06 04:55:07 -06:00
getMobility ( simulator , static_cast < int > ( subsetPerfID ) , mob , deferred_logger ) ;
2020-10-09 06:38:33 -05:00
2020-10-26 11:33:07 -05:00
const auto & fs = fluidState ( subsetPerfID ) ;
setToZero ( connPI ) ;
2020-10-09 06:38:33 -05:00
2020-11-06 14:13:10 -06:00
if ( this - > isInjector ( ) ) {
this - > computeConnLevelInjInd ( fs , preferred_phase , connPICalc ,
mob , connPI , deferred_logger ) ;
2020-10-09 06:38:33 -05:00
}
2020-11-06 14:17:53 -06:00
else { // Production or zero flow rate
2020-11-06 14:13:10 -06:00
this - > computeConnLevelProdInd ( fs , connPICalc , mob , connPI ) ;
2020-11-06 14:17:53 -06:00
}
2020-10-09 06:38:33 -05:00
2020-10-26 11:33:07 -05:00
addVector ( connPI , wellPI ) ;
2020-10-09 06:38:33 -05:00
+ + subsetPerfID ;
connPI + = np ;
}
2020-11-06 14:15:25 -06:00
2024-09-30 02:05:25 -05:00
// Sum with communication in case of distributed well.
const auto & comm = this - > parallel_well_info_ . communication ( ) ;
if ( comm . size ( ) > 1 ) {
comm . sum ( wellPI , np ) ;
}
2020-11-06 14:15:25 -06:00
assert ( static_cast < int > ( subsetPerfID ) = = this - > number_of_perforations_ & &
" Internal logic error in processing connections for PI/II " ) ;
2020-10-09 06:38:33 -05:00
}
2023-06-06 14:31:17 -05:00
template < typename TypeTag >
2024-02-20 15:32:18 -06:00
typename MultisegmentWell < TypeTag > : : Scalar
2023-06-06 14:31:17 -05:00
MultisegmentWell < TypeTag > : :
connectionDensity ( const int globalConnIdx ,
[[maybe_unused]] const int openConnIdx ) const
{
// Simple approximation: Mixture density at reservoir connection is
// mixture density at connection's segment.
const auto segNum = this - > wellEcl ( )
. getConnections ( ) [ globalConnIdx ] . segment ( ) ;
const auto segIdx = this - > wellEcl ( )
. getSegments ( ) . segmentNumberToIndex ( segNum ) ;
return this - > segments_ . density ( segIdx ) . value ( ) ;
}
2019-08-15 04:29:30 -05:00
template < typename TypeTag >
void
MultisegmentWell < TypeTag > : :
2020-10-29 06:57:29 -05:00
addWellContributions ( SparseMatrixAdapter & jacobian ) const
2019-08-15 04:29:30 -05:00
{
2022-11-11 14:41:24 -06:00
this - > linSys_ . extract ( jacobian ) ;
2019-08-15 04:29:30 -05:00
}
2022-04-28 02:32:07 -05:00
template < typename TypeTag >
void
MultisegmentWell < TypeTag > : :
addWellPressureEquations ( PressureMatrix & jacobian ,
const BVector & weights ,
const int pressureVarIndex ,
2022-11-11 14:41:24 -06:00
const bool use_well_weights ,
2024-02-17 11:13:46 -06:00
const WellState < Scalar > & well_state ) const
2022-04-28 02:32:07 -05:00
{
2022-05-05 07:21:28 -05:00
// Add the pressure contribution to the cpr system for the well
2022-11-11 14:41:24 -06:00
this - > linSys_ . extractCPRPressureMatrix ( jacobian ,
weights ,
pressureVarIndex ,
use_well_weights ,
* this ,
this - > SPres ,
well_state ) ;
2022-04-28 02:32:07 -05:00
}
2022-11-07 23:38:12 -06:00
2019-08-15 04:29:30 -05:00
2021-09-13 02:36:16 -05:00
template < typename TypeTag >
template < class Value >
void
MultisegmentWell < TypeTag > : :
computePerfRate ( const Value & pressure_cell ,
const Value & rs ,
const Value & rv ,
const std : : vector < Value > & b_perfcells ,
const std : : vector < Value > & mob_perfcells ,
2023-11-14 05:45:25 -06:00
const std : : vector < Scalar > & Tw ,
2021-09-13 02:36:16 -05:00
const int perf ,
const Value & segment_pressure ,
const Value & segment_density ,
const bool & allow_cf ,
const std : : vector < Value > & cmix_s ,
std : : vector < Value > & cq_s ,
Value & perf_press ,
2024-02-20 07:18:30 -06:00
PerforationRates < Scalar > & perf_rates ,
2021-09-13 02:36:16 -05:00
DeferredLogger & deferred_logger ) const
{
2024-10-15 07:11:34 -05:00
const int local_perf_index = this - > pw_info_ . globalToLocal ( perf ) ;
if ( local_perf_index < 0 ) // then the perforation is not on this process
return ;
2021-09-13 02:36:16 -05:00
// pressure difference between the segment and the perforation
2022-12-19 08:52:36 -06:00
const Value perf_seg_press_diff = this - > gravity ( ) * segment_density *
this - > segments_ . perforation_depth_diff ( perf ) ;
2021-09-13 02:36:16 -05:00
// pressure difference between the perforation and the grid cell
2024-10-15 07:11:34 -05:00
const Scalar cell_perf_press_diff = this - > cell_perforation_pressure_diffs_ [ local_perf_index ] ;
2021-09-13 02:36:16 -05:00
2023-02-17 05:17:41 -06:00
// perforation pressure is the wellbore pressure corrected to perforation depth
// (positive sign due to convention in segments_.perforation_depth_diff() )
perf_press = segment_pressure + perf_seg_press_diff ;
2023-03-30 08:49:27 -05:00
// cell pressure corrected to perforation depth
2023-02-17 05:17:41 -06:00
const Value cell_press_at_perf = pressure_cell - cell_perf_press_diff ;
2023-03-30 08:49:27 -05:00
2021-09-13 02:36:16 -05:00
// Pressure drawdown (also used to determine direction of flow)
2023-02-17 05:17:41 -06:00
const Value drawdown = cell_press_at_perf - perf_press ;
2021-09-13 02:36:16 -05:00
// producing perforations
2023-05-05 06:58:23 -05:00
if ( drawdown > 0.0 ) {
// Do nothing if crossflow is not allowed
2021-09-13 02:36:16 -05:00
if ( ! allow_cf & & this - > isInjector ( ) ) {
return ;
}
// compute component volumetric rates at standard conditions
for ( int comp_idx = 0 ; comp_idx < this - > numComponents ( ) ; + + comp_idx ) {
2023-11-14 05:45:25 -06:00
const Value cq_p = - Tw [ comp_idx ] * ( mob_perfcells [ comp_idx ] * drawdown ) ;
2021-09-13 02:36:16 -05:00
cq_s [ comp_idx ] = b_perfcells [ comp_idx ] * cq_p ;
}
2019-08-15 04:29:30 -05:00
2021-09-13 02:36:16 -05:00
if ( FluidSystem : : phaseIsActive ( FluidSystem : : oilPhaseIdx ) & & FluidSystem : : phaseIsActive ( FluidSystem : : gasPhaseIdx ) ) {
const unsigned oilCompIdx = Indices : : canonicalToActiveComponentIndex ( FluidSystem : : oilCompIdx ) ;
const unsigned gasCompIdx = Indices : : canonicalToActiveComponentIndex ( FluidSystem : : gasCompIdx ) ;
const Value cq_s_oil = cq_s [ oilCompIdx ] ;
const Value cq_s_gas = cq_s [ gasCompIdx ] ;
cq_s [ gasCompIdx ] + = rs * cq_s_oil ;
cq_s [ oilCompIdx ] + = rv * cq_s_gas ;
}
} else { // injecting perforations
// Do nothing if crossflow is not allowed
if ( ! allow_cf & & this - > isProducer ( ) ) {
return ;
}
// for injecting perforations, we use total mobility
Value total_mob = mob_perfcells [ 0 ] ;
for ( int comp_idx = 1 ; comp_idx < this - > numComponents ( ) ; + + comp_idx ) {
total_mob + = mob_perfcells [ comp_idx ] ;
}
// compute volume ratio between connection and at standard conditions
Value volume_ratio = 0.0 ;
if ( FluidSystem : : phaseIsActive ( FluidSystem : : waterPhaseIdx ) ) {
const unsigned waterCompIdx = Indices : : canonicalToActiveComponentIndex ( FluidSystem : : waterCompIdx ) ;
volume_ratio + = cmix_s [ waterCompIdx ] / b_perfcells [ waterCompIdx ] ;
}
if ( FluidSystem : : phaseIsActive ( FluidSystem : : oilPhaseIdx ) & & FluidSystem : : phaseIsActive ( FluidSystem : : gasPhaseIdx ) ) {
const unsigned oilCompIdx = Indices : : canonicalToActiveComponentIndex ( FluidSystem : : oilCompIdx ) ;
const unsigned gasCompIdx = Indices : : canonicalToActiveComponentIndex ( FluidSystem : : gasCompIdx ) ;
// Incorporate RS/RV factors if both oil and gas active
// TODO: not sure we use rs rv from the perforation cells when handling injecting perforations
// basically, for injecting perforations, the wellbore is the upstreaming side.
const Value d = 1.0 - rv * rs ;
if ( getValue ( d ) = = 0.0 ) {
2023-11-08 08:16:17 -06:00
OPM_DEFLOG_PROBLEM ( NumericalProblem ,
fmt : : format ( " Zero d value obtained for well {} "
" during flux calculation with rs {} and rv {} " ,
this - > name ( ) , rs , rv ) ,
deferred_logger ) ;
2021-09-13 02:36:16 -05:00
}
const Value tmp_oil = ( cmix_s [ oilCompIdx ] - rv * cmix_s [ gasCompIdx ] ) / d ;
volume_ratio + = tmp_oil / b_perfcells [ oilCompIdx ] ;
const Value tmp_gas = ( cmix_s [ gasCompIdx ] - rs * cmix_s [ oilCompIdx ] ) / d ;
volume_ratio + = tmp_gas / b_perfcells [ gasCompIdx ] ;
} else { // not having gas and oil at the same time
if ( FluidSystem : : phaseIsActive ( FluidSystem : : oilPhaseIdx ) ) {
const unsigned oilCompIdx = Indices : : canonicalToActiveComponentIndex ( FluidSystem : : oilCompIdx ) ;
volume_ratio + = cmix_s [ oilCompIdx ] / b_perfcells [ oilCompIdx ] ;
}
if ( FluidSystem : : phaseIsActive ( FluidSystem : : gasPhaseIdx ) ) {
const unsigned gasCompIdx = Indices : : canonicalToActiveComponentIndex ( FluidSystem : : gasCompIdx ) ;
volume_ratio + = cmix_s [ gasCompIdx ] / b_perfcells [ gasCompIdx ] ;
}
}
// injecting connections total volumerates at standard conditions
2023-11-14 05:45:25 -06:00
for ( int componentIdx = 0 ; componentIdx < this - > numComponents ( ) ; + + componentIdx ) {
2023-11-14 07:06:23 -06:00
const Value cqt_i = - Tw [ componentIdx ] * ( total_mob * drawdown ) ;
Value cqt_is = cqt_i / volume_ratio ;
cq_s [ componentIdx ] = cmix_s [ componentIdx ] * cqt_is ;
2021-09-13 02:36:16 -05:00
}
} // end for injection perforations
// calculating the perforation solution gas rate and solution oil rates
if ( this - > isProducer ( ) ) {
if ( FluidSystem : : phaseIsActive ( FluidSystem : : oilPhaseIdx ) & & FluidSystem : : phaseIsActive ( FluidSystem : : gasPhaseIdx ) ) {
const unsigned oilCompIdx = Indices : : canonicalToActiveComponentIndex ( FluidSystem : : oilCompIdx ) ;
const unsigned gasCompIdx = Indices : : canonicalToActiveComponentIndex ( FluidSystem : : gasCompIdx ) ;
// TODO: the formulations here remain to be tested with cases with strong crossflow through production wells
// s means standard condition, r means reservoir condition
// q_os = q_or * b_o + rv * q_gr * b_g
// q_gs = q_gr * g_g + rs * q_or * b_o
// d = 1.0 - rs * rv
// q_or = 1 / (b_o * d) * (q_os - rv * q_gs)
// q_gr = 1 / (b_g * d) * (q_gs - rs * q_os)
2024-02-20 15:32:18 -06:00
const Scalar d = 1.0 - getValue ( rv ) * getValue ( rs ) ;
2021-09-13 02:36:16 -05:00
// vaporized oil into gas
// rv * q_gr * b_g = rv * (q_gs - rs * q_os) / d
2023-05-04 06:02:17 -05:00
perf_rates . vap_oil = getValue ( rv ) * ( getValue ( cq_s [ gasCompIdx ] ) - getValue ( rs ) * getValue ( cq_s [ oilCompIdx ] ) ) / d ;
2021-09-13 02:36:16 -05:00
// dissolved of gas in oil
// rs * q_or * b_o = rs * (q_os - rv * q_gs) / d
2023-05-04 06:02:17 -05:00
perf_rates . dis_gas = getValue ( rs ) * ( getValue ( cq_s [ oilCompIdx ] ) - getValue ( rv ) * getValue ( cq_s [ gasCompIdx ] ) ) / d ;
2021-09-13 02:36:16 -05:00
}
}
}
2019-08-15 04:29:30 -05:00
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2023-05-07 14:56:09 -05:00
template < class Value >
2017-09-08 04:31:43 -05:00
void
MultisegmentWell < TypeTag > : :
2023-05-07 14:56:09 -05:00
computePerfRate ( const IntensiveQuantities & int_quants ,
const std : : vector < Value > & mob_perfcells ,
2023-11-14 05:45:25 -06:00
const std : : vector < Scalar > & Tw ,
2023-05-07 14:56:09 -05:00
const int seg ,
const int perf ,
const Value & segment_pressure ,
const bool & allow_cf ,
std : : vector < Value > & cq_s ,
Value & perf_press ,
2024-02-20 07:18:30 -06:00
PerforationRates < Scalar > & perf_rates ,
2023-05-07 14:56:09 -05:00
DeferredLogger & deferred_logger ) const
2019-02-07 07:43:17 -06:00
2017-09-08 04:31:43 -05:00
{
2023-05-07 14:56:09 -05:00
auto obtain = [ this ] ( const Eval & value )
{
if constexpr ( std : : is_same_v < Value , Scalar > ) {
2023-06-27 06:41:12 -05:00
static_cast < void > ( this ) ; // suppress clang warning
2023-05-07 14:56:09 -05:00
return getValue ( value ) ;
} else {
return this - > extendEval ( value ) ;
}
} ;
auto obtainN = [ ] ( const auto & value )
{
if constexpr ( std : : is_same_v < Value , Scalar > ) {
return getValue ( value ) ;
} else {
return value ;
2017-12-04 02:14:08 -06:00
}
2023-05-07 14:56:09 -05:00
} ;
2021-09-13 02:36:16 -05:00
const auto & fs = int_quants . fluidState ( ) ;
2023-05-07 14:56:09 -05:00
const Value pressure_cell = obtain ( this - > getPerfCellPressure ( fs ) ) ;
const Value rs = obtain ( fs . Rs ( ) ) ;
const Value rv = obtain ( fs . Rv ( ) ) ;
2021-09-13 02:36:16 -05:00
// not using number_of_phases_ because of solvent
2023-05-07 14:56:09 -05:00
std : : vector < Value > b_perfcells ( this - > num_components_ , 0.0 ) ;
2017-09-08 04:09:54 -05:00
2021-09-13 02:36:16 -05:00
for ( unsigned phaseIdx = 0 ; phaseIdx < FluidSystem : : numPhases ; + + phaseIdx ) {
if ( ! FluidSystem : : phaseIsActive ( phaseIdx ) ) {
continue ;
}
const unsigned compIdx = Indices : : canonicalToActiveComponentIndex ( FluidSystem : : solventComponentIndex ( phaseIdx ) ) ;
2023-05-07 14:56:09 -05:00
b_perfcells [ compIdx ] = obtain ( fs . invB ( phaseIdx ) ) ;
2021-09-13 02:36:16 -05:00
}
2023-05-07 14:56:09 -05:00
std : : vector < Value > cmix_s ( this - > numComponents ( ) , 0.0 ) ;
2021-09-13 02:36:16 -05:00
for ( int comp_idx = 0 ; comp_idx < this - > numComponents ( ) ; + + comp_idx ) {
2023-05-07 14:56:09 -05:00
cmix_s [ comp_idx ] = obtainN ( this - > primary_variables_ . surfaceVolumeFraction ( seg , comp_idx ) ) ;
2021-09-13 02:36:16 -05:00
}
this - > computePerfRate ( pressure_cell ,
rs ,
rv ,
b_perfcells ,
mob_perfcells ,
Tw ,
perf ,
segment_pressure ,
2023-05-07 14:56:09 -05:00
obtainN ( this - > segments_ . density ( seg ) ) ,
2021-09-13 02:36:16 -05:00
allow_cf ,
cmix_s ,
cq_s ,
perf_press ,
2023-05-04 06:02:17 -05:00
perf_rates ,
2021-09-13 02:36:16 -05:00
deferred_logger ) ;
}
2017-09-08 04:09:54 -05:00
2017-09-28 08:12:09 -05:00
template < typename TypeTag >
2017-09-08 04:09:54 -05:00
void
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
computeSegmentFluidProperties ( const Simulator & simulator , DeferredLogger & deferred_logger )
2017-09-08 04:09:54 -05:00
{
2017-09-28 03:04:04 -05:00
// TODO: the concept of phases and components are rather confusing in this function.
// needs to be addressed sooner or later.
2017-09-08 04:09:54 -05:00
// get the temperature for later use. It is only useful when we are not handling
// thermal related simulation
// basically, it is a single value for all the segments
2017-10-12 08:36:54 -05:00
2017-09-08 04:09:54 -05:00
EvalWell temperature ;
2020-07-02 06:44:01 -05:00
EvalWell saltConcentration ;
2017-09-08 04:09:54 -05:00
// not sure how to handle the pvt region related to segment
// for the current approach, we use the pvt region of the first perforated cell
// although there are some text indicating using the pvt region of the lowest
// perforated cell
// TODO: later to investigate how to handle the pvt region
int pvt_region_index ;
{
// using the first perforated cell
2021-09-06 05:58:16 -05:00
const int cell_idx = this - > well_cells_ [ 0 ] ;
2024-02-06 04:55:07 -06:00
const auto & intQuants = simulator . model ( ) . intensiveQuantities ( cell_idx , /*timeIdx=*/ 0 ) ;
2017-09-08 04:09:54 -05:00
const auto & fs = intQuants . fluidState ( ) ;
temperature . setValue ( fs . temperature ( FluidSystem : : oilPhaseIdx ) . value ( ) ) ;
2021-06-03 08:34:14 -05:00
saltConcentration = this - > extendEval ( fs . saltConcentration ( ) ) ;
2017-09-08 04:09:54 -05:00
pvt_region_index = fs . pvtRegionIndex ( ) ;
}
2022-12-19 07:20:55 -06:00
this - > segments_ . computeFluidProperties ( temperature ,
saltConcentration ,
this - > primary_variables_ ,
pvt_region_index ,
deferred_logger ) ;
2021-06-03 08:34:14 -05:00
}
2017-12-04 02:14:08 -06:00
2021-09-13 02:36:16 -05:00
template < typename TypeTag >
2023-05-07 14:56:09 -05:00
template < class Value >
2021-09-13 02:36:16 -05:00
void
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
getMobility ( const Simulator & simulator ,
2024-10-15 07:11:34 -05:00
const int local_perf_index ,
2023-05-12 03:02:05 -05:00
std : : vector < Value > & mob ,
DeferredLogger & deferred_logger ) const
2021-09-13 02:36:16 -05:00
{
2023-05-07 14:56:09 -05:00
auto obtain = [ this ] ( const Eval & value )
{
if constexpr ( std : : is_same_v < Value , Scalar > ) {
2023-06-27 06:41:12 -05:00
static_cast < void > ( this ) ; // suppress clang warning
2023-05-07 14:56:09 -05:00
return getValue ( value ) ;
} else {
2023-05-12 03:02:05 -05:00
return this - > extendEval ( value ) ;
2023-05-07 14:56:09 -05:00
}
} ;
2023-06-27 06:59:10 -05:00
2024-10-15 07:11:34 -05:00
WellInterface < TypeTag > : : getMobility ( simulator , local_perf_index , mob , obtain , deferred_logger ) ;
2023-10-09 04:14:48 -05:00
2023-06-21 04:48:57 -05:00
if ( this - > isInjector ( ) & & this - > well_ecl_ . getInjMultMode ( ) ! = Well : : InjMultMode : : NONE ) {
2024-10-15 07:11:34 -05:00
const auto perf_ecl_index = this - > perforationData ( ) [ local_perf_index ] . ecl_index ;
2023-06-27 06:59:10 -05:00
const Connection & con = this - > well_ecl_ . getConnections ( ) [ perf_ecl_index ] ;
const int seg = this - > segmentNumberToIndex ( con . segment ( ) ) ;
2023-06-21 04:48:57 -05:00
// from the reference results, it looks like MSW uses segment pressure instead of BHP here
// Note: this is against the documented definition.
// we can change this depending on what we want
2024-02-20 15:32:18 -06:00
const Scalar segment_pres = this - > primary_variables_ . getSegmentPressure ( seg ) . value ( ) ;
const Scalar perf_seg_press_diff = this - > gravity ( ) * this - > segments_ . density ( seg ) . value ( )
2024-10-15 07:11:34 -05:00
* this - > segments_ . perforation_depth_diff ( local_perf_index ) ;
2024-02-20 15:32:18 -06:00
const Scalar perf_press = segment_pres + perf_seg_press_diff ;
2024-10-15 07:11:34 -05:00
const Scalar multiplier = this - > getInjMult ( local_perf_index , segment_pres , perf_press , deferred_logger ) ;
2023-08-15 02:32:10 -05:00
for ( std : : size_t i = 0 ; i < mob . size ( ) ; + + i ) {
2023-06-21 04:48:57 -05:00
mob [ i ] * = multiplier ;
}
}
2021-09-13 02:36:16 -05:00
}
2017-09-11 10:31:42 -05:00
2019-11-01 04:04:44 -05:00
template < typename TypeTag >
2024-02-20 15:32:18 -06:00
typename MultisegmentWell < TypeTag > : : Scalar
2019-11-01 04:04:44 -05:00
MultisegmentWell < TypeTag > : :
2021-04-26 02:31:29 -05:00
getRefDensity ( ) const
2019-11-01 04:04:44 -05:00
{
2022-12-19 08:52:36 -06:00
return this - > segments_ . getRefDensity ( ) ;
2017-09-12 10:13:30 -05:00
}
2017-09-12 10:32:27 -05:00
2020-11-27 00:57:55 -06:00
template < typename TypeTag >
void
MultisegmentWell < TypeTag > : :
2024-02-17 11:13:46 -06:00
checkOperabilityUnderBHPLimit ( const WellState < Scalar > & /*well_state*/ ,
const Simulator & simulator ,
DeferredLogger & deferred_logger )
2020-11-27 00:57:55 -06:00
{
2024-02-06 04:55:07 -06:00
const auto & summaryState = simulator . vanguard ( ) . summaryState ( ) ;
2024-02-20 15:32:18 -06:00
const Scalar bhp_limit = WellBhpThpCalculator ( * this ) . mostStrictBhpFromBhpLimits ( summaryState ) ;
2020-11-27 00:57:55 -06:00
// Crude but works: default is one atmosphere.
// TODO: a better way to detect whether the BHP is defaulted or not
const bool bhp_limit_not_defaulted = bhp_limit > 1.5 * unit : : barsa ;
if ( bhp_limit_not_defaulted | | ! this - > wellHasTHPConstraints ( summaryState ) ) {
// if the BHP limit is not defaulted or the well does not have a THP limit
// we need to check the BHP limit
2024-02-20 15:32:18 -06:00
Scalar total_ipr_mass_rate = 0.0 ;
2022-06-07 04:30:35 -05:00
for ( unsigned phaseIdx = 0 ; phaseIdx < FluidSystem : : numPhases ; + + phaseIdx )
{
if ( ! FluidSystem : : phaseIsActive ( phaseIdx ) ) {
continue ;
}
const unsigned compIdx = Indices : : canonicalToActiveComponentIndex ( FluidSystem : : solventComponentIndex ( phaseIdx ) ) ;
2024-02-20 15:32:18 -06:00
const Scalar ipr_rate = this - > ipr_a_ [ compIdx ] - this - > ipr_b_ [ compIdx ] * bhp_limit ;
2022-06-07 04:30:35 -05:00
2024-02-20 15:32:18 -06:00
const Scalar rho = FluidSystem : : referenceDensity ( phaseIdx , Base : : pvtRegionIdx ( ) ) ;
2022-05-31 07:31:26 -05:00
total_ipr_mass_rate + = ipr_rate * rho ;
2020-11-30 02:23:16 -06:00
}
2022-05-31 07:31:26 -05:00
if ( ( this - > isProducer ( ) & & total_ipr_mass_rate < 0. ) | | ( this - > isInjector ( ) & & total_ipr_mass_rate > 0. ) ) {
2020-11-30 02:23:16 -06:00
this - > operability_status_ . operable_under_only_bhp_limit = false ;
2020-11-27 00:57:55 -06:00
}
// checking whether running under BHP limit will violate THP limit
if ( this - > operability_status_ . operable_under_only_bhp_limit & & this - > wellHasTHPConstraints ( summaryState ) ) {
// option 1: calculate well rates based on the BHP limit.
// option 2: stick with the above IPR curve
// we use IPR here
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > well_rates_bhp_limit ;
2024-02-06 04:55:07 -06:00
computeWellRatesWithBhp ( simulator , bhp_limit , well_rates_bhp_limit , deferred_logger ) ;
2017-09-12 10:32:27 -05:00
2024-02-20 15:32:18 -06:00
const Scalar thp_limit = this - > getTHPConstraint ( summaryState ) ;
const Scalar thp = WellBhpThpCalculator ( * this ) . calculateThpFromBhp ( well_rates_bhp_limit ,
2022-10-19 02:55:14 -05:00
bhp_limit ,
this - > getRefDensity ( ) ,
2023-12-01 06:21:34 -06:00
this - > wellEcl ( ) . alq_value ( summaryState ) ,
2023-05-09 16:44:05 -05:00
thp_limit ,
2022-10-19 02:55:14 -05:00
deferred_logger ) ;
2021-10-01 06:47:53 -05:00
if ( ( this - > isProducer ( ) & & thp < thp_limit ) | | ( this - > isInjector ( ) & & thp > thp_limit ) ) {
2020-11-27 00:57:55 -06:00
this - > operability_status_ . obey_thp_limit_under_bhp_limit = false ;
}
}
} else {
// defaulted BHP and there is a THP constraint
// default BHP limit is about 1 atm.
// when applied the hydrostatic pressure correction dp,
// most likely we get a negative value (bhp + dp)to search in the VFP table,
// which is not desirable.
// we assume we can operate under defaulted BHP limit and will violate the THP limit
// when operating under defaulted BHP limit.
this - > operability_status_ . operable_under_only_bhp_limit = true ;
this - > operability_status_ . obey_thp_limit_under_bhp_limit = false ;
}
}
template < typename TypeTag >
2018-11-15 07:31:41 -06:00
void
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
updateIPR ( const Simulator & simulator , DeferredLogger & deferred_logger ) const
2018-11-15 07:31:41 -06:00
{
2020-11-27 00:57:55 -06:00
// TODO: not handling solvent related here for now
2019-12-13 04:00:01 -06:00
2020-11-27 00:57:55 -06:00
// initialize all the values to be zero to begin with
2021-09-06 05:58:16 -05:00
std : : fill ( this - > ipr_a_ . begin ( ) , this - > ipr_a_ . end ( ) , 0. ) ;
2021-09-06 05:58:16 -05:00
std : : fill ( this - > ipr_b_ . begin ( ) , this - > ipr_b_ . end ( ) , 0. ) ;
2020-11-27 00:57:55 -06:00
2021-06-03 08:34:14 -05:00
const int nseg = this - > numberOfSegments ( ) ;
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > seg_dp ( nseg , 0.0 ) ;
2020-11-27 00:57:55 -06:00
for ( int seg = 0 ; seg < nseg ; + + seg ) {
// calculating the perforation rate for each perforation that belongs to this segment
2024-02-20 15:32:18 -06:00
const Scalar dp = this - > getSegmentDp ( seg ,
2023-01-12 07:23:22 -06:00
this - > segments_ . density ( seg ) . value ( ) ,
seg_dp ) ;
2022-06-07 07:14:59 -05:00
seg_dp [ seg ] = dp ;
2022-12-19 08:52:36 -06:00
for ( const int perf : this - > segments_ . perforations ( ) [ seg ] ) {
2024-10-15 07:11:34 -05:00
const int local_perf_index = this - > pw_info_ . globalToLocal ( perf ) ;
if ( local_perf_index < 0 ) // then the perforation is not on this process
continue ;
2022-12-19 06:52:21 -06:00
std : : vector < Scalar > mob ( this - > num_components_ , 0.0 ) ;
2020-11-27 00:57:55 -06:00
2022-12-19 06:52:21 -06:00
// TODO: maybe we should store the mobility somewhere, so that we only need to calculate it one per iteration
2024-10-15 07:11:34 -05:00
getMobility ( simulator , local_perf_index , mob , deferred_logger ) ;
2020-11-27 00:57:55 -06:00
2024-10-15 07:11:34 -05:00
const int cell_idx = this - > well_cells_ [ local_perf_index ] ;
2024-02-06 04:55:07 -06:00
const auto & int_quantities = simulator . model ( ) . intensiveQuantities ( cell_idx , /*timeIdx=*/ 0 ) ;
2022-12-19 06:52:21 -06:00
const auto & fs = int_quantities . fluidState ( ) ;
// pressure difference between the segment and the perforation
2024-02-20 15:32:18 -06:00
const Scalar perf_seg_press_diff = this - > segments_ . getPressureDiffSegPerf ( seg , perf ) ;
2022-12-19 06:52:21 -06:00
// pressure difference between the perforation and the grid cell
2024-10-15 07:11:34 -05:00
const Scalar cell_perf_press_diff = this - > cell_perforation_pressure_diffs_ [ local_perf_index ] ;
2024-02-20 15:32:18 -06:00
const Scalar pressure_cell = this - > getPerfCellPressure ( fs ) . value ( ) ;
2022-12-19 06:52:21 -06:00
// calculating the b for the connection
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > b_perf ( this - > num_components_ ) ;
2023-08-15 02:32:10 -05:00
for ( std : : size_t phase = 0 ; phase < FluidSystem : : numPhases ; + + phase ) {
2022-12-19 06:52:21 -06:00
if ( ! FluidSystem : : phaseIsActive ( phase ) ) {
continue ;
}
const unsigned comp_idx = Indices : : canonicalToActiveComponentIndex ( FluidSystem : : solventComponentIndex ( phase ) ) ;
b_perf [ comp_idx ] = fs . invB ( phase ) . value ( ) ;
2020-11-27 00:57:55 -06:00
}
2022-12-19 06:52:21 -06:00
// the pressure difference between the connection and BHP
2024-02-20 15:32:18 -06:00
const Scalar h_perf = cell_perf_press_diff + perf_seg_press_diff + dp ;
const Scalar pressure_diff = pressure_cell - h_perf ;
2020-11-27 00:57:55 -06:00
2022-12-19 06:52:21 -06:00
// do not take into consideration the crossflow here.
if ( ( this - > isProducer ( ) & & pressure_diff < 0. ) | | ( this - > isInjector ( ) & & pressure_diff > 0. ) ) {
deferred_logger . debug ( " CROSSFLOW_IPR " ,
" cross flow found when updateIPR for well " + this - > name ( ) ) ;
}
2020-11-27 00:57:55 -06:00
2022-12-19 06:52:21 -06:00
// the well index associated with the connection
2024-02-20 15:32:18 -06:00
const Scalar trans_mult = simulator . problem ( ) . template wellTransMultiplier < Scalar > ( int_quantities , cell_idx ) ;
2024-02-06 04:55:07 -06:00
const auto & wellstate_nupcol = simulator . problem ( ) . wellModel ( ) . nupcolWellState ( ) . well ( this - > index_of_well_ ) ;
2024-10-15 07:11:34 -05:00
const std : : vector < Scalar > tw_perf = this - > wellIndex ( local_perf_index , int_quantities , trans_mult , wellstate_nupcol ) ;
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > ipr_a_perf ( this - > ipr_a_ . size ( ) ) ;
std : : vector < Scalar > ipr_b_perf ( this - > ipr_b_ . size ( ) ) ;
2022-12-19 06:52:21 -06:00
for ( int comp_idx = 0 ; comp_idx < this - > num_components_ ; + + comp_idx ) {
2024-02-20 15:32:18 -06:00
const Scalar tw_mob = tw_perf [ comp_idx ] * mob [ comp_idx ] * b_perf [ comp_idx ] ;
2022-12-19 06:52:21 -06:00
ipr_a_perf [ comp_idx ] + = tw_mob * pressure_diff ;
ipr_b_perf [ comp_idx ] + = tw_mob ;
}
2020-11-27 00:57:55 -06:00
2022-12-19 06:52:21 -06:00
// we need to handle the rs and rv when both oil and gas are present
if ( FluidSystem : : phaseIsActive ( FluidSystem : : oilPhaseIdx ) & & FluidSystem : : phaseIsActive ( FluidSystem : : gasPhaseIdx ) ) {
const unsigned oil_comp_idx = Indices : : canonicalToActiveComponentIndex ( FluidSystem : : oilCompIdx ) ;
const unsigned gas_comp_idx = Indices : : canonicalToActiveComponentIndex ( FluidSystem : : gasCompIdx ) ;
2024-02-20 15:32:18 -06:00
const Scalar rs = ( fs . Rs ( ) ) . value ( ) ;
const Scalar rv = ( fs . Rv ( ) ) . value ( ) ;
2020-11-27 00:57:55 -06:00
2024-02-20 15:32:18 -06:00
const Scalar dis_gas_a = rs * ipr_a_perf [ oil_comp_idx ] ;
const Scalar vap_oil_a = rv * ipr_a_perf [ gas_comp_idx ] ;
2020-11-27 00:57:55 -06:00
2022-12-19 06:52:21 -06:00
ipr_a_perf [ gas_comp_idx ] + = dis_gas_a ;
ipr_a_perf [ oil_comp_idx ] + = vap_oil_a ;
2020-11-27 00:57:55 -06:00
2024-02-20 15:32:18 -06:00
const Scalar dis_gas_b = rs * ipr_b_perf [ oil_comp_idx ] ;
const Scalar vap_oil_b = rv * ipr_b_perf [ gas_comp_idx ] ;
2020-11-27 00:57:55 -06:00
2022-12-19 06:52:21 -06:00
ipr_b_perf [ gas_comp_idx ] + = dis_gas_b ;
ipr_b_perf [ oil_comp_idx ] + = vap_oil_b ;
}
2020-11-27 00:57:55 -06:00
2023-08-15 02:32:10 -05:00
for ( std : : size_t comp_idx = 0 ; comp_idx < ipr_a_perf . size ( ) ; + + comp_idx ) {
2022-12-19 06:52:21 -06:00
this - > ipr_a_ [ comp_idx ] + = ipr_a_perf [ comp_idx ] ;
this - > ipr_b_ [ comp_idx ] + = ipr_b_perf [ comp_idx ] ;
}
2020-11-27 00:57:55 -06:00
}
2019-12-13 04:00:01 -06:00
}
2024-09-30 02:05:25 -05:00
this - > parallel_well_info_ . communication ( ) . sum ( this - > ipr_a_ . data ( ) , this - > ipr_a_ . size ( ) ) ;
this - > parallel_well_info_ . communication ( ) . sum ( this - > ipr_b_ . data ( ) , this - > ipr_b_ . size ( ) ) ;
2020-11-27 00:57:55 -06:00
}
2019-12-13 04:00:01 -06:00
2023-10-26 10:28:05 -05:00
template < typename TypeTag >
void
MultisegmentWell < TypeTag > : :
2024-02-17 11:13:46 -06:00
updateIPRImplicit ( const Simulator & simulator ,
WellState < Scalar > & well_state ,
DeferredLogger & deferred_logger )
2023-10-26 10:28:05 -05:00
{
2023-11-07 05:31:26 -06:00
// Compute IPR based on *converged* well-equation:
// For a component rate r the derivative dr/dbhp is obtained by
2023-12-07 05:15:40 -06:00
// dr/dbhp = - (partial r/partial x) * inv(partial Eq/partial x) * (partial Eq/partial bhp_target)
// where Eq(x)=0 is the well equation setup with bhp control and primary variables x
2023-11-07 05:31:26 -06:00
// We shouldn't have zero rates at this stage, but check
bool zero_rates ;
auto rates = well_state . well ( this - > index_of_well_ ) . surface_rates ;
zero_rates = true ;
for ( std : : size_t p = 0 ; p < rates . size ( ) ; + + p ) {
zero_rates & = rates [ p ] = = 0.0 ;
}
auto & ws = well_state . well ( this - > index_of_well_ ) ;
if ( zero_rates ) {
const auto msg = fmt : : format ( " updateIPRImplicit: Well {} has zero rate, IPRs might be problematic " , this - > name ( ) ) ;
deferred_logger . debug ( msg ) ;
/*
2023-11-08 14:58:26 -06:00
// could revert to standard approach here:
2024-02-06 04:55:07 -06:00
updateIPR ( simulator , deferred_logger ) ;
2023-11-07 05:31:26 -06:00
for ( int comp_idx = 0 ; comp_idx < this - > num_components_ ; + + comp_idx ) {
2024-02-06 04:55:07 -06:00
const int idx = this - > modelCompIdxToFlowCompIdx ( comp_idx ) ;
2023-11-07 06:31:12 -06:00
ws . implicit_ipr_a [ idx ] = this - > ipr_a_ [ comp_idx ] ;
ws . implicit_ipr_b [ idx ] = this - > ipr_b_ [ comp_idx ] ;
2023-11-07 05:31:26 -06:00
}
return ;
*/
}
2024-02-06 04:55:07 -06:00
const auto & group_state = simulator . problem ( ) . wellModel ( ) . groupState ( ) ;
2023-11-07 05:31:26 -06:00
2023-11-07 06:31:12 -06:00
std : : fill ( ws . implicit_ipr_a . begin ( ) , ws . implicit_ipr_a . end ( ) , 0. ) ;
std : : fill ( ws . implicit_ipr_b . begin ( ) , ws . implicit_ipr_b . end ( ) , 0. ) ;
2023-11-07 05:31:26 -06:00
//WellState well_state_copy = well_state;
auto inj_controls = Well : : InjectionControls ( 0 ) ;
auto prod_controls = Well : : ProductionControls ( 0 ) ;
prod_controls . addControl ( Well : : ProducerCMode : : BHP ) ;
prod_controls . bhp_limit = well_state . well ( this - > index_of_well_ ) . bhp ;
// Set current control to bhp, and bhp value in state, modify bhp limit in control object.
const auto cmode = ws . production_cmode ;
ws . production_cmode = Well : : ProducerCMode : : BHP ;
2024-02-06 04:55:07 -06:00
const double dt = simulator . timeStepSize ( ) ;
assembleWellEqWithoutIteration ( simulator , dt , inj_controls , prod_controls , well_state , group_state , deferred_logger ) ;
2023-11-07 05:31:26 -06:00
BVectorWell rhs ( this - > numberOfSegments ( ) ) ;
rhs = 0.0 ;
rhs [ 0 ] [ SPres ] = - 1.0 ;
const BVectorWell x_well = this - > linSys_ . solve ( rhs ) ;
constexpr int num_eq = MSWEval : : numWellEq ;
for ( int comp_idx = 0 ; comp_idx < this - > num_components_ ; + + comp_idx ) {
const EvalWell comp_rate = this - > primary_variables_ . getQs ( comp_idx ) ;
2024-02-06 05:12:13 -06:00
const int idx = this - > modelCompIdxToFlowCompIdx ( comp_idx ) ;
2023-11-07 05:31:26 -06:00
for ( size_t pvIdx = 0 ; pvIdx < num_eq ; + + pvIdx ) {
2023-12-07 05:15:40 -06:00
// well primary variable derivatives in EvalWell start at position Indices::numEq
2023-11-07 06:31:12 -06:00
ws . implicit_ipr_b [ idx ] - = x_well [ 0 ] [ pvIdx ] * comp_rate . derivative ( pvIdx + Indices : : numEq ) ;
2023-11-07 05:31:26 -06:00
}
2023-11-07 06:31:12 -06:00
ws . implicit_ipr_a [ idx ] = ws . implicit_ipr_b [ idx ] * ws . bhp - comp_rate . value ( ) ;
2023-11-07 05:31:26 -06:00
}
// reset cmode
ws . production_cmode = cmode ;
2023-10-26 10:28:05 -05:00
}
2020-11-27 00:57:55 -06:00
template < typename TypeTag >
void
MultisegmentWell < TypeTag > : :
2024-02-17 11:13:46 -06:00
checkOperabilityUnderTHPLimit ( const Simulator & simulator ,
const WellState < Scalar > & well_state ,
DeferredLogger & deferred_logger )
2020-11-27 00:57:55 -06:00
{
2024-02-06 04:55:07 -06:00
const auto & summaryState = simulator . vanguard ( ) . summaryState ( ) ;
2022-02-07 04:28:35 -06:00
const auto obtain_bhp = this - > isProducer ( )
? computeBhpAtThpLimitProd (
2024-02-06 04:55:07 -06:00
well_state , simulator , summaryState , deferred_logger )
: computeBhpAtThpLimitInj ( simulator , summaryState , deferred_logger ) ;
2020-11-27 00:57:55 -06:00
if ( obtain_bhp ) {
this - > operability_status_ . can_obtain_bhp_with_thp_limit = true ;
2024-02-20 15:32:18 -06:00
const Scalar bhp_limit = WellBhpThpCalculator ( * this ) . mostStrictBhpFromBhpLimits ( summaryState ) ;
2020-11-27 00:57:55 -06:00
this - > operability_status_ . obey_bhp_limit_with_thp_limit = ( * obtain_bhp > = bhp_limit ) ;
2024-02-20 15:32:18 -06:00
const Scalar thp_limit = this - > getTHPConstraint ( summaryState ) ;
2021-10-01 06:47:53 -05:00
if ( this - > isProducer ( ) & & * obtain_bhp < thp_limit ) {
2020-11-27 00:57:55 -06:00
const std : : string msg = " obtained bhp " + std : : to_string ( unit : : convert : : to ( * obtain_bhp , unit : : barsa ) )
2021-06-03 08:34:14 -05:00
+ " bars is SMALLER than thp limit "
+ std : : to_string ( unit : : convert : : to ( thp_limit , unit : : barsa ) )
2021-09-06 05:58:16 -05:00
+ " bars as a producer for well " + this - > name ( ) ;
2021-06-03 08:34:14 -05:00
deferred_logger . debug ( msg ) ;
2017-09-12 10:32:27 -05:00
}
2021-10-01 06:47:53 -05:00
else if ( this - > isInjector ( ) & & * obtain_bhp > thp_limit ) {
const std : : string msg = " obtained bhp " + std : : to_string ( unit : : convert : : to ( * obtain_bhp , unit : : barsa ) )
+ " bars is LARGER than thp limit "
+ std : : to_string ( unit : : convert : : to ( thp_limit , unit : : barsa ) )
+ " bars as a injector for well " + this - > name ( ) ;
deferred_logger . debug ( msg ) ;
}
2021-06-03 08:34:14 -05:00
} else {
// Shutting wells that can not find bhp value from thp
// when under THP control
this - > operability_status_ . can_obtain_bhp_with_thp_limit = false ;
this - > operability_status_ . obey_bhp_limit_with_thp_limit = false ;
if ( ! this - > wellIsStopped ( ) ) {
2024-02-20 15:32:18 -06:00
const Scalar thp_limit = this - > getTHPConstraint ( summaryState ) ;
2021-06-03 08:34:14 -05:00
deferred_logger . debug ( " could not find bhp value at thp limit "
+ std : : to_string ( unit : : convert : : to ( thp_limit , unit : : barsa ) )
2021-09-06 05:58:16 -05:00
+ " bar for well " + this - > name ( ) + " , the well might need to be closed " ) ;
2017-09-12 10:32:27 -05:00
}
2019-11-29 10:23:26 -06:00
}
2017-09-12 10:32:27 -05:00
}
2017-09-26 03:07:06 -05:00
2021-06-03 08:34:14 -05:00
2017-10-10 07:30:23 -05:00
template < typename TypeTag >
2020-06-02 07:54:45 -05:00
bool
2017-10-10 07:30:23 -05:00
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
iterateWellEqWithControl ( const Simulator & simulator ,
2020-06-02 07:54:45 -05:00
const double dt ,
const Well : : InjectionControls & inj_controls ,
const Well : : ProductionControls & prod_controls ,
2024-02-17 11:13:46 -06:00
WellState < Scalar > & well_state ,
2024-02-17 11:13:46 -06:00
const GroupState < Scalar > & group_state ,
2021-05-05 04:22:44 -05:00
DeferredLogger & deferred_logger )
2017-10-10 07:30:23 -05:00
{
2021-09-29 09:01:16 -05:00
if ( ! this - > isOperableAndSolvable ( ) & & ! this - > wellIsStopped ( ) ) return true ;
2020-11-27 00:57:55 -06:00
2021-09-06 05:58:16 -05:00
const int max_iter_number = this - > param_ . max_inner_iter_ms_wells_ ;
2022-04-08 04:07:34 -05:00
{
// getWellFiniteResiduals returns false for nan/inf residuals
const auto & [ isFinite , residuals ] = this - > getFiniteWellResiduals ( Base : : B_avg_ , deferred_logger ) ;
if ( ! isFinite )
return false ;
}
2019-03-08 06:50:34 -06:00
std : : vector < std : : vector < Scalar > > residual_history ;
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > measure_history ;
2017-10-10 07:30:23 -05:00
int it = 0 ;
2019-03-08 06:50:34 -06:00
// relaxation factor
2024-02-20 15:32:18 -06:00
Scalar relaxation_factor = 1. ;
const Scalar min_relaxation_factor = 0.6 ;
2019-11-18 09:33:10 -06:00
bool converged = false ;
int stagnate_count = 0 ;
2020-05-12 05:05:30 -05:00
bool relax_convergence = false ;
2022-05-03 01:25:26 -05:00
this - > regularize_ = false ;
2019-11-18 09:33:10 -06:00
for ( ; it < max_iter_number ; + + it , + + debug_cost_counter_ ) {
2017-10-10 07:30:23 -05:00
2024-02-06 04:55:07 -06:00
assembleWellEqWithoutIteration ( simulator , dt , inj_controls , prod_controls , well_state , group_state , deferred_logger ) ;
2017-10-10 07:30:23 -05:00
2023-12-04 16:15:08 -06:00
BVectorWell dx_well ;
try {
dx_well = this - > linSys_ . solve ( ) ;
}
catch ( const NumericalProblem & exp ) {
// Add information about the well and log to deferred logger
// (Logging done inside of solve() method will only be seen if
// this is the process with rank zero)
deferred_logger . problem ( " In MultisegmentWell::iterateWellEqWithControl for well "
+ this - > name ( ) + " : " + exp . what ( ) ) ;
throw ;
}
2017-10-10 07:30:23 -05:00
2022-04-06 03:09:21 -05:00
if ( it > this - > param_ . strict_inner_iter_wells_ ) {
2020-05-12 05:05:30 -05:00
relax_convergence = true ;
2022-04-06 03:09:21 -05:00
this - > regularize_ = true ;
}
2017-10-10 07:30:23 -05:00
2024-03-19 06:50:34 -05:00
const auto report = getWellConvergence ( simulator , well_state , Base : : B_avg_ , deferred_logger , relax_convergence ) ;
2018-10-23 07:05:19 -05:00
if ( report . converged ( ) ) {
2019-11-18 09:33:10 -06:00
converged = true ;
2017-10-10 07:30:23 -05:00
break ;
}
2022-04-08 04:07:34 -05:00
{
// getFinteWellResiduals returns false for nan/inf residuals
const auto & [ isFinite , residuals ] = this - > getFiniteWellResiduals ( Base : : B_avg_ , deferred_logger ) ;
if ( ! isFinite )
return false ;
residual_history . push_back ( residuals ) ;
measure_history . push_back ( this - > getResidualMeasureValue ( well_state ,
2021-06-03 08:34:14 -05:00
residual_history [ it ] ,
2021-09-06 05:58:16 -05:00
this - > param_ . tolerance_wells_ ,
this - > param_ . tolerance_pressure_ms_wells_ ,
2021-06-03 08:34:14 -05:00
deferred_logger ) ) ;
2022-04-08 04:07:34 -05:00
}
2019-03-08 06:50:34 -06:00
bool is_oscillate = false ;
bool is_stagnate = false ;
2023-09-25 06:07:39 -05:00
this - > detectOscillations ( measure_history , is_oscillate , is_stagnate ) ;
2022-11-05 16:51:59 -05:00
// TODO: maybe we should have more sophisticated strategy to recover the relaxation factor,
2019-03-08 06:50:34 -06:00
// for example, to recover it to be bigger
if ( is_oscillate | | is_stagnate ) {
2019-11-18 09:33:10 -06:00
// HACK!
2020-05-12 05:05:30 -05:00
std : : ostringstream sstr ;
if ( relaxation_factor = = min_relaxation_factor ) {
2019-11-18 09:33:10 -06:00
// Still stagnating, terminate iterations if 5 iterations pass.
+ + stagnate_count ;
2020-05-12 05:05:30 -05:00
if ( stagnate_count = = 6 ) {
2021-09-06 05:58:16 -05:00
sstr < < " well " < < this - > name ( ) < < " observes severe stagnation and/or oscillation. We relax the tolerance and check for convergence. \n " ;
2024-03-19 06:50:34 -05:00
const auto reportStag = getWellConvergence ( simulator , well_state , Base : : B_avg_ , deferred_logger , true ) ;
2020-05-12 05:05:30 -05:00
if ( reportStag . converged ( ) ) {
converged = true ;
2021-09-06 05:58:16 -05:00
sstr < < " well " < < this - > name ( ) < < " manages to get converged with relaxed tolerances in " < < it < < " inner iterations " ;
2020-05-12 05:05:30 -05:00
deferred_logger . debug ( sstr . str ( ) ) ;
2020-06-02 07:54:45 -05:00
return converged ;
2020-05-12 05:05:30 -05:00
}
2019-11-18 09:33:10 -06:00
}
}
2019-03-08 06:50:34 -06:00
// a factor value to reduce the relaxation_factor
2024-02-20 15:32:18 -06:00
const Scalar reduction_mutliplier = 0.9 ;
2019-03-08 06:50:34 -06:00
relaxation_factor = std : : max ( relaxation_factor * reduction_mutliplier , min_relaxation_factor ) ;
// debug output
if ( is_stagnate ) {
2021-09-06 05:58:16 -05:00
sstr < < " well " < < this - > name ( ) < < " observes stagnation in inner iteration " < < it < < " \n " ;
2019-11-29 10:23:26 -06:00
2019-03-08 06:50:34 -06:00
}
if ( is_oscillate ) {
2021-09-06 05:58:16 -05:00
sstr < < " well " < < this - > name ( ) < < " observes oscillation in inner iteration " < < it < < " \n " ;
2019-03-08 06:50:34 -06:00
}
sstr < < " relaxation_factor is " < < relaxation_factor < < " now \n " ;
2022-04-06 03:09:21 -05:00
this - > regularize_ = true ;
2019-03-08 06:50:34 -06:00
deferred_logger . debug ( sstr . str ( ) ) ;
}
2024-03-19 06:50:34 -05:00
updateWellState ( simulator , dx_well , well_state , deferred_logger , relaxation_factor ) ;
2017-10-10 07:30:23 -05:00
initPrimaryVariablesEvaluation ( ) ;
}
2019-03-08 06:50:34 -06:00
// TODO: we should decide whether to keep the updated well_state, or recover to use the old well_state
2019-11-18 09:33:10 -06:00
if ( converged ) {
2019-03-08 06:50:34 -06:00
std : : ostringstream sstr ;
2021-09-06 05:58:16 -05:00
sstr < < " Well " < < this - > name ( ) < < " converged in " < < it < < " inner iterations. " ;
2020-05-12 05:05:30 -05:00
if ( relax_convergence )
2021-09-06 03:08:28 -05:00
sstr < < " (A relaxed tolerance was used after " < < this - > param_ . strict_inner_iter_wells_ < < " iterations) " ;
2019-03-08 06:50:34 -06:00
deferred_logger . debug ( sstr . str ( ) ) ;
} else {
std : : ostringstream sstr ;
2021-09-06 05:58:16 -05:00
sstr < < " Well " < < this - > name ( ) < < " did not converge in " < < it < < " inner iterations. " ;
2020-06-10 09:36:17 -05:00
# define EXTRA_DEBUG_MSW 0
# if EXTRA_DEBUG_MSW
2021-09-06 05:58:16 -05:00
sstr < < " ***** Outputting the residual history for well " < < this - > name ( ) < < " during inner iterations: " ;
2019-03-08 06:50:34 -06:00
for ( int i = 0 ; i < it ; + + i ) {
const auto & residual = residual_history [ i ] ;
sstr < < " residual at " < < i < < " th iteration " ;
for ( const auto & res : residual ) {
sstr < < " " < < res ;
}
sstr < < " " < < measure_history [ i ] < < " \n " ;
}
2020-06-10 09:36:17 -05:00
# endif
2022-11-14 04:36:48 -06:00
# undef EXTRA_DEBUG_MSW
2019-03-08 06:50:34 -06:00
deferred_logger . debug ( sstr . str ( ) ) ;
}
2020-06-02 07:54:45 -05:00
return converged ;
2017-10-10 07:30:23 -05:00
}
2023-08-16 05:50:06 -05:00
template < typename TypeTag >
bool
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
iterateWellEqWithSwitching ( const Simulator & simulator ,
2023-08-16 05:50:06 -05:00
const double dt ,
const Well : : InjectionControls & inj_controls ,
const Well : : ProductionControls & prod_controls ,
2024-02-17 11:13:46 -06:00
WellState < Scalar > & well_state ,
2024-02-17 11:13:46 -06:00
const GroupState < Scalar > & group_state ,
2023-10-26 10:28:05 -05:00
DeferredLogger & deferred_logger ,
2023-10-30 14:58:09 -05:00
const bool fixed_control /*false*/ ,
const bool fixed_status /*false*/ )
2023-08-16 05:50:06 -05:00
{
const int max_iter_number = this - > param_ . max_inner_iter_ms_wells_ ;
{
// getWellFiniteResiduals returns false for nan/inf residuals
const auto & [ isFinite , residuals ] = this - > getFiniteWellResiduals ( Base : : B_avg_ , deferred_logger ) ;
if ( ! isFinite )
return false ;
}
std : : vector < std : : vector < Scalar > > residual_history ;
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > measure_history ;
2023-08-16 05:50:06 -05:00
int it = 0 ;
// relaxation factor
2024-02-20 15:32:18 -06:00
Scalar relaxation_factor = 1. ;
const Scalar min_relaxation_factor = 0.6 ;
2023-08-16 05:50:06 -05:00
bool converged = false ;
2023-10-05 02:31:36 -05:00
[[maybe_unused]] int stagnate_count = 0 ;
2023-08-16 05:50:06 -05:00
bool relax_convergence = false ;
this - > regularize_ = false ;
2024-02-06 04:55:07 -06:00
const auto & summary_state = simulator . vanguard ( ) . summaryState ( ) ;
2017-10-10 07:30:23 -05:00
2023-12-07 05:15:40 -06:00
// Always take a few (more than one) iterations after a switch before allowing a new switch
// The optimal number here is subject to further investigation, but it has been observerved
// that unless this number is >1, we may get stuck in a cycle
2023-11-08 14:07:47 -06:00
const int min_its_after_switch = 3 ;
2023-08-16 05:50:06 -05:00
int its_since_last_switch = min_its_after_switch ;
int switch_count = 0 ;
2023-11-10 09:03:48 -06:00
// if we fail to solve eqs, we reset status/operability before leaving
2023-11-08 14:07:47 -06:00
const auto well_status_orig = this - > wellStatus_ ;
2023-11-10 09:03:48 -06:00
const auto operability_orig = this - > operability_status_ ;
2023-11-08 14:07:47 -06:00
auto well_status_cur = well_status_orig ;
2024-01-26 09:24:09 -06:00
// don't allow opening wells that are stopped from schedule or has a stopped well state
const bool allow_open = this - > well_ecl_ . getStatus ( ) = = WellStatus : : OPEN & &
well_state . well ( this - > index_of_well_ ) . status = = WellStatus : : OPEN ;
// don't allow switcing for wells under zero rate target or requested fixed status and control
2024-03-19 06:50:34 -05:00
const bool allow_switching = ! this - > wellUnderZeroRateTarget ( simulator , well_state , deferred_logger ) & &
2024-01-26 09:24:09 -06:00
( ! fixed_control | | ! fixed_status ) & & allow_open ;
2023-08-16 05:50:06 -05:00
bool changed = false ;
2023-10-09 04:14:48 -05:00
bool final_check = false ;
2023-11-10 09:03:48 -06:00
// well needs to be set operable or else solving/updating of re-opened wells is skipped
this - > operability_status_ . resetOperability ( ) ;
this - > operability_status_ . solvable = true ;
2023-08-16 05:50:06 -05:00
for ( ; it < max_iter_number ; + + it , + + debug_cost_counter_ ) {
its_since_last_switch + + ;
2023-11-08 14:07:47 -06:00
if ( allow_switching & & its_since_last_switch > = min_its_after_switch ) {
2024-02-20 15:32:18 -06:00
const Scalar wqTotal = this - > primary_variables_ . getWQTotal ( ) . value ( ) ;
changed = this - > updateWellControlAndStatusLocalIteration ( simulator , well_state , group_state ,
inj_controls , prod_controls , wqTotal ,
deferred_logger , fixed_control , fixed_status ) ;
2023-08-16 05:50:06 -05:00
if ( changed ) {
its_since_last_switch = 0 ;
switch_count + + ;
2023-11-08 14:07:47 -06:00
if ( well_status_cur ! = this - > wellStatus_ ) {
well_status_cur = this - > wellStatus_ ;
}
2023-08-16 05:50:06 -05:00
}
if ( ! changed & & final_check ) {
break ;
} else {
final_check = false ;
}
}
2024-02-06 04:55:07 -06:00
assembleWellEqWithoutIteration ( simulator , dt , inj_controls , prod_controls , well_state , group_state , deferred_logger ) ;
2023-08-16 05:50:06 -05:00
const BVectorWell dx_well = this - > linSys_ . solve ( ) ;
if ( it > this - > param_ . strict_inner_iter_wells_ ) {
relax_convergence = true ;
this - > regularize_ = true ;
}
2024-03-19 06:50:34 -05:00
const auto report = getWellConvergence ( simulator , well_state , Base : : B_avg_ , deferred_logger , relax_convergence ) ;
2023-08-16 05:50:06 -05:00
converged = report . converged ( ) ;
2024-09-30 02:06:50 -05:00
if ( this - > parallel_well_info_ . communication ( ) . size ( ) > 1 & & converged ! = this - > parallel_well_info_ . communication ( ) . min ( converged ) ) {
OPM_THROW ( std : : runtime_error , fmt : : format ( " Misalignment of the parallel simulation run in iterateWellEqWithSwitching - the well calculation succeeded on rank {} but failed on other ranks. " , this - > parallel_well_info_ . communication ( ) . rank ( ) ) ) ;
}
2023-08-16 05:50:06 -05:00
if ( converged ) {
// if equations are sufficiently linear they might converge in less than min_its_after_switch
// in this case, make sure all constraints are satisfied before returning
if ( switch_count > 0 & & its_since_last_switch < min_its_after_switch ) {
final_check = true ;
its_since_last_switch = min_its_after_switch ;
} else {
break ;
2023-09-25 06:07:39 -05:00
}
2023-08-16 05:50:06 -05:00
}
2023-09-25 06:07:39 -05:00
// getFinteWellResiduals returns false for nan/inf residuals
2023-08-16 05:50:06 -05:00
{
const auto & [ isFinite , residuals ] = this - > getFiniteWellResiduals ( Base : : B_avg_ , deferred_logger ) ;
if ( ! isFinite )
return false ;
residual_history . push_back ( residuals ) ;
}
2023-09-25 06:07:39 -05:00
if ( ! converged ) {
measure_history . push_back ( this - > getResidualMeasureValue ( well_state ,
residual_history [ it ] ,
this - > param_ . tolerance_wells_ ,
this - > param_ . tolerance_pressure_ms_wells_ ,
deferred_logger ) ) ;
bool is_oscillate = false ;
bool is_stagnate = false ;
this - > detectOscillations ( measure_history , is_oscillate , is_stagnate ) ;
// TODO: maybe we should have more sophisticated strategy to recover the relaxation factor,
// for example, to recover it to be bigger
if ( is_oscillate | | is_stagnate ) {
// HACK!
std : : string message ;
if ( relaxation_factor = = min_relaxation_factor ) {
+ + stagnate_count ;
if ( false ) { // this disables the usage of the relaxed tolerance
fmt : : format_to ( std : : back_inserter ( message ) , " Well {} observes severe stagnation and/or oscillation. "
" We relax the tolerance and check for convergence. \n " , this - > name ( ) ) ;
2024-03-19 06:50:34 -05:00
const auto reportStag = getWellConvergence ( simulator , well_state , Base : : B_avg_ ,
2023-09-25 06:07:39 -05:00
deferred_logger , true ) ;
if ( reportStag . converged ( ) ) {
converged = true ;
fmt : : format_to ( std : : back_inserter ( message ) , " Well {} manages to get converged with relaxed tolerances in {} inner iterations " , this - > name ( ) , it ) ;
deferred_logger . debug ( message ) ;
return converged ;
}
2023-08-16 05:50:06 -05:00
}
}
2023-09-25 06:07:39 -05:00
// a factor value to reduce the relaxation_factor
2024-02-20 15:32:18 -06:00
constexpr Scalar reduction_mutliplier = 0.9 ;
2023-09-25 06:07:39 -05:00
relaxation_factor = std : : max ( relaxation_factor * reduction_mutliplier , min_relaxation_factor ) ;
2023-08-16 05:50:06 -05:00
2023-09-25 06:07:39 -05:00
// debug output
if ( is_stagnate ) {
fmt : : format_to ( std : : back_inserter ( message ) , " well {} observes stagnation in inner iteration {} \n " , this - > name ( ) , it ) ;
}
if ( is_oscillate ) {
fmt : : format_to ( std : : back_inserter ( message ) , " well {} observes oscillation in inner iteration {} \n " , this - > name ( ) , it ) ;
}
fmt : : format_to ( std : : back_inserter ( message ) , " relaxation_factor is {} now \n " , relaxation_factor ) ;
2023-08-16 05:50:06 -05:00
2023-09-25 06:07:39 -05:00
this - > regularize_ = true ;
deferred_logger . debug ( message ) ;
2023-08-16 05:50:06 -05:00
}
}
2024-03-19 06:50:34 -05:00
updateWellState ( simulator , dx_well , well_state , deferred_logger , relaxation_factor ) ;
2023-08-16 05:50:06 -05:00
initPrimaryVariablesEvaluation ( ) ;
}
if ( converged ) {
if ( allow_switching ) {
// update operability if status change
const bool is_stopped = this - > wellIsStopped ( ) ;
if ( this - > wellHasTHPConstraints ( summary_state ) ) {
this - > operability_status_ . can_obtain_bhp_with_thp_limit = ! is_stopped ;
this - > operability_status_ . obey_thp_limit_under_bhp_limit = ! is_stopped ;
} else {
this - > operability_status_ . operable_under_only_bhp_limit = ! is_stopped ;
}
2023-09-25 06:07:39 -05:00
}
std : : string message = fmt : : format ( " Well {} converged in {} inner iterations ( "
" {} control/status switches). " , this - > name ( ) , it , switch_count ) ;
if ( relax_convergence ) {
message . append ( fmt : : format ( " (A relaxed tolerance was used after {} iterations) " ,
this - > param_ . strict_inner_iter_wells_ ) ) ;
}
deferred_logger . debug ( message ) ;
2023-08-16 05:50:06 -05:00
} else {
2023-11-08 14:07:47 -06:00
this - > wellStatus_ = well_status_orig ;
2023-11-10 09:03:48 -06:00
this - > operability_status_ = operability_orig ;
2023-12-07 05:15:40 -06:00
const std : : string message = fmt : : format ( " Well {} did not converge in {} inner iterations ( "
2023-09-25 06:07:39 -05:00
" {} control/status switches). " , this - > name ( ) , it , switch_count ) ;
deferred_logger . debug ( message ) ;
2024-03-02 09:21:53 -06:00
this - > primary_variables_ . outputLowLimitPressureSegments ( deferred_logger ) ;
2023-08-16 05:50:06 -05:00
}
return converged ;
}
2017-10-10 07:30:23 -05:00
template < typename TypeTag >
void
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
assembleWellEqWithoutIteration ( const Simulator & simulator ,
2017-10-10 07:30:23 -05:00
const double dt ,
2019-12-05 05:40:35 -06:00
const Well : : InjectionControls & inj_controls ,
const Well : : ProductionControls & prod_controls ,
2024-02-17 11:13:46 -06:00
WellState < Scalar > & well_state ,
2024-02-17 11:13:46 -06:00
const GroupState < Scalar > & group_state ,
2021-05-05 04:22:44 -05:00
DeferredLogger & deferred_logger )
2017-10-10 07:30:23 -05:00
{
2021-09-29 09:01:16 -05:00
if ( ! this - > isOperableAndSolvable ( ) & & ! this - > wellIsStopped ( ) ) return ;
2020-11-27 00:57:55 -06:00
2019-04-12 08:43:39 -05:00
// update the upwinding segments
2022-12-19 08:43:42 -06:00
this - > segments_ . updateUpwindingSegments ( this - > primary_variables_ ) ;
2019-04-12 08:43:39 -05:00
2020-06-08 06:53:13 -05:00
// calculate the fluid properties needed.
2024-02-06 04:55:07 -06:00
computeSegmentFluidProperties ( simulator , deferred_logger ) ;
2020-06-08 06:53:13 -05:00
2017-10-10 07:30:23 -05:00
// clear all entries
2022-11-11 14:41:24 -06:00
this - > linSys_ . clear ( ) ;
2020-08-18 10:54:15 -05:00
2021-08-05 03:14:02 -05:00
auto & ws = well_state . well ( this - > index_of_well_ ) ;
2023-05-04 06:24:23 -05:00
ws . phase_mixing_rates . fill ( 0.0 ) ;
2019-04-23 06:30:12 -05:00
2017-10-10 07:30:23 -05:00
// for the black oil cases, there will be four equations,
// the first three of them are the mass balance equations, the last one is the pressure equations.
//
// but for the top segment, the pressure equation will be the well control equation, and the other three will be the same.
2024-02-06 04:55:07 -06:00
const bool allow_cf = this - > getAllowCrossFlow ( ) | | openCrossFlowAvoidSingularity ( simulator ) ;
2017-10-10 07:30:23 -05:00
2021-06-03 08:34:14 -05:00
const int nseg = this - > numberOfSegments ( ) ;
2017-10-10 07:30:23 -05:00
for ( int seg = 0 ; seg < nseg ; + + seg ) {
2019-02-27 07:47:31 -06:00
// calculating the accumulation term
2022-11-05 16:51:59 -05:00
// TODO: without considering the efficiency factor for now
2017-10-10 07:30:23 -05:00
{
2024-02-06 04:55:07 -06:00
const EvalWell segment_surface_volume = getSegmentSurfaceVolume ( simulator , seg ) ;
2020-05-20 09:36:05 -05:00
// Add a regularization_factor to increase the accumulation term
2020-09-30 03:04:39 -05:00
// This will make the system less stiff and help convergence for
2020-05-20 09:36:05 -05:00
// difficult cases
2022-06-29 05:54:42 -05:00
const Scalar regularization_factor = this - > regularize_ ? this - > param_ . regularization_factor_wells_ : 1.0 ;
2017-10-10 07:30:23 -05:00
// for each component
2021-09-06 05:58:16 -05:00
for ( int comp_idx = 0 ; comp_idx < this - > num_components_ ; + + comp_idx ) {
2022-12-19 02:52:48 -06:00
const EvalWell accumulation_term = regularization_factor * ( segment_surface_volume * this - > primary_variables_ . surfaceVolumeFraction ( seg , comp_idx )
2019-04-08 04:40:37 -05:00
- segment_fluid_initial_ [ seg ] [ comp_idx ] ) / dt ;
2024-02-22 07:50:58 -06:00
MultisegmentWellAssemble ( * this ) .
2022-11-18 05:09:43 -06:00
assembleAccumulationTerm ( seg , comp_idx , accumulation_term , this - > linSys_ ) ;
2017-10-10 07:30:23 -05:00
}
}
2019-04-08 04:40:37 -05:00
// considering the contributions due to flowing out from the segment
{
2022-12-19 08:52:36 -06:00
const int seg_upwind = this - > segments_ . upwinding_segment ( seg ) ;
2021-09-06 05:58:16 -05:00
for ( int comp_idx = 0 ; comp_idx < this - > num_components_ ; + + comp_idx ) {
2022-12-19 02:52:48 -06:00
const EvalWell segment_rate =
this - > primary_variables_ . getSegmentRateUpwinding ( seg ,
seg_upwind ,
comp_idx ) *
this - > well_efficiency_factor_ ;
2024-02-22 07:50:58 -06:00
MultisegmentWellAssemble ( * this ) .
2022-11-18 05:09:43 -06:00
assembleOutflowTerm ( seg , seg_upwind , comp_idx , segment_rate , this - > linSys_ ) ;
2019-04-08 04:40:37 -05:00
}
}
2017-10-10 07:30:23 -05:00
// considering the contributions from the inlet segments
{
2022-12-19 08:52:36 -06:00
for ( const int inlet : this - > segments_ . inlets ( ) [ seg ] ) {
const int inlet_upwind = this - > segments_ . upwinding_segment ( inlet ) ;
2021-09-06 05:58:16 -05:00
for ( int comp_idx = 0 ; comp_idx < this - > num_components_ ; + + comp_idx ) {
2022-12-19 02:52:48 -06:00
const EvalWell inlet_rate =
this - > primary_variables_ . getSegmentRateUpwinding ( inlet ,
inlet_upwind ,
comp_idx ) *
this - > well_efficiency_factor_ ;
2024-02-22 07:50:58 -06:00
MultisegmentWellAssemble ( * this ) .
2022-11-18 05:09:43 -06:00
assembleInflowTerm ( seg , inlet , inlet_upwind , comp_idx , inlet_rate , this - > linSys_ ) ;
2017-10-10 07:30:23 -05:00
}
}
}
// calculating the perforation rate for each perforation that belongs to this segment
2022-12-19 02:52:48 -06:00
const EvalWell seg_pressure = this - > primary_variables_ . getSegmentPressure ( seg ) ;
2021-08-06 06:14:00 -05:00
auto & perf_data = ws . perf_data ;
2021-06-06 10:09:44 -05:00
auto & perf_rates = perf_data . phase_rates ;
2021-06-05 01:52:18 -05:00
auto & perf_press_state = perf_data . pressure ;
2022-12-19 08:52:36 -06:00
for ( const int perf : this - > segments_ . perforations ( ) [ seg ] ) {
2024-10-15 07:11:34 -05:00
const int local_perf_index = this - > pw_info_ . globalToLocal ( perf ) ;
if ( local_perf_index < 0 ) // then the perforation is not on this process
continue ;
const int cell_idx = this - > well_cells_ [ local_perf_index ] ;
2024-02-06 04:55:07 -06:00
const auto & int_quants = simulator . model ( ) . intensiveQuantities ( cell_idx , /*timeIdx=*/ 0 ) ;
2021-09-06 05:58:16 -05:00
std : : vector < EvalWell > mob ( this - > num_components_ , 0.0 ) ;
2024-10-15 07:11:34 -05:00
getMobility ( simulator , local_perf_index , mob , deferred_logger ) ;
2024-02-20 15:32:18 -06:00
const Scalar trans_mult = simulator . problem ( ) . template wellTransMultiplier < Scalar > ( int_quants , cell_idx ) ;
2024-02-06 04:55:07 -06:00
const auto & wellstate_nupcol = simulator . problem ( ) . wellModel ( ) . nupcolWellState ( ) . well ( this - > index_of_well_ ) ;
2024-10-15 07:11:34 -05:00
const std : : vector < Scalar > Tw = this - > wellIndex ( local_perf_index , int_quants , trans_mult , wellstate_nupcol ) ;
2021-09-06 05:58:16 -05:00
std : : vector < EvalWell > cq_s ( this - > num_components_ , 0.0 ) ;
2019-04-05 03:01:16 -05:00
EvalWell perf_press ;
2024-02-20 07:18:30 -06:00
PerforationRates < Scalar > perfRates ;
2023-05-07 14:56:09 -05:00
computePerfRate ( int_quants , mob , Tw , seg , perf , seg_pressure ,
allow_cf , cq_s , perf_press , perfRates , deferred_logger ) ;
2019-04-23 06:30:12 -05:00
// updating the solution gas rate and solution oil rate
2019-10-23 02:09:45 -05:00
if ( this - > isProducer ( ) ) {
2023-05-04 06:24:23 -05:00
ws . phase_mixing_rates [ ws . dissolved_gas ] + = perfRates . dis_gas ;
ws . phase_mixing_rates [ ws . vaporized_oil ] + = perfRates . vap_oil ;
2024-10-15 07:11:34 -05:00
perf_data . phase_mixing_rates [ local_perf_index ] [ ws . dissolved_gas ] = perfRates . dis_gas ;
perf_data . phase_mixing_rates [ local_perf_index ] [ ws . vaporized_oil ] = perfRates . vap_oil ;
2019-04-23 06:30:12 -05:00
}
2019-04-05 03:01:16 -05:00
// store the perf pressure and rates
2021-09-06 05:58:16 -05:00
for ( int comp_idx = 0 ; comp_idx < this - > num_components_ ; + + comp_idx ) {
2024-10-15 07:11:34 -05:00
perf_rates [ local_perf_index * this - > number_of_phases_ + this - > modelCompIdxToFlowCompIdx ( comp_idx ) ] = cq_s [ comp_idx ] . value ( ) ;
2019-04-05 03:01:16 -05:00
}
2024-10-15 07:11:34 -05:00
perf_press_state [ local_perf_index ] = perf_press . value ( ) ;
2017-10-10 07:30:23 -05:00
2021-09-06 05:58:16 -05:00
for ( int comp_idx = 0 ; comp_idx < this - > num_components_ ; + + comp_idx ) {
2017-10-10 07:30:23 -05:00
// the cq_s entering mass balance equations need to consider the efficiency factors.
2021-09-06 05:58:16 -05:00
const EvalWell cq_s_effective = cq_s [ comp_idx ] * this - > well_efficiency_factor_ ;
2017-10-10 07:30:23 -05:00
2024-10-15 07:11:34 -05:00
this - > connectionRates_ [ local_perf_index ] [ comp_idx ] = Base : : restrictEval ( cq_s_effective ) ;
2017-10-10 07:30:23 -05:00
2024-02-22 07:50:58 -06:00
MultisegmentWellAssemble ( * this ) .
2024-10-15 07:11:34 -05:00
assemblePerforationEq ( seg , local_perf_index , comp_idx , cq_s_effective , this - > linSys_ ) ;
2017-10-10 07:30:23 -05:00
}
}
2024-03-19 06:50:34 -05:00
// the fourth equation, the pressure drop equation
2017-10-10 07:30:23 -05:00
if ( seg = = 0 ) { // top segment, pressure equation is the control equation
2024-02-06 04:55:07 -06:00
const auto & summaryState = simulator . vanguard ( ) . summaryState ( ) ;
const Schedule & schedule = simulator . vanguard ( ) . schedule ( ) ;
2024-03-19 06:50:34 -05:00
const bool stopped_or_zero_target = this - > stoppedOrZeroRateTarget ( simulator , well_state , deferred_logger ) ;
2024-02-22 07:50:58 -06:00
MultisegmentWellAssemble ( * this ) .
2022-11-18 04:57:37 -06:00
assembleControlEq ( well_state ,
2021-06-03 08:34:14 -05:00
group_state ,
schedule ,
summaryState ,
inj_controls ,
prod_controls ,
getRefDensity ( ) ,
2022-12-19 04:49:00 -06:00
this - > primary_variables_ ,
2022-11-18 04:57:37 -06:00
this - > linSys_ ,
2024-03-19 06:50:34 -05:00
stopped_or_zero_target ,
2021-06-03 08:34:14 -05:00
deferred_logger ) ;
2017-10-10 07:30:23 -05:00
} else {
2024-02-06 04:55:07 -06:00
const UnitSystem & unit_system = simulator . vanguard ( ) . eclState ( ) . getDeckUnitSystem ( ) ;
const auto & summary_state = simulator . vanguard ( ) . summaryState ( ) ;
2023-09-15 08:21:19 -05:00
this - > assemblePressureEq ( seg , unit_system , well_state , summary_state , this - > param_ . use_average_density_ms_wells_ , deferred_logger ) ;
2017-10-10 07:30:23 -05:00
}
}
2022-11-11 14:41:24 -06:00
2024-09-30 02:05:25 -05:00
this - > parallel_well_info_ . communication ( ) . sum ( this - > ipr_a_ . data ( ) , this - > ipr_a_ . size ( ) ) ;
2022-11-11 14:41:24 -06:00
this - > linSys_ . createSolver ( ) ;
2017-10-10 07:30:23 -05:00
}
2019-11-24 04:11:50 -06:00
2021-01-08 08:27:43 -06:00
2019-11-07 07:01:21 -06:00
template < typename TypeTag >
bool
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
openCrossFlowAvoidSingularity ( const Simulator & simulator ) const
2019-11-07 07:01:21 -06:00
{
2024-02-06 04:55:07 -06:00
return ! this - > getAllowCrossFlow ( ) & & allDrawDownWrongDirection ( simulator ) ;
2019-11-07 07:01:21 -06:00
}
template < typename TypeTag >
bool
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
allDrawDownWrongDirection ( const Simulator & simulator ) const
2019-11-07 07:01:21 -06:00
{
bool all_drawdown_wrong_direction = true ;
2021-06-03 08:34:14 -05:00
const int nseg = this - > numberOfSegments ( ) ;
2017-10-10 07:30:23 -05:00
2019-11-07 07:01:21 -06:00
for ( int seg = 0 ; seg < nseg ; + + seg ) {
2022-12-19 02:52:48 -06:00
const EvalWell segment_pressure = this - > primary_variables_ . getSegmentPressure ( seg ) ;
2022-12-19 08:52:36 -06:00
for ( const int perf : this - > segments_ . perforations ( ) [ seg ] ) {
2024-10-15 07:11:34 -05:00
const int local_perf_index = this - > pw_info_ . globalToLocal ( perf ) ;
if ( local_perf_index < 0 ) // then the perforation is not on this process
continue ;
2018-11-22 08:44:09 -06:00
2024-10-15 07:11:34 -05:00
const int cell_idx = this - > well_cells_ [ local_perf_index ] ;
2024-02-06 04:55:07 -06:00
const auto & intQuants = simulator . model ( ) . intensiveQuantities ( cell_idx , /*timeIdx=*/ 0 ) ;
2019-11-07 07:01:21 -06:00
const auto & fs = intQuants . fluidState ( ) ;
// pressure difference between the segment and the perforation
2022-12-19 08:52:36 -06:00
const EvalWell perf_seg_press_diff = this - > segments_ . getPressureDiffSegPerf ( seg , perf ) ;
2019-11-07 07:01:21 -06:00
// pressure difference between the perforation and the grid cell
2024-10-15 07:11:34 -05:00
const Scalar cell_perf_press_diff = this - > cell_perforation_pressure_diffs_ [ local_perf_index ] ;
2019-11-07 07:01:21 -06:00
2024-02-20 15:32:18 -06:00
const Scalar pressure_cell = this - > getPerfCellPressure ( fs ) . value ( ) ;
const Scalar perf_press = pressure_cell - cell_perf_press_diff ;
2019-11-07 07:01:21 -06:00
// Pressure drawdown (also used to determine direction of flow)
// TODO: not 100% sure about the sign of the seg_perf_press_diff
const EvalWell drawdown = perf_press - ( segment_pressure + perf_seg_press_diff ) ;
// for now, if there is one perforation can produce/inject in the correct
// direction, we consider this well can still produce/inject.
// TODO: it can be more complicated than this to cause wrong-signed rates
2019-10-23 02:09:45 -05:00
if ( ( drawdown < 0. & & this - > isInjector ( ) ) | |
( drawdown > 0. & & this - > isProducer ( ) ) ) {
2019-11-07 07:01:21 -06:00
all_drawdown_wrong_direction = false ;
break ;
}
}
}
2024-09-30 02:05:25 -05:00
const auto & comm = this - > parallel_well_info_ . communication ( ) ;
if ( comm . size ( ) > 1 )
{
all_drawdown_wrong_direction =
( comm . min ( all_drawdown_wrong_direction ? 1 : 0 ) = = 1 ) ;
}
2018-11-22 08:44:09 -06:00
2019-11-07 07:01:21 -06:00
return all_drawdown_wrong_direction ;
}
2018-11-22 08:44:09 -06:00
2018-11-30 05:15:51 -06:00
template < typename TypeTag >
void
MultisegmentWell < TypeTag > : :
2024-02-17 11:13:46 -06:00
updateWaterThroughput ( const double /*dt*/ , WellState < Scalar > & /*well_state*/ ) const
2018-11-30 05:15:51 -06:00
{
}
2019-02-12 05:48:51 -06:00
template < typename TypeTag >
typename MultisegmentWell < TypeTag > : : EvalWell
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
getSegmentSurfaceVolume ( const Simulator & simulator , const int seg_idx ) const
2019-02-12 05:48:51 -06:00
{
EvalWell temperature ;
2020-07-02 06:44:01 -05:00
EvalWell saltConcentration ;
2019-02-12 05:48:51 -06:00
int pvt_region_index ;
{
// using the pvt region of first perforated cell
// TODO: it should be a member of the WellInterface, initialized properly
2021-09-06 05:58:16 -05:00
const int cell_idx = this - > well_cells_ [ 0 ] ;
2024-02-06 04:55:07 -06:00
const auto & intQuants = simulator . model ( ) . intensiveQuantities ( cell_idx , /*timeIdx=*/ 0 ) ;
2019-02-12 05:48:51 -06:00
const auto & fs = intQuants . fluidState ( ) ;
temperature . setValue ( fs . temperature ( FluidSystem : : oilPhaseIdx ) . value ( ) ) ;
2021-06-03 08:34:14 -05:00
saltConcentration = this - > extendEval ( fs . saltConcentration ( ) ) ;
2019-02-12 05:48:51 -06:00
pvt_region_index = fs . pvtRegionIndex ( ) ;
}
2022-12-19 07:20:55 -06:00
return this - > segments_ . getSurfaceVolume ( temperature ,
saltConcentration ,
this - > primary_variables_ ,
pvt_region_index ,
seg_idx ) ;
2019-04-12 08:43:39 -05:00
}
2019-11-18 09:33:10 -06:00
2022-02-07 04:28:35 -06:00
template < typename TypeTag >
2024-02-20 15:32:18 -06:00
std : : optional < typename MultisegmentWell < TypeTag > : : Scalar >
2022-02-07 04:28:35 -06:00
MultisegmentWell < TypeTag > : :
2024-02-17 11:13:46 -06:00
computeBhpAtThpLimitProd ( const WellState < Scalar > & well_state ,
2024-02-06 04:55:07 -06:00
const Simulator & simulator ,
2022-02-07 04:28:35 -06:00
const SummaryState & summary_state ,
DeferredLogger & deferred_logger ) const
{
return this - > MultisegmentWell < TypeTag > : : computeBhpAtThpLimitProdWithAlq (
2024-02-06 04:55:07 -06:00
simulator ,
2022-02-07 04:28:35 -06:00
summary_state ,
2022-10-19 05:23:02 -05:00
this - > getALQ ( well_state ) ,
2024-10-23 02:09:09 -05:00
deferred_logger ,
/*iterate_if_no_solution */ true ) ;
2022-02-07 04:28:35 -06:00
}
2019-11-18 09:33:10 -06:00
2019-11-29 10:23:26 -06:00
2019-11-18 09:33:10 -06:00
template < typename TypeTag >
2024-02-20 15:32:18 -06:00
std : : optional < typename MultisegmentWell < TypeTag > : : Scalar >
2019-11-18 09:33:10 -06:00
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
computeBhpAtThpLimitProdWithAlq ( const Simulator & simulator ,
2022-02-07 04:28:35 -06:00
const SummaryState & summary_state ,
2024-02-20 15:32:18 -06:00
const Scalar alq_value ,
2024-10-23 02:09:09 -05:00
DeferredLogger & deferred_logger ,
bool iterate_if_no_solution ) const
2019-11-18 09:33:10 -06:00
{
// Make the frates() function.
2024-02-20 15:32:18 -06:00
auto frates = [ this , & simulator , & deferred_logger ] ( const Scalar bhp ) {
2019-11-18 09:33:10 -06:00
// Not solving the well equations here, which means we are
// calculating at the current Fg/Fw values of the
// well. This does not matter unless the well is
// crossflowing, and then it is likely still a good
// approximation.
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > rates ( 3 ) ;
2024-02-06 04:55:07 -06:00
computeWellRatesWithBhp ( simulator , bhp , rates , deferred_logger ) ;
2019-11-18 09:33:10 -06:00
return rates ;
} ;
2022-10-19 05:04:02 -05:00
auto bhpAtLimit = WellBhpThpCalculator ( * this ) .
computeBhpAtThpLimitProd ( frates ,
2021-06-03 08:34:14 -05:00
summary_state ,
2024-02-06 04:55:07 -06:00
this - > maxPerfPress ( simulator ) ,
2022-10-19 05:04:02 -05:00
this - > getRefDensity ( ) ,
alq_value ,
2022-11-23 08:37:04 -06:00
this - > getTHPConstraint ( summary_state ) ,
2022-10-19 05:04:02 -05:00
deferred_logger ) ;
2021-09-13 03:58:50 -05:00
2022-10-19 05:04:02 -05:00
if ( bhpAtLimit )
2021-09-13 08:11:38 -05:00
return bhpAtLimit ;
2021-09-13 03:58:50 -05:00
2024-10-23 02:09:09 -05:00
if ( ! iterate_if_no_solution )
return std : : nullopt ;
2024-02-20 15:32:18 -06:00
auto fratesIter = [ this , & simulator , & deferred_logger ] ( const Scalar bhp ) {
2021-09-13 08:11:38 -05:00
// Solver the well iterations to see if we are
// able to get a solution with an update
// solution
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > rates ( 3 ) ;
2024-02-06 04:55:07 -06:00
computeWellRatesWithBhpIterations ( simulator , bhp , rates , deferred_logger ) ;
2021-09-13 03:58:50 -05:00
return rates ;
} ;
2022-10-19 05:04:02 -05:00
return WellBhpThpCalculator ( * this ) .
computeBhpAtThpLimitProd ( fratesIter ,
2021-09-13 03:58:50 -05:00
summary_state ,
2024-02-06 04:55:07 -06:00
this - > maxPerfPress ( simulator ) ,
2022-10-19 05:04:02 -05:00
this - > getRefDensity ( ) ,
alq_value ,
2022-11-23 08:37:04 -06:00
this - > getTHPConstraint ( summary_state ) ,
2022-10-19 05:04:02 -05:00
deferred_logger ) ;
2019-11-24 04:11:50 -06:00
}
2019-11-18 09:33:10 -06:00
template < typename TypeTag >
2024-02-20 15:32:18 -06:00
std : : optional < typename MultisegmentWell < TypeTag > : : Scalar >
2019-11-18 09:33:10 -06:00
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
computeBhpAtThpLimitInj ( const Simulator & simulator ,
2019-11-18 09:33:10 -06:00
const SummaryState & summary_state ,
DeferredLogger & deferred_logger ) const
{
// Make the frates() function.
2024-02-20 15:32:18 -06:00
auto frates = [ this , & simulator , & deferred_logger ] ( const Scalar bhp ) {
2019-11-18 09:33:10 -06:00
// Not solving the well equations here, which means we are
// calculating at the current Fg/Fw values of the
// well. This does not matter unless the well is
// crossflowing, and then it is likely still a good
// approximation.
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > rates ( 3 ) ;
2024-02-06 04:55:07 -06:00
computeWellRatesWithBhp ( simulator , bhp , rates , deferred_logger ) ;
2019-11-18 09:33:10 -06:00
return rates ;
} ;
2022-10-19 02:55:14 -05:00
auto bhpAtLimit = WellBhpThpCalculator ( * this ) .
2021-09-13 03:58:50 -05:00
computeBhpAtThpLimitInj ( frates ,
summary_state ,
2022-10-19 02:55:14 -05:00
this - > getRefDensity ( ) ,
0.05 ,
100 ,
false ,
2021-09-13 03:58:50 -05:00
deferred_logger ) ;
2022-10-19 02:55:14 -05:00
if ( bhpAtLimit )
2021-09-13 08:11:38 -05:00
return bhpAtLimit ;
2021-09-13 03:58:50 -05:00
2024-02-20 15:32:18 -06:00
auto fratesIter = [ this , & simulator , & deferred_logger ] ( const Scalar bhp ) {
2021-09-13 08:11:38 -05:00
// Solver the well iterations to see if we are
// able to get a solution with an update
// solution
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > rates ( 3 ) ;
2024-02-06 04:55:07 -06:00
computeWellRatesWithBhpIterations ( simulator , bhp , rates , deferred_logger ) ;
2021-09-13 03:58:50 -05:00
return rates ;
} ;
2022-10-19 02:55:14 -05:00
return WellBhpThpCalculator ( * this ) .
computeBhpAtThpLimitInj ( fratesIter ,
summary_state ,
this - > getRefDensity ( ) ,
0.05 ,
100 ,
false ,
deferred_logger ) ;
2019-11-18 09:33:10 -06:00
}
2019-11-24 04:11:50 -06:00
2019-11-18 09:33:10 -06:00
template < typename TypeTag >
2024-02-20 15:32:18 -06:00
typename MultisegmentWell < TypeTag > : : Scalar
2019-11-18 09:33:10 -06:00
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
maxPerfPress ( const Simulator & simulator ) const
2019-11-18 09:33:10 -06:00
{
2024-02-20 15:32:18 -06:00
Scalar max_pressure = 0.0 ;
2021-06-03 08:34:14 -05:00
const int nseg = this - > numberOfSegments ( ) ;
2019-11-18 09:33:10 -06:00
for ( int seg = 0 ; seg < nseg ; + + seg ) {
2022-12-19 08:52:36 -06:00
for ( const int perf : this - > segments_ . perforations ( ) [ seg ] ) {
2024-10-15 07:11:34 -05:00
const int local_perf_index = this - > pw_info_ . globalToLocal ( perf ) ;
if ( local_perf_index < 0 ) // then the perforation is not on this process
continue ;
const int cell_idx = this - > well_cells_ [ local_perf_index ] ;
2024-02-06 04:55:07 -06:00
const auto & int_quants = simulator . model ( ) . intensiveQuantities ( cell_idx , /*timeIdx=*/ 0 ) ;
2019-11-18 09:33:10 -06:00
const auto & fs = int_quants . fluidState ( ) ;
2024-02-20 15:32:18 -06:00
Scalar pressure_cell = this - > getPerfCellPressure ( fs ) . value ( ) ;
2019-11-18 09:33:10 -06:00
max_pressure = std : : max ( max_pressure , pressure_cell ) ;
}
}
return max_pressure ;
}
2019-11-24 04:11:50 -06:00
2020-10-15 07:15:05 -05:00
template < typename TypeTag >
2024-02-20 15:32:18 -06:00
std : : vector < typename MultisegmentWell < TypeTag > : : Scalar >
2020-10-15 07:15:05 -05:00
MultisegmentWell < TypeTag > : :
2024-02-06 04:55:07 -06:00
computeCurrentWellRates ( const Simulator & simulator ,
2020-10-15 07:15:05 -05:00
DeferredLogger & deferred_logger ) const
{
// Calculate the rates that follow from the current primary variables.
2021-09-13 02:36:16 -05:00
std : : vector < Scalar > well_q_s ( this - > num_components_ , 0.0 ) ;
2024-02-06 04:55:07 -06:00
const bool allow_cf = this - > getAllowCrossFlow ( ) | | openCrossFlowAvoidSingularity ( simulator ) ;
2021-06-03 08:34:14 -05:00
const int nseg = this - > numberOfSegments ( ) ;
2020-10-15 09:24:55 -05:00
for ( int seg = 0 ; seg < nseg ; + + seg ) {
// calculating the perforation rate for each perforation that belongs to this segment
2022-12-19 02:52:48 -06:00
const Scalar seg_pressure = getValue ( this - > primary_variables_ . getSegmentPressure ( seg ) ) ;
2022-12-19 08:52:36 -06:00
for ( const int perf : this - > segments_ . perforations ( ) [ seg ] ) {
2024-10-15 07:11:34 -05:00
const int local_perf_index = this - > pw_info_ . globalToLocal ( perf ) ;
if ( local_perf_index < 0 ) // then the perforation is not on this process
continue ;
const int cell_idx = this - > well_cells_ [ local_perf_index ] ;
2024-02-06 04:55:07 -06:00
const auto & int_quants = simulator . model ( ) . intensiveQuantities ( cell_idx , /*timeIdx=*/ 0 ) ;
2021-09-13 02:36:16 -05:00
std : : vector < Scalar > mob ( this - > num_components_ , 0.0 ) ;
2024-10-15 07:11:34 -05:00
getMobility ( simulator , local_perf_index , mob , deferred_logger ) ;
2024-02-20 15:32:18 -06:00
const Scalar trans_mult = simulator . problem ( ) . template wellTransMultiplier < Scalar > ( int_quants , cell_idx ) ;
2024-02-06 04:55:07 -06:00
const auto & wellstate_nupcol = simulator . problem ( ) . wellModel ( ) . nupcolWellState ( ) . well ( this - > index_of_well_ ) ;
2024-10-15 07:11:34 -05:00
const std : : vector < Scalar > Tw = this - > wellIndex ( local_perf_index , int_quants , trans_mult , wellstate_nupcol ) ;
2021-09-13 02:36:16 -05:00
std : : vector < Scalar > cq_s ( this - > num_components_ , 0.0 ) ;
2023-05-07 14:56:09 -05:00
Scalar perf_press = 0.0 ;
2024-02-20 07:18:30 -06:00
PerforationRates < Scalar > perf_rates ;
2023-05-07 14:56:09 -05:00
computePerfRate ( int_quants , mob , Tw , seg , perf , seg_pressure ,
allow_cf , cq_s , perf_press , perf_rates , deferred_logger ) ;
2021-09-06 05:58:16 -05:00
for ( int comp = 0 ; comp < this - > num_components_ ; + + comp ) {
2020-10-15 09:24:55 -05:00
well_q_s [ comp ] + = cq_s [ comp ] ;
}
2020-10-15 07:15:05 -05:00
}
}
2024-09-30 02:05:25 -05:00
const auto & comm = this - > parallel_well_info_ . communication ( ) ;
if ( comm . size ( ) > 1 )
{
comm . sum ( well_q_s . data ( ) , well_q_s . size ( ) ) ;
}
2021-09-13 02:36:16 -05:00
return well_q_s ;
2020-03-24 05:07:14 -05:00
}
2020-10-15 07:15:05 -05:00
2023-06-13 05:58:07 -05:00
template < typename TypeTag >
2024-02-20 15:32:18 -06:00
std : : vector < typename MultisegmentWell < TypeTag > : : Scalar >
2023-06-13 05:58:07 -05:00
MultisegmentWell < TypeTag > : :
getPrimaryVars ( ) const
{
const int num_seg = this - > numberOfSegments ( ) ;
2023-06-13 05:58:07 -05:00
constexpr int num_eq = MSWEval : : numWellEq ;
2024-02-20 15:32:18 -06:00
std : : vector < Scalar > retval ( num_seg * num_eq ) ;
2023-06-13 05:58:07 -05:00
for ( int ii = 0 ; ii < num_seg ; + + ii ) {
const auto & pv = this - > primary_variables_ . value ( ii ) ;
2023-06-13 05:58:07 -05:00
std : : copy ( pv . begin ( ) , pv . end ( ) , retval . begin ( ) + ii * num_eq ) ;
2023-06-13 05:58:07 -05:00
}
return retval ;
}
template < typename TypeTag >
int
MultisegmentWell < TypeTag > : :
2024-02-20 15:32:18 -06:00
setPrimaryVars ( typename std : : vector < Scalar > : : const_iterator it )
2023-06-13 05:58:07 -05:00
{
const int num_seg = this - > numberOfSegments ( ) ;
2023-06-13 05:58:07 -05:00
constexpr int num_eq = MSWEval : : numWellEq ;
2024-02-20 15:32:18 -06:00
std : : array < Scalar , num_eq > tmp ;
2023-06-13 05:58:07 -05:00
for ( int ii = 0 ; ii < num_seg ; + + ii ) {
2024-10-15 09:31:38 -05:00
const auto start = it + ii * num_eq ;
2023-06-13 05:58:07 -05:00
std : : copy ( start , start + num_eq , tmp . begin ( ) ) ;
2023-06-13 05:58:07 -05:00
this - > primary_variables_ . setValue ( ii , tmp ) ;
}
2023-06-13 05:58:07 -05:00
return num_seg * num_eq ;
2023-06-13 05:58:07 -05:00
}
2020-10-15 07:15:05 -05:00
} // namespace Opm