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