2019-02-07 07:43:17 -06:00
|
|
|
/*
|
2021-03-01 17:14:34 -06:00
|
|
|
Copyright 2016 - 2019 SINTEF Digital, Mathematics & Cybernetics.
|
|
|
|
Copyright 2016 - 2018 Equinor ASA.
|
|
|
|
Copyright 2017 Dr. Blatt - HPC-Simulation-Software & Services
|
|
|
|
Copyright 2016 - 2018 Norce AS
|
2017-02-13 09:45:06 -06:00
|
|
|
|
2021-03-01 17:14:34 -06:00
|
|
|
This file is part of the Open Porous Media project (OPM).
|
2017-02-13 09:45:06 -06:00
|
|
|
|
2021-03-01 17:14:34 -06:00
|
|
|
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.
|
2019-02-07 07:43:17 -06:00
|
|
|
|
2021-03-01 17:14:34 -06:00
|
|
|
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.
|
2019-02-07 07:43:17 -06:00
|
|
|
|
2021-03-01 17:14:34 -06:00
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
2019-02-07 07:43:17 -06:00
|
|
|
*/
|
|
|
|
|
2019-05-07 06:06:02 -05:00
|
|
|
#include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>
|
2019-11-25 03:34:50 -06:00
|
|
|
#include <opm/core/props/phaseUsageFromDeck.hpp>
|
2022-07-08 03:16:53 -05:00
|
|
|
#include <opm/grid/utility/cartesianToCompressed.hpp>
|
2021-12-14 01:30:15 -06:00
|
|
|
#include <opm/input/eclipse/Units/UnitSystem.hpp>
|
2020-10-29 17:30:09 -05:00
|
|
|
|
2021-05-31 07:31:56 -05:00
|
|
|
#include <opm/simulators/wells/VFPProperties.hpp>
|
|
|
|
|
2020-09-22 07:12:15 -05:00
|
|
|
#include <algorithm>
|
2020-10-29 17:16:31 -05:00
|
|
|
#include <utility>
|
2021-03-01 17:14:34 -06:00
|
|
|
|
2020-10-01 11:27:57 -05:00
|
|
|
#include <fmt/format.h>
|
2020-07-20 14:38:30 -05:00
|
|
|
|
2017-02-13 09:45:06 -06:00
|
|
|
namespace Opm {
|
2017-05-03 06:34:15 -05:00
|
|
|
template<typename TypeTag>
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2021-04-25 13:47:31 -05:00
|
|
|
BlackoilWellModel(Simulator& ebosSimulator, const PhaseUsage& phase_usage)
|
2021-05-19 07:51:14 -05:00
|
|
|
: BlackoilWellModelGeneric(ebosSimulator.vanguard().schedule(),
|
|
|
|
ebosSimulator.vanguard().summaryState(),
|
|
|
|
ebosSimulator.vanguard().eclState(),
|
|
|
|
phase_usage,
|
|
|
|
ebosSimulator.gridView().comm())
|
|
|
|
, ebosSimulator_(ebosSimulator)
|
2017-02-13 09:45:06 -06:00
|
|
|
{
|
2021-05-19 07:51:14 -05:00
|
|
|
terminal_output_ = ((ebosSimulator.gridView().comm().rank() == 0) &&
|
|
|
|
EWOMS_GET_PARAM(TypeTag, bool, EnableTerminalOutput));
|
2019-10-23 02:09:45 -05:00
|
|
|
|
2020-10-05 13:02:13 -05:00
|
|
|
local_num_cells_ = ebosSimulator_.gridView().size(0);
|
2021-03-01 17:14:34 -06:00
|
|
|
|
2020-10-05 05:43:55 -05:00
|
|
|
// Number of cells the global grid view
|
2020-10-05 13:02:13 -05:00
|
|
|
global_num_cells_ = ebosSimulator_.vanguard().globalNumCells();
|
2019-10-23 02:09:45 -05:00
|
|
|
|
2021-12-01 07:00:21 -06:00
|
|
|
{
|
|
|
|
auto& parallel_wells = ebosSimulator.vanguard().parallelWells();
|
2021-10-19 05:24:04 -05:00
|
|
|
|
2022-09-08 06:33:39 -05:00
|
|
|
this->parallel_well_info_.reserve(parallel_wells.size());
|
|
|
|
for( const auto& name_bool : parallel_wells) {
|
|
|
|
this->parallel_well_info_.emplace_back(name_bool, grid().comm());
|
|
|
|
}
|
2021-03-01 17:14:34 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
this->alternative_well_rate_init_ =
|
|
|
|
EWOMS_GET_PARAM(TypeTag, bool, AlternativeWellRateInit);
|
2018-08-16 04:51:36 -05:00
|
|
|
}
|
|
|
|
|
2021-04-25 13:47:31 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
BlackoilWellModel(Simulator& ebosSimulator) :
|
|
|
|
BlackoilWellModel(ebosSimulator, phaseUsageFromDeck(ebosSimulator.vanguard().eclState()))
|
|
|
|
{}
|
|
|
|
|
|
|
|
|
2018-08-16 04:51:36 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
2019-04-03 10:26:57 -05:00
|
|
|
init()
|
2018-08-16 04:51:36 -05:00
|
|
|
{
|
|
|
|
extractLegacyCellPvtRegionIndex_();
|
|
|
|
extractLegacyDepth_();
|
|
|
|
|
2017-11-08 06:57:36 -06:00
|
|
|
gravity_ = ebosSimulator_.problem().gravity()[2];
|
2017-02-13 09:45:06 -06:00
|
|
|
|
2018-03-23 06:56:19 -05:00
|
|
|
initial_step_ = true;
|
2018-08-16 04:51:36 -05:00
|
|
|
|
|
|
|
// add the eWoms auxiliary module for the wells to the list
|
|
|
|
ebosSimulator_.model().addAuxiliaryModule(this);
|
2018-11-15 07:37:01 -06:00
|
|
|
|
2020-10-05 13:02:13 -05:00
|
|
|
is_cell_perforated_.resize(local_num_cells_, false);
|
2018-08-16 04:51:36 -05:00
|
|
|
}
|
|
|
|
|
2021-06-07 07:35:34 -05:00
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
2022-04-12 01:44:52 -05:00
|
|
|
initWellContainer(const int reportStepIdx)
|
2021-06-07 07:35:34 -05:00
|
|
|
{
|
2022-04-12 01:44:52 -05:00
|
|
|
const uint64_t effective_events_mask = ScheduleEvents::WELL_STATUS_CHANGE
|
|
|
|
+ ScheduleEvents::NEW_WELL;
|
|
|
|
const auto& events = schedule()[reportStepIdx].wellgroup_events();
|
2021-06-07 07:35:34 -05:00
|
|
|
for (auto& wellPtr : this->well_container_) {
|
2022-04-12 01:44:52 -05:00
|
|
|
const bool well_opened_this_step = report_step_starts_ && events.hasEvent(wellPtr->name(), effective_events_mask);
|
2021-06-07 07:35:34 -05:00
|
|
|
wellPtr->init(&this->phase_usage_, this->depth_, this->gravity_,
|
2022-04-12 01:44:52 -05:00
|
|
|
this->local_num_cells_, this->B_avg_, well_opened_this_step);
|
2021-06-07 07:35:34 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2018-08-16 04:51:36 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
addNeighbors(std::vector<NeighborSet>& neighbors) const
|
|
|
|
{
|
|
|
|
if (!param_.matrix_add_well_contributions_) {
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
|
|
|
// Create cartesian to compressed mapping
|
2019-11-13 16:16:11 -06:00
|
|
|
const auto& schedule_wells = schedule().getWellsatEnd();
|
2018-08-16 04:51:36 -05:00
|
|
|
|
|
|
|
// initialize the additional cell connections introduced by wells.
|
2019-05-02 05:51:25 -05:00
|
|
|
for (const auto& well : schedule_wells)
|
2018-08-16 04:51:36 -05:00
|
|
|
{
|
|
|
|
std::vector<int> wellCells;
|
|
|
|
// All possible connections of the well
|
2019-05-02 05:51:25 -05:00
|
|
|
const auto& connectionSet = well.getConnections();
|
2018-08-16 04:51:36 -05:00
|
|
|
wellCells.reserve(connectionSet.size());
|
|
|
|
|
|
|
|
for ( size_t c=0; c < connectionSet.size(); c++ )
|
|
|
|
{
|
|
|
|
const auto& connection = connectionSet.get(c);
|
2022-09-05 05:42:54 -05:00
|
|
|
int compressed_idx = compressedIndexForInterior(connection.global_index());
|
2018-08-16 04:51:36 -05:00
|
|
|
if ( compressed_idx >= 0 ) { // Ignore connections in inactive/remote cells.
|
|
|
|
wellCells.push_back(compressed_idx);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
for (int cellIdx : wellCells) {
|
|
|
|
neighbors[cellIdx].insert(wellCells.begin(),
|
|
|
|
wellCells.end());
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
2019-10-02 05:48:29 -05:00
|
|
|
linearize(SparseMatrixAdapter& jacobian, GlobalEqVector& res)
|
2018-08-16 04:51:36 -05:00
|
|
|
{
|
2021-10-08 06:45:57 -05:00
|
|
|
if (!param_.matrix_add_well_contributions_)
|
|
|
|
{
|
|
|
|
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
|
|
|
{
|
|
|
|
// if the well contributions are not supposed to be included explicitly in
|
|
|
|
// the matrix, we only apply the vector part of the Schur complement here.
|
|
|
|
for (const auto& well: well_container_) {
|
|
|
|
// r = r - duneC_^T * invDuneD_ * resWell_
|
|
|
|
well->apply(res);
|
|
|
|
}
|
2019-02-25 03:51:30 -06:00
|
|
|
}
|
2021-10-08 06:45:57 -05:00
|
|
|
OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel::linearize failed: ",
|
|
|
|
ebosSimulator_.gridView().comm());
|
2019-02-25 03:51:30 -06:00
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
2018-08-16 04:51:36 -05:00
|
|
|
for (const auto& well: well_container_) {
|
2019-10-02 05:48:29 -05:00
|
|
|
well->addWellContributions(jacobian);
|
2018-08-16 04:51:36 -05:00
|
|
|
|
|
|
|
// applying the well residual to reservoir residuals
|
|
|
|
// r = r - duneC_^T * invDuneD_ * resWell_
|
|
|
|
well->apply(res);
|
|
|
|
}
|
2017-11-08 06:57:36 -06:00
|
|
|
}
|
2017-02-13 09:45:06 -06:00
|
|
|
|
|
|
|
|
2017-05-03 06:34:15 -05:00
|
|
|
template<typename TypeTag>
|
2017-02-13 09:45:06 -06:00
|
|
|
void
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2017-11-08 06:57:36 -06:00
|
|
|
beginReportStep(const int timeStepIdx)
|
2017-02-13 09:45:06 -06:00
|
|
|
{
|
2021-05-05 04:22:44 -05:00
|
|
|
DeferredLogger local_deferredLogger;
|
2021-09-15 13:06:26 -05:00
|
|
|
|
2020-01-29 01:41:41 -06:00
|
|
|
report_step_starts_ = true;
|
2019-02-07 07:43:17 -06:00
|
|
|
|
2018-02-01 09:27:42 -06:00
|
|
|
const Grid& grid = ebosSimulator_.vanguard().grid();
|
2019-05-29 00:44:23 -05:00
|
|
|
const auto& summaryState = ebosSimulator_.vanguard().summaryState();
|
2020-12-09 04:02:04 -06:00
|
|
|
// Make wells_ecl_ contain only this partition's wells.
|
2021-04-28 03:22:29 -05:00
|
|
|
wells_ecl_ = getLocalWells(timeStepIdx);
|
2021-09-27 05:00:39 -05:00
|
|
|
this->local_parallel_well_info_ = createLocalParallelWellInfo(wells_ecl_);
|
2020-10-09 06:38:33 -05:00
|
|
|
|
2021-09-15 13:06:26 -05:00
|
|
|
// at least initializeWellState might be throw
|
|
|
|
// exception in opm-material (UniformTabulated2DFunction.hpp)
|
|
|
|
// playing it safe by extending the scope a bit.
|
2021-09-20 03:56:11 -05:00
|
|
|
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
2021-09-15 13:06:26 -05:00
|
|
|
{
|
2017-11-16 04:42:34 -06:00
|
|
|
|
2021-09-15 13:06:26 -05:00
|
|
|
// The well state initialize bhp with the cell pressure in the top cell.
|
|
|
|
// We must therefore provide it with updated cell pressures
|
|
|
|
this->initializeWellPerfData();
|
|
|
|
this->initializeWellState(timeStepIdx, summaryState);
|
2017-11-08 06:57:36 -06:00
|
|
|
|
2021-09-15 13:06:26 -05:00
|
|
|
// handling MS well related
|
|
|
|
if (param_.use_multisegment_well_&& anyMSWellOpenLocal()) { // if we use MultisegmentWell model
|
|
|
|
this->wellState().initWellStateMSWell(wells_ecl_, &this->prevWellState());
|
|
|
|
}
|
2019-08-07 07:13:11 -05:00
|
|
|
|
2021-09-15 13:06:26 -05:00
|
|
|
const Group& fieldGroup = schedule().getGroup("FIELD", timeStepIdx);
|
|
|
|
WellGroupHelpers::setCmodeGroup(fieldGroup, schedule(), summaryState, timeStepIdx, this->wellState(), this->groupState());
|
2017-11-08 06:57:36 -06:00
|
|
|
|
2021-09-15 13:06:26 -05:00
|
|
|
// Compute reservoir volumes for RESV controls.
|
|
|
|
rateConverter_.reset(new RateConverterType (phase_usage_,
|
|
|
|
std::vector<int>(local_num_cells_, 0)));
|
|
|
|
rateConverter_->template defineState<ElementContext>(ebosSimulator_);
|
2021-09-01 08:04:54 -05:00
|
|
|
|
2021-09-15 13:06:26 -05:00
|
|
|
// Compute regional average pressures used by gpmaint
|
|
|
|
if (schedule_[timeStepIdx].has_gpmaint()) {
|
|
|
|
const auto& fp = this->eclState_.fieldProps();
|
|
|
|
const auto& fipnum = fp.get_int("FIPNUM");
|
|
|
|
regionalAveragePressureCalculator_.reset(new AverageRegionalPressureType (phase_usage_,fipnum));
|
2021-01-14 12:22:34 -06:00
|
|
|
}
|
2020-10-29 17:30:09 -05:00
|
|
|
|
2021-09-15 13:06:26 -05:00
|
|
|
{
|
|
|
|
const auto& sched_state = this->schedule()[timeStepIdx];
|
|
|
|
// update VFP properties
|
2022-02-22 08:25:10 -06:00
|
|
|
vfp_properties_.reset(new VFPProperties( sched_state.vfpinj(), sched_state.vfpprod(), this->prevWellState()));
|
2021-09-15 13:06:26 -05:00
|
|
|
this->initializeWellProdIndCalculators();
|
|
|
|
if (sched_state.events().hasEvent(ScheduleEvents::Events::WELL_PRODUCTIVITY_INDEX)) {
|
|
|
|
this->runWellPIScaling(timeStepIdx, local_deferredLogger);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
2021-09-20 03:56:11 -05:00
|
|
|
OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger, "beginReportStep() failed: ",
|
2021-05-25 05:57:11 -05:00
|
|
|
terminal_output_, grid.comm());
|
2021-03-26 16:10:16 -05:00
|
|
|
// Store the current well state, to be able to recover in the case of failed iterations
|
2021-04-22 10:31:21 -05:00
|
|
|
this->commitWGState();
|
2018-06-06 08:17:59 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// called at the beginning of a time step
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
2021-03-01 08:17:52 -06:00
|
|
|
beginTimeStep()
|
|
|
|
{
|
2020-05-15 04:21:32 -05:00
|
|
|
updatePerforationIntensiveQuantities();
|
2021-03-18 08:49:52 -05:00
|
|
|
updateAverageFormationFactor();
|
2021-05-05 04:22:44 -05:00
|
|
|
DeferredLogger local_deferredLogger;
|
2019-02-07 07:43:17 -06:00
|
|
|
|
2021-04-22 10:31:21 -05:00
|
|
|
this->resetWGState();
|
2018-08-16 04:51:36 -05:00
|
|
|
const int reportStepIdx = ebosSimulator_.episodeIndex();
|
2021-05-19 07:51:14 -05:00
|
|
|
updateAndCommunicateGroupData(reportStepIdx,
|
|
|
|
ebosSimulator_.model().newtonMethod().numIterations());
|
|
|
|
this->wellState().gliftTimeStepInit();
|
2018-08-16 04:51:36 -05:00
|
|
|
const double simulationTime = ebosSimulator_.time();
|
2021-09-20 03:56:11 -05:00
|
|
|
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
|
|
|
{
|
2019-02-07 07:43:17 -06:00
|
|
|
// test wells
|
|
|
|
wellTesting(reportStepIdx, simulationTime, local_deferredLogger);
|
2018-06-06 08:17:59 -05:00
|
|
|
|
2019-02-07 07:43:17 -06:00
|
|
|
// create the well container
|
2021-06-07 07:49:41 -05:00
|
|
|
createWellContainer(reportStepIdx);
|
2017-11-08 06:57:36 -06:00
|
|
|
|
2022-09-13 03:11:10 -05:00
|
|
|
// Wells are active if they are active wells on at least one process.
|
|
|
|
const Grid& grid = ebosSimulator_.vanguard().grid();
|
|
|
|
wells_active_ = !this->well_container_.empty();
|
|
|
|
wells_active_ = grid.comm().max(wells_active_);
|
|
|
|
|
2019-02-07 07:43:17 -06:00
|
|
|
// do the initialization for all the wells
|
|
|
|
// TODO: to see whether we can postpone of the intialization of the well containers to
|
|
|
|
// optimize the usage of the following several member variables
|
2022-04-12 01:44:52 -05:00
|
|
|
this->initWellContainer(reportStepIdx);
|
2017-08-07 04:35:59 -05:00
|
|
|
|
2019-02-07 07:43:17 -06:00
|
|
|
// update the updated cell flag
|
|
|
|
std::fill(is_cell_perforated_.begin(), is_cell_perforated_.end(), false);
|
|
|
|
for (auto& well : well_container_) {
|
|
|
|
well->updatePerforatedCell(is_cell_perforated_);
|
|
|
|
}
|
2018-11-15 07:37:01 -06:00
|
|
|
|
2019-02-07 07:43:17 -06:00
|
|
|
// calculate the efficiency factors for each well
|
2019-08-07 07:13:11 -05:00
|
|
|
calculateEfficiencyFactors(reportStepIdx);
|
2017-08-07 04:35:59 -05:00
|
|
|
|
2021-01-28 07:33:21 -06:00
|
|
|
if constexpr (has_polymer_)
|
2019-02-07 07:43:17 -06:00
|
|
|
{
|
2020-08-27 02:13:30 -05:00
|
|
|
if (PolymerModule::hasPlyshlog() || getPropValue<TypeTag, Properties::EnablePolymerMW>() ) {
|
2021-01-28 05:52:19 -06:00
|
|
|
setRepRadiusPerfLength();
|
2019-02-07 07:43:17 -06:00
|
|
|
}
|
2017-11-08 06:57:36 -06:00
|
|
|
}
|
2021-05-25 05:57:11 -05:00
|
|
|
|
2021-03-18 02:29:23 -05:00
|
|
|
}
|
2021-09-20 03:56:11 -05:00
|
|
|
OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger, "beginTimeStep() failed: ",
|
2021-05-25 05:57:11 -05:00
|
|
|
terminal_output_, ebosSimulator_.vanguard().grid().comm());
|
2019-02-07 07:43:17 -06:00
|
|
|
|
2018-04-07 14:41:34 -05:00
|
|
|
for (auto& well : well_container_) {
|
|
|
|
well->setVFPProperties(vfp_properties_.get());
|
2021-06-14 03:38:45 -05:00
|
|
|
well->setGuideRate(&guideRate_);
|
2018-04-07 14:41:34 -05:00
|
|
|
}
|
2017-08-07 04:35:59 -05:00
|
|
|
|
2018-10-31 09:32:50 -05:00
|
|
|
// Close completions due to economical reasons
|
2018-06-06 08:17:59 -05:00
|
|
|
for (auto& well : well_container_) {
|
2021-09-21 02:30:02 -05:00
|
|
|
well->closeCompletions(wellTestState());
|
2018-06-06 08:17:59 -05:00
|
|
|
}
|
2019-09-30 05:49:36 -05:00
|
|
|
|
2020-10-15 10:56:11 -05:00
|
|
|
if (alternative_well_rate_init_) {
|
2020-12-07 03:05:46 -06:00
|
|
|
// Update the well rates of well_state_, if only single-phase rates, to
|
|
|
|
// have proper multi-phase rates proportional to rates at bhp zero.
|
|
|
|
// This is done only for producers, as injectors will only have a single
|
|
|
|
// nonzero phase anyway.
|
2020-10-15 10:56:11 -05:00
|
|
|
for (auto& well : well_container_) {
|
2020-12-07 03:05:46 -06:00
|
|
|
if (well->isProducer()) {
|
2021-03-26 16:10:16 -05:00
|
|
|
well->updateWellStateRates(ebosSimulator_, this->wellState(), local_deferredLogger);
|
2020-12-07 03:05:46 -06:00
|
|
|
}
|
2020-10-15 10:56:11 -05:00
|
|
|
}
|
2020-05-15 04:21:32 -05:00
|
|
|
}
|
|
|
|
|
2022-03-23 06:29:40 -05:00
|
|
|
// calculate the well potentials
|
|
|
|
try {
|
|
|
|
updateWellPotentials(reportStepIdx,
|
|
|
|
/*onlyAfterEvent*/true,
|
|
|
|
ebosSimulator_.vanguard().summaryConfig(),
|
|
|
|
local_deferredLogger);
|
|
|
|
} catch ( std::runtime_error& e ) {
|
|
|
|
const std::string msg = "A zero well potential is returned for output purposes. ";
|
|
|
|
local_deferredLogger.warning("WELL_POTENTIAL_CALCULATION_FAILED", msg);
|
|
|
|
}
|
|
|
|
|
2021-03-19 05:09:14 -05:00
|
|
|
//update guide rates
|
2019-12-12 02:22:37 -06:00
|
|
|
const auto& comm = ebosSimulator_.vanguard().grid().comm();
|
2021-03-19 05:09:14 -05:00
|
|
|
const auto& summaryState = ebosSimulator_.vanguard().summaryState();
|
|
|
|
std::vector<double> pot(numPhases(), 0.0);
|
|
|
|
const Group& fieldGroup = schedule().getGroup("FIELD", reportStepIdx);
|
2021-06-28 16:12:18 -05:00
|
|
|
WellGroupHelpers::updateGuideRates(fieldGroup, schedule(), summaryState, this->phase_usage_, reportStepIdx, simulationTime,
|
|
|
|
this->wellState(), this->groupState(), comm, &this->guideRate_, pot, local_deferredLogger);
|
2021-09-20 03:56:11 -05:00
|
|
|
std::string exc_msg;
|
|
|
|
auto exc_type = ExceptionType::NONE;
|
2021-09-21 07:03:35 -05:00
|
|
|
// update gpmaint targets
|
|
|
|
if (schedule_[reportStepIdx].has_gpmaint()) {
|
|
|
|
regionalAveragePressureCalculator_->template defineState<ElementContext>(ebosSimulator_);
|
|
|
|
const double dt = ebosSimulator_.timeStepSize();
|
|
|
|
WellGroupHelpers::updateGpMaintTargetForGroups(fieldGroup,
|
|
|
|
schedule_, *regionalAveragePressureCalculator_, reportStepIdx, dt, this->wellState(), this->groupState());
|
|
|
|
}
|
2021-01-19 07:43:32 -06:00
|
|
|
try {
|
2021-03-09 06:37:03 -06:00
|
|
|
// Compute initial well solution for new wells and injectors that change injection type i.e. WAG.
|
2021-01-19 07:43:32 -06:00
|
|
|
for (auto& well : well_container_) {
|
|
|
|
const uint64_t effective_events_mask = ScheduleEvents::WELL_STATUS_CHANGE
|
2021-03-09 06:37:03 -06:00
|
|
|
+ ScheduleEvents::INJECTION_TYPE_CHANGED
|
2021-04-16 06:38:56 -05:00
|
|
|
+ ScheduleEvents::WELL_SWITCHED_INJECTOR_PRODUCER
|
2021-03-09 06:37:03 -06:00
|
|
|
+ ScheduleEvents::NEW_WELL;
|
2021-01-19 07:43:32 -06:00
|
|
|
|
|
|
|
const auto& events = schedule()[reportStepIdx].wellgroup_events();
|
|
|
|
const bool event = report_step_starts_ && events.hasEvent(well->name(), effective_events_mask);
|
2021-09-21 03:32:56 -05:00
|
|
|
const bool dyn_status_change = this->wellState().well(well->name()).status
|
|
|
|
!= this->prevWellState().well(well->name()).status;
|
|
|
|
|
|
|
|
if (event || dyn_status_change) {
|
2021-01-19 07:43:32 -06:00
|
|
|
try {
|
2021-06-10 08:09:05 -05:00
|
|
|
well->updateWellStateWithTarget(ebosSimulator_, this->groupState(), this->wellState(), local_deferredLogger);
|
2021-03-26 16:10:16 -05:00
|
|
|
well->calculateExplicitQuantities(ebosSimulator_, this->wellState(), local_deferredLogger);
|
2021-04-22 10:31:21 -05:00
|
|
|
well->solveWellEquation(ebosSimulator_, this->wellState(), this->groupState(), local_deferredLogger);
|
2021-03-18 02:29:23 -05:00
|
|
|
} catch (const std::exception& e) {
|
2021-01-19 07:43:32 -06:00
|
|
|
const std::string msg = "Compute initial well solution for new well " + well->name() + " failed. Continue with zero initial rates";
|
|
|
|
local_deferredLogger.warning("WELL_INITIAL_SOLVE_FAILED", msg);
|
|
|
|
}
|
|
|
|
}
|
2021-01-04 07:00:59 -06:00
|
|
|
}
|
2021-03-18 02:29:23 -05:00
|
|
|
}
|
2021-09-20 03:56:11 -05:00
|
|
|
// Catch clauses for all errors setting exc_type and exc_msg
|
|
|
|
OPM_PARALLEL_CATCH_CLAUSE(exc_type, exc_msg);
|
2021-03-18 02:29:23 -05:00
|
|
|
|
|
|
|
if (exc_type != ExceptionType::NONE) {
|
2021-01-19 07:43:32 -06:00
|
|
|
const std::string msg = "Compute initial well solution for new wells failed. Continue with zero initial rates";
|
|
|
|
local_deferredLogger.warning("WELL_INITIAL_SOLVE_FAILED", msg);
|
2021-01-04 07:00:59 -06:00
|
|
|
}
|
|
|
|
|
2020-09-30 03:04:39 -05:00
|
|
|
logAndCheckForExceptionsAndThrow(local_deferredLogger,
|
2021-05-25 05:57:11 -05:00
|
|
|
exc_type, "beginTimeStep() failed: " + exc_msg, terminal_output_, comm);
|
2020-09-30 03:04:39 -05:00
|
|
|
|
2017-08-10 08:27:05 -05:00
|
|
|
}
|
|
|
|
|
2017-11-08 08:48:30 -06:00
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
2021-03-01 08:17:52 -06:00
|
|
|
BlackoilWellModel<TypeTag>::wellTesting(const int timeStepIdx,
|
|
|
|
const double simulationTime,
|
2021-05-05 04:22:44 -05:00
|
|
|
DeferredLogger& deferred_logger)
|
2021-03-01 08:17:52 -06:00
|
|
|
{
|
2021-01-12 07:08:41 -06:00
|
|
|
const auto& wtest_config = schedule()[timeStepIdx].wtest_config();
|
2021-09-29 08:05:23 -05:00
|
|
|
if (!wtest_config.empty()) { // there is a WTEST request
|
2021-08-24 02:00:06 -05:00
|
|
|
const std::vector<std::string> wellsForTesting = wellTestState()
|
2021-10-08 08:23:03 -05:00
|
|
|
.test_wells(wtest_config, simulationTime);
|
2021-08-24 02:00:06 -05:00
|
|
|
for (const std::string& well_name : wellsForTesting) {
|
2021-11-26 03:39:07 -06:00
|
|
|
|
|
|
|
const Well& wellEcl = schedule().getWell(well_name, timeStepIdx);
|
|
|
|
if (wellEcl.getStatus() == Well::Status::SHUT)
|
2021-10-08 08:17:43 -05:00
|
|
|
continue;
|
2018-06-06 08:17:59 -05:00
|
|
|
|
2021-10-08 08:17:43 -05:00
|
|
|
WellInterfacePtr well = createWellForWellTest(well_name, timeStepIdx, deferred_logger);
|
2019-02-03 01:13:11 -06:00
|
|
|
// some preparation before the well can be used
|
2022-04-12 01:44:52 -05:00
|
|
|
well->init(&phase_usage_, depth_, gravity_, local_num_cells_, B_avg_, true);
|
2021-11-26 03:39:07 -06:00
|
|
|
|
2019-08-07 07:13:11 -05:00
|
|
|
double well_efficiency_factor = wellEcl.getEfficiencyFactor();
|
2021-03-01 08:17:52 -06:00
|
|
|
WellGroupHelpers::accumulateGroupEfficiencyFactor(schedule().getGroup(wellEcl.groupName(), timeStepIdx),
|
|
|
|
schedule(), timeStepIdx, well_efficiency_factor);
|
|
|
|
|
2019-02-03 01:13:11 -06:00
|
|
|
well->setWellEfficiencyFactor(well_efficiency_factor);
|
|
|
|
well->setVFPProperties(vfp_properties_.get());
|
2021-06-14 03:38:45 -05:00
|
|
|
well->setGuideRate(&guideRate_);
|
2017-11-08 08:48:30 -06:00
|
|
|
|
2021-08-24 02:00:06 -05:00
|
|
|
well->wellTesting(ebosSimulator_, simulationTime, this->wellState(), this->groupState(), wellTestState(), deferred_logger);
|
2019-02-03 01:13:11 -06:00
|
|
|
}
|
2019-01-18 07:04:30 -06:00
|
|
|
}
|
2017-11-08 08:48:30 -06:00
|
|
|
}
|
|
|
|
|
2021-03-01 08:17:52 -06:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2017-11-08 08:48:30 -06:00
|
|
|
// called at the end of a report step
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
2021-03-01 08:17:52 -06:00
|
|
|
endReportStep()
|
|
|
|
{
|
2020-11-27 11:26:59 -06:00
|
|
|
// Clear the communication data structures for above values.
|
2021-09-27 05:00:39 -05:00
|
|
|
for (auto&& pinfo : this->local_parallel_well_info_)
|
2020-11-27 11:26:59 -06:00
|
|
|
{
|
2021-09-27 05:00:39 -05:00
|
|
|
pinfo.get().clear();
|
2020-11-27 11:26:59 -06:00
|
|
|
}
|
2017-11-08 08:48:30 -06:00
|
|
|
}
|
|
|
|
|
2021-03-01 08:17:52 -06:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2017-11-08 08:48:30 -06:00
|
|
|
// called at the end of a report step
|
|
|
|
template<typename TypeTag>
|
2020-05-07 09:13:39 -05:00
|
|
|
const SimulatorReportSingle&
|
2017-11-08 08:48:30 -06:00
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
lastReport() const {return last_report_; }
|
|
|
|
|
2021-03-01 08:17:52 -06:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2017-11-08 08:48:30 -06:00
|
|
|
// called at the end of a time step
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
2021-03-01 08:17:52 -06:00
|
|
|
timeStepSucceeded(const double& simulationTime, const double dt)
|
|
|
|
{
|
2021-03-01 08:43:43 -06:00
|
|
|
this->closed_this_step_.clear();
|
2019-02-03 01:13:11 -06:00
|
|
|
|
2020-01-29 01:41:41 -06:00
|
|
|
// time step is finished and we are not any more at the beginning of an report step
|
|
|
|
report_step_starts_ = false;
|
2020-04-24 08:25:38 -05:00
|
|
|
const int reportStepIdx = ebosSimulator_.episodeIndex();
|
2020-01-29 01:41:41 -06:00
|
|
|
|
2021-05-05 04:22:44 -05:00
|
|
|
DeferredLogger local_deferredLogger;
|
2018-02-14 06:34:35 -06:00
|
|
|
for (const auto& well : well_container_) {
|
2020-08-27 02:13:30 -05:00
|
|
|
if (getPropValue<TypeTag, Properties::EnablePolymerMW>() && well->isInjector()) {
|
2021-03-26 16:10:16 -05:00
|
|
|
well->updateWaterThroughput(dt, this->wellState());
|
2018-11-30 05:15:51 -06:00
|
|
|
}
|
2018-02-14 06:34:35 -06:00
|
|
|
}
|
2021-10-08 02:47:22 -05:00
|
|
|
// report well switching
|
|
|
|
for (const auto& well : well_container_) {
|
|
|
|
well->reportWellSwitching(this->wellState().well(well->indexOfWell()), local_deferredLogger);
|
|
|
|
}
|
2018-11-07 07:53:43 -06:00
|
|
|
|
2019-08-07 07:13:11 -05:00
|
|
|
// update the rate converter with current averages pressures etc in
|
|
|
|
rateConverter_->template defineState<ElementContext>(ebosSimulator_);
|
|
|
|
|
|
|
|
// calculate the well potentials
|
2019-02-07 07:43:17 -06:00
|
|
|
try {
|
2021-06-07 06:04:29 -05:00
|
|
|
updateWellPotentials(reportStepIdx,
|
|
|
|
/*onlyAfterEvent*/false,
|
|
|
|
ebosSimulator_.vanguard().summaryConfig(),
|
|
|
|
local_deferredLogger);
|
2019-02-07 07:43:17 -06:00
|
|
|
} catch ( std::runtime_error& e ) {
|
2018-11-07 07:53:43 -06:00
|
|
|
const std::string msg = "A zero well potential is returned for output purposes. ";
|
2019-02-03 01:13:11 -06:00
|
|
|
local_deferredLogger.warning("WELL_POTENTIAL_CALCULATION_FAILED", msg);
|
2018-11-07 07:53:43 -06:00
|
|
|
}
|
2020-05-04 08:56:34 -05:00
|
|
|
|
2021-09-21 02:30:02 -05:00
|
|
|
updateWellTestState(simulationTime, wellTestState());
|
2021-08-09 06:21:10 -05:00
|
|
|
|
2020-05-04 08:56:34 -05:00
|
|
|
// check group sales limits at the end of the timestep
|
2021-09-21 07:03:35 -05:00
|
|
|
const Group& fieldGroup = schedule_.getGroup("FIELD", reportStepIdx);
|
2021-05-19 07:51:14 -05:00
|
|
|
checkGconsaleLimits(fieldGroup, this->wellState(),
|
|
|
|
ebosSimulator_.episodeIndex(), local_deferredLogger);
|
2020-05-04 08:56:34 -05:00
|
|
|
|
2020-10-09 06:38:33 -05:00
|
|
|
this->calculateProductivityIndexValues(local_deferredLogger);
|
|
|
|
|
2021-04-22 10:31:21 -05:00
|
|
|
this->commitWGState();
|
2021-05-25 05:57:11 -05:00
|
|
|
|
|
|
|
const Opm::Parallel::Communication& comm = grid().comm();
|
|
|
|
DeferredLogger global_deferredLogger = gatherDeferredLogger(local_deferredLogger, comm);
|
2019-02-03 01:13:11 -06:00
|
|
|
if (terminal_output_) {
|
|
|
|
global_deferredLogger.logMessages();
|
|
|
|
}
|
2021-03-08 08:11:50 -06:00
|
|
|
|
|
|
|
//reporting output temperatures
|
|
|
|
this->computeWellTemperature();
|
2017-11-08 08:48:30 -06:00
|
|
|
}
|
2017-08-10 08:27:05 -05:00
|
|
|
|
2022-08-10 03:01:54 -05:00
|
|
|
|
2022-06-17 04:15:15 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
computeTotalRatesForDof(RateVector& rate,
|
2022-08-10 03:01:54 -05:00
|
|
|
unsigned elemIdx) const
|
2022-06-17 04:15:15 -05:00
|
|
|
{
|
|
|
|
rate = 0;
|
2022-08-10 03:01:54 -05:00
|
|
|
|
2022-06-17 04:15:15 -05:00
|
|
|
if (!is_cell_perforated_[elemIdx])
|
|
|
|
return;
|
|
|
|
|
|
|
|
for (const auto& well : well_container_)
|
|
|
|
well->addCellRates(rate, elemIdx);
|
|
|
|
}
|
2018-08-16 04:51:36 -05:00
|
|
|
|
2022-08-10 03:01:54 -05:00
|
|
|
|
2018-08-16 04:51:36 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
template <class Context>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
computeTotalRatesForDof(RateVector& rate,
|
|
|
|
const Context& context,
|
|
|
|
unsigned spaceIdx,
|
|
|
|
unsigned timeIdx) const
|
|
|
|
{
|
|
|
|
rate = 0;
|
|
|
|
int elemIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
|
2018-11-15 07:37:01 -06:00
|
|
|
|
|
|
|
if (!is_cell_perforated_[elemIdx])
|
|
|
|
return;
|
|
|
|
|
2018-08-16 04:51:36 -05:00
|
|
|
for (const auto& well : well_container_)
|
|
|
|
well->addCellRates(rate, elemIdx);
|
|
|
|
}
|
|
|
|
|
2018-11-14 06:18:48 -06:00
|
|
|
|
|
|
|
|
2020-10-19 17:16:52 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
initializeWellState(const int timeStepIdx,
|
|
|
|
const SummaryState& summaryState)
|
|
|
|
{
|
|
|
|
std::vector<double> cellPressures(this->local_num_cells_, 0.0);
|
|
|
|
ElementContext elemCtx(ebosSimulator_);
|
|
|
|
|
|
|
|
const auto& gridView = ebosSimulator_.vanguard().gridView();
|
|
|
|
const auto& elemEndIt = gridView.template end</*codim=*/0>();
|
2021-09-20 04:12:27 -05:00
|
|
|
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
|
|
|
|
2020-10-19 17:16:52 -05:00
|
|
|
for (auto elemIt = gridView.template begin</*codim=*/0>();
|
|
|
|
elemIt != elemEndIt;
|
|
|
|
++elemIt)
|
|
|
|
{
|
|
|
|
if (elemIt->partitionType() != Dune::InteriorEntity) {
|
|
|
|
continue;
|
|
|
|
}
|
|
|
|
|
|
|
|
elemCtx.updatePrimaryStencil(*elemIt);
|
|
|
|
elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
|
|
|
|
|
|
|
|
const auto& fs = elemCtx.intensiveQuantities(/*spaceIdx=*/0, /*timeIdx=*/0).fluidState();
|
|
|
|
// copy of get perfpressure in Standard well except for value
|
|
|
|
double& perf_pressure = cellPressures[elemCtx.globalSpaceIndex(/*spaceIdx=*/0, /*timeIdx=*/0)];
|
|
|
|
if (Indices::oilEnabled) {
|
|
|
|
perf_pressure = fs.pressure(FluidSystem::oilPhaseIdx).value();
|
|
|
|
} else if (Indices::waterEnabled) {
|
|
|
|
perf_pressure = fs.pressure(FluidSystem::waterPhaseIdx).value();
|
|
|
|
} else {
|
|
|
|
perf_pressure = fs.pressure(FluidSystem::gasPhaseIdx).value();
|
|
|
|
}
|
|
|
|
}
|
2021-05-25 05:57:11 -05:00
|
|
|
OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel::initializeWellState() failed: ", ebosSimulator_.vanguard().grid().comm());
|
2020-10-19 17:16:52 -05:00
|
|
|
|
2021-03-26 16:10:16 -05:00
|
|
|
this->wellState().init(cellPressures, schedule(), wells_ecl_, local_parallel_well_info_, timeStepIdx,
|
2021-04-25 13:47:31 -05:00
|
|
|
&this->prevWellState(), well_perf_data_,
|
2021-04-28 03:22:29 -05:00
|
|
|
summaryState);
|
2020-10-19 17:16:52 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2017-06-15 10:19:49 -05:00
|
|
|
template<typename TypeTag>
|
2021-06-07 07:49:41 -05:00
|
|
|
void
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2019-10-23 02:09:45 -05:00
|
|
|
createWellContainer(const int time_step)
|
2017-06-15 10:19:49 -05:00
|
|
|
{
|
2021-05-05 04:22:44 -05:00
|
|
|
DeferredLogger local_deferredLogger;
|
2019-12-13 04:08:36 -06:00
|
|
|
|
2019-11-01 09:11:21 -05:00
|
|
|
const int nw = numLocalWells();
|
2017-06-15 10:19:49 -05:00
|
|
|
|
2021-06-07 07:49:41 -05:00
|
|
|
well_container_.clear();
|
|
|
|
|
2017-08-21 03:23:42 -05:00
|
|
|
if (nw > 0) {
|
2021-06-07 07:49:41 -05:00
|
|
|
well_container_.reserve(nw);
|
2021-03-01 08:17:52 -06:00
|
|
|
|
2017-06-15 10:19:49 -05:00
|
|
|
for (int w = 0; w < nw; ++w) {
|
2019-10-23 02:09:45 -05:00
|
|
|
const Well& well_ecl = wells_ecl_[w];
|
2021-06-02 04:02:46 -05:00
|
|
|
|
|
|
|
if (well_ecl.getConnections().empty()) {
|
|
|
|
// No connections in this well. Nothing to do.
|
|
|
|
continue;
|
|
|
|
}
|
|
|
|
|
2019-10-23 02:09:45 -05:00
|
|
|
const std::string& well_name = well_ecl.name();
|
2021-03-01 08:43:43 -06:00
|
|
|
const auto well_status = this->schedule()
|
|
|
|
.getWell(well_name, time_step).getStatus();
|
|
|
|
|
|
|
|
if ((well_ecl.getStatus() == Well::Status::SHUT) ||
|
|
|
|
(well_status == Well::Status::SHUT))
|
|
|
|
{
|
|
|
|
// Due to ACTIONX the well might have been closed behind our back.
|
|
|
|
if (well_ecl.getStatus() != Well::Status::SHUT) {
|
|
|
|
this->closed_this_step_.insert(well_name);
|
2021-03-26 16:10:16 -05:00
|
|
|
this->wellState().shutWell(w);
|
2021-03-01 08:43:43 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
continue;
|
|
|
|
}
|
2017-06-15 10:19:49 -05:00
|
|
|
|
2019-07-31 09:15:41 -05:00
|
|
|
// A new WCON keywords can re-open a well that was closed/shut due to Physical limit
|
2021-10-07 08:18:53 -05:00
|
|
|
if (this->wellTestState().well_is_closed(well_name)) {
|
2019-07-31 09:15:41 -05:00
|
|
|
// TODO: more checking here, to make sure this standard more specific and complete
|
|
|
|
// maybe there is some WCON keywords will not open the well
|
2021-08-04 03:08:31 -05:00
|
|
|
auto& events = this->wellState().well(w).events;
|
2021-05-20 06:32:18 -05:00
|
|
|
if (events.hasEvent(WellState::event_mask)) {
|
2021-09-21 02:30:02 -05:00
|
|
|
if (wellTestState().lastTestTime(well_name) == ebosSimulator_.time()) {
|
2019-07-31 09:15:41 -05:00
|
|
|
// The well was shut this timestep, we are most likely retrying
|
|
|
|
// a timestep without the well in question, after it caused
|
|
|
|
// repeated timestep cuts. It should therefore not be opened,
|
|
|
|
// even if it was new or received new targets this report step.
|
2021-05-20 06:32:18 -05:00
|
|
|
events.clearEvent(WellState::event_mask);
|
2019-07-31 09:15:41 -05:00
|
|
|
} else {
|
2021-10-07 08:30:10 -05:00
|
|
|
wellTestState().open_well(well_name);
|
2021-10-07 07:53:06 -05:00
|
|
|
wellTestState().open_completions(well_name);
|
2018-12-14 03:04:59 -06:00
|
|
|
}
|
2018-11-17 16:36:31 -06:00
|
|
|
}
|
2019-07-31 09:15:41 -05:00
|
|
|
}
|
2018-11-17 16:30:27 -06:00
|
|
|
|
2019-07-31 09:15:41 -05:00
|
|
|
// TODO: should we do this for all kinds of closing reasons?
|
2021-09-21 02:30:02 -05:00
|
|
|
// something like wellTestState().hasWell(well_name)?
|
2019-08-07 07:13:11 -05:00
|
|
|
bool wellIsStopped = false;
|
2021-10-07 08:18:53 -05:00
|
|
|
if (wellTestState().well_is_closed(well_name))
|
2021-03-01 08:17:52 -06:00
|
|
|
{
|
|
|
|
if (well_ecl.getAutomaticShutIn()) {
|
2019-07-31 09:15:41 -05:00
|
|
|
// shut wells are not added to the well container
|
2021-03-26 16:10:16 -05:00
|
|
|
this->wellState().shutWell(w);
|
2019-07-31 09:15:41 -05:00
|
|
|
continue;
|
|
|
|
} else {
|
2022-06-07 03:01:02 -05:00
|
|
|
if (!well_ecl.getAllowCrossFlow()) {
|
|
|
|
// stopped wells where cross flow is not allowed
|
|
|
|
// are not added to the well container
|
|
|
|
this->wellState().shutWell(w);
|
|
|
|
continue;
|
|
|
|
}
|
2019-08-07 07:13:11 -05:00
|
|
|
// stopped wells are added to the container but marked as stopped
|
2021-03-26 16:10:16 -05:00
|
|
|
this->wellState().stopWell(w);
|
2019-08-07 07:13:11 -05:00
|
|
|
wellIsStopped = true;
|
2018-06-06 08:17:59 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2019-12-13 04:08:36 -06:00
|
|
|
// If a production well disallows crossflow and its
|
|
|
|
// (prediction type) rate control is zero, then it is effectively shut.
|
|
|
|
if (!well_ecl.getAllowCrossFlow() && well_ecl.isProducer() && well_ecl.predictionMode()) {
|
|
|
|
const auto& summaryState = ebosSimulator_.vanguard().summaryState();
|
2021-03-01 08:17:52 -06:00
|
|
|
const auto prod_controls = well_ecl.productionControls(summaryState);
|
|
|
|
|
|
|
|
auto is_zero = [](const double x)
|
|
|
|
{
|
|
|
|
return std::isfinite(x) && !std::isnormal(x);
|
|
|
|
};
|
|
|
|
|
2019-12-13 04:08:36 -06:00
|
|
|
bool zero_rate_control = false;
|
|
|
|
switch (prod_controls.cmode) {
|
|
|
|
case Well::ProducerCMode::ORAT:
|
2021-03-01 08:17:52 -06:00
|
|
|
zero_rate_control = is_zero(prod_controls.oil_rate);
|
2019-12-13 04:08:36 -06:00
|
|
|
break;
|
2021-03-01 08:17:52 -06:00
|
|
|
|
2019-12-13 04:08:36 -06:00
|
|
|
case Well::ProducerCMode::WRAT:
|
2021-03-01 08:17:52 -06:00
|
|
|
zero_rate_control = is_zero(prod_controls.water_rate);
|
2019-12-13 04:08:36 -06:00
|
|
|
break;
|
2021-03-01 08:17:52 -06:00
|
|
|
|
2019-12-13 04:08:36 -06:00
|
|
|
case Well::ProducerCMode::GRAT:
|
2021-03-01 08:17:52 -06:00
|
|
|
zero_rate_control = is_zero(prod_controls.gas_rate);
|
2019-12-13 04:08:36 -06:00
|
|
|
break;
|
2021-03-01 08:17:52 -06:00
|
|
|
|
2019-12-13 04:08:36 -06:00
|
|
|
case Well::ProducerCMode::LRAT:
|
2021-03-01 08:17:52 -06:00
|
|
|
zero_rate_control = is_zero(prod_controls.liquid_rate);
|
2019-12-13 04:08:36 -06:00
|
|
|
break;
|
2021-03-01 08:17:52 -06:00
|
|
|
|
2019-12-13 04:08:36 -06:00
|
|
|
case Well::ProducerCMode::RESV:
|
2021-03-01 08:17:52 -06:00
|
|
|
zero_rate_control = is_zero(prod_controls.resv_rate);
|
2019-12-13 04:08:36 -06:00
|
|
|
break;
|
|
|
|
default:
|
|
|
|
// Might still have zero rate controls, but is pressure controlled.
|
|
|
|
zero_rate_control = false;
|
2021-03-01 08:17:52 -06:00
|
|
|
break;
|
2019-12-13 04:08:36 -06:00
|
|
|
}
|
2021-03-01 08:17:52 -06:00
|
|
|
|
2019-12-13 04:08:36 -06:00
|
|
|
if (zero_rate_control) {
|
|
|
|
// Treat as shut, do not add to container.
|
|
|
|
local_deferredLogger.info(" Well shut due to zero rate control and disallowing crossflow: " + well_ecl.name());
|
2021-03-26 16:10:16 -05:00
|
|
|
this->wellState().shutWell(w);
|
2019-12-13 04:08:36 -06:00
|
|
|
continue;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2019-12-06 06:12:36 -06:00
|
|
|
if (well_status == Well::Status::STOP) {
|
2021-03-26 16:10:16 -05:00
|
|
|
this->wellState().stopWell(w);
|
2019-12-06 06:12:36 -06:00
|
|
|
wellIsStopped = true;
|
|
|
|
}
|
|
|
|
|
2021-06-07 07:49:41 -05:00
|
|
|
well_container_.emplace_back(this->createWellPointer(w, time_step));
|
2020-12-01 11:04:46 -06:00
|
|
|
|
2019-08-07 07:13:11 -05:00
|
|
|
if (wellIsStopped)
|
2021-06-07 07:49:41 -05:00
|
|
|
well_container_.back()->stopWell();
|
2017-06-15 10:19:49 -05:00
|
|
|
}
|
|
|
|
}
|
2019-02-07 07:43:17 -06:00
|
|
|
|
2019-12-13 04:08:36 -06:00
|
|
|
// Collect log messages and print.
|
2021-05-25 05:57:11 -05:00
|
|
|
|
|
|
|
const Opm::Parallel::Communication& comm = grid().comm();
|
|
|
|
DeferredLogger global_deferredLogger = gatherDeferredLogger(local_deferredLogger, comm);
|
2019-12-13 04:08:36 -06:00
|
|
|
if (terminal_output_) {
|
|
|
|
global_deferredLogger.logMessages();
|
|
|
|
}
|
|
|
|
|
2021-06-07 05:15:36 -05:00
|
|
|
well_container_generic_.clear();
|
2021-06-07 07:49:41 -05:00
|
|
|
for (auto& w : well_container_)
|
2021-06-07 05:15:36 -05:00
|
|
|
well_container_generic_.push_back(w.get());
|
2021-03-01 18:11:19 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2020-12-01 11:04:46 -06:00
|
|
|
template <typename TypeTag>
|
|
|
|
typename BlackoilWellModel<TypeTag>::WellInterfacePtr
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
createWellPointer(const int wellID, const int time_step) const
|
|
|
|
{
|
|
|
|
const auto is_multiseg = this->wells_ecl_[wellID].isMultiSegment();
|
|
|
|
|
|
|
|
if (! (this->param_.use_multisegment_well_ && is_multiseg)) {
|
|
|
|
return this->template createTypedWellPointer<StandardWell<TypeTag>>(wellID, time_step);
|
|
|
|
}
|
|
|
|
else {
|
|
|
|
return this->template createTypedWellPointer<MultisegmentWell<TypeTag>>(wellID, time_step);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <typename TypeTag>
|
|
|
|
template <typename WellType>
|
|
|
|
std::unique_ptr<WellType>
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
createTypedWellPointer(const int wellID, const int time_step) const
|
|
|
|
{
|
|
|
|
// Use the pvtRegionIdx from the top cell
|
|
|
|
const auto& perf_data = this->well_perf_data_[wellID];
|
|
|
|
|
2020-11-27 11:22:23 -06:00
|
|
|
// Cater for case where local part might have no perforations.
|
2021-03-01 08:17:52 -06:00
|
|
|
const auto pvtreg = perf_data.empty()
|
|
|
|
? 0 : pvt_region_idx_[perf_data.front().cell_index];
|
|
|
|
|
2021-09-27 05:00:39 -05:00
|
|
|
const auto& parallel_well_info = this->local_parallel_well_info_[wellID].get();
|
2021-03-01 08:17:52 -06:00
|
|
|
const auto global_pvtreg = parallel_well_info.broadcastFirstPerforationValue(pvtreg);
|
2020-11-27 11:22:23 -06:00
|
|
|
|
2020-12-01 11:04:46 -06:00
|
|
|
return std::make_unique<WellType>(this->wells_ecl_[wellID],
|
2020-11-27 11:22:23 -06:00
|
|
|
parallel_well_info,
|
2020-12-01 11:04:46 -06:00
|
|
|
time_step,
|
|
|
|
this->param_,
|
|
|
|
*this->rateConverter_,
|
2020-11-27 11:22:23 -06:00
|
|
|
global_pvtreg,
|
2020-12-01 11:04:46 -06:00
|
|
|
this->numComponents(),
|
|
|
|
this->numPhases(),
|
|
|
|
wellID,
|
|
|
|
perf_data);
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2018-10-31 08:56:56 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
typename BlackoilWellModel<TypeTag>::WellInterfacePtr
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
createWellForWellTest(const std::string& well_name,
|
2019-02-07 07:43:17 -06:00
|
|
|
const int report_step,
|
2021-05-05 04:22:44 -05:00
|
|
|
DeferredLogger& deferred_logger) const
|
2018-10-31 08:56:56 -05:00
|
|
|
{
|
|
|
|
// Finding the location of the well in wells_ecl
|
|
|
|
const int nw_wells_ecl = wells_ecl_.size();
|
|
|
|
int index_well_ecl = 0;
|
|
|
|
for (; index_well_ecl < nw_wells_ecl; ++index_well_ecl) {
|
2019-05-02 05:51:25 -05:00
|
|
|
if (well_name == wells_ecl_[index_well_ecl].name()) {
|
2018-10-31 08:56:56 -05:00
|
|
|
break;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
// It should be able to find in wells_ecl.
|
|
|
|
if (index_well_ecl == nw_wells_ecl) {
|
2019-02-07 07:43:17 -06:00
|
|
|
OPM_DEFLOG_THROW(std::logic_error, "Could not find well " << well_name << " in wells_ecl ", deferred_logger);
|
2018-10-31 08:56:56 -05:00
|
|
|
}
|
|
|
|
|
2020-12-01 11:04:46 -06:00
|
|
|
return this->createWellPointer(index_well_ecl, report_step);
|
2018-10-31 08:56:56 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2017-05-03 06:34:15 -05:00
|
|
|
template<typename TypeTag>
|
2017-11-08 06:57:36 -06:00
|
|
|
void
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2017-11-08 06:57:36 -06:00
|
|
|
assemble(const int iterationIdx,
|
2018-06-29 02:12:14 -05:00
|
|
|
const double dt)
|
2017-02-13 09:45:06 -06:00
|
|
|
{
|
2017-11-08 06:57:36 -06:00
|
|
|
|
2021-05-05 04:22:44 -05:00
|
|
|
DeferredLogger local_deferredLogger;
|
2020-09-30 03:04:39 -05:00
|
|
|
if (this->glift_debug) {
|
2020-10-01 11:27:57 -05:00
|
|
|
const std::string msg = fmt::format(
|
|
|
|
"assemble() : iteration {}" , iterationIdx);
|
|
|
|
gliftDebug(msg, local_deferredLogger);
|
2020-09-30 03:04:39 -05:00
|
|
|
}
|
2020-05-07 09:13:39 -05:00
|
|
|
last_report_ = SimulatorReportSingle();
|
2020-09-07 08:05:02 -05:00
|
|
|
Dune::Timer perfTimer;
|
|
|
|
perfTimer.start();
|
2017-11-08 06:57:36 -06:00
|
|
|
|
2017-11-08 08:48:30 -06:00
|
|
|
if ( ! wellsActive() ) {
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
2019-02-03 01:13:11 -06:00
|
|
|
|
2017-11-08 06:57:36 -06:00
|
|
|
updatePerforationIntensiveQuantities();
|
2017-02-13 09:45:06 -06:00
|
|
|
|
2021-09-23 13:15:19 -05:00
|
|
|
if (iterationIdx == 0) {
|
|
|
|
// try-catch is needed here as updateWellControls
|
|
|
|
// contains global communication and has either to
|
|
|
|
// be reached by all processes or all need to abort
|
|
|
|
// before.
|
|
|
|
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
|
|
|
{
|
2019-02-07 07:43:17 -06:00
|
|
|
calculateExplicitQuantities(local_deferredLogger);
|
2019-05-24 09:45:27 -05:00
|
|
|
prepareTimeStep(local_deferredLogger);
|
2019-02-07 07:43:17 -06:00
|
|
|
}
|
2021-09-23 13:15:19 -05:00
|
|
|
OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger, "assemble() failed (It=0): ",
|
2021-05-25 05:57:11 -05:00
|
|
|
terminal_output_, grid().comm());
|
2021-09-23 13:15:19 -05:00
|
|
|
}
|
|
|
|
updateWellControls(local_deferredLogger, /* check group controls */ true);
|
2019-11-06 09:16:19 -06:00
|
|
|
|
2022-02-08 08:46:08 -06:00
|
|
|
bool alq_updated = false;
|
2021-09-23 13:15:19 -05:00
|
|
|
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
|
|
|
{
|
2019-02-07 07:43:17 -06:00
|
|
|
// Set the well primary variables based on the value of well solutions
|
|
|
|
initPrimaryVariablesEvaluation();
|
2017-02-13 09:45:06 -06:00
|
|
|
|
2022-02-08 08:46:08 -06:00
|
|
|
alq_updated = maybeDoGasLiftOptimize(local_deferredLogger);
|
2021-03-18 08:49:52 -05:00
|
|
|
assembleWellEq(dt, local_deferredLogger);
|
2021-03-18 02:29:23 -05:00
|
|
|
}
|
2021-09-20 03:56:11 -05:00
|
|
|
OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger, "assemble() failed: ",
|
2021-05-25 05:57:11 -05:00
|
|
|
terminal_output_, grid().comm());
|
2021-11-18 05:57:16 -06:00
|
|
|
|
|
|
|
//update guide rates
|
2022-01-17 06:01:42 -06:00
|
|
|
const int reportStepIdx = ebosSimulator_.episodeIndex();
|
2022-02-08 08:46:08 -06:00
|
|
|
if (alq_updated || guideRateUpdateIsNeeded(reportStepIdx)) {
|
2021-11-18 05:57:16 -06:00
|
|
|
const double simulationTime = ebosSimulator_.time();
|
|
|
|
const auto& comm = ebosSimulator_.vanguard().grid().comm();
|
|
|
|
const auto& summaryState = ebosSimulator_.vanguard().summaryState();
|
|
|
|
std::vector<double> pot(numPhases(), 0.0);
|
|
|
|
const Group& fieldGroup = schedule().getGroup("FIELD", reportStepIdx);
|
|
|
|
WellGroupHelpers::updateGuideRates(fieldGroup, schedule(), summaryState, this->phase_usage_, reportStepIdx, simulationTime,
|
|
|
|
this->wellState(), this->groupState(), comm, &this->guideRate_, pot, local_deferredLogger);
|
|
|
|
}
|
2017-11-08 06:57:36 -06:00
|
|
|
last_report_.converged = true;
|
2020-09-07 08:05:02 -05:00
|
|
|
last_report_.assemble_time_well += perfTimer.stop();
|
2017-02-13 09:45:06 -06:00
|
|
|
}
|
|
|
|
|
2021-03-25 16:40:12 -05:00
|
|
|
template<typename TypeTag>
|
2022-02-08 08:46:08 -06:00
|
|
|
bool
|
2021-03-25 16:40:12 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2021-05-05 04:22:44 -05:00
|
|
|
maybeDoGasLiftOptimize(DeferredLogger& deferred_logger)
|
2021-03-03 06:59:26 -06:00
|
|
|
{
|
2022-01-19 03:14:47 -06:00
|
|
|
bool do_glift_optimization = false;
|
2022-02-08 08:46:08 -06:00
|
|
|
int num_wells_changed = 0;
|
2022-01-19 03:14:47 -06:00
|
|
|
const double simulation_time = ebosSimulator_.time();
|
|
|
|
const double min_wait = ebosSimulator_.vanguard().schedule().glo(ebosSimulator_.episodeIndex()).min_wait();
|
|
|
|
// We only optimize if a min_wait time has past.
|
|
|
|
// If all_newton is true we still want to optimize several times pr timestep
|
|
|
|
// i.e. we also optimize if check simulation_time == last_glift_opt_time_
|
|
|
|
// that is when the last_glift_opt_time is already updated with the current time step
|
|
|
|
if ( simulation_time == last_glift_opt_time_ || simulation_time >= (last_glift_opt_time_ + min_wait)) {
|
|
|
|
do_glift_optimization = true;
|
|
|
|
last_glift_opt_time_ = simulation_time;
|
|
|
|
}
|
|
|
|
|
2021-09-27 07:04:50 -05:00
|
|
|
if (do_glift_optimization) {
|
2021-05-27 01:31:49 -05:00
|
|
|
GLiftOptWells glift_wells;
|
|
|
|
GLiftProdWells prod_wells;
|
|
|
|
GLiftWellStateMap state_map;
|
|
|
|
// NOTE: To make GasLiftGroupInfo (see below) independent of the TypeTag
|
|
|
|
// associated with *this (i.e. BlackoilWellModel<TypeTag>) we observe
|
|
|
|
// that GasLiftGroupInfo's only dependence on *this is that it needs to
|
|
|
|
// access the eclipse Wells in the well container (the eclipse Wells
|
|
|
|
// themselves are independent of the TypeTag).
|
|
|
|
// Hence, we extract them from the well container such that we can pass
|
|
|
|
// them to the GasLiftGroupInfo constructor.
|
|
|
|
GLiftEclWells ecl_well_map;
|
|
|
|
initGliftEclWellMap(ecl_well_map);
|
|
|
|
GasLiftGroupInfo group_info {
|
|
|
|
ecl_well_map,
|
|
|
|
ebosSimulator_.vanguard().schedule(),
|
|
|
|
ebosSimulator_.vanguard().summaryState(),
|
|
|
|
ebosSimulator_.episodeIndex(),
|
|
|
|
ebosSimulator_.model().newtonMethod().numIterations(),
|
|
|
|
phase_usage_,
|
|
|
|
deferred_logger,
|
2021-06-22 02:52:22 -05:00
|
|
|
this->wellState(),
|
Improve debugging tools in gaslift code.
Introduces a gaslift debugging variable in ALQState in WellState. This
variable will persist between timesteps in contrast to when debugging
variables are defined in GasLiftSingleWell, GasLiftGroupState, or GasLiftStage2.
Currently only an integer variable debug_counter is added to ALQState,
which can be used as follows: First debugging is switched on globally
for BlackOilWellModel, GasLiftSingleWell, GasLiftGroupState, and
GasLiftStage2 by setting glift_debug to a true value in BlackOilWellModelGeneric.
Then, the following debugging code can be added to e.g. one of
GasLiftSingleWell, GasLiftGroupState, or GasLiftStage2 :
auto count = debugUpdateGlobalCounter_();
if (count == some_integer) {
displayDebugMessage_("stop here");
}
Here, the integer "some_integer" is determined typically by looking at
the debugging output of a previous run. This can be done since the
call to debugUpdateGlobalCounter_() will print out the current value
of the counter and then increment the counter by one. And it will be
easy to recognize these values in the debug ouput. If you find a place
in the output that looks suspect, just take a note of the counter
value in the output around that point and insert the value for
"some_integer", then after recompiling the code with the desired value
for "some_integer", it is now easy to set a breakpoint in GDB at the
line
displayDebugMessage_("stop here").
shown in the above snippet. This should improve the ability to quickly
to set a breakpoint in GDB around at a given time and point in the simulation.
2022-01-23 13:37:26 -06:00
|
|
|
ebosSimulator_.vanguard().grid().comm(),
|
|
|
|
this->glift_debug
|
2021-05-27 01:31:49 -05:00
|
|
|
};
|
2021-06-22 02:52:22 -05:00
|
|
|
group_info.initialize();
|
2021-05-27 01:31:49 -05:00
|
|
|
gasLiftOptimizationStage1(
|
|
|
|
deferred_logger, prod_wells, glift_wells, group_info, state_map);
|
|
|
|
gasLiftOptimizationStage2(
|
|
|
|
deferred_logger, prod_wells, glift_wells, state_map,
|
|
|
|
ebosSimulator_.episodeIndex());
|
|
|
|
if (this->glift_debug) gliftDebugShowALQ(deferred_logger);
|
2022-02-08 08:46:08 -06:00
|
|
|
num_wells_changed = glift_wells.size();
|
2021-05-27 01:31:49 -05:00
|
|
|
}
|
2022-03-24 06:42:46 -05:00
|
|
|
num_wells_changed = this->comm_.sum(num_wells_changed);
|
2022-02-08 08:46:08 -06:00
|
|
|
return num_wells_changed > 0;
|
2021-05-27 01:31:49 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
gasLiftOptimizationStage1(DeferredLogger& deferred_logger,
|
|
|
|
GLiftProdWells &prod_wells, GLiftOptWells &glift_wells,
|
|
|
|
GasLiftGroupInfo &group_info, GLiftWellStateMap &state_map)
|
|
|
|
{
|
|
|
|
auto comm = ebosSimulator_.vanguard().grid().comm();
|
2021-06-17 03:44:32 -05:00
|
|
|
int num_procs = comm.size();
|
2021-05-27 01:31:49 -05:00
|
|
|
// NOTE: Gas lift optimization stage 1 seems to be difficult
|
|
|
|
// to do in parallel since the wells are optimized on different
|
|
|
|
// processes and each process needs to know the current ALQ allocated
|
|
|
|
// to each group it is a memeber of in order to check group limits and avoid
|
|
|
|
// allocating more ALQ than necessary. (Surplus ALQ is removed in
|
|
|
|
// stage 2). In stage1, as each well is adding ALQ, the current group ALQ needs
|
|
|
|
// to be communicated to the other processes. But there is no common
|
|
|
|
// synchronization point that all process will reach in the
|
|
|
|
// runOptimizeLoop_() in GasLiftSingleWell.cpp.
|
|
|
|
//
|
|
|
|
// TODO: Maybe a better solution could be invented by distributing
|
|
|
|
// wells according to certain parent groups. Then updated group rates
|
|
|
|
// might not have to be communicated to the other processors.
|
|
|
|
|
|
|
|
// Currently, the best option seems to be to run this part sequentially
|
|
|
|
// (not in parallel).
|
|
|
|
//
|
|
|
|
// TODO: The simplest approach seems to be if a) one process could take
|
|
|
|
// ownership of all the wells (the union of all the wells in the
|
|
|
|
// well_container_ of each process) then this process could do the
|
|
|
|
// optimization, while the other processes could wait for it to
|
|
|
|
// finish (e.g. comm.barrier()), or alternatively, b) if all
|
|
|
|
// processes could take ownership of all the wells. Then there
|
|
|
|
// would be no need for synchronization here..
|
|
|
|
//
|
2021-06-17 03:44:32 -05:00
|
|
|
for (int i = 0; i< num_procs; i++) {
|
2021-05-27 01:31:49 -05:00
|
|
|
int num_rates_to_sync = 0; // communication variable
|
|
|
|
GLiftSyncGroups groups_to_sync;
|
2021-06-17 03:44:32 -05:00
|
|
|
if (comm.rank() == i) {
|
2021-05-27 01:31:49 -05:00
|
|
|
// Run stage1: Optimize single wells while also checking group limits
|
|
|
|
for (const auto& well : well_container_) {
|
|
|
|
// NOTE: Only the wells in "group_info" needs to be optimized
|
|
|
|
if (group_info.hasWell(well->name())) {
|
2022-02-07 04:28:35 -06:00
|
|
|
gasLiftOptimizationStage1SingleWell(
|
|
|
|
well.get(), deferred_logger, prod_wells, glift_wells,
|
|
|
|
group_info, state_map, groups_to_sync
|
Improve debugging tools in gaslift code.
Introduces a gaslift debugging variable in ALQState in WellState. This
variable will persist between timesteps in contrast to when debugging
variables are defined in GasLiftSingleWell, GasLiftGroupState, or GasLiftStage2.
Currently only an integer variable debug_counter is added to ALQState,
which can be used as follows: First debugging is switched on globally
for BlackOilWellModel, GasLiftSingleWell, GasLiftGroupState, and
GasLiftStage2 by setting glift_debug to a true value in BlackOilWellModelGeneric.
Then, the following debugging code can be added to e.g. one of
GasLiftSingleWell, GasLiftGroupState, or GasLiftStage2 :
auto count = debugUpdateGlobalCounter_();
if (count == some_integer) {
displayDebugMessage_("stop here");
}
Here, the integer "some_integer" is determined typically by looking at
the debugging output of a previous run. This can be done since the
call to debugUpdateGlobalCounter_() will print out the current value
of the counter and then increment the counter by one. And it will be
easy to recognize these values in the debug ouput. If you find a place
in the output that looks suspect, just take a note of the counter
value in the output around that point and insert the value for
"some_integer", then after recompiling the code with the desired value
for "some_integer", it is now easy to set a breakpoint in GDB at the
line
displayDebugMessage_("stop here").
shown in the above snippet. This should improve the ability to quickly
to set a breakpoint in GDB around at a given time and point in the simulation.
2022-01-23 13:37:26 -06:00
|
|
|
);
|
2021-05-27 01:31:49 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
num_rates_to_sync = groups_to_sync.size();
|
|
|
|
}
|
|
|
|
// Since "group_info" is not used in stage2, there is no need to
|
|
|
|
// communicate rates if this is the last iteration...
|
|
|
|
if (i == (num_procs - 1))
|
|
|
|
break;
|
|
|
|
num_rates_to_sync = comm.sum(num_rates_to_sync);
|
|
|
|
if (num_rates_to_sync > 0) {
|
|
|
|
std::vector<int> group_indexes;
|
|
|
|
group_indexes.reserve(num_rates_to_sync);
|
|
|
|
std::vector<double> group_alq_rates;
|
|
|
|
group_alq_rates.reserve(num_rates_to_sync);
|
|
|
|
std::vector<double> group_oil_rates;
|
|
|
|
group_oil_rates.reserve(num_rates_to_sync);
|
|
|
|
std::vector<double> group_gas_rates;
|
|
|
|
group_gas_rates.reserve(num_rates_to_sync);
|
2021-11-04 07:12:05 -05:00
|
|
|
std::vector<double> group_water_rates;
|
|
|
|
group_water_rates.reserve(num_rates_to_sync);
|
2021-05-27 01:31:49 -05:00
|
|
|
if (comm.rank() == i) {
|
|
|
|
for (auto idx : groups_to_sync) {
|
2021-11-04 07:12:05 -05:00
|
|
|
auto [oil_rate, gas_rate, water_rate, alq] = group_info.getRates(idx);
|
2021-05-27 01:31:49 -05:00
|
|
|
group_indexes.push_back(idx);
|
|
|
|
group_oil_rates.push_back(oil_rate);
|
|
|
|
group_gas_rates.push_back(gas_rate);
|
2021-11-04 07:12:05 -05:00
|
|
|
group_water_rates.push_back(water_rate);
|
2021-05-27 01:31:49 -05:00
|
|
|
group_alq_rates.push_back(alq);
|
|
|
|
}
|
2021-11-09 05:48:48 -06:00
|
|
|
} else {
|
|
|
|
group_indexes.resize(num_rates_to_sync);
|
|
|
|
group_oil_rates.resize(num_rates_to_sync);
|
|
|
|
group_gas_rates.resize(num_rates_to_sync);
|
|
|
|
group_water_rates.resize(num_rates_to_sync);
|
|
|
|
group_alq_rates.resize(num_rates_to_sync);
|
2021-05-27 01:31:49 -05:00
|
|
|
}
|
|
|
|
// TODO: We only need to broadcast to processors with index
|
|
|
|
// j > i since we do not use the "group_info" in stage 2. In
|
|
|
|
// this case we should use comm.send() and comm.receive()
|
|
|
|
// instead of comm.broadcast() to communicate with specific
|
|
|
|
// processes, and these processes only need to receive the
|
|
|
|
// data if they are going to check the group rates in stage1
|
|
|
|
// Another similar idea is to only communicate the rates to
|
|
|
|
// process j = i + 1
|
2021-11-09 05:48:48 -06:00
|
|
|
Mpi::broadcast(comm, i, group_indexes, group_oil_rates,
|
|
|
|
group_gas_rates, group_water_rates, group_alq_rates);
|
2021-05-27 01:31:49 -05:00
|
|
|
if (comm.rank() != i) {
|
|
|
|
for (int j=0; j<num_rates_to_sync; j++) {
|
|
|
|
group_info.updateRate(group_indexes[j],
|
2021-11-04 07:12:05 -05:00
|
|
|
group_oil_rates[j], group_gas_rates[j], group_water_rates[j], group_alq_rates[j]);
|
2021-05-27 01:31:49 -05:00
|
|
|
}
|
|
|
|
}
|
Improve debugging tools in gaslift code.
Introduces a gaslift debugging variable in ALQState in WellState. This
variable will persist between timesteps in contrast to when debugging
variables are defined in GasLiftSingleWell, GasLiftGroupState, or GasLiftStage2.
Currently only an integer variable debug_counter is added to ALQState,
which can be used as follows: First debugging is switched on globally
for BlackOilWellModel, GasLiftSingleWell, GasLiftGroupState, and
GasLiftStage2 by setting glift_debug to a true value in BlackOilWellModelGeneric.
Then, the following debugging code can be added to e.g. one of
GasLiftSingleWell, GasLiftGroupState, or GasLiftStage2 :
auto count = debugUpdateGlobalCounter_();
if (count == some_integer) {
displayDebugMessage_("stop here");
}
Here, the integer "some_integer" is determined typically by looking at
the debugging output of a previous run. This can be done since the
call to debugUpdateGlobalCounter_() will print out the current value
of the counter and then increment the counter by one. And it will be
easy to recognize these values in the debug ouput. If you find a place
in the output that looks suspect, just take a note of the counter
value in the output around that point and insert the value for
"some_integer", then after recompiling the code with the desired value
for "some_integer", it is now easy to set a breakpoint in GDB at the
line
displayDebugMessage_("stop here").
shown in the above snippet. This should improve the ability to quickly
to set a breakpoint in GDB around at a given time and point in the simulation.
2022-01-23 13:37:26 -06:00
|
|
|
if (this->glift_debug) {
|
|
|
|
int counter = 0;
|
|
|
|
if (comm.rank() == i) {
|
|
|
|
counter = this->wellState().gliftGetDebugCounter();
|
|
|
|
}
|
|
|
|
counter = comm.sum(counter);
|
|
|
|
if (comm.rank() != i) {
|
|
|
|
this->wellState().gliftSetDebugCounter(counter);
|
|
|
|
}
|
|
|
|
}
|
2021-05-27 01:31:49 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2022-02-07 04:28:35 -06:00
|
|
|
// NOTE: this method cannot be const since it passes this->wellState()
|
|
|
|
// (see below) to the GasLiftSingleWell constructor which accepts WellState
|
|
|
|
// as a non-const reference..
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
gasLiftOptimizationStage1SingleWell(WellInterface<TypeTag> *well,
|
|
|
|
DeferredLogger& deferred_logger,
|
|
|
|
GLiftProdWells &prod_wells, GLiftOptWells &glift_wells,
|
|
|
|
GasLiftGroupInfo &group_info, GLiftWellStateMap &state_map,
|
|
|
|
GLiftSyncGroups& sync_groups)
|
|
|
|
{
|
|
|
|
const auto& summary_state = ebosSimulator_.vanguard().summaryState();
|
|
|
|
std::unique_ptr<GasLiftSingleWell> glift
|
|
|
|
= std::make_unique<GasLiftSingleWell>(
|
|
|
|
*well, ebosSimulator_, summary_state,
|
|
|
|
deferred_logger, this->wellState(), this->groupState(),
|
2022-03-24 06:42:46 -05:00
|
|
|
group_info, sync_groups, this->comm_, this->glift_debug);
|
2022-02-07 04:28:35 -06:00
|
|
|
auto state = glift->runOptimize(
|
|
|
|
ebosSimulator_.model().newtonMethod().numIterations());
|
|
|
|
if (state) {
|
|
|
|
state_map.insert({well->name(), std::move(state)});
|
|
|
|
glift_wells.insert({well->name(), std::move(glift)});
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
prod_wells.insert({well->name(), well});
|
|
|
|
}
|
|
|
|
|
|
|
|
|
2021-05-27 01:31:49 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
initGliftEclWellMap(GLiftEclWells &ecl_well_map)
|
|
|
|
{
|
|
|
|
for ( const auto& well: well_container_ ) {
|
|
|
|
ecl_well_map.try_emplace(
|
|
|
|
well->name(), &(well->wellEcl()), well->indexOfWell());
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2021-03-03 06:59:26 -06:00
|
|
|
|
2017-05-03 06:34:15 -05:00
|
|
|
template<typename TypeTag>
|
2017-02-13 09:45:06 -06:00
|
|
|
void
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2021-05-05 04:22:44 -05:00
|
|
|
assembleWellEq(const double dt, DeferredLogger& deferred_logger)
|
2017-02-13 09:45:06 -06:00
|
|
|
{
|
2018-06-06 08:17:59 -05:00
|
|
|
for (auto& well : well_container_) {
|
2021-04-22 10:31:21 -05:00
|
|
|
well->assembleWellEq(ebosSimulator_, dt, this->wellState(), this->groupState(), deferred_logger);
|
2017-07-21 07:21:17 -05:00
|
|
|
}
|
2017-02-13 10:07:34 -06:00
|
|
|
}
|
|
|
|
|
2018-11-05 08:54:48 -06:00
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
apply( BVector& r) const
|
|
|
|
{
|
|
|
|
if ( ! localWellsActive() ) {
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
|
|
|
for (auto& well : well_container_) {
|
|
|
|
well->apply(r);
|
|
|
|
}
|
|
|
|
}
|
2017-02-13 10:07:34 -06:00
|
|
|
|
|
|
|
|
2017-07-21 07:21:17 -05:00
|
|
|
// Ax = A x - C D^-1 B x
|
2017-05-03 06:34:15 -05:00
|
|
|
template<typename TypeTag>
|
2017-02-13 10:07:34 -06:00
|
|
|
void
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2017-06-07 07:23:43 -05:00
|
|
|
apply(const BVector& x, BVector& Ax) const
|
2017-02-13 10:07:34 -06:00
|
|
|
{
|
2017-07-21 08:30:34 -05:00
|
|
|
// TODO: do we still need localWellsActive()?
|
2018-02-26 08:47:25 -06:00
|
|
|
if ( ! localWellsActive() ) {
|
2017-02-13 10:07:34 -06:00
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
2017-07-21 07:21:17 -05:00
|
|
|
for (auto& well : well_container_) {
|
2018-03-02 13:47:04 -06:00
|
|
|
well->apply(x, Ax);
|
2017-07-21 07:21:17 -05:00
|
|
|
}
|
2017-02-13 10:07:34 -06:00
|
|
|
}
|
|
|
|
|
2020-06-25 11:44:49 -05:00
|
|
|
#if HAVE_CUDA || HAVE_OPENCL
|
2020-03-13 08:21:59 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
getWellContributions(WellContributions& wellContribs) const
|
|
|
|
{
|
2020-05-15 09:00:09 -05:00
|
|
|
// prepare for StandardWells
|
2021-09-06 05:47:21 -05:00
|
|
|
wellContribs.setBlockSize(StandardWell<TypeTag>::Indices::numEq, StandardWell<TypeTag>::numStaticWellEq);
|
2020-09-24 14:34:46 -05:00
|
|
|
|
2020-03-13 08:21:59 -05:00
|
|
|
for(unsigned int i = 0; i < well_container_.size(); i++){
|
|
|
|
auto& well = well_container_[i];
|
|
|
|
std::shared_ptr<StandardWell<TypeTag> > derived = std::dynamic_pointer_cast<StandardWell<TypeTag> >(well);
|
2020-05-15 09:00:09 -05:00
|
|
|
if (derived) {
|
|
|
|
unsigned int numBlocks;
|
|
|
|
derived->getNumBlocks(numBlocks);
|
|
|
|
wellContribs.addNumBlocks(numBlocks);
|
|
|
|
}
|
2020-03-13 08:21:59 -05:00
|
|
|
}
|
2020-05-15 09:00:09 -05:00
|
|
|
|
|
|
|
// allocate memory for data from StandardWells
|
2020-03-18 09:08:48 -05:00
|
|
|
wellContribs.alloc();
|
2020-05-15 09:00:09 -05:00
|
|
|
|
2020-03-13 08:21:59 -05:00
|
|
|
for(unsigned int i = 0; i < well_container_.size(); i++){
|
|
|
|
auto& well = well_container_[i];
|
2020-05-15 09:00:09 -05:00
|
|
|
// maybe WellInterface could implement addWellContribution()
|
|
|
|
auto derived_std = std::dynamic_pointer_cast<StandardWell<TypeTag> >(well);
|
|
|
|
if (derived_std) {
|
|
|
|
derived_std->addWellContribution(wellContribs);
|
2020-03-19 10:06:49 -05:00
|
|
|
} else {
|
2020-05-15 09:00:09 -05:00
|
|
|
auto derived_ms = std::dynamic_pointer_cast<MultisegmentWell<TypeTag> >(well);
|
|
|
|
if (derived_ms) {
|
|
|
|
derived_ms->addWellContribution(wellContribs);
|
|
|
|
} else {
|
|
|
|
OpmLog::warning("Warning unknown type of well");
|
|
|
|
}
|
2020-03-19 10:06:49 -05:00
|
|
|
}
|
2020-03-13 08:21:59 -05:00
|
|
|
}
|
|
|
|
}
|
2020-03-18 11:48:28 -05:00
|
|
|
#endif
|
2017-02-13 10:07:34 -06:00
|
|
|
|
2017-07-21 07:21:17 -05:00
|
|
|
// Ax = Ax - alpha * C D^-1 B x
|
2017-05-03 06:34:15 -05:00
|
|
|
template<typename TypeTag>
|
2017-02-13 10:07:34 -06:00
|
|
|
void
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2017-06-07 07:23:43 -05:00
|
|
|
applyScaleAdd(const Scalar alpha, const BVector& x, BVector& Ax) const
|
2017-02-13 10:07:34 -06:00
|
|
|
{
|
2018-02-26 08:47:25 -06:00
|
|
|
if ( ! localWellsActive() ) {
|
2017-02-13 10:07:34 -06:00
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
|
|
|
if( scaleAddRes_.size() != Ax.size() ) {
|
|
|
|
scaleAddRes_.resize( Ax.size() );
|
|
|
|
}
|
|
|
|
|
|
|
|
scaleAddRes_ = 0.0;
|
2017-07-21 07:21:17 -05:00
|
|
|
// scaleAddRes_ = - C D^-1 B x
|
2017-02-13 10:07:34 -06:00
|
|
|
apply( x, scaleAddRes_ );
|
2017-07-21 07:21:17 -05:00
|
|
|
// Ax = Ax + alpha * scaleAddRes_
|
2017-02-13 10:07:34 -06:00
|
|
|
Ax.axpy( alpha, scaleAddRes_ );
|
|
|
|
}
|
|
|
|
|
2022-05-05 07:21:28 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
addWellContributions(SparseMatrixAdapter& jacobian) const
|
|
|
|
{
|
|
|
|
for ( const auto& well: well_container_ ) {
|
|
|
|
well->addWellContributions(jacobian);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
addWellPressureEquations(PressureMatrix& jacobian, const BVector& weights,const bool use_well_weights) const
|
|
|
|
{
|
|
|
|
int nw = this->numLocalWellsEnd();
|
|
|
|
int rdofs = local_num_cells_;
|
2022-06-02 06:33:23 -05:00
|
|
|
for ( int i = 0; i < nw; i++ ){
|
2022-05-05 07:21:28 -05:00
|
|
|
int wdof = rdofs + i;
|
|
|
|
jacobian[wdof][wdof] = 1.0;// better scaling ?
|
|
|
|
}
|
2022-06-10 04:08:24 -05:00
|
|
|
|
2022-06-02 06:33:23 -05:00
|
|
|
for ( const auto& well : well_container_ ) {
|
2022-05-05 07:21:28 -05:00
|
|
|
well->addWellPressureEquations(jacobian, weights, pressureVarIndex, use_well_weights, this->wellState());
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
int
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
numLocalWellsEnd() const
|
|
|
|
{
|
|
|
|
auto w = schedule().getWellsatEnd();
|
|
|
|
w.erase(std::remove_if(w.begin(), w.end(), not_on_process_), w.end());
|
|
|
|
return w.size();
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
std::vector<std::vector<int>>
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
getMaxWellConnections() const
|
|
|
|
{
|
|
|
|
std::vector<std::vector<int>> wells;
|
2017-02-13 10:07:34 -06:00
|
|
|
|
2022-05-05 07:21:28 -05:00
|
|
|
auto schedule_wells = schedule().getWellsatEnd();
|
|
|
|
schedule_wells.erase(std::remove_if(schedule_wells.begin(), schedule_wells.end(), not_on_process_), schedule_wells.end());
|
|
|
|
wells.reserve(schedule_wells.size());
|
|
|
|
|
|
|
|
// initialize the additional cell connections introduced by wells.
|
|
|
|
for ( const auto& well : schedule_wells )
|
|
|
|
{
|
|
|
|
std::vector<int> compressed_well_perforations;
|
|
|
|
// All possible completions of the well
|
|
|
|
const auto& completionSet = well.getConnections();
|
|
|
|
compressed_well_perforations.reserve(completionSet.size());
|
|
|
|
|
2022-07-08 03:16:53 -05:00
|
|
|
for (const auto& connection: well.getConnections())
|
2022-05-05 07:21:28 -05:00
|
|
|
{
|
2022-08-17 07:13:13 -05:00
|
|
|
const int compressed_idx = compressedIndexForInterior(connection.global_index());
|
2022-05-05 07:21:28 -05:00
|
|
|
if ( compressed_idx >= 0 ) // Ignore completions in inactive/remote cells.
|
|
|
|
{
|
|
|
|
compressed_well_perforations.push_back(compressed_idx);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2022-06-02 06:33:23 -05:00
|
|
|
// also include wells with no perforations in case
|
|
|
|
std::sort(compressed_well_perforations.begin(),
|
|
|
|
compressed_well_perforations.end());
|
2022-05-05 07:21:28 -05:00
|
|
|
|
2022-06-02 06:33:23 -05:00
|
|
|
wells.push_back(compressed_well_perforations);
|
2022-05-05 07:21:28 -05:00
|
|
|
}
|
|
|
|
return wells;
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
addWellPressureEquationsStruct(PressureMatrix& jacobian) const
|
|
|
|
{
|
|
|
|
int nw = this->numLocalWellsEnd();
|
|
|
|
int rdofs = local_num_cells_;
|
|
|
|
for(int i=0; i < nw; i++){
|
|
|
|
int wdof = rdofs + i;
|
|
|
|
jacobian.entry(wdof,wdof) = 1.0;// better scaling ?
|
|
|
|
}
|
|
|
|
std::vector<std::vector<int>> wellconnections = getMaxWellConnections();
|
|
|
|
for(int i=0; i < nw; i++){
|
|
|
|
const auto& perfcells = wellconnections[i];
|
|
|
|
for(int perfcell : perfcells){
|
|
|
|
int wdof = rdofs + i;
|
|
|
|
jacobian.entry(wdof,perfcell) = 0.0;
|
|
|
|
jacobian.entry(perfcell, wdof) = 0.0;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
int
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
numLocalNonshutWells() const
|
|
|
|
{
|
|
|
|
return well_container_.size();
|
|
|
|
}
|
|
|
|
|
2017-05-03 06:34:15 -05:00
|
|
|
template<typename TypeTag>
|
2017-02-13 10:07:34 -06:00
|
|
|
void
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2017-11-08 06:57:36 -06:00
|
|
|
recoverWellSolutionAndUpdateWellState(const BVector& x)
|
2017-02-13 10:07:34 -06:00
|
|
|
{
|
2021-05-25 05:57:11 -05:00
|
|
|
|
2021-05-05 04:22:44 -05:00
|
|
|
DeferredLogger local_deferredLogger;
|
2021-09-20 03:56:11 -05:00
|
|
|
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
|
|
|
{
|
2019-02-07 07:43:17 -06:00
|
|
|
if (localWellsActive()) {
|
|
|
|
for (auto& well : well_container_) {
|
2021-03-26 16:10:16 -05:00
|
|
|
well->recoverWellSolutionAndUpdateWellState(x, this->wellState(), local_deferredLogger);
|
2019-02-07 07:43:17 -06:00
|
|
|
}
|
|
|
|
}
|
2021-05-25 05:57:11 -05:00
|
|
|
|
2021-03-18 02:29:23 -05:00
|
|
|
}
|
2021-09-20 03:56:11 -05:00
|
|
|
OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger,
|
|
|
|
"recoverWellSolutionAndUpdateWellState() failed: ",
|
2021-05-25 05:57:11 -05:00
|
|
|
terminal_output_, ebosSimulator_.vanguard().grid().comm());
|
|
|
|
|
2017-02-13 10:07:34 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2017-05-03 06:34:15 -05:00
|
|
|
template<typename TypeTag>
|
2017-02-13 10:07:34 -06:00
|
|
|
void
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2017-08-21 08:41:25 -05:00
|
|
|
initPrimaryVariablesEvaluation() const
|
2017-02-13 10:07:34 -06:00
|
|
|
{
|
2017-07-25 03:15:27 -05:00
|
|
|
for (auto& well : well_container_) {
|
2017-08-21 08:41:25 -05:00
|
|
|
well->initPrimaryVariablesEvaluation();
|
2017-06-19 05:43:08 -05:00
|
|
|
}
|
2017-02-13 10:07:34 -06:00
|
|
|
}
|
|
|
|
|
2017-02-14 08:06:57 -06:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2017-02-14 04:34:03 -06:00
|
|
|
|
2017-05-03 06:34:15 -05:00
|
|
|
template<typename TypeTag>
|
2018-11-13 07:02:55 -06:00
|
|
|
ConvergenceReport
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2020-04-01 02:55:29 -05:00
|
|
|
getWellConvergence(const std::vector<Scalar>& B_avg, bool checkGroupConvergence) const
|
2017-02-14 04:34:03 -06:00
|
|
|
{
|
2019-02-07 07:43:17 -06:00
|
|
|
|
2021-05-05 04:22:44 -05:00
|
|
|
DeferredLogger local_deferredLogger;
|
2018-11-13 07:02:55 -06:00
|
|
|
// Get global (from all processes) convergence report.
|
|
|
|
ConvergenceReport local_report;
|
2021-09-02 06:47:38 -05:00
|
|
|
const int iterationIdx = ebosSimulator_.model().newtonMethod().numIterations();
|
2017-06-28 04:15:04 -05:00
|
|
|
for (const auto& well : well_container_) {
|
2022-04-19 07:13:48 -05:00
|
|
|
if (well->isOperableAndSolvable() || well->wellIsStopped()) {
|
2021-09-06 03:08:28 -05:00
|
|
|
local_report += well->getWellConvergence(this->wellState(), B_avg, local_deferredLogger, iterationIdx > param_.strict_outer_iter_wells_ );
|
2022-04-05 07:42:27 -05:00
|
|
|
} else {
|
|
|
|
ConvergenceReport report;
|
|
|
|
using CR = ConvergenceReport;
|
|
|
|
report.setWellFailed({CR::WellFailure::Type::Unsolvable, CR::Severity::Normal, -1, well->name()});
|
|
|
|
local_report += report;
|
2018-11-17 16:14:51 -06:00
|
|
|
}
|
2017-08-22 07:49:30 -05:00
|
|
|
}
|
2021-05-25 05:57:11 -05:00
|
|
|
|
|
|
|
const Opm::Parallel::Communication comm = grid().comm();
|
|
|
|
DeferredLogger global_deferredLogger = gatherDeferredLogger(local_deferredLogger, comm);
|
2019-02-07 07:43:17 -06:00
|
|
|
if (terminal_output_) {
|
|
|
|
global_deferredLogger.logMessages();
|
|
|
|
}
|
|
|
|
|
2021-05-25 05:57:11 -05:00
|
|
|
ConvergenceReport report = gatherConvergenceReport(local_report, comm);
|
2017-08-22 07:49:30 -05:00
|
|
|
|
2018-11-13 07:02:55 -06:00
|
|
|
// Log debug messages for NaN or too large residuals.
|
2019-02-03 01:13:11 -06:00
|
|
|
if (terminal_output_) {
|
|
|
|
for (const auto& f : report.wellFailures()) {
|
|
|
|
if (f.severity() == ConvergenceReport::Severity::NotANumber) {
|
|
|
|
OpmLog::debug("NaN residual found with phase " + std::to_string(f.phase()) + " for well " + f.wellName());
|
|
|
|
} else if (f.severity() == ConvergenceReport::Severity::TooLarge) {
|
|
|
|
OpmLog::debug("Too large residual found with phase " + std::to_string(f.phase()) + " for well " + f.wellName());
|
|
|
|
}
|
2017-03-24 06:12:06 -05:00
|
|
|
}
|
|
|
|
}
|
2020-04-08 03:41:20 -05:00
|
|
|
|
2020-04-01 02:55:29 -05:00
|
|
|
if (checkGroupConvergence) {
|
|
|
|
const int reportStepIdx = ebosSimulator_.episodeIndex();
|
|
|
|
const Group& fieldGroup = schedule().getGroup("FIELD", reportStepIdx);
|
2021-05-19 07:51:14 -05:00
|
|
|
bool violated = checkGroupConstraints(fieldGroup,
|
|
|
|
ebosSimulator_.episodeIndex(),
|
|
|
|
global_deferredLogger);
|
2020-04-01 02:55:29 -05:00
|
|
|
report.setGroupConverged(!violated);
|
|
|
|
}
|
2018-11-13 07:02:55 -06:00
|
|
|
return report;
|
2017-02-14 04:34:03 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2017-05-03 06:34:15 -05:00
|
|
|
template<typename TypeTag>
|
2017-02-14 04:34:03 -06:00
|
|
|
void
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2021-05-05 04:22:44 -05:00
|
|
|
calculateExplicitQuantities(DeferredLogger& deferred_logger) const
|
2017-02-14 04:34:03 -06:00
|
|
|
{
|
2021-09-29 09:01:16 -05:00
|
|
|
// TODO: checking isOperableAndSolvable() ?
|
2018-11-17 16:14:51 -06:00
|
|
|
for (auto& well : well_container_) {
|
2021-03-26 16:10:16 -05:00
|
|
|
well->calculateExplicitQuantities(ebosSimulator_, this->wellState(), deferred_logger);
|
2018-11-17 16:14:51 -06:00
|
|
|
}
|
2017-02-14 04:34:03 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2017-05-03 06:34:15 -05:00
|
|
|
template<typename TypeTag>
|
2017-02-14 04:34:03 -06:00
|
|
|
void
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2021-05-05 04:22:44 -05:00
|
|
|
updateWellControls(DeferredLogger& deferred_logger, const bool checkGroupControls)
|
2017-02-14 04:34:03 -06:00
|
|
|
{
|
2019-01-21 01:26:28 -06:00
|
|
|
// Even if there are no wells active locally, we cannot
|
|
|
|
// return as the DeferredLogger uses global communication.
|
|
|
|
// For no well active globally we simply return.
|
2017-04-12 10:37:34 -05:00
|
|
|
if( !wellsActive() ) return ;
|
2017-02-14 04:34:03 -06:00
|
|
|
|
2021-05-19 07:51:14 -05:00
|
|
|
const int episodeIdx = ebosSimulator_.episodeIndex();
|
|
|
|
const int iterationIdx = ebosSimulator_.model().newtonMethod().numIterations();
|
2021-08-27 00:57:16 -05:00
|
|
|
const auto& comm = ebosSimulator_.vanguard().grid().comm();
|
2021-05-19 07:51:14 -05:00
|
|
|
updateAndCommunicateGroupData(episodeIdx, iterationIdx);
|
2020-02-10 08:16:09 -06:00
|
|
|
|
2021-06-07 05:20:49 -05:00
|
|
|
updateNetworkPressures(episodeIdx);
|
2020-05-15 04:21:32 -05:00
|
|
|
|
2020-02-10 08:16:09 -06:00
|
|
|
std::set<std::string> switched_wells;
|
|
|
|
std::set<std::string> switched_groups;
|
|
|
|
|
|
|
|
if (checkGroupControls) {
|
|
|
|
// Check group individual constraints.
|
2021-08-27 00:57:16 -05:00
|
|
|
bool changed_individual = updateGroupIndividualControls(deferred_logger, switched_groups,
|
|
|
|
episodeIdx, iterationIdx);
|
2020-02-10 08:16:09 -06:00
|
|
|
|
2021-08-27 00:57:16 -05:00
|
|
|
if (changed_individual)
|
|
|
|
updateAndCommunicate(episodeIdx, iterationIdx, deferred_logger);
|
2020-02-10 08:16:09 -06:00
|
|
|
|
|
|
|
// Check group's constraints from higher levels.
|
2021-08-27 00:57:16 -05:00
|
|
|
bool changed_higher = updateGroupHigherControls(deferred_logger,
|
|
|
|
episodeIdx,
|
|
|
|
switched_groups);
|
2020-02-10 08:16:09 -06:00
|
|
|
|
2021-08-27 00:57:16 -05:00
|
|
|
if (changed_higher)
|
|
|
|
updateAndCommunicate(episodeIdx, iterationIdx, deferred_logger);
|
2019-11-22 05:29:47 -06:00
|
|
|
|
2020-02-10 08:16:09 -06:00
|
|
|
// Check wells' group constraints and communicate.
|
2021-08-27 00:57:16 -05:00
|
|
|
bool changed_well_group = false;
|
2020-02-10 08:16:09 -06:00
|
|
|
for (const auto& well : well_container_) {
|
|
|
|
const auto mode = WellInterface<TypeTag>::IndividualOrGroup::Group;
|
2021-08-27 00:57:16 -05:00
|
|
|
const bool changed_well = well->updateWellControl(ebosSimulator_, mode, this->wellState(), this->groupState(), deferred_logger);
|
|
|
|
if (changed_well) {
|
2020-02-10 08:16:09 -06:00
|
|
|
switched_wells.insert(well->name());
|
2021-08-27 00:57:16 -05:00
|
|
|
changed_well_group = changed_well || changed_well_group;
|
2020-02-10 08:16:09 -06:00
|
|
|
}
|
|
|
|
}
|
2021-08-27 00:57:16 -05:00
|
|
|
changed_well_group = comm.sum(changed_well_group);
|
|
|
|
if (changed_well_group)
|
|
|
|
updateAndCommunicate(episodeIdx, iterationIdx, deferred_logger);
|
2019-08-07 07:13:11 -05:00
|
|
|
}
|
|
|
|
|
2020-02-10 08:16:09 -06:00
|
|
|
// Check individual well constraints and communicate.
|
2021-08-27 00:57:16 -05:00
|
|
|
bool changed_well_individual = false;
|
2017-06-28 04:15:04 -05:00
|
|
|
for (const auto& well : well_container_) {
|
2020-02-10 08:16:09 -06:00
|
|
|
if (switched_wells.count(well->name())) {
|
|
|
|
continue;
|
|
|
|
}
|
|
|
|
const auto mode = WellInterface<TypeTag>::IndividualOrGroup::Individual;
|
2021-08-27 00:57:16 -05:00
|
|
|
const bool changed_well = well->updateWellControl(ebosSimulator_, mode, this->wellState(), this->groupState(), deferred_logger);
|
|
|
|
if (changed_well) {
|
|
|
|
changed_well_individual = changed_well || changed_well_individual;
|
|
|
|
}
|
2019-01-18 07:04:30 -06:00
|
|
|
}
|
2021-08-27 00:57:16 -05:00
|
|
|
changed_well_individual = comm.sum(changed_well_individual);
|
|
|
|
if (changed_well_individual)
|
|
|
|
updateAndCommunicate(episodeIdx, iterationIdx, deferred_logger);
|
2021-09-24 05:09:02 -05:00
|
|
|
|
|
|
|
// update wsolvent fraction for REIN wells
|
|
|
|
const Group& fieldGroup = schedule().getGroup("FIELD", episodeIdx);
|
|
|
|
updateWsolvent(fieldGroup, episodeIdx, this->nupcolWellState());
|
2017-02-14 04:34:03 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
|
2021-08-25 07:06:56 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
updateAndCommunicate(const int reportStepIdx,
|
|
|
|
const int iterationIdx,
|
|
|
|
DeferredLogger& deferred_logger)
|
|
|
|
{
|
|
|
|
updateAndCommunicateGroupData(reportStepIdx, iterationIdx);
|
|
|
|
// if a well or group change control it affects all wells that are under the same group
|
|
|
|
for (const auto& well : well_container_) {
|
|
|
|
well->updateWellStateWithTarget(ebosSimulator_, this->groupState(), this->wellState(), deferred_logger);
|
|
|
|
}
|
2021-08-27 00:57:16 -05:00
|
|
|
updateAndCommunicateGroupData(reportStepIdx, iterationIdx);
|
2021-08-25 07:06:56 -05:00
|
|
|
}
|
2017-02-14 06:39:53 -06:00
|
|
|
|
2017-05-03 06:34:15 -05:00
|
|
|
template<typename TypeTag>
|
2017-02-14 06:39:53 -06:00
|
|
|
void
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2018-06-28 06:47:10 -05:00
|
|
|
updateWellTestState(const double& simulationTime, WellTestState& wellTestState) const
|
2017-02-14 06:39:53 -06:00
|
|
|
{
|
2021-05-05 04:22:44 -05:00
|
|
|
DeferredLogger local_deferredLogger;
|
2017-07-26 04:01:26 -05:00
|
|
|
for (const auto& well : well_container_) {
|
2021-09-20 04:16:32 -05:00
|
|
|
const auto& wname = well->name();
|
2021-10-07 08:18:53 -05:00
|
|
|
const auto wasClosed = wellTestState.well_is_closed(wname);
|
2021-09-21 04:34:19 -05:00
|
|
|
well->checkWellOperability(ebosSimulator_, this->wellState(), local_deferredLogger);
|
2021-09-20 04:16:32 -05:00
|
|
|
well->updateWellTestState(this->wellState().well(wname), simulationTime, /*writeMessageToOPMLog=*/ true, wellTestState, local_deferredLogger);
|
2021-03-01 08:43:43 -06:00
|
|
|
|
2021-10-07 08:18:53 -05:00
|
|
|
if (!wasClosed && wellTestState.well_is_closed(wname)) {
|
2021-09-20 04:16:32 -05:00
|
|
|
this->closed_this_step_.insert(wname);
|
2021-03-01 08:43:43 -06:00
|
|
|
}
|
2019-02-03 01:13:11 -06:00
|
|
|
}
|
2021-05-25 05:57:11 -05:00
|
|
|
|
|
|
|
const Opm::Parallel::Communication comm = grid().comm();
|
|
|
|
DeferredLogger global_deferredLogger = gatherDeferredLogger(local_deferredLogger, comm);
|
2019-02-03 01:13:11 -06:00
|
|
|
if (terminal_output_) {
|
|
|
|
global_deferredLogger.logMessages();
|
2017-07-26 04:01:26 -05:00
|
|
|
}
|
2017-02-14 06:39:53 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
|
2017-05-03 06:34:15 -05:00
|
|
|
template<typename TypeTag>
|
2017-02-14 06:39:53 -06:00
|
|
|
void
|
2021-06-07 07:09:36 -05:00
|
|
|
BlackoilWellModel<TypeTag>::computePotentials(const std::size_t widx,
|
|
|
|
const WellState& well_state_copy,
|
|
|
|
std::string& exc_msg,
|
|
|
|
ExceptionType::ExcEnum& exc_type,
|
|
|
|
DeferredLogger& deferred_logger)
|
2017-02-14 06:39:53 -06:00
|
|
|
{
|
2017-11-08 06:57:36 -06:00
|
|
|
const int np = numPhases();
|
2021-06-07 07:09:36 -05:00
|
|
|
std::vector<double> potentials;
|
|
|
|
const auto& well= well_container_[widx];
|
|
|
|
try {
|
|
|
|
well->computeWellPotentials(ebosSimulator_, well_state_copy, potentials, deferred_logger);
|
|
|
|
}
|
2021-09-20 03:56:11 -05:00
|
|
|
// catch all possible exception and store type and message.
|
|
|
|
OPM_PARALLEL_CATCH_CLAUSE(exc_type, exc_msg);
|
2021-06-07 07:09:36 -05:00
|
|
|
// Store it in the well state
|
|
|
|
// potentials is resized and set to zero in the beginning of well->ComputeWellPotentials
|
|
|
|
// and updated only if sucessfull. i.e. the potentials are zero for exceptions
|
2021-08-05 03:57:15 -05:00
|
|
|
auto& ws = this->wellState().well(well->indexOfWell());
|
2021-06-07 07:09:36 -05:00
|
|
|
for (int p = 0; p < np; ++p) {
|
2021-11-18 05:57:16 -06:00
|
|
|
// make sure the potentials are positive
|
|
|
|
ws.well_potentials[p] = std::max(0.0, potentials[p]);
|
2019-02-07 07:43:17 -06:00
|
|
|
}
|
2017-02-14 08:06:57 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
2021-03-01 18:11:19 -06:00
|
|
|
template <typename TypeTag>
|
2020-10-09 06:38:33 -05:00
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
calculateProductivityIndexValues(DeferredLogger& deferred_logger)
|
|
|
|
{
|
2021-03-01 18:11:19 -06:00
|
|
|
for (const auto& wellPtr : this->well_container_) {
|
|
|
|
this->calculateProductivityIndexValues(wellPtr.get(), deferred_logger);
|
2020-10-09 06:38:33 -05:00
|
|
|
}
|
2021-03-01 18:11:19 -06:00
|
|
|
}
|
2020-10-09 06:38:33 -05:00
|
|
|
|
2021-03-01 18:11:19 -06:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
calculateProductivityIndexValuesShutWells(const int reportStepIdx,
|
|
|
|
DeferredLogger& deferred_logger)
|
|
|
|
{
|
|
|
|
// For the purpose of computing PI/II values, it is sufficient to
|
|
|
|
// construct StandardWell instances only. We don't need to form
|
|
|
|
// well objects that honour the 'isMultisegment()' flag of the
|
|
|
|
// corresponding "this->wells_ecl_[shutWell]".
|
|
|
|
|
|
|
|
for (const auto& shutWell : this->local_shut_wells_) {
|
2021-06-02 04:02:46 -05:00
|
|
|
if (this->wells_ecl_[shutWell].getConnections().empty()) {
|
|
|
|
// No connections in this well. Nothing to do.
|
|
|
|
continue;
|
|
|
|
}
|
|
|
|
|
2021-03-01 18:11:19 -06:00
|
|
|
auto wellPtr = this->template createTypedWellPointer
|
|
|
|
<StandardWell<TypeTag>>(shutWell, reportStepIdx);
|
|
|
|
|
|
|
|
wellPtr->init(&this->phase_usage_, this->depth_, this->gravity_,
|
2022-04-12 01:44:52 -05:00
|
|
|
this->local_num_cells_, this->B_avg_, true);
|
2021-03-01 18:11:19 -06:00
|
|
|
|
|
|
|
this->calculateProductivityIndexValues(wellPtr.get(), deferred_logger);
|
2020-10-09 06:38:33 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2021-03-01 18:11:19 -06:00
|
|
|
template <typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
calculateProductivityIndexValues(const WellInterface<TypeTag>* wellPtr,
|
|
|
|
DeferredLogger& deferred_logger)
|
|
|
|
{
|
|
|
|
wellPtr->updateProductivityIndex(this->ebosSimulator_,
|
|
|
|
this->prod_index_calc_[wellPtr->indexOfWell()],
|
|
|
|
this->wellState(),
|
|
|
|
deferred_logger);
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2017-05-03 06:34:15 -05:00
|
|
|
template<typename TypeTag>
|
2017-03-16 10:39:05 -05:00
|
|
|
void
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2021-05-05 04:22:44 -05:00
|
|
|
prepareTimeStep(DeferredLogger& deferred_logger)
|
2017-03-16 10:39:05 -05:00
|
|
|
{
|
2021-09-23 08:22:48 -05:00
|
|
|
for (const auto& well : well_container_) {
|
2021-09-24 02:46:22 -05:00
|
|
|
auto& events = this->wellState().well(well->indexOfWell()).events;
|
|
|
|
if (events.hasEvent(WellState::event_mask)) {
|
|
|
|
well->updateWellStateWithTarget(ebosSimulator_, this->groupState(), this->wellState(), deferred_logger);
|
2022-04-27 06:36:52 -05:00
|
|
|
well->updatePrimaryVariables(this->wellState(), deferred_logger);
|
|
|
|
well->initPrimaryVariablesEvaluation();
|
2021-09-24 02:46:22 -05:00
|
|
|
// There is no new well control change input within a report step,
|
|
|
|
// so next time step, the well does not consider to have effective events anymore.
|
|
|
|
events.clearEvent(WellState::event_mask);
|
2021-09-23 08:22:48 -05:00
|
|
|
}
|
2021-09-24 02:46:22 -05:00
|
|
|
// solve the well equation initially to improve the initial solution of the well model
|
2021-11-18 05:57:16 -06:00
|
|
|
if (param_.solve_welleq_initially_ && well->isOperableAndSolvable()) {
|
|
|
|
try {
|
|
|
|
well->solveWellEquation(ebosSimulator_, this->wellState(), this->groupState(), deferred_logger);
|
|
|
|
} catch (const std::exception& e) {
|
|
|
|
const std::string msg = "Compute initial well solution for " + well->name() + " initially failed. Continue with the privious rates";
|
|
|
|
deferred_logger.warning("WELL_INITIAL_SOLVE_FAILED", msg);
|
|
|
|
}
|
2021-09-24 02:46:22 -05:00
|
|
|
}
|
|
|
|
}
|
2021-09-23 08:22:48 -05:00
|
|
|
updatePrimaryVariables(deferred_logger);
|
2017-09-04 06:54:41 -05:00
|
|
|
}
|
|
|
|
|
2017-06-07 02:29:31 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2021-03-18 08:49:52 -05:00
|
|
|
updateAverageFormationFactor()
|
2017-06-23 05:24:50 -05:00
|
|
|
{
|
2021-03-18 08:49:52 -05:00
|
|
|
std::vector< Scalar > B_avg(numComponents(), Scalar() );
|
2018-02-01 09:27:42 -06:00
|
|
|
const auto& grid = ebosSimulator_.vanguard().grid();
|
2017-06-23 05:24:50 -05:00
|
|
|
const auto& gridView = grid.leafGridView();
|
2017-11-08 06:57:36 -06:00
|
|
|
ElementContext elemCtx(ebosSimulator_);
|
2017-06-23 05:24:50 -05:00
|
|
|
const auto& elemEndIt = gridView.template end</*codim=*/0, Dune::Interior_Partition>();
|
2021-09-20 04:12:27 -05:00
|
|
|
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
2017-06-23 05:24:50 -05:00
|
|
|
|
|
|
|
for (auto elemIt = gridView.template begin</*codim=*/0, Dune::Interior_Partition>();
|
|
|
|
elemIt != elemEndIt; ++elemIt)
|
|
|
|
{
|
|
|
|
elemCtx.updatePrimaryStencil(*elemIt);
|
|
|
|
elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
|
|
|
|
|
|
|
|
const auto& intQuants = elemCtx.intensiveQuantities(/*spaceIdx=*/0, /*timeIdx=*/0);
|
|
|
|
const auto& fs = intQuants.fluidState();
|
|
|
|
|
2017-12-04 02:14:08 -06:00
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
|
2017-06-23 05:24:50 -05:00
|
|
|
{
|
2017-12-04 02:14:08 -06:00
|
|
|
if (!FluidSystem::phaseIsActive(phaseIdx)) {
|
|
|
|
continue;
|
|
|
|
}
|
|
|
|
|
|
|
|
const unsigned compIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
|
|
|
|
auto& B = B_avg[ compIdx ];
|
2017-06-23 05:24:50 -05:00
|
|
|
|
2017-12-04 02:14:08 -06:00
|
|
|
B += 1 / fs.invB(phaseIdx).value();
|
2017-06-23 05:24:50 -05:00
|
|
|
}
|
2021-01-28 07:33:21 -06:00
|
|
|
if constexpr (has_solvent_) {
|
2017-06-27 08:16:22 -05:00
|
|
|
auto& B = B_avg[solventSaturationIdx];
|
2017-06-23 05:24:50 -05:00
|
|
|
B += 1 / intQuants.solventInverseFormationVolumeFactor().value();
|
|
|
|
}
|
|
|
|
}
|
2021-05-25 05:57:11 -05:00
|
|
|
OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel::updateAverageFormationFactor() failed: ", grid.comm())
|
2017-06-23 05:24:50 -05:00
|
|
|
|
|
|
|
// compute global average
|
|
|
|
grid.comm().sum(B_avg.data(), B_avg.size());
|
|
|
|
for(auto& bval: B_avg)
|
|
|
|
{
|
2020-10-05 13:02:13 -05:00
|
|
|
bval/=global_num_cells_;
|
2017-06-23 05:24:50 -05:00
|
|
|
}
|
2021-03-18 08:49:52 -05:00
|
|
|
B_avg_ = B_avg;
|
2017-06-23 05:24:50 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
|
2017-08-08 03:44:10 -05:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
2017-09-26 03:52:05 -05:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2021-05-05 04:22:44 -05:00
|
|
|
updatePrimaryVariables(DeferredLogger& deferred_logger)
|
2017-08-08 03:44:10 -05:00
|
|
|
{
|
|
|
|
for (const auto& well : well_container_) {
|
2021-03-26 16:10:16 -05:00
|
|
|
well->updatePrimaryVariables(this->wellState(), deferred_logger);
|
2017-11-08 06:57:36 -06:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2017-11-08 08:48:30 -06:00
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::extractLegacyCellPvtRegionIndex_()
|
|
|
|
{
|
2018-02-01 09:27:42 -06:00
|
|
|
const auto& grid = ebosSimulator_.vanguard().grid();
|
2017-11-08 08:48:30 -06:00
|
|
|
const auto& eclProblem = ebosSimulator_.problem();
|
|
|
|
const unsigned numCells = grid.size(/*codim=*/0);
|
|
|
|
|
|
|
|
pvt_region_idx_.resize(numCells);
|
|
|
|
for (unsigned cellIdx = 0; cellIdx < numCells; ++cellIdx) {
|
|
|
|
pvt_region_idx_[cellIdx] =
|
|
|
|
eclProblem.pvtRegionIndex(cellIdx);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
// The number of components in the model.
|
|
|
|
template<typename TypeTag>
|
|
|
|
int
|
|
|
|
BlackoilWellModel<TypeTag>::numComponents() const
|
|
|
|
{
|
2019-11-01 08:58:10 -05:00
|
|
|
if (wellsActive() && numPhases() < 3) {
|
2019-10-09 08:24:23 -05:00
|
|
|
return numPhases();
|
2017-11-08 08:48:30 -06:00
|
|
|
}
|
|
|
|
int numComp = FluidSystem::numComponents;
|
2021-01-28 07:33:21 -06:00
|
|
|
if constexpr (has_solvent_) {
|
2017-11-08 08:48:30 -06:00
|
|
|
numComp ++;
|
|
|
|
}
|
|
|
|
|
|
|
|
return numComp;
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::extractLegacyDepth_()
|
|
|
|
{
|
2021-12-01 07:00:21 -06:00
|
|
|
const auto& eclProblem = ebosSimulator_.problem();
|
|
|
|
depth_.resize(local_num_cells_);
|
|
|
|
for (unsigned cellIdx = 0; cellIdx < local_num_cells_; ++cellIdx) {
|
|
|
|
depth_[cellIdx] = eclProblem.dofCenterDepth(cellIdx);
|
2017-11-08 08:48:30 -06:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
updatePerforationIntensiveQuantities() {
|
|
|
|
ElementContext elemCtx(ebosSimulator_);
|
|
|
|
const auto& gridView = ebosSimulator_.gridView();
|
|
|
|
const auto& elemEndIt = gridView.template end</*codim=*/0, Dune::Interior_Partition>();
|
2021-09-20 04:12:27 -05:00
|
|
|
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
2017-11-08 08:48:30 -06:00
|
|
|
for (auto elemIt = gridView.template begin</*codim=*/0, Dune::Interior_Partition>();
|
|
|
|
elemIt != elemEndIt;
|
|
|
|
++elemIt)
|
|
|
|
{
|
2018-11-16 09:02:47 -06:00
|
|
|
|
2017-11-08 08:48:30 -06:00
|
|
|
elemCtx.updatePrimaryStencil(*elemIt);
|
2018-11-16 09:02:47 -06:00
|
|
|
int elemIdx = elemCtx.globalSpaceIndex(0, 0);
|
|
|
|
|
|
|
|
if (!is_cell_perforated_[elemIdx]) {
|
|
|
|
continue;
|
|
|
|
}
|
2017-11-08 08:48:30 -06:00
|
|
|
elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
|
|
|
|
}
|
2021-05-25 05:57:11 -05:00
|
|
|
OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel::updatePerforationIntensiveQuantities() failed: ", ebosSimulator_.vanguard().grid().comm());
|
2017-11-08 08:48:30 -06:00
|
|
|
}
|
|
|
|
|
2017-11-08 06:57:36 -06:00
|
|
|
|
2019-09-23 08:15:55 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
typename BlackoilWellModel<TypeTag>::WellInterfacePtr
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
getWell(const std::string& well_name) const
|
|
|
|
{
|
|
|
|
// finding the iterator of the well in wells_ecl
|
|
|
|
auto well = std::find_if(well_container_.begin(),
|
|
|
|
well_container_.end(),
|
|
|
|
[&well_name](const WellInterfacePtr& elem)->bool {
|
|
|
|
return elem->name() == well_name;
|
|
|
|
});
|
|
|
|
|
|
|
|
assert(well != well_container_.end());
|
|
|
|
|
|
|
|
return *well;
|
|
|
|
}
|
|
|
|
|
2021-08-19 05:15:35 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
bool
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
hasWell(const std::string& well_name) const
|
|
|
|
{
|
|
|
|
return std::any_of(well_container_.begin(), well_container_.end(),
|
|
|
|
[&well_name](const WellInterfacePtr& elem) -> bool
|
|
|
|
{
|
|
|
|
return elem->name() == well_name;
|
|
|
|
});
|
|
|
|
}
|
|
|
|
|
2019-09-23 08:15:55 -05:00
|
|
|
|
2020-02-10 08:16:09 -06:00
|
|
|
|
|
|
|
|
2021-05-19 07:51:14 -05:00
|
|
|
template <typename TypeTag>
|
|
|
|
int
|
2020-02-10 08:16:09 -06:00
|
|
|
BlackoilWellModel<TypeTag>::
|
2021-05-19 07:51:14 -05:00
|
|
|
reportStepIndex() const
|
2020-02-10 08:16:09 -06:00
|
|
|
{
|
2021-05-19 07:51:14 -05:00
|
|
|
return std::max(this->ebosSimulator_.episodeIndex(), 0);
|
2020-03-24 03:24:45 -05:00
|
|
|
}
|
2020-04-08 03:41:20 -05:00
|
|
|
|
|
|
|
|
2019-08-07 07:13:11 -05:00
|
|
|
|
2020-03-24 03:24:45 -05:00
|
|
|
|
2019-08-07 07:13:11 -05:00
|
|
|
|
2020-04-24 08:25:38 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
2021-05-19 07:51:14 -05:00
|
|
|
calcRates(const int fipnum,
|
|
|
|
const int pvtreg,
|
|
|
|
std::vector<double>& resv_coeff)
|
2020-04-24 08:25:38 -05:00
|
|
|
{
|
2021-05-19 07:51:14 -05:00
|
|
|
rateConverter_->calcCoeff(fipnum, pvtreg, resv_coeff);
|
2020-03-24 03:24:45 -05:00
|
|
|
}
|
2020-05-04 08:56:34 -05:00
|
|
|
|
2021-06-08 08:29:15 -05:00
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
calcInjRates(const int fipnum,
|
|
|
|
const int pvtreg,
|
|
|
|
std::vector<double>& resv_coeff)
|
|
|
|
{
|
|
|
|
rateConverter_->calcInjCoeff(fipnum, pvtreg, resv_coeff);
|
|
|
|
}
|
|
|
|
|
2021-03-08 08:11:50 -06:00
|
|
|
|
|
|
|
template <typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
|
|
|
computeWellTemperature()
|
|
|
|
{
|
|
|
|
if (!has_energy_)
|
|
|
|
return;
|
|
|
|
|
|
|
|
int np = numPhases();
|
|
|
|
double cellInternalEnergy;
|
|
|
|
double cellBinv;
|
|
|
|
double cellDensity;
|
|
|
|
double perfPhaseRate;
|
|
|
|
const int nw = numLocalWells();
|
|
|
|
for (auto wellID = 0*nw; wellID < nw; ++wellID) {
|
|
|
|
const Well& well = wells_ecl_[wellID];
|
|
|
|
if (well.isInjector())
|
|
|
|
continue;
|
|
|
|
|
|
|
|
int connpos = 0;
|
|
|
|
for (int i = 0; i < wellID; ++i) {
|
|
|
|
connpos += well_perf_data_[i].size();
|
|
|
|
}
|
|
|
|
connpos *= np;
|
2022-08-29 06:41:32 -05:00
|
|
|
std::array<double,2> weighted{0.0,0.0};
|
|
|
|
auto& [weighted_temperature, total_weight] = weighted;
|
2021-03-08 08:11:50 -06:00
|
|
|
|
2021-09-27 05:00:39 -05:00
|
|
|
auto& well_info = local_parallel_well_info_[wellID].get();
|
2021-03-08 08:11:50 -06:00
|
|
|
const int num_perf_this_well = well_info.communication().sum(well_perf_data_[wellID].size());
|
2021-08-06 06:14:00 -05:00
|
|
|
auto& ws = this->wellState().well(wellID);
|
|
|
|
auto& perf_data = ws.perf_data;
|
2021-06-06 10:09:44 -05:00
|
|
|
auto& perf_phase_rate = perf_data.phase_rates;
|
2021-03-08 08:11:50 -06:00
|
|
|
|
|
|
|
for (int perf = 0; perf < num_perf_this_well; ++perf) {
|
|
|
|
const int cell_idx = well_perf_data_[wellID][perf].cell_index;
|
|
|
|
const auto& intQuants = *(ebosSimulator_.model().cachedIntensiveQuantities(cell_idx, /*timeIdx=*/0));
|
|
|
|
const auto& fs = intQuants.fluidState();
|
|
|
|
|
2021-12-06 03:26:54 -06:00
|
|
|
// we on only have one temperature pr cell any phaseIdx will do
|
2021-03-08 08:11:50 -06:00
|
|
|
double cellTemperatures = fs.temperature(/*phaseIdx*/0).value();
|
2021-12-06 03:26:54 -06:00
|
|
|
|
2021-03-08 08:11:50 -06:00
|
|
|
double weight_factor = 0.0;
|
2021-12-06 03:26:54 -06:00
|
|
|
for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
|
|
|
|
{
|
|
|
|
if (!FluidSystem::phaseIsActive(phaseIdx)) {
|
|
|
|
continue;
|
|
|
|
}
|
2021-03-08 08:11:50 -06:00
|
|
|
cellInternalEnergy = fs.enthalpy(phaseIdx).value() - fs.pressure(phaseIdx).value() / fs.density(phaseIdx).value();
|
|
|
|
cellBinv = fs.invB(phaseIdx).value();
|
|
|
|
cellDensity = fs.density(phaseIdx).value();
|
2021-05-20 09:16:12 -05:00
|
|
|
perfPhaseRate = perf_phase_rate[ perf*np + phaseIdx ];
|
2021-03-08 08:11:50 -06:00
|
|
|
weight_factor += cellDensity * perfPhaseRate/cellBinv * cellInternalEnergy/cellTemperatures;
|
|
|
|
}
|
|
|
|
total_weight += weight_factor;
|
|
|
|
weighted_temperature += weight_factor * cellTemperatures;
|
|
|
|
}
|
2022-08-29 06:41:32 -05:00
|
|
|
well_info.communication().sum(weighted.data(), 2);
|
2021-08-03 14:25:03 -05:00
|
|
|
this->wellState().well(wellID).temperature = weighted_temperature/total_weight;
|
2021-03-08 08:11:50 -06:00
|
|
|
}
|
|
|
|
}
|
2021-06-03 05:58:39 -05:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilWellModel<TypeTag>::
|
2021-06-11 07:58:38 -05:00
|
|
|
assignWellTracerRates(data::Wells& wsrpt) const
|
2021-06-03 05:58:39 -05:00
|
|
|
{
|
2021-06-11 07:58:38 -05:00
|
|
|
const auto & wellTracerRates = ebosSimulator_.problem().tracerModel().getWellTracerRates();
|
2021-06-03 05:58:39 -05:00
|
|
|
|
|
|
|
if (wellTracerRates.empty())
|
|
|
|
return; // no tracers
|
|
|
|
|
|
|
|
for (const auto& wTR : wellTracerRates) {
|
|
|
|
std::string wellName = wTR.first.first;
|
|
|
|
auto xwPos = wsrpt.find(wellName);
|
|
|
|
if (xwPos == wsrpt.end()) { // No well results.
|
|
|
|
continue;
|
|
|
|
}
|
|
|
|
std::string tracerName = wTR.first.second;
|
|
|
|
double rate = wTR.second;
|
|
|
|
xwPos->second.rates.set(data::Rates::opt::tracer, rate, tracerName);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2017-02-13 09:45:06 -06:00
|
|
|
} // namespace Opm
|