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