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