From b7fd9e42f43cf7aca3e83aae778eba284d8d5268 Mon Sep 17 00:00:00 2001 From: =?UTF-8?q?B=C3=A5rd=20Skaflestad?= Date: Tue, 1 Dec 2020 18:04:46 +0100 Subject: [PATCH] Well Model: Add Helper for Well Pointer Creation This commit adds a new helper function, WellInterfacePtr createWellPointer(wellID, reportStep) const which is responsible for creating appropriately typed derived well pointers depending on well types (multi-segment vs. standard). This, in turn, allows us to centralise this logic and use the same factory function both when creating the 'well_container_' and when forming the well-test objects. Finally, this helper will become useful for calculating PI/II values of shut/stopped wells in the context of WELPI. --- opm/simulators/wells/BlackoilWellModel.hpp | 18 +++- .../wells/BlackoilWellModel_impl.hpp | 101 ++++++++---------- 2 files changed, 63 insertions(+), 56 deletions(-) diff --git a/opm/simulators/wells/BlackoilWellModel.hpp b/opm/simulators/wells/BlackoilWellModel.hpp index 7da8f2dec..43b079a0c 100644 --- a/opm/simulators/wells/BlackoilWellModel.hpp +++ b/opm/simulators/wells/BlackoilWellModel.hpp @@ -31,8 +31,15 @@ #include #include -#include +#include +#include +#include +#include #include +#include +#include + +#include #include @@ -289,6 +296,15 @@ namespace Opm { // create the well container std::vector createWellContainer(const int time_step); + WellInterfacePtr + createWellPointer(const int wellID, + const int time_step) const; + + template + std::unique_ptr + createTypedWellPointer(const int wellID, + const int time_step) const; + WellInterfacePtr createWellForWellTest(const std::string& well_name, const int report_step, Opm::DeferredLogger& deferred_logger) const; WellState well_state_; diff --git a/opm/simulators/wells/BlackoilWellModel_impl.hpp b/opm/simulators/wells/BlackoilWellModel_impl.hpp index 52e3cf1be..32cda46d1 100644 --- a/opm/simulators/wells/BlackoilWellModel_impl.hpp +++ b/opm/simulators/wells/BlackoilWellModel_impl.hpp @@ -752,33 +752,8 @@ namespace Opm { wellIsStopped = true; } - // Use the pvtRegionIdx from the top cell - const int well_cell_top = well_perf_data_[w][0].cell_index; - const int pvtreg = pvt_region_idx_[well_cell_top]; + well_container.emplace_back(this->createWellPointer(w, time_step)); - if (!well_ecl.isMultiSegment() || !param_.use_multisegment_well_) { - well_container.emplace_back(new StandardWell(well_ecl, - time_step, - param_, - *rateConverter_, - pvtreg, - numComponents(), - numPhases(), - w, - well_state_.firstPerfIndex()[w], - well_perf_data_[w])); - } else { - well_container.emplace_back(new MultisegmentWell(well_ecl, - time_step, - param_, - *rateConverter_, - pvtreg, - numComponents(), - numPhases(), - w, - well_state_.firstPerfIndex()[w], - well_perf_data_[w])); - } if (wellIsStopped) well_container.back()->stopWell(); } @@ -797,6 +772,50 @@ namespace Opm { + template + typename BlackoilWellModel::WellInterfacePtr + BlackoilWellModel:: + 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>(wellID, time_step); + } + else { + return this->template createTypedWellPointer>(wellID, time_step); + } + } + + + + + + template + template + std::unique_ptr + BlackoilWellModel:: + 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]; + + return std::make_unique(this->wells_ecl_[wellID], + time_step, + this->param_, + *this->rateConverter_, + this->pvt_region_idx_[perf_data.front().cell_index], + this->numComponents(), + this->numPhases(), + wellID, + this->well_state_.firstPerfIndex()[wellID], + perf_data); + } + + + + + template typename BlackoilWellModel::WellInterfacePtr BlackoilWellModel:: @@ -817,35 +836,7 @@ namespace Opm { OPM_DEFLOG_THROW(std::logic_error, "Could not find well " << well_name << " in wells_ecl ", deferred_logger); } - const Well& well_ecl = wells_ecl_[index_well_ecl]; - - // Use the pvtRegionIdx from the top cell - const int well_cell_top = well_perf_data_[index_well_ecl][0].cell_index; - const int pvtreg = pvt_region_idx_[well_cell_top]; - - if (!well_ecl.isMultiSegment() || !param_.use_multisegment_well_) { - return WellInterfacePtr(new StandardWell(well_ecl, - report_step, - param_, - *rateConverter_, - pvtreg, - numComponents(), - numPhases(), - index_well_ecl, - well_state_.firstPerfIndex()[index_well_ecl], - well_perf_data_[index_well_ecl])); - } else { - return WellInterfacePtr(new MultisegmentWell(well_ecl, - report_step, - param_, - *rateConverter_, - pvtreg, - numComponents(), - numPhases(), - index_well_ecl, - well_state_.firstPerfIndex()[index_well_ecl], - well_perf_data_[index_well_ecl])); - } + return this->createWellPointer(index_well_ecl, report_step); }