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