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
54#endif
55
58
59#if COMPILE_GPU_BRIDGE
61#endif
62
63#include <algorithm>
64#include <cassert>
65#include <iomanip>
66#include <utility>
67#include <optional>
68
69#include <fmt/format.h>
70
71namespace Opm {
72 template<typename TypeTag>
75 : WellConnectionModule(*this, simulator.gridView().comm())
76 , BlackoilWellModelGeneric<Scalar, IndexTraits>(simulator.vanguard().schedule(),
77 gaslift_,
78 network_,
79 simulator.vanguard().summaryState(),
80 simulator.vanguard().eclState(),
81 FluidSystem::phaseUsage(),
82 simulator.gridView().comm())
83 , simulator_(simulator)
84 , guide_rate_handler_{
85 *this,
86 simulator.vanguard().schedule(),
87 simulator.vanguard().summaryState(),
88 simulator.vanguard().grid().comm()
89 }
90 , gaslift_(this->terminal_output_)
91 , network_(*this)
92 {
93 local_num_cells_ = simulator_.gridView().size(0);
94
95 // Number of cells the global grid view
96 global_num_cells_ = simulator_.vanguard().globalNumCells();
97
98 {
99 auto& parallel_wells = simulator.vanguard().parallelWells();
100
101 this->parallel_well_info_.reserve(parallel_wells.size());
102 for( const auto& name_bool : parallel_wells) {
103 this->parallel_well_info_.emplace_back(name_bool, grid().comm());
104 }
105 }
106
108 Parameters::Get<Parameters::AlternativeWellRateInit>();
109
110 using SourceDataSpan =
111 typename PAvgDynamicSourceData<Scalar>::template SourceDataSpan<Scalar>;
112
114 [this](const std::size_t globalIndex)
115 { return this->compressedIndexForInterior(globalIndex); },
116 [this](const int localCell, SourceDataSpan sourceTerms)
117 {
118 using Item = typename SourceDataSpan::Item;
119
120 const auto* intQuants = this->simulator_.model()
121 .cachedIntensiveQuantities(localCell, /*timeIndex = */0);
122 const auto& fs = intQuants->fluidState();
123
124 sourceTerms
125 .set(Item::PoreVol, intQuants->porosity().value() *
126 this->simulator_.model().dofTotalVolume(localCell))
127 .set(Item::Depth, this->depth_[localCell]);
128
129 constexpr auto io = FluidSystem::oilPhaseIdx;
130 constexpr auto ig = FluidSystem::gasPhaseIdx;
131 constexpr auto iw = FluidSystem::waterPhaseIdx;
132
133 // Ideally, these would be 'constexpr'.
134 const auto haveOil = FluidSystem::phaseIsActive(io);
135 const auto haveGas = FluidSystem::phaseIsActive(ig);
136 const auto haveWat = FluidSystem::phaseIsActive(iw);
137
138 auto weightedPhaseDensity = [&fs](const auto ip)
139 {
140 return fs.saturation(ip).value() * fs.density(ip).value();
141 };
142
143 if (haveOil) { sourceTerms.set(Item::Pressure, fs.pressure(io).value()); }
144 else if (haveGas) { sourceTerms.set(Item::Pressure, fs.pressure(ig).value()); }
145 else { sourceTerms.set(Item::Pressure, fs.pressure(iw).value()); }
146
147 // Strictly speaking, assumes SUM(s[p]) == 1.
148 auto rho = 0.0;
149 if (haveOil) { rho += weightedPhaseDensity(io); }
150 if (haveGas) { rho += weightedPhaseDensity(ig); }
151 if (haveWat) { rho += weightedPhaseDensity(iw); }
152
153 sourceTerms.set(Item::MixtureDensity, rho);
154 }
155 );
156 }
157
158 template<typename TypeTag>
159 void
161 init()
162 {
163 extractLegacyCellPvtRegionIndex_();
164 extractLegacyDepth_();
165
166 gravity_ = simulator_.problem().gravity()[2];
167
168 this->initial_step_ = true;
169
170 // add the eWoms auxiliary module for the wells to the list
171 simulator_.model().addAuxiliaryModule(this);
172
173 is_cell_perforated_.resize(local_num_cells_, false);
174 }
175
176
177 template<typename TypeTag>
178 void
180 initWellContainer(const int reportStepIdx)
181 {
182 const uint64_t effective_events_mask = ScheduleEvents::WELL_STATUS_CHANGE
183 + ScheduleEvents::NEW_WELL;
184 const auto& events = this->schedule()[reportStepIdx].wellgroup_events();
185 for (auto& wellPtr : this->well_container_) {
186 const bool well_opened_this_step = this->report_step_starts_ &&
187 events.hasEvent(wellPtr->name(),
188 effective_events_mask);
189 wellPtr->init(this->depth_, this->gravity_,
190 this->B_avg_, well_opened_this_step);
191 }
192 }
193
194 template<typename TypeTag>
195 void
197 beginReportStep(const int timeStepIdx)
198 {
199 DeferredLogger local_deferredLogger{};
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_LOG(local_deferredLogger,
233 "beginReportStep() failed: ",
234 this->terminal_output_, comm)
235
236 // Store the current well and group states in order to recover in
237 // the case of failed iterations
238 this->commitWGState();
239
240 this->wellStructureChangedDynamically_ = false;
241 }
242
243
244
245
246
247 template <typename TypeTag>
248 void
250 initializeLocalWellStructure(const int reportStepIdx,
251 const bool enableWellPIScaling)
252 {
253 DeferredLogger local_deferredLogger{};
254
255 const auto& comm = this->simulator_.vanguard().grid().comm();
256
257 // Wells_ecl_ holds this rank's wells, both open and stopped/shut.
258 this->wells_ecl_ = this->getLocalWells(reportStepIdx);
259 this->local_parallel_well_info_ =
260 this->createLocalParallelWellInfo(this->wells_ecl_);
261
262 // At least initializeWellState() might be throw an exception in
263 // UniformTabulated2DFunction. Playing it safe by extending the
264 // scope a bit.
266 {
267 this->initializeWellPerfData();
268 this->initializeWellState(reportStepIdx);
269 this->wbp_.initializeWBPCalculationService();
270
271 if (this->param_.use_multisegment_well_ && this->anyMSWellOpenLocal()) {
272 this->wellState().initWellStateMSWell(this->wells_ecl_, &this->prevWellState());
273 }
274
275 this->initializeWellProdIndCalculators();
276
277 if (enableWellPIScaling && this->schedule()[reportStepIdx].events()
278 .hasEvent(ScheduleEvents::Events::WELL_PRODUCTIVITY_INDEX))
279 {
280 this->runWellPIScaling(reportStepIdx, local_deferredLogger);
281 }
282 }
283 OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger,
284 "Failed to initialize local well structure: ",
285 this->terminal_output_, comm)
286 }
287
288
289
290
291
292 template <typename TypeTag>
293 void
295 initializeGroupStructure(const int reportStepIdx)
296 {
297 DeferredLogger local_deferredLogger{};
298
299 const auto& comm = this->simulator_.vanguard().grid().comm();
300
302 {
303 const auto& fieldGroup =
304 this->schedule().getGroup("FIELD", reportStepIdx);
305
306 this->groupStateHelper().setCmodeGroup(fieldGroup);
307
308 // Define per region average pressure calculators for use by
309 // pressure maintenance groups (GPMAINT keyword).
310 if (this->schedule()[reportStepIdx].has_gpmaint()) {
311 this->groupStateHelper().setRegionAveragePressureCalculator(
312 fieldGroup,
313 this->eclState_.fieldProps(),
314 this->regionalAveragePressureCalculator_
315 );
316 }
317 }
318 OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger,
319 "Failed to initialize group structure: ",
320 this->terminal_output_, comm)
321 }
322
323
324
325
326
327 // called at the beginning of a time step
328 template<typename TypeTag>
329 void
332 {
333 OPM_TIMEBLOCK(beginTimeStep);
334
335 this->updateAverageFormationFactor();
336
337 DeferredLogger local_deferredLogger;
338
339#ifdef RESERVOIR_COUPLING_ENABLED
340 auto rescoup_logger_guard = this->setupRescoupScopedLogger(local_deferredLogger);
341#endif
342 this->switched_prod_groups_.clear();
343 this->switched_inj_groups_.clear();
344
345 if (this->wellStructureChangedDynamically_) {
346 // Something altered the well structure/topology. Possibly
347 // WELSPECS/COMPDAT and/or WELOPEN run from an ACTIONX block.
348 // Reconstruct the local wells to account for the new well
349 // structure.
350 const auto reportStepIdx =
351 this->simulator_.episodeIndex();
352
353 // Disable WELPI scaling when well structure is updated in the
354 // middle of a report step.
355 const auto enableWellPIScaling = false;
356
357 this->initializeLocalWellStructure(reportStepIdx, enableWellPIScaling);
358 this->initializeGroupStructure(reportStepIdx);
359
360 this->commitWGState();
361
362 // Reset topology flag to signal that we've handled this
363 // structure change. That way we don't end up here in
364 // subsequent calls to beginTimeStep() unless there's a new
365 // dynamic change to the well structure during a report step.
366 this->wellStructureChangedDynamically_ = false;
367 }
368
369 this->resetWGState();
370 const int reportStepIdx = simulator_.episodeIndex();
371
372 this->wellState().updateWellsDefaultALQ(this->schedule(), reportStepIdx, this->summaryState());
373 this->wellState().gliftTimeStepInit();
374
375 const double simulationTime = simulator_.time();
377 {
378 // test wells
379 wellTesting(reportStepIdx, simulationTime, local_deferredLogger);
380
381 // create the well container
382 createWellContainer(reportStepIdx);
383
384#ifdef RESERVOIR_COUPLING_ENABLED
385 // Receive all slave group data early to ensure it's available for any calculations
386 if (this->isReservoirCouplingMaster()) {
387 this->receiveSlaveGroupData();
388 }
389#endif
390
391 // we need to update the group data after the well is created
392 // to make sure we get the correct mapping.
393 this->updateAndCommunicateGroupData(reportStepIdx,
394 simulator_.model().newtonMethod().numIterations(),
395 param_.nupcol_group_rate_tolerance_, /*update_wellgrouptarget*/ false,
396 local_deferredLogger);
397
398 // Wells are active if they are active wells on at least one process.
399 const Grid& grid = simulator_.vanguard().grid();
400 this->wells_active_ = grid.comm().max(!this->well_container_.empty());
401
402 // do the initialization for all the wells
403 // TODO: to see whether we can postpone of the intialization of the well containers to
404 // optimize the usage of the following several member variables
405 this->initWellContainer(reportStepIdx);
406
407 // update the updated cell flag
408 std::fill(is_cell_perforated_.begin(), is_cell_perforated_.end(), false);
409 for (auto& well : well_container_) {
410 well->updatePerforatedCell(is_cell_perforated_);
411 }
412
413 // calculate the efficiency factors for each well
414 this->calculateEfficiencyFactors(reportStepIdx);
415
416 if constexpr (has_polymer_)
417 {
418 if (PolymerModule::hasPlyshlog() || getPropValue<TypeTag, Properties::EnablePolymerMW>() ) {
419 this->setRepRadiusPerfLength();
420 }
421 }
422
423 }
424
425 OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger, "beginTimeStep() failed: ",
426 this->terminal_output_, simulator_.vanguard().grid().comm());
427
428 for (auto& well : well_container_) {
429 well->setVFPProperties(this->vfp_properties_.get());
430 well->setGuideRate(&this->guideRate_);
431 }
432
433 this->updateFiltrationModelsPreStep(local_deferredLogger);
434
435 // Close completions due to economic reasons
436 for (auto& well : well_container_) {
437 well->closeCompletions(this->wellTestState());
438 }
439
440 // we need the inj_multiplier from the previous time step
441 this->initInjMult();
442
443 if (alternative_well_rate_init_) {
444 // Update the well rates of well_state_, if only single-phase rates, to
445 // have proper multi-phase rates proportional to rates at bhp zero.
446 // This is done only for producers, as injectors will only have a single
447 // nonzero phase anyway.
448 for (const auto& well : well_container_) {
449 if (well->isProducer() && !well->wellIsStopped()) {
450 well->initializeProducerWellState(simulator_, this->wellState(), local_deferredLogger);
451 }
452 }
453 }
454
455 for (const auto& well : well_container_) {
456 if (well->isVFPActive(local_deferredLogger)){
457 well->setPrevSurfaceRates(this->wellState(), this->prevWellState());
458 }
459 }
460 try {
461 this->updateWellPotentials(reportStepIdx,
462 /*onlyAfterEvent*/true,
463 simulator_.vanguard().summaryConfig(),
464 local_deferredLogger);
465 } catch ( std::runtime_error& e ) {
466 const std::string msg = "A zero well potential is returned for output purposes. ";
467 local_deferredLogger.warning("WELL_POTENTIAL_CALCULATION_FAILED", msg);
468 }
469 this->guide_rate_handler_.setLogger(&local_deferredLogger);
470 //update guide rates
471 this->guide_rate_handler_.updateGuideRates(
472 reportStepIdx, simulationTime, this->wellState(), this->groupState()
473 );
474#ifdef RESERVOIR_COUPLING_ENABLED
475 if (this->isReservoirCouplingSlave()) {
476 this->sendSlaveGroupDataToMaster();
477 }
478#endif
479 std::string exc_msg;
480 auto exc_type = ExceptionType::NONE;
481 // update gpmaint targets
482 if (this->schedule_[reportStepIdx].has_gpmaint()) {
483 for (const auto& calculator : regionalAveragePressureCalculator_) {
484 calculator.second->template defineState<ElementContext>(simulator_);
485 }
486 const double dt = simulator_.timeStepSize();
487 const Group& fieldGroup = this->schedule().getGroup("FIELD", reportStepIdx);
488 this->groupStateHelper().updateGpMaintTargetForGroups(fieldGroup,
489 regionalAveragePressureCalculator_,
490 dt);
491 }
492
493 this->updateAndCommunicateGroupData(reportStepIdx,
494 simulator_.model().newtonMethod().numIterations(),
495 param_.nupcol_group_rate_tolerance_,
496 /*update_wellgrouptarget*/ true,
497 local_deferredLogger);
498 try {
499 // Compute initial well solution for new wells and injectors that change injection type i.e. WAG.
500 for (auto& well : well_container_) {
501 const uint64_t effective_events_mask = ScheduleEvents::WELL_STATUS_CHANGE
502 + ScheduleEvents::INJECTION_TYPE_CHANGED
503 + ScheduleEvents::WELL_SWITCHED_INJECTOR_PRODUCER
504 + ScheduleEvents::NEW_WELL;
505
506 const auto& events = this->schedule()[reportStepIdx].wellgroup_events();
507 const bool event = this->report_step_starts_ && events.hasEvent(well->name(), effective_events_mask);
508 const bool dyn_status_change = this->wellState().well(well->name()).status
509 != this->prevWellState().well(well->name()).status;
510
511 if (event || dyn_status_change) {
512 try {
513 well->scaleSegmentRatesAndPressure(this->wellState());
514 well->calculateExplicitQuantities(simulator_, this->groupStateHelper(), local_deferredLogger);
515 well->updateWellStateWithTarget(simulator_, this->groupStateHelper(), this->wellState(), local_deferredLogger);
516 well->updatePrimaryVariables(this->groupStateHelper(), local_deferredLogger);
517 well->solveWellEquation(
518 simulator_, this->groupStateHelper(), this->wellState(), local_deferredLogger
519 );
520 } catch (const std::exception& e) {
521 const std::string msg = "Compute initial well solution for new well " + well->name() + " failed. Continue with zero initial rates";
522 local_deferredLogger.warning("WELL_INITIAL_SOLVE_FAILED", msg);
523 }
524 }
525 }
526 }
527 // Catch clauses for all errors setting exc_type and exc_msg
528 OPM_PARALLEL_CATCH_CLAUSE(exc_type, exc_msg);
529
530 if (exc_type != ExceptionType::NONE) {
531 const std::string msg = "Compute initial well solution for new wells failed. Continue with zero initial rates";
532 local_deferredLogger.warning("WELL_INITIAL_SOLVE_FAILED", msg);
533 }
534
535 const auto& comm = simulator_.vanguard().grid().comm();
536 logAndCheckForExceptionsAndThrow(local_deferredLogger,
537 exc_type, "beginTimeStep() failed: " + exc_msg, this->terminal_output_, comm);
538
539 }
540
541#ifdef RESERVOIR_COUPLING_ENABLED
542 // Automatically manages the lifecycle of the DeferredLogger pointer
543 // in the reservoir coupling logger. Ensures the logger is properly
544 // cleared when it goes out of scope, preventing dangling pointer issues:
545 //
546 // - The ScopedLoggerGuard constructor sets the logger pointer
547 // - When the guard goes out of scope, the destructor clears the pointer
548 // - Move semantics transfer ownership safely when returning from this function
549 // - The moved-from guard is "nullified" and its destructor does nothing
550 // - Only the final guard in the caller will clear the logger
551 template<typename TypeTag>
552 std::optional<ReservoirCoupling::ScopedLoggerGuard>
555 if (this->isReservoirCouplingMaster()) {
557 this->reservoirCouplingMaster().getLogger(),
558 &local_logger
559 };
560 } else if (this->isReservoirCouplingSlave()) {
561 return ReservoirCoupling::ScopedLoggerGuard{
562 this->reservoirCouplingSlave().getLogger(),
563 &local_logger
564 };
565 }
566 return std::nullopt;
567 }
568
569 template<typename TypeTag>
570 void
571 BlackoilWellModel<TypeTag>::
572 receiveSlaveGroupData()
573 {
574 assert(this->isReservoirCouplingMaster());
575 RescoupReceiveSlaveGroupData<Scalar, IndexTraits> slave_group_data_receiver{
576 this->groupStateHelper(),
577 };
578 slave_group_data_receiver.receiveSlaveGroupData();
579 }
580
581 template<typename TypeTag>
582 void
583 BlackoilWellModel<TypeTag>::
584 sendSlaveGroupDataToMaster()
585 {
586 assert(this->isReservoirCouplingSlave());
587 RescoupSendSlaveGroupData<Scalar, IndexTraits> slave_group_data_sender{this->groupStateHelper()};
588 slave_group_data_sender.sendSlaveGroupDataToMaster();
589 }
590#endif // RESERVOIR_COUPLING_ENABLED
591
592 template<typename TypeTag>
593 void
595 const double simulationTime,
596 DeferredLogger& deferred_logger)
597 {
598 for (const std::string& well_name : this->getWellsForTesting(timeStepIdx, simulationTime)) {
599 const Well& wellEcl = this->schedule().getWell(well_name, timeStepIdx);
600 if (wellEcl.getStatus() == Well::Status::SHUT)
601 continue;
602
603 WellInterfacePtr well = createWellForWellTest(well_name, timeStepIdx, deferred_logger);
604 // some preparation before the well can be used
605 well->init(depth_, gravity_, B_avg_, true);
606
607 Scalar well_efficiency_factor = wellEcl.getEfficiencyFactor() *
608 this->wellState().getGlobalEfficiencyScalingFactor(well_name);
609 this->groupStateHelper().accumulateGroupEfficiencyFactor(
610 this->schedule().getGroup(wellEcl.groupName(), timeStepIdx),
611 well_efficiency_factor
612 );
613
614 well->setWellEfficiencyFactor(well_efficiency_factor);
615 well->setVFPProperties(this->vfp_properties_.get());
616 well->setGuideRate(&this->guideRate_);
617
618 // initialize rates/previous rates to prevent zero fractions in vfp-interpolation
619 if (well->isProducer() && alternative_well_rate_init_) {
620 well->initializeProducerWellState(simulator_, this->wellState(), deferred_logger);
621 }
622 if (well->isVFPActive(deferred_logger)) {
623 well->setPrevSurfaceRates(this->wellState(), this->prevWellState());
624 }
625
626 const auto& network = this->schedule()[timeStepIdx].network();
627 if (network.active()) {
628 this->network_.initializeWell(*well);
629 }
630 try {
631 using GLiftEclWells = typename GasLiftGroupInfo<Scalar, IndexTraits>::GLiftEclWells;
632 GLiftEclWells ecl_well_map;
633 gaslift_.initGliftEclWellMap(well_container_, ecl_well_map);
634 well->wellTesting(simulator_,
635 simulationTime,
636 this->groupStateHelper(),
637 this->wellState(),
638 this->wellTestState(),
639 ecl_well_map,
640 this->well_open_times_,
641 deferred_logger);
642 } catch (const std::exception& e) {
643 const std::string msg = fmt::format("Exception during testing of well: {}. The well will not open.\n Exception message: {}", wellEcl.name(), e.what());
644 deferred_logger.warning("WELL_TESTING_FAILED", msg);
645 }
646 }
647 }
648
649
650
651
652
653 // called at the end of a report step
654 template<typename TypeTag>
655 void
658 {
659 // Clear the communication data structures for above values.
660 for (auto&& pinfo : this->local_parallel_well_info_)
661 {
662 pinfo.get().clear();
663 }
664 }
665
666
667
668
669
670 // called at the end of a report step
671 template<typename TypeTag>
674 lastReport() const {return last_report_; }
675
676
677
678
679
680 // called at the end of a time step
681 template<typename TypeTag>
682 void
684 timeStepSucceeded(const double simulationTime, const double dt)
685 {
686 this->closed_this_step_.clear();
687
688 // time step is finished and we are not any more at the beginning of an report step
689 this->report_step_starts_ = false;
690 const int reportStepIdx = simulator_.episodeIndex();
691
692 DeferredLogger local_deferredLogger;
693 for (const auto& well : well_container_) {
694 if (getPropValue<TypeTag, Properties::EnablePolymerMW>() && well->isInjector()) {
695 well->updateWaterThroughput(dt, this->wellState());
696 }
697 }
698 // update connection transmissibility factor and d factor (if applicable) in the wellstate
699 for (const auto& well : well_container_) {
700 well->updateConnectionTransmissibilityFactor(simulator_, this->wellState().well(well->indexOfWell()));
701 well->updateConnectionDFactor(simulator_, this->wellState().well(well->indexOfWell()));
702 }
703
704 if (Indices::waterEnabled) {
705 this->updateFiltrationModelsPostStep(dt, FluidSystem::waterPhaseIdx, local_deferredLogger);
706 }
707
708 // WINJMULT: At the end of the time step, update the inj_multiplier saved in WellState for later use
709 this->updateInjMult(local_deferredLogger);
710
711 // report well switching
712 for (const auto& well : well_container_) {
713 well->reportWellSwitching(this->wellState().well(well->indexOfWell()), local_deferredLogger);
714 }
715 // report group switching
716 if (this->terminal_output_) {
717 this->reportGroupSwitching(local_deferredLogger);
718 }
719
720 // update the rate converter with current averages pressures etc in
721 rateConverter_->template defineState<ElementContext>(simulator_);
722
723 // calculate the well potentials
724 try {
725 this->updateWellPotentials(reportStepIdx,
726 /*onlyAfterEvent*/false,
727 simulator_.vanguard().summaryConfig(),
728 local_deferredLogger);
729 } catch ( std::runtime_error& e ) {
730 const std::string msg = "A zero well potential is returned for output purposes. ";
731 local_deferredLogger.warning("WELL_POTENTIAL_CALCULATION_FAILED", msg);
732 }
733
734 updateWellTestState(simulationTime, this->wellTestState());
735
736 // check group sales limits at the end of the timestep
737 const Group& fieldGroup = this->schedule_.getGroup("FIELD", reportStepIdx);
738 this->checkGEconLimits(fieldGroup, simulationTime,
739 simulator_.episodeIndex(), local_deferredLogger);
740 this->checkGconsaleLimits(fieldGroup, this->wellState(),
741 simulator_.episodeIndex(), local_deferredLogger);
742
743 this->calculateProductivityIndexValues(local_deferredLogger);
744
745 this->commitWGState();
746
747 const Opm::Parallel::Communication& comm = grid().comm();
748 DeferredLogger global_deferredLogger = gatherDeferredLogger(local_deferredLogger, comm);
749 if (this->terminal_output_) {
750 global_deferredLogger.logMessages();
751 }
752
753 //reporting output temperatures
754 this->computeWellTemperature();
755 }
756
757
758 template<typename TypeTag>
759 void
762 unsigned elemIdx) const
763 {
764 rate = 0;
765
766 if (!is_cell_perforated_[elemIdx] || cellRates_.count(elemIdx) == 0) {
767 return;
768 }
769
770 rate = cellRates_.at(elemIdx);
771 }
772
773
774 template<typename TypeTag>
775 template <class Context>
776 void
779 const Context& context,
780 unsigned spaceIdx,
781 unsigned timeIdx) const
782 {
783 rate = 0;
784 int elemIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
785
786 if (!is_cell_perforated_[elemIdx] || cellRates_.count(elemIdx) == 0) {
787 return;
788 }
789
790 rate = cellRates_.at(elemIdx);
791 }
792
793
794
795 template<typename TypeTag>
796 void
798 initializeWellState(const int timeStepIdx)
799 {
800 const auto pressIx = []()
801 {
802 if (Indices::oilEnabled) { return FluidSystem::oilPhaseIdx; }
803 if (Indices::waterEnabled) { return FluidSystem::waterPhaseIdx; }
804
805 return FluidSystem::gasPhaseIdx;
806 }();
807
808 auto cellPressures = std::vector<Scalar>(this->local_num_cells_, Scalar{0});
809 auto cellTemperatures = std::vector<Scalar>(this->local_num_cells_, Scalar{0});
810
811 auto elemCtx = ElementContext { this->simulator_ };
812 const auto& gridView = this->simulator_.vanguard().gridView();
813
815 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
816 elemCtx.updatePrimaryStencil(elem);
817 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
818
819 const auto ix = elemCtx.globalSpaceIndex(/*spaceIdx=*/0, /*timeIdx=*/0);
820 const auto& fs = elemCtx.intensiveQuantities(/*spaceIdx=*/0, /*timeIdx=*/0).fluidState();
821
822 cellPressures[ix] = fs.pressure(pressIx).value();
823 cellTemperatures[ix] = fs.temperature(0).value();
824 }
825 OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel::initializeWellState() failed: ",
826 this->simulator_.vanguard().grid().comm());
827
828 this->wellState().init(cellPressures, cellTemperatures, this->schedule(), this->wells_ecl_,
829 this->local_parallel_well_info_, timeStepIdx,
830 &this->prevWellState(), this->well_perf_data_,
831 this->summaryState(), simulator_.vanguard().enableDistributedWells());
832 }
833
834
835
836
837
838 template<typename TypeTag>
839 void
841 createWellContainer(const int report_step)
842 {
843 DeferredLogger local_deferredLogger;
844
845 const int nw = this->numLocalWells();
846
847 well_container_.clear();
848
849 if (nw > 0) {
850 well_container_.reserve(nw);
851
852 const auto& wmatcher = this->schedule().wellMatcher(report_step);
853 const auto& wcycle = this->schedule()[report_step].wcycle.get();
854
855 // First loop and check for status changes. This is necessary
856 // as wcycle needs the updated open/close times.
857 std::for_each(this->wells_ecl_.begin(), this->wells_ecl_.end(),
858 [this, &wg_events = this->report_step_start_events_](const auto& well_ecl)
859 {
860 if (!well_ecl.hasConnections()) {
861 // No connections in this well. Nothing to do.
862 return;
863 }
864
865 constexpr auto events_mask = ScheduleEvents::WELL_STATUS_CHANGE |
866 ScheduleEvents::REQUEST_OPEN_WELL |
867 ScheduleEvents::REQUEST_SHUT_WELL;
868 const bool well_event =
869 this->report_step_starts_ &&
870 wg_events.hasEvent(well_ecl.name(), events_mask);
871 // WCYCLE is suspendended by explicit SHUT events by the user.
872 // and restarted after explicit OPEN events.
873 // Note: OPEN or SHUT event does not necessary mean the well
874 // actually opened or shut at this point as the simulator could
875 // have done this by operabilty checks and well testing. This
876 // may need further testing and imply code changes to cope with
877 // these corner cases.
878 if (well_event) {
879 if (well_ecl.getStatus() == WellStatus::OPEN) {
880 this->well_open_times_.insert_or_assign(well_ecl.name(),
881 this->simulator_.time());
882 this->well_close_times_.erase(well_ecl.name());
883 } else if (well_ecl.getStatus() == WellStatus::SHUT) {
884 this->well_close_times_.insert_or_assign(well_ecl.name(),
885 this->simulator_.time());
886 this->well_open_times_.erase(well_ecl.name());
887 }
888 }
889 });
890
891 // Grab wcycle states. This needs to run before the schedule gets processed
892 const auto cycle_states = wcycle.wellStatus(this->simulator_.time(),
893 wmatcher,
894 this->well_open_times_,
895 this->well_close_times_);
896
897 for (int w = 0; w < nw; ++w) {
898 const Well& well_ecl = this->wells_ecl_[w];
899
900 if (!well_ecl.hasConnections()) {
901 // No connections in this well. Nothing to do.
902 continue;
903 }
904
905 const std::string& well_name = well_ecl.name();
906 const auto well_status = this->schedule()
907 .getWell(well_name, report_step).getStatus();
908
909 const bool shut_event = this->wellState().well(w).events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE)
910 && well_status == Well::Status::SHUT;
911 const bool open_event = this->wellState().well(w).events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE)
912 && well_status == Well::Status::OPEN;
913 const auto& ws = this->wellState().well(well_name);
914
915 if (shut_event && ws.status != Well::Status::SHUT) {
916 this->closed_this_step_.insert(well_name);
917 this->wellState().shutWell(w);
918 } else if (open_event && ws.status != Well::Status::OPEN) {
919 this->wellState().openWell(w);
920 }
921
922 // A new WCON keywords can re-open a well that was closed/shut due to Physical limit
923 if (this->wellTestState().well_is_closed(well_name)) {
924 // The well was shut this timestep, we are most likely retrying
925 // a timestep without the well in question, after it caused
926 // repeated timestep cuts. It should therefore not be opened,
927 // even if it was new or received new targets this report step.
928 const bool closed_this_step = (this->wellTestState().lastTestTime(well_name) == simulator_.time());
929 // TODO: more checking here, to make sure this standard more specific and complete
930 // maybe there is some WCON keywords will not open the well
931 auto& events = this->wellState().well(w).events;
932 if (events.hasEvent(ScheduleEvents::REQUEST_OPEN_WELL)) {
933 if (!closed_this_step) {
934 this->wellTestState().open_well(well_name);
935 this->wellTestState().open_completions(well_name);
936 this->well_open_times_.insert_or_assign(well_name,
937 this->simulator_.time());
938 this->well_close_times_.erase(well_name);
939 }
940 events.clearEvent(ScheduleEvents::REQUEST_OPEN_WELL);
941 }
942 }
943
944 // TODO: should we do this for all kinds of closing reasons?
945 // something like wellTestState().hasWell(well_name)?
946 if (this->wellTestState().well_is_closed(well_name))
947 {
948 if (well_ecl.getAutomaticShutIn()) {
949 // shut wells are not added to the well container
950 this->wellState().shutWell(w);
951 this->well_close_times_.erase(well_name);
952 this->well_open_times_.erase(well_name);
953 continue;
954 } else {
955 if (!well_ecl.getAllowCrossFlow()) {
956 // stopped wells where cross flow is not allowed
957 // are not added to the well container
958 this->wellState().shutWell(w);
959 this->well_close_times_.erase(well_name);
960 this->well_open_times_.erase(well_name);
961 continue;
962 }
963 // stopped wells are added to the container but marked as stopped
964 this->wellState().stopWell(w);
965 }
966 }
967
968 // shut wells with zero rante constraints and disallowing
969 if (!well_ecl.getAllowCrossFlow()) {
970 const bool any_zero_rate_constraint = well_ecl.isProducer()
971 ? well_ecl.productionControls(this->summaryState_).anyZeroRateConstraint()
972 : well_ecl.injectionControls(this->summaryState_).anyZeroRateConstraint();
973 if (any_zero_rate_constraint) {
974 // Treat as shut, do not add to container.
975 local_deferredLogger.debug(fmt::format(" Well {} gets shut due to having zero rate constraint and disallowing crossflow ", well_ecl.name()) );
976 this->wellState().shutWell(w);
977 this->well_close_times_.erase(well_name);
978 this->well_open_times_.erase(well_name);
979 continue;
980 }
981 }
982
983 if (!wcycle.empty()) {
984 const auto it = cycle_states.find(well_name);
985 if (it != cycle_states.end()) {
986 if (!it->second || well_status == Well::Status::SHUT) {
987 // If well is shut in schedule we keep it shut
988 if (well_status == Well::Status::SHUT) {
989 this->well_open_times_.erase(well_name);
990 this->well_close_times_.erase(well_name);
991 }
992 this->wellState().shutWell(w);
993 continue;
994 } else {
995 this->wellState().openWell(w);
996 }
997 }
998 }
999
1000 // We dont add SHUT wells to the container
1001 if (ws.status == Well::Status::SHUT) {
1002 continue;
1003 }
1004
1005 well_container_.emplace_back(this->createWellPointer(w, report_step));
1006
1007 if (ws.status == Well::Status::STOP) {
1008 well_container_.back()->stopWell();
1009 this->well_close_times_.erase(well_name);
1010 this->well_open_times_.erase(well_name);
1011 }
1012 }
1013
1014 if (!wcycle.empty()) {
1015 const auto schedule_open =
1016 [&wg_events = this->report_step_start_events_](const std::string& name)
1017 {
1018 return wg_events.hasEvent(name, ScheduleEvents::REQUEST_OPEN_WELL);
1019 };
1020 for (const auto& [wname, wscale] : wcycle.efficiencyScale(this->simulator_.time(),
1021 this->simulator_.timeStepSize(),
1022 wmatcher,
1023 this->well_open_times_,
1024 schedule_open))
1025 {
1026 this->wellState().updateEfficiencyScalingFactor(wname, wscale);
1027 this->schedule_.add_event(ScheduleEvents::WELLGROUP_EFFICIENCY_UPDATE, report_step);
1028 }
1029 }
1030 }
1031
1032 // Collect log messages and print.
1033
1034 const Opm::Parallel::Communication& comm = grid().comm();
1035 DeferredLogger global_deferredLogger = gatherDeferredLogger(local_deferredLogger, comm);
1036 if (this->terminal_output_) {
1037 global_deferredLogger.logMessages();
1038 }
1039
1040 this->well_container_generic_.clear();
1041 for (auto& w : well_container_) {
1042 this->well_container_generic_.push_back(w.get());
1043 }
1044
1045 this->network_.initialize(report_step);
1046
1047 this->wbp_.registerOpenWellsForWBPCalculation();
1048 }
1049
1050
1051
1052
1053
1054 template <typename TypeTag>
1055 typename BlackoilWellModel<TypeTag>::WellInterfacePtr
1057 createWellPointer(const int wellID, const int report_step) const
1058 {
1059 const auto is_multiseg = this->wells_ecl_[wellID].isMultiSegment();
1060
1061 if (! (this->param_.use_multisegment_well_ && is_multiseg)) {
1062 return this->template createTypedWellPointer<StandardWell<TypeTag>>(wellID, report_step);
1063 }
1064 else {
1065 return this->template createTypedWellPointer<MultisegmentWell<TypeTag>>(wellID, report_step);
1066 }
1067 }
1068
1069
1070
1071
1072
1073 template <typename TypeTag>
1074 template <typename WellType>
1075 std::unique_ptr<WellType>
1077 createTypedWellPointer(const int wellID, const int time_step) const
1078 {
1079 // Use the pvtRegionIdx from the top cell
1080 const auto& perf_data = this->well_perf_data_[wellID];
1081
1082 // Cater for case where local part might have no perforations.
1083 const auto pvtreg = perf_data.empty()
1084 ? 0 : this->pvt_region_idx_[perf_data.front().cell_index];
1085
1086 const auto& parallel_well_info = this->local_parallel_well_info_[wellID].get();
1087 const auto global_pvtreg = parallel_well_info.broadcastFirstPerforationValue(pvtreg);
1088
1089 return std::make_unique<WellType>(this->wells_ecl_[wellID],
1090 parallel_well_info,
1091 time_step,
1092 this->param_,
1093 *this->rateConverter_,
1094 global_pvtreg,
1095 this->numConservationQuantities(),
1096 this->numPhases(),
1097 wellID,
1098 perf_data);
1099 }
1100
1101
1102
1103
1104
1105 template<typename TypeTag>
1108 createWellForWellTest(const std::string& well_name,
1109 const int report_step,
1110 DeferredLogger& deferred_logger) const
1111 {
1112 // Finding the location of the well in wells_ecl
1113 const auto it = std::find_if(this->wells_ecl_.begin(),
1114 this->wells_ecl_.end(),
1115 [&well_name](const auto& w)
1116 { return well_name == w.name(); });
1117 // It should be able to find in wells_ecl.
1118 if (it == this->wells_ecl_.end()) {
1119 OPM_DEFLOG_THROW(std::logic_error,
1120 fmt::format("Could not find well {} in wells_ecl ", well_name),
1121 deferred_logger);
1122 }
1123
1124 const int pos = static_cast<int>(std::distance(this->wells_ecl_.begin(), it));
1125 return this->createWellPointer(pos, report_step);
1126 }
1127
1128
1129
1130 template<typename TypeTag>
1131 void
1133 assemble(const int iterationIdx,
1134 const double dt)
1135 {
1136 OPM_TIMEFUNCTION();
1137 DeferredLogger local_deferredLogger;
1138
1139 this->guide_rate_handler_.setLogger(&local_deferredLogger);
1141 if (gaslift_.terminalOutput()) {
1142 const std::string msg =
1143 fmt::format("assemble() : iteration {}" , iterationIdx);
1144 gaslift_.gliftDebug(msg, local_deferredLogger);
1145 }
1146 }
1147 last_report_ = SimulatorReportSingle();
1148 Dune::Timer perfTimer;
1149 perfTimer.start();
1150 this->closed_offending_wells_.clear();
1151
1152 {
1153 const int episodeIdx = simulator_.episodeIndex();
1154 const auto& network = this->schedule()[episodeIdx].network();
1155 if (!this->wellsActive() && !network.active()) {
1156 return;
1157 }
1158 }
1159
1160 if (iterationIdx == 0 && this->wellsActive()) {
1161 OPM_TIMEBLOCK(firstIterationAssmble);
1162 // try-catch is needed here as updateWellControls
1163 // contains global communication and has either to
1164 // be reached by all processes or all need to abort
1165 // before.
1167 {
1168 calculateExplicitQuantities(local_deferredLogger);
1169 prepareTimeStep(local_deferredLogger);
1170 }
1171 OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger,
1172 "assemble() failed (It=0): ",
1173 this->terminal_output_, grid().comm());
1174 }
1175
1176 const bool well_group_control_changed = updateWellControlsAndNetwork(false, dt, local_deferredLogger);
1177
1178 // even when there is no wells active, the network nodal pressure still need to be updated through updateWellControlsAndNetwork()
1179 // but there is no need to assemble the well equations
1180 if ( ! this->wellsActive() ) {
1181 return;
1182 }
1183
1184 assembleWellEqWithoutIteration(dt, local_deferredLogger);
1185 // Pre-compute cell rates to we don't have to do this for every cell during linearization...
1186 updateCellRates();
1187
1188 // if group or well control changes we don't consider the
1189 // case converged
1190 last_report_.well_group_control_changed = well_group_control_changed;
1191 last_report_.assemble_time_well += perfTimer.stop();
1192 }
1193
1194
1195
1196
1197 template<typename TypeTag>
1198 bool
1200 updateWellControlsAndNetwork(const bool mandatory_network_balance,
1201 const double dt,
1202 DeferredLogger& local_deferredLogger)
1203 {
1204 OPM_TIMEFUNCTION();
1205 // not necessarily that we always need to update once of the network solutions
1206 bool do_network_update = true;
1207 bool well_group_control_changed = false;
1208 Scalar network_imbalance = 0.0;
1209 // after certain number of the iterations, we use relaxed tolerance for the network update
1210 const std::size_t iteration_to_relax = param_.network_max_strict_outer_iterations_;
1211 // after certain number of the iterations, we terminate
1212 const std::size_t max_iteration = param_.network_max_outer_iterations_;
1213 std::size_t network_update_iteration = 0;
1214 network_needs_more_balancing_force_another_newton_iteration_ = false;
1215 while (do_network_update) {
1216 if (network_update_iteration >= max_iteration ) {
1217 // only output to terminal if we at the last newton iterations where we try to balance the network.
1218 const int episodeIdx = simulator_.episodeIndex();
1219 const int iterationIdx = simulator_.model().newtonMethod().numIterations();
1220 if (this->network_.shouldBalance(episodeIdx, iterationIdx + 1)) {
1221 if (this->terminal_output_) {
1222 const std::string msg = fmt::format("Maximum of {:d} network iterations has been used and we stop the update, \n"
1223 "and try again after the next Newton iteration (imbalance = {:.2e} bar, ctrl_change = {})",
1224 max_iteration, network_imbalance*1.0e-5, well_group_control_changed);
1225 local_deferredLogger.debug(msg);
1226 }
1227 // To avoid stopping the newton iterations too early, before the network is converged,
1228 // we need to report it
1229 network_needs_more_balancing_force_another_newton_iteration_ = true;
1230 } else {
1231 if (this->terminal_output_) {
1232 const std::string msg = fmt::format("Maximum of {:d} network iterations has been used and we stop the update. \n"
1233 "The simulator will continue with unconverged network results (imbalance = {:.2e} bar, ctrl_change = {})",
1234 max_iteration, network_imbalance*1.0e-5, well_group_control_changed);
1235 local_deferredLogger.info(msg);
1236 }
1237 }
1238 break;
1239 }
1240 if (this->terminal_output_ && (network_update_iteration == iteration_to_relax) ) {
1241 local_deferredLogger.debug("We begin using relaxed tolerance for network update now after " + std::to_string(iteration_to_relax) + " iterations ");
1242 }
1243 const bool relax_network_balance = network_update_iteration >= iteration_to_relax;
1244 // Never optimize gas lift in last iteration, to allow network convergence (unless max_iter < 2)
1245 const bool optimize_gas_lift = ( (network_update_iteration + 1) < std::max(max_iteration, static_cast<std::size_t>(2)) );
1246 std::tie(well_group_control_changed, do_network_update, network_imbalance) =
1247 updateWellControlsAndNetworkIteration(mandatory_network_balance, relax_network_balance, optimize_gas_lift, dt,local_deferredLogger);
1248 ++network_update_iteration;
1249 }
1250 return well_group_control_changed;
1251 }
1252
1253
1254
1255
1256 template<typename TypeTag>
1257 std::tuple<bool, bool, typename BlackoilWellModel<TypeTag>::Scalar>
1259 updateWellControlsAndNetworkIteration(const bool mandatory_network_balance,
1260 const bool relax_network_tolerance,
1261 const bool optimize_gas_lift,
1262 const double dt,
1263 DeferredLogger& local_deferredLogger)
1264 {
1265 OPM_TIMEFUNCTION();
1266 const int iterationIdx = simulator_.model().newtonMethod().numIterations();
1267 const int reportStepIdx = simulator_.episodeIndex();
1268 this->updateAndCommunicateGroupData(reportStepIdx, iterationIdx,
1269 param_.nupcol_group_rate_tolerance_, /*update_wellgrouptarget*/ true, local_deferredLogger);
1270 // We need to call updateWellControls before we update the network as
1271 // network updates are only done on thp controlled wells.
1272 // Note that well controls are allowed to change during updateNetwork
1273 // and in prepareWellsBeforeAssembling during well solves.
1274 bool well_group_control_changed = updateWellControls(local_deferredLogger);
1275 const auto [more_inner_network_update, network_imbalance] =
1276 this->network_.update(mandatory_network_balance,
1277 local_deferredLogger,
1278 relax_network_tolerance);
1279
1280 bool alq_updated = false;
1282 {
1283 if (optimize_gas_lift) {
1284 // we need to update the potentials if the thp limit as been modified by
1285 // the network balancing
1286 const bool updatePotentials = (this->network_.shouldBalance(reportStepIdx, iterationIdx) ||
1287 mandatory_network_balance);
1288 alq_updated = gaslift_.maybeDoGasLiftOptimize(simulator_,
1289 well_container_,
1290 this->network_.nodePressures(),
1291 updatePotentials,
1292 this->wellState(),
1293 this->groupState(),
1294 local_deferredLogger);
1295 }
1296 prepareWellsBeforeAssembling(dt, local_deferredLogger);
1297 }
1298 OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger,
1299 "updateWellControlsAndNetworkIteration() failed: ",
1300 this->terminal_output_, grid().comm());
1301
1302 // update guide rates
1303 if (alq_updated || BlackoilWellModelGuideRates(*this).
1304 guideRateUpdateIsNeeded(reportStepIdx)) {
1305 const double simulationTime = simulator_.time();
1306 // NOTE: For reservoir coupling: Slave group potentials are only communicated
1307 // at the start of the time step, see beginTimeStep(). Here, we assume those
1308 // potentials remain unchanged during the time step when updating guide rates below.
1309 this->guide_rate_handler_.updateGuideRates(
1310 reportStepIdx, simulationTime, this->wellState(), this->groupState()
1311 );
1312 }
1313 // we need to re-iterate the network when the well group controls changed or gaslift/alq is changed or
1314 // the inner iterations are did not converge
1315 const bool more_network_update = this->network_.shouldBalance(reportStepIdx, iterationIdx) &&
1316 (more_inner_network_update || well_group_control_changed || alq_updated);
1317 return {well_group_control_changed, more_network_update, network_imbalance};
1318 }
1319
1320 template<typename TypeTag>
1321 void
1323 assembleWellEq(const double dt, DeferredLogger& deferred_logger)
1324 {
1325 OPM_TIMEFUNCTION();
1326 for (auto& well : well_container_) {
1327 well->assembleWellEq(simulator_, dt, this->groupStateHelper(), deferred_logger);
1328 }
1329 }
1330
1331
1332 template<typename TypeTag>
1333 void
1335 prepareWellsBeforeAssembling(const double dt, DeferredLogger& deferred_logger)
1336 {
1337 OPM_TIMEFUNCTION();
1338 for (auto& well : well_container_) {
1339 well->prepareWellBeforeAssembling(
1340 simulator_, dt, this->groupStateHelper(), this->wellState(), deferred_logger
1341 );
1342 }
1343 }
1344
1345
1346 template<typename TypeTag>
1347 void
1349 assembleWellEqWithoutIteration(const double dt, DeferredLogger& deferred_logger)
1350 {
1351 OPM_TIMEFUNCTION();
1352 // We make sure that all processes throw in case there is an exception
1353 // on one of them (WetGasPvt::saturationPressure might throw if not converged)
1355
1356 for (auto& well: well_container_) {
1357 well->assembleWellEqWithoutIteration(simulator_, this->groupStateHelper(), dt, this->wellState(),
1358 deferred_logger,
1359 /*solving_with_zero_rate=*/false);
1360 }
1361 OPM_END_PARALLEL_TRY_CATCH_LOG(deferred_logger, "BlackoilWellModel::assembleWellEqWithoutIteration failed: ",
1362 this->terminal_output_, grid().comm());
1363
1364 }
1365
1366 template<typename TypeTag>
1367 void
1370 {
1371 // Pre-compute cell rates for all wells
1372 cellRates_.clear();
1373 for (const auto& well : well_container_) {
1374 well->addCellRates(cellRates_);
1375 }
1376 }
1377
1378 template<typename TypeTag>
1379 void
1381 updateCellRatesForDomain(int domainIndex, const std::map<std::string, int>& well_domain_map)
1382 {
1383 // Pre-compute cell rates only for wells in the specified domain
1384 cellRates_.clear();
1385 for (const auto& well : well_container_) {
1386 const auto it = well_domain_map.find(well->name());
1387 if (it != well_domain_map.end() && it->second == domainIndex) {
1388 well->addCellRates(cellRates_);
1389 }
1390 }
1391 }
1392
1393#if COMPILE_GPU_BRIDGE
1394 template<typename TypeTag>
1395 void
1398 {
1399 // prepare for StandardWells
1401
1402 for(unsigned int i = 0; i < well_container_.size(); i++){
1403 auto& well = well_container_[i];
1404 auto derived = dynamic_cast<StandardWell<TypeTag>*>(well.get());
1405 if (derived) {
1406 wellContribs.addNumBlocks(derived->linSys().getNumBlocks());
1407 }
1408 }
1409
1410 // allocate memory for data from StandardWells
1411 wellContribs.alloc();
1412
1413 for(unsigned int i = 0; i < well_container_.size(); i++){
1414 auto& well = well_container_[i];
1415 // maybe WellInterface could implement addWellContribution()
1416 auto derived_std = dynamic_cast<StandardWell<TypeTag>*>(well.get());
1417 if (derived_std) {
1418 derived_std->linSys().extract(derived_std->numStaticWellEq, wellContribs);
1419 } else {
1420 auto derived_ms = dynamic_cast<MultisegmentWell<TypeTag>*>(well.get());
1421 if (derived_ms) {
1422 derived_ms->linSys().extract(wellContribs);
1423 } else {
1424 OpmLog::warning("Warning unknown type of well");
1425 }
1426 }
1427 }
1428 }
1429#endif
1430
1431 template<typename TypeTag>
1432 void
1434 addWellContributions(SparseMatrixAdapter& jacobian) const
1435 {
1436 for ( const auto& well: well_container_ ) {
1437 well->addWellContributions(jacobian);
1438 }
1439 }
1440
1441 template<typename TypeTag>
1442 void
1445 const BVector& weights,
1446 const bool use_well_weights) const
1447 {
1448 int nw = this->numLocalWellsEnd();
1449 int rdofs = local_num_cells_;
1450 for ( int i = 0; i < nw; i++ ) {
1451 int wdof = rdofs + i;
1452 jacobian[wdof][wdof] = 1.0;// better scaling ?
1453 }
1454
1455 for (const auto& well : well_container_) {
1456 well->addWellPressureEquations(jacobian,
1457 weights,
1458 pressureVarIndex,
1459 use_well_weights,
1460 this->wellState());
1461 }
1462 }
1463
1464 template <typename TypeTag>
1466 addReservoirSourceTerms(GlobalEqVector& residual,
1467 const std::vector<typename SparseMatrixAdapter::MatrixBlock*>& diagMatAddress) const
1468 {
1469 // NB this loop may write multiple times to the same element
1470 // if a cell is perforated by more than one well, so it should
1471 // not be OpenMP-parallelized.
1472 for (const auto& well : well_container_) {
1473 if (!well->isOperableAndSolvable() && !well->wellIsStopped()) {
1474 continue;
1475 }
1476 const auto& cells = well->cells();
1477 const auto& rates = well->connectionRates();
1478 for (unsigned perfIdx = 0; perfIdx < rates.size(); ++perfIdx) {
1479 unsigned cellIdx = cells[perfIdx];
1480 auto rate = rates[perfIdx];
1481 rate *= -1.0;
1482 VectorBlockType res(0.0);
1483 using MatrixBlockType = typename SparseMatrixAdapter::MatrixBlock;
1484 MatrixBlockType bMat(0.0);
1485 simulator_.model().linearizer().setResAndJacobi(res, bMat, rate);
1486 residual[cellIdx] += res;
1487 *diagMatAddress[cellIdx] += bMat;
1488 }
1489 }
1490 }
1491
1492
1493 template<typename TypeTag>
1494 void
1497 {
1498 int nw = this->numLocalWellsEnd();
1499 int rdofs = local_num_cells_;
1500 for (int i = 0; i < nw; ++i) {
1501 int wdof = rdofs + i;
1502 jacobian.entry(wdof,wdof) = 1.0;// better scaling ?
1503 }
1504 const auto wellconnections = this->getMaxWellConnections();
1505 for (int i = 0; i < nw; ++i) {
1506 const auto& perfcells = wellconnections[i];
1507 for (int perfcell : perfcells) {
1508 int wdof = rdofs + i;
1509 jacobian.entry(wdof, perfcell) = 0.0;
1510 jacobian.entry(perfcell, wdof) = 0.0;
1511 }
1512 }
1513 }
1514
1515
1516 template<typename TypeTag>
1517 void
1520 {
1521 DeferredLogger local_deferredLogger;
1523 {
1524 for (const auto& well : well_container_) {
1525 const auto& cells = well->cells();
1526 x_local_.resize(cells.size());
1527
1528 for (size_t i = 0; i < cells.size(); ++i) {
1529 x_local_[i] = x[cells[i]];
1530 }
1531 well->recoverWellSolutionAndUpdateWellState(simulator_, x_local_,
1532 this->groupStateHelper(), this->wellState(), local_deferredLogger);
1533 }
1534 }
1535 OPM_END_PARALLEL_TRY_CATCH_LOG(local_deferredLogger,
1536 "recoverWellSolutionAndUpdateWellState() failed: ",
1537 this->terminal_output_, simulator_.vanguard().grid().comm());
1538 }
1539
1540
1541 template<typename TypeTag>
1542 void
1544 recoverWellSolutionAndUpdateWellStateDomain(const BVector& x, const int domainIdx)
1545 {
1546 if (!nldd_) {
1547 OPM_THROW(std::logic_error, "Attempt to call NLDD method without a NLDD solver");
1548 }
1549
1550 return nldd_->recoverWellSolutionAndUpdateWellState(x, domainIdx);
1551 }
1552
1553
1554 template<typename TypeTag>
1557 getWellConvergence(const std::vector<Scalar>& B_avg, bool checkWellGroupControlsAndNetwork) const
1558 {
1559
1560 DeferredLogger local_deferredLogger;
1561 // Get global (from all processes) convergence report.
1562 ConvergenceReport local_report;
1563 const int iterationIdx = simulator_.model().newtonMethod().numIterations();
1564 for (const auto& well : well_container_) {
1565 if (well->isOperableAndSolvable() || well->wellIsStopped()) {
1566 local_report += well->getWellConvergence(
1567 this->groupStateHelper(), B_avg, local_deferredLogger,
1568 iterationIdx > param_.strict_outer_iter_wells_);
1569 } else {
1570 ConvergenceReport report;
1571 using CR = ConvergenceReport;
1572 report.setWellFailed({CR::WellFailure::Type::Unsolvable, CR::Severity::Normal, -1, well->name()});
1573 local_report += report;
1574 }
1575 }
1576
1577 const Opm::Parallel::Communication comm = grid().comm();
1578 DeferredLogger global_deferredLogger = gatherDeferredLogger(local_deferredLogger, comm);
1579 ConvergenceReport report = gatherConvergenceReport(local_report, comm);
1580
1581 if (checkWellGroupControlsAndNetwork) {
1582 // the well_group_control_changed info is already communicated
1583 report.setWellGroupTargetsViolated(this->lastReport().well_group_control_changed);
1584 report.setNetworkNotYetBalancedForceAnotherNewtonIteration(network_needs_more_balancing_force_another_newton_iteration_);
1585 }
1586
1587 if (this->terminal_output_) {
1588 global_deferredLogger.logMessages();
1589
1590 // Log debug messages for NaN or too large residuals.
1591 for (const auto& f : report.wellFailures()) {
1592 if (f.severity() == ConvergenceReport::Severity::NotANumber) {
1593 OpmLog::debug("NaN residual found with phase " + std::to_string(f.phase()) + " for well " + f.wellName());
1594 } else if (f.severity() == ConvergenceReport::Severity::TooLarge) {
1595 OpmLog::debug("Too large residual found with phase " + std::to_string(f.phase()) + " for well " + f.wellName());
1596 }
1597 }
1598 }
1599 return report;
1600 }
1601
1602
1603
1604
1605
1606 template<typename TypeTag>
1607 void
1609 calculateExplicitQuantities(DeferredLogger& deferred_logger) const
1610 {
1611 // TODO: checking isOperableAndSolvable() ?
1612 for (auto& well : well_container_) {
1613 well->calculateExplicitQuantities(simulator_, this->groupStateHelper(), deferred_logger);
1614 }
1615 }
1616
1617
1618
1619
1620
1621 template<typename TypeTag>
1622 bool
1624 updateWellControls(DeferredLogger& deferred_logger)
1625 {
1626 OPM_TIMEFUNCTION();
1627 if (!this->wellsActive()) {
1628 return false;
1629 }
1630 const int episodeIdx = simulator_.episodeIndex();
1631 const int iterationIdx = simulator_.model().newtonMethod().numIterations();
1632 const auto& comm = simulator_.vanguard().grid().comm();
1633 size_t iter = 0;
1634 bool changed_well_group = false;
1635 const Group& fieldGroup = this->schedule().getGroup("FIELD", episodeIdx);
1636 // Check group individual constraints.
1637 // iterate a few times to make sure all constraints are honored
1638 const std::size_t max_iter = param_.well_group_constraints_max_iterations_;
1639 while(!changed_well_group && iter < max_iter) {
1640 changed_well_group = updateGroupControls(fieldGroup, deferred_logger, episodeIdx, iterationIdx);
1641
1642 // Check wells' group constraints and communicate.
1643 bool changed_well_to_group = false;
1644 {
1645 OPM_TIMEBLOCK(UpdateWellControls);
1646 // For MS Wells a linear solve is performed below and the matrix might be singular.
1647 // We need to communicate the exception thrown to the others and rethrow.
1649 for (const auto& well : well_container_) {
1651 const bool changed_well = well->updateWellControl(
1652 simulator_, mode, this->groupStateHelper(), this->wellState(), deferred_logger
1653 );
1654 if (changed_well) {
1655 changed_well_to_group = changed_well || changed_well_to_group;
1656 }
1657 }
1658 OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel: updating well controls failed: ",
1659 simulator_.gridView().comm());
1660 }
1661
1662 changed_well_to_group = comm.sum(static_cast<int>(changed_well_to_group));
1663 if (changed_well_to_group) {
1664 updateAndCommunicate(episodeIdx, iterationIdx, deferred_logger);
1665 changed_well_group = true;
1666 }
1667
1668 // Check individual well constraints and communicate.
1669 bool changed_well_individual = false;
1670 {
1671 // For MS Wells a linear solve is performed below and the matrix might be singular.
1672 // We need to communicate the exception thrown to the others and rethrow.
1674 for (const auto& well : well_container_) {
1676 const bool changed_well = well->updateWellControl(
1677 simulator_, mode, this->groupStateHelper(), this->wellState(), deferred_logger
1678 );
1679 if (changed_well) {
1680 changed_well_individual = changed_well || changed_well_individual;
1681 }
1682 }
1683 OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel: updating well controls failed: ",
1684 simulator_.gridView().comm());
1685 }
1686
1687 changed_well_individual = comm.sum(static_cast<int>(changed_well_individual));
1688 if (changed_well_individual) {
1689 updateAndCommunicate(episodeIdx, iterationIdx, deferred_logger);
1690 changed_well_group = true;
1691 }
1692 iter++;
1693 }
1694
1695 // update wsolvent fraction for REIN wells
1696 this->updateWsolvent(fieldGroup, episodeIdx, this->nupcolWellState());
1697
1698 return changed_well_group;
1699 }
1700
1701
1702 template<typename TypeTag>
1703 void
1705 updateAndCommunicate(const int reportStepIdx,
1706 const int iterationIdx,
1707 DeferredLogger& deferred_logger)
1708 {
1709 this->updateAndCommunicateGroupData(reportStepIdx,
1710 iterationIdx,
1711 param_.nupcol_group_rate_tolerance_,
1712 /*update_wellgrouptarget*/ true,
1713 deferred_logger);
1714
1715 // updateWellStateWithTarget might throw for multisegment wells hence we
1716 // have a parallel try catch here to thrown on all processes.
1718 // if a well or group change control it affects all wells that are under the same group
1719 for (const auto& well : well_container_) {
1720 // We only want to update wells under group-control here
1721 const auto& ws = this->wellState().well(well->indexOfWell());
1722 if (ws.production_cmode == Well::ProducerCMode::GRUP ||
1723 ws.injection_cmode == Well::InjectorCMode::GRUP)
1724 {
1725 well->updateWellStateWithTarget(
1726 simulator_, this->groupStateHelper(), this->wellState(), deferred_logger
1727 );
1728 }
1729 }
1730 OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel::updateAndCommunicate failed: ",
1731 simulator_.gridView().comm())
1732 this->updateAndCommunicateGroupData(reportStepIdx,
1733 iterationIdx,
1734 param_.nupcol_group_rate_tolerance_,
1735 /*update_wellgrouptarget*/ true,
1736 deferred_logger);
1737 }
1738
1739 template<typename TypeTag>
1740 bool
1742 updateGroupControls(const Group& group,
1743 DeferredLogger& deferred_logger,
1744 const int reportStepIdx,
1745 const int iterationIdx)
1746 {
1747 OPM_TIMEFUNCTION();
1748 bool changed = false;
1749 // restrict the number of group switches but only after nupcol iterations.
1750 const int nupcol = this->schedule()[reportStepIdx].nupcol();
1751 const int max_number_of_group_switches = param_.max_number_of_group_switches_;
1752 const bool update_group_switching_log = iterationIdx >= nupcol;
1753 const bool changed_hc = this->checkGroupHigherConstraints(group, deferred_logger, reportStepIdx, max_number_of_group_switches, update_group_switching_log);
1754 if (changed_hc) {
1755 changed = true;
1756 updateAndCommunicate(reportStepIdx, iterationIdx, deferred_logger);
1757 }
1758
1759 bool changed_individual =
1761 updateGroupIndividualControl(group,
1762 reportStepIdx,
1763 max_number_of_group_switches,
1764 update_group_switching_log,
1765 this->switched_inj_groups_,
1766 this->switched_prod_groups_,
1767 this->closed_offending_wells_,
1768 this->groupState(),
1769 this->wellState(),
1770 deferred_logger);
1771
1772 if (changed_individual) {
1773 changed = true;
1774 updateAndCommunicate(reportStepIdx, iterationIdx, deferred_logger);
1775 }
1776 // call recursively down the group hierarchy
1777 for (const std::string& groupName : group.groups()) {
1778 bool changed_this = updateGroupControls( this->schedule().getGroup(groupName, reportStepIdx), deferred_logger, reportStepIdx,iterationIdx);
1779 changed = changed || changed_this;
1780 }
1781 return changed;
1782 }
1783
1784 template<typename TypeTag>
1785 void
1787 updateWellTestState(const double simulationTime, WellTestState& wellTestState)
1788 {
1789 OPM_TIMEFUNCTION();
1790 DeferredLogger local_deferredLogger;
1791 for (const auto& well : well_container_) {
1792 const auto& wname = well->name();
1793 const auto wasClosed = wellTestState.well_is_closed(wname);
1794 well->checkWellOperability(simulator_,
1795 this->wellState(),
1796 this->groupStateHelper(),
1797 local_deferredLogger);
1798 const bool under_zero_target =
1799 well->wellUnderZeroGroupRateTarget(this->groupStateHelper(), local_deferredLogger);
1800 well->updateWellTestState(this->wellState().well(wname),
1801 simulationTime,
1802 /*writeMessageToOPMLog=*/ true,
1803 under_zero_target,
1804 wellTestState,
1805 local_deferredLogger);
1806
1807 if (!wasClosed && wellTestState.well_is_closed(wname)) {
1808 this->closed_this_step_.insert(wname);
1809
1810 // maybe open a new well
1811 const WellEconProductionLimits& econ_production_limits = well->wellEcl().getEconLimits();
1812 if (econ_production_limits.validFollowonWell()) {
1813 const auto episode_idx = simulator_.episodeIndex();
1814 const auto follow_on_well = econ_production_limits.followonWell();
1815 if (!this->schedule().hasWell(follow_on_well, episode_idx)) {
1816 const auto msg = fmt::format("Well {} was closed. But the given follow on well {} does not exist."
1817 "The simulator continues without opening a follow on well.",
1818 wname, follow_on_well);
1819 local_deferredLogger.warning(msg);
1820 }
1821 auto& ws = this->wellState().well(follow_on_well);
1822 const bool success = ws.updateStatus(WellStatus::OPEN);
1823 if (success) {
1824 const auto msg = fmt::format("Well {} was closed. The follow on well {} opens instead.", wname, follow_on_well);
1825 local_deferredLogger.info(msg);
1826 } else {
1827 const auto msg = fmt::format("Well {} was closed. The follow on well {} is already open.", wname, follow_on_well);
1828 local_deferredLogger.warning(msg);
1829 }
1830 }
1831
1832 }
1833 }
1834
1835 for (const auto& [group_name, to] : this->closed_offending_wells_) {
1836 if (this->hasOpenLocalWell(to.second) &&
1837 !this->wasDynamicallyShutThisTimeStep(to.second))
1838 {
1839 wellTestState.close_well(to.second,
1840 WellTestConfig::Reason::GROUP,
1841 simulationTime);
1842 this->updateClosedWellsThisStep(to.second);
1843 const std::string msg =
1844 fmt::format("Procedure on exceeding {} limit is WELL for group {}. "
1845 "Well {} is {}.",
1846 to.first,
1847 group_name,
1848 to.second,
1849 "shut");
1850 local_deferredLogger.info(msg);
1851 }
1852 }
1853
1854 const Opm::Parallel::Communication comm = grid().comm();
1855 DeferredLogger global_deferredLogger =
1856 gatherDeferredLogger(local_deferredLogger, comm);
1857
1858 if (this->terminal_output_) {
1859 global_deferredLogger.logMessages();
1860 }
1861 }
1862
1863
1864 template<typename TypeTag>
1865 void
1867 const WellState<Scalar, IndexTraits>& well_state_copy,
1868 std::string& exc_msg,
1869 ExceptionType::ExcEnum& exc_type,
1870 DeferredLogger& deferred_logger)
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, deferred_logger);
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(), deferred_logger
1976 );
1977 well->updatePrimaryVariables(this->groupStateHelper(), deferred_logger);
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(), deferred_logger
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(deferred_logger);
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(), deferred_logger);
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)
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)
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:295
void init()
Definition: BlackoilWellModel_impl.hpp:161
const Simulator & simulator() const
Definition: BlackoilWellModel.hpp:371
std::vector< Scalar > depth_
Definition: BlackoilWellModel.hpp:478
std::size_t global_num_cells_
Definition: BlackoilWellModel.hpp:474
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilWellModel.hpp:107
void initWellContainer(const int reportStepIdx) override
Definition: BlackoilWellModel_impl.hpp:180
void beginReportStep(const int time_step)
Definition: BlackoilWellModel_impl.hpp:197
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:798
const Grid & grid() const
Definition: BlackoilWellModel.hpp:368
const SimulatorReportSingle & lastReport() const
Definition: BlackoilWellModel_impl.hpp:674
void addWellContributions(SparseMatrixAdapter &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1434
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: BlackoilWellModel.hpp:290
bool empty() const
Definition: BlackoilWellModel.hpp:335
void computeTotalRatesForDof(RateVector &rate, unsigned globalIdx) const
Definition: BlackoilWellModel_impl.hpp:761
void beginTimeStep()
Definition: BlackoilWellModel_impl.hpp:331
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: BlackoilWellModel.hpp:108
void calculateExplicitQuantities(DeferredLogger &deferred_logger) const
Definition: BlackoilWellModel_impl.hpp:1609
void updatePrimaryVariables(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:2060
void initializeLocalWellStructure(const int reportStepIdx, const bool enableWellPIScaling)
Definition: BlackoilWellModel_impl.hpp:250
Dune::BlockVector< VectorBlockType > BVector
Definition: BlackoilWellModel.hpp:130
BlackoilWellModel(Simulator &simulator)
Definition: BlackoilWellModel_impl.hpp:74
void wellTesting(const int timeStepIdx, const double simulationTime, DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:594
typename FluidSystem::IndexTraitsType IndexTraits
Definition: BlackoilWellModel.hpp:114
std::size_t local_num_cells_
Definition: BlackoilWellModel.hpp:476
bool alternative_well_rate_init_
Definition: BlackoilWellModel.hpp:479
void timeStepSucceeded(const double simulationTime, const double dt)
Definition: BlackoilWellModel_impl.hpp:684
Simulator & simulator_
Definition: BlackoilWellModel.hpp:448
void createWellContainer(const int report_step) override
Definition: BlackoilWellModel_impl.hpp:841
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:343
void endReportStep()
Definition: BlackoilWellModel_impl.hpp:657
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:75
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
bool updateWellControl(const Simulator &simulator, const IndividualOrGroup iog, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:189
virtual void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const =0
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
Opm::DeferredLogger gatherDeferredLogger(const Opm::DeferredLogger &local_deferredlogger, Parallel::Communication communicator)
Create a global log combining local logs.
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