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