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#if HAVE_MPI
27#define RESERVOIR_COUPLING_ENABLED
28#endif
29#ifdef RESERVOIR_COUPLING_ENABLED
32#endif
33
34#include <dune/common/fmatrix.hh>
35#include <dune/istl/bcrsmatrix.hh>
36#include <dune/istl/matrixmatrix.hh>
37
38#include <opm/common/OpmLog/OpmLog.hpp>
39
40#include <opm/input/eclipse/Schedule/Group/Group.hpp>
41#include <opm/input/eclipse/Schedule/Group/GuideRate.hpp>
42#include <opm/input/eclipse/Schedule/Schedule.hpp>
43#include <opm/input/eclipse/Schedule/Well/WellTestState.hpp>
44
45#include <opm/material/densead/Math.hpp>
46
49
51
54
56
81
82#include <cstddef>
83#include <map>
84#include <memory>
85#include <optional>
86#include <string>
87#include <tuple>
88#include <vector>
89
90namespace Opm {
91
92template<class Scalar> class BlackoilWellModelNldd;
93template<class T, template <typename, typename...> class Storage> class SparseTable;
94
95#if COMPILE_GPU_BRIDGE
96template<class Scalar> class WellContributions;
97#endif
98
100 template<typename TypeTag>
101 class BlackoilWellModel : public WellConnectionAuxiliaryModule<TypeTag, BlackoilWellModel<TypeTag>>
102 , public BlackoilWellModelGeneric<GetPropType<TypeTag, Properties::Scalar>,
103 typename GetPropType<TypeTag, Properties::FluidSystem>::IndexTraitsType>
104 {
105 public:
106 // --------- Types ---------
118
120 using IndexTraits = typename FluidSystem::IndexTraitsType;
122
124
125 static const int numEq = Indices::numEq;
126 static const int solventSaturationIdx = Indices::solventSaturationIdx;
127 static constexpr bool has_solvent_ = getPropValue<TypeTag, Properties::EnableSolvent>();
128 static constexpr bool has_polymer_ = getPropValue<TypeTag, Properties::EnablePolymer>();
129 static constexpr EnergyModules energyModuleType_ = getPropValue<TypeTag, Properties::EnergyModuleType>();
130 static constexpr bool has_energy_ = (energyModuleType_ == EnergyModules::FullyImplicitThermal);
131 static constexpr bool has_micp_ = Indices::enableMICP;
132
133 // TODO: where we should put these types, WellInterface or Well Model?
134 // or there is some other strategy, like TypeTag
135 using VectorBlockType = Dune::FieldVector<Scalar, numEq>;
136 using BVector = Dune::BlockVector<VectorBlockType>;
137
140
141 // For the conversion between the surface volume rate and reservoir voidage rate
142 using RateConverterType = RateConverter::
143 SurfaceToReservoirVoidage<FluidSystem, std::vector<int> >;
144
145 // For computing average pressured used by gpmaint
146 using AverageRegionalPressureType = RegionAverageCalculator::
147 AverageRegionalPressure<FluidSystem, std::vector<int> >;
148
150
151 void init();
152 void initWellContainer(const int reportStepIdx) override;
153
155 {
156 OPM_TIMEBLOCK(beginEpsiode);
157 beginReportStep(simulator_.episodeIndex());
158 }
159
160 void beginTimeStep();
161
163 {
164 OPM_TIMEBLOCK(beginIteration);
165 assemble(simulator_.model().newtonMethod().numIterations(),
166 simulator_.timeStepSize());
167 }
168
170 { }
171
173 {
174 OPM_TIMEBLOCK(endTimeStep);
175 timeStepSucceeded(simulator_.time(), simulator_.timeStepSize());
176 }
177
179 {
181 }
182
184 unsigned globalIdx) const;
185
186 template <class Context>
188 const Context& context,
189 unsigned spaceIdx,
190 unsigned timeIdx) const;
191
192
193 using WellInterfacePtr = std::unique_ptr<WellInterface<TypeTag>>;
194
196 void initFromRestartFile(const RestartValue& restartValues)
197 {
198 initFromRestartFile(restartValues,
199 this->simulator_.vanguard().transferWTestState(),
200 grid().size(0),
202 this->simulator_.vanguard().enableDistributedWells());
203 }
204
206 void prepareDeserialize(const int report_step)
207 {
208 prepareDeserialize(report_step, grid().size(0),
210 this->simulator_.vanguard().enableDistributedWells());
211 }
212
213 data::Wells wellData() const
214 {
215 auto wsrpt = this->wellState()
216 .report(this->simulator_.vanguard().globalCell().data(),
217 [this](const int well_index)
218 {
219 return this->wasDynamicallyShutThisTimeStep(well_index);
220 });
221
223 .assignWellGuideRates(wsrpt, this->reportStepIndex());
224
225 this->assignWellTracerRates(wsrpt);
226
227 if (const auto& rspec = eclState().runspec();
228 rspec.co2Storage() || rspec.h2Storage())
229 {
230 // The gas reference density (surface condition) is the
231 // same for all PVT regions in CO2STORE/H2STORE runs so,
232 // for simplicity, we use region zero (0) here.
233
234 this->assignMassGasRate(wsrpt, FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, 0));
235 }
236
237 this->assignWellTargets(wsrpt);
238
239 this->assignDynamicWellStatus(wsrpt);
240
241 // Assigning (a subset of the) property values in shut
242 // connections should be the last step of wellData().
243 this->assignShutConnections(wsrpt, this->reportStepIndex());
244
245 return wsrpt;
246 }
247
248 data::WellBlockAveragePressures wellBlockAveragePressures() const
249 {
251 }
252
253#if COMPILE_GPU_BRIDGE
254 // accumulate the contributions of all Wells in the WellContributions object
255 void getWellContributions(WellContributions<Scalar>& x) const;
256#endif
257
258 // Check if well equations is converged.
259 ConvergenceReport getWellConvergence(const std::vector<Scalar>& B_avg, const bool checkWellGroupControlsAndNetwork = false) const;
260
261 const SimulatorReportSingle& lastReport() const;
262
263 void addWellContributions(SparseMatrixAdapter& jacobian) const;
264
265 // add source from wells to the reservoir matrix
267 const std::vector<typename SparseMatrixAdapter::MatrixBlock*>& diagMatAddress) const;
268
269 // called at the beginning of a report step
270 void beginReportStep(const int time_step);
271
272 // it should be able to go to prepareTimeStep(), however, the updateWellControls()
273 // makes it a little more difficult. unless we introduce if (iterationIdx != 0) to avoid doing the above function
274 // twice at the beginning of the time step
277 void calculateExplicitQuantities(DeferredLogger& deferred_logger) const;
278 // some preparation work, mostly related to group control and RESV,
279 // at the beginning of each time step (Not report step)
280 void prepareTimeStep(DeferredLogger& deferred_logger);
281
282 bool
283 updateWellControls(DeferredLogger& deferred_logger);
284
285 void updateAndCommunicate(const int reportStepIdx,
286 const int iterationIdx,
287 DeferredLogger& deferred_logger);
288
289 bool updateGroupControls(const Group& group,
290 DeferredLogger& deferred_logger,
291 const int reportStepIdx,
292 const int iterationIdx);
293
294 const WellInterface<TypeTag>& getWell(const std::string& well_name) const;
295
296 using PressureMatrix = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, 1, 1>>;
297
299 const BVector& weights,
300 const bool use_well_weights) const;
303 const BVector& weights,
304 const bool use_well_weights,
305 const int domainIndex) const
306 {
307 if (!nldd_) {
308 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
309 }
310 return nldd_->addWellPressureEquations(jacobian,
311 weights,
312 use_well_weights,
313 domainIndex);
314 }
315
317 const std::vector<WellInterfacePtr>& localNonshutWells() const
318 {
319 return well_container_;
320 }
321
323 {
324 if (!nldd_) {
325 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
326 }
327 return nldd_->well_local_cells();
328 }
329
330 const std::map<std::string, int>& well_domain() const
331 {
332 if (!nldd_) {
333 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
334 }
335
336 return nldd_->well_domain();
337 }
338
339 auto begin() const { return well_container_.begin(); }
340 auto end() const { return well_container_.end(); }
341 bool empty() const { return well_container_.empty(); }
342
345
348
349 int compressedIndexForInterior(int cartesian_cell_idx) const override
350 {
351 return simulator_.vanguard().compressedIndexForInterior(cartesian_cell_idx);
352 }
353
354 int compressedIndexForInteriorLGR(const std::string& lgr_tag, const Connection& conn) const override
355 {
356 return simulator_.vanguard().compressedIndexForInteriorLGR(lgr_tag, conn);
357 }
358
359 // using the solution x to recover the solution xw for wells and applying
360 // xw to update Well State
362
363 // using the solution x to recover the solution xw for wells and applying
364 // xw to update Well State
366 const int domainIdx);
367 // Update cellRates_ with contributions from all wells
368 void updateCellRates();
369
370 // Update cellRates_ with contributions from wells in a specific domain
371 void updateCellRatesForDomain(int domainIndex,
372 const std::map<std::string, int>& well_domain_map);
373
374 const Grid& grid() const
375 { return simulator_.vanguard().grid(); }
376
377 const Simulator& simulator() const
378 { return simulator_; }
379
381 { nldd_ = mod; }
382
383#ifdef RESERVOIR_COUPLING_ENABLED
385 return *(this->simulator_.reservoirCouplingMaster());
386 }
388 return *(this->simulator_.reservoirCouplingSlave());
389 }
391 return this->simulator_.reservoirCouplingMaster() != nullptr;
392 }
394 return this->simulator_.reservoirCouplingSlave() != nullptr;
395 }
397 {
398 this->guide_rate_handler_.setReservoirCouplingMaster(master);
400 }
402 {
403 this->guide_rate_handler_.setReservoirCouplingSlave(slave);
405 }
406
409
412
421 std::optional<ReservoirCoupling::ScopedLoggerGuard>
423
424 #endif
425
426 bool updateWellControlsAndNetwork(const bool mandatory_network_balance,
427 const double dt,
428 DeferredLogger& local_deferredLogger);
429
430 // TODO: finding a better naming
431 void assembleWellEqWithoutIteration(const double dt, DeferredLogger& deferred_logger);
432
433 const std::vector<Scalar>& B_avg() const
434 { return B_avg_; }
435
436 const ModelParameters& param() const
437 { return param_; }
438
439 protected:
441
442 // a vector of all the wells.
443 std::vector<WellInterfacePtr> well_container_{};
444
445 std::vector<bool> is_cell_perforated_{};
446
447 void initializeWellState(const int timeStepIdx);
448
449 // create the well container
450 void createWellContainer(const int report_step) override;
451
453 createWellPointer(const int wellID,
454 const int report_step) const;
455
456 template <typename WellType>
457 std::unique_ptr<WellType>
458 createTypedWellPointer(const int wellID,
459 const int time_step) const;
460
461 WellInterfacePtr createWellForWellTest(const std::string& well_name,
462 const int report_step,
463 DeferredLogger& deferred_logger) const;
464
466 std::size_t global_num_cells_{};
467 // the number of the cells in the local grid
468 std::size_t local_num_cells_{};
470 std::vector<Scalar> depth_{};
472 std::unique_ptr<RateConverterType> rateConverter_{};
473 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>> regionalAveragePressureCalculator_{};
474
477
478 // A flag to tell the convergence report whether we need to take another newton step
480
481 std::vector<Scalar> B_avg_{};
482
483 const EquilGrid& equilGrid() const
484 { return simulator_.vanguard().equilGrid(); }
485
486 const EclipseState& eclState() const
487 { return simulator_.vanguard().eclState(); }
488
489 // compute the well fluxes and assemble them in to the reservoir equations as source terms
490 // and in the well equations.
491 void assemble(const int iterationIdx,
492 const double dt);
493
494 // well controls and network pressures affect each other and are solved in an iterative manner.
495 // the function handles one iteration of updating well controls and network pressures.
496 // it is possible to decouple the update of well controls and network pressures further.
497 // the returned two booleans are {continue_due_to_network, well_group_control_changed}, respectively
498 std::tuple<bool, bool, Scalar> updateWellControlsAndNetworkIteration(const bool mandatory_network_balance,
499 const bool relax_network_tolerance,
500 const bool optimize_gas_lift,
501 const double dt,
502 DeferredLogger& local_deferredLogger);
503
512 void initializeLocalWellStructure(const int reportStepIdx,
513 const bool enableWellPIScaling);
514
518 void initializeGroupStructure(const int reportStepIdx);
519
520 // called at the end of a time step
521 void timeStepSucceeded(const double simulationTime, const double dt);
522
523 // called at the end of a report step
524 void endReportStep();
525
526 // setting the well_solutions_ based on well_state.
527 void updatePrimaryVariables(DeferredLogger& deferred_logger);
528
530
531 void computePotentials(const std::size_t widx,
532 const WellState<Scalar, IndexTraits>& well_state_copy,
533 std::string& exc_msg,
534 ExceptionType::ExcEnum& exc_type,
535 DeferredLogger& deferred_logger) override;
536
537 const std::vector<Scalar>& wellPerfEfficiencyFactors() const;
538
539 void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger& deferred_logger) override;
540 void calculateProductivityIndexValues(DeferredLogger& deferred_logger) override;
542 DeferredLogger& deferred_logger);
543
544 // The number of conservation quantities.
545 int numConservationQuantities() const;
546
547 int reportStepIndex() const;
548
549 void assembleWellEq(const double dt, DeferredLogger& deferred_logger);
550
551 void prepareWellsBeforeAssembling(const double dt, DeferredLogger& deferred_logger);
552
554
555 void extractLegacyDepth_();
556
558 void updateWellTestState(const double simulationTime, WellTestState& wellTestState);
559
560 void wellTesting(const int timeStepIdx, const double simulationTime, DeferredLogger& deferred_logger);
561
562 void calcResvCoeff(const int fipnum,
563 const int pvtreg,
564 const std::vector<Scalar>& production_rates,
565 std::vector<Scalar>& resv_coeff) override;
566
567 void calcInjResvCoeff(const int fipnum,
568 const int pvtreg,
569 std::vector<Scalar>& resv_coeff) override;
570
572
573 private:
576 BlackoilWellModelNldd<TypeTag>* nldd_ = nullptr;
577
578 // These members are used to avoid reallocation in specific functions
579 // instead of using local variables.
580 // Their state is not relevant between function calls, so they can
581 // (and must) be mutable, as the functions using them are const.
582 mutable BVector x_local_;
583
584 // Store cell rates after assembling to avoid iterating all wells and connections for every element
585 std::map<int, RateVector> cellRates_;
586
587 void assignWellTracerRates(data::Wells& wsrpt) const;
588 };
589
590} // namespace Opm
591
593
594#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
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:104
void initializeGroupStructure(const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:295
bool network_needs_more_balancing_force_another_newton_iteration_
Definition: BlackoilWellModel.hpp:479
ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster()
Definition: BlackoilWellModel.hpp:384
void endTimeStep()
Definition: BlackoilWellModel.hpp:172
void prepareTimeStep(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1962
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:1259
bool updateGroupControls(const Group &group, DeferredLogger &deferred_logger, const int reportStepIdx, const int iterationIdx)
Definition: BlackoilWellModel_impl.hpp:1742
WellInterfacePtr createWellPointer(const int wellID, const int report_step) const
Definition: BlackoilWellModel_impl.hpp:1057
GuideRateHandler< Scalar, IndexTraits > guide_rate_handler_
Definition: BlackoilWellModel.hpp:476
const EclipseState & eclState() const
Definition: BlackoilWellModel.hpp:486
void calcInjResvCoeff(const int fipnum, const int pvtreg, std::vector< Scalar > &resv_coeff) override
Definition: BlackoilWellModel_impl.hpp:2155
const std::vector< WellInterfacePtr > & localNonshutWells() const
Get list of local nonshut wells.
Definition: BlackoilWellModel.hpp:317
void prepareDeserialize(const int report_step)
Definition: BlackoilWellModel.hpp:206
void init()
Definition: BlackoilWellModel_impl.hpp:161
void setNlddAdapter(BlackoilWellModelNldd< TypeTag > *mod)
Definition: BlackoilWellModel.hpp:380
void calcResvCoeff(const int fipnum, const int pvtreg, const std::vector< Scalar > &production_rates, std::vector< Scalar > &resv_coeff) override
Definition: BlackoilWellModel_impl.hpp:2144
const Simulator & simulator() const
Definition: BlackoilWellModel.hpp:377
std::vector< Scalar > depth_
Definition: BlackoilWellModel.hpp:470
static constexpr EnergyModules energyModuleType_
Definition: BlackoilWellModel.hpp:129
void updateAndCommunicate(const int reportStepIdx, const int iterationIdx, DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1705
std::size_t global_num_cells_
Definition: BlackoilWellModel.hpp:466
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilWellModel.hpp:113
const std::vector< Scalar > & B_avg() const
Definition: BlackoilWellModel.hpp:433
SimulatorReportSingle last_report_
Definition: BlackoilWellModel.hpp:475
static const int solventSaturationIdx
Definition: BlackoilWellModel.hpp:126
Scalar gravity_
Definition: BlackoilWellModel.hpp:469
void initWellContainer(const int reportStepIdx) override
Definition: BlackoilWellModel_impl.hpp:180
void beginReportStep(const int time_step)
Definition: BlackoilWellModel_impl.hpp:197
const WellInterface< TypeTag > & getWell(const std::string &well_name) const
Definition: BlackoilWellModel_impl.hpp:2115
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: BlackoilWellModel.hpp:109
void beginIteration()
Definition: BlackoilWellModel.hpp:162
bool isReservoirCouplingMaster() const
Definition: BlackoilWellModel.hpp:390
Dune::FieldVector< Scalar, numEq > VectorBlockType
Definition: BlackoilWellModel.hpp:135
void initFromRestartFile(const RestartValue &restartValues)
Definition: BlackoilWellModel.hpp:196
const ModelParameters param_
Definition: BlackoilWellModel.hpp:465
void beginEpisode()
Definition: BlackoilWellModel.hpp:154
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: BlackoilWellModel.hpp:110
GetPropType< TypeTag, Properties::Grid > Grid
Definition: BlackoilWellModel.hpp:107
GetPropType< TypeTag, Properties::EquilGrid > EquilGrid
Definition: BlackoilWellModel.hpp:108
int numConservationQuantities() const
Definition: BlackoilWellModel_impl.hpp:2087
void assembleWellEq(const double dt, DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1323
const std::map< std::string, int > & well_domain() const
Definition: BlackoilWellModel.hpp:330
bool updateWellControls(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1624
void endIteration()
Definition: BlackoilWellModel.hpp:169
static constexpr bool has_polymer_
Definition: BlackoilWellModel.hpp:128
int reportStepIndex() const
Definition: BlackoilWellModel_impl.hpp:2132
void calculateProductivityIndexValues(DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1904
void extractLegacyDepth_()
Definition: BlackoilWellModel_impl.hpp:2103
void extractLegacyCellPvtRegionIndex_()
Definition: BlackoilWellModel_impl.hpp:2071
void recoverWellSolutionAndUpdateWellStateDomain(const BVector &x, const int domainIdx)
Definition: BlackoilWellModel_impl.hpp:1544
void updateAverageFormationFactor()
Definition: BlackoilWellModel_impl.hpp:2014
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilWellModel.hpp:112
bool isReservoirCouplingSlave() const
Definition: BlackoilWellModel.hpp:393
void initializeWellState(const int timeStepIdx)
Definition: BlackoilWellModel_impl.hpp:798
void assemble(const int iterationIdx, const double dt)
Definition: BlackoilWellModel_impl.hpp:1133
const Grid & grid() const
Definition: BlackoilWellModel.hpp:374
void computeWellTemperature()
Definition: BlackoilWellModel_impl.hpp:2166
static const int numEq
Definition: BlackoilWellModel.hpp:125
void addWellPressureEquations(PressureMatrix &jacobian, const BVector &weights, const bool use_well_weights) const
Definition: BlackoilWellModel_impl.hpp:1444
const SimulatorReportSingle & lastReport() const
Definition: BlackoilWellModel_impl.hpp:674
bool updateWellControlsAndNetwork(const bool mandatory_network_balance, const double dt, DeferredLogger &local_deferredLogger)
Definition: BlackoilWellModel_impl.hpp:1200
void endEpisode()
Definition: BlackoilWellModel.hpp:178
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: BlackoilWellModel.hpp:116
void addWellContributions(SparseMatrixAdapter &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1434
const EquilGrid & equilGrid() const
Definition: BlackoilWellModel.hpp:483
WellInterfacePtr createWellForWellTest(const std::string &well_name, const int report_step, DeferredLogger &deferred_logger) const
Definition: BlackoilWellModel_impl.hpp:1108
GetPropType< TypeTag, Properties::Indices > Indices
Definition: BlackoilWellModel.hpp:111
std::optional< ReservoirCoupling::ScopedLoggerGuard > setupRescoupScopedLogger(DeferredLogger &local_logger)
Setup RAII guard for reservoir coupling logger.
int compressedIndexForInteriorLGR(const std::string &lgr_tag, const Connection &conn) const override
Definition: BlackoilWellModel.hpp:354
void assembleWellEqWithoutIteration(const double dt, DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1349
void sendSlaveGroupDataToMaster()
Send comprehensive slave group data to master.
std::vector< WellInterfacePtr > well_container_
Definition: BlackoilWellModel.hpp:443
void receiveSlaveGroupData()
Receive comprehensive slave group data from slaves.
static constexpr std::size_t pressureVarIndex
Definition: BlackoilWellModel.hpp:123
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: BlackoilWellModel.hpp:296
bool empty() const
Definition: BlackoilWellModel.hpp:341
void computeTotalRatesForDof(RateVector &rate, unsigned globalIdx) const
Definition: BlackoilWellModel_impl.hpp:761
bool addMatrixContributions() const
Definition: BlackoilWellModel.hpp:343
void beginTimeStep()
Definition: BlackoilWellModel_impl.hpp:331
std::vector< Scalar > B_avg_
Definition: BlackoilWellModel.hpp:481
const std::vector< Scalar > & wellPerfEfficiencyFactors() const
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: BlackoilWellModel.hpp:114
void calculateExplicitQuantities(DeferredLogger &deferred_logger) const
Definition: BlackoilWellModel_impl.hpp:1609
std::map< std::string, std::unique_ptr< AverageRegionalPressureType > > regionalAveragePressureCalculator_
Definition: BlackoilWellModel.hpp:473
void updatePrimaryVariables(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:2062
GetPropType< TypeTag, Properties::GlobalEqVector > GlobalEqVector
Definition: BlackoilWellModel.hpp:115
void initializeLocalWellStructure(const int reportStepIdx, const bool enableWellPIScaling)
Definition: BlackoilWellModel_impl.hpp:250
static constexpr bool has_solvent_
Definition: BlackoilWellModel.hpp:127
void setReservoirCouplingMaster(ReservoirCouplingMaster< Scalar > *master)
Definition: BlackoilWellModel.hpp:396
Dune::BlockVector< VectorBlockType > BVector
Definition: BlackoilWellModel.hpp:136
static constexpr bool has_micp_
Definition: BlackoilWellModel.hpp:131
std::unique_ptr< RateConverterType > rateConverter_
Definition: BlackoilWellModel.hpp:472
BlackoilWellModel(Simulator &simulator)
Definition: BlackoilWellModel_impl.hpp:74
const SparseTable< int > & well_local_cells() const
Definition: BlackoilWellModel.hpp:322
void wellTesting(const int timeStepIdx, const double simulationTime, DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:594
std::vector< bool > is_cell_perforated_
Definition: BlackoilWellModel.hpp:445
ConvergenceReport getWellConvergence(const std::vector< Scalar > &B_avg, const bool checkWellGroupControlsAndNetwork=false) const
Definition: BlackoilWellModel_impl.hpp:1557
int numStrictIterations() const
Definition: BlackoilWellModel.hpp:346
typename FluidSystem::IndexTraitsType IndexTraits
Definition: BlackoilWellModel.hpp:120
void updateCellRatesForDomain(int domainIndex, const std::map< std::string, int > &well_domain_map)
Definition: BlackoilWellModel_impl.hpp:1381
data::WellBlockAveragePressures wellBlockAveragePressures() const
Definition: BlackoilWellModel.hpp:248
void updateCellRates()
Definition: BlackoilWellModel_impl.hpp:1369
auto begin() const
Definition: BlackoilWellModel.hpp:339
std::size_t local_num_cells_
Definition: BlackoilWellModel.hpp:468
auto end() const
Definition: BlackoilWellModel.hpp:340
bool alternative_well_rate_init_
Definition: BlackoilWellModel.hpp:471
void timeStepSucceeded(const double simulationTime, const double dt)
Definition: BlackoilWellModel_impl.hpp:684
std::unique_ptr< WellType > createTypedWellPointer(const int wellID, const int time_step) const
Definition: BlackoilWellModel_impl.hpp:1077
Simulator & simulator_
Definition: BlackoilWellModel.hpp:440
void setReservoirCouplingSlave(ReservoirCouplingSlave< Scalar > *slave)
Definition: BlackoilWellModel.hpp:401
void createWellContainer(const int report_step) override
Definition: BlackoilWellModel_impl.hpp:841
std::unique_ptr< WellInterface< TypeTag > > WellInterfacePtr
Definition: BlackoilWellModel.hpp:193
ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave()
Definition: BlackoilWellModel.hpp:387
data::Wells wellData() const
Definition: BlackoilWellModel.hpp:213
const ModelParameters & param() const
Definition: BlackoilWellModel.hpp:436
void updateWellTestState(const double simulationTime, WellTestState &wellTestState)
upate the wellTestState related to economic limits
Definition: BlackoilWellModel_impl.hpp:1787
void addWellPressureEquationsDomain(PressureMatrix &jacobian, const BVector &weights, const bool use_well_weights, const int domainIndex) const
Definition: BlackoilWellModel.hpp:302
void addReservoirSourceTerms(GlobalEqVector &residual, const std::vector< typename SparseMatrixAdapter::MatrixBlock * > &diagMatAddress) const
Definition: BlackoilWellModel_impl.hpp:1466
int compressedIndexForInterior(int cartesian_cell_idx) const override
get compressed index for interior cells (-1, otherwise
Definition: BlackoilWellModel.hpp:349
static constexpr bool has_energy_
Definition: BlackoilWellModel.hpp:130
void recoverWellSolutionAndUpdateWellState(const BVector &x)
Definition: BlackoilWellModel_impl.hpp:1519
void addWellPressureEquationsStruct(PressureMatrix &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1496
void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1918
void computePotentials(const std::size_t widx, const WellState< Scalar, IndexTraits > &well_state_copy, std::string &exc_msg, ExceptionType::ExcEnum &exc_type, DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1868
void endReportStep()
Definition: BlackoilWellModel_impl.hpp:657
void prepareWellsBeforeAssembling(const double dt, DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1335
data::WellBlockAveragePressures computeWellBlockAveragePressures(const Scalar gravity) const
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
Definition: GroupStateHelper.hpp:59
void setReservoirCouplingMaster(ReservoirCouplingMaster< Scalar > *reservoir_coupling_master)
Definition: GroupStateHelper.hpp:245
void setReservoirCouplingSlave(ReservoirCouplingSlave< Scalar > *reservoir_coupling_slave)
Definition: GroupStateHelper.hpp:250
Handles computation and reporting of guide rates for wells and groups.
Definition: GuideRateHandler.hpp:53
Definition: RateConverter.hpp:70
Definition: RegionAverageCalculator.hpp:60
Definition: ReservoirCouplingMaster.hpp:38
Definition: ReservoirCouplingSlave.hpp:40
Definition: BlackoilWellModel.hpp:93
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:182
bool matrix_add_well_contributions_
Whether to add influences of wells between cells to the matrix and preconditioner matrix.
Definition: BlackoilModelParameters.hpp:288
bool use_multisegment_well_
Definition: BlackoilModelParameters.hpp:282
int strict_outer_iter_wells_
Newton iteration where wells are stricly convergent.
Definition: BlackoilModelParameters.hpp:238
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34