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