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 template<typename TypeTag>
646 void
647 BlackoilWellModel<TypeTag>::
648 rescoupSyncSummaryData()
649 {
650 // Reservoir coupling: exchange production data between slaves and master.
651 //
652 // Master side: after its first substep, the master blocks here until all
653 // slaves have completed the sync step and sent their production data.
654 // This ensures evalSummaryState() (called next in endTimeStep) and all
655 // subsequent master substeps have correct slave production rates.
656 //
657 // Slave side: on the last substep of the sync step, the slave sends its
658 // production data to the master. The master is already waiting at this
659 // point (blocked on MPI_Recv from its first substep's timeStepSucceeded).
660 if (this->isReservoirCouplingMaster()) {
661 if (this->reservoirCouplingMaster().needsSlaveDataReceive()) {
662 this->receiveSlaveGroupData();
663 this->reservoirCouplingMaster().setNeedsSlaveDataReceive(false);
664 }
665 }
666 if (this->isReservoirCouplingSlave()) {
667 if (this->reservoirCouplingSlave().isLastSubstepOfSyncTimestep()) {
668 this->sendSlaveGroupDataToMaster();
669 }
670 }
671 }
672
673#endif // RESERVOIR_COUPLING_ENABLED
674
675 template<typename TypeTag>
676 void
678 const double simulationTime,
679 DeferredLogger& deferred_logger)
680 {
681 for (const std::string& well_name : this->getWellsForTesting(timeStepIdx, simulationTime)) {
682 const Well& wellEcl = this->schedule().getWell(well_name, timeStepIdx);
683 if (wellEcl.getStatus() == Well::Status::SHUT)
684 continue;
685
686 WellInterfacePtr well = createWellForWellTest(well_name, timeStepIdx, deferred_logger);
687 // some preparation before the well can be used
688 well->init(depth_, gravity_, B_avg_, true);
689
690 Scalar well_efficiency_factor = wellEcl.getEfficiencyFactor() *
691 this->wellState().getGlobalEfficiencyScalingFactor(well_name);
692 this->groupStateHelper().accumulateGroupEfficiencyFactor(
693 this->schedule().getGroup(wellEcl.groupName(), timeStepIdx),
694 well_efficiency_factor
695 );
696
697 well->setWellEfficiencyFactor(well_efficiency_factor);
698 well->setVFPProperties(this->vfp_properties_.get());
699 well->setGuideRate(&this->guideRate_);
700
701 // initialize rates/previous rates to prevent zero fractions in vfp-interpolation
702 if (well->isProducer() && alternative_well_rate_init_) {
703 well->initializeProducerWellState(simulator_, this->wellState(), deferred_logger);
704 }
705 if (well->isVFPActive(deferred_logger)) {
706 well->setPrevSurfaceRates(this->wellState(), this->prevWellState());
707 }
708
709 const auto& network = this->schedule()[timeStepIdx].network();
710 if (network.active()) {
711 this->network_.initializeWell(*well);
712 }
713 try {
714 using GLiftEclWells = typename GasLiftGroupInfo<Scalar, IndexTraits>::GLiftEclWells;
715 GLiftEclWells ecl_well_map;
716 gaslift_.initGliftEclWellMap(well_container_, ecl_well_map);
717 well->wellTesting(simulator_,
718 simulationTime,
719 this->groupStateHelper(),
720 this->wellState(),
721 this->wellTestState(),
722 ecl_well_map,
723 this->well_open_times_);
724 } catch (const std::exception& e) {
725 const std::string msg =
726 fmt::format(fmt::runtime("Exception during testing of well: {}. The well will not open.\n"
727 "Exception message: {}"), wellEcl.name(), e.what());
728 deferred_logger.warning("WELL_TESTING_FAILED", msg);
729 }
730 }
731 }
732
733 // called at the end of a report step
734 template<typename TypeTag>
735 void
738 {
739 // Clear the communication data structures for above values.
740 for (auto&& pinfo : this->local_parallel_well_info_)
741 {
742 pinfo.get().clear();
743 }
744 }
745
746
747
748
749
750 // called at the end of a report step
751 template<typename TypeTag>
754 lastReport() const {return last_report_; }
755
756
757
758
759
760 // called at the end of a time step
761 template<typename TypeTag>
762 void
764 timeStepSucceeded(const double simulationTime, const double dt)
765 {
766 this->closed_this_step_.clear();
767
768 // time step is finished and we are not any more at the beginning of an report step
769 this->report_step_starts_ = false;
770 const int reportStepIdx = simulator_.episodeIndex();
771
772 auto logger_guard = this->groupStateHelper().pushLogger();
773 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
774 for (const auto& well : well_container_) {
775 if (getPropValue<TypeTag, Properties::EnablePolymerMW>() && well->isInjector()) {
776 well->updateWaterThroughput(dt, this->wellState());
777 }
778 }
779 // update connection transmissibility factor and d factor (if applicable) in the wellstate
780 for (const auto& well : well_container_) {
781 well->updateConnectionTransmissibilityFactor(simulator_, this->wellState().well(well->indexOfWell()));
782 well->updateConnectionDFactor(simulator_, this->wellState().well(well->indexOfWell()));
783 }
784
785 if (Indices::waterEnabled) {
786 this->updateFiltrationModelsPostStep(dt, FluidSystem::waterPhaseIdx, local_deferredLogger);
787 }
788
789 // WINJMULT: At the end of the time step, update the inj_multiplier saved in WellState for later use
790 this->updateInjMult(local_deferredLogger);
791
792 // report well switching
793 for (const auto& well : well_container_) {
794 well->reportWellSwitching(this->wellState().well(well->indexOfWell()), local_deferredLogger);
795 }
796 // report group switching
797 if (this->terminal_output_) {
798 this->reportGroupSwitching(local_deferredLogger);
799 }
800
801 // update the rate converter with current averages pressures etc in
802 rateConverter_->template defineState<ElementContext>(simulator_);
803
804 // calculate the well potentials
805 try {
806 this->updateWellPotentials(reportStepIdx,
807 /*onlyAfterEvent*/false,
808 simulator_.vanguard().summaryConfig(),
809 local_deferredLogger);
810 } catch ( std::runtime_error& e ) {
811 const std::string msg = "A zero well potential is returned for output purposes. ";
812 local_deferredLogger.warning("WELL_POTENTIAL_CALCULATION_FAILED", msg);
813 }
814
815 updateWellTestState(simulationTime, this->wellTestState());
816
817 // check group sales limits at the end of the timestep
818 const Group& fieldGroup = this->schedule_.getGroup("FIELD", reportStepIdx);
819 this->checkGEconLimits(fieldGroup, simulationTime,
820 simulator_.episodeIndex(), local_deferredLogger);
821 this->checkGconsaleLimits(fieldGroup, this->wellState(),
822 simulator_.episodeIndex(), local_deferredLogger);
823
824 this->calculateProductivityIndexValues(local_deferredLogger);
825
826 this->groupStateHelper().updateNONEProductionGroups();
827
828#ifdef RESERVOIR_COUPLING_ENABLED
829 this->rescoupSyncSummaryData();
830#endif
831 this->commitWGState();
832
833 //reporting output temperatures
834 this->computeWellTemperature();
835 }
836
837
838 template<typename TypeTag>
839 void
842 unsigned elemIdx) const
843 {
844 rate = 0;
845
846 if (!is_cell_perforated_[elemIdx] || cellRates_.count(elemIdx) == 0) {
847 return;
848 }
849
850 rate = cellRates_.at(elemIdx);
851 }
852
853
854 template<typename TypeTag>
855 template <class Context>
856 void
859 const Context& context,
860 unsigned spaceIdx,
861 unsigned timeIdx) const
862 {
863 rate = 0;
864 int elemIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
865
866 if (!is_cell_perforated_[elemIdx] || cellRates_.count(elemIdx) == 0) {
867 return;
868 }
869
870 rate = cellRates_.at(elemIdx);
871 }
872
873
874
875 template<typename TypeTag>
876 void
878 initializeWellState(const int timeStepIdx)
879 {
880 const auto pressIx = []()
881 {
882 if (Indices::oilEnabled) { return FluidSystem::oilPhaseIdx; }
883 if (Indices::waterEnabled) { return FluidSystem::waterPhaseIdx; }
884
885 return FluidSystem::gasPhaseIdx;
886 }();
887
888 auto cellPressures = std::vector<Scalar>(this->local_num_cells_, Scalar{0});
889 auto cellTemperatures = std::vector<Scalar>(this->local_num_cells_, Scalar{0});
890
891 auto elemCtx = ElementContext { this->simulator_ };
892 const auto& gridView = this->simulator_.vanguard().gridView();
893
895 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
896 elemCtx.updatePrimaryStencil(elem);
897 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
898
899 const auto ix = elemCtx.globalSpaceIndex(/*spaceIdx=*/0, /*timeIdx=*/0);
900 const auto& fs = elemCtx.intensiveQuantities(/*spaceIdx=*/0, /*timeIdx=*/0).fluidState();
901
902 cellPressures[ix] = fs.pressure(pressIx).value();
903 cellTemperatures[ix] = fs.temperature(0).value();
904 }
905 OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel::initializeWellState() failed: ",
906 this->simulator_.vanguard().grid().comm());
907
908 this->wellState().init(cellPressures, cellTemperatures, this->schedule(), this->wells_ecl_,
909 this->local_parallel_well_info_, timeStepIdx,
910 &this->prevWellState(), this->well_perf_data_,
911 this->summaryState(), simulator_.vanguard().enableDistributedWells());
912 }
913
914
915
916
917
918 template<typename TypeTag>
919 void
921 createWellContainer(const int report_step)
922 {
923 auto logger_guard = this->groupStateHelper().pushLogger();
924 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
925
926 const int nw = this->numLocalWells();
927
928 well_container_.clear();
929
930 if (nw > 0) {
931 well_container_.reserve(nw);
932
933 const auto& wmatcher = this->schedule().wellMatcher(report_step);
934 const auto& wcycle = this->schedule()[report_step].wcycle.get();
935
936 // First loop and check for status changes. This is necessary
937 // as wcycle needs the updated open/close times.
938 std::ranges::for_each(this->wells_ecl_,
939 [this, &wg_events = this->report_step_start_events_](const auto& well_ecl)
940 {
941 if (!well_ecl.hasConnections()) {
942 // No connections in this well. Nothing to do.
943 return;
944 }
945
946 constexpr auto events_mask = ScheduleEvents::WELL_STATUS_CHANGE |
947 ScheduleEvents::REQUEST_OPEN_WELL |
948 ScheduleEvents::REQUEST_SHUT_WELL;
949 const bool well_event =
950 this->report_step_starts_ &&
951 wg_events.hasEvent(well_ecl.name(), events_mask);
952 // WCYCLE is suspendended by explicit SHUT events by the user.
953 // and restarted after explicit OPEN events.
954 // Note: OPEN or SHUT event does not necessary mean the well
955 // actually opened or shut at this point as the simulator could
956 // have done this by operabilty checks and well testing. This
957 // may need further testing and imply code changes to cope with
958 // these corner cases.
959 if (well_event) {
960 if (well_ecl.getStatus() == WellStatus::OPEN) {
961 this->well_open_times_.insert_or_assign(well_ecl.name(),
962 this->simulator_.time());
963 this->well_close_times_.erase(well_ecl.name());
964 } else if (well_ecl.getStatus() == WellStatus::SHUT) {
965 this->well_close_times_.insert_or_assign(well_ecl.name(),
966 this->simulator_.time());
967 this->well_open_times_.erase(well_ecl.name());
968 }
969 }
970 });
971
972 // Grab wcycle states. This needs to run before the schedule gets processed
973 const auto cycle_states = wcycle.wellStatus(this->simulator_.time(),
974 wmatcher,
975 this->well_open_times_,
976 this->well_close_times_);
977
978 for (int w = 0; w < nw; ++w) {
979 const Well& well_ecl = this->wells_ecl_[w];
980
981 if (!well_ecl.hasConnections()) {
982 // No connections in this well. Nothing to do.
983 continue;
984 }
985
986 const std::string& well_name = well_ecl.name();
987 const auto well_status = this->schedule()
988 .getWell(well_name, report_step).getStatus();
989
990 const bool shut_event = this->wellState().well(w).events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE)
991 && well_status == Well::Status::SHUT;
992 const bool open_event = this->wellState().well(w).events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE)
993 && well_status == Well::Status::OPEN;
994 const auto& ws = this->wellState().well(well_name);
995
996 if (shut_event && ws.status != Well::Status::SHUT) {
997 this->closed_this_step_.insert(well_name);
998 this->wellState().shutWell(w);
999 } else if (open_event && ws.status != Well::Status::OPEN) {
1000 this->wellState().openWell(w);
1001 }
1002
1003 // A new WCON keywords can re-open a well that was closed/shut due to Physical limit
1004 if (this->wellTestState().well_is_closed(well_name)) {
1005 // The well was shut this timestep, we are most likely retrying
1006 // a timestep without the well in question, after it caused
1007 // repeated timestep cuts. It should therefore not be opened,
1008 // even if it was new or received new targets this report step.
1009 const bool closed_this_step = (this->wellTestState().lastTestTime(well_name) == simulator_.time());
1010 // TODO: more checking here, to make sure this standard more specific and complete
1011 // maybe there is some WCON keywords will not open the well
1012 auto& events = this->wellState().well(w).events;
1013 if (events.hasEvent(ScheduleEvents::REQUEST_OPEN_WELL)) {
1014 if (!closed_this_step) {
1015 this->wellTestState().open_well(well_name);
1016 this->wellTestState().open_completions(well_name);
1017 this->well_open_times_.insert_or_assign(well_name,
1018 this->simulator_.time());
1019 this->well_close_times_.erase(well_name);
1020 }
1021 events.clearEvent(ScheduleEvents::REQUEST_OPEN_WELL);
1022 }
1023 }
1024
1025 // TODO: should we do this for all kinds of closing reasons?
1026 // something like wellTestState().hasWell(well_name)?
1027 if (this->wellTestState().well_is_closed(well_name))
1028 {
1029 if (well_ecl.getAutomaticShutIn()) {
1030 // shut wells are not added to the well container
1031 this->wellState().shutWell(w);
1032 this->well_close_times_.erase(well_name);
1033 this->well_open_times_.erase(well_name);
1034 continue;
1035 } else {
1036 if (!well_ecl.getAllowCrossFlow()) {
1037 // stopped wells where cross flow is not allowed
1038 // are not added to the well container
1039 this->wellState().shutWell(w);
1040 this->well_close_times_.erase(well_name);
1041 this->well_open_times_.erase(well_name);
1042 continue;
1043 }
1044 // stopped wells are added to the container but marked as stopped
1045 this->wellState().stopWell(w);
1046 }
1047 }
1048
1049 // shut wells with zero rante constraints and disallowing
1050 if (!well_ecl.getAllowCrossFlow()) {
1051 const bool any_zero_rate_constraint = well_ecl.isProducer()
1052 ? well_ecl.productionControls(this->summaryState_).anyZeroRateConstraint()
1053 : well_ecl.injectionControls(this->summaryState_).anyZeroRateConstraint();
1054 if (any_zero_rate_constraint) {
1055 // Treat as shut, do not add to container.
1056 local_deferredLogger.debug(fmt::format(fmt::runtime(" Well {} gets shut due to having zero rate constraint and disallowing crossflow "), well_ecl.name()));
1057 this->wellState().shutWell(w);
1058 this->well_close_times_.erase(well_name);
1059 this->well_open_times_.erase(well_name);
1060 continue;
1061 }
1062 }
1063
1064 if (!wcycle.empty()) {
1065 const auto it = cycle_states.find(well_name);
1066 if (it != cycle_states.end()) {
1067 if (!it->second || well_status == Well::Status::SHUT) {
1068 // If well is shut in schedule we keep it shut
1069 if (well_status == Well::Status::SHUT) {
1070 this->well_open_times_.erase(well_name);
1071 this->well_close_times_.erase(well_name);
1072 }
1073 this->wellState().shutWell(w);
1074 continue;
1075 } else {
1076 this->wellState().openWell(w);
1077 }
1078 }
1079 }
1080
1081 // We dont add SHUT wells to the container
1082 if (ws.status == Well::Status::SHUT) {
1083 continue;
1084 }
1085
1086 well_container_.emplace_back(this->createWellPointer(w, report_step));
1087
1088 if (ws.status == Well::Status::STOP) {
1089 well_container_.back()->stopWell();
1090 this->well_close_times_.erase(well_name);
1091 this->well_open_times_.erase(well_name);
1092 }
1093 }
1094
1095 if (!wcycle.empty()) {
1096 const auto schedule_open =
1097 [&wg_events = this->report_step_start_events_](const std::string& name)
1098 {
1099 return wg_events.hasEvent(name, ScheduleEvents::REQUEST_OPEN_WELL);
1100 };
1101 for (const auto& [wname, wscale] : wcycle.efficiencyScale(this->simulator_.time(),
1102 this->simulator_.timeStepSize(),
1103 wmatcher,
1104 this->well_open_times_,
1105 schedule_open))
1106 {
1107 this->wellState().updateEfficiencyScalingFactor(wname, wscale);
1108 this->schedule_.add_event(ScheduleEvents::WELLGROUP_EFFICIENCY_UPDATE, report_step);
1109 }
1110 }
1111 }
1112
1113 this->well_container_generic_.clear();
1114 for (auto& w : well_container_) {
1115 this->well_container_generic_.push_back(w.get());
1116 }
1117
1118 this->network_.initialize(report_step);
1119
1120 this->wbp_.registerOpenWellsForWBPCalculation();
1121 }
1122
1123
1124
1125
1126
1127 template <typename TypeTag>
1130 createWellPointer(const int wellID, const int report_step) const
1131 {
1132 const auto is_multiseg = this->wells_ecl_[wellID].isMultiSegment();
1133
1134 if (! (this->param_.use_multisegment_well_ && is_multiseg)) {
1135 return this->template createTypedWellPointer<StandardWell<TypeTag>>(wellID, report_step);
1136 }
1137 else {
1138 return this->template createTypedWellPointer<MultisegmentWell<TypeTag>>(wellID, report_step);
1139 }
1140 }
1141
1142
1143
1144
1145
1146 template <typename TypeTag>
1147 template <typename WellType>
1148 std::unique_ptr<WellType>
1150 createTypedWellPointer(const int wellID, const int time_step) const
1151 {
1152 // Use the pvtRegionIdx from the top cell
1153 const auto& perf_data = this->well_perf_data_[wellID];
1154
1155 // Cater for case where local part might have no perforations.
1156 const auto pvtreg = perf_data.empty()
1157 ? 0 : this->pvt_region_idx_[perf_data.front().cell_index];
1158
1159 const auto& parallel_well_info = this->local_parallel_well_info_[wellID].get();
1160 const auto global_pvtreg = parallel_well_info.broadcastFirstPerforationValue(pvtreg);
1161
1162 return std::make_unique<WellType>(this->wells_ecl_[wellID],
1163 parallel_well_info,
1164 time_step,
1165 this->param_,
1166 *this->rateConverter_,
1167 global_pvtreg,
1168 this->numConservationQuantities(),
1169 this->numPhases(),
1170 wellID,
1171 perf_data);
1172 }
1173
1174
1175
1176
1177
1178 template<typename TypeTag>
1181 createWellForWellTest(const std::string& well_name,
1182 const int report_step,
1183 DeferredLogger& deferred_logger) const
1184 {
1185 // Finding the location of the well in wells_ecl
1186 const auto it =
1187 std::ranges::find_if(this->wells_ecl_,
1188 [&well_name](const auto& w)
1189 { return well_name == w.name(); });
1190 // It should be able to find in wells_ecl.
1191 if (it == this->wells_ecl_.end()) {
1192 OPM_DEFLOG_THROW(std::logic_error,
1193 fmt::format(fmt::runtime("Could not find well {} in wells_ecl"), well_name),
1194 deferred_logger);
1195 }
1196
1197 const int pos = static_cast<int>(std::distance(this->wells_ecl_.begin(), it));
1198 return this->createWellPointer(pos, report_step);
1199 }
1200
1201
1202
1203 template<typename TypeTag>
1204 void
1206 assemble(const double dt)
1207 {
1208 OPM_TIMEFUNCTION();
1209 auto logger_guard = this->groupStateHelper().pushLogger();
1210 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
1211
1212 const auto& iterCtx = simulator_.problem().iterationContext();
1213
1215 if (gaslift_.terminalOutput()) {
1216 const std::string msg =
1217 fmt::format(fmt::runtime("assemble() : iteration {}"), iterCtx.iteration());
1218 gaslift_.gliftDebug(msg, local_deferredLogger);
1219 }
1220 }
1221 last_report_ = SimulatorReportSingle();
1222 Dune::Timer perfTimer;
1223 perfTimer.start();
1224 this->closed_offending_wells_.clear();
1225
1226 {
1227 const int episodeIdx = simulator_.episodeIndex();
1228 const auto& network = this->schedule()[episodeIdx].network();
1229 if (!this->wellsActive() && !network.active()) {
1230 return;
1231 }
1232 }
1233
1234 // Timestep initialization: should run once at the start of each timestep.
1235 if (iterCtx.needsTimestepInit() && this->wellsActive()) {
1236 OPM_TIMEBLOCK(firstIterationAssemble);
1237 // try-catch is needed here as updateWellControls
1238 // contains global communication and has either to
1239 // be reached by all processes or all need to abort
1240 // before.
1242 {
1243 calculateExplicitQuantities();
1244 prepareTimeStep(local_deferredLogger);
1245 }
1246 OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger,
1247 "assemble() failed during well initialization: ",
1248 this->terminal_output_, grid().comm());
1249 }
1250
1251 const bool well_group_control_changed = updateWellControlsAndNetwork(false, dt, local_deferredLogger);
1252
1253 // even when there is no wells active, the network nodal pressure still need to be updated through updateWellControlsAndNetwork()
1254 // but there is no need to assemble the well equations
1255 if ( ! this->wellsActive() ) {
1256 return;
1257 }
1258
1259 assembleWellEqWithoutIteration(dt);
1260 // Pre-compute cell rates to we don't have to do this for every cell during linearization...
1261 updateCellRates();
1262
1263 // if group or well control changes we don't consider the
1264 // case converged
1265 last_report_.well_group_control_changed = well_group_control_changed;
1266 last_report_.assemble_time_well += perfTimer.stop();
1267 }
1268
1269
1270
1271
1272 template<typename TypeTag>
1273 bool
1275 updateWellControlsAndNetwork(const bool mandatory_network_balance,
1276 const double dt,
1277 DeferredLogger& local_deferredLogger)
1278 {
1279 OPM_TIMEFUNCTION();
1280 // not necessarily that we always need to update once of the network solutions
1281 bool do_network_update = true;
1282 bool well_group_control_changed = false;
1283 Scalar network_imbalance = 0.0;
1284 // after certain number of the iterations, we use relaxed tolerance for the network update
1285 const std::size_t iteration_to_relax = param_.network_max_strict_outer_iterations_;
1286 // after certain number of the iterations, we terminate
1287 const std::size_t max_iteration = param_.network_max_outer_iterations_;
1288 std::size_t network_update_iteration = 0;
1289 network_needs_more_balancing_force_another_newton_iteration_ = false;
1290 while (do_network_update) {
1291 if (network_update_iteration >= max_iteration ) {
1292 // only output to terminal if we at the last newton iterations where we try to balance the network.
1293 const int episodeIdx = simulator_.episodeIndex();
1294 const auto& iterCtx = simulator_.problem().iterationContext();
1295 if (this->network_.willBalanceOnNextIteration(episodeIdx, iterCtx)) {
1296 if (this->terminal_output_) {
1297 const std::string msg = fmt::format("Maximum of {:d} network iterations has been used and we stop the update, \n"
1298 "and try again after the next Newton iteration (imbalance = {:.2e} bar)",
1299 max_iteration, network_imbalance*1.0e-5);
1300 local_deferredLogger.debug(msg);
1301 }
1302 // To avoid stopping the newton iterations too early, before the network is converged,
1303 // we need to report it
1304 network_needs_more_balancing_force_another_newton_iteration_ = true;
1305 } else {
1306 if (this->terminal_output_) {
1307 const std::string msg = fmt::format("Maximum of {:d} network iterations has been used and we stop the update. \n"
1308 "The simulator will continue with unconverged network results (imbalance = {:.2e} bar)",
1309 max_iteration, network_imbalance*1.0e-5);
1310 local_deferredLogger.info(msg);
1311 }
1312 }
1313 break;
1314 }
1315 if (this->terminal_output_ && (network_update_iteration == iteration_to_relax) ) {
1316 local_deferredLogger.debug("We begin using relaxed tolerance for network update now after " + std::to_string(iteration_to_relax) + " iterations ");
1317 }
1318 const bool relax_network_balance = network_update_iteration >= iteration_to_relax;
1319 // Never optimize gas lift in last iteration, to allow network convergence (unless max_iter < 2)
1320 const bool optimize_gas_lift = ( (network_update_iteration + 1) < std::max(max_iteration, static_cast<std::size_t>(2)) );
1321 std::tie(well_group_control_changed, do_network_update, network_imbalance) =
1322 updateWellControlsAndNetworkIteration(mandatory_network_balance, relax_network_balance, optimize_gas_lift, dt,local_deferredLogger);
1323 ++network_update_iteration;
1324 }
1325 return well_group_control_changed;
1326 }
1327
1328
1329
1330
1331 template<typename TypeTag>
1332 std::tuple<bool, bool, typename BlackoilWellModel<TypeTag>::Scalar>
1334 updateWellControlsAndNetworkIteration(const bool mandatory_network_balance,
1335 const bool relax_network_tolerance,
1336 const bool optimize_gas_lift,
1337 const double dt,
1338 DeferredLogger& local_deferredLogger)
1339 {
1340 OPM_TIMEFUNCTION();
1341 const auto& iterCtx = simulator_.problem().iterationContext();
1342 const int reportStepIdx = simulator_.episodeIndex();
1343 this->updateAndCommunicateGroupData(reportStepIdx, iterCtx,
1344 param_.nupcol_group_rate_tolerance_, /*update_wellgrouptarget*/ true);
1345 // We need to call updateWellControls before we update the network as
1346 // network updates are only done on thp controlled wells.
1347 // Note that well controls are allowed to change during updateNetwork
1348 // and in prepareWellsBeforeAssembling during well solves.
1349 bool well_group_control_changed = updateWellControls(local_deferredLogger);
1350 const auto [more_inner_network_update, network_imbalance] =
1351 this->network_.update(mandatory_network_balance,
1352 local_deferredLogger,
1353 relax_network_tolerance);
1354
1355 bool alq_updated = false;
1357 {
1358 if (optimize_gas_lift) {
1359 // we need to update the potentials if the thp limit as been modified by
1360 // the network balancing
1361 const bool updatePotentials = (this->network_.shouldBalance(reportStepIdx, iterCtx) ||
1362 mandatory_network_balance);
1363 alq_updated = gaslift_.maybeDoGasLiftOptimize(simulator_,
1364 well_container_,
1365 this->network_.nodePressures(),
1366 updatePotentials,
1367 this->wellState(),
1368 this->groupState(),
1369 local_deferredLogger);
1370 }
1371 prepareWellsBeforeAssembling(dt);
1372 }
1373 OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger,
1374 "updateWellControlsAndNetworkIteration() failed: ",
1375 this->terminal_output_, grid().comm());
1376
1377 // update guide rates
1378 if (alq_updated || BlackoilWellModelGuideRates(*this).
1379 guideRateUpdateIsNeeded(reportStepIdx)) {
1380 const double simulationTime = simulator_.time();
1381 // NOTE: For reservoir coupling: Slave group potentials are only communicated
1382 // at the start of the time step, see beginTimeStep(). Here, we assume those
1383 // potentials remain unchanged during the time step when updating guide rates below.
1384 this->guide_rate_handler_.updateGuideRates(
1385 reportStepIdx, simulationTime, this->wellState(), this->groupState()
1386 );
1387 }
1388 // we need to re-iterate the network when the well group controls changed or gaslift/alq is changed or
1389 // the inner iterations are did not converge
1390 const bool more_network_update = this->network_.shouldBalance(reportStepIdx, iterCtx) &&
1391 (more_inner_network_update || alq_updated);
1392 return {well_group_control_changed, more_network_update, network_imbalance};
1393 }
1394
1395 template<typename TypeTag>
1396 void
1398 assembleWellEq(const double dt)
1399 {
1400 OPM_TIMEFUNCTION();
1401 for (auto& well : well_container_) {
1402 well->assembleWellEq(simulator_, dt, this->groupStateHelper(), this->wellState());
1403 }
1404 }
1405
1406
1407 template<typename TypeTag>
1408 void
1410 prepareWellsBeforeAssembling(const double dt)
1411 {
1412 OPM_TIMEFUNCTION();
1413 for (auto& well : well_container_) {
1414 well->prepareWellBeforeAssembling(
1415 simulator_, dt, this->groupStateHelper(), this->wellState()
1416 );
1417 }
1418 }
1419
1420
1421 template<typename TypeTag>
1422 void
1424 assembleWellEqWithoutIteration(const double dt)
1425 {
1426 OPM_TIMEFUNCTION();
1427 auto& deferred_logger = this->groupStateHelper().deferredLogger();
1428 // We make sure that all processes throw in case there is an exception
1429 // on one of them (WetGasPvt::saturationPressure might throw if not converged)
1431
1432 for (auto& well: well_container_) {
1433 well->assembleWellEqWithoutIteration(simulator_, this->groupStateHelper(), dt, this->wellState(),
1434 /*solving_with_zero_rate=*/false);
1435 }
1436 OPM_END_PARALLEL_TRY_CATCH_LOG(deferred_logger, "BlackoilWellModel::assembleWellEqWithoutIteration failed: ",
1437 this->terminal_output_, grid().comm());
1438
1439 }
1440
1441 template<typename TypeTag>
1442 void
1445 {
1446 // Pre-compute cell rates for all wells
1447 cellRates_.clear();
1448 for (const auto& well : well_container_) {
1449 well->addCellRates(cellRates_);
1450 }
1451 }
1452
1453 template<typename TypeTag>
1454 void
1456 updateCellRatesForDomain(int domainIndex, const std::map<std::string, int>& well_domain_map)
1457 {
1458 // Pre-compute cell rates only for wells in the specified domain
1459 cellRates_.clear();
1460 for (const auto& well : well_container_) {
1461 const auto it = well_domain_map.find(well->name());
1462 if (it != well_domain_map.end() && it->second == domainIndex) {
1463 well->addCellRates(cellRates_);
1464 }
1465 }
1466 }
1467
1468#if COMPILE_GPU_BRIDGE
1469 template<typename TypeTag>
1470 void
1473 {
1474 // prepare for StandardWells
1476
1477 for(unsigned int i = 0; i < well_container_.size(); i++){
1478 auto& well = well_container_[i];
1479 auto derived = dynamic_cast<StandardWell<TypeTag>*>(well.get());
1480 if (derived) {
1481 wellContribs.addNumBlocks(derived->linSys().getNumBlocks());
1482 }
1483 }
1484
1485 // allocate memory for data from StandardWells
1486 wellContribs.alloc();
1487
1488 for(unsigned int i = 0; i < well_container_.size(); i++){
1489 auto& well = well_container_[i];
1490 // maybe WellInterface could implement addWellContribution()
1491 auto derived_std = dynamic_cast<StandardWell<TypeTag>*>(well.get());
1492 if (derived_std) {
1493 derived_std->linSys().extract(derived_std->numStaticWellEq, wellContribs);
1494 } else {
1495 auto derived_ms = dynamic_cast<MultisegmentWell<TypeTag>*>(well.get());
1496 if (derived_ms) {
1497 derived_ms->linSys().extract(wellContribs);
1498 } else {
1499 OpmLog::warning("Warning unknown type of well");
1500 }
1501 }
1502 }
1503 }
1504#endif
1505
1506 template<typename TypeTag>
1507 void
1509 addWellContributions(SparseMatrixAdapter& jacobian) const
1510 {
1511 for ( const auto& well: well_container_ ) {
1512 well->addWellContributions(jacobian);
1513 }
1514 }
1515
1516 template<typename TypeTag>
1517 void
1520 const BVector& weights,
1521 const bool use_well_weights) const
1522 {
1523 int nw = this->numLocalWellsEnd();
1524 int rdofs = local_num_cells_;
1525 for ( int i = 0; i < nw; i++ ) {
1526 int wdof = rdofs + i;
1527 jacobian[wdof][wdof] = 1.0;// better scaling ?
1528 }
1529
1530 for (const auto& well : well_container_) {
1531 well->addWellPressureEquations(jacobian,
1532 weights,
1533 pressureVarIndex,
1534 use_well_weights,
1535 this->wellState());
1536 }
1537 }
1538
1539 template <typename TypeTag>
1541 addReservoirSourceTerms(GlobalEqVector& residual,
1542 const std::vector<typename SparseMatrixAdapter::MatrixBlock*>& diagMatAddress) const
1543 {
1544 // NB this loop may write multiple times to the same element
1545 // if a cell is perforated by more than one well, so it should
1546 // not be OpenMP-parallelized.
1547 for (const auto& well : well_container_) {
1548 if (!well->isOperableAndSolvable() && !well->wellIsStopped()) {
1549 continue;
1550 }
1551 const auto& cells = well->cells();
1552 const auto& rates = well->connectionRates();
1553 for (unsigned perfIdx = 0; perfIdx < rates.size(); ++perfIdx) {
1554 unsigned cellIdx = cells[perfIdx];
1555 auto rate = rates[perfIdx];
1556 rate *= -1.0;
1557 VectorBlockType res(0.0);
1558 using MatrixBlockType = typename SparseMatrixAdapter::MatrixBlock;
1559 MatrixBlockType bMat(0.0);
1560 simulator_.model().linearizer().setResAndJacobi(res, bMat, rate);
1561 residual[cellIdx] += res;
1562 *diagMatAddress[cellIdx] += bMat;
1563 }
1564 }
1565 }
1566
1567
1568 template<typename TypeTag>
1569 void
1572 {
1573 int nw = this->numLocalWellsEnd();
1574 int rdofs = local_num_cells_;
1575 for (int i = 0; i < nw; ++i) {
1576 int wdof = rdofs + i;
1577 jacobian.entry(wdof,wdof) = 1.0;// better scaling ?
1578 }
1579 const auto wellconnections = this->getMaxWellConnections();
1580 for (int i = 0; i < nw; ++i) {
1581 const auto& perfcells = wellconnections[i];
1582 for (int perfcell : perfcells) {
1583 int wdof = rdofs + i;
1584 jacobian.entry(wdof, perfcell) = 0.0;
1585 jacobian.entry(perfcell, wdof) = 0.0;
1586 }
1587 }
1588 }
1589
1590
1591 template<typename TypeTag>
1592 void
1595 {
1596 auto loggerGuard = this->groupStateHelper().pushLogger();
1598 {
1599 for (const auto& well : well_container_) {
1600 const auto& cells = well->cells();
1601 x_local_.resize(cells.size());
1602
1603 for (size_t i = 0; i < cells.size(); ++i) {
1604 x_local_[i] = x[cells[i]];
1605 }
1606 well->recoverWellSolutionAndUpdateWellState(simulator_, x_local_,
1607 this->groupStateHelper(), this->wellState());
1608 }
1609 }
1610 OPM_END_PARALLEL_TRY_CATCH("recoverWellSolutionAndUpdateWellState() failed: ",
1611 simulator_.vanguard().grid().comm());
1612 }
1613
1614
1615 template<typename TypeTag>
1616 void
1618 recoverWellSolutionAndUpdateWellStateDomain(const BVector& x, const int domainIdx)
1619 {
1620 if (!nldd_) {
1621 OPM_THROW(std::logic_error, "Attempt to call NLDD method without a NLDD solver");
1622 }
1623
1624 return nldd_->recoverWellSolutionAndUpdateWellState(x, domainIdx);
1625 }
1626
1627
1628 template<typename TypeTag>
1631 getWellConvergence(const std::vector<Scalar>& B_avg, bool checkWellGroupControlsAndNetwork) const
1632 {
1633 // Get global (from all processes) convergence report.
1634 ConvergenceReport local_report;
1635 const auto& iterCtx = simulator_.problem().iterationContext();
1636 const bool relaxTolerance = iterCtx.shouldRelax(param_.strict_outer_iter_wells_ + 1);
1637 {
1638 auto logger_guard = this->groupStateHelper().pushLogger();
1639 for (const auto& well : well_container_) {
1640 if (well->isOperableAndSolvable() || well->wellIsStopped()) {
1641 local_report += well->getWellConvergence(
1642 this->groupStateHelper(), B_avg,
1643 relaxTolerance);
1644 } else {
1645 ConvergenceReport report;
1646 using CR = ConvergenceReport;
1647 report.setWellFailed({CR::WellFailure::Type::Unsolvable, CR::Severity::Normal, -1, well->name()});
1648 local_report += report;
1649 }
1650 }
1651 } // logger_guard goes out of scope here, before the OpmLog::debug() calls below
1652
1653 const Opm::Parallel::Communication comm = grid().comm();
1654 ConvergenceReport report = gatherConvergenceReport(local_report, comm);
1655
1656 if (checkWellGroupControlsAndNetwork) {
1657 // the well_group_control_changed info is already communicated
1658 report.setWellGroupTargetsViolated(this->lastReport().well_group_control_changed);
1659 report.setNetworkNotYetBalancedForceAnotherNewtonIteration(network_needs_more_balancing_force_another_newton_iteration_);
1660 }
1661
1662 if (this->terminal_output_) {
1663 // Log debug messages for NaN or too large residuals.
1664 for (const auto& f : report.wellFailures()) {
1665 if (f.severity() == ConvergenceReport::Severity::NotANumber) {
1666 OpmLog::debug("NaN residual found with phase " + std::to_string(f.phase()) + " for well " + f.wellName());
1667 } else if (f.severity() == ConvergenceReport::Severity::TooLarge) {
1668 OpmLog::debug("Too large residual found with phase " + std::to_string(f.phase()) + " for well " + f.wellName());
1669 }
1670 }
1671 }
1672 return report;
1673 }
1674
1675
1676
1677
1678
1679 template<typename TypeTag>
1680 void
1683 {
1684 // TODO: checking isOperableAndSolvable() ?
1685 for (auto& well : well_container_) {
1686 well->calculateExplicitQuantities(simulator_, this->groupStateHelper());
1687 }
1688 }
1689
1690
1691
1692
1693
1694 template<typename TypeTag>
1695 bool
1697 updateWellControls(DeferredLogger& deferred_logger)
1698 {
1699 OPM_TIMEFUNCTION();
1700 if (!this->wellsActive()) {
1701 return false;
1702 }
1703 const int episodeIdx = simulator_.episodeIndex();
1704 const auto& comm = simulator_.vanguard().grid().comm();
1705 size_t iter = 0;
1706 bool changed_well_group = false;
1707 const Group& fieldGroup = this->schedule().getGroup("FIELD", episodeIdx);
1708 // Check group individual constraints.
1709 // iterate a few times to make sure all constraints are honored
1710 const std::size_t max_iter = param_.well_group_constraints_max_iterations_;
1711 while(!changed_well_group && iter < max_iter) {
1712 changed_well_group = updateGroupControls(fieldGroup, deferred_logger, episodeIdx);
1713
1714 // Check wells' group constraints and communicate.
1715 bool changed_well_to_group = false;
1716 {
1717 OPM_TIMEBLOCK(UpdateWellControls);
1718 // For MS Wells a linear solve is performed below and the matrix might be singular.
1719 // We need to communicate the exception thrown to the others and rethrow.
1721 for (const auto& well : well_container_) {
1723 const bool changed_well = well->updateWellControl(
1724 simulator_, mode, this->groupStateHelper(), this->wellState()
1725 );
1726 if (changed_well) {
1727 changed_well_to_group = changed_well || changed_well_to_group;
1728 }
1729 }
1730 OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel: updating well controls failed: ",
1731 simulator_.gridView().comm());
1732 }
1733
1734 changed_well_to_group = comm.sum(static_cast<int>(changed_well_to_group));
1735 if (changed_well_to_group) {
1736 updateAndCommunicate(episodeIdx);
1737 changed_well_group = true;
1738 }
1739
1740 // Check individual well constraints and communicate.
1741 bool changed_well_individual = false;
1742 {
1743 // For MS Wells a linear solve is performed below and the matrix might be singular.
1744 // We need to communicate the exception thrown to the others and rethrow.
1746 for (const auto& well : well_container_) {
1748 const bool changed_well = well->updateWellControl(
1749 simulator_, mode, this->groupStateHelper(), this->wellState()
1750 );
1751 if (changed_well) {
1752 changed_well_individual = changed_well || changed_well_individual;
1753 }
1754 }
1755 OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel: updating well controls failed: ",
1756 simulator_.gridView().comm());
1757 }
1758
1759 changed_well_individual = comm.sum(static_cast<int>(changed_well_individual));
1760 if (changed_well_individual) {
1761 updateAndCommunicate(episodeIdx);
1762 changed_well_group = true;
1763 }
1764 iter++;
1765 }
1766
1767 // update wsolvent fraction for REIN wells
1768 this->updateWsolvent(fieldGroup, episodeIdx, this->nupcolWellState());
1769
1770 return changed_well_group;
1771 }
1772
1773
1774 template<typename TypeTag>
1775 void
1777 updateAndCommunicate(const int reportStepIdx)
1778 {
1779 const auto& iterCtx = simulator_.problem().iterationContext();
1780 this->updateAndCommunicateGroupData(reportStepIdx,
1781 iterCtx,
1782 param_.nupcol_group_rate_tolerance_,
1783 /*update_wellgrouptarget*/ true);
1784
1785 // updateWellStateWithTarget might throw for multisegment wells hence we
1786 // have a parallel try catch here to thrown on all processes.
1788 // if a well or group change control it affects all wells that are under the same group
1789 for (const auto& well : well_container_) {
1790 // We only want to update wells under group-control here
1791 const auto& ws = this->wellState().well(well->indexOfWell());
1792 if (ws.production_cmode == Well::ProducerCMode::GRUP ||
1793 ws.injection_cmode == Well::InjectorCMode::GRUP)
1794 {
1795 well->updateWellStateWithTarget(
1796 simulator_, this->groupStateHelper(), this->wellState()
1797 );
1798 }
1799 }
1800 OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel::updateAndCommunicate failed: ",
1801 simulator_.gridView().comm())
1802 this->updateAndCommunicateGroupData(reportStepIdx,
1803 iterCtx,
1804 param_.nupcol_group_rate_tolerance_,
1805 /*update_wellgrouptarget*/ true);
1806 }
1807
1808 template<typename TypeTag>
1809 bool
1811 updateGroupControls(const Group& group,
1812 DeferredLogger& deferred_logger,
1813 const int reportStepIdx)
1814 {
1815 OPM_TIMEFUNCTION();
1816 const auto& iterCtx = simulator_.problem().iterationContext();
1817 bool changed = false;
1818 // restrict the number of group switches but only after nupcol iterations.
1819 const int nupcol = this->schedule()[reportStepIdx].nupcol();
1820 const int max_number_of_group_switches = param_.max_number_of_group_switches_;
1821 const bool update_group_switching_log = !iterCtx.withinNupcol(nupcol);
1822 const bool changed_hc = this->checkGroupHigherConstraints(group, deferred_logger, reportStepIdx, max_number_of_group_switches, update_group_switching_log);
1823 if (changed_hc) {
1824 changed = true;
1825 updateAndCommunicate(reportStepIdx);
1826 }
1827
1828 bool changed_individual =
1830 updateGroupIndividualControl(group,
1831 reportStepIdx,
1832 max_number_of_group_switches,
1833 update_group_switching_log,
1834 this->switched_inj_groups_,
1835 this->switched_prod_groups_,
1836 this->closed_offending_wells_,
1837 this->groupState(),
1838 this->wellState(),
1839 deferred_logger);
1840
1841 if (changed_individual) {
1842 changed = true;
1843 updateAndCommunicate(reportStepIdx);
1844 }
1845 // call recursively down the group hierarchy
1846 for (const std::string& groupName : group.groups()) {
1847 bool changed_this = updateGroupControls(this->schedule().getGroup(groupName, reportStepIdx), deferred_logger, reportStepIdx);
1848 changed = changed || changed_this;
1849 }
1850 return changed;
1851 }
1852
1853 template<typename TypeTag>
1854 void
1856 updateWellTestState(const double simulationTime, WellTestState& wellTestState)
1857 {
1858 OPM_TIMEFUNCTION();
1859 auto logger_guard = this->groupStateHelper().pushLogger();
1860 auto& local_deferredLogger = this->groupStateHelper().deferredLogger();
1861 for (const auto& well : well_container_) {
1862 const auto& wname = well->name();
1863 const auto wasClosed = wellTestState.well_is_closed(wname);
1864 well->checkWellOperability(simulator_,
1865 this->wellState(),
1866 this->groupStateHelper());
1867 const bool under_zero_target =
1868 well->wellUnderZeroGroupRateTarget(this->groupStateHelper());
1869 well->updateWellTestState(this->wellState().well(wname),
1870 simulationTime,
1871 /*writeMessageToOPMLog=*/ true,
1872 under_zero_target,
1873 wellTestState,
1874 local_deferredLogger);
1875
1876 if (!wasClosed && wellTestState.well_is_closed(wname)) {
1877 this->closed_this_step_.insert(wname);
1878
1879 // maybe open a new well
1880 const WellEconProductionLimits& econ_production_limits = well->wellEcl().getEconLimits();
1881 if (econ_production_limits.validFollowonWell()) {
1882 const auto episode_idx = simulator_.episodeIndex();
1883 const auto follow_on_well = econ_production_limits.followonWell();
1884 if (!this->schedule().hasWell(follow_on_well, episode_idx)) {
1885 const auto msg = fmt::format("Well {} was closed. But the given follow on well {} does not exist."
1886 "The simulator continues without opening a follow on well.",
1887 wname, follow_on_well);
1888 local_deferredLogger.warning(msg);
1889 }
1890 auto& ws = this->wellState().well(follow_on_well);
1891 const bool success = ws.updateStatus(WellStatus::OPEN);
1892 if (success) {
1893 const auto msg = fmt::format("Well {} was closed. The follow on well {} opens instead.", wname, follow_on_well);
1894 local_deferredLogger.info(msg);
1895 } else {
1896 const auto msg = fmt::format("Well {} was closed. The follow on well {} is already open.", wname, follow_on_well);
1897 local_deferredLogger.warning(msg);
1898 }
1899 }
1900
1901 }
1902 }
1903
1904 for (const auto& [group_name, to] : this->closed_offending_wells_) {
1905 if (this->hasOpenLocalWell(to.second) &&
1906 !this->wasDynamicallyShutThisTimeStep(to.second))
1907 {
1908 wellTestState.close_well(to.second,
1909 WellTestConfig::Reason::GROUP,
1910 simulationTime);
1911 this->updateClosedWellsThisStep(to.second);
1912 const std::string msg =
1913 fmt::format("Procedure on exceeding {} limit is WELL for group {}. "
1914 "Well {} is {}.",
1915 to.first,
1916 group_name,
1917 to.second,
1918 "shut");
1919 local_deferredLogger.info(msg);
1920 }
1921 }
1922 }
1923
1924
1925 template<typename TypeTag>
1926 void
1928 const WellState<Scalar, IndexTraits>& well_state_copy,
1929 std::string& exc_msg,
1930 ExceptionType::ExcEnum& exc_type)
1931 {
1932 OPM_TIMEFUNCTION();
1933 const int np = this->numPhases();
1934 std::vector<Scalar> potentials;
1935 const auto& well = well_container_[widx];
1936 std::string cur_exc_msg;
1937 auto cur_exc_type = ExceptionType::NONE;
1938 try {
1939 well->computeWellPotentials(simulator_, well_state_copy, this->groupStateHelper(), potentials);
1940 }
1941 // catch all possible exception and store type and message.
1942 OPM_PARALLEL_CATCH_CLAUSE(cur_exc_type, cur_exc_msg);
1943 if (cur_exc_type != ExceptionType::NONE) {
1944 exc_msg += fmt::format("\nFor well {}: {}", well->name(), cur_exc_msg);
1945 }
1946 exc_type = std::max(exc_type, cur_exc_type);
1947 // Store it in the well state
1948 // potentials is resized and set to zero in the beginning of well->ComputeWellPotentials
1949 // and updated only if sucessfull. i.e. the potentials are zero for exceptions
1950 auto& ws = this->wellState().well(well->indexOfWell());
1951 for (int p = 0; p < np; ++p) {
1952 // make sure the potentials are positive
1953 ws.well_potentials[p] = std::max(Scalar{0.0}, potentials[p]);
1954 }
1955 }
1956
1957
1958
1959 template <typename TypeTag>
1960 void
1963 {
1964 for (const auto& wellPtr : this->well_container_) {
1965 this->calculateProductivityIndexValues(wellPtr.get(), deferred_logger);
1966 }
1967 }
1968
1969
1970
1971
1972
1973 template <typename TypeTag>
1974 void
1976 calculateProductivityIndexValuesShutWells(const int reportStepIdx,
1977 DeferredLogger& deferred_logger)
1978 {
1979 // For the purpose of computing PI/II values, it is sufficient to
1980 // construct StandardWell instances only. We don't need to form
1981 // well objects that honour the 'isMultisegment()' flag of the
1982 // corresponding "this->wells_ecl_[shutWell]".
1983
1984 for (const auto& shutWell : this->local_shut_wells_) {
1985 if (!this->wells_ecl_[shutWell].hasConnections()) {
1986 // No connections in this well. Nothing to do.
1987 continue;
1988 }
1989
1990 auto wellPtr = this->template createTypedWellPointer
1991 <StandardWell<TypeTag>>(shutWell, reportStepIdx);
1992
1993 wellPtr->init(this->depth_, this->gravity_, this->B_avg_, true);
1994
1995 this->calculateProductivityIndexValues(wellPtr.get(), deferred_logger);
1996 }
1997 }
1998
1999
2000
2001
2002
2003 template <typename TypeTag>
2004 void
2007 DeferredLogger& deferred_logger)
2008 {
2009 wellPtr->updateProductivityIndex(this->simulator_,
2010 this->prod_index_calc_[wellPtr->indexOfWell()],
2011 this->wellState(),
2012 deferred_logger);
2013 }
2014
2015
2016
2017 template<typename TypeTag>
2018 void
2020 prepareTimeStep(DeferredLogger& deferred_logger)
2021 {
2022 // Check if there is a network with active prediction wells at this time step.
2023 const auto episodeIdx = simulator_.episodeIndex();
2024 this->network_.updateActiveState(episodeIdx);
2025
2026 // Rebalance the network initially if any wells in the network have status changes
2027 // (Need to check this before clearing events)
2028 const bool do_prestep_network_rebalance =
2029 param_.pre_solve_network_ && this->network_.needPreStepRebalance(episodeIdx);
2030
2031 for (const auto& well : well_container_) {
2032 auto& events = this->wellState().well(well->indexOfWell()).events;
2033 if (events.hasEvent(WellState<Scalar, IndexTraits>::event_mask)) {
2034 well->updateWellStateWithTarget(
2035 simulator_, this->groupStateHelper(), this->wellState()
2036 );
2037 well->updatePrimaryVariables(this->groupStateHelper());
2038 // There is no new well control change input within a report step,
2039 // so next time step, the well does not consider to have effective events anymore.
2041 }
2042 // these events only work for the first time step within the report step
2043 if (events.hasEvent(ScheduleEvents::REQUEST_OPEN_WELL)) {
2044 events.clearEvent(ScheduleEvents::REQUEST_OPEN_WELL);
2045 }
2046 // solve the well equation initially to improve the initial solution of the well model
2047 if (param_.solve_welleq_initially_ && well->isOperableAndSolvable()) {
2048 try {
2049 well->solveWellEquation(
2050 simulator_, this->groupStateHelper(), this->wellState()
2051 );
2052 } catch (const std::exception& e) {
2053 const std::string msg = "Compute initial well solution for " + well->name() + " initially failed. Continue with the previous rates";
2054 deferred_logger.warning("WELL_INITIAL_SOLVE_FAILED", msg);
2055 }
2056 }
2057 // If we're using local well solves that include control switches, they also update
2058 // operability, so reset before main iterations begin
2059 well->resetWellOperability();
2060 }
2061 updatePrimaryVariables();
2062
2063 // Actually do the pre-step network rebalance, using the updated well states and initial solutions
2064 if (do_prestep_network_rebalance) {
2065 network_.doPreStepRebalance(deferred_logger);
2066 }
2067 }
2068
2069 template<typename TypeTag>
2070 void
2073 {
2074 std::vector< Scalar > B_avg(numConservationQuantities(), Scalar() );
2075 const auto& grid = simulator_.vanguard().grid();
2076 const auto& gridView = grid.leafGridView();
2077 ElementContext elemCtx(simulator_);
2078
2080 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
2081 elemCtx.updatePrimaryStencil(elem);
2082 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
2083
2084 const auto& intQuants = elemCtx.intensiveQuantities(/*spaceIdx=*/0, /*timeIdx=*/0);
2085 const auto& fs = intQuants.fluidState();
2086
2087 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
2088 {
2089 if (!FluidSystem::phaseIsActive(phaseIdx)) {
2090 continue;
2091 }
2092
2093 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
2094 auto& B = B_avg[ compIdx ];
2095
2096 B += 1 / fs.invB(phaseIdx).value();
2097 }
2098 if constexpr (has_solvent_) {
2099 auto& B = B_avg[solventSaturationIdx];
2100 B += 1 / intQuants.solventInverseFormationVolumeFactor().value();
2101 }
2102 }
2103 OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel::updateAverageFormationFactor() failed: ", grid.comm())
2104
2105 // compute global average
2106 grid.comm().sum(B_avg.data(), B_avg.size());
2107 B_avg_.resize(B_avg.size());
2108 std::ranges::transform(B_avg, B_avg_.begin(),
2109 [gcells = global_num_cells_](const auto bval)
2110 { return bval / gcells; });
2111 }
2112
2113
2114
2115
2116
2117 template<typename TypeTag>
2118 void
2121 {
2122 for (const auto& well : well_container_) {
2123 well->updatePrimaryVariables(this->groupStateHelper());
2124 }
2125 }
2126
2127 template<typename TypeTag>
2128 void
2130 {
2131 const auto& grid = simulator_.vanguard().grid();
2132 const auto& eclProblem = simulator_.problem();
2133 const unsigned numCells = grid.size(/*codim=*/0);
2134
2135 this->pvt_region_idx_.resize(numCells);
2136 for (unsigned cellIdx = 0; cellIdx < numCells; ++cellIdx) {
2137 this->pvt_region_idx_[cellIdx] =
2138 eclProblem.pvtRegionIndex(cellIdx);
2139 }
2140 }
2141
2142 // The number of components in the model.
2143 template<typename TypeTag>
2144 int
2146 {
2147 // The numPhases() functions returns 1-3, depending on which
2148 // of the (oil, water, gas) phases are active. For each of those phases,
2149 // if the phase is active the corresponding component is present and
2150 // conserved.
2151 // Apart from (oil, water, gas), in the current well model only solvent
2152 // is explicitly modelled as a conserved quantity (polymer, energy, salt
2153 // etc. are not), unlike the reservoir part where all such quantities are
2154 // conserved. This function must therefore be updated when/if we add
2155 // more conserved quantities in the well model.
2156 return this->numPhases() + has_solvent_;
2157 }
2158
2159 template<typename TypeTag>
2160 void
2162 {
2163 const auto& eclProblem = simulator_.problem();
2164 depth_.resize(local_num_cells_);
2165 for (unsigned cellIdx = 0; cellIdx < local_num_cells_; ++cellIdx) {
2166 depth_[cellIdx] = eclProblem.dofCenterDepth(cellIdx);
2167 }
2168 }
2169
2170 template<typename TypeTag>
2173 getWell(const std::string& well_name) const
2174 {
2175 // finding the iterator of the well in wells_ecl
2176 const auto well =
2177 std::ranges::find_if(well_container_,
2178 [&well_name](const WellInterfacePtr& elem) -> bool
2179 { return elem->name() == well_name; });
2180
2181 assert(well != well_container_.end());
2182
2183 return **well;
2184 }
2185
2186 template <typename TypeTag>
2187 int
2189 reportStepIndex() const
2190 {
2191 return std::max(this->simulator_.episodeIndex(), 0);
2192 }
2193
2194
2195
2196
2197
2198 template<typename TypeTag>
2199 void
2201 calcResvCoeff(const int fipnum,
2202 const int pvtreg,
2203 const std::vector<Scalar>& production_rates,
2204 std::vector<Scalar>& resv_coeff) const
2205 {
2206 rateConverter_->calcCoeff(fipnum, pvtreg, production_rates, resv_coeff);
2207 }
2208
2209 template<typename TypeTag>
2210 void
2212 calcInjResvCoeff(const int fipnum,
2213 const int pvtreg,
2214 std::vector<Scalar>& resv_coeff) const
2215 {
2216 rateConverter_->calcInjCoeff(fipnum, pvtreg, resv_coeff);
2217 }
2218
2219
2220 template <typename TypeTag>
2221 void
2224 {
2225 if constexpr (energyModuleType_ == EnergyModules::FullyImplicitThermal) {
2226 const int np = this->numPhases();
2227 const int nw = this->numLocalWells();
2228 for (auto wellID = 0*nw; wellID < nw; ++wellID) {
2229 const Well& well = this->wells_ecl_[wellID];
2230 auto& ws = this->wellState().well(wellID);
2231 if (well.isInjector()) {
2232 if (ws.status != WellStatus::STOP) {
2233 this->wellState().well(wellID).temperature = well.inj_temperature();
2234 continue;
2235 }
2236 }
2237 std::array<Scalar,2> weighted{0.0,0.0};
2238 auto& [weighted_temperature, total_weight] = weighted;
2239 const auto& well_info = this->local_parallel_well_info_[wellID].get();
2240 using int_type = decltype(this->well_perf_data_[wellID].size());
2241 for (int_type perf = 0, end_perf = this->well_perf_data_[wellID].size(); perf < end_perf; ++perf) {
2242 const int cell_idx = this->well_perf_data_[wellID][perf].cell_index;
2243 const auto& intQuants = simulator_.model().intensiveQuantities(cell_idx, /*timeIdx=*/0);
2244 const auto& fs = intQuants.fluidState();
2245 Scalar weight_factor = computeTemperatureWeightFactor(perf, np, fs, ws);
2246 total_weight += weight_factor;
2247 weighted_temperature += weight_factor * fs.temperature(/*phaseIdx*/0).value();
2248 }
2249 well_info.communication().sum(weighted.data(), 2);
2250 this->wellState().well(wellID).temperature = weighted_temperature / total_weight;
2251 }
2252 }
2253 }
2254
2255
2256 template <typename TypeTag>
2258 assignWellTracerRates(data::Wells& wsrpt) const
2259 {
2260 const auto reportStepIdx = static_cast<unsigned int>(this->reportStepIndex());
2261 const auto& trMod = this->simulator_.problem().tracerModel();
2262
2263 BlackoilWellModelGeneric<Scalar, IndexTraits>::assignWellTracerRates(wsrpt, trMod.getWellTracerRates(), reportStepIdx);
2264 BlackoilWellModelGeneric<Scalar, IndexTraits>::assignWellTracerRates(wsrpt, trMod.getWellFreeTracerRates(), reportStepIdx);
2265 BlackoilWellModelGeneric<Scalar, IndexTraits>::assignWellTracerRates(wsrpt, trMod.getWellSolTracerRates(), reportStepIdx);
2266
2267 this->assignMswTracerRates(wsrpt, trMod.getMswTracerRates(), reportStepIdx);
2268 }
2269
2270 template <typename TypeTag>
2271 void BlackoilWellModel<TypeTag>::
2272 assignWellSpeciesRates(data::Wells& wsrpt) const
2273 {
2274 const auto reportStepIdx = static_cast<unsigned int>(this->reportStepIndex());
2275 const auto& geochemMod = this->simulator_.problem().geochemistryModel();
2276
2277 BlackoilWellModelGeneric<Scalar, IndexTraits>::assignWellTracerRates(wsrpt, geochemMod.getWellSpeciesRates(), reportStepIdx);
2278
2279 this->assignMswTracerRates(wsrpt, geochemMod.getMswSpeciesRates(), reportStepIdx);
2280 }
2281
2282 template <typename TypeTag>
2283 [[nodiscard]] auto BlackoilWellModel<TypeTag>::rsConstInfo() const
2284 -> typename WellState<Scalar,IndexTraits>::RsConstInfo
2285 {
2286 if (! FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) ||
2287 ! FluidSystem::enableConstantRs())
2288 {
2289 return {};
2290 }
2291
2292 const auto& rsConstTables = this->eclState_
2293 .getTableManager().getRsconstTables();
2294
2295 if (rsConstTables.empty() ||
2296 (rsConstTables[0].numRows() != std::size_t{1}))
2297 {
2298 return {};
2299 }
2300
2301 const auto rsConst = rsConstTables[0].getColumn(0).front();
2302
2303 return { true, static_cast<Scalar>(rsConst) };
2304 }
2305
2306} // namespace Opm
2307
2308#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:2201
void prepareTimeStep(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:2020
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:1334
WellInterfacePtr createWellPointer(const int wellID, const int report_step) const
Definition: BlackoilWellModel_impl.hpp:1130
void prepareWellsBeforeAssembling(const double dt)
Definition: BlackoilWellModel_impl.hpp:1410
void init()
Definition: BlackoilWellModel_impl.hpp:166
const Simulator & simulator() const
Definition: BlackoilWellModel.hpp:370
std::vector< Scalar > depth_
Definition: BlackoilWellModel.hpp:517
std::size_t global_num_cells_
Definition: BlackoilWellModel.hpp:513
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:2173
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:2145
bool updateWellControls(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1697
int reportStepIndex() const
Definition: BlackoilWellModel_impl.hpp:2189
void calculateProductivityIndexValues(DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1962
void extractLegacyDepth_()
Definition: BlackoilWellModel_impl.hpp:2161
void extractLegacyCellPvtRegionIndex_()
Definition: BlackoilWellModel_impl.hpp:2129
void recoverWellSolutionAndUpdateWellStateDomain(const BVector &x, const int domainIdx)
Definition: BlackoilWellModel_impl.hpp:1618
void updateAverageFormationFactor()
Definition: BlackoilWellModel_impl.hpp:2072
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilWellModel.hpp:106
void initializeWellState(const int timeStepIdx)
Definition: BlackoilWellModel_impl.hpp:878
const Grid & grid() const
Definition: BlackoilWellModel.hpp:367
void updatePrimaryVariables()
Definition: BlackoilWellModel_impl.hpp:2120
void computeWellTemperature()
Definition: BlackoilWellModel_impl.hpp:2223
void addWellPressureEquations(PressureMatrix &jacobian, const BVector &weights, const bool use_well_weights) const
Definition: BlackoilWellModel_impl.hpp:1519
const SimulatorReportSingle & lastReport() const
Definition: BlackoilWellModel_impl.hpp:754
bool updateWellControlsAndNetwork(const bool mandatory_network_balance, const double dt, DeferredLogger &local_deferredLogger)
Definition: BlackoilWellModel_impl.hpp:1275
void addWellContributions(SparseMatrixAdapter &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1509
void assembleWellEq(const double dt)
Definition: BlackoilWellModel_impl.hpp:1398
WellInterfacePtr createWellForWellTest(const std::string &well_name, const int report_step, DeferredLogger &deferred_logger) const
Definition: BlackoilWellModel_impl.hpp:1181
void calculateExplicitQuantities() const
Definition: BlackoilWellModel_impl.hpp:1682
void updateAndCommunicate(const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:1777
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: BlackoilWellModel.hpp:289
void computeTotalRatesForDof(RateVector &rate, unsigned globalIdx) const
Definition: BlackoilWellModel_impl.hpp:841
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:1811
void calcInjResvCoeff(const int fipnum, const int pvtreg, std::vector< Scalar > &resv_coeff) const override
Definition: BlackoilWellModel_impl.hpp:2212
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:677
ConvergenceReport getWellConvergence(const std::vector< Scalar > &B_avg, const bool checkWellGroupControlsAndNetwork=false) const
Definition: BlackoilWellModel_impl.hpp:1631
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:1456
void assembleWellEqWithoutIteration(const double dt)
Definition: BlackoilWellModel_impl.hpp:1424
void updateCellRates()
Definition: BlackoilWellModel_impl.hpp:1444
void assemble(const double dt)
Definition: BlackoilWellModel_impl.hpp:1206
std::size_t local_num_cells_
Definition: BlackoilWellModel.hpp:515
bool alternative_well_rate_init_
Definition: BlackoilWellModel.hpp:518
void timeStepSucceeded(const double simulationTime, const double dt)
Definition: BlackoilWellModel_impl.hpp:764
std::unique_ptr< WellType > createTypedWellPointer(const int wellID, const int time_step) const
Definition: BlackoilWellModel_impl.hpp:1150
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:1927
Simulator & simulator_
Definition: BlackoilWellModel.hpp:487
void createWellContainer(const int report_step) override
Definition: BlackoilWellModel_impl.hpp:921
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:1856
void addReservoirSourceTerms(GlobalEqVector &residual, const std::vector< typename SparseMatrixAdapter::MatrixBlock * > &diagMatAddress) const
Definition: BlackoilWellModel_impl.hpp:1541
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:1594
void addWellPressureEquationsStruct(PressureMatrix &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1571
void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1976
void endReportStep()
Definition: BlackoilWellModel_impl.hpp:737
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