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