2021-05-19 07:51:14 -05:00
|
|
|
/*
|
|
|
|
Copyright 2016 SINTEF ICT, Applied Mathematics.
|
|
|
|
Copyright 2016 - 2017 Statoil ASA.
|
|
|
|
Copyright 2017 Dr. Blatt - HPC-Simulation-Software & Services
|
|
|
|
Copyright 2016 - 2018 IRIS AS
|
|
|
|
|
|
|
|
This file is part of the Open Porous Media project (OPM).
|
|
|
|
|
|
|
|
OPM is free software: you can redistribute it and/or modify
|
|
|
|
it under the terms of the GNU General Public License as published by
|
|
|
|
the Free Software Foundation, either version 3 of the License, or
|
|
|
|
(at your option) any later version.
|
|
|
|
|
|
|
|
OPM is distributed in the hope that it will be useful,
|
|
|
|
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
GNU General Public License for more details.
|
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
*/
|
|
|
|
|
|
|
|
#ifndef OPM_BLACKOILWELLMODEL_GENERIC_HEADER_INCLUDED
|
|
|
|
#define OPM_BLACKOILWELLMODEL_GENERIC_HEADER_INCLUDED
|
|
|
|
|
|
|
|
#include <opm/output/data/GuideRateValue.hpp>
|
2021-10-07 06:21:35 -05:00
|
|
|
|
2021-12-14 01:30:15 -06:00
|
|
|
#include <opm/input/eclipse/Schedule/Well/WellTestState.hpp>
|
|
|
|
#include <opm/input/eclipse/Schedule/Group/GuideRate.hpp>
|
2021-05-19 07:51:14 -05:00
|
|
|
|
2021-06-07 06:04:29 -05:00
|
|
|
#include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>
|
2021-10-07 06:21:35 -05:00
|
|
|
|
2021-05-19 07:51:14 -05:00
|
|
|
#include <opm/simulators/wells/PerforationData.hpp>
|
|
|
|
#include <opm/simulators/wells/WellProdIndexCalculator.hpp>
|
2021-06-07 05:15:36 -05:00
|
|
|
#include <opm/simulators/wells/WGState.hpp>
|
2021-05-19 07:51:14 -05:00
|
|
|
|
2021-10-07 06:21:35 -05:00
|
|
|
#include <functional>
|
|
|
|
#include <map>
|
|
|
|
#include <memory>
|
|
|
|
#include <string>
|
|
|
|
#include <unordered_map>
|
|
|
|
#include <unordered_set>
|
|
|
|
#include <vector>
|
|
|
|
|
2021-05-19 07:51:14 -05:00
|
|
|
namespace Opm {
|
2021-10-07 06:21:35 -05:00
|
|
|
class DeferredLogger;
|
|
|
|
class EclipseState;
|
|
|
|
class GasLiftSingleWellGeneric;
|
|
|
|
class GasLiftWellState;
|
|
|
|
class Group;
|
|
|
|
class GuideRateConfig;
|
|
|
|
class ParallelWellInfo;
|
|
|
|
class RestartValue;
|
|
|
|
class Schedule;
|
|
|
|
class SummaryConfig;
|
|
|
|
class VFPProperties;
|
|
|
|
class WellInterfaceGeneric;
|
|
|
|
class WellState;
|
|
|
|
} // namespace Opm
|
2021-05-19 07:51:14 -05:00
|
|
|
|
2021-10-07 06:21:35 -05:00
|
|
|
namespace Opm { namespace data {
|
|
|
|
struct GroupData;
|
|
|
|
struct GroupGuideRates;
|
|
|
|
class GroupAndNetworkValues;
|
|
|
|
struct NodeData;
|
|
|
|
}} // namespace Opm::data
|
|
|
|
|
|
|
|
namespace Opm {
|
2021-05-19 07:51:14 -05:00
|
|
|
|
|
|
|
/// Class for handling the blackoil well model.
|
|
|
|
class BlackoilWellModelGeneric
|
|
|
|
{
|
|
|
|
public:
|
|
|
|
// --------- Types ---------
|
2021-10-07 06:21:35 -05:00
|
|
|
using GLiftOptWells = std::map<std::string, std::unique_ptr<GasLiftSingleWellGeneric>>;
|
|
|
|
using GLiftProdWells = std::map<std::string, const WellInterfaceGeneric*>;
|
|
|
|
using GLiftWellStateMap = std::map<std::string, std::unique_ptr<GasLiftWellState>>;
|
2021-05-19 07:51:14 -05:00
|
|
|
|
2021-06-07 08:01:10 -05:00
|
|
|
BlackoilWellModelGeneric(Schedule& schedule,
|
2021-05-19 07:51:14 -05:00
|
|
|
const SummaryState& summaryState,
|
|
|
|
const EclipseState& eclState,
|
|
|
|
const PhaseUsage& phase_usage,
|
2021-05-25 05:57:11 -05:00
|
|
|
const Parallel::Communication& comm);
|
2021-05-19 07:51:14 -05:00
|
|
|
|
|
|
|
virtual ~BlackoilWellModelGeneric() = default;
|
|
|
|
|
|
|
|
int numLocalWells() const;
|
|
|
|
int numPhases() const;
|
|
|
|
|
|
|
|
/// return true if wells are available in the reservoir
|
|
|
|
bool wellsActive() const;
|
|
|
|
bool hasWell(const std::string& wname);
|
|
|
|
/// return true if wells are available on this process
|
|
|
|
bool localWellsActive() const;
|
|
|
|
// whether there exists any multisegment well open on this process
|
|
|
|
bool anyMSWellOpenLocal() const;
|
|
|
|
|
|
|
|
const Well& getWellEcl(const std::string& well_name) const;
|
|
|
|
std::vector<Well> getLocalWells(const int timeStepIdx) const;
|
|
|
|
const Schedule& schedule() const { return schedule_; }
|
|
|
|
const PhaseUsage& phaseUsage() const { return phase_usage_; }
|
|
|
|
const GroupState& groupState() const { return this->active_wgstate_.group_state; }
|
|
|
|
|
|
|
|
/*
|
|
|
|
Immutable version of the currently active wellstate.
|
|
|
|
*/
|
|
|
|
const WellState& wellState() const
|
|
|
|
{
|
|
|
|
return this->active_wgstate_.well_state;
|
|
|
|
}
|
|
|
|
|
|
|
|
/*
|
|
|
|
Mutable version of the currently active wellstate.
|
|
|
|
*/
|
|
|
|
WellState& wellState()
|
|
|
|
{
|
|
|
|
return this->active_wgstate_.well_state;
|
|
|
|
}
|
|
|
|
|
2021-06-23 08:46:33 -05:00
|
|
|
GroupState& groupState() { return this->active_wgstate_.group_state; }
|
|
|
|
|
2021-09-21 02:30:02 -05:00
|
|
|
WellTestState& wellTestState() { return this->active_wgstate_.well_test_state; }
|
|
|
|
|
|
|
|
const WellTestState& wellTestState() const { return this->active_wgstate_.well_test_state; }
|
|
|
|
|
2021-06-23 08:46:33 -05:00
|
|
|
|
2021-05-19 07:51:14 -05:00
|
|
|
double wellPI(const int well_index) const;
|
|
|
|
double wellPI(const std::string& well_name) const;
|
|
|
|
|
|
|
|
void updateEclWells(const int timeStepIdx,
|
2021-11-01 02:49:46 -05:00
|
|
|
const std::unordered_set<std::string>& wells,
|
|
|
|
const SummaryState& st);
|
2021-05-19 07:51:14 -05:00
|
|
|
|
|
|
|
|
|
|
|
void loadRestartData(const data::Wells& rst_wells,
|
|
|
|
const data::GroupAndNetworkValues& grpNwrkValues,
|
|
|
|
const PhaseUsage& phases,
|
|
|
|
const bool handle_ms_well,
|
|
|
|
WellState& well_state);
|
|
|
|
|
|
|
|
void initFromRestartFile(const RestartValue& restartValues,
|
2021-10-11 07:38:41 -05:00
|
|
|
WellTestState wtestState,
|
2021-05-19 07:51:14 -05:00
|
|
|
const size_t numCells,
|
|
|
|
bool handle_ms_well);
|
|
|
|
|
|
|
|
void setWellsActive(const bool wells_active);
|
|
|
|
|
|
|
|
/*
|
|
|
|
Will assign the internal member last_valid_well_state_ to the
|
|
|
|
current value of the this->active_well_state_. The state stored
|
|
|
|
with storeWellState() can then subsequently be recovered with the
|
|
|
|
resetWellState() method.
|
|
|
|
*/
|
|
|
|
void commitWGState()
|
|
|
|
{
|
|
|
|
this->last_valid_wgstate_ = this->active_wgstate_;
|
|
|
|
}
|
|
|
|
|
|
|
|
data::GroupAndNetworkValues groupAndNetworkData(const int reportStepIdx) const;
|
|
|
|
|
2021-06-07 05:20:49 -05:00
|
|
|
/// Return true if any well has a THP constraint.
|
|
|
|
bool hasTHPConstraints() const;
|
|
|
|
|
2021-10-14 06:14:38 -05:00
|
|
|
/// Shut down any single well
|
2021-06-07 05:20:49 -05:00
|
|
|
/// Returns true if the well was actually found and shut.
|
2021-10-14 06:14:38 -05:00
|
|
|
bool forceShutWellByName(const std::string& wellname,
|
|
|
|
const double simulation_time);
|
2021-06-07 05:20:49 -05:00
|
|
|
|
2021-05-19 07:51:14 -05:00
|
|
|
protected:
|
|
|
|
|
|
|
|
/*
|
|
|
|
The dynamic state of the well model is maintained with an instance
|
|
|
|
of the WellState class. Currently we have
|
|
|
|
three different wellstate instances:
|
|
|
|
|
|
|
|
1. The currently active wellstate is in the active_well_state_
|
|
|
|
member. That is the state which is mutated by the simulator.
|
|
|
|
|
|
|
|
2. In the case timestep fails to converge and we must go back and
|
|
|
|
try again with a smaller timestep we need to recover the last
|
|
|
|
valid wellstate. This is maintained with the
|
|
|
|
last_valid_well_state_ member and the functions
|
|
|
|
commitWellState() and resetWellState().
|
|
|
|
|
|
|
|
3. For the NUPCOL functionality we should either use the
|
|
|
|
currently active wellstate or a wellstate frozen at max
|
|
|
|
nupcol iterations. This is handled with the member
|
|
|
|
nupcol_well_state_ and the initNupcolWellState() function.
|
|
|
|
*/
|
|
|
|
|
|
|
|
|
|
|
|
/*
|
|
|
|
Will return the last good wellstate. This is typcially used when
|
|
|
|
initializing a new report step where the Schedule object might
|
|
|
|
have introduced new wells. The wellstate returned by
|
|
|
|
prevWellState() must have been stored with the commitWellState()
|
|
|
|
function first.
|
|
|
|
*/
|
|
|
|
const WellState& prevWellState() const
|
|
|
|
{
|
|
|
|
return this->last_valid_wgstate_.well_state;
|
|
|
|
}
|
|
|
|
|
|
|
|
const WGState& prevWGState() const
|
|
|
|
{
|
|
|
|
return this->last_valid_wgstate_;
|
|
|
|
}
|
|
|
|
/*
|
|
|
|
Will return the currently active nupcolWellState; must initialize
|
|
|
|
the internal nupcol wellstate with initNupcolWellState() first.
|
|
|
|
*/
|
|
|
|
const WellState& nupcolWellState() const
|
|
|
|
{
|
|
|
|
return this->nupcol_wgstate_.well_state;
|
|
|
|
}
|
|
|
|
|
|
|
|
/*
|
|
|
|
Will store a copy of the input argument well_state in the
|
|
|
|
last_valid_well_state_ member, that state can then be recovered
|
|
|
|
with a subsequent call to resetWellState().
|
|
|
|
*/
|
|
|
|
void commitWGState(WGState wgstate)
|
|
|
|
{
|
|
|
|
this->last_valid_wgstate_ = std::move(wgstate);
|
|
|
|
}
|
|
|
|
|
|
|
|
/*
|
|
|
|
Will update the internal variable active_well_state_ to whatever
|
|
|
|
was stored in the last_valid_well_state_ member. This function
|
|
|
|
works in pair with commitWellState() which should be called first.
|
|
|
|
*/
|
|
|
|
void resetWGState()
|
|
|
|
{
|
|
|
|
this->active_wgstate_ = this->last_valid_wgstate_;
|
|
|
|
}
|
|
|
|
|
|
|
|
/*
|
|
|
|
Will store the current active wellstate in the nupcol_well_state_
|
|
|
|
member. This can then be subsequently retrieved with accessor
|
|
|
|
nupcolWellState().
|
|
|
|
*/
|
|
|
|
void updateNupcolWGState()
|
|
|
|
{
|
|
|
|
this->nupcol_wgstate_ = this->active_wgstate_;
|
|
|
|
}
|
|
|
|
|
|
|
|
/// \brief Create the parallel well information
|
|
|
|
/// \param localWells The local wells from ECL schedule
|
2021-09-27 05:00:39 -05:00
|
|
|
std::vector<std::reference_wrapper<ParallelWellInfo>> createLocalParallelWellInfo(const std::vector<Well>& wells);
|
2021-05-19 07:51:14 -05:00
|
|
|
|
|
|
|
void initializeWellProdIndCalculators();
|
|
|
|
void initializeWellPerfData();
|
|
|
|
|
|
|
|
bool wasDynamicallyShutThisTimeStep(const int well_index) const;
|
|
|
|
|
2021-06-07 05:20:49 -05:00
|
|
|
void updateNetworkPressures(const int reportStepIdx);
|
|
|
|
|
2021-05-19 07:51:14 -05:00
|
|
|
void updateWsolvent(const Group& group,
|
|
|
|
const int reportStepIdx,
|
|
|
|
const WellState& wellState);
|
|
|
|
void setWsolvent(const Group& group,
|
|
|
|
const int reportStepIdx,
|
|
|
|
double wsolvent);
|
|
|
|
virtual void calcRates(const int fipnum,
|
|
|
|
const int pvtreg,
|
|
|
|
std::vector<double>& resv_coeff) = 0;
|
2021-06-08 08:29:15 -05:00
|
|
|
virtual void calcInjRates(const int fipnum,
|
|
|
|
const int pvtreg,
|
|
|
|
std::vector<double>& resv_coeff) = 0;
|
2021-05-19 07:51:14 -05:00
|
|
|
|
|
|
|
data::GuideRateValue getGuideRateValues(const Group& group) const;
|
|
|
|
data::GuideRateValue getGuideRateValues(const Well& well) const;
|
|
|
|
data::GuideRateValue getGuideRateInjectionGroupValues(const Group& group) const;
|
|
|
|
void getGuideRateValues(const GuideRate::RateVector& qs,
|
|
|
|
const bool is_inj,
|
|
|
|
const std::string& wgname,
|
|
|
|
data::GuideRateValue& grval) const;
|
|
|
|
|
2021-10-25 06:06:02 -05:00
|
|
|
void assignWellGuideRates(data::Wells& wsrpt,
|
|
|
|
const int reportStepIdx) const;
|
2021-05-19 07:51:14 -05:00
|
|
|
void assignShutConnections(data::Wells& wsrpt,
|
|
|
|
const int reportStepIndex) const;
|
|
|
|
void assignGroupControl(const Group& group,
|
|
|
|
data::GroupData& gdata) const;
|
|
|
|
void assignGroupGuideRates(const Group& group,
|
|
|
|
const std::unordered_map<std::string, data::GroupGuideRates>& groupGuideRates,
|
|
|
|
data::GroupData& gdata) const;
|
|
|
|
void assignGroupValues(const int reportStepIdx,
|
|
|
|
std::map<std::string, data::GroupData>& gvalues) const;
|
|
|
|
void assignNodeValues(std::map<std::string, data::NodeData>& nodevalues) const;
|
|
|
|
|
2021-10-07 06:21:35 -05:00
|
|
|
void loadRestartConnectionData(const std::vector<data::Rates::opt>& phs,
|
|
|
|
const data::Well& rst_well,
|
|
|
|
const std::vector<PerforationData>& old_perf_data,
|
|
|
|
SingleWellState& ws);
|
|
|
|
|
|
|
|
void loadRestartSegmentData(const std::string& well_name,
|
|
|
|
const std::vector<data::Rates::opt>& phs,
|
|
|
|
const data::Well& rst_well,
|
|
|
|
SingleWellState& ws);
|
|
|
|
|
|
|
|
void loadRestartWellData(const std::string& well_name,
|
|
|
|
const bool handle_ms_well,
|
|
|
|
const std::vector<data::Rates::opt>& phs,
|
|
|
|
const data::Well& rst_well,
|
|
|
|
const std::vector<PerforationData>& old_perf_data,
|
|
|
|
SingleWellState& ws);
|
|
|
|
|
|
|
|
void loadRestartGroupData(const std::string& group,
|
|
|
|
const data::GroupData& value);
|
|
|
|
|
2021-10-07 06:41:36 -05:00
|
|
|
void loadRestartGuideRates(const int report_step,
|
|
|
|
const GuideRateModel::Target target,
|
|
|
|
const data::Wells& rst_wells);
|
|
|
|
|
|
|
|
void loadRestartGuideRates(const int report_step,
|
|
|
|
const GuideRateConfig& config,
|
|
|
|
const std::map<std::string, data::GroupData>& rst_groups);
|
|
|
|
|
2021-05-19 07:51:14 -05:00
|
|
|
std::unordered_map<std::string, data::GroupGuideRates>
|
|
|
|
calculateAllGroupGuiderates(const int reportStepIdx) const;
|
|
|
|
|
2021-06-07 05:20:49 -05:00
|
|
|
void calculateEfficiencyFactors(const int reportStepIdx);
|
|
|
|
|
2021-05-19 07:51:14 -05:00
|
|
|
bool checkGroupConstraints(const Group& group,
|
|
|
|
const int reportStepIdx,
|
|
|
|
DeferredLogger& deferred_logger) const;
|
|
|
|
|
2021-06-09 08:15:12 -05:00
|
|
|
std::pair<Group::InjectionCMode, double> checkGroupInjectionConstraints(const Group& group,
|
2021-05-19 07:51:14 -05:00
|
|
|
const int reportStepIdx,
|
|
|
|
const Phase& phase) const;
|
2021-06-09 08:15:12 -05:00
|
|
|
std::pair<Group::ProductionCMode, double> checkGroupProductionConstraints(const Group& group,
|
2021-05-19 07:51:14 -05:00
|
|
|
const int reportStepIdx,
|
|
|
|
DeferredLogger& deferred_logger) const;
|
|
|
|
|
|
|
|
void checkGconsaleLimits(const Group& group,
|
|
|
|
WellState& well_state,
|
|
|
|
const int reportStepIdx,
|
|
|
|
DeferredLogger& deferred_logger);
|
|
|
|
|
2021-08-27 00:57:16 -05:00
|
|
|
bool checkGroupHigherConstraints(const Group& group,
|
2021-05-19 07:51:14 -05:00
|
|
|
DeferredLogger& deferred_logger,
|
|
|
|
const int reportStepIdx,
|
|
|
|
std::set<std::string>& switched_groups);
|
|
|
|
|
2021-08-27 00:57:16 -05:00
|
|
|
bool updateGroupIndividualControl(const Group& group,
|
2021-05-19 07:51:14 -05:00
|
|
|
DeferredLogger& deferred_logger,
|
|
|
|
const int reportStepIdx,
|
|
|
|
std::set<std::string>& switched_groups);
|
|
|
|
|
2021-08-27 00:57:16 -05:00
|
|
|
bool updateGroupIndividualControls(DeferredLogger& deferred_logger,
|
2021-05-19 07:51:14 -05:00
|
|
|
std::set<std::string>& switched_groups,
|
|
|
|
const int reportStepIdx,
|
|
|
|
const int iterationIdx);
|
|
|
|
|
2021-08-27 00:57:16 -05:00
|
|
|
bool updateGroupHigherControls(DeferredLogger& deferred_logger,
|
2021-05-19 07:51:14 -05:00
|
|
|
const int reportStepIdx,
|
|
|
|
std::set<std::string>& switched_groups);
|
|
|
|
|
|
|
|
void actionOnBrokenConstraints(const Group& group,
|
|
|
|
const Group::ExceedAction& exceed_action,
|
|
|
|
const Group::ProductionCMode& newControl,
|
|
|
|
DeferredLogger& deferred_logger);
|
|
|
|
void actionOnBrokenConstraints(const Group& group,
|
|
|
|
const Group::InjectionCMode& newControl,
|
|
|
|
const Phase& controlPhase,
|
|
|
|
DeferredLogger& deferred_logger);
|
|
|
|
|
|
|
|
void updateAndCommunicateGroupData(const int reportStepIdx,
|
|
|
|
const int iterationIdx);
|
|
|
|
|
2021-06-07 05:20:49 -05:00
|
|
|
void inferLocalShutWells();
|
|
|
|
|
2021-06-07 06:04:29 -05:00
|
|
|
void setRepRadiusPerfLength();
|
|
|
|
|
2021-06-07 06:04:29 -05:00
|
|
|
void gliftDebug(const std::string& msg,
|
|
|
|
DeferredLogger& deferred_logger) const;
|
|
|
|
|
|
|
|
void gliftDebugShowALQ(DeferredLogger& deferred_logger);
|
|
|
|
|
|
|
|
void gasLiftOptimizationStage2(DeferredLogger& deferred_logger,
|
|
|
|
GLiftProdWells& prod_wells,
|
|
|
|
GLiftOptWells& glift_wells,
|
|
|
|
GLiftWellStateMap& map,
|
|
|
|
const int episodeIndex);
|
|
|
|
|
2021-06-07 06:04:29 -05:00
|
|
|
virtual void computePotentials(const std::size_t widx,
|
|
|
|
const WellState& well_state_copy,
|
|
|
|
std::string& exc_msg,
|
|
|
|
ExceptionType::ExcEnum& exc_type,
|
|
|
|
DeferredLogger& deferred_logger) = 0;
|
|
|
|
|
|
|
|
// Calculating well potentials for each well
|
|
|
|
void updateWellPotentials(const int reportStepIdx,
|
|
|
|
const bool onlyAfterEvent,
|
|
|
|
const SummaryConfig& summaryConfig,
|
|
|
|
DeferredLogger& deferred_logger);
|
|
|
|
|
2022-01-17 06:01:42 -06:00
|
|
|
bool guideRateUpdateIsNeeded(const int reportStepIdx) const;
|
2021-11-18 05:57:16 -06:00
|
|
|
|
2021-06-07 08:01:10 -05:00
|
|
|
// create the well container
|
|
|
|
virtual void createWellContainer(const int time_step) = 0;
|
2022-04-12 01:44:52 -05:00
|
|
|
virtual void initWellContainer(const int reportStepIdx) = 0;
|
2021-06-07 06:04:29 -05:00
|
|
|
|
2021-06-07 08:01:10 -05:00
|
|
|
virtual void calculateProductivityIndexValuesShutWells(const int reportStepIdx,
|
|
|
|
DeferredLogger& deferred_logger) = 0;
|
|
|
|
virtual void calculateProductivityIndexValues(DeferredLogger& deferred_logger) = 0;
|
|
|
|
|
|
|
|
void runWellPIScaling(const int timeStepIdx,
|
|
|
|
DeferredLogger& local_deferredLogger);
|
|
|
|
|
2021-10-19 09:50:42 -05:00
|
|
|
/// \brief get compressed index for interior cells (-1, otherwise
|
|
|
|
virtual int compressedIndexForInterior(int cartesian_cell_idx) const = 0;
|
2021-09-27 04:28:27 -05:00
|
|
|
|
|
|
|
|
2021-06-07 08:01:10 -05:00
|
|
|
Schedule& schedule_;
|
2021-05-19 07:51:14 -05:00
|
|
|
const SummaryState& summaryState_;
|
|
|
|
const EclipseState& eclState_;
|
2021-05-25 05:57:11 -05:00
|
|
|
const Parallel::Communication& comm_;
|
2021-05-19 07:51:14 -05:00
|
|
|
|
|
|
|
PhaseUsage phase_usage_;
|
|
|
|
bool terminal_output_{false};
|
|
|
|
bool wells_active_{false};
|
|
|
|
bool initial_step_{};
|
2021-06-07 06:04:29 -05:00
|
|
|
bool report_step_starts_{};
|
2021-05-19 07:51:14 -05:00
|
|
|
|
2021-06-07 08:01:10 -05:00
|
|
|
std::optional<int> last_run_wellpi_{};
|
|
|
|
|
2021-05-19 07:51:14 -05:00
|
|
|
std::vector<Well> wells_ecl_;
|
|
|
|
std::vector<std::vector<PerforationData>> well_perf_data_;
|
|
|
|
std::function<bool(const Well&)> not_on_process_{};
|
|
|
|
|
2021-06-07 05:15:36 -05:00
|
|
|
// a vector of all the wells.
|
|
|
|
std::vector<WellInterfaceGeneric*> well_container_generic_{};
|
|
|
|
|
2021-06-07 05:20:49 -05:00
|
|
|
std::vector<int> local_shut_wells_{};
|
|
|
|
|
2021-05-19 07:51:14 -05:00
|
|
|
std::vector<ParallelWellInfo> parallel_well_info_;
|
2021-09-27 05:00:39 -05:00
|
|
|
std::vector<std::reference_wrapper<ParallelWellInfo>> local_parallel_well_info_;
|
2021-05-19 07:51:14 -05:00
|
|
|
|
|
|
|
std::vector<WellProdIndexCalculator> prod_index_calc_;
|
|
|
|
|
|
|
|
std::vector<int> pvt_region_idx_;
|
|
|
|
|
|
|
|
mutable std::unordered_set<std::string> closed_this_step_;
|
|
|
|
|
2021-06-14 03:38:45 -05:00
|
|
|
GuideRate guideRate_;
|
2021-06-07 05:20:49 -05:00
|
|
|
std::unique_ptr<VFPProperties> vfp_properties_{};
|
2021-05-19 07:51:14 -05:00
|
|
|
std::map<std::string, double> node_pressures_; // Storing network pressures for output.
|
|
|
|
|
|
|
|
/*
|
|
|
|
The various wellState members should be accessed and modified
|
|
|
|
through the accessor functions wellState(), prevWellState(),
|
|
|
|
commitWellState(), resetWellState(), nupcolWellState() and
|
|
|
|
updateNupcolWellState().
|
|
|
|
*/
|
|
|
|
WGState active_wgstate_;
|
|
|
|
WGState last_valid_wgstate_;
|
|
|
|
WGState nupcol_wgstate_;
|
2021-06-07 06:04:29 -05:00
|
|
|
|
2021-06-07 06:04:29 -05:00
|
|
|
bool glift_debug = false;
|
|
|
|
|
2022-01-19 03:14:47 -06:00
|
|
|
double last_glift_opt_time_ = -1.0;
|
|
|
|
|
2021-10-07 06:21:35 -05:00
|
|
|
private:
|
2021-06-07 06:04:29 -05:00
|
|
|
WellInterfaceGeneric* getGenWell(const std::string& well_name);
|
2021-05-19 07:51:14 -05:00
|
|
|
};
|
|
|
|
|
|
|
|
|
|
|
|
} // namespace Opm
|
|
|
|
|
|
|
|
#endif
|