BlackoilWellModel_impl.hpp
Go to the documentation of this file.
1/*
2 Copyright 2016 - 2019 SINTEF Digital, Mathematics & Cybernetics.
3 Copyright 2016 - 2018 Equinor ASA.
4 Copyright 2017 Dr. Blatt - HPC-Simulation-Software & Services
5 Copyright 2016 - 2018 Norce AS
6
7 This file is part of the Open Porous Media project (OPM).
8
9 OPM is free software: you can redistribute it and/or modify
10 it under the terms of the GNU General Public License as published by
11 the Free Software Foundation, either version 3 of the License, or
12 (at your option) any later version.
13
14 OPM is distributed in the hope that it will be useful,
15 but WITHOUT ANY WARRANTY; without even the implied warranty of
16 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 GNU General Public License for more details.
18
19 You should have received a copy of the GNU General Public License
20 along with OPM. If not, see <http://www.gnu.org/licenses/>.
21*/
22
23#ifndef OPM_BLACKOILWELLMODEL_IMPL_HEADER_INCLUDED
24#define OPM_BLACKOILWELLMODEL_IMPL_HEADER_INCLUDED
25
26// Improve IDE experience
27#ifndef OPM_BLACKOILWELLMODEL_HEADER_INCLUDED
28#include <config.h>
30#endif
31
32#include <opm/grid/utility/cartesianToCompressed.hpp>
33
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>
40
41#include <opm/input/eclipse/Units/UnitSystem.hpp>
42
50
51#ifdef RESERVOIR_COUPLING_ENABLED
56#endif
57
59#if HAVE_MPI
61#endif
62
63#if COMPILE_GPU_BRIDGE
65#endif
66
67#include <algorithm>
68#include <cassert>
69#include <cstddef>
70#include <iomanip>
71#include <optional>
72#include <utility>
73
74#include <fmt/format.h>
75
76namespace Opm {
77 template<typename TypeTag>
80 : WellConnectionModule(*this, simulator.gridView().comm())
81 , BlackoilWellModelGeneric<Scalar, IndexTraits>(simulator.vanguard().schedule(),
82 gaslift_,
83 network_,
84 simulator.vanguard().summaryState(),
85 simulator.vanguard().eclState(),
86 FluidSystem::phaseUsage(),
87 simulator.gridView().comm())
88 , simulator_(simulator)
89 , guide_rate_handler_{
90 *this,
91 simulator.vanguard().schedule(),
92 simulator.vanguard().summaryState(),
93 simulator.vanguard().grid().comm()
94 }
95 , gaslift_(this->terminal_output_)
96 , network_(*this)
97 {
98 local_num_cells_ = simulator_.gridView().size(0);
99
100 // Number of cells the global grid view
101 global_num_cells_ = simulator_.vanguard().globalNumCells();
102
103 {
104 auto& parallel_wells = simulator.vanguard().parallelWells();
105
106 this->parallel_well_info_.reserve(parallel_wells.size());
107 for( const auto& name_bool : parallel_wells) {
108 this->parallel_well_info_.emplace_back(name_bool, grid().comm());
109 }
110 }
111
113 Parameters::Get<Parameters::AlternativeWellRateInit>();
114
115 using SourceDataSpan =
116 typename PAvgDynamicSourceData<Scalar>::template SourceDataSpan<Scalar>;
117
119 [this](const std::size_t globalIndex)
120 { return this->compressedIndexForInterior(globalIndex); },
121 [this](const int localCell, SourceDataSpan sourceTerms)
122 {
123 using Item = typename SourceDataSpan::Item;
124
125 const auto* intQuants = this->simulator_.model()
126 .cachedIntensiveQuantities(localCell, /*timeIndex = */0);
127 const auto& fs = intQuants->fluidState();
128
129 sourceTerms
130 .set(Item::PoreVol, intQuants->porosity().value() *
131 this->simulator_.model().dofTotalVolume(localCell))
132 .set(Item::Depth, this->depth_[localCell]);
133
134 constexpr auto io = FluidSystem::oilPhaseIdx;
135 constexpr auto ig = FluidSystem::gasPhaseIdx;
136 constexpr auto iw = FluidSystem::waterPhaseIdx;
137
138 // Ideally, these would be 'constexpr'.
139 const auto haveOil = FluidSystem::phaseIsActive(io);
140 const auto haveGas = FluidSystem::phaseIsActive(ig);
141 const auto haveWat = FluidSystem::phaseIsActive(iw);
142
143 auto weightedPhaseDensity = [&fs](const auto ip)
144 {
145 return fs.saturation(ip).value() * fs.density(ip).value();
146 };
147
148 if (haveOil) { sourceTerms.set(Item::Pressure, fs.pressure(io).value()); }
149 else if (haveGas) { sourceTerms.set(Item::Pressure, fs.pressure(ig).value()); }
150 else { sourceTerms.set(Item::Pressure, fs.pressure(iw).value()); }
151
152 // Strictly speaking, assumes SUM(s[p]) == 1.
153 auto rho = 0.0;
154 if (haveOil) { rho += weightedPhaseDensity(io); }
155 if (haveGas) { rho += weightedPhaseDensity(ig); }
156 if (haveWat) { rho += weightedPhaseDensity(iw); }
157
158 sourceTerms.set(Item::MixtureDensity, rho);
159 }
160 );
161 }
162
163 template<typename TypeTag>
164 void
166 init()
167 {
168 extractLegacyCellPvtRegionIndex_();
169 extractLegacyDepth_();
170
171 gravity_ = simulator_.problem().gravity()[2];
172
173 this->initial_step_ = true;
174
175 // add the eWoms auxiliary module for the wells to the list
176 simulator_.model().addAuxiliaryModule(this);
177
178 is_cell_perforated_.resize(local_num_cells_, false);
179 }
180
181
182 template<typename TypeTag>
183 void
185 initWellContainer(const int reportStepIdx)
186 {
187 const uint64_t effective_events_mask = ScheduleEvents::WELL_STATUS_CHANGE
188 + ScheduleEvents::NEW_WELL;
189 const auto& events = this->schedule()[reportStepIdx].wellgroup_events();
190 for (auto& wellPtr : this->well_container_) {
191 const bool well_opened_this_step = this->report_step_starts_ &&
192 events.hasEvent(wellPtr->name(),
193 effective_events_mask);
194 wellPtr->init(this->depth_, this->gravity_,
195 this->B_avg_, well_opened_this_step);
196 }
197 }
198
199 template<typename TypeTag>
200 void
202 beginReportStep(const int timeStepIdx)
203 {
204 this->groupStateHelper().setReportStep(timeStepIdx);
205 this->report_step_starts_ = true;
206 this->report_step_start_events_ = this->schedule()[timeStepIdx].wellgroup_events();
207
208 this->rateConverter_ = std::make_unique<RateConverterType>
209 (std::vector<int>(this->local_num_cells_, 0));
210
211 {
212 // WELPI scaling runs at start of report step.
213 const auto enableWellPIScaling = true;
214 this->initializeLocalWellStructure(timeStepIdx, enableWellPIScaling);
215 }
216
217 this->initializeGroupStructure(timeStepIdx);
218
219 const auto& comm = this->simulator_.vanguard().grid().comm();
220
222 {
223 // Create facility for calculating reservoir voidage volumes for
224 // purpose of RESV controls.
225 this->rateConverter_->template defineState<ElementContext>(this->simulator_);
226
227 // Update VFP properties.
228 {
229 const auto& sched_state = this->schedule()[timeStepIdx];
230
231 this->vfp_properties_ = std::make_unique<VFPProperties<Scalar, IndexTraits>>
232 (sched_state.vfpinj(), sched_state.vfpprod(), this->wellState());
233 }
234 }
235 OPM_END_PARALLEL_TRY_CATCH("beginReportStep() failed: ", comm)
236
237 // Store the current well and group states in order to recover in
238 // the case of failed iterations
239 this->commitWGState();
240
241 this->wellStructureChangedDynamically_ = false;
242 }
243
244
245
246
247
248 template <typename TypeTag>
249 void
251 initializeLocalWellStructure(const int reportStepIdx,
252 const bool enableWellPIScaling)
253 {
254 auto logger_guard = this->groupStateHelper().pushLogger();
255 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
256
257 const auto& comm = this->simulator_.vanguard().grid().comm();
258
259 // Wells_ecl_ holds this rank's wells, both open and stopped/shut.
260 this->wells_ecl_ = this->getLocalWells(reportStepIdx);
261 this->local_parallel_well_info_ =
262 this->createLocalParallelWellInfo(this->wells_ecl_);
263
264 // At least initializeWellState() might be throw an exception in
265 // UniformTabulated2DFunction. Playing it safe by extending the
266 // scope a bit.
268 {
269 this->initializeWellPerfData();
270 this->initializeWellState(reportStepIdx);
271 this->wbp_.initializeWBPCalculationService();
272
273 if (this->param_.use_multisegment_well_ && this->anyMSWellOpenLocal()) {
274 this->wellState().initWellStateMSWell(this->wells_ecl_, &this->prevWellState());
275 }
276
277 this->initializeWellProdIndCalculators();
278
279 if (enableWellPIScaling && this->schedule()[reportStepIdx].events()
280 .hasEvent(ScheduleEvents::Events::WELL_PRODUCTIVITY_INDEX))
281 {
282 this->runWellPIScaling(reportStepIdx, local_deferredLogger);
283 }
284 }
285 OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger,
286 "Failed to initialize local well structure: ",
287 this->terminal_output_, comm)
288 }
289
290
291
292
293
294 template <typename TypeTag>
295 void
297 initializeGroupStructure(const int reportStepIdx)
298 {
299 const auto& comm = this->simulator_.vanguard().grid().comm();
300
302 {
303 const auto& fieldGroup =
304 this->schedule().getGroup("FIELD", reportStepIdx);
305
306 this->groupStateHelper().setCmodeGroup(fieldGroup);
307
308 // Define per region average pressure calculators for use by
309 // pressure maintenance groups (GPMAINT keyword).
310 if (this->schedule()[reportStepIdx].has_gpmaint()) {
311 this->groupStateHelper().setRegionAveragePressureCalculator(
312 fieldGroup,
313 this->eclState_.fieldProps(),
314 this->regionalAveragePressureCalculator_
315 );
316 }
317 }
318 OPM_END_PARALLEL_TRY_CATCH("Failed to initialize group structure: ", comm)
319 }
320
321
322
323
324
325 // called at the beginning of a time step
326 template<typename TypeTag>
327 void
330 {
331 OPM_TIMEBLOCK(beginTimeStep);
332
333 this->updateAverageFormationFactor();
334
335 auto logger_guard = this->groupStateHelper().pushLogger();
336 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
337
338#ifdef RESERVOIR_COUPLING_ENABLED
339 auto rescoup_logger_guard = this->setupRescoupScopedLogger(local_deferredLogger);
340#endif
341
342 this->switched_prod_groups_.clear();
343 this->switched_inj_groups_.clear();
344
345 if (this->wellStructureChangedDynamically_) {
346 // Something altered the well structure/topology. Possibly
347 // WELSPECS/COMPDAT and/or WELOPEN run from an ACTIONX block.
348 // Reconstruct the local wells to account for the new well
349 // structure.
350 const auto reportStepIdx =
351 this->simulator_.episodeIndex();
352
353 // Disable WELPI scaling when well structure is updated in the
354 // middle of a report step.
355 const auto enableWellPIScaling = false;
356
357 this->initializeLocalWellStructure(reportStepIdx, enableWellPIScaling);
358 this->initializeGroupStructure(reportStepIdx);
359
360 this->commitWGState();
361
362 // Reset topology flag to signal that we've handled this
363 // structure change. That way we don't end up here in
364 // subsequent calls to beginTimeStep() unless there's a new
365 // dynamic change to the well structure during a report step.
366 this->wellStructureChangedDynamically_ = false;
367 }
368
369 this->resetWGState();
370 const int reportStepIdx = simulator_.episodeIndex();
371
372 this->wellState().updateWellsDefaultALQ(this->schedule(), reportStepIdx, this->summaryState());
373 this->wellState().gliftTimeStepInit();
374
375 const double simulationTime = simulator_.time();
376 const auto& iterCtx = simulator_.problem().iterationContext();
378 {
379 // test wells
380 wellTesting(reportStepIdx, simulationTime, local_deferredLogger);
381
382 // create the well container
383 createWellContainer(reportStepIdx);
384
385#ifdef RESERVOIR_COUPLING_ENABLED
386 if (this->isReservoirCouplingMaster()) {
387 if (this->reservoirCouplingMaster().isFirstSubstepOfSyncTimestep()) {
388 this->receiveSlaveGroupData();
389 }
390 }
391#endif
392
393 // we need to update the group data after the well is created
394 // to make sure we get the correct mapping.
395 this->updateAndCommunicateGroupData(reportStepIdx,
396 iterCtx,
397 param_.nupcol_group_rate_tolerance_, /*update_wellgrouptarget*/ false);
398
399 // Wells are active if they are active wells on at least one process.
400 const Grid& grid = simulator_.vanguard().grid();
401 this->wells_active_ = grid.comm().max(!this->well_container_.empty());
402
403 // do the initialization for all the wells
404 // TODO: to see whether we can postpone of the intialization of the well containers to
405 // optimize the usage of the following several member variables
406 this->initWellContainer(reportStepIdx);
407
408 // update the updated cell flag
409 std::fill(is_cell_perforated_.begin(), is_cell_perforated_.end(), false);
410 for (auto& well : well_container_) {
411 well->updatePerforatedCell(is_cell_perforated_);
412 }
413
414 // calculate the efficiency factors for each well
415 this->calculateEfficiencyFactors(reportStepIdx);
416
417 if constexpr (has_polymer_)
418 {
419 if (PolymerModule::hasPlyshlog() || getPropValue<TypeTag, Properties::EnablePolymerMW>() ) {
420 this->setRepRadiusPerfLength();
421 }
422 }
423
424 }
425
426 OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger, "beginTimeStep() failed: ",
427 this->terminal_output_, simulator_.vanguard().grid().comm());
428
429 for (auto& well : well_container_) {
430 well->setVFPProperties(this->vfp_properties_.get());
431 well->setGuideRate(&this->guideRate_);
432 }
433
434 this->updateFiltrationModelsPreStep(local_deferredLogger);
435
436 // Close completions due to economic reasons
437 for (auto& well : well_container_) {
438 well->closeCompletions(this->wellTestState());
439 }
440
441 // we need the inj_multiplier from the previous time step
442 this->initInjMult();
443
444 if (alternative_well_rate_init_) {
445 // Update the well rates of well_state_, if only single-phase rates, to
446 // have proper multi-phase rates proportional to rates at bhp zero.
447 // This is done only for producers, as injectors will only have a single
448 // nonzero phase anyway.
449 for (const auto& well : well_container_) {
450 if (well->isProducer() && !well->wellIsStopped()) {
451 well->initializeProducerWellState(simulator_, this->wellState(), local_deferredLogger);
452 }
453 }
454 }
455
456 for (const auto& well : well_container_) {
457 if (well->isVFPActive(local_deferredLogger)){
458 well->setPrevSurfaceRates(this->wellState(), this->prevWellState());
459 }
460 }
461 try {
462 this->updateWellPotentials(reportStepIdx,
463 /*onlyAfterEvent*/true,
464 simulator_.vanguard().summaryConfig(),
465 local_deferredLogger);
466 } catch ( std::runtime_error& e ) {
467 const std::string msg = "A zero well potential is returned for output purposes. ";
468 local_deferredLogger.warning("WELL_POTENTIAL_CALCULATION_FAILED", msg);
469 }
470 //update guide rates
471 this->guide_rate_handler_.updateGuideRates(
472 reportStepIdx, simulationTime, this->wellState(), this->groupState()
473 );
474 bool slave_needs_well_solution = false;
475#ifdef RESERVOIR_COUPLING_ENABLED
476 if (this->isReservoirCouplingSlave()) {
477 if (this->reservoirCouplingSlave().isFirstSubstepOfSyncTimestep()) {
478 this->sendSlaveGroupDataToMaster();
479 this->receiveGroupConstraintsFromMaster();
480 this->groupStateHelper().updateSlaveGroupCmodesFromMaster();
481 slave_needs_well_solution = true;
482 }
483 }
484#endif
485 std::string exc_msg;
486 auto exc_type = ExceptionType::NONE;
487 // update gpmaint targets
488 if (this->schedule_[reportStepIdx].has_gpmaint()) {
489 for (const auto& calculator : regionalAveragePressureCalculator_) {
490 calculator.second->template defineState<ElementContext>(simulator_);
491 }
492 const double dt = simulator_.timeStepSize();
493 const Group& fieldGroup = this->schedule().getGroup("FIELD", reportStepIdx);
494 try {
495 this->groupStateHelper().updateGpMaintTargetForGroups(fieldGroup,
496 regionalAveragePressureCalculator_,
497 dt);
498 }
499 OPM_PARALLEL_CATCH_CLAUSE(exc_type, exc_msg);
500 }
501
502 this->updateAndCommunicateGroupData(reportStepIdx,
503 iterCtx,
504 param_.nupcol_group_rate_tolerance_,
505 /*update_wellgrouptarget*/ true);
506 try {
507 // Compute initial well solution for new wells and injectors that change injection type i.e. WAG.
508 for (auto& well : well_container_) {
509 const uint64_t effective_events_mask = ScheduleEvents::WELL_STATUS_CHANGE
510 + ScheduleEvents::INJECTION_TYPE_CHANGED
511 + ScheduleEvents::WELL_SWITCHED_INJECTOR_PRODUCER
512 + ScheduleEvents::NEW_WELL;
513
514 const auto& events = this->schedule()[reportStepIdx].wellgroup_events();
515 const bool event = this->report_step_starts_ && events.hasEvent(well->name(), effective_events_mask);
516 const bool dyn_status_change = this->wellState().well(well->name()).status
517 != this->prevWellState().well(well->name()).status;
518
519 if (event || dyn_status_change || slave_needs_well_solution) {
520 try {
521 well->scaleSegmentRatesAndPressure(this->wellState());
522 well->calculateExplicitQuantities(simulator_, this->groupStateHelper());
523 well->updateWellStateWithTarget(simulator_, this->groupStateHelper(), this->wellState());
524 well->updatePrimaryVariables(this->groupStateHelper());
525 well->solveWellEquation(
526 simulator_, this->groupStateHelper(), this->wellState()
527 );
528 } catch (const std::exception& e) {
529 const std::string msg = "Compute initial well solution for new well " + well->name() + " failed. Continue with zero initial rates";
530 local_deferredLogger.warning("WELL_INITIAL_SOLVE_FAILED", msg);
531 }
532 }
533 }
534 }
535 // Catch clauses for all errors setting exc_type and exc_msg
536 OPM_PARALLEL_CATCH_CLAUSE(exc_type, exc_msg);
537
538#ifdef RESERVOIR_COUPLING_ENABLED
539 if (this->isReservoirCouplingSlave()) {
540 if (slave_needs_well_solution) {
541 this->updateAndCommunicateGroupData(reportStepIdx,
542 iterCtx,
543 param_.nupcol_group_rate_tolerance_,
544 /*update_wellgrouptarget*/ false);
545 this->sendSlaveGroupDataToMaster();
546 }
547 }
548#endif
549
550#ifdef RESERVOIR_COUPLING_ENABLED
551 if (this->isReservoirCouplingMaster()) {
552 if (this->reservoirCouplingMaster().isFirstSubstepOfSyncTimestep()) {
553 this->sendMasterGroupConstraintsToSlaves();
554 this->receiveSlaveGroupData();
555 }
556 }
557#endif
558
559 if (exc_type != ExceptionType::NONE) {
560 const std::string msg = "Compute initial well solution for new wells failed. Continue with zero initial rates";
561 local_deferredLogger.warning("WELL_INITIAL_SOLVE_FAILED", msg);
562 }
563
564 const auto& comm = simulator_.vanguard().grid().comm();
565 logAndCheckForExceptionsAndThrow(local_deferredLogger,
566 exc_type, "beginTimeStep() failed: " + exc_msg, this->terminal_output_, comm);
567
568 }
569
570#ifdef RESERVOIR_COUPLING_ENABLED
571 // Automatically manages the lifecycle of the DeferredLogger pointer
572 // in the reservoir coupling logger. Ensures the logger is properly
573 // cleared when it goes out of scope, preventing dangling pointer issues:
574 //
575 // - The ScopedLoggerGuard constructor sets the logger pointer
576 // - When the guard goes out of scope, the destructor clears the pointer
577 // - Move semantics transfer ownership safely when returning from this function
578 // - The moved-from guard is "nullified" and its destructor does nothing
579 // - Only the final guard in the caller will clear the logger
580 template<typename TypeTag>
581 std::optional<ReservoirCoupling::ScopedLoggerGuard>
584 if (this->isReservoirCouplingMaster()) {
586 this->reservoirCouplingMaster().logger(),
587 &local_logger
588 };
589 } else if (this->isReservoirCouplingSlave()) {
590 return ReservoirCoupling::ScopedLoggerGuard{
591 this->reservoirCouplingSlave().logger(),
592 &local_logger
593 };
594 }
595 return std::nullopt;
596 }
597
598 template<typename TypeTag>
599 void
600 BlackoilWellModel<TypeTag>::
601 receiveSlaveGroupData()
602 {
603 assert(this->isReservoirCouplingMaster());
604 RescoupReceiveSlaveGroupData<Scalar, IndexTraits> slave_group_data_receiver{
605 this->groupStateHelper(),
606 };
607 slave_group_data_receiver.receiveSlaveGroupData();
608 }
609
610 template<typename TypeTag>
611 void
612 BlackoilWellModel<TypeTag>::
613 sendSlaveGroupDataToMaster()
614 {
615 assert(this->isReservoirCouplingSlave());
616 RescoupSendSlaveGroupData<Scalar, IndexTraits> slave_group_data_sender{this->groupStateHelper()};
617 slave_group_data_sender.sendSlaveGroupDataToMaster();
618 }
619
620 template<typename TypeTag>
621 void
622 BlackoilWellModel<TypeTag>::
623 sendMasterGroupConstraintsToSlaves()
624 {
625 // This function is called by the master process to send the group constraints to the slaves.
626 RescoupConstraintsCalculator<Scalar, IndexTraits> constraints_calculator{
627 this->guide_rate_handler_,
628 this->groupStateHelper()
629 };
630 constraints_calculator.calculateMasterGroupConstraintsAndSendToSlaves();
631 }
632
633 template<typename TypeTag>
634 void
635 BlackoilWellModel<TypeTag>::
636 receiveGroupConstraintsFromMaster()
637 {
638 RescoupReceiveGroupConstraints<Scalar, IndexTraits> constraint_receiver{
639 this->guide_rate_handler_,
640 this->groupStateHelper()
641 };
642 constraint_receiver.receiveGroupConstraintsFromMaster();
643 }
644
645#endif // RESERVOIR_COUPLING_ENABLED
646
647 template<typename TypeTag>
648 void
650 const double simulationTime,
651 DeferredLogger& deferred_logger)
652 {
653 for (const std::string& well_name : this->getWellsForTesting(timeStepIdx, simulationTime)) {
654 const Well& wellEcl = this->schedule().getWell(well_name, timeStepIdx);
655 if (wellEcl.getStatus() == Well::Status::SHUT)
656 continue;
657
658 WellInterfacePtr well = createWellForWellTest(well_name, timeStepIdx, deferred_logger);
659 // some preparation before the well can be used
660 well->init(depth_, gravity_, B_avg_, true);
661
662 Scalar well_efficiency_factor = wellEcl.getEfficiencyFactor() *
663 this->wellState().getGlobalEfficiencyScalingFactor(well_name);
664 this->groupStateHelper().accumulateGroupEfficiencyFactor(
665 this->schedule().getGroup(wellEcl.groupName(), timeStepIdx),
666 well_efficiency_factor
667 );
668
669 well->setWellEfficiencyFactor(well_efficiency_factor);
670 well->setVFPProperties(this->vfp_properties_.get());
671 well->setGuideRate(&this->guideRate_);
672
673 // initialize rates/previous rates to prevent zero fractions in vfp-interpolation
674 if (well->isProducer() && alternative_well_rate_init_) {
675 well->initializeProducerWellState(simulator_, this->wellState(), deferred_logger);
676 }
677 if (well->isVFPActive(deferred_logger)) {
678 well->setPrevSurfaceRates(this->wellState(), this->prevWellState());
679 }
680
681 const auto& network = this->schedule()[timeStepIdx].network();
682 if (network.active()) {
683 this->network_.initializeWell(*well);
684 }
685 try {
686 using GLiftEclWells = typename GasLiftGroupInfo<Scalar, IndexTraits>::GLiftEclWells;
687 GLiftEclWells ecl_well_map;
688 gaslift_.initGliftEclWellMap(well_container_, ecl_well_map);
689 well->wellTesting(simulator_,
690 simulationTime,
691 this->groupStateHelper(),
692 this->wellState(),
693 this->wellTestState(),
694 ecl_well_map,
695 this->well_open_times_);
696 } catch (const std::exception& e) {
697 const std::string msg =
698 fmt::format(fmt::runtime("Exception during testing of well: {}. The well will not open.\n"
699 "Exception message: {}"), wellEcl.name(), e.what());
700 deferred_logger.warning("WELL_TESTING_FAILED", msg);
701 }
702 }
703 }
704
705 // called at the end of a report step
706 template<typename TypeTag>
707 void
710 {
711 // Clear the communication data structures for above values.
712 for (auto&& pinfo : this->local_parallel_well_info_)
713 {
714 pinfo.get().clear();
715 }
716 }
717
718
719
720
721
722 // called at the end of a report step
723 template<typename TypeTag>
726 lastReport() const {return last_report_; }
727
728
729
730
731
732 // called at the end of a time step
733 template<typename TypeTag>
734 void
736 timeStepSucceeded(const double simulationTime, const double dt)
737 {
738 this->closed_this_step_.clear();
739
740 // time step is finished and we are not any more at the beginning of an report step
741 this->report_step_starts_ = false;
742 const int reportStepIdx = simulator_.episodeIndex();
743
744 auto logger_guard = this->groupStateHelper().pushLogger();
745 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
746 for (const auto& well : well_container_) {
747 if (getPropValue<TypeTag, Properties::EnablePolymerMW>() && well->isInjector()) {
748 well->updateWaterThroughput(dt, this->wellState());
749 }
750 }
751 // update connection transmissibility factor and d factor (if applicable) in the wellstate
752 for (const auto& well : well_container_) {
753 well->updateConnectionTransmissibilityFactor(simulator_, this->wellState().well(well->indexOfWell()));
754 well->updateConnectionDFactor(simulator_, this->wellState().well(well->indexOfWell()));
755 }
756
757 if (Indices::waterEnabled) {
758 this->updateFiltrationModelsPostStep(dt, FluidSystem::waterPhaseIdx, local_deferredLogger);
759 }
760
761 // WINJMULT: At the end of the time step, update the inj_multiplier saved in WellState for later use
762 this->updateInjMult(local_deferredLogger);
763
764 // report well switching
765 for (const auto& well : well_container_) {
766 well->reportWellSwitching(this->wellState().well(well->indexOfWell()), local_deferredLogger);
767 }
768 // report group switching
769 if (this->terminal_output_) {
770 this->reportGroupSwitching(local_deferredLogger);
771 }
772
773 // update the rate converter with current averages pressures etc in
774 rateConverter_->template defineState<ElementContext>(simulator_);
775
776 // calculate the well potentials
777 try {
778 this->updateWellPotentials(reportStepIdx,
779 /*onlyAfterEvent*/false,
780 simulator_.vanguard().summaryConfig(),
781 local_deferredLogger);
782 } catch ( std::runtime_error& e ) {
783 const std::string msg = "A zero well potential is returned for output purposes. ";
784 local_deferredLogger.warning("WELL_POTENTIAL_CALCULATION_FAILED", msg);
785 }
786
787 updateWellTestState(simulationTime, this->wellTestState());
788
789 // check group sales limits at the end of the timestep
790 const Group& fieldGroup = this->schedule_.getGroup("FIELD", reportStepIdx);
791 this->checkGEconLimits(fieldGroup, simulationTime,
792 simulator_.episodeIndex(), local_deferredLogger);
793 this->checkGconsaleLimits(fieldGroup, this->wellState(),
794 simulator_.episodeIndex(), local_deferredLogger);
795
796 this->calculateProductivityIndexValues(local_deferredLogger);
797
798 this->groupStateHelper().updateNONEProductionGroups();
799
800 this->commitWGState();
801
802 //reporting output temperatures
803 this->computeWellTemperature();
804 }
805
806
807 template<typename TypeTag>
808 void
811 unsigned elemIdx) const
812 {
813 rate = 0;
814
815 if (!is_cell_perforated_[elemIdx] || cellRates_.count(elemIdx) == 0) {
816 return;
817 }
818
819 rate = cellRates_.at(elemIdx);
820 }
821
822
823 template<typename TypeTag>
824 template <class Context>
825 void
828 const Context& context,
829 unsigned spaceIdx,
830 unsigned timeIdx) const
831 {
832 rate = 0;
833 int elemIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
834
835 if (!is_cell_perforated_[elemIdx] || cellRates_.count(elemIdx) == 0) {
836 return;
837 }
838
839 rate = cellRates_.at(elemIdx);
840 }
841
842
843
844 template<typename TypeTag>
845 void
847 initializeWellState(const int timeStepIdx)
848 {
849 const auto pressIx = []()
850 {
851 if (Indices::oilEnabled) { return FluidSystem::oilPhaseIdx; }
852 if (Indices::waterEnabled) { return FluidSystem::waterPhaseIdx; }
853
854 return FluidSystem::gasPhaseIdx;
855 }();
856
857 auto cellPressures = std::vector<Scalar>(this->local_num_cells_, Scalar{0});
858 auto cellTemperatures = std::vector<Scalar>(this->local_num_cells_, Scalar{0});
859
860 auto elemCtx = ElementContext { this->simulator_ };
861 const auto& gridView = this->simulator_.vanguard().gridView();
862
864 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
865 elemCtx.updatePrimaryStencil(elem);
866 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
867
868 const auto ix = elemCtx.globalSpaceIndex(/*spaceIdx=*/0, /*timeIdx=*/0);
869 const auto& fs = elemCtx.intensiveQuantities(/*spaceIdx=*/0, /*timeIdx=*/0).fluidState();
870
871 cellPressures[ix] = fs.pressure(pressIx).value();
872 cellTemperatures[ix] = fs.temperature(0).value();
873 }
874 OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel::initializeWellState() failed: ",
875 this->simulator_.vanguard().grid().comm());
876
877 this->wellState().init(cellPressures, cellTemperatures, this->schedule(), this->wells_ecl_,
878 this->local_parallel_well_info_, timeStepIdx,
879 &this->prevWellState(), this->well_perf_data_,
880 this->summaryState(), simulator_.vanguard().enableDistributedWells());
881 }
882
883
884
885
886
887 template<typename TypeTag>
888 void
890 createWellContainer(const int report_step)
891 {
892 auto logger_guard = this->groupStateHelper().pushLogger();
893 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
894
895 const int nw = this->numLocalWells();
896
897 well_container_.clear();
898
899 if (nw > 0) {
900 well_container_.reserve(nw);
901
902 const auto& wmatcher = this->schedule().wellMatcher(report_step);
903 const auto& wcycle = this->schedule()[report_step].wcycle.get();
904
905 // First loop and check for status changes. This is necessary
906 // as wcycle needs the updated open/close times.
907 std::ranges::for_each(this->wells_ecl_,
908 [this, &wg_events = this->report_step_start_events_](const auto& well_ecl)
909 {
910 if (!well_ecl.hasConnections()) {
911 // No connections in this well. Nothing to do.
912 return;
913 }
914
915 constexpr auto events_mask = ScheduleEvents::WELL_STATUS_CHANGE |
916 ScheduleEvents::REQUEST_OPEN_WELL |
917 ScheduleEvents::REQUEST_SHUT_WELL;
918 const bool well_event =
919 this->report_step_starts_ &&
920 wg_events.hasEvent(well_ecl.name(), events_mask);
921 // WCYCLE is suspendended by explicit SHUT events by the user.
922 // and restarted after explicit OPEN events.
923 // Note: OPEN or SHUT event does not necessary mean the well
924 // actually opened or shut at this point as the simulator could
925 // have done this by operabilty checks and well testing. This
926 // may need further testing and imply code changes to cope with
927 // these corner cases.
928 if (well_event) {
929 if (well_ecl.getStatus() == WellStatus::OPEN) {
930 this->well_open_times_.insert_or_assign(well_ecl.name(),
931 this->simulator_.time());
932 this->well_close_times_.erase(well_ecl.name());
933 } else if (well_ecl.getStatus() == WellStatus::SHUT) {
934 this->well_close_times_.insert_or_assign(well_ecl.name(),
935 this->simulator_.time());
936 this->well_open_times_.erase(well_ecl.name());
937 }
938 }
939 });
940
941 // Grab wcycle states. This needs to run before the schedule gets processed
942 const auto cycle_states = wcycle.wellStatus(this->simulator_.time(),
943 wmatcher,
944 this->well_open_times_,
945 this->well_close_times_);
946
947 for (int w = 0; w < nw; ++w) {
948 const Well& well_ecl = this->wells_ecl_[w];
949
950 if (!well_ecl.hasConnections()) {
951 // No connections in this well. Nothing to do.
952 continue;
953 }
954
955 const std::string& well_name = well_ecl.name();
956 const auto well_status = this->schedule()
957 .getWell(well_name, report_step).getStatus();
958
959 const bool shut_event = this->wellState().well(w).events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE)
960 && well_status == Well::Status::SHUT;
961 const bool open_event = this->wellState().well(w).events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE)
962 && well_status == Well::Status::OPEN;
963 const auto& ws = this->wellState().well(well_name);
964
965 if (shut_event && ws.status != Well::Status::SHUT) {
966 this->closed_this_step_.insert(well_name);
967 this->wellState().shutWell(w);
968 } else if (open_event && ws.status != Well::Status::OPEN) {
969 this->wellState().openWell(w);
970 }
971
972 // A new WCON keywords can re-open a well that was closed/shut due to Physical limit
973 if (this->wellTestState().well_is_closed(well_name)) {
974 // The well was shut this timestep, we are most likely retrying
975 // a timestep without the well in question, after it caused
976 // repeated timestep cuts. It should therefore not be opened,
977 // even if it was new or received new targets this report step.
978 const bool closed_this_step = (this->wellTestState().lastTestTime(well_name) == simulator_.time());
979 // TODO: more checking here, to make sure this standard more specific and complete
980 // maybe there is some WCON keywords will not open the well
981 auto& events = this->wellState().well(w).events;
982 if (events.hasEvent(ScheduleEvents::REQUEST_OPEN_WELL)) {
983 if (!closed_this_step) {
984 this->wellTestState().open_well(well_name);
985 this->wellTestState().open_completions(well_name);
986 this->well_open_times_.insert_or_assign(well_name,
987 this->simulator_.time());
988 this->well_close_times_.erase(well_name);
989 }
990 events.clearEvent(ScheduleEvents::REQUEST_OPEN_WELL);
991 }
992 }
993
994 // TODO: should we do this for all kinds of closing reasons?
995 // something like wellTestState().hasWell(well_name)?
996 if (this->wellTestState().well_is_closed(well_name))
997 {
998 if (well_ecl.getAutomaticShutIn()) {
999 // shut wells are not added to the well container
1000 this->wellState().shutWell(w);
1001 this->well_close_times_.erase(well_name);
1002 this->well_open_times_.erase(well_name);
1003 continue;
1004 } else {
1005 if (!well_ecl.getAllowCrossFlow()) {
1006 // stopped wells where cross flow is not allowed
1007 // are not added to the well container
1008 this->wellState().shutWell(w);
1009 this->well_close_times_.erase(well_name);
1010 this->well_open_times_.erase(well_name);
1011 continue;
1012 }
1013 // stopped wells are added to the container but marked as stopped
1014 this->wellState().stopWell(w);
1015 }
1016 }
1017
1018 // shut wells with zero rante constraints and disallowing
1019 if (!well_ecl.getAllowCrossFlow()) {
1020 const bool any_zero_rate_constraint = well_ecl.isProducer()
1021 ? well_ecl.productionControls(this->summaryState_).anyZeroRateConstraint()
1022 : well_ecl.injectionControls(this->summaryState_).anyZeroRateConstraint();
1023 if (any_zero_rate_constraint) {
1024 // Treat as shut, do not add to container.
1025 local_deferredLogger.debug(fmt::format(fmt::runtime(" Well {} gets shut due to having zero rate constraint and disallowing crossflow "), well_ecl.name()));
1026 this->wellState().shutWell(w);
1027 this->well_close_times_.erase(well_name);
1028 this->well_open_times_.erase(well_name);
1029 continue;
1030 }
1031 }
1032
1033 if (!wcycle.empty()) {
1034 const auto it = cycle_states.find(well_name);
1035 if (it != cycle_states.end()) {
1036 if (!it->second || well_status == Well::Status::SHUT) {
1037 // If well is shut in schedule we keep it shut
1038 if (well_status == Well::Status::SHUT) {
1039 this->well_open_times_.erase(well_name);
1040 this->well_close_times_.erase(well_name);
1041 }
1042 this->wellState().shutWell(w);
1043 continue;
1044 } else {
1045 this->wellState().openWell(w);
1046 }
1047 }
1048 }
1049
1050 // We dont add SHUT wells to the container
1051 if (ws.status == Well::Status::SHUT) {
1052 continue;
1053 }
1054
1055 well_container_.emplace_back(this->createWellPointer(w, report_step));
1056
1057 if (ws.status == Well::Status::STOP) {
1058 well_container_.back()->stopWell();
1059 this->well_close_times_.erase(well_name);
1060 this->well_open_times_.erase(well_name);
1061 }
1062 }
1063
1064 if (!wcycle.empty()) {
1065 const auto schedule_open =
1066 [&wg_events = this->report_step_start_events_](const std::string& name)
1067 {
1068 return wg_events.hasEvent(name, ScheduleEvents::REQUEST_OPEN_WELL);
1069 };
1070 for (const auto& [wname, wscale] : wcycle.efficiencyScale(this->simulator_.time(),
1071 this->simulator_.timeStepSize(),
1072 wmatcher,
1073 this->well_open_times_,
1074 schedule_open))
1075 {
1076 this->wellState().updateEfficiencyScalingFactor(wname, wscale);
1077 this->schedule_.add_event(ScheduleEvents::WELLGROUP_EFFICIENCY_UPDATE, report_step);
1078 }
1079 }
1080 }
1081
1082 this->well_container_generic_.clear();
1083 for (auto& w : well_container_) {
1084 this->well_container_generic_.push_back(w.get());
1085 }
1086
1087 this->network_.initialize(report_step);
1088
1089 this->wbp_.registerOpenWellsForWBPCalculation();
1090 }
1091
1092
1093
1094
1095
1096 template <typename TypeTag>
1099 createWellPointer(const int wellID, const int report_step) const
1100 {
1101 const auto is_multiseg = this->wells_ecl_[wellID].isMultiSegment();
1102
1103 if (! (this->param_.use_multisegment_well_ && is_multiseg)) {
1104 return this->template createTypedWellPointer<StandardWell<TypeTag>>(wellID, report_step);
1105 }
1106 else {
1107 return this->template createTypedWellPointer<MultisegmentWell<TypeTag>>(wellID, report_step);
1108 }
1109 }
1110
1111
1112
1113
1114
1115 template <typename TypeTag>
1116 template <typename WellType>
1117 std::unique_ptr<WellType>
1119 createTypedWellPointer(const int wellID, const int time_step) const
1120 {
1121 // Use the pvtRegionIdx from the top cell
1122 const auto& perf_data = this->well_perf_data_[wellID];
1123
1124 // Cater for case where local part might have no perforations.
1125 const auto pvtreg = perf_data.empty()
1126 ? 0 : this->pvt_region_idx_[perf_data.front().cell_index];
1127
1128 const auto& parallel_well_info = this->local_parallel_well_info_[wellID].get();
1129 const auto global_pvtreg = parallel_well_info.broadcastFirstPerforationValue(pvtreg);
1130
1131 return std::make_unique<WellType>(this->wells_ecl_[wellID],
1132 parallel_well_info,
1133 time_step,
1134 this->param_,
1135 *this->rateConverter_,
1136 global_pvtreg,
1137 this->numConservationQuantities(),
1138 this->numPhases(),
1139 wellID,
1140 perf_data);
1141 }
1142
1143
1144
1145
1146
1147 template<typename TypeTag>
1150 createWellForWellTest(const std::string& well_name,
1151 const int report_step,
1152 DeferredLogger& deferred_logger) const
1153 {
1154 // Finding the location of the well in wells_ecl
1155 const auto it =
1156 std::ranges::find_if(this->wells_ecl_,
1157 [&well_name](const auto& w)
1158 { return well_name == w.name(); });
1159 // It should be able to find in wells_ecl.
1160 if (it == this->wells_ecl_.end()) {
1161 OPM_DEFLOG_THROW(std::logic_error,
1162 fmt::format(fmt::runtime("Could not find well {} in wells_ecl"), well_name),
1163 deferred_logger);
1164 }
1165
1166 const int pos = static_cast<int>(std::distance(this->wells_ecl_.begin(), it));
1167 return this->createWellPointer(pos, report_step);
1168 }
1169
1170
1171
1172 template<typename TypeTag>
1173 void
1175 assemble(const double dt)
1176 {
1177 OPM_TIMEFUNCTION();
1178 auto logger_guard = this->groupStateHelper().pushLogger();
1179 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
1180
1181 const auto& iterCtx = simulator_.problem().iterationContext();
1182
1184 if (gaslift_.terminalOutput()) {
1185 const std::string msg =
1186 fmt::format(fmt::runtime("assemble() : iteration {}"), iterCtx.iteration());
1187 gaslift_.gliftDebug(msg, local_deferredLogger);
1188 }
1189 }
1190 last_report_ = SimulatorReportSingle();
1191 Dune::Timer perfTimer;
1192 perfTimer.start();
1193 this->closed_offending_wells_.clear();
1194
1195 {
1196 const int episodeIdx = simulator_.episodeIndex();
1197 const auto& network = this->schedule()[episodeIdx].network();
1198 if (!this->wellsActive() && !network.active()) {
1199 return;
1200 }
1201 }
1202
1203 // Timestep initialization: should run once at the start of each timestep.
1204 if (iterCtx.needsTimestepInit() && this->wellsActive()) {
1205 OPM_TIMEBLOCK(firstIterationAssemble);
1206 // try-catch is needed here as updateWellControls
1207 // contains global communication and has either to
1208 // be reached by all processes or all need to abort
1209 // before.
1211 {
1212 calculateExplicitQuantities();
1213 prepareTimeStep(local_deferredLogger);
1214 }
1215 OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger,
1216 "assemble() failed during well initialization: ",
1217 this->terminal_output_, grid().comm());
1218 }
1219
1220 const bool well_group_control_changed = updateWellControlsAndNetwork(false, dt, local_deferredLogger);
1221
1222 // even when there is no wells active, the network nodal pressure still need to be updated through updateWellControlsAndNetwork()
1223 // but there is no need to assemble the well equations
1224 if ( ! this->wellsActive() ) {
1225 return;
1226 }
1227
1228 assembleWellEqWithoutIteration(dt);
1229 // Pre-compute cell rates to we don't have to do this for every cell during linearization...
1230 updateCellRates();
1231
1232 // if group or well control changes we don't consider the
1233 // case converged
1234 last_report_.well_group_control_changed = well_group_control_changed;
1235 last_report_.assemble_time_well += perfTimer.stop();
1236 }
1237
1238
1239
1240
1241 template<typename TypeTag>
1242 bool
1244 updateWellControlsAndNetwork(const bool mandatory_network_balance,
1245 const double dt,
1246 DeferredLogger& local_deferredLogger)
1247 {
1248 OPM_TIMEFUNCTION();
1249 // not necessarily that we always need to update once of the network solutions
1250 bool do_network_update = true;
1251 bool well_group_control_changed = false;
1252 Scalar network_imbalance = 0.0;
1253 // after certain number of the iterations, we use relaxed tolerance for the network update
1254 const std::size_t iteration_to_relax = param_.network_max_strict_outer_iterations_;
1255 // after certain number of the iterations, we terminate
1256 const std::size_t max_iteration = param_.network_max_outer_iterations_;
1257 std::size_t network_update_iteration = 0;
1258 network_needs_more_balancing_force_another_newton_iteration_ = false;
1259 while (do_network_update) {
1260 if (network_update_iteration >= max_iteration ) {
1261 // only output to terminal if we at the last newton iterations where we try to balance the network.
1262 const int episodeIdx = simulator_.episodeIndex();
1263 const auto& iterCtx = simulator_.problem().iterationContext();
1264 if (this->network_.willBalanceOnNextIteration(episodeIdx, iterCtx)) {
1265 if (this->terminal_output_) {
1266 const std::string msg = fmt::format("Maximum of {:d} network iterations has been used and we stop the update, \n"
1267 "and try again after the next Newton iteration (imbalance = {:.2e} bar)",
1268 max_iteration, network_imbalance*1.0e-5);
1269 local_deferredLogger.debug(msg);
1270 }
1271 // To avoid stopping the newton iterations too early, before the network is converged,
1272 // we need to report it
1273 network_needs_more_balancing_force_another_newton_iteration_ = true;
1274 } else {
1275 if (this->terminal_output_) {
1276 const std::string msg = fmt::format("Maximum of {:d} network iterations has been used and we stop the update. \n"
1277 "The simulator will continue with unconverged network results (imbalance = {:.2e} bar)",
1278 max_iteration, network_imbalance*1.0e-5);
1279 local_deferredLogger.info(msg);
1280 }
1281 }
1282 break;
1283 }
1284 if (this->terminal_output_ && (network_update_iteration == iteration_to_relax) ) {
1285 local_deferredLogger.debug("We begin using relaxed tolerance for network update now after " + std::to_string(iteration_to_relax) + " iterations ");
1286 }
1287 const bool relax_network_balance = network_update_iteration >= iteration_to_relax;
1288 // Never optimize gas lift in last iteration, to allow network convergence (unless max_iter < 2)
1289 const bool optimize_gas_lift = ( (network_update_iteration + 1) < std::max(max_iteration, static_cast<std::size_t>(2)) );
1290 std::tie(well_group_control_changed, do_network_update, network_imbalance) =
1291 updateWellControlsAndNetworkIteration(mandatory_network_balance, relax_network_balance, optimize_gas_lift, dt,local_deferredLogger);
1292 ++network_update_iteration;
1293 }
1294 return well_group_control_changed;
1295 }
1296
1297
1298
1299
1300 template<typename TypeTag>
1301 std::tuple<bool, bool, typename BlackoilWellModel<TypeTag>::Scalar>
1303 updateWellControlsAndNetworkIteration(const bool mandatory_network_balance,
1304 const bool relax_network_tolerance,
1305 const bool optimize_gas_lift,
1306 const double dt,
1307 DeferredLogger& local_deferredLogger)
1308 {
1309 OPM_TIMEFUNCTION();
1310 const auto& iterCtx = simulator_.problem().iterationContext();
1311 const int reportStepIdx = simulator_.episodeIndex();
1312 this->updateAndCommunicateGroupData(reportStepIdx, iterCtx,
1313 param_.nupcol_group_rate_tolerance_, /*update_wellgrouptarget*/ true);
1314 // We need to call updateWellControls before we update the network as
1315 // network updates are only done on thp controlled wells.
1316 // Note that well controls are allowed to change during updateNetwork
1317 // and in prepareWellsBeforeAssembling during well solves.
1318 bool well_group_control_changed = updateWellControls(local_deferredLogger);
1319 const auto [more_inner_network_update, network_imbalance] =
1320 this->network_.update(mandatory_network_balance,
1321 local_deferredLogger,
1322 relax_network_tolerance);
1323
1324 bool alq_updated = false;
1326 {
1327 if (optimize_gas_lift) {
1328 // we need to update the potentials if the thp limit as been modified by
1329 // the network balancing
1330 const bool updatePotentials = (this->network_.shouldBalance(reportStepIdx, iterCtx) ||
1331 mandatory_network_balance);
1332 alq_updated = gaslift_.maybeDoGasLiftOptimize(simulator_,
1333 well_container_,
1334 this->network_.nodePressures(),
1335 updatePotentials,
1336 this->wellState(),
1337 this->groupState(),
1338 local_deferredLogger);
1339 }
1340 prepareWellsBeforeAssembling(dt);
1341 }
1342 OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger,
1343 "updateWellControlsAndNetworkIteration() failed: ",
1344 this->terminal_output_, grid().comm());
1345
1346 // update guide rates
1347 if (alq_updated || BlackoilWellModelGuideRates(*this).
1348 guideRateUpdateIsNeeded(reportStepIdx)) {
1349 const double simulationTime = simulator_.time();
1350 // NOTE: For reservoir coupling: Slave group potentials are only communicated
1351 // at the start of the time step, see beginTimeStep(). Here, we assume those
1352 // potentials remain unchanged during the time step when updating guide rates below.
1353 this->guide_rate_handler_.updateGuideRates(
1354 reportStepIdx, simulationTime, this->wellState(), this->groupState()
1355 );
1356 }
1357 // we need to re-iterate the network when the well group controls changed or gaslift/alq is changed or
1358 // the inner iterations are did not converge
1359 const bool more_network_update = this->network_.shouldBalance(reportStepIdx, iterCtx) &&
1360 (more_inner_network_update || alq_updated);
1361 return {well_group_control_changed, more_network_update, network_imbalance};
1362 }
1363
1364 template<typename TypeTag>
1365 void
1367 assembleWellEq(const double dt)
1368 {
1369 OPM_TIMEFUNCTION();
1370 for (auto& well : well_container_) {
1371 well->assembleWellEq(simulator_, dt, this->groupStateHelper(), this->wellState());
1372 }
1373 }
1374
1375
1376 template<typename TypeTag>
1377 void
1379 prepareWellsBeforeAssembling(const double dt)
1380 {
1381 OPM_TIMEFUNCTION();
1382 for (auto& well : well_container_) {
1383 well->prepareWellBeforeAssembling(
1384 simulator_, dt, this->groupStateHelper(), this->wellState()
1385 );
1386 }
1387 }
1388
1389
1390 template<typename TypeTag>
1391 void
1393 assembleWellEqWithoutIteration(const double dt)
1394 {
1395 OPM_TIMEFUNCTION();
1396 auto& deferred_logger = this->groupStateHelper().deferredLogger();
1397 // We make sure that all processes throw in case there is an exception
1398 // on one of them (WetGasPvt::saturationPressure might throw if not converged)
1400
1401 for (auto& well: well_container_) {
1402 well->assembleWellEqWithoutIteration(simulator_, this->groupStateHelper(), dt, this->wellState(),
1403 /*solving_with_zero_rate=*/false);
1404 }
1405 OPM_END_PARALLEL_TRY_CATCH_LOG(deferred_logger, "BlackoilWellModel::assembleWellEqWithoutIteration failed: ",
1406 this->terminal_output_, grid().comm());
1407
1408 }
1409
1410 template<typename TypeTag>
1411 void
1414 {
1415 // Pre-compute cell rates for all wells
1416 cellRates_.clear();
1417 for (const auto& well : well_container_) {
1418 well->addCellRates(cellRates_);
1419 }
1420 }
1421
1422 template<typename TypeTag>
1423 void
1425 updateCellRatesForDomain(int domainIndex, const std::map<std::string, int>& well_domain_map)
1426 {
1427 // Pre-compute cell rates only for wells in the specified domain
1428 cellRates_.clear();
1429 for (const auto& well : well_container_) {
1430 const auto it = well_domain_map.find(well->name());
1431 if (it != well_domain_map.end() && it->second == domainIndex) {
1432 well->addCellRates(cellRates_);
1433 }
1434 }
1435 }
1436
1437#if COMPILE_GPU_BRIDGE
1438 template<typename TypeTag>
1439 void
1442 {
1443 // prepare for StandardWells
1445
1446 for(unsigned int i = 0; i < well_container_.size(); i++){
1447 auto& well = well_container_[i];
1448 auto derived = dynamic_cast<StandardWell<TypeTag>*>(well.get());
1449 if (derived) {
1450 wellContribs.addNumBlocks(derived->linSys().getNumBlocks());
1451 }
1452 }
1453
1454 // allocate memory for data from StandardWells
1455 wellContribs.alloc();
1456
1457 for(unsigned int i = 0; i < well_container_.size(); i++){
1458 auto& well = well_container_[i];
1459 // maybe WellInterface could implement addWellContribution()
1460 auto derived_std = dynamic_cast<StandardWell<TypeTag>*>(well.get());
1461 if (derived_std) {
1462 derived_std->linSys().extract(derived_std->numStaticWellEq, wellContribs);
1463 } else {
1464 auto derived_ms = dynamic_cast<MultisegmentWell<TypeTag>*>(well.get());
1465 if (derived_ms) {
1466 derived_ms->linSys().extract(wellContribs);
1467 } else {
1468 OpmLog::warning("Warning unknown type of well");
1469 }
1470 }
1471 }
1472 }
1473#endif
1474
1475 template<typename TypeTag>
1476 void
1478 addWellContributions(SparseMatrixAdapter& jacobian) const
1479 {
1480 for ( const auto& well: well_container_ ) {
1481 well->addWellContributions(jacobian);
1482 }
1483 }
1484
1485 template<typename TypeTag>
1486 void
1489 const BVector& weights,
1490 const bool use_well_weights) const
1491 {
1492 int nw = this->numLocalWellsEnd();
1493 int rdofs = local_num_cells_;
1494 for ( int i = 0; i < nw; i++ ) {
1495 int wdof = rdofs + i;
1496 jacobian[wdof][wdof] = 1.0;// better scaling ?
1497 }
1498
1499 for (const auto& well : well_container_) {
1500 well->addWellPressureEquations(jacobian,
1501 weights,
1502 pressureVarIndex,
1503 use_well_weights,
1504 this->wellState());
1505 }
1506 }
1507
1508 template <typename TypeTag>
1510 addReservoirSourceTerms(GlobalEqVector& residual,
1511 const std::vector<typename SparseMatrixAdapter::MatrixBlock*>& diagMatAddress) const
1512 {
1513 // NB this loop may write multiple times to the same element
1514 // if a cell is perforated by more than one well, so it should
1515 // not be OpenMP-parallelized.
1516 for (const auto& well : well_container_) {
1517 if (!well->isOperableAndSolvable() && !well->wellIsStopped()) {
1518 continue;
1519 }
1520 const auto& cells = well->cells();
1521 const auto& rates = well->connectionRates();
1522 for (unsigned perfIdx = 0; perfIdx < rates.size(); ++perfIdx) {
1523 unsigned cellIdx = cells[perfIdx];
1524 auto rate = rates[perfIdx];
1525 rate *= -1.0;
1526 VectorBlockType res(0.0);
1527 using MatrixBlockType = typename SparseMatrixAdapter::MatrixBlock;
1528 MatrixBlockType bMat(0.0);
1529 simulator_.model().linearizer().setResAndJacobi(res, bMat, rate);
1530 residual[cellIdx] += res;
1531 *diagMatAddress[cellIdx] += bMat;
1532 }
1533 }
1534 }
1535
1536
1537 template<typename TypeTag>
1538 void
1541 {
1542 int nw = this->numLocalWellsEnd();
1543 int rdofs = local_num_cells_;
1544 for (int i = 0; i < nw; ++i) {
1545 int wdof = rdofs + i;
1546 jacobian.entry(wdof,wdof) = 1.0;// better scaling ?
1547 }
1548 const auto wellconnections = this->getMaxWellConnections();
1549 for (int i = 0; i < nw; ++i) {
1550 const auto& perfcells = wellconnections[i];
1551 for (int perfcell : perfcells) {
1552 int wdof = rdofs + i;
1553 jacobian.entry(wdof, perfcell) = 0.0;
1554 jacobian.entry(perfcell, wdof) = 0.0;
1555 }
1556 }
1557 }
1558
1559
1560 template<typename TypeTag>
1561 void
1564 {
1565 auto loggerGuard = this->groupStateHelper().pushLogger();
1567 {
1568 for (const auto& well : well_container_) {
1569 const auto& cells = well->cells();
1570 x_local_.resize(cells.size());
1571
1572 for (size_t i = 0; i < cells.size(); ++i) {
1573 x_local_[i] = x[cells[i]];
1574 }
1575 well->recoverWellSolutionAndUpdateWellState(simulator_, x_local_,
1576 this->groupStateHelper(), this->wellState());
1577 }
1578 }
1579 OPM_END_PARALLEL_TRY_CATCH("recoverWellSolutionAndUpdateWellState() failed: ",
1580 simulator_.vanguard().grid().comm());
1581 }
1582
1583
1584 template<typename TypeTag>
1585 void
1587 recoverWellSolutionAndUpdateWellStateDomain(const BVector& x, const int domainIdx)
1588 {
1589 if (!nldd_) {
1590 OPM_THROW(std::logic_error, "Attempt to call NLDD method without a NLDD solver");
1591 }
1592
1593 return nldd_->recoverWellSolutionAndUpdateWellState(x, domainIdx);
1594 }
1595
1596
1597 template<typename TypeTag>
1600 getWellConvergence(const std::vector<Scalar>& B_avg, bool checkWellGroupControlsAndNetwork) const
1601 {
1602 // Get global (from all processes) convergence report.
1603 ConvergenceReport local_report;
1604 const auto& iterCtx = simulator_.problem().iterationContext();
1605 const bool relaxTolerance = iterCtx.shouldRelax(param_.strict_outer_iter_wells_ + 1);
1606 {
1607 auto logger_guard = this->groupStateHelper().pushLogger();
1608 for (const auto& well : well_container_) {
1609 if (well->isOperableAndSolvable() || well->wellIsStopped()) {
1610 local_report += well->getWellConvergence(
1611 this->groupStateHelper(), B_avg,
1612 relaxTolerance);
1613 } else {
1614 ConvergenceReport report;
1615 using CR = ConvergenceReport;
1616 report.setWellFailed({CR::WellFailure::Type::Unsolvable, CR::Severity::Normal, -1, well->name()});
1617 local_report += report;
1618 }
1619 }
1620 } // logger_guard goes out of scope here, before the OpmLog::debug() calls below
1621
1622 const Opm::Parallel::Communication comm = grid().comm();
1623 ConvergenceReport report = gatherConvergenceReport(local_report, comm);
1624
1625 if (checkWellGroupControlsAndNetwork) {
1626 // the well_group_control_changed info is already communicated
1627 report.setWellGroupTargetsViolated(this->lastReport().well_group_control_changed);
1628 report.setNetworkNotYetBalancedForceAnotherNewtonIteration(network_needs_more_balancing_force_another_newton_iteration_);
1629 }
1630
1631 if (this->terminal_output_) {
1632 // Log debug messages for NaN or too large residuals.
1633 for (const auto& f : report.wellFailures()) {
1634 if (f.severity() == ConvergenceReport::Severity::NotANumber) {
1635 OpmLog::debug("NaN residual found with phase " + std::to_string(f.phase()) + " for well " + f.wellName());
1636 } else if (f.severity() == ConvergenceReport::Severity::TooLarge) {
1637 OpmLog::debug("Too large residual found with phase " + std::to_string(f.phase()) + " for well " + f.wellName());
1638 }
1639 }
1640 }
1641 return report;
1642 }
1643
1644
1645
1646
1647
1648 template<typename TypeTag>
1649 void
1652 {
1653 // TODO: checking isOperableAndSolvable() ?
1654 for (auto& well : well_container_) {
1655 well->calculateExplicitQuantities(simulator_, this->groupStateHelper());
1656 }
1657 }
1658
1659
1660
1661
1662
1663 template<typename TypeTag>
1664 bool
1666 updateWellControls(DeferredLogger& deferred_logger)
1667 {
1668 OPM_TIMEFUNCTION();
1669 if (!this->wellsActive()) {
1670 return false;
1671 }
1672 const int episodeIdx = simulator_.episodeIndex();
1673 const auto& comm = simulator_.vanguard().grid().comm();
1674 size_t iter = 0;
1675 bool changed_well_group = false;
1676 const Group& fieldGroup = this->schedule().getGroup("FIELD", episodeIdx);
1677 // Check group individual constraints.
1678 // iterate a few times to make sure all constraints are honored
1679 const std::size_t max_iter = param_.well_group_constraints_max_iterations_;
1680 while(!changed_well_group && iter < max_iter) {
1681 changed_well_group = updateGroupControls(fieldGroup, deferred_logger, episodeIdx);
1682
1683 // Check wells' group constraints and communicate.
1684 bool changed_well_to_group = false;
1685 {
1686 OPM_TIMEBLOCK(UpdateWellControls);
1687 // For MS Wells a linear solve is performed below and the matrix might be singular.
1688 // We need to communicate the exception thrown to the others and rethrow.
1690 for (const auto& well : well_container_) {
1692 const bool changed_well = well->updateWellControl(
1693 simulator_, mode, this->groupStateHelper(), this->wellState()
1694 );
1695 if (changed_well) {
1696 changed_well_to_group = changed_well || changed_well_to_group;
1697 }
1698 }
1699 OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel: updating well controls failed: ",
1700 simulator_.gridView().comm());
1701 }
1702
1703 changed_well_to_group = comm.sum(static_cast<int>(changed_well_to_group));
1704 if (changed_well_to_group) {
1705 updateAndCommunicate(episodeIdx);
1706 changed_well_group = true;
1707 }
1708
1709 // Check individual well constraints and communicate.
1710 bool changed_well_individual = false;
1711 {
1712 // For MS Wells a linear solve is performed below and the matrix might be singular.
1713 // We need to communicate the exception thrown to the others and rethrow.
1715 for (const auto& well : well_container_) {
1717 const bool changed_well = well->updateWellControl(
1718 simulator_, mode, this->groupStateHelper(), this->wellState()
1719 );
1720 if (changed_well) {
1721 changed_well_individual = changed_well || changed_well_individual;
1722 }
1723 }
1724 OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel: updating well controls failed: ",
1725 simulator_.gridView().comm());
1726 }
1727
1728 changed_well_individual = comm.sum(static_cast<int>(changed_well_individual));
1729 if (changed_well_individual) {
1730 updateAndCommunicate(episodeIdx);
1731 changed_well_group = true;
1732 }
1733 iter++;
1734 }
1735
1736 // update wsolvent fraction for REIN wells
1737 this->updateWsolvent(fieldGroup, episodeIdx, this->nupcolWellState());
1738
1739 return changed_well_group;
1740 }
1741
1742
1743 template<typename TypeTag>
1744 void
1746 updateAndCommunicate(const int reportStepIdx)
1747 {
1748 const auto& iterCtx = simulator_.problem().iterationContext();
1749 this->updateAndCommunicateGroupData(reportStepIdx,
1750 iterCtx,
1751 param_.nupcol_group_rate_tolerance_,
1752 /*update_wellgrouptarget*/ true);
1753
1754 // updateWellStateWithTarget might throw for multisegment wells hence we
1755 // have a parallel try catch here to thrown on all processes.
1757 // if a well or group change control it affects all wells that are under the same group
1758 for (const auto& well : well_container_) {
1759 // We only want to update wells under group-control here
1760 const auto& ws = this->wellState().well(well->indexOfWell());
1761 if (ws.production_cmode == Well::ProducerCMode::GRUP ||
1762 ws.injection_cmode == Well::InjectorCMode::GRUP)
1763 {
1764 well->updateWellStateWithTarget(
1765 simulator_, this->groupStateHelper(), this->wellState()
1766 );
1767 }
1768 }
1769 OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel::updateAndCommunicate failed: ",
1770 simulator_.gridView().comm())
1771 this->updateAndCommunicateGroupData(reportStepIdx,
1772 iterCtx,
1773 param_.nupcol_group_rate_tolerance_,
1774 /*update_wellgrouptarget*/ true);
1775 }
1776
1777 template<typename TypeTag>
1778 bool
1780 updateGroupControls(const Group& group,
1781 DeferredLogger& deferred_logger,
1782 const int reportStepIdx)
1783 {
1784 OPM_TIMEFUNCTION();
1785 const auto& iterCtx = simulator_.problem().iterationContext();
1786 bool changed = false;
1787 // restrict the number of group switches but only after nupcol iterations.
1788 const int nupcol = this->schedule()[reportStepIdx].nupcol();
1789 const int max_number_of_group_switches = param_.max_number_of_group_switches_;
1790 const bool update_group_switching_log = !iterCtx.withinNupcol(nupcol);
1791 const bool changed_hc = this->checkGroupHigherConstraints(group, deferred_logger, reportStepIdx, max_number_of_group_switches, update_group_switching_log);
1792 if (changed_hc) {
1793 changed = true;
1794 updateAndCommunicate(reportStepIdx);
1795 }
1796
1797 bool changed_individual =
1799 updateGroupIndividualControl(group,
1800 reportStepIdx,
1801 max_number_of_group_switches,
1802 update_group_switching_log,
1803 this->switched_inj_groups_,
1804 this->switched_prod_groups_,
1805 this->closed_offending_wells_,
1806 this->groupState(),
1807 this->wellState(),
1808 deferred_logger);
1809
1810 if (changed_individual) {
1811 changed = true;
1812 updateAndCommunicate(reportStepIdx);
1813 }
1814 // call recursively down the group hierarchy
1815 for (const std::string& groupName : group.groups()) {
1816 bool changed_this = updateGroupControls(this->schedule().getGroup(groupName, reportStepIdx), deferred_logger, reportStepIdx);
1817 changed = changed || changed_this;
1818 }
1819 return changed;
1820 }
1821
1822 template<typename TypeTag>
1823 void
1825 updateWellTestState(const double simulationTime, WellTestState& wellTestState)
1826 {
1827 OPM_TIMEFUNCTION();
1828 auto logger_guard = this->groupStateHelper().pushLogger();
1829 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
1830 for (const auto& well : well_container_) {
1831 const auto& wname = well->name();
1832 const auto wasClosed = wellTestState.well_is_closed(wname);
1833 well->checkWellOperability(simulator_,
1834 this->wellState(),
1835 this->groupStateHelper());
1836 const bool under_zero_target =
1837 well->wellUnderZeroGroupRateTarget(this->groupStateHelper());
1838 well->updateWellTestState(this->wellState().well(wname),
1839 simulationTime,
1840 /*writeMessageToOPMLog=*/ true,
1841 under_zero_target,
1842 wellTestState,
1843 local_deferredLogger);
1844
1845 if (!wasClosed && wellTestState.well_is_closed(wname)) {
1846 this->closed_this_step_.insert(wname);
1847
1848 // maybe open a new well
1849 const WellEconProductionLimits& econ_production_limits = well->wellEcl().getEconLimits();
1850 if (econ_production_limits.validFollowonWell()) {
1851 const auto episode_idx = simulator_.episodeIndex();
1852 const auto follow_on_well = econ_production_limits.followonWell();
1853 if (!this->schedule().hasWell(follow_on_well, episode_idx)) {
1854 const auto msg = fmt::format("Well {} was closed. But the given follow on well {} does not exist."
1855 "The simulator continues without opening a follow on well.",
1856 wname, follow_on_well);
1857 local_deferredLogger.warning(msg);
1858 }
1859 auto& ws = this->wellState().well(follow_on_well);
1860 const bool success = ws.updateStatus(WellStatus::OPEN);
1861 if (success) {
1862 const auto msg = fmt::format("Well {} was closed. The follow on well {} opens instead.", wname, follow_on_well);
1863 local_deferredLogger.info(msg);
1864 } else {
1865 const auto msg = fmt::format("Well {} was closed. The follow on well {} is already open.", wname, follow_on_well);
1866 local_deferredLogger.warning(msg);
1867 }
1868 }
1869
1870 }
1871 }
1872
1873 for (const auto& [group_name, to] : this->closed_offending_wells_) {
1874 if (this->hasOpenLocalWell(to.second) &&
1875 !this->wasDynamicallyShutThisTimeStep(to.second))
1876 {
1877 wellTestState.close_well(to.second,
1878 WellTestConfig::Reason::GROUP,
1879 simulationTime);
1880 this->updateClosedWellsThisStep(to.second);
1881 const std::string msg =
1882 fmt::format("Procedure on exceeding {} limit is WELL for group {}. "
1883 "Well {} is {}.",
1884 to.first,
1885 group_name,
1886 to.second,
1887 "shut");
1888 local_deferredLogger.info(msg);
1889 }
1890 }
1891 }
1892
1893
1894 template<typename TypeTag>
1895 void
1897 const WellState<Scalar, IndexTraits>& well_state_copy,
1898 std::string& exc_msg,
1899 ExceptionType::ExcEnum& exc_type)
1900 {
1901 OPM_TIMEFUNCTION();
1902 const int np = this->numPhases();
1903 std::vector<Scalar> potentials;
1904 const auto& well = well_container_[widx];
1905 std::string cur_exc_msg;
1906 auto cur_exc_type = ExceptionType::NONE;
1907 try {
1908 well->computeWellPotentials(simulator_, well_state_copy, this->groupStateHelper(), potentials);
1909 }
1910 // catch all possible exception and store type and message.
1911 OPM_PARALLEL_CATCH_CLAUSE(cur_exc_type, cur_exc_msg);
1912 if (cur_exc_type != ExceptionType::NONE) {
1913 exc_msg += fmt::format("\nFor well {}: {}", well->name(), cur_exc_msg);
1914 }
1915 exc_type = std::max(exc_type, cur_exc_type);
1916 // Store it in the well state
1917 // potentials is resized and set to zero in the beginning of well->ComputeWellPotentials
1918 // and updated only if sucessfull. i.e. the potentials are zero for exceptions
1919 auto& ws = this->wellState().well(well->indexOfWell());
1920 for (int p = 0; p < np; ++p) {
1921 // make sure the potentials are positive
1922 ws.well_potentials[p] = std::max(Scalar{0.0}, potentials[p]);
1923 }
1924 }
1925
1926
1927
1928 template <typename TypeTag>
1929 void
1932 {
1933 for (const auto& wellPtr : this->well_container_) {
1934 this->calculateProductivityIndexValues(wellPtr.get(), deferred_logger);
1935 }
1936 }
1937
1938
1939
1940
1941
1942 template <typename TypeTag>
1943 void
1945 calculateProductivityIndexValuesShutWells(const int reportStepIdx,
1946 DeferredLogger& deferred_logger)
1947 {
1948 // For the purpose of computing PI/II values, it is sufficient to
1949 // construct StandardWell instances only. We don't need to form
1950 // well objects that honour the 'isMultisegment()' flag of the
1951 // corresponding "this->wells_ecl_[shutWell]".
1952
1953 for (const auto& shutWell : this->local_shut_wells_) {
1954 if (!this->wells_ecl_[shutWell].hasConnections()) {
1955 // No connections in this well. Nothing to do.
1956 continue;
1957 }
1958
1959 auto wellPtr = this->template createTypedWellPointer
1960 <StandardWell<TypeTag>>(shutWell, reportStepIdx);
1961
1962 wellPtr->init(this->depth_, this->gravity_, this->B_avg_, true);
1963
1964 this->calculateProductivityIndexValues(wellPtr.get(), deferred_logger);
1965 }
1966 }
1967
1968
1969
1970
1971
1972 template <typename TypeTag>
1973 void
1976 DeferredLogger& deferred_logger)
1977 {
1978 wellPtr->updateProductivityIndex(this->simulator_,
1979 this->prod_index_calc_[wellPtr->indexOfWell()],
1980 this->wellState(),
1981 deferred_logger);
1982 }
1983
1984
1985
1986 template<typename TypeTag>
1987 void
1989 prepareTimeStep(DeferredLogger& deferred_logger)
1990 {
1991 // Check if there is a network with active prediction wells at this time step.
1992 const auto episodeIdx = simulator_.episodeIndex();
1993 this->network_.updateActiveState(episodeIdx);
1994
1995 // Rebalance the network initially if any wells in the network have status changes
1996 // (Need to check this before clearing events)
1997 const bool do_prestep_network_rebalance =
1998 param_.pre_solve_network_ && this->network_.needPreStepRebalance(episodeIdx);
1999
2000 for (const auto& well : well_container_) {
2001 auto& events = this->wellState().well(well->indexOfWell()).events;
2002 if (events.hasEvent(WellState<Scalar, IndexTraits>::event_mask)) {
2003 well->updateWellStateWithTarget(
2004 simulator_, this->groupStateHelper(), this->wellState()
2005 );
2006 well->updatePrimaryVariables(this->groupStateHelper());
2007 // There is no new well control change input within a report step,
2008 // so next time step, the well does not consider to have effective events anymore.
2010 }
2011 // these events only work for the first time step within the report step
2012 if (events.hasEvent(ScheduleEvents::REQUEST_OPEN_WELL)) {
2013 events.clearEvent(ScheduleEvents::REQUEST_OPEN_WELL);
2014 }
2015 // solve the well equation initially to improve the initial solution of the well model
2016 if (param_.solve_welleq_initially_ && well->isOperableAndSolvable()) {
2017 try {
2018 well->solveWellEquation(
2019 simulator_, this->groupStateHelper(), this->wellState()
2020 );
2021 } catch (const std::exception& e) {
2022 const std::string msg = "Compute initial well solution for " + well->name() + " initially failed. Continue with the previous rates";
2023 deferred_logger.warning("WELL_INITIAL_SOLVE_FAILED", msg);
2024 }
2025 }
2026 // If we're using local well solves that include control switches, they also update
2027 // operability, so reset before main iterations begin
2028 well->resetWellOperability();
2029 }
2030 updatePrimaryVariables();
2031
2032 // Actually do the pre-step network rebalance, using the updated well states and initial solutions
2033 if (do_prestep_network_rebalance) {
2034 network_.doPreStepRebalance(deferred_logger);
2035 }
2036 }
2037
2038 template<typename TypeTag>
2039 void
2042 {
2043 std::vector< Scalar > B_avg(numConservationQuantities(), Scalar() );
2044 const auto& grid = simulator_.vanguard().grid();
2045 const auto& gridView = grid.leafGridView();
2046 ElementContext elemCtx(simulator_);
2047
2049 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
2050 elemCtx.updatePrimaryStencil(elem);
2051 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
2052
2053 const auto& intQuants = elemCtx.intensiveQuantities(/*spaceIdx=*/0, /*timeIdx=*/0);
2054 const auto& fs = intQuants.fluidState();
2055
2056 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
2057 {
2058 if (!FluidSystem::phaseIsActive(phaseIdx)) {
2059 continue;
2060 }
2061
2062 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
2063 auto& B = B_avg[ compIdx ];
2064
2065 B += 1 / fs.invB(phaseIdx).value();
2066 }
2067 if constexpr (has_solvent_) {
2068 auto& B = B_avg[solventSaturationIdx];
2069 B += 1 / intQuants.solventInverseFormationVolumeFactor().value();
2070 }
2071 }
2072 OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel::updateAverageFormationFactor() failed: ", grid.comm())
2073
2074 // compute global average
2075 grid.comm().sum(B_avg.data(), B_avg.size());
2076 B_avg_.resize(B_avg.size());
2077 std::ranges::transform(B_avg, B_avg_.begin(),
2078 [gcells = global_num_cells_](const auto bval)
2079 { return bval / gcells; });
2080 }
2081
2082
2083
2084
2085
2086 template<typename TypeTag>
2087 void
2090 {
2091 for (const auto& well : well_container_) {
2092 well->updatePrimaryVariables(this->groupStateHelper());
2093 }
2094 }
2095
2096 template<typename TypeTag>
2097 void
2099 {
2100 const auto& grid = simulator_.vanguard().grid();
2101 const auto& eclProblem = simulator_.problem();
2102 const unsigned numCells = grid.size(/*codim=*/0);
2103
2104 this->pvt_region_idx_.resize(numCells);
2105 for (unsigned cellIdx = 0; cellIdx < numCells; ++cellIdx) {
2106 this->pvt_region_idx_[cellIdx] =
2107 eclProblem.pvtRegionIndex(cellIdx);
2108 }
2109 }
2110
2111 // The number of components in the model.
2112 template<typename TypeTag>
2113 int
2115 {
2116 // The numPhases() functions returns 1-3, depending on which
2117 // of the (oil, water, gas) phases are active. For each of those phases,
2118 // if the phase is active the corresponding component is present and
2119 // conserved.
2120 // Apart from (oil, water, gas), in the current well model only solvent
2121 // is explicitly modelled as a conserved quantity (polymer, energy, salt
2122 // etc. are not), unlike the reservoir part where all such quantities are
2123 // conserved. This function must therefore be updated when/if we add
2124 // more conserved quantities in the well model.
2125 return this->numPhases() + has_solvent_;
2126 }
2127
2128 template<typename TypeTag>
2129 void
2131 {
2132 const auto& eclProblem = simulator_.problem();
2133 depth_.resize(local_num_cells_);
2134 for (unsigned cellIdx = 0; cellIdx < local_num_cells_; ++cellIdx) {
2135 depth_[cellIdx] = eclProblem.dofCenterDepth(cellIdx);
2136 }
2137 }
2138
2139 template<typename TypeTag>
2142 getWell(const std::string& well_name) const
2143 {
2144 // finding the iterator of the well in wells_ecl
2145 const auto well =
2146 std::ranges::find_if(well_container_,
2147 [&well_name](const WellInterfacePtr& elem) -> bool
2148 { return elem->name() == well_name; });
2149
2150 assert(well != well_container_.end());
2151
2152 return **well;
2153 }
2154
2155 template <typename TypeTag>
2156 int
2158 reportStepIndex() const
2159 {
2160 return std::max(this->simulator_.episodeIndex(), 0);
2161 }
2162
2163
2164
2165
2166
2167 template<typename TypeTag>
2168 void
2170 calcResvCoeff(const int fipnum,
2171 const int pvtreg,
2172 const std::vector<Scalar>& production_rates,
2173 std::vector<Scalar>& resv_coeff) const
2174 {
2175 rateConverter_->calcCoeff(fipnum, pvtreg, production_rates, resv_coeff);
2176 }
2177
2178 template<typename TypeTag>
2179 void
2181 calcInjResvCoeff(const int fipnum,
2182 const int pvtreg,
2183 std::vector<Scalar>& resv_coeff) const
2184 {
2185 rateConverter_->calcInjCoeff(fipnum, pvtreg, resv_coeff);
2186 }
2187
2188
2189 template <typename TypeTag>
2190 void
2193 {
2194 if constexpr (energyModuleType_ == EnergyModules::FullyImplicitThermal) {
2195 const int np = this->numPhases();
2196 const int nw = this->numLocalWells();
2197 for (auto wellID = 0*nw; wellID < nw; ++wellID) {
2198 const Well& well = this->wells_ecl_[wellID];
2199 auto& ws = this->wellState().well(wellID);
2200 if (well.isInjector()) {
2201 if (ws.status != WellStatus::STOP) {
2202 this->wellState().well(wellID).temperature = well.inj_temperature();
2203 continue;
2204 }
2205 }
2206 std::array<Scalar,2> weighted{0.0,0.0};
2207 auto& [weighted_temperature, total_weight] = weighted;
2208 const auto& well_info = this->local_parallel_well_info_[wellID].get();
2209 using int_type = decltype(this->well_perf_data_[wellID].size());
2210 for (int_type perf = 0, end_perf = this->well_perf_data_[wellID].size(); perf < end_perf; ++perf) {
2211 const int cell_idx = this->well_perf_data_[wellID][perf].cell_index;
2212 const auto& intQuants = simulator_.model().intensiveQuantities(cell_idx, /*timeIdx=*/0);
2213 const auto& fs = intQuants.fluidState();
2214 Scalar weight_factor = computeTemperatureWeightFactor(perf, np, fs, ws);
2215 total_weight += weight_factor;
2216 weighted_temperature += weight_factor * fs.temperature(/*phaseIdx*/0).value();
2217 }
2218 well_info.communication().sum(weighted.data(), 2);
2219 this->wellState().well(wellID).temperature = weighted_temperature / total_weight;
2220 }
2221 }
2222 }
2223
2224
2225 template <typename TypeTag>
2227 assignWellTracerRates(data::Wells& wsrpt) const
2228 {
2229 const auto reportStepIdx = static_cast<unsigned int>(this->reportStepIndex());
2230 const auto& trMod = this->simulator_.problem().tracerModel();
2231
2232 BlackoilWellModelGeneric<Scalar, IndexTraits>::assignWellTracerRates(wsrpt, trMod.getWellTracerRates(), reportStepIdx);
2233 BlackoilWellModelGeneric<Scalar, IndexTraits>::assignWellTracerRates(wsrpt, trMod.getWellFreeTracerRates(), reportStepIdx);
2234 BlackoilWellModelGeneric<Scalar, IndexTraits>::assignWellTracerRates(wsrpt, trMod.getWellSolTracerRates(), reportStepIdx);
2235
2236 this->assignMswTracerRates(wsrpt, trMod.getMswTracerRates(), reportStepIdx);
2237 }
2238
2239 template <typename TypeTag>
2240 void BlackoilWellModel<TypeTag>::
2241 assignWellSpeciesRates(data::Wells& wsrpt) const
2242 {
2243 const auto reportStepIdx = static_cast<unsigned int>(this->reportStepIndex());
2244 const auto& geochemMod = this->simulator_.problem().geochemistryModel();
2245
2246 BlackoilWellModelGeneric<Scalar, IndexTraits>::assignWellTracerRates(wsrpt, geochemMod.getWellSpeciesRates(), reportStepIdx);
2247
2248 this->assignMswTracerRates(wsrpt, geochemMod.getMswSpeciesRates(), reportStepIdx);
2249 }
2250
2251 template <typename TypeTag>
2252 [[nodiscard]] auto BlackoilWellModel<TypeTag>::rsConstInfo() const
2253 -> typename WellState<Scalar,IndexTraits>::RsConstInfo
2254 {
2255 if (! FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) ||
2256 ! FluidSystem::enableConstantRs())
2257 {
2258 return {};
2259 }
2260
2261 const auto& rsConstTables = this->eclState_
2262 .getTableManager().getRsconstTables();
2263
2264 if (rsConstTables.empty() ||
2265 (rsConstTables[0].numRows() != std::size_t{1}))
2266 {
2267 return {};
2268 }
2269
2270 const auto rsConst = rsConstTables[0].getColumn(0).front();
2271
2272 return { true, static_cast<Scalar>(rsConst) };
2273 }
2274
2275} // namespace Opm
2276
2277#endif // OPM_BLACKOILWELLMODEL_IMPL_HEADER_INCLUDED
#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:97
BlackoilWellModelWBP< GetPropType< TypeTag, Properties::Scalar >, GetPropType< TypeTag, Properties::FluidSystem >::IndexTraitsType > wbp_
Definition: BlackoilWellModelGeneric.hpp:518
std::vector< ParallelWellInfo< GetPropType< TypeTag, Properties::Scalar > > > parallel_well_info_
Definition: BlackoilWellModelGeneric.hpp:545
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:98
void initializeGroupStructure(const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:297
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:2170
void prepareTimeStep(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1989
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:1303
WellInterfacePtr createWellPointer(const int wellID, const int report_step) const
Definition: BlackoilWellModel_impl.hpp:1099
void prepareWellsBeforeAssembling(const double dt)
Definition: BlackoilWellModel_impl.hpp:1379
void init()
Definition: BlackoilWellModel_impl.hpp:166
const Simulator & simulator() const
Definition: BlackoilWellModel.hpp:370
std::vector< Scalar > depth_
Definition: BlackoilWellModel.hpp:516
std::size_t global_num_cells_
Definition: BlackoilWellModel.hpp:512
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilWellModel.hpp:107
void initWellContainer(const int reportStepIdx) override
Definition: BlackoilWellModel_impl.hpp:185
void beginReportStep(const int time_step)
Definition: BlackoilWellModel_impl.hpp:202
const WellInterface< TypeTag > & getWell(const std::string &well_name) const
Definition: BlackoilWellModel_impl.hpp:2142
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: BlackoilWellModel.hpp:103
Dune::FieldVector< Scalar, numEq > VectorBlockType
Definition: BlackoilWellModel.hpp:130
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: BlackoilWellModel.hpp:104
GetPropType< TypeTag, Properties::Grid > Grid
Definition: BlackoilWellModel.hpp:101
int numConservationQuantities() const
Definition: BlackoilWellModel_impl.hpp:2114
bool updateWellControls(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1666
int reportStepIndex() const
Definition: BlackoilWellModel_impl.hpp:2158
void calculateProductivityIndexValues(DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1931
void extractLegacyDepth_()
Definition: BlackoilWellModel_impl.hpp:2130
void extractLegacyCellPvtRegionIndex_()
Definition: BlackoilWellModel_impl.hpp:2098
void recoverWellSolutionAndUpdateWellStateDomain(const BVector &x, const int domainIdx)
Definition: BlackoilWellModel_impl.hpp:1587
void updateAverageFormationFactor()
Definition: BlackoilWellModel_impl.hpp:2041
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilWellModel.hpp:106
void initializeWellState(const int timeStepIdx)
Definition: BlackoilWellModel_impl.hpp:847
const Grid & grid() const
Definition: BlackoilWellModel.hpp:367
void updatePrimaryVariables()
Definition: BlackoilWellModel_impl.hpp:2089
void computeWellTemperature()
Definition: BlackoilWellModel_impl.hpp:2192
void addWellPressureEquations(PressureMatrix &jacobian, const BVector &weights, const bool use_well_weights) const
Definition: BlackoilWellModel_impl.hpp:1488
const SimulatorReportSingle & lastReport() const
Definition: BlackoilWellModel_impl.hpp:726
bool updateWellControlsAndNetwork(const bool mandatory_network_balance, const double dt, DeferredLogger &local_deferredLogger)
Definition: BlackoilWellModel_impl.hpp:1244
void addWellContributions(SparseMatrixAdapter &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1478
void assembleWellEq(const double dt)
Definition: BlackoilWellModel_impl.hpp:1367
WellInterfacePtr createWellForWellTest(const std::string &well_name, const int report_step, DeferredLogger &deferred_logger) const
Definition: BlackoilWellModel_impl.hpp:1150
void calculateExplicitQuantities() const
Definition: BlackoilWellModel_impl.hpp:1651
void updateAndCommunicate(const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:1746
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: BlackoilWellModel.hpp:289
void computeTotalRatesForDof(RateVector &rate, unsigned globalIdx) const
Definition: BlackoilWellModel_impl.hpp:810
void beginTimeStep()
Definition: BlackoilWellModel_impl.hpp:329
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: BlackoilWellModel.hpp:108
bool updateGroupControls(const Group &group, DeferredLogger &deferred_logger, const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:1780
void calcInjResvCoeff(const int fipnum, const int pvtreg, std::vector< Scalar > &resv_coeff) const override
Definition: BlackoilWellModel_impl.hpp:2181
void initializeLocalWellStructure(const int reportStepIdx, const bool enableWellPIScaling)
Definition: BlackoilWellModel_impl.hpp:251
Dune::BlockVector< VectorBlockType > BVector
Definition: BlackoilWellModel.hpp:131
BlackoilWellModel(Simulator &simulator)
Definition: BlackoilWellModel_impl.hpp:79
void wellTesting(const int timeStepIdx, const double simulationTime, DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:649
ConvergenceReport getWellConvergence(const std::vector< Scalar > &B_avg, const bool checkWellGroupControlsAndNetwork=false) const
Definition: BlackoilWellModel_impl.hpp:1600
typename FluidSystem::IndexTraitsType IndexTraits
Definition: BlackoilWellModel.hpp:114
void updateCellRatesForDomain(int domainIndex, const std::map< std::string, int > &well_domain_map)
Definition: BlackoilWellModel_impl.hpp:1425
void assembleWellEqWithoutIteration(const double dt)
Definition: BlackoilWellModel_impl.hpp:1393
void updateCellRates()
Definition: BlackoilWellModel_impl.hpp:1413
void assemble(const double dt)
Definition: BlackoilWellModel_impl.hpp:1175
std::size_t local_num_cells_
Definition: BlackoilWellModel.hpp:514
bool alternative_well_rate_init_
Definition: BlackoilWellModel.hpp:517
void timeStepSucceeded(const double simulationTime, const double dt)
Definition: BlackoilWellModel_impl.hpp:736
std::unique_ptr< WellType > createTypedWellPointer(const int wellID, const int time_step) const
Definition: BlackoilWellModel_impl.hpp:1119
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:1896
Simulator & simulator_
Definition: BlackoilWellModel.hpp:486
void createWellContainer(const int report_step) override
Definition: BlackoilWellModel_impl.hpp:890
std::unique_ptr< WellInterface< TypeTag > > WellInterfacePtr
Definition: BlackoilWellModel.hpp:187
void updateWellTestState(const double simulationTime, WellTestState &wellTestState)
upate the wellTestState related to economic limits
Definition: BlackoilWellModel_impl.hpp:1825
void addReservoirSourceTerms(GlobalEqVector &residual, const std::vector< typename SparseMatrixAdapter::MatrixBlock * > &diagMatAddress) const
Definition: BlackoilWellModel_impl.hpp:1510
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:1563
void addWellPressureEquationsStruct(PressureMatrix &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1540
void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1945
void endReportStep()
Definition: BlackoilWellModel_impl.hpp:709
void initializeSources(typename ParallelWBPCalculation< Scalar >::GlobalToLocal index, typename ParallelWBPCalculation< Scalar >::Evaluator eval)
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
Guard for managing DeferredLogger lifecycle in ReservoirCoupling.
Definition: ReservoirCoupling.hpp:88
Definition: StandardWell.hpp:60
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:77
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)
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34