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