23#ifndef OPM_BLACKOILWELLMODEL_IMPL_HEADER_INCLUDED
24#define OPM_BLACKOILWELLMODEL_IMPL_HEADER_INCLUDED
27#ifndef OPM_BLACKOILWELLMODEL_HEADER_INCLUDED
32#include <opm/grid/utility/cartesianToCompressed.hpp>
34#include <opm/input/eclipse/Schedule/Network/Balance.hpp>
35#include <opm/input/eclipse/Schedule/Network/ExtNetwork.hpp>
36#include <opm/input/eclipse/Schedule/Well/PAvgDynamicSourceData.hpp>
37#include <opm/input/eclipse/Schedule/Well/WellMatcher.hpp>
38#include <opm/input/eclipse/Schedule/Well/WellTestConfig.hpp>
39#include <opm/input/eclipse/Schedule/Well/WellEconProductionLimits.hpp>
41#include <opm/input/eclipse/Units/UnitSystem.hpp>
67#include <fmt/format.h>
70 template<
typename TypeTag>
77 simulator.vanguard().summaryState(),
78 simulator.vanguard().eclState(),
80 simulator.gridView().comm(),
82 , simulator_(simulator)
83 , guide_rate_handler_{
85 simulator.vanguard().schedule(),
86 simulator.vanguard().summaryState(),
87 simulator.vanguard().grid().comm()
89 , gaslift_(this->terminal_output_)
92 , rescoupHelper_(*this)
101 auto& parallel_wells =
simulator.vanguard().parallelWells();
104 for(
const auto& name_bool : parallel_wells) {
110 Parameters::Get<Parameters::AlternativeWellRateInit>();
112 using SourceDataSpan =
113 typename PAvgDynamicSourceData<Scalar>::template SourceDataSpan<Scalar>;
115 this->
wbp_.initializeSources(
116 [
this](
const std::size_t globalIndex)
118 [
this](
const int localCell, SourceDataSpan sourceTerms)
120 using Item =
typename SourceDataSpan::Item;
122 const auto* intQuants = this->
simulator_.model()
123 .cachedIntensiveQuantities(localCell, 0);
124 const auto& fs = intQuants->fluidState();
127 .set(Item::PoreVol, intQuants->porosity().value() *
128 this->
simulator_.model().dofTotalVolume(localCell))
129 .set(Item::Depth, this->
depth_[localCell]);
131 constexpr auto io = FluidSystem::oilPhaseIdx;
132 constexpr auto ig = FluidSystem::gasPhaseIdx;
133 constexpr auto iw = FluidSystem::waterPhaseIdx;
136 const auto haveOil = FluidSystem::phaseIsActive(io);
137 const auto haveGas = FluidSystem::phaseIsActive(ig);
138 const auto haveWat = FluidSystem::phaseIsActive(iw);
140 auto weightedPhaseDensity = [&fs](
const auto ip)
142 return fs.saturation(ip).value() * fs.density(ip).value();
145 if (haveOil) { sourceTerms.set(Item::Pressure, fs.pressure(io).value()); }
146 else if (haveGas) { sourceTerms.set(Item::Pressure, fs.pressure(ig).value()); }
147 else { sourceTerms.set(Item::Pressure, fs.pressure(iw).value()); }
151 if (haveOil) { rho += weightedPhaseDensity(io); }
152 if (haveGas) { rho += weightedPhaseDensity(ig); }
153 if (haveWat) { rho += weightedPhaseDensity(iw); }
155 sourceTerms.set(Item::MixtureDensity, rho);
160 template<
typename TypeTag>
165 extractLegacyCellPvtRegionIndex_();
166 extractLegacyDepth_();
168 gravity_ = simulator_.problem().gravity()[2];
170 this->initial_step_ =
true;
173 simulator_.model().addAuxiliaryModule(
this);
175 is_cell_perforated_.resize(local_num_cells_,
false);
179 template<
typename TypeTag>
184 const uint64_t effective_events_mask = ScheduleEvents::WELL_STATUS_CHANGE
185 + ScheduleEvents::NEW_WELL;
186 const auto& events = this->schedule()[reportStepIdx].wellgroup_events();
187 for (
auto& wellPtr : this->well_container_) {
188 const bool well_opened_this_step = this->report_step_starts_ &&
189 events.hasEvent(wellPtr->name(),
190 effective_events_mask);
191 wellPtr->init(this->depth_, this->gravity_,
192 this->B_avg_, well_opened_this_step);
196 template<
typename TypeTag>
201 this->groupStateHelper().setReportStep(timeStepIdx);
202 this->report_step_starts_ =
true;
203 this->report_step_start_events_ = this->schedule()[timeStepIdx].wellgroup_events();
205 this->rateConverter_ = std::make_unique<RateConverterType>
206 (std::vector<int>(this->local_num_cells_, 0));
210 const auto enableWellPIScaling =
true;
211 this->initializeLocalWellStructure(timeStepIdx, enableWellPIScaling);
214 this->initializeGroupStructure(timeStepIdx);
216 const auto& comm = this->simulator_.vanguard().grid().comm();
222 this->rateConverter_->template defineState<ElementContext>(this->simulator_);
226 const auto& sched_state = this->schedule()[timeStepIdx];
228 this->vfp_properties_ = std::make_unique<VFPProperties<Scalar, IndexTraits>>
229 (sched_state.vfpinj(), sched_state.vfpprod(), this->wellState());
236 this->commitWGState();
238 this->wellStructureChangedDynamically_ =
false;
245 template <
typename TypeTag>
249 const bool enableWellPIScaling)
251 auto logger_guard = this->groupStateHelper().pushLogger();
252 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
254 const auto& comm = this->simulator_.vanguard().grid().comm();
257 this->wells_ecl_ = this->getLocalWells(reportStepIdx);
258 this->local_parallel_well_info_ =
259 this->createLocalParallelWellInfo(this->wells_ecl_);
266 this->initializeWellPerfData();
267 this->initializeWellState(reportStepIdx);
268 this->wbp_.initializeWBPCalculationService();
270 if (this->param_.use_multisegment_well_ && this->anyMSWellOpenLocal()) {
271 this->wellState().initWellStateMSWell(this->wells_ecl_, &this->prevWellState());
274 this->initializeWellProdIndCalculators();
276 if (enableWellPIScaling && this->schedule()[reportStepIdx].events()
277 .hasEvent(ScheduleEvents::Events::WELL_PRODUCTIVITY_INDEX))
279 this->runWellPIScaling(reportStepIdx, local_deferredLogger);
283 "Failed to initialize local well structure: ",
284 this->terminal_output_, comm)
291 template <
typename TypeTag>
296 const auto& comm = this->simulator_.vanguard().grid().comm();
300 const auto& fieldGroup =
301 this->schedule().getGroup(
"FIELD", reportStepIdx);
303 this->groupStateHelper().setCmodeGroup(fieldGroup);
307 if (this->schedule()[reportStepIdx].has_gpmaint()) {
308 this->groupStateHelper().setRegionAveragePressureCalculator(
310 this->eclState_.fieldProps(),
311 this->regionalAveragePressureCalculator_
323 template<
typename TypeTag>
328 OPM_TIMEBLOCK(beginTimeStep);
330 this->updateAverageFormationFactor();
332 auto logger_guard = this->groupStateHelper().pushLogger();
333 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
335#ifdef RESERVOIR_COUPLING_ENABLED
336 auto rescoup_logger_guard = this->rescoupHelper_.setupScopedLogger(local_deferredLogger);
339 this->switched_prod_groups_.clear();
340 this->switched_inj_groups_.clear();
342 if (this->wellStructureChangedDynamically_) {
347 const auto reportStepIdx =
348 this->simulator_.episodeIndex();
352 const auto enableWellPIScaling =
false;
354 this->initializeLocalWellStructure(reportStepIdx, enableWellPIScaling);
355 this->initializeGroupStructure(reportStepIdx);
357 this->commitWGState();
363 this->wellStructureChangedDynamically_ =
false;
366 this->resetWGState();
367 const int reportStepIdx = simulator_.episodeIndex();
369 this->wellState().updateWellsDefaultALQ(this->schedule(), reportStepIdx, this->summaryState());
370 this->wellState().gliftTimeStepInit();
372 const double simulationTime = simulator_.time();
376 wellTesting(reportStepIdx, simulationTime, local_deferredLogger);
379 createWellContainer(reportStepIdx);
381#ifdef RESERVOIR_COUPLING_ENABLED
382 if (this->isReservoirCouplingMaster()) {
383 if (this->reservoirCouplingMaster().isFirstSubstepOfSyncTimestep()) {
384 this->rescoupHelper_.receiveSlaveGroupData();
391 this->updateAndCommunicateGroupData(reportStepIdx,
false);
394 const Grid& grid = simulator_.vanguard().grid();
395 this->wells_active_ = grid.comm().max(!this->well_container_.empty());
400 this->initWellContainer(reportStepIdx);
403 std::fill(is_cell_perforated_.begin(), is_cell_perforated_.end(),
false);
404 for (
auto& well : well_container_) {
405 well->updatePerforatedCell(is_cell_perforated_);
409 this->calculateEfficiencyFactors(reportStepIdx);
411 if constexpr (has_polymer_)
413 if (PolymerModule::hasPlyshlog() || getPropValue<TypeTag, Properties::EnablePolymerMW>() ) {
414 this->setRepRadiusPerfLength();
421 this->terminal_output_, simulator_.vanguard().grid().comm());
423 for (
auto& well : well_container_) {
424 well->setVFPProperties(this->vfp_properties_.get());
425 well->setGuideRate(&this->guideRate_);
428 this->updateFiltrationModelsPreStep(local_deferredLogger);
431 for (
auto& well : well_container_) {
432 well->closeCompletions(this->wellTestState());
438 if (alternative_well_rate_init_) {
443 for (
const auto& well : well_container_) {
444 if (well->isProducer() && !well->wellIsStopped()) {
445 well->initializeProducerWellState(simulator_, this->wellState(), local_deferredLogger);
450 for (
const auto& well : well_container_) {
451 if (well->isVFPActive(local_deferredLogger)){
452 well->setPrevSurfaceRates(this->wellState(), this->prevWellState());
456 this->updateWellPotentials(reportStepIdx,
458 simulator_.vanguard().summaryConfig(),
459 local_deferredLogger);
460 }
catch ( std::runtime_error& e ) {
461 const std::string msg =
"A zero well potential is returned for output purposes. ";
462 local_deferredLogger.warning(
"WELL_POTENTIAL_CALCULATION_FAILED", msg);
465 this->guide_rate_handler_.updateGuideRates(
466 reportStepIdx, simulationTime, this->wellState(), this->groupState()
468 bool slave_needs_well_solution =
false;
469#ifdef RESERVOIR_COUPLING_ENABLED
470 if (this->isReservoirCouplingSlave()) {
471 if (this->reservoirCouplingSlave().isFirstSubstepOfSyncTimestep()) {
472 this->rescoupHelper_.sendSlaveGroupDataToMaster();
473 this->rescoupHelper_.receiveGroupConstraintsFromMaster();
474 this->rescoupHelper_.receiveCoupledNetworkActiveStatus();
475 this->groupStateHelper().updateSlaveGroupCmodesFromMaster();
476 this->reservoirCouplingSlave().markSlaveGroupsInSchedule(
477 this->schedule_, reportStepIdx);
478 slave_needs_well_solution =
true;
485 if (this->schedule_[reportStepIdx].has_gpmaint()) {
486 for (
const auto& calculator : regionalAveragePressureCalculator_) {
487 calculator.second->template defineState<ElementContext>(simulator_);
489 const double dt = simulator_.timeStepSize();
490 const Group& fieldGroup = this->schedule().getGroup(
"FIELD", reportStepIdx);
492 this->groupStateHelper().updateGpMaintTargetForGroups(fieldGroup,
493 regionalAveragePressureCalculator_,
499 this->updateAndCommunicateGroupData(reportStepIdx,
true);
502 for (
auto& well : well_container_) {
503 const uint64_t effective_events_mask = ScheduleEvents::WELL_STATUS_CHANGE
504 + ScheduleEvents::INJECTION_TYPE_CHANGED
505 + ScheduleEvents::WELL_SWITCHED_INJECTOR_PRODUCER
506 + ScheduleEvents::NEW_WELL;
508 const auto& events = this->schedule()[reportStepIdx].wellgroup_events();
509 const bool event = this->report_step_starts_ && events.hasEvent(well->name(), effective_events_mask);
510 const bool dyn_status_change = this->wellState().well(well->name()).status
511 != this->prevWellState().well(well->name()).status;
513 if (event || dyn_status_change || slave_needs_well_solution) {
515 well->scaleSegmentRatesAndPressure(this->wellState());
516 well->calculateExplicitQuantities(simulator_, this->groupStateHelper());
517 well->updateWellStateWithTarget(simulator_, this->groupStateHelper(), this->wellState());
518 well->updatePrimaryVariables(this->groupStateHelper());
519 well->solveWellEquation(
520 simulator_, this->groupStateHelper(), this->wellState()
522 }
catch (
const std::exception& e) {
523 const std::string msg =
"Compute initial well solution for new well " + well->name() +
" failed. Continue with zero initial rates";
524 local_deferredLogger.warning(
"WELL_INITIAL_SOLVE_FAILED", msg);
532#ifdef RESERVOIR_COUPLING_ENABLED
533 if (slave_needs_well_solution) {
535 this->updateAndCommunicateGroupData(reportStepIdx,
false);
536 this->rescoupHelper_.sendSlaveGroupDataToMaster();
538 else if (this->isReservoirCouplingMaster()) {
539 if (this->reservoirCouplingMaster().isFirstSubstepOfSyncTimestep()) {
540 this->rescoupHelper_.sendMasterGroupConstraintsToSlaves();
541 this->rescoupHelper_.sendCoupledNetworkActiveStatus();
542 this->rescoupHelper_.receiveSlaveGroupData();
548 const std::string msg =
"Compute initial well solution for new wells failed. Continue with zero initial rates";
549 local_deferredLogger.warning(
"WELL_INITIAL_SOLVE_FAILED", msg);
552 const auto& comm = simulator_.vanguard().grid().comm();
554 exc_type,
"beginTimeStep() failed: " + exc_msg, this->terminal_output_, comm);
558 template<
typename TypeTag>
561 const double simulationTime,
564 for (
const std::string& well_name : this->getWellsForTesting(timeStepIdx, simulationTime)) {
565 const Well& wellEcl = this->schedule().getWell(well_name, timeStepIdx);
566 if (wellEcl.getStatus() == Well::Status::SHUT)
569 WellInterfacePtr well = createWellForWellTest(well_name, timeStepIdx, deferred_logger);
571 well->init(depth_, gravity_, B_avg_,
true);
573 Scalar well_efficiency_factor = wellEcl.getEfficiencyFactor() *
574 this->wellState().getGlobalEfficiencyScalingFactor(well_name);
575 this->groupStateHelper().accumulateGroupEfficiencyFactor(
576 this->schedule().getGroup(wellEcl.groupName(), timeStepIdx),
577 well_efficiency_factor
580 well->setWellEfficiencyFactor(well_efficiency_factor);
581 well->setVFPProperties(this->vfp_properties_.get());
582 well->setGuideRate(&this->guideRate_);
585 if (well->isProducer() && alternative_well_rate_init_) {
586 well->initializeProducerWellState(simulator_, this->wellState(), deferred_logger);
588 if (well->isVFPActive(deferred_logger)) {
589 well->setPrevSurfaceRates(this->wellState(), this->prevWellState());
592 const auto& network = this->schedule()[timeStepIdx].network();
593 if (network.active()) {
594 this->network_.initializeWell(*well);
598 GLiftEclWells ecl_well_map;
599 gaslift_.initGliftEclWellMap(well_container_, ecl_well_map);
600 well->wellTesting(simulator_,
602 this->groupStateHelper(),
604 this->wellTestState(),
606 this->well_open_times_);
607 }
catch (
const std::exception& e) {
608 const std::string msg =
609 fmt::format(fmt::runtime(
"Exception during testing of well: {}. The well will not open.\n"
610 "Exception message: {}"), wellEcl.name(), e.what());
611 deferred_logger.
warning(
"WELL_TESTING_FAILED", msg);
617 template<
typename TypeTag>
623 for (
auto&& pinfo : this->local_parallel_well_info_)
634 template<
typename TypeTag>
644 template<
typename TypeTag>
649 this->closed_this_step_.clear();
652 this->report_step_starts_ =
false;
653 const int reportStepIdx = simulator_.episodeIndex();
655 auto logger_guard = this->groupStateHelper().pushLogger();
656 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
657 for (
const auto& well : well_container_) {
658 if (getPropValue<TypeTag, Properties::EnablePolymerMW>() && well->isInjector()) {
659 well->updateWaterThroughput(dt, this->wellState());
663 for (
const auto& well : well_container_) {
664 well->updateConnectionTransmissibilityFactor(simulator_, this->wellState().well(well->indexOfWell()));
665 well->updateConnectionDFactor(simulator_, this->wellState().well(well->indexOfWell()));
668 if (Indices::waterEnabled) {
669 this->updateFiltrationModelsPostStep(dt, FluidSystem::waterPhaseIdx, local_deferredLogger);
673 this->updateInjMult(local_deferredLogger);
676 for (
const auto& well : well_container_) {
677 well->reportWellSwitching(this->wellState().well(well->indexOfWell()), local_deferredLogger);
680 if (this->terminal_output_) {
681 this->reportGroupSwitching(local_deferredLogger);
685 rateConverter_->template defineState<ElementContext>(simulator_);
689 this->updateWellPotentials(reportStepIdx,
691 simulator_.vanguard().summaryConfig(),
692 local_deferredLogger);
693 }
catch ( std::runtime_error& e ) {
694 const std::string msg =
"A zero well potential is returned for output purposes. ";
695 local_deferredLogger.warning(
"WELL_POTENTIAL_CALCULATION_FAILED", msg);
698 updateWellTestState(simulationTime, this->wellTestState());
701 const Group& fieldGroup = this->schedule_.getGroup(
"FIELD", reportStepIdx);
702 this->checkGEconLimits(fieldGroup, simulationTime,
703 simulator_.episodeIndex(), local_deferredLogger);
704 this->checkGconsaleLimits(fieldGroup, this->wellState(),
705 simulator_.episodeIndex(), local_deferredLogger);
707 this->calculateProductivityIndexValues(local_deferredLogger);
709 this->groupStateHelper().updateNONEProductionGroups();
711#ifdef RESERVOIR_COUPLING_ENABLED
712 this->rescoupHelper_.rescoupSyncSummaryData();
714 this->commitWGState();
717 this->computeWellTemperature();
721 template<
typename TypeTag>
725 unsigned elemIdx)
const
729 if (!is_cell_perforated_[elemIdx] || cellRates_.count(elemIdx) == 0) {
733 rate = cellRates_.at(elemIdx);
737 template<
typename TypeTag>
738 template <
class Context>
742 const Context& context,
744 unsigned timeIdx)
const
747 int elemIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
749 if (!is_cell_perforated_[elemIdx] || cellRates_.count(elemIdx) == 0) {
753 rate = cellRates_.at(elemIdx);
758 template<
typename TypeTag>
763 const auto pressIx = []()
765 if (Indices::oilEnabled) {
return FluidSystem::oilPhaseIdx; }
766 if (Indices::waterEnabled) {
return FluidSystem::waterPhaseIdx; }
768 return FluidSystem::gasPhaseIdx;
771 auto cellPressures = std::vector<Scalar>(this->local_num_cells_,
Scalar{0});
772 auto cellTemperatures = std::vector<Scalar>(this->local_num_cells_,
Scalar{0});
775 const auto& gridView = this->simulator_.vanguard().gridView();
778 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
779 elemCtx.updatePrimaryStencil(elem);
780 elemCtx.updatePrimaryIntensiveQuantities(0);
782 const auto ix = elemCtx.globalSpaceIndex(0, 0);
783 const auto& fs = elemCtx.intensiveQuantities(0, 0).fluidState();
785 cellPressures[ix] = fs.pressure(pressIx).value();
786 cellTemperatures[ix] = fs.temperature(0).value();
789 this->simulator_.vanguard().grid().comm());
791 this->wellState().init(cellPressures, cellTemperatures, this->schedule(), this->wells_ecl_,
792 this->local_parallel_well_info_, timeStepIdx,
793 &this->prevWellState(), this->well_perf_data_,
794 this->summaryState(), simulator_.vanguard().enableDistributedWells());
801 template<
typename TypeTag>
806 auto logger_guard = this->groupStateHelper().pushLogger();
807 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
809 const int nw = this->numLocalWells();
811 well_container_.clear();
814 well_container_.reserve(nw);
816 const auto& wmatcher = this->schedule().wellMatcher(report_step);
817 const auto& wcycle = this->schedule()[report_step].wcycle.get();
821 std::ranges::for_each(this->wells_ecl_,
822 [
this, &wg_events = this->report_step_start_events_](
const auto& well_ecl)
824 if (!well_ecl.hasConnections()) {
829 constexpr auto events_mask = ScheduleEvents::WELL_STATUS_CHANGE |
830 ScheduleEvents::REQUEST_OPEN_WELL |
831 ScheduleEvents::REQUEST_SHUT_WELL;
832 const bool well_event =
833 this->report_step_starts_ &&
834 wg_events.hasEvent(well_ecl.name(), events_mask);
843 if (well_ecl.getStatus() == WellStatus::OPEN) {
844 this->well_open_times_.insert_or_assign(well_ecl.name(),
845 this->simulator_.time());
846 this->well_close_times_.erase(well_ecl.name());
847 }
else if (well_ecl.getStatus() == WellStatus::SHUT) {
848 this->well_close_times_.insert_or_assign(well_ecl.name(),
849 this->simulator_.time());
850 this->well_open_times_.erase(well_ecl.name());
856 const auto cycle_states = wcycle.wellStatus(this->simulator_.time(),
858 this->well_open_times_,
859 this->well_close_times_);
861 for (
int w = 0; w < nw; ++w) {
862 const Well& well_ecl = this->wells_ecl_[w];
864 if (!well_ecl.hasConnections()) {
869 const std::string& well_name = well_ecl.name();
870 const auto well_status = this->schedule()
871 .getWell(well_name, report_step).getStatus();
873 const bool shut_event = this->wellState().well(w).events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE)
874 && well_status == Well::Status::SHUT;
875 const bool open_event = this->wellState().well(w).events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE)
876 && well_status == Well::Status::OPEN;
877 const auto& ws = this->wellState().well(well_name);
879 if (shut_event && ws.status != Well::Status::SHUT) {
880 this->closed_this_step_.insert(well_name);
881 this->wellState().shutWell(w);
882 }
else if (open_event && ws.status != Well::Status::OPEN) {
883 this->wellState().openWell(w);
887 if (this->wellTestState().well_is_closed(well_name)) {
892 const bool closed_this_step = (this->wellTestState().lastTestTime(well_name) == simulator_.time());
895 auto& events = this->wellState().well(w).events;
896 if (events.hasEvent(ScheduleEvents::REQUEST_OPEN_WELL)) {
897 if (!closed_this_step) {
898 this->wellTestState().open_well(well_name);
899 this->wellTestState().open_completions(well_name);
900 this->well_open_times_.insert_or_assign(well_name,
901 this->simulator_.time());
902 this->well_close_times_.erase(well_name);
904 events.clearEvent(ScheduleEvents::REQUEST_OPEN_WELL);
914 auto& well_test_state = this->wellTestState();
915 const auto& sched_state = this->schedule()[report_step];
916 const bool may_reopen_completions =
917 ws.status == Well::Status::OPEN &&
918 sched_state.events().hasEvent(ScheduleEvents::COMPLETION_CHANGE) &&
919 well_test_state.num_closed_completions() > 0;
921 if (may_reopen_completions) {
922 const auto& completion_events = sched_state.wellcompletion_events();
923 for (
const auto& connection : well_ecl.getConnections()) {
924 const int complnum = connection.complnum();
928 if (!well_test_state.completion_is_closed(well_name, complnum) ||
929 !completion_events.hasEvent(well_name, complnum, ScheduleEvents::REQUEST_OPEN_COMPLETION)) {
934 const bool closed_this_step =
935 (well_test_state.lastCompletionCloseTime(well_name, complnum) == simulator_.time());
936 if (closed_this_step) {
940 well_test_state.open_completion(well_name, complnum);
941 local_deferredLogger.info(
942 fmt::format(
"Completion {} - block ({}, {}, {}) for well {} "
943 "is reopened due to an explicit WELOPEN/COMPDAT "
946 connection.getI() + 1,
947 connection.getJ() + 1,
948 connection.getK() + 1,
955 if (this->wellTestState().well_is_closed(well_name))
957 if (well_ecl.getAutomaticShutIn() ||
958 !well_ecl.getAllowCrossFlow() ||
959 this->allConnectionsClosed(well_ecl))
961 this->wellState().shutWell(w);
962 this->well_close_times_.erase(well_name);
963 this->well_open_times_.erase(well_name);
967 this->wellState().stopWell(w);
971 if (!well_ecl.getAllowCrossFlow()) {
972 const bool any_zero_rate_constraint = well_ecl.isProducer()
973 ? well_ecl.productionControls(this->summaryState_).anyZeroRateConstraint()
974 : well_ecl.injectionControls(this->summaryState_).anyZeroRateConstraint();
975 if (any_zero_rate_constraint) {
977 local_deferredLogger.debug(fmt::format(fmt::runtime(
" Well {} gets shut due to having zero rate constraint and disallowing crossflow "), well_ecl.name()));
978 this->wellState().shutWell(w);
979 this->well_close_times_.erase(well_name);
980 this->well_open_times_.erase(well_name);
985 if (!wcycle.empty()) {
986 const auto it = cycle_states.find(well_name);
987 if (it != cycle_states.end()) {
988 if (!it->second || well_status == Well::Status::SHUT) {
990 if (well_status == Well::Status::SHUT) {
991 this->well_open_times_.erase(well_name);
992 this->well_close_times_.erase(well_name);
994 this->wellState().shutWell(w);
997 this->wellState().openWell(w);
1003 if (ws.status == Well::Status::SHUT) {
1007 well_container_.emplace_back(this->createWellPointer(w, report_step));
1009 if (ws.status == Well::Status::STOP) {
1010 well_container_.back()->stopWell();
1011 this->well_close_times_.erase(well_name);
1012 this->well_open_times_.erase(well_name);
1016 if (!wcycle.empty()) {
1017 const auto schedule_open =
1018 [&wg_events = this->report_step_start_events_](
const std::string& name)
1020 return wg_events.hasEvent(name, ScheduleEvents::REQUEST_OPEN_WELL);
1022 for (
const auto& [wname, wscale] : wcycle.efficiencyScale(this->simulator_.time(),
1023 this->simulator_.timeStepSize(),
1025 this->well_open_times_,
1028 this->wellState().updateEfficiencyScalingFactor(wname, wscale);
1029 this->schedule_.add_event(ScheduleEvents::WELLGROUP_EFFICIENCY_UPDATE, report_step);
1034 this->well_container_generic_.clear();
1035 for (
auto& w : well_container_) {
1036 this->well_container_generic_.push_back(w.get());
1039 this->network_.initialize(report_step);
1041 this->wbp_.registerOpenWellsForWBPCalculation();
1048 template <
typename TypeTag>
1053 const auto is_multiseg = this->wells_ecl_[wellID].isMultiSegment();
1055 if (! (this->param_.use_multisegment_well_ && is_multiseg)) {
1056 return this->
template createTypedWellPointer<StandardWell<TypeTag>>(wellID, report_step);
1059 return this->
template createTypedWellPointer<MultisegmentWell<TypeTag>>(wellID, report_step);
1067 template <
typename TypeTag>
1068 template <
typename WellType>
1069 std::unique_ptr<WellType>
1074 const auto& perf_data = this->well_perf_data_[wellID];
1077 const auto pvtreg = perf_data.empty()
1078 ? 0 : this->pvt_region_idx_[perf_data.front().cell_index];
1080 const auto& parallel_well_info = this->local_parallel_well_info_[wellID].get();
1081 const auto global_pvtreg = parallel_well_info.broadcastFirstPerforationValue(pvtreg);
1083 return std::make_unique<WellType>(this->wells_ecl_[wellID],
1087 *this->rateConverter_,
1089 this->numConservationQuantities(),
1099 template<
typename TypeTag>
1103 const int report_step,
1108 std::ranges::find_if(this->wells_ecl_,
1109 [&well_name](
const auto& w)
1110 {
return well_name == w.name(); });
1112 if (it == this->wells_ecl_.end()) {
1114 fmt::format(fmt::runtime(
"Could not find well {} in wells_ecl"), well_name),
1118 const int pos =
static_cast<int>(std::distance(this->wells_ecl_.begin(), it));
1119 return this->createWellPointer(pos, report_step);
1124 template<
typename TypeTag>
1130 auto logger_guard = this->groupStateHelper().pushLogger();
1131 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
1133 const auto& iterCtx = simulator_.problem().iterationContext();
1136 if (gaslift_.terminalOutput()) {
1137 const std::string msg =
1138 fmt::format(fmt::runtime(
"assemble() : iteration {}"), iterCtx.iteration());
1139 gaslift_.gliftDebug(msg, local_deferredLogger);
1143 Dune::Timer perfTimer;
1145 this->closed_offending_wells_.clear();
1148 if (iterCtx.needsTimestepInit()) {
1149 this->updateNetworkActiveState_();
1151 const int episodeIdx = simulator_.episodeIndex();
1152 const auto& network = this->schedule()[episodeIdx].network();
1153 if (!this->wellsActive() && !network.active()) {
1159 if (iterCtx.needsTimestepInit() && this->wellsActive()) {
1160 OPM_TIMEBLOCK(firstIterationAssemble);
1167 calculateExplicitQuantities();
1168 prepareTimeStep(local_deferredLogger);
1171 "assemble() failed during well initialization: ",
1172 this->terminal_output_, grid().comm());
1175 const bool well_group_control_changed = updateWellControlsAndNetwork(
1178 local_deferredLogger);
1182 if ( ! this->wellsActive() ) {
1186 assembleWellEqWithoutIteration(dt);
1192 last_report_.well_group_control_changed = well_group_control_changed;
1193 last_report_.assemble_time_well += perfTimer.stop();
1199 template<
typename TypeTag>
1208 bool do_network_update =
true;
1209 bool well_group_control_changed =
false;
1210 Scalar network_imbalance = 0.0;
1212 const std::size_t iteration_to_relax = param_.network_max_strict_outer_iterations_;
1214 const std::size_t max_iteration = param_.network_max_outer_iterations_;
1215 std::size_t network_update_iteration = 0;
1216 network_needs_more_balancing_force_another_newton_iteration_ =
false;
1217 while (do_network_update) {
1218 if (!this->isRescoupSlaveCoupledNetworkIteration_()
1219 && network_update_iteration >= max_iteration ) {
1221 const int episodeIdx = simulator_.episodeIndex();
1222 if (this->network_.willBalanceOnNextIteration(episodeIdx)) {
1223 if (this->terminal_output_) {
1224 const std::string msg = fmt::format(
"Maximum of {:d} network iterations has been used and we stop the update, \n"
1225 "and try again after the next Newton iteration (imbalance = {:.2e} bar)",
1226 max_iteration, network_imbalance*1.0e-5);
1227 local_deferredLogger.
debug(msg);
1231 network_needs_more_balancing_force_another_newton_iteration_ =
true;
1233 if (this->terminal_output_) {
1234 const std::string msg = fmt::format(
"Maximum of {:d} network iterations has been used and we stop the update. \n"
1235 "The simulator will continue with unconverged network results (imbalance = {:.2e} bar)",
1236 max_iteration, network_imbalance*1.0e-5);
1237 local_deferredLogger.
info(msg);
1242 if (this->terminal_output_ && (network_update_iteration == iteration_to_relax) ) {
1243 local_deferredLogger.
debug(
"We begin using relaxed tolerance for network update now after " +
std::to_string(iteration_to_relax) +
" iterations ");
1245 const bool relax_network_balance = network_update_iteration >= iteration_to_relax;
1247 const bool optimize_gas_lift = ( (network_update_iteration + 1) < std::max(max_iteration,
static_cast<std::size_t
>(2)) );
1248 std::tie(well_group_control_changed, do_network_update, network_imbalance) =
1249 updateWellControlsAndNetworkIteration(mandatory_network_balance, relax_network_balance, optimize_gas_lift, dt,local_deferredLogger);
1250 ++network_update_iteration;
1252 if (this->isRescoupMasterCoupledNetworkIteration_()) {
1253 this->sendSlaveNetworkLoopTerminationSignal_();
1255 return well_group_control_changed;
1261 template<
typename TypeTag>
1262 std::tuple<bool, bool, typename BlackoilWellModel<TypeTag>::Scalar>
1265 const bool relax_network_tolerance,
1266 const bool optimize_gas_lift,
1271 const int reportStepIdx = simulator_.episodeIndex();
1273#ifdef RESERVOIR_COUPLING_ENABLED
1274 if (this->isRescoupSlaveCoupledNetworkIteration_()) {
1275 this->rescoupHelper_.receiveMasterGroupNodePressuresFromMaster();
1278 this->updateAndCommunicateGroupData(reportStepIdx,
true);
1283 bool well_group_control_changed = updateWellControls(local_deferredLogger);
1284 const auto [more_inner_network_update, network_imbalance] =
1285 this->network_.update(mandatory_network_balance,
1286 local_deferredLogger,
1287 relax_network_tolerance);
1288#ifdef RESERVOIR_COUPLING_ENABLED
1289 if (this->isReservoirCouplingMaster()) {
1290 this->rescoupHelper_.maybeExchangeNetworkOuterIterationWithSlaves(more_inner_network_update);
1294 bool alq_updated =
false;
1297 if (optimize_gas_lift) {
1300 const bool updatePotentials = (this->network_.shouldBalance(reportStepIdx) ||
1301 mandatory_network_balance);
1302 alq_updated = gaslift_.maybeDoGasLiftOptimize(simulator_,
1304 this->network_.nodePressures(),
1308 local_deferredLogger);
1310 prepareWellsBeforeAssembling(dt);
1313 "updateWellControlsAndNetworkIteration() failed: ",
1314 this->terminal_output_, grid().comm());
1318 guideRateUpdateIsNeeded(reportStepIdx)) {
1319 const double simulationTime = simulator_.time();
1323 this->guide_rate_handler_.updateGuideRates(
1324 reportStepIdx, simulationTime, this->wellState(), this->groupState()
1329 bool more_network_update = this->network_.shouldBalance(reportStepIdx) &&
1330 (more_inner_network_update || alq_updated);
1332 if (this->isRescoupSlaveOnSyncStepFirstSubstep_()
1333 && this->isRescoupSlaveConnectedToMasterNetwork_()) {
1340 const bool more_cross_rescoup_update =
1341 this->maybeSendSlaveGroupFlowToMaster_(reportStepIdx);
1342 more_network_update = more_network_update || more_cross_rescoup_update;
1345 return {well_group_control_changed, more_network_update, network_imbalance};
1348 template<
typename TypeTag>
1354 for (
auto& well : well_container_) {
1355 well->assembleWellEq(simulator_, dt, this->groupStateHelper(), this->wellState());
1360 template<
typename TypeTag>
1366 for (
auto& well : well_container_) {
1367 well->prepareWellBeforeAssembling(
1368 simulator_, dt, this->groupStateHelper(), this->wellState()
1374 template<
typename TypeTag>
1380 auto& deferred_logger = this->groupStateHelper().deferredLogger();
1385 for (
auto& well: well_container_) {
1386 well->assembleWellEqWithoutIteration(simulator_, this->groupStateHelper(), dt, this->wellState(),
1390 this->terminal_output_, grid().comm());
1394 template<
typename TypeTag>
1401 for (
const auto& well : well_container_) {
1402 well->addCellRates(cellRates_);
1406 template<
typename TypeTag>
1413 for (
const auto& well : well_container_) {
1414 const auto it = well_domain_map.find(well->name());
1415 if (it != well_domain_map.end() && it->second == domainIndex) {
1416 well->addCellRates(cellRates_);
1421#if COMPILE_GPU_BRIDGE
1422 template<
typename TypeTag>
1430 for(
unsigned int i = 0; i < well_container_.size(); i++){
1431 auto& well = well_container_[i];
1434 wellContribs.
addNumBlocks(derived->linSys().getNumBlocks());
1439 wellContribs.
alloc();
1441 for(
unsigned int i = 0; i < well_container_.size(); i++){
1442 auto& well = well_container_[i];
1444 auto derived_std =
dynamic_cast<StandardWell<TypeTag>*
>(well.get());
1446 derived_std->linSys().extract(derived_std->numStaticWellEq, wellContribs);
1448 auto derived_ms =
dynamic_cast<MultisegmentWell<TypeTag>*
>(well.get());
1450 derived_ms->linSys().extract(wellContribs);
1452 OpmLog::warning(
"Warning unknown type of well");
1459 template<
typename TypeTag>
1464 for (
const auto& well: well_container_ ) {
1465 well->addWellContributions(jacobian);
1469 template<
typename TypeTag>
1474 const bool use_well_weights)
const
1476 int nw = this->numLocalWellsEnd();
1477 int rdofs = local_num_cells_;
1478 for (
int i = 0; i < nw; i++ ) {
1479 int wdof = rdofs + i;
1480 jacobian[wdof][wdof] = 1.0;
1483 for (
const auto& well : well_container_) {
1484 well->addWellPressureEquations(jacobian,
1492 template <
typename TypeTag>
1495 const std::vector<typename SparseMatrixAdapter::MatrixBlock*>& diagMatAddress)
const
1500 for (
const auto& well : well_container_) {
1501 if (!well->isOperableAndSolvable() && !well->wellIsStopped()) {
1504 const auto& cells = well->cells();
1505 const auto& rates = well->connectionRates();
1506 for (
unsigned perfIdx = 0; perfIdx < rates.size(); ++perfIdx) {
1507 unsigned cellIdx = cells[perfIdx];
1508 auto rate = rates[perfIdx];
1511 using MatrixBlockType =
typename SparseMatrixAdapter::MatrixBlock;
1512 MatrixBlockType bMat(0.0);
1513 simulator_.model().linearizer().setResAndJacobi(res, bMat, rate);
1514 residual[cellIdx] += res;
1515 *diagMatAddress[cellIdx] += bMat;
1521 template<
typename TypeTag>
1526 int nw = this->numLocalWellsEnd();
1527 int rdofs = local_num_cells_;
1528 const auto wellconnections = this->getMaxWellConnections();
1529 for (
int i = 0; i < nw; ++i) {
1530 int wdof = rdofs + i;
1531 jacobian.entry(wdof,wdof) = 0.0;
1532 const auto& perfcells = wellconnections[i];
1533 for (
int perfcell : perfcells) {
1534 jacobian.entry(wdof, perfcell) = 0.0;
1535 jacobian.entry(perfcell, wdof) = 0.0;
1541 template<
typename TypeTag>
1546 auto loggerGuard = this->groupStateHelper().pushLogger();
1549 for (
const auto& well : well_container_) {
1550 const auto& cells = well->cells();
1551 x_local_.resize(cells.size());
1553 for (
size_t i = 0; i < cells.size(); ++i) {
1554 x_local_[i] = x[cells[i]];
1556 well->recoverWellSolutionAndUpdateWellState(simulator_, x_local_,
1557 this->groupStateHelper(), this->wellState());
1561 simulator_.vanguard().grid().comm());
1565 template<
typename TypeTag>
1571 OPM_THROW(std::logic_error,
"Attempt to call NLDD method without a NLDD solver");
1574 return nldd_->recoverWellSolutionAndUpdateWellState(x, domainIdx);
1578 template<
typename TypeTag>
1581 getWellConvergence(
const std::vector<Scalar>& B_avg,
bool checkWellGroupControlsAndNetwork)
const
1585 const auto& iterCtx = simulator_.problem().iterationContext();
1586 const bool relaxTolerance = iterCtx.shouldRelax(param_.strict_outer_iter_wells_ + 1);
1588 auto logger_guard = this->groupStateHelper().pushLogger();
1589 for (
const auto& well : well_container_) {
1590 if (well->isOperableAndSolvable() || well->wellIsStopped()) {
1591 local_report += well->getWellConvergence(
1592 this->groupStateHelper(), B_avg,
1597 report.
setWellFailed({CR::WellFailure::Type::Unsolvable, CR::Severity::Normal, -1, well->name()});
1598 local_report += report;
1606 if (checkWellGroupControlsAndNetwork) {
1612 if (this->terminal_output_) {
1616 OpmLog::debug(
"NaN residual found with phase " +
std::to_string(f.phase()) +
" for well " + f.wellName());
1618 OpmLog::debug(
"Too large residual found with phase " +
std::to_string(f.phase()) +
" for well " + f.wellName());
1629 template<
typename TypeTag>
1635 for (
auto& well : well_container_) {
1636 well->calculateExplicitQuantities(simulator_, this->groupStateHelper());
1644 template<
typename TypeTag>
1650 if (!this->wellsActive()) {
1653 const int episodeIdx = simulator_.episodeIndex();
1654 const auto& comm = simulator_.vanguard().grid().comm();
1656 bool changed_well_group =
false;
1657 const Group& fieldGroup = this->schedule().getGroup(
"FIELD", episodeIdx);
1660 const std::size_t max_iter = param_.well_group_constraints_max_iterations_;
1661 while(!changed_well_group && iter < max_iter) {
1662 changed_well_group = updateGroupControls(fieldGroup, deferred_logger, episodeIdx);
1665 bool changed_well_to_group =
false;
1667 OPM_TIMEBLOCK(UpdateWellControls);
1671 for (
const auto& well : well_container_) {
1674 simulator_, mode, this->groupStateHelper(), this->wellState()
1677 changed_well_to_group = changed_well || changed_well_to_group;
1681 simulator_.gridView().comm());
1684 changed_well_to_group = comm.sum(
static_cast<int>(changed_well_to_group));
1685 if (changed_well_to_group) {
1686 updateAndCommunicate(episodeIdx);
1687 changed_well_group =
true;
1691 bool changed_well_individual =
false;
1696 for (
const auto& well : well_container_) {
1699 simulator_, mode, this->groupStateHelper(), this->wellState()
1702 changed_well_individual = changed_well || changed_well_individual;
1706 simulator_.gridView().comm());
1709 changed_well_individual = comm.sum(
static_cast<int>(changed_well_individual));
1710 if (changed_well_individual) {
1711 updateAndCommunicate(episodeIdx);
1712 changed_well_group =
true;
1718 this->updateWsolvent(fieldGroup, episodeIdx, this->nupcolWellState());
1720 return changed_well_group;
1724 template<
typename TypeTag>
1729 this->updateAndCommunicateGroupData(reportStepIdx,
true);
1735 for (
const auto& well : well_container_) {
1737 const auto& ws = this->wellState().well(well->indexOfWell());
1738 if (ws.production_cmode == Well::ProducerCMode::GRUP ||
1739 ws.injection_cmode == Well::InjectorCMode::GRUP)
1741 well->updateWellStateWithTarget(
1742 simulator_, this->groupStateHelper(), this->wellState()
1747 simulator_.gridView().comm())
1748 this->updateAndCommunicateGroupData(reportStepIdx,
true);
1751 template<
typename TypeTag>
1756 const int reportStepIdx)
1765 if (this->isReservoirCouplingMasterGroup(group.name())) {
1768 const auto& iterCtx = simulator_.problem().iterationContext();
1769 bool changed =
false;
1771 const int nupcol = this->schedule()[reportStepIdx].nupcol();
1772 const bool update_group_switching_log = !iterCtx.withinNupcol(nupcol);
1773 const bool changed_hc = this->checkGroupHigherConstraints(
1774 group, deferred_logger, reportStepIdx, update_group_switching_log);
1777 updateAndCommunicate(reportStepIdx);
1780 bool changed_individual =
1782 updateGroupIndividualControl(group,
1784 param_.max_number_of_group_switches_,
1785 update_group_switching_log,
1786 this->switched_inj_groups_,
1787 this->switched_prod_groups_,
1788 this->closed_offending_wells_,
1793 if (changed_individual) {
1795 updateAndCommunicate(reportStepIdx);
1798 for (
const std::string& groupName : group.groups()) {
1799 bool changed_this = updateGroupControls(this->schedule().getGroup(groupName, reportStepIdx), deferred_logger, reportStepIdx);
1800 changed = changed || changed_this;
1805 template<
typename TypeTag>
1811 auto logger_guard = this->groupStateHelper().pushLogger();
1812 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
1813 for (
const auto& well : well_container_) {
1814 const auto& wname = well->name();
1815 const auto wasClosed = wellTestState.well_is_closed(wname);
1816 well->checkWellOperability(simulator_,
1818 this->groupStateHelper());
1819 const bool under_zero_target =
1820 well->wellUnderZeroGroupRateTarget(this->groupStateHelper());
1821 well->updateWellTestState(this->wellState().well(wname),
1826 local_deferredLogger);
1828 if (!wasClosed && wellTestState.well_is_closed(wname)) {
1829 this->closed_this_step_.insert(wname);
1832 const WellEconProductionLimits& econ_production_limits = well->wellEcl().getEconLimits();
1833 if (econ_production_limits.validFollowonWell()) {
1834 const auto episode_idx = simulator_.episodeIndex();
1835 const auto follow_on_well = econ_production_limits.followonWell();
1836 if (!this->schedule().hasWell(follow_on_well, episode_idx)) {
1837 const auto msg = fmt::format(
"Well {} was closed. But the given follow on well {} does not exist."
1838 "The simulator continues without opening a follow on well.",
1839 wname, follow_on_well);
1840 local_deferredLogger.warning(msg);
1842 auto& ws = this->wellState().well(follow_on_well);
1843 const bool success = ws.updateStatus(WellStatus::OPEN);
1845 const auto msg = fmt::format(
"Well {} was closed. The follow on well {} opens instead.", wname, follow_on_well);
1846 local_deferredLogger.info(msg);
1848 const auto msg = fmt::format(
"Well {} was closed. The follow on well {} is already open.", wname, follow_on_well);
1849 local_deferredLogger.warning(msg);
1856 for (
const auto& [group_name, to] : this->closed_offending_wells_) {
1857 if (this->hasOpenLocalWell(to.second) &&
1858 !this->wasDynamicallyShutThisTimeStep(to.second))
1860 wellTestState.close_well(to.second,
1861 WellTestConfig::Reason::GROUP,
1863 this->updateClosedWellsThisStep(to.second);
1864 const std::string msg =
1865 fmt::format(
"Procedure on exceeding {} limit is WELL for group {}. "
1871 local_deferredLogger.info(msg);
1877 template<
typename TypeTag>
1881 std::string& exc_msg,
1885 const int np = this->numPhases();
1886 std::vector<Scalar> potentials;
1887 const auto& well = well_container_[widx];
1888 std::string cur_exc_msg;
1891 well->computeWellPotentials(simulator_, well_state_copy, this->groupStateHelper(), potentials);
1896 exc_msg += fmt::format(
"\nFor well {}: {}", well->name(), cur_exc_msg);
1898 exc_type = std::max(exc_type, cur_exc_type);
1902 auto& ws = this->wellState().well(well->indexOfWell());
1903 for (
int p = 0; p < np; ++p) {
1905 ws.well_potentials[p] = std::max(
Scalar{0.0}, potentials[p]);
1911 template <
typename TypeTag>
1916 for (
const auto& wellPtr : this->well_container_) {
1917 this->calculateProductivityIndexValues(wellPtr.get(), deferred_logger);
1925 template <
typename TypeTag>
1936 for (
const auto& shutWell : this->local_shut_wells_) {
1937 if (!this->wells_ecl_[shutWell].hasConnections()) {
1942 auto wellPtr = this->
template createTypedWellPointer
1945 wellPtr->
init(this->depth_, this->gravity_, this->B_avg_,
true);
1947 this->calculateProductivityIndexValues(wellPtr.get(), deferred_logger);
1955 template <
typename TypeTag>
1969 template<
typename TypeTag>
1974 const auto episodeIdx = simulator_.episodeIndex();
1977 const bool do_prestep_network_rebalance = this->shouldDoPreStepNetworkRebalance_(episodeIdx);
1979 for (
const auto& well : well_container_) {
1980 auto& events = this->wellState().well(well->indexOfWell()).events;
1982 well->updateWellStateWithTarget(
1983 simulator_, this->groupStateHelper(), this->wellState()
1985 well->updatePrimaryVariables(this->groupStateHelper());
1991 if (events.hasEvent(ScheduleEvents::REQUEST_OPEN_WELL)) {
1992 events.clearEvent(ScheduleEvents::REQUEST_OPEN_WELL);
1995 if (param_.solve_welleq_initially_ && well->isOperableAndSolvable()) {
1997 well->solveWellEquation(
1998 simulator_, this->groupStateHelper(), this->wellState()
2000 }
catch (
const std::exception& e) {
2001 const std::string msg =
"Compute initial well solution for " + well->name() +
" initially failed. Continue with the previous rates";
2002 deferred_logger.
warning(
"WELL_INITIAL_SOLVE_FAILED", msg);
2007 well->resetWellOperability();
2009 updatePrimaryVariables();
2012 if (do_prestep_network_rebalance) {
2013 network_.doPreStepRebalance(deferred_logger);
2017 template<
typename TypeTag>
2022 std::vector< Scalar > B_avg(numConservationQuantities(),
Scalar() );
2023 const auto& grid = simulator_.vanguard().grid();
2024 const auto& gridView = grid.leafGridView();
2028 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
2029 elemCtx.updatePrimaryStencil(elem);
2030 elemCtx.updatePrimaryIntensiveQuantities(0);
2032 const auto& intQuants = elemCtx.intensiveQuantities(0, 0);
2033 const auto& fs = intQuants.fluidState();
2035 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
2037 if (!FluidSystem::phaseIsActive(phaseIdx)) {
2041 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
2042 auto& B = B_avg[ compIdx ];
2044 B += 1 / fs.invB(phaseIdx).value();
2046 if constexpr (has_solvent_) {
2047 auto& B = B_avg[solventSaturationIdx];
2048 B += 1 / intQuants.solventInverseFormationVolumeFactor().value();
2054 grid.comm().sum(B_avg.data(), B_avg.size());
2055 B_avg_.resize(B_avg.size());
2056 std::ranges::transform(B_avg, B_avg_.begin(),
2057 [gcells = global_num_cells_](
const auto bval)
2058 { return bval / gcells; });
2065 template<
typename TypeTag>
2070 for (
const auto& well : well_container_) {
2071 well->updatePrimaryVariables(this->groupStateHelper());
2075 template<
typename TypeTag>
2079 const auto& grid = simulator_.vanguard().grid();
2080 const auto& eclProblem = simulator_.problem();
2081 const unsigned numCells = grid.size(0);
2083 this->pvt_region_idx_.resize(numCells);
2084 for (
unsigned cellIdx = 0; cellIdx < numCells; ++cellIdx) {
2085 this->pvt_region_idx_[cellIdx] =
2086 eclProblem.pvtRegionIndex(cellIdx);
2091 template<
typename TypeTag>
2104 return this->numPhases() + has_solvent_;
2107 template<
typename TypeTag>
2111 const auto& eclProblem = simulator_.problem();
2112 depth_.resize(local_num_cells_);
2113 for (
unsigned cellIdx = 0; cellIdx < local_num_cells_; ++cellIdx) {
2114 depth_[cellIdx] = eclProblem.dofCenterDepth(cellIdx);
2118 template<
typename TypeTag>
2121 getWell(
const std::string& well_name)
const
2125 std::ranges::find_if(well_container_,
2127 {
return elem->name() == well_name; });
2129 assert(well != well_container_.end());
2134 template <
typename TypeTag>
2139 return std::max(this->simulator_.episodeIndex(), 0);
2146 template<
typename TypeTag>
2151 const std::vector<Scalar>& production_rates,
2152 std::vector<Scalar>& resv_coeff)
const
2154 rateConverter_->calcCoeff(fipnum, pvtreg, production_rates, resv_coeff);
2157 template<
typename TypeTag>
2162 std::vector<Scalar>& resv_coeff)
const
2164 rateConverter_->calcInjCoeff(fipnum, pvtreg, resv_coeff);
2168 template <
typename TypeTag>
2173 if constexpr (energyModuleType_ == EnergyModules::FullyImplicitThermal) {
2174 const int np = this->numPhases();
2175 const int nw = this->numLocalWells();
2176 for (
auto wellID = 0*nw; wellID < nw; ++wellID) {
2177 const Well& well = this->wells_ecl_[wellID];
2178 auto& ws = this->wellState().well(wellID);
2179 if (well.isInjector()) {
2180 if (ws.status != WellStatus::STOP) {
2181 this->wellState().well(wellID).temperature = well.inj_temperature();
2185 std::array<Scalar,2> weighted{0.0,0.0};
2186 auto& [weighted_temperature, total_weight] = weighted;
2187 const auto& well_info = this->local_parallel_well_info_[wellID].get();
2188 using int_type =
decltype(this->well_perf_data_[wellID].size());
2189 for (int_type perf = 0, end_perf = this->well_perf_data_[wellID].size(); perf < end_perf; ++perf) {
2190 const int cell_idx = this->well_perf_data_[wellID][perf].cell_index;
2191 const auto& intQuants = simulator_.model().intensiveQuantities(cell_idx, 0);
2192 const auto& fs = intQuants.fluidState();
2193 Scalar weight_factor = computeTemperatureWeightFactor(perf, np, fs, ws);
2194 total_weight += weight_factor;
2195 weighted_temperature += weight_factor * fs.temperature(0).value();
2197 well_info.communication().sum(weighted.data(), 2);
2198 this->wellState().well(wellID).temperature = weighted_temperature / total_weight;
2204 template <
typename TypeTag>
2206 -> typename
WellState<Scalar,IndexTraits>::RsConstInfo
2208 if (! FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) ||
2209 ! FluidSystem::enableConstantRs())
2214 const auto& rsConstTables = this->eclState_
2215 .getTableManager().getRsconstTables();
2217 if (rsConstTables.empty() ||
2218 (rsConstTables[0].numRows() != std::size_t{1}))
2223 const auto rsConst = rsConstTables[0].getColumn(0).front();
2225 return {
true,
static_cast<Scalar
>(rsConst) };
2231 template <
typename TypeTag>
2232 void BlackoilWellModel<TypeTag>::
2233 assignWellTracerRates_(data::Wells& wsrpt)
const
2235 const auto reportStepIdx =
static_cast<unsigned int>(this->reportStepIndex());
2236 const auto& trMod = this->simulator_.problem().tracerModel();
2242 this->assignMswTracerRates(wsrpt, trMod.getMswTracerRates(), reportStepIdx);
2245 template <
typename TypeTag>
2246 void BlackoilWellModel<TypeTag>::
2247 assignWellSpeciesRates_(data::Wells& wsrpt)
const
2249 const auto reportStepIdx =
static_cast<unsigned int>(this->reportStepIndex());
2250 const auto& geochemMod = this->simulator_.problem().geochemistryModel();
2254 this->assignMswTracerRates(wsrpt, geochemMod.getMswSpeciesRates(), reportStepIdx);
2257 template <
typename TypeTag>
2258 bool BlackoilWellModel<TypeTag>::isRescoupMasterCoupledNetworkIteration_()
const
2260#ifdef RESERVOIR_COUPLING_ENABLED
2261 return this->rescoupHelper_.masterIsInCoupledNetworkIteration();
2267 template <
typename TypeTag>
2268 bool BlackoilWellModel<TypeTag>::isRescoupSlaveCoupledNetworkIteration_()
const
2279#ifdef RESERVOIR_COUPLING_ENABLED
2280 return this->isRescoupSlaveOnSyncStepFirstSubstep_()
2281 && this->reservoirCouplingSlave().connectedToMasterCoupledNetwork()
2282 && !this->reservoirCouplingSlave().lastReceivedMasterGroupNodePressuresIsFinal();
2288 template <
typename TypeTag>
2289 bool BlackoilWellModel<TypeTag>::isRescoupSlaveOnSyncStepFirstSubstep_()
const
2296#ifdef RESERVOIR_COUPLING_ENABLED
2297 return this->isReservoirCouplingSlave()
2298 && this->reservoirCouplingSlave().isFirstSubstepOfSyncTimestep();
2304 template <
typename TypeTag>
2305 bool BlackoilWellModel<TypeTag>::isRescoupSlaveConnectedToMasterNetwork_()
const
2307#ifdef RESERVOIR_COUPLING_ENABLED
2308 return this->isReservoirCouplingSlave()
2309 && this->reservoirCouplingSlave().connectedToMasterCoupledNetwork();
2315 template <
typename TypeTag>
2316 bool BlackoilWellModel<TypeTag>::
2317 maybeSendSlaveGroupFlowToMaster_([[maybe_unused]]
const int reportStepIdx)
2319#ifdef RESERVOIR_COUPLING_ENABLED
2326 assert(this->isReservoirCouplingSlave());
2327 const bool is_final =
2328 this->reservoirCouplingSlave().lastReceivedMasterGroupNodePressuresIsFinal();
2330 this->updateAndCommunicateGroupData(reportStepIdx,
false);
2331 this->rescoupHelper_.sendSlaveGroupDataToMaster();
2340 template <
typename TypeTag>
2341 void BlackoilWellModel<TypeTag>::
2342 sendSlaveNetworkLoopTerminationSignal_()
2344#ifdef RESERVOIR_COUPLING_ENABLED
2349 assert(this->isReservoirCouplingMaster());
2350 this->rescoupHelper_.sendMasterGroupNodePressuresToSlaves(
true);
2354 template <
typename TypeTag>
2355 bool BlackoilWellModel<TypeTag>::
2356 shouldDoPreStepNetworkRebalance_(
const int episodeIdx)
const
2369 return param_.pre_solve_network_
2370 && this->network_.needPreStepRebalance(episodeIdx)
2371 && !this->isRescoupCoupledNetworkParticipant_();
2374 template <
typename TypeTag>
2375 bool BlackoilWellModel<TypeTag>::isRescoupCoupledNetworkParticipant_()
const
2377#ifdef RESERVOIR_COUPLING_ENABLED
2378 if (this->isReservoirCouplingMaster()) {
2379 return this->rescoupHelper_.masterNetworkHasMasterGroupLeaves();
2381 if (this->isReservoirCouplingSlave()) {
2382 return this->reservoirCouplingSlave().connectedToMasterCoupledNetwork();
2390 template <
typename TypeTag>
2391 void BlackoilWellModel<TypeTag>::
2392 updateNetworkActiveState_()
2403 const int episodeIdx = this->simulator_.episodeIndex();
2404 this->network_.updateActiveState(episodeIdx);
#define OPM_END_PARALLEL_TRY_CATCH_LOG(obptc_logger, obptc_prefix, obptc_output, comm)
Catch exception, log, and throw in a parallel try-catch clause.
Definition: DeferredLoggingErrorHelpers.hpp:202
#define OPM_DEFLOG_THROW(Exception, message, deferred_logger)
Definition: DeferredLoggingErrorHelpers.hpp:45
#define OPM_END_PARALLEL_TRY_CATCH(prefix, comm)
Catch exception and throw in a parallel try-catch clause.
Definition: DeferredLoggingErrorHelpers.hpp:192
#define OPM_PARALLEL_CATCH_CLAUSE(obptc_exc_type, obptc_exc_msg)
Inserts catch classes for the parallel try-catch.
Definition: DeferredLoggingErrorHelpers.hpp:166
#define OPM_BEGIN_PARALLEL_TRY_CATCH()
Macro to setup the try of a parallel try-catch.
Definition: DeferredLoggingErrorHelpers.hpp:158
void logAndCheckForExceptionsAndThrow(Opm::DeferredLogger &deferred_logger, Opm::ExceptionType::ExcEnum exc_type, const std::string &message, const bool terminal_output, Opm::Parallel::Communication comm)
Definition: DeferredLoggingErrorHelpers.hpp:111
Class for handling constraints for the blackoil well model.
Definition: BlackoilWellModelConstraints.hpp:42
Class for handling the gaslift in the blackoil well model.
Definition: BlackoilWellModelGasLift.hpp:96
Class for handling the blackoil well model.
Definition: BlackoilWellModelGeneric.hpp:98
BlackoilWellModelWBP< Scalar, IndexTraits > wbp_
Definition: BlackoilWellModelGeneric.hpp:523
const Parallel::Communication & comm() const
Definition: BlackoilWellModelGeneric.hpp:227
std::vector< ParallelWellInfo< Scalar > > parallel_well_info_
Definition: BlackoilWellModelGeneric.hpp:550
void assignWellTracerRates(data::Wells &wsrpt, const WellTracerRates &wellTracerRates, const unsigned reportStep) const
Class for handling the guide rates in the blackoil well model.
Definition: BlackoilWellModelGuideRates.hpp:47
Class for handling the blackoil well model.
Definition: BlackoilWellModel.hpp:99
void initializeGroupStructure(const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:294
void calcResvCoeff(const int fipnum, const int pvtreg, const std::vector< Scalar > &production_rates, std::vector< Scalar > &resv_coeff) const override
Definition: BlackoilWellModel_impl.hpp:2149
void prepareTimeStep(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1972
std::tuple< bool, bool, Scalar > updateWellControlsAndNetworkIteration(const bool mandatory_network_balance, const bool relax_network_tolerance, const bool optimize_gas_lift, const double dt, DeferredLogger &local_deferredLogger)
Definition: BlackoilWellModel_impl.hpp:1264
WellInterfacePtr createWellPointer(const int wellID, const int report_step) const
Definition: BlackoilWellModel_impl.hpp:1051
void prepareWellsBeforeAssembling(const double dt)
Definition: BlackoilWellModel_impl.hpp:1363
void init()
Definition: BlackoilWellModel_impl.hpp:163
const Simulator & simulator() const
Definition: BlackoilWellModel.hpp:370
std::vector< Scalar > depth_
Definition: BlackoilWellModel.hpp:538
std::size_t global_num_cells_
Definition: BlackoilWellModel.hpp:534
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilWellModel.hpp:108
void initWellContainer(const int reportStepIdx) override
Definition: BlackoilWellModel_impl.hpp:182
void beginReportStep(const int time_step)
Definition: BlackoilWellModel_impl.hpp:199
const WellInterface< TypeTag > & getWell(const std::string &well_name) const
Definition: BlackoilWellModel_impl.hpp:2121
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: BlackoilWellModel.hpp:104
Dune::FieldVector< Scalar, numEq > VectorBlockType
Definition: BlackoilWellModel.hpp:132
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: BlackoilWellModel.hpp:105
GetPropType< TypeTag, Properties::Grid > Grid
Definition: BlackoilWellModel.hpp:102
int numConservationQuantities() const
Definition: BlackoilWellModel_impl.hpp:2093
bool updateWellControls(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1647
int reportStepIndex() const
Definition: BlackoilWellModel_impl.hpp:2137
void calculateProductivityIndexValues(DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1914
void extractLegacyDepth_()
Definition: BlackoilWellModel_impl.hpp:2109
void extractLegacyCellPvtRegionIndex_()
Definition: BlackoilWellModel_impl.hpp:2077
void recoverWellSolutionAndUpdateWellStateDomain(const BVector &x, const int domainIdx)
Definition: BlackoilWellModel_impl.hpp:1568
void updateAverageFormationFactor()
Definition: BlackoilWellModel_impl.hpp:2020
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilWellModel.hpp:107
BlackoilWellModel(Simulator &simulator, const NewtonIterationContext &iter_ctx)
Definition: BlackoilWellModel_impl.hpp:72
void initializeWellState(const int timeStepIdx)
Definition: BlackoilWellModel_impl.hpp:761
const Grid & grid() const
Definition: BlackoilWellModel.hpp:367
void updatePrimaryVariables()
Definition: BlackoilWellModel_impl.hpp:2068
void computeWellTemperature()
Definition: BlackoilWellModel_impl.hpp:2171
void addWellPressureEquations(PressureMatrix &jacobian, const BVector &weights, const bool use_well_weights) const
Definition: BlackoilWellModel_impl.hpp:1472
const SimulatorReportSingle & lastReport() const
Definition: BlackoilWellModel_impl.hpp:637
bool updateWellControlsAndNetwork(const bool mandatory_network_balance, const double dt, DeferredLogger &local_deferredLogger)
Definition: BlackoilWellModel_impl.hpp:1202
void addWellContributions(SparseMatrixAdapter &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1462
void assembleWellEq(const double dt)
Definition: BlackoilWellModel_impl.hpp:1351
WellInterfacePtr createWellForWellTest(const std::string &well_name, const int report_step, DeferredLogger &deferred_logger) const
Definition: BlackoilWellModel_impl.hpp:1102
void calculateExplicitQuantities() const
Definition: BlackoilWellModel_impl.hpp:1632
void updateAndCommunicate(const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:1727
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: BlackoilWellModel.hpp:289
void computeTotalRatesForDof(RateVector &rate, unsigned globalIdx) const
Definition: BlackoilWellModel_impl.hpp:724
void beginTimeStep()
Definition: BlackoilWellModel_impl.hpp:326
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: BlackoilWellModel.hpp:109
bool updateGroupControls(const Group &group, DeferredLogger &deferred_logger, const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:1754
void calcInjResvCoeff(const int fipnum, const int pvtreg, std::vector< Scalar > &resv_coeff) const override
Definition: BlackoilWellModel_impl.hpp:2160
void initializeLocalWellStructure(const int reportStepIdx, const bool enableWellPIScaling)
Definition: BlackoilWellModel_impl.hpp:248
Dune::BlockVector< VectorBlockType > BVector
Definition: BlackoilWellModel.hpp:133
void wellTesting(const int timeStepIdx, const double simulationTime, DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:560
ConvergenceReport getWellConvergence(const std::vector< Scalar > &B_avg, const bool checkWellGroupControlsAndNetwork=false) const
Definition: BlackoilWellModel_impl.hpp:1581
typename FluidSystem::IndexTraitsType IndexTraits
Definition: BlackoilWellModel.hpp:115
void updateCellRatesForDomain(int domainIndex, const std::map< std::string, int > &well_domain_map)
Definition: BlackoilWellModel_impl.hpp:1409
void assembleWellEqWithoutIteration(const double dt)
Definition: BlackoilWellModel_impl.hpp:1377
void updateCellRates()
Definition: BlackoilWellModel_impl.hpp:1397
void assemble(const double dt)
Definition: BlackoilWellModel_impl.hpp:1127
std::size_t local_num_cells_
Definition: BlackoilWellModel.hpp:536
bool alternative_well_rate_init_
Definition: BlackoilWellModel.hpp:539
void timeStepSucceeded(const double simulationTime, const double dt)
Definition: BlackoilWellModel_impl.hpp:647
std::unique_ptr< WellType > createTypedWellPointer(const int wellID, const int time_step) const
Definition: BlackoilWellModel_impl.hpp:1071
void computePotentials(const std::size_t widx, const WellState< Scalar, IndexTraits > &well_state_copy, std::string &exc_msg, ExceptionType::ExcEnum &exc_type) override
Definition: BlackoilWellModel_impl.hpp:1879
Simulator & simulator_
Definition: BlackoilWellModel.hpp:508
void createWellContainer(const int report_step) override
Definition: BlackoilWellModel_impl.hpp:804
std::unique_ptr< WellInterface< TypeTag > > WellInterfacePtr
Definition: BlackoilWellModel.hpp:189
void updateWellTestState(const double simulationTime, WellTestState &wellTestState)
upate the wellTestState related to economic limits
Definition: BlackoilWellModel_impl.hpp:1808
void addReservoirSourceTerms(GlobalEqVector &residual, const std::vector< typename SparseMatrixAdapter::MatrixBlock * > &diagMatAddress) const
Definition: BlackoilWellModel_impl.hpp:1494
int compressedIndexForInterior(int cartesian_cell_idx) const override
get compressed index for interior cells (-1, otherwise
Definition: BlackoilWellModel.hpp:342
void recoverWellSolutionAndUpdateWellState(const BVector &x)
Definition: BlackoilWellModel_impl.hpp:1544
void addWellPressureEquationsStruct(PressureMatrix &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1524
void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1928
void endReportStep()
Definition: BlackoilWellModel_impl.hpp:620
Definition: ConvergenceReport.hpp:38
void setWellFailed(const WellFailure &wf)
Definition: ConvergenceReport.hpp:272
void setWellGroupTargetsViolated(const bool wellGroupTargetsViolated)
Definition: ConvergenceReport.hpp:290
const std::vector< WellFailure > & wellFailures() const
Definition: ConvergenceReport.hpp:380
void setNetworkNotYetBalancedForceAnotherNewtonIteration(const bool network_needs_more_balancing_force_another_newton_iteration)
Definition: ConvergenceReport.hpp:295
Definition: DeferredLogger.hpp:57
void info(const std::string &tag, const std::string &message)
void warning(const std::string &tag, const std::string &message)
void debug(const std::string &tag, const std::string &message)
std::map< std::string, std::pair< const Well *, int > > GLiftEclWells
Definition: GasLiftGroupInfo.hpp:65
Definition: StandardWell.hpp:55
virtual void init(const std::vector< Scalar > &depth_arg, const Scalar gravity_arg, const std::vector< Scalar > &B_avg, const bool changed_to_open_this_step) override
Definition: StandardWell_impl.hpp:77
Definition: WellContributions.hpp:51
void alloc()
Allocate memory for the StandardWells.
void setBlockSize(unsigned int dim, unsigned int dim_wells)
void addNumBlocks(unsigned int numBlocks)
int indexOfWell() const
Index of well in the wells struct and wellState.
Definition: WellInterface.hpp:78
virtual void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const =0
bool updateWellControl(const Simulator &simulator, const IndividualOrGroup iog, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:190
Definition: WellState.hpp:66
ExcEnum
Definition: DeferredLogger.hpp:45
@ NONE
Definition: DeferredLogger.hpp:46
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilbioeffectsmodules.hh:45
ConvergenceReport gatherConvergenceReport(const ConvergenceReport &local_report, Parallel::Communication communicator)
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
#define RESERVOIR_COUPLING_ENABLED
Definition: simulator.hh:32
Context for iteration-dependent decisions in the Newton solver.
Definition: NewtonIterationContext.hpp:43
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34