BlackoilWellModel.hpp
Go to the documentation of this file.
1/*
2 Copyright 2016 SINTEF ICT, Applied Mathematics.
3 Copyright 2016 - 2017 Statoil ASA.
4 Copyright 2017 Dr. Blatt - HPC-Simulation-Software & Services
5 Copyright 2016 - 2018 IRIS 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_HEADER_INCLUDED
24#define OPM_BLACKOILWELLMODEL_HEADER_INCLUDED
25
26#include <dune/common/fmatrix.hh>
27
28#include <dune/istl/bcrsmatrix.hh>
29#include <dune/istl/matrixmatrix.hh>
30
31#include <opm/common/OpmLog/OpmLog.hpp>
32
33#include <opm/input/eclipse/Schedule/Group/Group.hpp>
34#include <opm/input/eclipse/Schedule/Group/GuideRate.hpp>
35#include <opm/input/eclipse/Schedule/Schedule.hpp>
36#include <opm/input/eclipse/Schedule/Well/WellTestState.hpp>
37
38#include <opm/material/densead/Math.hpp>
39
42
44
47
49
76
77#include <cstddef>
78#include <map>
79#include <memory>
80#include <optional>
81#include <string>
82#include <tuple>
83#include <vector>
84
85namespace Opm {
86
87template<class Scalar> class BlackoilWellModelNldd;
88template<class T, template <typename, typename...> class Storage> class SparseTable;
89
90#if COMPILE_GPU_BRIDGE
91template<class Scalar> class WellContributions;
92#endif
93
95 template<typename TypeTag>
96 class BlackoilWellModel : public WellConnectionAuxiliaryModule<TypeTag, BlackoilWellModel<TypeTag>>
97 , public BlackoilWellModelGeneric<GetPropType<TypeTag, Properties::Scalar>,
98 typename GetPropType<TypeTag, Properties::FluidSystem>::IndexTraitsType>
99 {
100 public:
101 // --------- Types ---------
113
115 using IndexTraits = typename FluidSystem::IndexTraitsType;
117
119
120 static const int numEq = Indices::numEq;
121 static const int solventSaturationIdx = Indices::solventSaturationIdx;
122 static constexpr bool has_solvent_ = getPropValue<TypeTag, Properties::EnableSolvent>();
123 static constexpr bool has_polymer_ = getPropValue<TypeTag, Properties::EnablePolymer>();
124 static constexpr EnergyModules energyModuleType_ = getPropValue<TypeTag, Properties::EnergyModuleType>();
125 static constexpr bool has_energy_ = (energyModuleType_ == EnergyModules::FullyImplicitThermal);
126 static constexpr bool has_micp_ = Indices::enableMICP;
127 static constexpr bool has_geochem_ = getPropValue<TypeTag, Properties::EnableGeochemistry>();
128 static constexpr bool has_bioeffects_ = getPropValue<TypeTag, Properties::EnableBioeffects>();
129
130 // TODO: where we should put these types, WellInterface or Well Model?
131 // or there is some other strategy, like TypeTag
132 using VectorBlockType = Dune::FieldVector<Scalar, numEq>;
133 using BVector = Dune::BlockVector<VectorBlockType>;
134
135 using PolymerModule = BlackOilPolymerModule<TypeTag, has_polymer_>;
136 using BioeffectsModule = BlackOilBioeffectsModule<TypeTag, has_bioeffects_>;
137
138 // For the conversion between the surface volume rate and reservoir voidage rate
139 using RateConverterType = RateConverter::
140 SurfaceToReservoirVoidage<FluidSystem, std::vector<int> >;
141
142 // For computing average pressured used by gpmaint
143 using AverageRegionalPressureType = RegionAverageCalculator::
144 AverageRegionalPressure<FluidSystem, std::vector<int> >;
145
147
148 void init();
149 void initWellContainer(const int reportStepIdx) override;
150
152 {
153 OPM_TIMEBLOCK(beginEpsiode);
154 beginReportStep(simulator_.episodeIndex());
155 }
156
157 void beginTimeStep();
158
160 {
161 OPM_TIMEBLOCK(beginIteration);
162 assemble(simulator_.timeStepSize());
163 }
164
166 { }
167
169 {
170 OPM_TIMEBLOCK(endTimeStep);
171 timeStepSucceeded(simulator_.time(), simulator_.timeStepSize());
172 }
173
175 {
177 }
178
180 unsigned globalIdx) const;
181
182 template <class Context>
184 const Context& context,
185 unsigned spaceIdx,
186 unsigned timeIdx) const;
187
188
189 using WellInterfacePtr = std::unique_ptr<WellInterface<TypeTag>>;
190
192 void initFromRestartFile(const RestartValue& restartValues)
193 {
194 initFromRestartFile(restartValues,
195 this->simulator_.vanguard().transferWTestState(),
196 grid().size(0),
197 this->simulator_.vanguard().enableDistributedWells());
198 }
199
201 void prepareDeserialize(const int report_step)
202 {
203 prepareDeserialize(report_step, grid().size(0),
204 this->simulator_.vanguard().enableDistributedWells());
205 }
206
207 data::Wells wellData() const
208 {
209 auto wsrpt = this->wellState()
210 .report(this->simulator_.vanguard().globalCell().data(),
211 [this](const int well_index)
212 { return this->wasDynamicallyShutThisTimeStep(well_index); },
213 this->rsConstInfo());
214
216 .assignWellGuideRates(wsrpt, this->reportStepIndex());
217
218 this->assignWellTracerRates_(wsrpt);
219
220 if constexpr (has_geochem_) {
221 this->assignWellSpeciesRates_(wsrpt);
222 }
223
224 if (const auto& rspec = eclState().runspec();
225 rspec.co2Storage() || rspec.h2Storage())
226 {
227 // The gas reference density (surface condition) is the
228 // same for all PVT regions in CO2STORE/H2STORE runs so,
229 // for simplicity, we use region zero (0) here.
230
231 this->assignMassGasRate(wsrpt, FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, 0));
232 }
233
234 this->assignWellTargets(wsrpt);
235
236 this->assignDynamicWellStatus(wsrpt);
237
238 // Assigning (a subset of the) property values in shut
239 // connections should be the last step of wellData().
240 this->assignShutConnections(wsrpt, this->reportStepIndex());
241
242 return wsrpt;
243 }
244
245 data::WellBlockAveragePressures wellBlockAveragePressures() const
246 {
247 return this->wbp_.computeWellBlockAveragePressures(this->gravity_);
248 }
249
250#if COMPILE_GPU_BRIDGE
251 // accumulate the contributions of all Wells in the WellContributions object
252 void getWellContributions(WellContributions<Scalar>& x) const;
253#endif
254
255 // Check if well equations is converged.
256 ConvergenceReport getWellConvergence(const std::vector<Scalar>& B_avg, const bool checkWellGroupControlsAndNetwork = false) const;
257
258 const SimulatorReportSingle& lastReport() const;
259
260 void addWellContributions(SparseMatrixAdapter& jacobian) const;
261
262 // add source from wells to the reservoir matrix
264 const std::vector<typename SparseMatrixAdapter::MatrixBlock*>& diagMatAddress) const;
265
266 // called at the beginning of a report step
267 void beginReportStep(const int time_step);
268
272 void calculateExplicitQuantities() const;
273
276 void prepareTimeStep(DeferredLogger& deferred_logger);
277
278 bool
279 updateWellControls(DeferredLogger& deferred_logger);
280
281 void updateAndCommunicate(const int reportStepIdx);
282
283 bool updateGroupControls(const Group& group,
284 DeferredLogger& deferred_logger,
285 const int reportStepIdx);
286
287 const WellInterface<TypeTag>& getWell(const std::string& well_name) const;
288
289 using PressureMatrix = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, 1, 1>>;
290
292 const BVector& weights,
293 const bool use_well_weights) const;
296 const BVector& weights,
297 const bool use_well_weights,
298 const int domainIndex) const
299 {
300 if (!nldd_) {
301 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
302 }
303 return nldd_->addWellPressureEquations(jacobian,
304 weights,
305 use_well_weights,
306 domainIndex);
307 }
308
310 const std::vector<WellInterfacePtr>& localNonshutWells() const
311 {
312 return well_container_;
313 }
314
316 {
317 if (!nldd_) {
318 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
319 }
320 return nldd_->well_local_cells();
321 }
322
323 const std::map<std::string, int>& well_domain() const
324 {
325 if (!nldd_) {
326 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
327 }
328
329 return nldd_->well_domain();
330 }
331
332 auto begin() const { return well_container_.begin(); }
333 auto end() const { return well_container_.end(); }
334 bool empty() const { return well_container_.empty(); }
335
338
341
342 int compressedIndexForInterior(int cartesian_cell_idx) const override
343 {
344 return simulator_.vanguard().compressedIndexForInterior(cartesian_cell_idx);
345 }
346
347 int compressedIndexForInteriorLGR(const std::string& lgr_tag, const Connection& conn) const override
348 {
349 return simulator_.vanguard().compressedIndexForInteriorLGR(lgr_tag, conn);
350 }
351
352 // using the solution x to recover the solution xw for wells and applying
353 // xw to update Well State
355
356 // using the solution x to recover the solution xw for wells and applying
357 // xw to update Well State
359 const int domainIdx);
360 // Update cellRates_ with contributions from all wells
361 void updateCellRates();
362
363 // Update cellRates_ with contributions from wells in a specific domain
364 void updateCellRatesForDomain(int domainIndex,
365 const std::map<std::string, int>& well_domain_map);
366
367 const Grid& grid() const
368 { return simulator_.vanguard().grid(); }
369
370 const Simulator& simulator() const
371 { return simulator_; }
373 { return simulator_; }
374
376 const BlackoilWellModelNetwork<TypeTag>& network() const { return network_; }
377
378 std::vector<WellInterfacePtr>& wellContainer() { return well_container_; }
379 const std::vector<WellInterfacePtr>& wellContainer() const { return well_container_; }
380
382 { nldd_ = mod; }
383
385 { return guide_rate_handler_; }
386
388 { return guide_rate_handler_; }
389
390 // === Reservoir Coupling ===
391
395
403 void updateGuideRates(const int report_step_idx,
404 const double sim_time)
405 {
406 this->guide_rate_handler_.updateGuideRates(
407 report_step_idx, sim_time, this->wellState(), this->groupState()
408 );
409 }
410
412 bool isReservoirCouplingMaster() const { return rescoup_.isMaster(); }
413
415 bool isReservoirCouplingSlave() const { return rescoup_.isSlave(); }
416
422 bool isReservoirCouplingMasterGroup(const std::string& group_name) const {
423 return rescoup_.isMasterGroup(group_name);
424 }
425
429 return rescoup_.master();
430 }
432 return rescoup_.master();
433 }
434
438 return rescoup_.slave();
439 }
441 return rescoup_.slave();
442 }
443
444#ifdef RESERVOIR_COUPLING_ENABLED
445 void setReservoirCouplingMaster(ReservoirCouplingMaster<Scalar>* master)
446 {
447 rescoup_.setMaster(master);
448 this->guide_rate_handler_.setReservoirCouplingMaster(master);
449 this->groupStateHelper().setReservoirCouplingMaster(master);
450 }
451 void setReservoirCouplingSlave(ReservoirCouplingSlave<Scalar>* slave)
452 {
453 rescoup_.setSlave(slave);
454 this->guide_rate_handler_.setReservoirCouplingSlave(slave);
455 this->groupStateHelper().setReservoirCouplingSlave(slave);
456 }
457#endif
458
459 bool updateWellControlsAndNetwork(const bool mandatory_network_balance,
460 const double dt,
461 DeferredLogger& local_deferredLogger);
462
463 // TODO: finding a better naming
464 void assembleWellEqWithoutIteration(const double dt);
465
466 const std::vector<Scalar>& B_avg() const
467 { return B_avg_; }
468
469 const ModelParameters& param() const
470 { return param_; }
471
477
478#ifdef RESERVOIR_COUPLING_ENABLED
482 BlackoilWellModelRescoup<TypeTag>& rescoupHelper() { return rescoupHelper_; }
483#endif
484
485
486 template<class FluidState, class SingleWellState>
487 static Scalar computeTemperatureWeightFactor(const int perf_index, const int np, const FluidState& fs, const SingleWellState& ws)
488 {
489 const auto& perf_phase_rate = ws.perf_data.phase_rates;
490 // we only have one temperature pr cell any phaseIdx will do
491 Scalar cellTemperatures = fs.temperature(/*phaseIdx*/0).value();
492 Scalar weight_factor = 0.0;
493 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
494 if (!FluidSystem::phaseIsActive(phaseIdx)) {
495 continue;
496 }
497 Scalar cellInternalEnergy = fs.enthalpy(phaseIdx).value() -
498 fs.pressure(phaseIdx).value() / fs.density(phaseIdx).value();
499 Scalar cellBinv = fs.invB(phaseIdx).value();
500 Scalar cellDensity = fs.density(phaseIdx).value();
501 Scalar perfPhaseRate = perf_phase_rate[perf_index*np + phaseIdx];
502 weight_factor += cellDensity * (perfPhaseRate / cellBinv) * (cellInternalEnergy / cellTemperatures);
503 }
504 return (std::abs(weight_factor) + 1e-13);
505 }
506
507 protected:
509
510 // a vector of all the wells.
511 std::vector<WellInterfacePtr> well_container_{};
512
513 std::vector<bool> is_cell_perforated_{};
514
515 void initializeWellState(const int timeStepIdx);
516
517 // create the well container
518 void createWellContainer(const int report_step) override;
519
521 createWellPointer(const int wellID,
522 const int report_step) const;
523
524 template <typename WellType>
525 std::unique_ptr<WellType>
526 createTypedWellPointer(const int wellID,
527 const int time_step) const;
528
529 WellInterfacePtr createWellForWellTest(const std::string& well_name,
530 const int report_step,
531 DeferredLogger& deferred_logger) const;
532
534 std::size_t global_num_cells_{};
535 // the number of the cells in the local grid
536 std::size_t local_num_cells_{};
538 std::vector<Scalar> depth_{};
540 std::unique_ptr<RateConverterType> rateConverter_{};
541 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>> regionalAveragePressureCalculator_{};
542
546
547 // A flag to tell the convergence report whether we need to take another newton step
549
550 std::vector<Scalar> B_avg_{};
551
552 const EquilGrid& equilGrid() const
553 { return simulator_.vanguard().equilGrid(); }
554
555 const EclipseState& eclState() const
556 { return simulator_.vanguard().eclState(); }
557
558 // compute the well fluxes and assemble them in to the reservoir equations as source terms
559 // and in the well equations.
560 void assemble(const double dt);
561
562 // well controls and network pressures affect each other and are solved in an iterative manner.
563 // the function handles one iteration of updating well controls and network pressures.
564 // it is possible to decouple the update of well controls and network pressures further.
565 // the returned two booleans are {continue_due_to_network, well_group_control_changed}, respectively
566 std::tuple<bool, bool, Scalar> updateWellControlsAndNetworkIteration(const bool mandatory_network_balance,
567 const bool relax_network_tolerance,
568 const bool optimize_gas_lift,
569 const double dt,
570 DeferredLogger& local_deferredLogger);
571
580 void initializeLocalWellStructure(const int reportStepIdx,
581 const bool enableWellPIScaling);
582
586 void initializeGroupStructure(const int reportStepIdx);
587
588 // called at the end of a time step
589 void timeStepSucceeded(const double simulationTime, const double dt);
590
591 // called at the end of a report step
592 void endReportStep();
593
594 // setting the well_solutions_ based on well_state.
596
598
599 void computePotentials(const std::size_t widx,
600 const WellState<Scalar, IndexTraits>& well_state_copy,
601 std::string& exc_msg,
602 ExceptionType::ExcEnum& exc_type) override;
603
604 const std::vector<Scalar>& wellPerfEfficiencyFactors() const;
605
606 void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger& deferred_logger) override;
607 void calculateProductivityIndexValues(DeferredLogger& deferred_logger) override;
609 DeferredLogger& deferred_logger);
610
611 // The number of conservation quantities.
612 int numConservationQuantities() const;
613
614 int reportStepIndex() const;
615
616 void assembleWellEq(const double dt);
617
618 void prepareWellsBeforeAssembling(const double dt);
619
621
622 void extractLegacyDepth_();
623
625 void updateWellTestState(const double simulationTime, WellTestState& wellTestState);
626
627 void wellTesting(const int timeStepIdx, const double simulationTime, DeferredLogger& deferred_logger);
628
629 void calcResvCoeff(const int fipnum,
630 const int pvtreg,
631 const std::vector<Scalar>& production_rates,
632 std::vector<Scalar>& resv_coeff) const override;
633
634 void calcInjResvCoeff(const int fipnum,
635 const int pvtreg,
636 std::vector<Scalar>& resv_coeff) const override;
637
639
640 private:
641 // Private helper methods (alphabetical order)
642 // --------------------------------------------
643
644 void assignWellSpeciesRates_(data::Wells& wsrpt) const;
645 void assignWellTracerRates_(data::Wells& wsrpt) const;
646
655 bool isRescoupCoupledNetworkParticipant_() const;
656
664 bool isRescoupMasterCoupledNetworkIteration_() const;
665
673 bool isRescoupSlaveCoupledNetworkIteration_() const;
674
681 bool isRescoupSlaveOnSyncStepFirstSubstep_() const;
682
691 bool isRescoupSlaveConnectedToMasterNetwork_() const;
692
699 bool maybeSendSlaveGroupFlowToMaster_(const int reportStepIdx);
700
705 void sendSlaveNetworkLoopTerminationSignal_();
706
709 bool shouldDoPreStepNetworkRebalance_(const int episodeIdx) const;
710
712 void updateNetworkActiveState_();
713
716#ifdef RESERVOIR_COUPLING_ENABLED
717 BlackoilWellModelRescoup<TypeTag> rescoupHelper_;
718#endif
719 BlackoilWellModelNldd<TypeTag>* nldd_ = nullptr;
720
721 // These members are used to avoid reallocation in specific functions
722 // instead of using local variables.
723 // Their state is not relevant between function calls, so they can
724 // (and must) be mutable, as the functions using them are const.
725 mutable BVector x_local_;
726
727 // Store cell rates after assembling to avoid iterating all wells and connections for every element
728 std::map<int, RateVector> cellRates_;
729
730 [[nodiscard]] auto rsConstInfo() const
731 -> typename WellState<Scalar,IndexTraits>::RsConstInfo;
732 };
733
734} // namespace Opm
735
736#include "BlackoilWellModel_impl.hpp"
737
738#endif // OPM_BLACKOILWELLMODEL_HEADER_INCLUDED
Class for handling the gaslift in the blackoil well model.
Definition: BlackoilWellModelGasLift.hpp:96
Class for handling the blackoil well model.
Definition: BlackoilWellModelGeneric.hpp:98
const WellState< Scalar, IndexTraits > & wellState() const
Definition: BlackoilWellModelGeneric.hpp:157
Class for handling the guide rates in the blackoil well model.
Definition: BlackoilWellModelGuideRates.hpp:47
void assignWellGuideRates(data::Wells &wsrpt, const int reportStepIdx) const
Assign well guide rates.
Class for handling the blackoil well network model.
Definition: BlackoilWellModelNetwork.hpp:47
Class for handling the blackoil well model in a NLDD solver.
Definition: BlackoilWellModelNldd.hpp:80
Class for handling the blackoil well model.
Definition: BlackoilWellModel.hpp:99
ReservoirCoupling::Proxy< Scalar > & rescoup()
Get the reservoir coupling proxy.
Definition: BlackoilWellModel.hpp:393
void initializeGroupStructure(const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:294
bool network_needs_more_balancing_force_another_newton_iteration_
Definition: BlackoilWellModel.hpp:548
void calcResvCoeff(const int fipnum, const int pvtreg, const std::vector< Scalar > &production_rates, std::vector< Scalar > &resv_coeff) const override
Definition: BlackoilWellModel_impl.hpp:2104
ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster()
Get reference to reservoir coupling master.
Definition: BlackoilWellModel.hpp:428
void endTimeStep()
Definition: BlackoilWellModel.hpp:168
void prepareTimeStep(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1927
std::tuple< bool, bool, Scalar > updateWellControlsAndNetworkIteration(const bool mandatory_network_balance, const bool relax_network_tolerance, const bool optimize_gas_lift, const double dt, DeferredLogger &local_deferredLogger)
Definition: BlackoilWellModel_impl.hpp:1219
static constexpr bool has_bioeffects_
Definition: BlackoilWellModel.hpp:128
WellInterfacePtr createWellPointer(const int wellID, const int report_step) const
Definition: BlackoilWellModel_impl.hpp:1006
GuideRateHandler< Scalar, IndexTraits > guide_rate_handler_
Definition: BlackoilWellModel.hpp:544
const EclipseState & eclState() const
Definition: BlackoilWellModel.hpp:555
void prepareWellsBeforeAssembling(const double dt)
Definition: BlackoilWellModel_impl.hpp:1318
const std::vector< WellInterfacePtr > & localNonshutWells() const
Get list of local nonshut wells.
Definition: BlackoilWellModel.hpp:310
void prepareDeserialize(const int report_step)
Definition: BlackoilWellModel.hpp:201
void init()
Definition: BlackoilWellModel_impl.hpp:163
void setNlddAdapter(BlackoilWellModelNldd< TypeTag > *mod)
Definition: BlackoilWellModel.hpp:381
const Simulator & simulator() const
Definition: BlackoilWellModel.hpp:370
std::vector< Scalar > depth_
Definition: BlackoilWellModel.hpp:538
static constexpr bool has_geochem_
Definition: BlackoilWellModel.hpp:127
static constexpr EnergyModules energyModuleType_
Definition: BlackoilWellModel.hpp:124
std::size_t global_num_cells_
Definition: BlackoilWellModel.hpp:534
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilWellModel.hpp:108
const std::vector< Scalar > & B_avg() const
Definition: BlackoilWellModel.hpp:466
SimulatorReportSingle last_report_
Definition: BlackoilWellModel.hpp:543
static const int solventSaturationIdx
Definition: BlackoilWellModel.hpp:121
Scalar gravity_
Definition: BlackoilWellModel.hpp:537
void initWellContainer(const int reportStepIdx) override
Definition: BlackoilWellModel_impl.hpp:182
void beginReportStep(const int time_step)
Definition: BlackoilWellModel_impl.hpp:199
const WellInterface< TypeTag > & getWell(const std::string &well_name) const
Definition: BlackoilWellModel_impl.hpp:2076
bool useTightRcNetworkCoupling() const
True when tight (per-sub-iteration) reservoir-coupling network coupling is in effect....
Definition: BlackoilWellModel.hpp:475
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: BlackoilWellModel.hpp:104
void beginIteration()
Definition: BlackoilWellModel.hpp:159
bool isReservoirCouplingMaster() const
Check if this process is a reservoir coupling master.
Definition: BlackoilWellModel.hpp:412
Dune::FieldVector< Scalar, numEq > VectorBlockType
Definition: BlackoilWellModel.hpp:132
const BlackoilWellModelNetwork< TypeTag > & network() const
Definition: BlackoilWellModel.hpp:376
void initFromRestartFile(const RestartValue &restartValues)
Definition: BlackoilWellModel.hpp:192
const ModelParameters param_
Definition: BlackoilWellModel.hpp:533
void beginEpisode()
Definition: BlackoilWellModel.hpp:151
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: BlackoilWellModel.hpp:105
GetPropType< TypeTag, Properties::Grid > Grid
Definition: BlackoilWellModel.hpp:102
GetPropType< TypeTag, Properties::EquilGrid > EquilGrid
Definition: BlackoilWellModel.hpp:103
int numConservationQuantities() const
Definition: BlackoilWellModel_impl.hpp:2048
const std::map< std::string, int > & well_domain() const
Definition: BlackoilWellModel.hpp:323
const ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave() const
Definition: BlackoilWellModel.hpp:440
bool updateWellControls(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1602
void endIteration()
Definition: BlackoilWellModel.hpp:165
static constexpr bool has_polymer_
Definition: BlackoilWellModel.hpp:123
int reportStepIndex() const
Definition: BlackoilWellModel_impl.hpp:2092
void calculateProductivityIndexValues(DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1869
void extractLegacyDepth_()
Definition: BlackoilWellModel_impl.hpp:2064
void extractLegacyCellPvtRegionIndex_()
Definition: BlackoilWellModel_impl.hpp:2032
std::vector< WellInterfacePtr > & wellContainer()
Definition: BlackoilWellModel.hpp:378
void recoverWellSolutionAndUpdateWellStateDomain(const BVector &x, const int domainIdx)
Definition: BlackoilWellModel_impl.hpp:1523
void updateAverageFormationFactor()
Definition: BlackoilWellModel_impl.hpp:1975
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilWellModel.hpp:107
BlackoilWellModel(Simulator &simulator, const NewtonIterationContext &iter_ctx)
Definition: BlackoilWellModel_impl.hpp:72
bool isReservoirCouplingSlave() const
Check if this process is a reservoir coupling slave.
Definition: BlackoilWellModel.hpp:415
const ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster() const
Definition: BlackoilWellModel.hpp:431
void initializeWellState(const int timeStepIdx)
Definition: BlackoilWellModel_impl.hpp:761
const Grid & grid() const
Definition: BlackoilWellModel.hpp:367
void updatePrimaryVariables()
Definition: BlackoilWellModel_impl.hpp:2023
const ReservoirCoupling::Proxy< Scalar > & rescoup() const
Definition: BlackoilWellModel.hpp:394
BlackOilPolymerModule< TypeTag, has_polymer_ > PolymerModule
Definition: BlackoilWellModel.hpp:135
void computeWellTemperature()
Definition: BlackoilWellModel_impl.hpp:2126
static const int numEq
Definition: BlackoilWellModel.hpp:120
void addWellPressureEquations(PressureMatrix &jacobian, const BVector &weights, const bool use_well_weights) const
Definition: BlackoilWellModel_impl.hpp:1427
const SimulatorReportSingle & lastReport() const
Definition: BlackoilWellModel_impl.hpp:637
bool updateWellControlsAndNetwork(const bool mandatory_network_balance, const double dt, DeferredLogger &local_deferredLogger)
Definition: BlackoilWellModel_impl.hpp:1157
void endEpisode()
Definition: BlackoilWellModel.hpp:174
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: BlackoilWellModel.hpp:111
void addWellContributions(SparseMatrixAdapter &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1417
const EquilGrid & equilGrid() const
Definition: BlackoilWellModel.hpp:552
void assembleWellEq(const double dt)
Definition: BlackoilWellModel_impl.hpp:1306
WellInterfacePtr createWellForWellTest(const std::string &well_name, const int report_step, DeferredLogger &deferred_logger) const
Definition: BlackoilWellModel_impl.hpp:1057
GetPropType< TypeTag, Properties::Indices > Indices
Definition: BlackoilWellModel.hpp:106
int compressedIndexForInteriorLGR(const std::string &lgr_tag, const Connection &conn) const override
Definition: BlackoilWellModel.hpp:347
void calculateExplicitQuantities() const
Definition: BlackoilWellModel_impl.hpp:1587
BlackoilWellModelNetwork< TypeTag > & network()
Definition: BlackoilWellModel.hpp:375
std::vector< WellInterfacePtr > well_container_
Definition: BlackoilWellModel.hpp:511
static constexpr std::size_t pressureVarIndex
Definition: BlackoilWellModel.hpp:118
void updateAndCommunicate(const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:1682
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:724
bool addMatrixContributions() const
Definition: BlackoilWellModel.hpp:336
void beginTimeStep()
Definition: BlackoilWellModel_impl.hpp:326
std::vector< Scalar > B_avg_
Definition: BlackoilWellModel.hpp:550
const std::vector< Scalar > & wellPerfEfficiencyFactors() const
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: BlackoilWellModel.hpp:109
bool updateGroupControls(const Group &group, DeferredLogger &deferred_logger, const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:1709
std::map< std::string, std::unique_ptr< AverageRegionalPressureType > > regionalAveragePressureCalculator_
Definition: BlackoilWellModel.hpp:541
void calcInjResvCoeff(const int fipnum, const int pvtreg, std::vector< Scalar > &resv_coeff) const override
Definition: BlackoilWellModel_impl.hpp:2115
static Scalar computeTemperatureWeightFactor(const int perf_index, const int np, const FluidState &fs, const SingleWellState &ws)
Definition: BlackoilWellModel.hpp:487
GetPropType< TypeTag, Properties::GlobalEqVector > GlobalEqVector
Definition: BlackoilWellModel.hpp:110
void initializeLocalWellStructure(const int reportStepIdx, const bool enableWellPIScaling)
Definition: BlackoilWellModel_impl.hpp:248
GuideRateHandler< Scalar, IndexTraits > & guideRateHandler()
Definition: BlackoilWellModel.hpp:384
static constexpr bool has_solvent_
Definition: BlackoilWellModel.hpp:122
Dune::BlockVector< VectorBlockType > BVector
Definition: BlackoilWellModel.hpp:133
static constexpr bool has_micp_
Definition: BlackoilWellModel.hpp:126
std::unique_ptr< RateConverterType > rateConverter_
Definition: BlackoilWellModel.hpp:540
const SparseTable< int > & well_local_cells() const
Definition: BlackoilWellModel.hpp:315
void wellTesting(const int timeStepIdx, const double simulationTime, DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:560
std::vector< bool > is_cell_perforated_
Definition: BlackoilWellModel.hpp:513
ConvergenceReport getWellConvergence(const std::vector< Scalar > &B_avg, const bool checkWellGroupControlsAndNetwork=false) const
Definition: BlackoilWellModel_impl.hpp:1536
int numStrictIterations() const
Definition: BlackoilWellModel.hpp:339
typename FluidSystem::IndexTraitsType IndexTraits
Definition: BlackoilWellModel.hpp:115
void updateCellRatesForDomain(int domainIndex, const std::map< std::string, int > &well_domain_map)
Definition: BlackoilWellModel_impl.hpp:1364
data::WellBlockAveragePressures wellBlockAveragePressures() const
Definition: BlackoilWellModel.hpp:245
void assembleWellEqWithoutIteration(const double dt)
Definition: BlackoilWellModel_impl.hpp:1332
void updateCellRates()
Definition: BlackoilWellModel_impl.hpp:1352
void assemble(const double dt)
Definition: BlackoilWellModel_impl.hpp:1082
auto begin() const
Definition: BlackoilWellModel.hpp:332
std::size_t local_num_cells_
Definition: BlackoilWellModel.hpp:536
auto end() const
Definition: BlackoilWellModel.hpp:333
BlackOilBioeffectsModule< TypeTag, has_bioeffects_ > BioeffectsModule
Definition: BlackoilWellModel.hpp:136
bool alternative_well_rate_init_
Definition: BlackoilWellModel.hpp:539
Simulator & simulator()
Definition: BlackoilWellModel.hpp:372
void timeStepSucceeded(const double simulationTime, const double dt)
Definition: BlackoilWellModel_impl.hpp:647
std::unique_ptr< WellType > createTypedWellPointer(const int wellID, const int time_step) const
Definition: BlackoilWellModel_impl.hpp:1026
void computePotentials(const std::size_t widx, const WellState< Scalar, IndexTraits > &well_state_copy, std::string &exc_msg, ExceptionType::ExcEnum &exc_type) override
Definition: BlackoilWellModel_impl.hpp:1834
Simulator & simulator_
Definition: BlackoilWellModel.hpp:508
void updateGuideRates(const int report_step_idx, const double sim_time)
Update guide rates for all wells and groups.
Definition: BlackoilWellModel.hpp:403
ReservoirCoupling::Proxy< Scalar > rescoup_
Definition: BlackoilWellModel.hpp:545
void createWellContainer(const int report_step) override
Definition: BlackoilWellModel_impl.hpp:804
std::unique_ptr< WellInterface< TypeTag > > WellInterfacePtr
Definition: BlackoilWellModel.hpp:189
ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave()
Get reference to reservoir coupling slave.
Definition: BlackoilWellModel.hpp:437
data::Wells wellData() const
Definition: BlackoilWellModel.hpp:207
const ModelParameters & param() const
Definition: BlackoilWellModel.hpp:469
void updateWellTestState(const double simulationTime, WellTestState &wellTestState)
upate the wellTestState related to economic limits
Definition: BlackoilWellModel_impl.hpp:1763
void addWellPressureEquationsDomain(PressureMatrix &jacobian, const BVector &weights, const bool use_well_weights, const int domainIndex) const
Definition: BlackoilWellModel.hpp:295
bool isReservoirCouplingMasterGroup(const std::string &group_name) const
Check if a group is a reservoir coupling master group.
Definition: BlackoilWellModel.hpp:422
const GuideRateHandler< Scalar, IndexTraits > & guideRateHandler() const
Definition: BlackoilWellModel.hpp:387
void addReservoirSourceTerms(GlobalEqVector &residual, const std::vector< typename SparseMatrixAdapter::MatrixBlock * > &diagMatAddress) const
Definition: BlackoilWellModel_impl.hpp:1449
int compressedIndexForInterior(int cartesian_cell_idx) const override
get compressed index for interior cells (-1, otherwise
Definition: BlackoilWellModel.hpp:342
static constexpr bool has_energy_
Definition: BlackoilWellModel.hpp:125
void recoverWellSolutionAndUpdateWellState(const BVector &x)
Definition: BlackoilWellModel_impl.hpp:1499
void addWellPressureEquationsStruct(PressureMatrix &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1479
void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1883
const std::vector< WellInterfacePtr > & wellContainer() const
Definition: BlackoilWellModel.hpp:379
void endReportStep()
Definition: BlackoilWellModel_impl.hpp:620
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
Definition: GroupStateHelper.hpp:56
Handles computation and reporting of guide rates for wells and groups.
Definition: GuideRateHandler.hpp:46
Definition: RateConverter.hpp:70
Definition: RegionAverageCalculator.hpp:60
Thin proxy for reservoir coupling master/slave pointers.
Definition: RescoupProxy.hpp:54
Definition: ReservoirCouplingMaster.hpp:38
Definition: ReservoirCouplingSlave.hpp:40
Definition: SingleWellState.hpp:44
PerfData< Scalar > perf_data
Definition: SingleWellState.hpp:156
Definition: BlackoilWellModel.hpp:88
Definition: WellConnectionAuxiliaryModule.hpp:39
Definition: WellContributions.hpp:51
Definition: WellInterface.hpp:77
Definition: WellState.hpp:66
ExcEnum
Definition: DeferredLogger.hpp:45
Definition: blackoilbioeffectsmodules.hh:45
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
Solver parameters for the NonlinearSystemBlackOilReservoir.
Definition: BlackoilModelParameters.hpp:201
bool matrix_add_well_contributions_
Whether to add influences of wells between cells to the matrix and preconditioner matrix.
Definition: BlackoilModelParameters.hpp:315
int strict_outer_iter_wells_
Newton iteration where wells are stricly convergent.
Definition: BlackoilModelParameters.hpp:265
bool rc_network_loose_coupling_
Definition: BlackoilModelParameters.hpp:357
Context for iteration-dependent decisions in the Newton solver.
Definition: NewtonIterationContext.hpp:43
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34