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
80
81#include <cstddef>
82#include <map>
83#include <memory>
84#include <optional>
85#include <string>
86#include <tuple>
87#include <vector>
88
89namespace Opm {
90
91template<class Scalar> class BlackoilWellModelNldd;
92template<class T, template <typename, typename...> class Storage> class SparseTable;
93
94#if COMPILE_GPU_BRIDGE
95template<class Scalar> class WellContributions;
96#endif
97
99 template<typename TypeTag>
100 class BlackoilWellModel : public WellConnectionAuxiliaryModule<TypeTag, BlackoilWellModel<TypeTag>>
101 , public BlackoilWellModelGeneric<GetPropType<TypeTag, Properties::Scalar>,
102 typename GetPropType<TypeTag, Properties::FluidSystem>::IndexTraitsType>
103 {
104 public:
105 // --------- Types ---------
117
119 using IndexTraits = typename FluidSystem::IndexTraitsType;
121
123
124 static const int numEq = Indices::numEq;
125 static const int solventSaturationIdx = Indices::solventSaturationIdx;
126 static constexpr bool has_solvent_ = getPropValue<TypeTag, Properties::EnableSolvent>();
127 static constexpr bool has_polymer_ = getPropValue<TypeTag, Properties::EnablePolymer>();
128 static constexpr bool has_energy_ = (getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal);
129 static constexpr bool has_micp_ = Indices::enableMICP;
130
131 // TODO: where we should put these types, WellInterface or Well Model?
132 // or there is some other strategy, like TypeTag
133 using VectorBlockType = Dune::FieldVector<Scalar, numEq>;
134 using BVector = Dune::BlockVector<VectorBlockType>;
135
138
139 // For the conversion between the surface volume rate and reservoir voidage rate
140 using RateConverterType = RateConverter::
141 SurfaceToReservoirVoidage<FluidSystem, std::vector<int> >;
142
143 // For computing average pressured used by gpmaint
144 using AverageRegionalPressureType = RegionAverageCalculator::
145 AverageRegionalPressure<FluidSystem, std::vector<int> >;
146
148
149 void init();
150 void initWellContainer(const int reportStepIdx) override;
151
153 {
154 OPM_TIMEBLOCK(beginEpsiode);
155 beginReportStep(simulator_.episodeIndex());
156 }
157
158 void beginTimeStep();
159
161 {
162 OPM_TIMEBLOCK(beginIteration);
163 assemble(simulator_.model().newtonMethod().numIterations(),
164 simulator_.timeStepSize());
165 }
166
168 { }
169
171 {
172 OPM_TIMEBLOCK(endTimeStep);
173 timeStepSucceeded(simulator_.time(), simulator_.timeStepSize());
174 }
175
177 {
179 }
180
182 unsigned globalIdx) const;
183
184 template <class Context>
186 const Context& context,
187 unsigned spaceIdx,
188 unsigned timeIdx) const;
189
190
191 using WellInterfacePtr = std::unique_ptr<WellInterface<TypeTag>>;
192
194 void initFromRestartFile(const RestartValue& restartValues)
195 {
196 initFromRestartFile(restartValues,
197 this->simulator_.vanguard().transferWTestState(),
198 grid().size(0),
200 this->simulator_.vanguard().enableDistributedWells());
201 }
202
204 void prepareDeserialize(const int report_step)
205 {
206 prepareDeserialize(report_step, grid().size(0),
208 this->simulator_.vanguard().enableDistributedWells());
209 }
210
211 data::Wells wellData() const
212 {
213 auto wsrpt = this->wellState()
214 .report(this->simulator_.vanguard().globalCell().data(),
215 [this](const int well_index)
216 {
217 return this->wasDynamicallyShutThisTimeStep(well_index);
218 });
219
221 .assignWellGuideRates(wsrpt, this->reportStepIndex());
222
223 this->assignWellTracerRates(wsrpt);
224
225 if (const auto& rspec = eclState().runspec();
226 rspec.co2Storage() || rspec.h2Storage())
227 {
228 // The gas reference density (surface condition) is the
229 // same for all PVT regions in CO2STORE/H2STORE runs so,
230 // for simplicity, we use region zero (0) here.
231
232 this->assignMassGasRate(wsrpt, FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, 0));
233 }
234
235 this->assignWellTargets(wsrpt);
236
237 this->assignDynamicWellStatus(wsrpt);
238
239 // Assigning (a subset of the) property values in shut
240 // connections should be the last step of wellData().
241 this->assignShutConnections(wsrpt, this->reportStepIndex());
242
243 return wsrpt;
244 }
245
246 data::WellBlockAveragePressures wellBlockAveragePressures() const
247 {
249 }
250
251#if COMPILE_GPU_BRIDGE
252 // accumulate the contributions of all Wells in the WellContributions object
253 void getWellContributions(WellContributions<Scalar>& x) const;
254#endif
255
256 // Check if well equations is converged.
257 ConvergenceReport getWellConvergence(const std::vector<Scalar>& B_avg, const bool checkWellGroupControlsAndNetwork = false) const;
258
259 const SimulatorReportSingle& lastReport() const;
260
261 void addWellContributions(SparseMatrixAdapter& jacobian) const;
262
263 // add source from wells to the reservoir matrix
265 const std::vector<typename SparseMatrixAdapter::MatrixBlock*>& diagMatAddress) const;
266
267 // called at the beginning of a report step
268 void beginReportStep(const int time_step);
269
270 // it should be able to go to prepareTimeStep(), however, the updateWellControls()
271 // makes it a little more difficult. unless we introduce if (iterationIdx != 0) to avoid doing the above function
272 // twice at the beginning of the time step
275 void calculateExplicitQuantities(DeferredLogger& deferred_logger) const;
276 // some preparation work, mostly related to group control and RESV,
277 // at the beginning of each time step (Not report step)
278 void prepareTimeStep(DeferredLogger& deferred_logger);
279
280 bool
281 updateWellControls(DeferredLogger& deferred_logger);
282
283 std::tuple<bool, Scalar>
284 updateNetworks(const bool mandatory_network_balance, DeferredLogger& deferred_logger, const bool relax_network_tolerance = false);
285
286
287 void updateAndCommunicate(const int reportStepIdx,
288 const int iterationIdx,
289 DeferredLogger& deferred_logger);
290
291 bool updateGroupControls(const Group& group,
292 DeferredLogger& deferred_logger,
293 const int reportStepIdx,
294 const int iterationIdx);
295
296 const WellInterface<TypeTag>& getWell(const std::string& well_name) const;
297
298 using PressureMatrix = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, 1, 1>>;
299
301 const BVector& weights,
302 const bool use_well_weights) const;
305 const BVector& weights,
306 const bool use_well_weights,
307 const int domainIndex) const
308 {
309 if (!nldd_) {
310 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
311 }
312 return nldd_->addWellPressureEquations(jacobian,
313 weights,
314 use_well_weights,
315 domainIndex);
316 }
317
319 const std::vector<WellInterfacePtr>& localNonshutWells() const
320 {
321 return well_container_;
322 }
323
325 {
326 if (!nldd_) {
327 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
328 }
329 return nldd_->well_local_cells();
330 }
331
332 const std::map<std::string, int>& well_domain() const
333 {
334 if (!nldd_) {
335 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
336 }
337
338 return nldd_->well_domain();
339 }
340
341 auto begin() const { return well_container_.begin(); }
342 auto end() const { return well_container_.end(); }
343 bool empty() const { return well_container_.empty(); }
344
347
350
351 int compressedIndexForInterior(int cartesian_cell_idx) const override
352 {
353 return simulator_.vanguard().compressedIndexForInterior(cartesian_cell_idx);
354 }
355
356 int compressedIndexForInteriorLGR(const std::string& lgr_tag, const Connection& conn) const override
357 {
358 return simulator_.vanguard().compressedIndexForInteriorLGR(lgr_tag, conn);
359 }
360
361 // using the solution x to recover the solution xw for wells and applying
362 // xw to update Well State
364
365 // using the solution x to recover the solution xw for wells and applying
366 // xw to update Well State
368 const int domainIdx);
369 // Update cellRates_ with contributions from all wells
370 void updateCellRates();
371
372 // Update cellRates_ with contributions from wells in a specific domain
373 void updateCellRatesForDomain(int domainIndex,
374 const std::map<std::string, int>& well_domain_map);
375
376 const Grid& grid() const
377 { return simulator_.vanguard().grid(); }
378
379 const Simulator& simulator() const
380 { return simulator_; }
381
383 { nldd_ = mod; }
384
385#ifdef RESERVOIR_COUPLING_ENABLED
387 return *(this->simulator_.reservoirCouplingMaster());
388 }
390 return *(this->simulator_.reservoirCouplingSlave());
391 }
393 return this->simulator_.reservoirCouplingMaster() != nullptr;
394 }
396 return this->simulator_.reservoirCouplingSlave() != nullptr;
397 }
399 {
400 this->guide_rate_handler_.setReservoirCouplingMaster(master);
401 this->wgHelper().setReservoirCouplingMaster(master);
402 }
404 {
405 this->guide_rate_handler_.setReservoirCouplingSlave(slave);
406 this->wgHelper().setReservoirCouplingSlave(slave);
407 }
408
411
414
423 std::optional<ReservoirCoupling::ScopedLoggerGuard>
425
426 #endif
427 protected:
429
430 // a vector of all the wells.
431 std::vector<WellInterfacePtr> well_container_{};
432
433 std::vector<bool> is_cell_perforated_{};
434
435 void initializeWellState(const int timeStepIdx);
436
437 // create the well container
438 void createWellContainer(const int report_step) override;
439
441 createWellPointer(const int wellID,
442 const int report_step) const;
443
444 template <typename WellType>
445 std::unique_ptr<WellType>
446 createTypedWellPointer(const int wellID,
447 const int time_step) const;
448
449 WellInterfacePtr createWellForWellTest(const std::string& well_name,
450 const int report_step,
451 DeferredLogger& deferred_logger) const;
452
454 std::size_t global_num_cells_{};
455 // the number of the cells in the local grid
456 std::size_t local_num_cells_{};
458 std::vector<Scalar> depth_{};
460 std::map<std::string, Scalar> well_group_thp_calc_;
461 std::unique_ptr<RateConverterType> rateConverter_{};
462 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>> regionalAveragePressureCalculator_{};
463
466
467 // A flag to tell the convergence report whether we need to take another newton step
469
470 // Pre-step network solve at static reservoir conditions (group and well states might be updated)
471 void doPreStepNetworkRebalance(DeferredLogger& deferred_logger);
472
473 std::vector<Scalar> B_avg_{};
474
475 const EquilGrid& equilGrid() const
476 { return simulator_.vanguard().equilGrid(); }
477
478 const EclipseState& eclState() const
479 { return simulator_.vanguard().eclState(); }
480
481 // compute the well fluxes and assemble them in to the reservoir equations as source terms
482 // and in the well equations.
483 void assemble(const int iterationIdx,
484 const double dt);
485
486 // well controls and network pressures affect each other and are solved in an iterative manner.
487 // the function handles one iteration of updating well controls and network pressures.
488 // it is possible to decouple the update of well controls and network pressures further.
489 // the returned two booleans are {continue_due_to_network, well_group_control_changed}, respectively
490 std::tuple<bool, bool, Scalar> updateWellControlsAndNetworkIteration(const bool mandatory_network_balance,
491 const bool relax_network_tolerance,
492 const bool optimize_gas_lift,
493 const double dt,
494 DeferredLogger& local_deferredLogger);
495
496 bool updateWellControlsAndNetwork(const bool mandatory_network_balance,
497 const double dt,
498 DeferredLogger& local_deferredLogger);
499
500 bool computeWellGroupThp(const double dt, DeferredLogger& local_deferredLogger);
501
510 void initializeLocalWellStructure(const int reportStepIdx,
511 const bool enableWellPIScaling);
512
516 void initializeGroupStructure(const int reportStepIdx);
517
518 // called at the end of a time step
519 void timeStepSucceeded(const double simulationTime, const double dt);
520
521 // called at the end of a report step
522 void endReportStep();
523
524 // setting the well_solutions_ based on well_state.
525 void updatePrimaryVariables(DeferredLogger& deferred_logger);
526
528
529 void computePotentials(const std::size_t widx,
530 const WellState<Scalar, IndexTraits>& well_state_copy,
531 std::string& exc_msg,
532 ExceptionType::ExcEnum& exc_type,
533 DeferredLogger& deferred_logger) override;
534
535 const std::vector<Scalar>& wellPerfEfficiencyFactors() const;
536
537 void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger& deferred_logger) override;
538 void calculateProductivityIndexValues(DeferredLogger& deferred_logger) override;
540 DeferredLogger& deferred_logger);
541
542 // The number of conservation quantities.
543 int numConservationQuantities() const;
544
545 int reportStepIndex() const;
546
547 void assembleWellEq(const double dt, DeferredLogger& deferred_logger);
548
549 void prepareWellsBeforeAssembling(const double dt, DeferredLogger& deferred_logger);
550
551 // TODO: finding a better naming
552 void assembleWellEqWithoutIteration(const double dt, DeferredLogger& deferred_logger);
553
555
556 void extractLegacyDepth_();
557
559 void updateWellTestState(const double simulationTime, WellTestState& wellTestState);
560
561 void wellTesting(const int timeStepIdx, const double simulationTime, DeferredLogger& deferred_logger);
562
563 void calcResvCoeff(const int fipnum,
564 const int pvtreg,
565 const std::vector<Scalar>& production_rates,
566 std::vector<Scalar>& resv_coeff) override;
567
568 void calcInjResvCoeff(const int fipnum,
569 const int pvtreg,
570 std::vector<Scalar>& resv_coeff) override;
571
573
574 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:95
BlackoilWellModelWBP< GetPropType< TypeTag, Properties::Scalar >, GetPropType< TypeTag, Properties::FluidSystem >::IndexTraitsType > wbp_
Definition: BlackoilWellModelGeneric.hpp:530
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 model in a NLDD solver.
Definition: BlackoilWellModelNldd.hpp:80
Class for handling the blackoil well model.
Definition: BlackoilWellModel.hpp:103
void initializeGroupStructure(const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:297
bool network_needs_more_balancing_force_another_newton_iteration_
Definition: BlackoilWellModel.hpp:468
ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster()
Definition: BlackoilWellModel.hpp:386
void endTimeStep()
Definition: BlackoilWellModel.hpp:170
void prepareTimeStep(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:2263
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:1315
bool updateGroupControls(const Group &group, DeferredLogger &deferred_logger, const int reportStepIdx, const int iterationIdx)
Definition: BlackoilWellModel_impl.hpp:2043
WellInterfacePtr createWellPointer(const int wellID, const int report_step) const
Definition: BlackoilWellModel_impl.hpp:1083
GuideRateHandler< Scalar, IndexTraits > guide_rate_handler_
Definition: BlackoilWellModel.hpp:465
const EclipseState & eclState() const
Definition: BlackoilWellModel.hpp:478
void calcInjResvCoeff(const int fipnum, const int pvtreg, std::vector< Scalar > &resv_coeff) override
Definition: BlackoilWellModel_impl.hpp:2453
const std::vector< WellInterfacePtr > & localNonshutWells() const
Get list of local nonshut wells.
Definition: BlackoilWellModel.hpp:319
void prepareDeserialize(const int report_step)
Definition: BlackoilWellModel.hpp:204
void init()
Definition: BlackoilWellModel_impl.hpp:163
void setNlddAdapter(BlackoilWellModelNldd< TypeTag > *mod)
Definition: BlackoilWellModel.hpp:382
void calcResvCoeff(const int fipnum, const int pvtreg, const std::vector< Scalar > &production_rates, std::vector< Scalar > &resv_coeff) override
Definition: BlackoilWellModel_impl.hpp:2442
const Simulator & simulator() const
Definition: BlackoilWellModel.hpp:379
std::vector< Scalar > depth_
Definition: BlackoilWellModel.hpp:458
void updateAndCommunicate(const int reportStepIdx, const int iterationIdx, DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:2006
std::size_t global_num_cells_
Definition: BlackoilWellModel.hpp:454
std::map< std::string, Scalar > well_group_thp_calc_
Definition: BlackoilWellModel.hpp:460
void doPreStepNetworkRebalance(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1159
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilWellModel.hpp:112
SimulatorReportSingle last_report_
Definition: BlackoilWellModel.hpp:464
static const int solventSaturationIdx
Definition: BlackoilWellModel.hpp:125
Scalar gravity_
Definition: BlackoilWellModel.hpp:457
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:2413
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: BlackoilWellModel.hpp:108
void beginIteration()
Definition: BlackoilWellModel.hpp:160
bool isReservoirCouplingMaster() const
Definition: BlackoilWellModel.hpp:392
Dune::FieldVector< Scalar, numEq > VectorBlockType
Definition: BlackoilWellModel.hpp:133
void initFromRestartFile(const RestartValue &restartValues)
Definition: BlackoilWellModel.hpp:194
const ModelParameters param_
Definition: BlackoilWellModel.hpp:453
void beginEpisode()
Definition: BlackoilWellModel.hpp:152
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: BlackoilWellModel.hpp:109
GetPropType< TypeTag, Properties::Grid > Grid
Definition: BlackoilWellModel.hpp:106
GetPropType< TypeTag, Properties::EquilGrid > EquilGrid
Definition: BlackoilWellModel.hpp:107
int numConservationQuantities() const
Definition: BlackoilWellModel_impl.hpp:2385
void assembleWellEq(const double dt, DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1568
const std::map< std::string, int > & well_domain() const
Definition: BlackoilWellModel.hpp:332
bool updateWellControls(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1868
void endIteration()
Definition: BlackoilWellModel.hpp:167
static constexpr bool has_polymer_
Definition: BlackoilWellModel.hpp:127
int reportStepIndex() const
Definition: BlackoilWellModel_impl.hpp:2430
void calculateProductivityIndexValues(DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:2205
void extractLegacyDepth_()
Definition: BlackoilWellModel_impl.hpp:2401
void extractLegacyCellPvtRegionIndex_()
Definition: BlackoilWellModel_impl.hpp:2369
void recoverWellSolutionAndUpdateWellStateDomain(const BVector &x, const int domainIdx)
Definition: BlackoilWellModel_impl.hpp:1788
void updateAverageFormationFactor()
Definition: BlackoilWellModel_impl.hpp:2312
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilWellModel.hpp:111
bool isReservoirCouplingSlave() const
Definition: BlackoilWellModel.hpp:395
void initializeWellState(const int timeStepIdx)
Definition: BlackoilWellModel_impl.hpp:809
void assemble(const int iterationIdx, const double dt)
Definition: BlackoilWellModel_impl.hpp:1189
const Grid & grid() const
Definition: BlackoilWellModel.hpp:376
void computeWellTemperature()
Definition: BlackoilWellModel_impl.hpp:2464
static const int numEq
Definition: BlackoilWellModel.hpp:124
void addWellPressureEquations(PressureMatrix &jacobian, const BVector &weights, const bool use_well_weights) const
Definition: BlackoilWellModel_impl.hpp:1688
const SimulatorReportSingle & lastReport() const
Definition: BlackoilWellModel_impl.hpp:685
bool updateWellControlsAndNetwork(const bool mandatory_network_balance, const double dt, DeferredLogger &local_deferredLogger)
Definition: BlackoilWellModel_impl.hpp:1256
void endEpisode()
Definition: BlackoilWellModel.hpp:176
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: BlackoilWellModel.hpp:115
void addWellContributions(SparseMatrixAdapter &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1678
const EquilGrid & equilGrid() const
Definition: BlackoilWellModel.hpp:475
std::tuple< bool, Scalar > updateNetworks(const bool mandatory_network_balance, DeferredLogger &deferred_logger, const bool relax_network_tolerance=false)
Definition: BlackoilWellModel_impl.hpp:1949
WellInterfacePtr createWellForWellTest(const std::string &well_name, const int report_step, DeferredLogger &deferred_logger) const
Definition: BlackoilWellModel_impl.hpp:1134
GetPropType< TypeTag, Properties::Indices > Indices
Definition: BlackoilWellModel.hpp:110
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:356
void assembleWellEqWithoutIteration(const double dt, DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1594
void sendSlaveGroupDataToMaster()
Send comprehensive slave group data to master.
std::vector< WellInterfacePtr > well_container_
Definition: BlackoilWellModel.hpp:431
void receiveSlaveGroupData()
Receive comprehensive slave group data from slaves.
static constexpr std::size_t pressureVarIndex
Definition: BlackoilWellModel.hpp:122
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: BlackoilWellModel.hpp:298
bool empty() const
Definition: BlackoilWellModel.hpp:343
void computeTotalRatesForDof(RateVector &rate, unsigned globalIdx) const
Definition: BlackoilWellModel_impl.hpp:772
bool addMatrixContributions() const
Definition: BlackoilWellModel.hpp:345
void beginTimeStep()
Definition: BlackoilWellModel_impl.hpp:333
std::vector< Scalar > B_avg_
Definition: BlackoilWellModel.hpp:473
const std::vector< Scalar > & wellPerfEfficiencyFactors() const
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: BlackoilWellModel.hpp:113
void calculateExplicitQuantities(DeferredLogger &deferred_logger) const
Definition: BlackoilWellModel_impl.hpp:1853
std::map< std::string, std::unique_ptr< AverageRegionalPressureType > > regionalAveragePressureCalculator_
Definition: BlackoilWellModel.hpp:462
void updatePrimaryVariables(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:2360
bool computeWellGroupThp(const double dt, DeferredLogger &local_deferredLogger)
Definition: BlackoilWellModel_impl.hpp:1382
GetPropType< TypeTag, Properties::GlobalEqVector > GlobalEqVector
Definition: BlackoilWellModel.hpp:114
void initializeLocalWellStructure(const int reportStepIdx, const bool enableWellPIScaling)
Definition: BlackoilWellModel_impl.hpp:252
static constexpr bool has_solvent_
Definition: BlackoilWellModel.hpp:126
void setReservoirCouplingMaster(ReservoirCouplingMaster< Scalar > *master)
Definition: BlackoilWellModel.hpp:398
Dune::BlockVector< VectorBlockType > BVector
Definition: BlackoilWellModel.hpp:134
static constexpr bool has_micp_
Definition: BlackoilWellModel.hpp:129
std::unique_ptr< RateConverterType > rateConverter_
Definition: BlackoilWellModel.hpp:461
BlackoilWellModel(Simulator &simulator)
Definition: BlackoilWellModel_impl.hpp:78
const SparseTable< int > & well_local_cells() const
Definition: BlackoilWellModel.hpp:324
void wellTesting(const int timeStepIdx, const double simulationTime, DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:597
std::vector< bool > is_cell_perforated_
Definition: BlackoilWellModel.hpp:433
ConvergenceReport getWellConvergence(const std::vector< Scalar > &B_avg, const bool checkWellGroupControlsAndNetwork=false) const
Definition: BlackoilWellModel_impl.hpp:1801
int numStrictIterations() const
Definition: BlackoilWellModel.hpp:348
typename FluidSystem::IndexTraitsType IndexTraits
Definition: BlackoilWellModel.hpp:119
void updateCellRatesForDomain(int domainIndex, const std::map< std::string, int > &well_domain_map)
Definition: BlackoilWellModel_impl.hpp:1625
data::WellBlockAveragePressures wellBlockAveragePressures() const
Definition: BlackoilWellModel.hpp:246
void updateCellRates()
Definition: BlackoilWellModel_impl.hpp:1613
auto begin() const
Definition: BlackoilWellModel.hpp:341
std::size_t local_num_cells_
Definition: BlackoilWellModel.hpp:456
auto end() const
Definition: BlackoilWellModel.hpp:342
bool alternative_well_rate_init_
Definition: BlackoilWellModel.hpp:459
void timeStepSucceeded(const double simulationTime, const double dt)
Definition: BlackoilWellModel_impl.hpp:695
std::unique_ptr< WellType > createTypedWellPointer(const int wellID, const int time_step) const
Definition: BlackoilWellModel_impl.hpp:1103
Simulator & simulator_
Definition: BlackoilWellModel.hpp:428
void setReservoirCouplingSlave(ReservoirCouplingSlave< Scalar > *slave)
Definition: BlackoilWellModel.hpp:403
void createWellContainer(const int report_step) override
Definition: BlackoilWellModel_impl.hpp:852
std::unique_ptr< WellInterface< TypeTag > > WellInterfacePtr
Definition: BlackoilWellModel.hpp:191
ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave()
Definition: BlackoilWellModel.hpp:389
data::Wells wellData() const
Definition: BlackoilWellModel.hpp:211
void updateWellTestState(const double simulationTime, WellTestState &wellTestState)
upate the wellTestState related to economic limits
Definition: BlackoilWellModel_impl.hpp:2088
void addWellPressureEquationsDomain(PressureMatrix &jacobian, const BVector &weights, const bool use_well_weights, const int domainIndex) const
Definition: BlackoilWellModel.hpp:304
void addReservoirSourceTerms(GlobalEqVector &residual, const std::vector< typename SparseMatrixAdapter::MatrixBlock * > &diagMatAddress) const
Definition: BlackoilWellModel_impl.hpp:1710
int compressedIndexForInterior(int cartesian_cell_idx) const override
get compressed index for interior cells (-1, otherwise
Definition: BlackoilWellModel.hpp:351
static constexpr bool has_energy_
Definition: BlackoilWellModel.hpp:128
void recoverWellSolutionAndUpdateWellState(const BVector &x)
Definition: BlackoilWellModel_impl.hpp:1763
void addWellPressureEquationsStruct(PressureMatrix &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1740
void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:2219
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:2169
void endReportStep()
Definition: BlackoilWellModel_impl.hpp:668
void prepareWellsBeforeAssembling(const double dt, DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1580
data::WellBlockAveragePressures computeWellBlockAveragePressures(const Scalar gravity) const
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
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:92
Definition: WellConnectionAuxiliaryModule.hpp:39
Definition: WellContributions.hpp:51
Definition: WellGroupHelper.hpp:59
void setReservoirCouplingSlave(ReservoirCouplingSlave< Scalar > *reservoir_coupling_slave)
Definition: WellGroupHelper.hpp:254
void setReservoirCouplingMaster(ReservoirCouplingMaster< Scalar > *reservoir_coupling_master)
Definition: WellGroupHelper.hpp:249
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:180
bool matrix_add_well_contributions_
Whether to add influences of wells between cells to the matrix and preconditioner matrix.
Definition: BlackoilModelParameters.hpp:286
bool use_multisegment_well_
Definition: BlackoilModelParameters.hpp:280
int strict_outer_iter_wells_
Newton iteration where wells are stricly convergent.
Definition: BlackoilModelParameters.hpp:236
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34