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