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