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
27
28#include <dune/common/fmatrix.hh>
29#include <dune/istl/bcrsmatrix.hh>
30#include <dune/istl/matrixmatrix.hh>
31
32#include <opm/common/OpmLog/OpmLog.hpp>
33
34#include <opm/input/eclipse/Schedule/Group/Group.hpp>
35#include <opm/input/eclipse/Schedule/Group/GuideRate.hpp>
36#include <opm/input/eclipse/Schedule/Schedule.hpp>
37#include <opm/input/eclipse/Schedule/Well/WellTestState.hpp>
38
39#include <opm/material/densead/Math.hpp>
40
43
45
48
50
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
127 // TODO: where we should put these types, WellInterface or Well Model?
128 // or there is some other strategy, like TypeTag
129 using VectorBlockType = Dune::FieldVector<Scalar, numEq>;
130 using BVector = Dune::BlockVector<VectorBlockType>;
131
134
135 // For the conversion between the surface volume rate and reservoir voidage rate
136 using RateConverterType = RateConverter::
137 SurfaceToReservoirVoidage<FluidSystem, std::vector<int> >;
138
139 // For computing average pressured used by gpmaint
140 using AverageRegionalPressureType = RegionAverageCalculator::
141 AverageRegionalPressure<FluidSystem, std::vector<int> >;
142
144
145 void init();
146 void initWellContainer(const int reportStepIdx) override;
147
149 {
150 OPM_TIMEBLOCK(beginEpsiode);
151 beginReportStep(simulator_.episodeIndex());
152 }
153
154 void beginTimeStep();
155
157 {
158 OPM_TIMEBLOCK(beginIteration);
159 assemble(simulator_.model().newtonMethod().numIterations(),
160 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 {
213 return this->wasDynamicallyShutThisTimeStep(well_index);
214 });
215
217 .assignWellGuideRates(wsrpt, this->reportStepIndex());
218
219 this->assignWellTracerRates(wsrpt);
220
221 if (const auto& rspec = eclState().runspec();
222 rspec.co2Storage() || rspec.h2Storage())
223 {
224 // The gas reference density (surface condition) is the
225 // same for all PVT regions in CO2STORE/H2STORE runs so,
226 // for simplicity, we use region zero (0) here.
227
228 this->assignMassGasRate(wsrpt, FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, 0));
229 }
230
231 this->assignWellTargets(wsrpt);
232
233 this->assignDynamicWellStatus(wsrpt);
234
235 // Assigning (a subset of the) property values in shut
236 // connections should be the last step of wellData().
237 this->assignShutConnections(wsrpt, this->reportStepIndex());
238
239 return wsrpt;
240 }
241
242 data::WellBlockAveragePressures wellBlockAveragePressures() const
243 {
245 }
246
247#if COMPILE_GPU_BRIDGE
248 // accumulate the contributions of all Wells in the WellContributions object
249 void getWellContributions(WellContributions<Scalar>& x) const;
250#endif
251
252 // Check if well equations is converged.
253 ConvergenceReport getWellConvergence(const std::vector<Scalar>& B_avg, const bool checkWellGroupControlsAndNetwork = false) const;
254
255 const SimulatorReportSingle& lastReport() const;
256
257 void addWellContributions(SparseMatrixAdapter& jacobian) const;
258
259 // add source from wells to the reservoir matrix
261 const std::vector<typename SparseMatrixAdapter::MatrixBlock*>& diagMatAddress) const;
262
263 // called at the beginning of a report step
264 void beginReportStep(const int time_step);
265
266 // it should be able to go to prepareTimeStep(), however, the updateWellControls()
267 // makes it a little more difficult. unless we introduce if (iterationIdx != 0) to avoid doing the above function
268 // twice at the beginning of the time step
271 void calculateExplicitQuantities() const;
272 // some preparation work, mostly related to group control and RESV,
273 // at the beginning of each time step (Not report step)
274 void prepareTimeStep(DeferredLogger& deferred_logger);
275
276 bool
277 updateWellControls(DeferredLogger& deferred_logger);
278
279 void updateAndCommunicate(const int reportStepIdx,
280 const int iterationIdx);
281
282 bool updateGroupControls(const Group& group,
283 DeferredLogger& deferred_logger,
284 const int reportStepIdx,
285 const int iterationIdx);
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 receiveGroupTargetsFromMaster(const int reportStepIdx);
436 void sendMasterGroupTargetsToSlaves();
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 protected:
465
466 // a vector of all the wells.
467 std::vector<WellInterfacePtr> well_container_{};
468
469 std::vector<bool> is_cell_perforated_{};
470
471 void initializeWellState(const int timeStepIdx);
472
473 // create the well container
474 void createWellContainer(const int report_step) override;
475
477 createWellPointer(const int wellID,
478 const int report_step) const;
479
480 template <typename WellType>
481 std::unique_ptr<WellType>
482 createTypedWellPointer(const int wellID,
483 const int time_step) const;
484
485 WellInterfacePtr createWellForWellTest(const std::string& well_name,
486 const int report_step,
487 DeferredLogger& deferred_logger) const;
488
490 std::size_t global_num_cells_{};
491 // the number of the cells in the local grid
492 std::size_t local_num_cells_{};
494 std::vector<Scalar> depth_{};
496 std::unique_ptr<RateConverterType> rateConverter_{};
497 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>> regionalAveragePressureCalculator_{};
498
502
503 // A flag to tell the convergence report whether we need to take another newton step
505
506 std::vector<Scalar> B_avg_{};
507
508 const EquilGrid& equilGrid() const
509 { return simulator_.vanguard().equilGrid(); }
510
511 const EclipseState& eclState() const
512 { return simulator_.vanguard().eclState(); }
513
514 // compute the well fluxes and assemble them in to the reservoir equations as source terms
515 // and in the well equations.
516 void assemble(const int iterationIdx,
517 const double dt);
518
519 // well controls and network pressures affect each other and are solved in an iterative manner.
520 // the function handles one iteration of updating well controls and network pressures.
521 // it is possible to decouple the update of well controls and network pressures further.
522 // the returned two booleans are {continue_due_to_network, well_group_control_changed}, respectively
523 std::tuple<bool, bool, Scalar> updateWellControlsAndNetworkIteration(const bool mandatory_network_balance,
524 const bool relax_network_tolerance,
525 const bool optimize_gas_lift,
526 const double dt,
527 DeferredLogger& local_deferredLogger);
528
537 void initializeLocalWellStructure(const int reportStepIdx,
538 const bool enableWellPIScaling);
539
543 void initializeGroupStructure(const int reportStepIdx);
544
545 // called at the end of a time step
546 void timeStepSucceeded(const double simulationTime, const double dt);
547
548 // called at the end of a report step
549 void endReportStep();
550
551 // setting the well_solutions_ based on well_state.
553
555
556 void computePotentials(const std::size_t widx,
557 const WellState<Scalar, IndexTraits>& well_state_copy,
558 std::string& exc_msg,
559 ExceptionType::ExcEnum& exc_type) override;
560
561 const std::vector<Scalar>& wellPerfEfficiencyFactors() const;
562
563 void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger& deferred_logger) override;
564 void calculateProductivityIndexValues(DeferredLogger& deferred_logger) override;
566 DeferredLogger& deferred_logger);
567
568 // The number of conservation quantities.
569 int numConservationQuantities() const;
570
571 int reportStepIndex() const;
572
573 void assembleWellEq(const double dt);
574
575 void prepareWellsBeforeAssembling(const double dt);
576
578
579 void extractLegacyDepth_();
580
582 void updateWellTestState(const double simulationTime, WellTestState& wellTestState);
583
584 void wellTesting(const int timeStepIdx, const double simulationTime, DeferredLogger& deferred_logger);
585
586 void calcResvCoeff(const int fipnum,
587 const int pvtreg,
588 const std::vector<Scalar>& production_rates,
589 std::vector<Scalar>& resv_coeff) const override;
590
591 void calcInjResvCoeff(const int fipnum,
592 const int pvtreg,
593 std::vector<Scalar>& resv_coeff) const override;
594
596
597 private:
600 BlackoilWellModelNldd<TypeTag>* nldd_ = nullptr;
601
602 // These members are used to avoid reallocation in specific functions
603 // instead of using local variables.
604 // Their state is not relevant between function calls, so they can
605 // (and must) be mutable, as the functions using them are const.
606 mutable BVector x_local_;
607
608 // Store cell rates after assembling to avoid iterating all wells and connections for every element
609 std::map<int, RateVector> cellRates_;
610
611 void assignWellTracerRates(data::Wells& wsrpt) const;
612 };
613
614} // namespace Opm
615
617
618#endif // OPM_BLACKOILWELLMODEL_HEADER_INCLUDED
Contains the high level supplements required to extend the black oil model by bioeffects.
Definition: blackoilbioeffectsmodules.hh:93
Contains the high level supplements required to extend the black oil model by polymer.
Definition: blackoilpolymermodules.hh:64
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
const GroupState< GetPropType< TypeTag, Properties::Scalar > > & groupState() const
Definition: BlackoilWellModelGeneric.hpp:141
BlackoilWellModelWBP< GetPropType< TypeTag, Properties::Scalar >, GetPropType< TypeTag, Properties::FluidSystem >::IndexTraitsType > wbp_
Definition: BlackoilWellModelGeneric.hpp:517
const WellState< GetPropType< TypeTag, Properties::Scalar >, GetPropType< TypeTag, Properties::FluidSystem >::IndexTraitsType > & wellState() const
Definition: BlackoilWellModelGeneric.hpp:151
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:294
bool network_needs_more_balancing_force_another_newton_iteration_
Definition: BlackoilWellModel.hpp:504
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:2142
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:1960
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:1274
bool updateGroupControls(const Group &group, DeferredLogger &deferred_logger, const int reportStepIdx, const int iterationIdx)
Definition: BlackoilWellModel_impl.hpp:1751
WellInterfacePtr createWellPointer(const int wellID, const int report_step) const
Definition: BlackoilWellModel_impl.hpp:1072
GuideRateHandler< Scalar, IndexTraits > guide_rate_handler_
Definition: BlackoilWellModel.hpp:500
const EclipseState & eclState() const
Definition: BlackoilWellModel.hpp:511
void prepareWellsBeforeAssembling(const double dt)
Definition: BlackoilWellModel_impl.hpp:1350
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:163
void setNlddAdapter(BlackoilWellModelNldd< TypeTag > *mod)
Definition: BlackoilWellModel.hpp:373
const Simulator & simulator() const
Definition: BlackoilWellModel.hpp:370
std::vector< Scalar > depth_
Definition: BlackoilWellModel.hpp:494
static constexpr EnergyModules energyModuleType_
Definition: BlackoilWellModel.hpp:123
std::size_t global_num_cells_
Definition: BlackoilWellModel.hpp:490
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:499
static const int solventSaturationIdx
Definition: BlackoilWellModel.hpp:120
Scalar gravity_
Definition: BlackoilWellModel.hpp:493
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:2113
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: BlackoilWellModel.hpp:103
void beginIteration()
Definition: BlackoilWellModel.hpp:156
bool isReservoirCouplingMaster() const
Check if this process is a reservoir coupling master.
Definition: BlackoilWellModel.hpp:398
Dune::FieldVector< Scalar, numEq > VectorBlockType
Definition: BlackoilWellModel.hpp:129
void initFromRestartFile(const RestartValue &restartValues)
Definition: BlackoilWellModel.hpp:190
const ModelParameters param_
Definition: BlackoilWellModel.hpp:489
void beginEpisode()
Definition: BlackoilWellModel.hpp:148
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:2085
const std::map< std::string, int > & well_domain() const
Definition: BlackoilWellModel.hpp:323
bool updateWellControls(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1636
void endIteration()
Definition: BlackoilWellModel.hpp:163
static constexpr bool has_polymer_
Definition: BlackoilWellModel.hpp:122
int reportStepIndex() const
Definition: BlackoilWellModel_impl.hpp:2130
void calculateProductivityIndexValues(DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1902
void extractLegacyDepth_()
Definition: BlackoilWellModel_impl.hpp:2101
void extractLegacyCellPvtRegionIndex_()
Definition: BlackoilWellModel_impl.hpp:2069
void recoverWellSolutionAndUpdateWellStateDomain(const BVector &x, const int domainIdx)
Definition: BlackoilWellModel_impl.hpp:1558
void updateAverageFormationFactor()
Definition: BlackoilWellModel_impl.hpp:2012
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:820
void assemble(const int iterationIdx, const double dt)
Definition: BlackoilWellModel_impl.hpp:1148
const Grid & grid() const
Definition: BlackoilWellModel.hpp:367
void updatePrimaryVariables()
Definition: BlackoilWellModel_impl.hpp:2060
const ReservoirCoupling::Proxy< Scalar > & rescoup() const
Definition: BlackoilWellModel.hpp:380
void computeWellTemperature()
Definition: BlackoilWellModel_impl.hpp:2164
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:1459
const SimulatorReportSingle & lastReport() const
Definition: BlackoilWellModel_impl.hpp:701
bool updateWellControlsAndNetwork(const bool mandatory_network_balance, const double dt, DeferredLogger &local_deferredLogger)
Definition: BlackoilWellModel_impl.hpp:1215
void endEpisode()
Definition: BlackoilWellModel.hpp:172
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: BlackoilWellModel.hpp:110
void addWellContributions(SparseMatrixAdapter &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1449
const EquilGrid & equilGrid() const
Definition: BlackoilWellModel.hpp:508
void assembleWellEq(const double dt)
Definition: BlackoilWellModel_impl.hpp:1338
WellInterfacePtr createWellForWellTest(const std::string &well_name, const int report_step, DeferredLogger &deferred_logger) const
Definition: BlackoilWellModel_impl.hpp:1123
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:1621
std::vector< WellInterfacePtr > well_container_
Definition: BlackoilWellModel.hpp:467
static constexpr std::size_t pressureVarIndex
Definition: BlackoilWellModel.hpp:117
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: BlackoilWellModel.hpp:289
bool empty() const
Definition: BlackoilWellModel.hpp:334
void computeTotalRatesForDof(RateVector &rate, unsigned globalIdx) const
Definition: BlackoilWellModel_impl.hpp:783
bool addMatrixContributions() const
Definition: BlackoilWellModel.hpp:336
void beginTimeStep()
Definition: BlackoilWellModel_impl.hpp:326
void updateAndCommunicate(const int reportStepIdx, const int iterationIdx)
Definition: BlackoilWellModel_impl.hpp:1717
std::vector< Scalar > B_avg_
Definition: BlackoilWellModel.hpp:506
const std::vector< Scalar > & wellPerfEfficiencyFactors() const
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: BlackoilWellModel.hpp:108
std::map< std::string, std::unique_ptr< AverageRegionalPressureType > > regionalAveragePressureCalculator_
Definition: BlackoilWellModel.hpp:497
void calcInjResvCoeff(const int fipnum, const int pvtreg, std::vector< Scalar > &resv_coeff) const override
Definition: BlackoilWellModel_impl.hpp:2153
GetPropType< TypeTag, Properties::GlobalEqVector > GlobalEqVector
Definition: BlackoilWellModel.hpp:109
void initializeLocalWellStructure(const int reportStepIdx, const bool enableWellPIScaling)
Definition: BlackoilWellModel_impl.hpp:248
static constexpr bool has_solvent_
Definition: BlackoilWellModel.hpp:121
Dune::BlockVector< VectorBlockType > BVector
Definition: BlackoilWellModel.hpp:130
static constexpr bool has_micp_
Definition: BlackoilWellModel.hpp:125
std::unique_ptr< RateConverterType > rateConverter_
Definition: BlackoilWellModel.hpp:496
BlackoilWellModel(Simulator &simulator)
Definition: BlackoilWellModel_impl.hpp:76
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:626
std::vector< bool > is_cell_perforated_
Definition: BlackoilWellModel.hpp:469
ConvergenceReport getWellConvergence(const std::vector< Scalar > &B_avg, const bool checkWellGroupControlsAndNetwork=false) const
Definition: BlackoilWellModel_impl.hpp:1571
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:1396
data::WellBlockAveragePressures wellBlockAveragePressures() const
Definition: BlackoilWellModel.hpp:242
void assembleWellEqWithoutIteration(const double dt)
Definition: BlackoilWellModel_impl.hpp:1364
void updateCellRates()
Definition: BlackoilWellModel_impl.hpp:1384
auto begin() const
Definition: BlackoilWellModel.hpp:332
std::size_t local_num_cells_
Definition: BlackoilWellModel.hpp:492
auto end() const
Definition: BlackoilWellModel.hpp:333
bool alternative_well_rate_init_
Definition: BlackoilWellModel.hpp:495
void timeStepSucceeded(const double simulationTime, const double dt)
Definition: BlackoilWellModel_impl.hpp:711
std::unique_ptr< WellType > createTypedWellPointer(const int wellID, const int time_step) const
Definition: BlackoilWellModel_impl.hpp:1092
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:1867
Simulator & simulator_
Definition: BlackoilWellModel.hpp:464
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:501
void createWellContainer(const int report_step) override
Definition: BlackoilWellModel_impl.hpp:863
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:1796
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:1481
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:1534
void addWellPressureEquationsStruct(PressureMatrix &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1511
void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1916
void endReportStep()
Definition: BlackoilWellModel_impl.hpp:684
data::WellBlockAveragePressures computeWellBlockAveragePressures(const Scalar gravity) const
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
Definition: GroupStateHelper.hpp:54
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: 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
ExcEnum
Definition: DeferredLogger.hpp:45
Definition: blackoilbioeffectsmodules.hh:43
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