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