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
75
76#include <cstddef>
77#include <map>
78#include <memory>
79#include <optional>
80#include <string>
81#include <tuple>
82#include <vector>
83
84namespace Opm {
85
86template<class Scalar> class BlackoilWellModelNldd;
87template<class T, template <typename, typename...> class Storage> class SparseTable;
88
89#if COMPILE_GPU_BRIDGE
90template<class Scalar> class WellContributions;
91#endif
92
94 template<typename TypeTag>
95 class BlackoilWellModel : public WellConnectionAuxiliaryModule<TypeTag, BlackoilWellModel<TypeTag>>
96 , public BlackoilWellModelGeneric<GetPropType<TypeTag, Properties::Scalar>,
97 typename GetPropType<TypeTag, Properties::FluidSystem>::IndexTraitsType>
98 {
99 public:
100 // --------- Types ---------
112
114 using IndexTraits = typename FluidSystem::IndexTraitsType;
116
118
119 static const int numEq = Indices::numEq;
120 static const int solventSaturationIdx = Indices::solventSaturationIdx;
121 static constexpr bool has_solvent_ = getPropValue<TypeTag, Properties::EnableSolvent>();
122 static constexpr bool has_polymer_ = getPropValue<TypeTag, Properties::EnablePolymer>();
123 static constexpr EnergyModules energyModuleType_ = getPropValue<TypeTag, Properties::EnergyModuleType>();
124 static constexpr bool has_energy_ = (energyModuleType_ == EnergyModules::FullyImplicitThermal);
125 static constexpr bool has_micp_ = Indices::enableMICP;
126 static constexpr bool has_geochem_ = getPropValue<TypeTag, Properties::EnableGeochemistry>();
127
128 // TODO: where we should put these types, WellInterface or Well Model?
129 // or there is some other strategy, like TypeTag
130 using VectorBlockType = Dune::FieldVector<Scalar, numEq>;
131 using BVector = Dune::BlockVector<VectorBlockType>;
132
135
136 // For the conversion between the surface volume rate and reservoir voidage rate
137 using RateConverterType = RateConverter::
138 SurfaceToReservoirVoidage<FluidSystem, std::vector<int> >;
139
140 // For computing average pressured used by gpmaint
141 using AverageRegionalPressureType = RegionAverageCalculator::
142 AverageRegionalPressure<FluidSystem, std::vector<int> >;
143
145
146 void init();
147 void initWellContainer(const int reportStepIdx) override;
148
150 {
151 OPM_TIMEBLOCK(beginEpsiode);
152 beginReportStep(simulator_.episodeIndex());
153 }
154
155 void beginTimeStep();
156
158 {
159 OPM_TIMEBLOCK(beginIteration);
160 assemble(simulator_.timeStepSize());
161 }
162
164 { }
165
167 {
168 OPM_TIMEBLOCK(endTimeStep);
169 timeStepSucceeded(simulator_.time(), simulator_.timeStepSize());
170 }
171
173 {
175 }
176
178 unsigned globalIdx) const;
179
180 template <class Context>
182 const Context& context,
183 unsigned spaceIdx,
184 unsigned timeIdx) const;
185
186
187 using WellInterfacePtr = std::unique_ptr<WellInterface<TypeTag>>;
188
190 void initFromRestartFile(const RestartValue& restartValues)
191 {
192 initFromRestartFile(restartValues,
193 this->simulator_.vanguard().transferWTestState(),
194 grid().size(0),
196 this->simulator_.vanguard().enableDistributedWells());
197 }
198
200 void prepareDeserialize(const int report_step)
201 {
202 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 {
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_; }
372
374 { nldd_ = mod; }
375
376 // === Reservoir Coupling ===
377
381
389 void updateGuideRates(const int report_step_idx,
390 const double sim_time)
391 {
392 this->guide_rate_handler_.updateGuideRates(
393 report_step_idx, sim_time, this->wellState(), this->groupState()
394 );
395 }
396
398 bool isReservoirCouplingMaster() const { return rescoup_.isMaster(); }
399
401 bool isReservoirCouplingSlave() const { return rescoup_.isSlave(); }
402
406 return rescoup_.master();
407 }
408
412 return rescoup_.slave();
413 }
414
415#ifdef RESERVOIR_COUPLING_ENABLED
416 void setReservoirCouplingMaster(ReservoirCouplingMaster<Scalar>* master)
417 {
418 rescoup_.setMaster(master);
419 this->guide_rate_handler_.setReservoirCouplingMaster(master);
420 this->groupStateHelper().setReservoirCouplingMaster(master);
421 }
422 void setReservoirCouplingSlave(ReservoirCouplingSlave<Scalar>* slave)
423 {
424 rescoup_.setSlave(slave);
425 this->guide_rate_handler_.setReservoirCouplingSlave(slave);
426 this->groupStateHelper().setReservoirCouplingSlave(slave);
427 }
428
430 void sendSlaveGroupDataToMaster();
431
433 void receiveSlaveGroupData();
434
435 void receiveGroupConstraintsFromMaster();
436 void sendMasterGroupConstraintsToSlaves();
437
446 std::optional<ReservoirCoupling::ScopedLoggerGuard>
447 setupRescoupScopedLogger(DeferredLogger& local_logger);
448#endif
449
450 bool updateWellControlsAndNetwork(const bool mandatory_network_balance,
451 const double dt,
452 DeferredLogger& local_deferredLogger);
453
454 // TODO: finding a better naming
455 void assembleWellEqWithoutIteration(const double dt);
456
457 const std::vector<Scalar>& B_avg() const
458 { return B_avg_; }
459
460 const ModelParameters& param() const
461 { return param_; }
462
463
464 template<class FluidState, class SingleWellState>
465 static Scalar computeTemperatureWeightFactor(const int perf_index, const int np, const FluidState& fs, const SingleWellState& ws)
466 {
467 const auto& perf_phase_rate = ws.perf_data.phase_rates;
468 // we only have one temperature pr cell any phaseIdx will do
469 Scalar cellTemperatures = fs.temperature(/*phaseIdx*/0).value();
470 Scalar weight_factor = 0.0;
471 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
472 if (!FluidSystem::phaseIsActive(phaseIdx)) {
473 continue;
474 }
475 Scalar cellInternalEnergy = fs.enthalpy(phaseIdx).value() -
476 fs.pressure(phaseIdx).value() / fs.density(phaseIdx).value();
477 Scalar cellBinv = fs.invB(phaseIdx).value();
478 Scalar cellDensity = fs.density(phaseIdx).value();
479 Scalar perfPhaseRate = perf_phase_rate[perf_index*np + phaseIdx];
480 weight_factor += cellDensity * (perfPhaseRate / cellBinv) * (cellInternalEnergy / cellTemperatures);
481 }
482 return (std::abs(weight_factor) + 1e-13);
483 }
484
485 protected:
487
488 // a vector of all the wells.
489 std::vector<WellInterfacePtr> well_container_{};
490
491 std::vector<bool> is_cell_perforated_{};
492
493 void initializeWellState(const int timeStepIdx);
494
495 // create the well container
496 void createWellContainer(const int report_step) override;
497
499 createWellPointer(const int wellID,
500 const int report_step) const;
501
502 template <typename WellType>
503 std::unique_ptr<WellType>
504 createTypedWellPointer(const int wellID,
505 const int time_step) const;
506
507 WellInterfacePtr createWellForWellTest(const std::string& well_name,
508 const int report_step,
509 DeferredLogger& deferred_logger) const;
510
512 std::size_t global_num_cells_{};
513 // the number of the cells in the local grid
514 std::size_t local_num_cells_{};
516 std::vector<Scalar> depth_{};
518 std::unique_ptr<RateConverterType> rateConverter_{};
519 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>> regionalAveragePressureCalculator_{};
520
524
525 // A flag to tell the convergence report whether we need to take another newton step
527
528 std::vector<Scalar> B_avg_{};
529
530 const EquilGrid& equilGrid() const
531 { return simulator_.vanguard().equilGrid(); }
532
533 const EclipseState& eclState() const
534 { return simulator_.vanguard().eclState(); }
535
536 // compute the well fluxes and assemble them in to the reservoir equations as source terms
537 // and in the well equations.
538 void assemble(const double dt);
539
540 // well controls and network pressures affect each other and are solved in an iterative manner.
541 // the function handles one iteration of updating well controls and network pressures.
542 // it is possible to decouple the update of well controls and network pressures further.
543 // the returned two booleans are {continue_due_to_network, well_group_control_changed}, respectively
544 std::tuple<bool, bool, Scalar> updateWellControlsAndNetworkIteration(const bool mandatory_network_balance,
545 const bool relax_network_tolerance,
546 const bool optimize_gas_lift,
547 const double dt,
548 DeferredLogger& local_deferredLogger);
549
558 void initializeLocalWellStructure(const int reportStepIdx,
559 const bool enableWellPIScaling);
560
564 void initializeGroupStructure(const int reportStepIdx);
565
566 // called at the end of a time step
567 void timeStepSucceeded(const double simulationTime, const double dt);
568
569 // called at the end of a report step
570 void endReportStep();
571
572 // setting the well_solutions_ based on well_state.
574
576
577 void computePotentials(const std::size_t widx,
578 const WellState<Scalar, IndexTraits>& well_state_copy,
579 std::string& exc_msg,
580 ExceptionType::ExcEnum& exc_type) override;
581
582 const std::vector<Scalar>& wellPerfEfficiencyFactors() const;
583
584 void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger& deferred_logger) override;
585 void calculateProductivityIndexValues(DeferredLogger& deferred_logger) override;
587 DeferredLogger& deferred_logger);
588
589 // The number of conservation quantities.
590 int numConservationQuantities() const;
591
592 int reportStepIndex() const;
593
594 void assembleWellEq(const double dt);
595
596 void prepareWellsBeforeAssembling(const double dt);
597
599
600 void extractLegacyDepth_();
601
603 void updateWellTestState(const double simulationTime, WellTestState& wellTestState);
604
605 void wellTesting(const int timeStepIdx, const double simulationTime, DeferredLogger& deferred_logger);
606
607 void calcResvCoeff(const int fipnum,
608 const int pvtreg,
609 const std::vector<Scalar>& production_rates,
610 std::vector<Scalar>& resv_coeff) const override;
611
612 void calcInjResvCoeff(const int fipnum,
613 const int pvtreg,
614 std::vector<Scalar>& resv_coeff) const override;
615
617
618 private:
621 BlackoilWellModelNldd<TypeTag>* nldd_ = nullptr;
622
623 // These members are used to avoid reallocation in specific functions
624 // instead of using local variables.
625 // Their state is not relevant between function calls, so they can
626 // (and must) be mutable, as the functions using them are const.
627 mutable BVector x_local_;
628
629 // Store cell rates after assembling to avoid iterating all wells and connections for every element
630 std::map<int, RateVector> cellRates_;
631
632 void assignWellTracerRates(data::Wells& wsrpt) const;
633 void assignWellSpeciesRates(data::Wells& wsrpt) const;
634
635 [[nodiscard]] auto rsConstInfo() const
636 -> typename WellState<Scalar,IndexTraits>::RsConstInfo;
637 };
638
639} // namespace Opm
640
641#include "BlackoilWellModel_impl.hpp"
642
643#endif // OPM_BLACKOILWELLMODEL_HEADER_INCLUDED
Contains the high level supplements required to extend the black oil model by bioeffects.
Definition: blackoilbioeffectsmodules.hh:95
Contains the high level supplements required to extend the black oil model by polymer.
Definition: blackoilpolymermodules.hh:65
Class for handling the gaslift in the blackoil well model.
Definition: BlackoilWellModelGasLift.hpp:96
Class for handling the blackoil well model.
Definition: BlackoilWellModelGeneric.hpp:97
const GroupState< GetPropType< TypeTag, Properties::Scalar > > & groupState() const
Definition: BlackoilWellModelGeneric.hpp:142
BlackoilWellModelWBP< GetPropType< TypeTag, Properties::Scalar >, GetPropType< TypeTag, Properties::FluidSystem >::IndexTraitsType > wbp_
Definition: BlackoilWellModelGeneric.hpp:518
const WellState< GetPropType< TypeTag, Properties::Scalar >, GetPropType< TypeTag, Properties::FluidSystem >::IndexTraitsType > & wellState() const
Definition: BlackoilWellModelGeneric.hpp:152
void assignMassGasRate(data::Wells &wsrpt, const GetPropType< TypeTag, Properties::Scalar > gasDensity) const
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:98
ReservoirCoupling::Proxy< Scalar > & rescoup()
Get the reservoir coupling proxy.
Definition: BlackoilWellModel.hpp:379
void initializeGroupStructure(const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:297
bool network_needs_more_balancing_force_another_newton_iteration_
Definition: BlackoilWellModel.hpp:526
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:2170
ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster()
Get reference to reservoir coupling master.
Definition: BlackoilWellModel.hpp:405
void endTimeStep()
Definition: BlackoilWellModel.hpp:166
void prepareTimeStep(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1989
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:1303
WellInterfacePtr createWellPointer(const int wellID, const int report_step) const
Definition: BlackoilWellModel_impl.hpp:1099
GuideRateHandler< Scalar, IndexTraits > guide_rate_handler_
Definition: BlackoilWellModel.hpp:522
const EclipseState & eclState() const
Definition: BlackoilWellModel.hpp:533
void prepareWellsBeforeAssembling(const double dt)
Definition: BlackoilWellModel_impl.hpp:1379
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:200
void init()
Definition: BlackoilWellModel_impl.hpp:166
void setNlddAdapter(BlackoilWellModelNldd< TypeTag > *mod)
Definition: BlackoilWellModel.hpp:373
const Simulator & simulator() const
Definition: BlackoilWellModel.hpp:370
std::vector< Scalar > depth_
Definition: BlackoilWellModel.hpp:516
static constexpr bool has_geochem_
Definition: BlackoilWellModel.hpp:126
static constexpr EnergyModules energyModuleType_
Definition: BlackoilWellModel.hpp:123
std::size_t global_num_cells_
Definition: BlackoilWellModel.hpp:512
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilWellModel.hpp:107
const std::vector< Scalar > & B_avg() const
Definition: BlackoilWellModel.hpp:457
SimulatorReportSingle last_report_
Definition: BlackoilWellModel.hpp:521
static const int solventSaturationIdx
Definition: BlackoilWellModel.hpp:120
Scalar gravity_
Definition: BlackoilWellModel.hpp:515
void initWellContainer(const int reportStepIdx) override
Definition: BlackoilWellModel_impl.hpp:185
void beginReportStep(const int time_step)
Definition: BlackoilWellModel_impl.hpp:202
const WellInterface< TypeTag > & getWell(const std::string &well_name) const
Definition: BlackoilWellModel_impl.hpp:2142
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: BlackoilWellModel.hpp:103
void beginIteration()
Definition: BlackoilWellModel.hpp:157
bool isReservoirCouplingMaster() const
Check if this process is a reservoir coupling master.
Definition: BlackoilWellModel.hpp:398
Dune::FieldVector< Scalar, numEq > VectorBlockType
Definition: BlackoilWellModel.hpp:130
void initFromRestartFile(const RestartValue &restartValues)
Definition: BlackoilWellModel.hpp:190
const ModelParameters param_
Definition: BlackoilWellModel.hpp:511
void beginEpisode()
Definition: BlackoilWellModel.hpp:149
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: BlackoilWellModel.hpp:104
GetPropType< TypeTag, Properties::Grid > Grid
Definition: BlackoilWellModel.hpp:101
GetPropType< TypeTag, Properties::EquilGrid > EquilGrid
Definition: BlackoilWellModel.hpp:102
int numConservationQuantities() const
Definition: BlackoilWellModel_impl.hpp:2114
const std::map< std::string, int > & well_domain() const
Definition: BlackoilWellModel.hpp:323
bool updateWellControls(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1666
void endIteration()
Definition: BlackoilWellModel.hpp:163
static constexpr bool has_polymer_
Definition: BlackoilWellModel.hpp:122
int reportStepIndex() const
Definition: BlackoilWellModel_impl.hpp:2158
void calculateProductivityIndexValues(DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1931
void extractLegacyDepth_()
Definition: BlackoilWellModel_impl.hpp:2130
void extractLegacyCellPvtRegionIndex_()
Definition: BlackoilWellModel_impl.hpp:2098
void recoverWellSolutionAndUpdateWellStateDomain(const BVector &x, const int domainIdx)
Definition: BlackoilWellModel_impl.hpp:1587
void updateAverageFormationFactor()
Definition: BlackoilWellModel_impl.hpp:2041
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilWellModel.hpp:106
bool isReservoirCouplingSlave() const
Check if this process is a reservoir coupling slave.
Definition: BlackoilWellModel.hpp:401
void initializeWellState(const int timeStepIdx)
Definition: BlackoilWellModel_impl.hpp:847
const Grid & grid() const
Definition: BlackoilWellModel.hpp:367
void updatePrimaryVariables()
Definition: BlackoilWellModel_impl.hpp:2089
const ReservoirCoupling::Proxy< Scalar > & rescoup() const
Definition: BlackoilWellModel.hpp:380
void computeWellTemperature()
Definition: BlackoilWellModel_impl.hpp:2192
static const int numEq
Definition: BlackoilWellModel.hpp:119
void addWellPressureEquations(PressureMatrix &jacobian, const BVector &weights, const bool use_well_weights) const
Definition: BlackoilWellModel_impl.hpp:1488
const SimulatorReportSingle & lastReport() const
Definition: BlackoilWellModel_impl.hpp:726
bool updateWellControlsAndNetwork(const bool mandatory_network_balance, const double dt, DeferredLogger &local_deferredLogger)
Definition: BlackoilWellModel_impl.hpp:1244
void endEpisode()
Definition: BlackoilWellModel.hpp:172
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: BlackoilWellModel.hpp:110
void addWellContributions(SparseMatrixAdapter &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1478
const EquilGrid & equilGrid() const
Definition: BlackoilWellModel.hpp:530
void assembleWellEq(const double dt)
Definition: BlackoilWellModel_impl.hpp:1367
WellInterfacePtr createWellForWellTest(const std::string &well_name, const int report_step, DeferredLogger &deferred_logger) const
Definition: BlackoilWellModel_impl.hpp:1150
GetPropType< TypeTag, Properties::Indices > Indices
Definition: BlackoilWellModel.hpp:105
int compressedIndexForInteriorLGR(const std::string &lgr_tag, const Connection &conn) const override
Definition: BlackoilWellModel.hpp:347
void calculateExplicitQuantities() const
Definition: BlackoilWellModel_impl.hpp:1651
std::vector< WellInterfacePtr > well_container_
Definition: BlackoilWellModel.hpp:489
static constexpr std::size_t pressureVarIndex
Definition: BlackoilWellModel.hpp:117
void updateAndCommunicate(const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:1746
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:810
bool addMatrixContributions() const
Definition: BlackoilWellModel.hpp:336
void beginTimeStep()
Definition: BlackoilWellModel_impl.hpp:329
std::vector< Scalar > B_avg_
Definition: BlackoilWellModel.hpp:528
const std::vector< Scalar > & wellPerfEfficiencyFactors() const
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: BlackoilWellModel.hpp:108
bool updateGroupControls(const Group &group, DeferredLogger &deferred_logger, const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:1780
std::map< std::string, std::unique_ptr< AverageRegionalPressureType > > regionalAveragePressureCalculator_
Definition: BlackoilWellModel.hpp:519
void calcInjResvCoeff(const int fipnum, const int pvtreg, std::vector< Scalar > &resv_coeff) const override
Definition: BlackoilWellModel_impl.hpp:2181
static Scalar computeTemperatureWeightFactor(const int perf_index, const int np, const FluidState &fs, const SingleWellState &ws)
Definition: BlackoilWellModel.hpp:465
GetPropType< TypeTag, Properties::GlobalEqVector > GlobalEqVector
Definition: BlackoilWellModel.hpp:109
void initializeLocalWellStructure(const int reportStepIdx, const bool enableWellPIScaling)
Definition: BlackoilWellModel_impl.hpp:251
static constexpr bool has_solvent_
Definition: BlackoilWellModel.hpp:121
Dune::BlockVector< VectorBlockType > BVector
Definition: BlackoilWellModel.hpp:131
static constexpr bool has_micp_
Definition: BlackoilWellModel.hpp:125
std::unique_ptr< RateConverterType > rateConverter_
Definition: BlackoilWellModel.hpp:518
BlackoilWellModel(Simulator &simulator)
Definition: BlackoilWellModel_impl.hpp:79
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:649
std::vector< bool > is_cell_perforated_
Definition: BlackoilWellModel.hpp:491
ConvergenceReport getWellConvergence(const std::vector< Scalar > &B_avg, const bool checkWellGroupControlsAndNetwork=false) const
Definition: BlackoilWellModel_impl.hpp:1600
int numStrictIterations() const
Definition: BlackoilWellModel.hpp:339
typename FluidSystem::IndexTraitsType IndexTraits
Definition: BlackoilWellModel.hpp:114
void updateCellRatesForDomain(int domainIndex, const std::map< std::string, int > &well_domain_map)
Definition: BlackoilWellModel_impl.hpp:1425
data::WellBlockAveragePressures wellBlockAveragePressures() const
Definition: BlackoilWellModel.hpp:245
void assembleWellEqWithoutIteration(const double dt)
Definition: BlackoilWellModel_impl.hpp:1393
void updateCellRates()
Definition: BlackoilWellModel_impl.hpp:1413
void assemble(const double dt)
Definition: BlackoilWellModel_impl.hpp:1175
auto begin() const
Definition: BlackoilWellModel.hpp:332
std::size_t local_num_cells_
Definition: BlackoilWellModel.hpp:514
auto end() const
Definition: BlackoilWellModel.hpp:333
bool alternative_well_rate_init_
Definition: BlackoilWellModel.hpp:517
void timeStepSucceeded(const double simulationTime, const double dt)
Definition: BlackoilWellModel_impl.hpp:736
std::unique_ptr< WellType > createTypedWellPointer(const int wellID, const int time_step) const
Definition: BlackoilWellModel_impl.hpp:1119
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:1896
Simulator & simulator_
Definition: BlackoilWellModel.hpp:486
void updateGuideRates(const int report_step_idx, const double sim_time)
Update guide rates for all wells and groups.
Definition: BlackoilWellModel.hpp:389
ReservoirCoupling::Proxy< Scalar > rescoup_
Definition: BlackoilWellModel.hpp:523
void createWellContainer(const int report_step) override
Definition: BlackoilWellModel_impl.hpp:890
std::unique_ptr< WellInterface< TypeTag > > WellInterfacePtr
Definition: BlackoilWellModel.hpp:187
ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave()
Get reference to reservoir coupling slave.
Definition: BlackoilWellModel.hpp:411
data::Wells wellData() const
Definition: BlackoilWellModel.hpp:207
const ModelParameters & param() const
Definition: BlackoilWellModel.hpp:460
void updateWellTestState(const double simulationTime, WellTestState &wellTestState)
upate the wellTestState related to economic limits
Definition: BlackoilWellModel_impl.hpp:1825
void addWellPressureEquationsDomain(PressureMatrix &jacobian, const BVector &weights, const bool use_well_weights, const int domainIndex) const
Definition: BlackoilWellModel.hpp:295
void addReservoirSourceTerms(GlobalEqVector &residual, const std::vector< typename SparseMatrixAdapter::MatrixBlock * > &diagMatAddress) const
Definition: BlackoilWellModel_impl.hpp:1510
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:124
void recoverWellSolutionAndUpdateWellState(const BVector &x)
Definition: BlackoilWellModel_impl.hpp:1563
void addWellPressureEquationsStruct(PressureMatrix &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1540
void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1945
void endReportStep()
Definition: BlackoilWellModel_impl.hpp:709
data::WellBlockAveragePressures computeWellBlockAveragePressures(const Scalar gravity) const
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
Definition: GroupStateHelper.hpp:57
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:43
PerfData< Scalar > perf_data
Definition: SingleWellState.hpp:143
Definition: BlackoilWellModel.hpp:87
Definition: WellConnectionAuxiliaryModule.hpp:39
Definition: WellContributions.hpp:51
Definition: WellInterface.hpp:77
Definition: WellState.hpp:66
data::Wells report(const int *globalCellIdxMap, const std::function< bool(const int)> &wasDynamicallyClosed, const RsConstInfo &rsConst=RsConstInfo{}) const
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 BlackoilModel.
Definition: BlackoilModelParameters.hpp:194
bool matrix_add_well_contributions_
Whether to add influences of wells between cells to the matrix and preconditioner matrix.
Definition: BlackoilModelParameters.hpp:308
bool use_multisegment_well_
Definition: BlackoilModelParameters.hpp:302
int strict_outer_iter_wells_
Newton iteration where wells are stricly convergent.
Definition: BlackoilModelParameters.hpp:258
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34