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_.timeStepSize());
160 }
161
163 { }
164
166 {
167 OPM_TIMEBLOCK(endTimeStep);
168 timeStepSucceeded(simulator_.time(), simulator_.timeStepSize());
169 }
170
172 {
174 }
175
177 unsigned globalIdx) const;
178
179 template <class Context>
181 const Context& context,
182 unsigned spaceIdx,
183 unsigned timeIdx) const;
184
185
186 using WellInterfacePtr = std::unique_ptr<WellInterface<TypeTag>>;
187
189 void initFromRestartFile(const RestartValue& restartValues)
190 {
191 initFromRestartFile(restartValues,
192 this->simulator_.vanguard().transferWTestState(),
193 grid().size(0),
195 this->simulator_.vanguard().enableDistributedWells());
196 }
197
199 void prepareDeserialize(const int report_step)
200 {
201 prepareDeserialize(report_step, grid().size(0),
203 this->simulator_.vanguard().enableDistributedWells());
204 }
205
206 data::Wells wellData() const
207 {
208 auto wsrpt = this->wellState()
209 .report(this->simulator_.vanguard().globalCell().data(),
210 [this](const int well_index)
211 {
212 return this->wasDynamicallyShutThisTimeStep(well_index);
213 });
214
216 .assignWellGuideRates(wsrpt, this->reportStepIndex());
217
218 this->assignWellTracerRates(wsrpt);
219
220 if (const auto& rspec = eclState().runspec();
221 rspec.co2Storage() || rspec.h2Storage())
222 {
223 // The gas reference density (surface condition) is the
224 // same for all PVT regions in CO2STORE/H2STORE runs so,
225 // for simplicity, we use region zero (0) here.
226
227 this->assignMassGasRate(wsrpt, FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, 0));
228 }
229
230 this->assignWellTargets(wsrpt);
231
232 this->assignDynamicWellStatus(wsrpt);
233
234 // Assigning (a subset of the) property values in shut
235 // connections should be the last step of wellData().
236 this->assignShutConnections(wsrpt, this->reportStepIndex());
237
238 return wsrpt;
239 }
240
241 data::WellBlockAveragePressures wellBlockAveragePressures() const
242 {
244 }
245
246#if COMPILE_GPU_BRIDGE
247 // accumulate the contributions of all Wells in the WellContributions object
248 void getWellContributions(WellContributions<Scalar>& x) const;
249#endif
250
251 // Check if well equations is converged.
252 ConvergenceReport getWellConvergence(const std::vector<Scalar>& B_avg, const bool checkWellGroupControlsAndNetwork = false) const;
253
254 const SimulatorReportSingle& lastReport() const;
255
256 void addWellContributions(SparseMatrixAdapter& jacobian) const;
257
258 // add source from wells to the reservoir matrix
260 const std::vector<typename SparseMatrixAdapter::MatrixBlock*>& diagMatAddress) const;
261
262 // called at the beginning of a report step
263 void beginReportStep(const int time_step);
264
268 void calculateExplicitQuantities() const;
269
272 void prepareTimeStep(DeferredLogger& deferred_logger);
273
274 bool
275 updateWellControls(DeferredLogger& deferred_logger);
276
277 void updateAndCommunicate(const int reportStepIdx);
278
279 bool updateGroupControls(const Group& group,
280 DeferredLogger& deferred_logger,
281 const int reportStepIdx);
282
283 const WellInterface<TypeTag>& getWell(const std::string& well_name) const;
284
285 using PressureMatrix = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, 1, 1>>;
286
288 const BVector& weights,
289 const bool use_well_weights) const;
292 const BVector& weights,
293 const bool use_well_weights,
294 const int domainIndex) const
295 {
296 if (!nldd_) {
297 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
298 }
299 return nldd_->addWellPressureEquations(jacobian,
300 weights,
301 use_well_weights,
302 domainIndex);
303 }
304
306 const std::vector<WellInterfacePtr>& localNonshutWells() const
307 {
308 return well_container_;
309 }
310
312 {
313 if (!nldd_) {
314 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
315 }
316 return nldd_->well_local_cells();
317 }
318
319 const std::map<std::string, int>& well_domain() const
320 {
321 if (!nldd_) {
322 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
323 }
324
325 return nldd_->well_domain();
326 }
327
328 auto begin() const { return well_container_.begin(); }
329 auto end() const { return well_container_.end(); }
330 bool empty() const { return well_container_.empty(); }
331
334
337
338 int compressedIndexForInterior(int cartesian_cell_idx) const override
339 {
340 return simulator_.vanguard().compressedIndexForInterior(cartesian_cell_idx);
341 }
342
343 int compressedIndexForInteriorLGR(const std::string& lgr_tag, const Connection& conn) const override
344 {
345 return simulator_.vanguard().compressedIndexForInteriorLGR(lgr_tag, conn);
346 }
347
348 // using the solution x to recover the solution xw for wells and applying
349 // xw to update Well State
351
352 // using the solution x to recover the solution xw for wells and applying
353 // xw to update Well State
355 const int domainIdx);
356 // Update cellRates_ with contributions from all wells
357 void updateCellRates();
358
359 // Update cellRates_ with contributions from wells in a specific domain
360 void updateCellRatesForDomain(int domainIndex,
361 const std::map<std::string, int>& well_domain_map);
362
363 const Grid& grid() const
364 { return simulator_.vanguard().grid(); }
365
366 const Simulator& simulator() const
367 { return simulator_; }
368
370 { nldd_ = mod; }
371
372 // === Reservoir Coupling ===
373
377
385 void updateGuideRates(const int report_step_idx,
386 const double sim_time)
387 {
388 this->guide_rate_handler_.updateGuideRates(
389 report_step_idx, sim_time, this->wellState(), this->groupState()
390 );
391 }
392
394 bool isReservoirCouplingMaster() const { return rescoup_.isMaster(); }
395
397 bool isReservoirCouplingSlave() const { return rescoup_.isSlave(); }
398
402 return rescoup_.master();
403 }
404
408 return rescoup_.slave();
409 }
410
411#ifdef RESERVOIR_COUPLING_ENABLED
412 void setReservoirCouplingMaster(ReservoirCouplingMaster<Scalar>* master)
413 {
414 rescoup_.setMaster(master);
415 this->guide_rate_handler_.setReservoirCouplingMaster(master);
416 this->groupStateHelper().setReservoirCouplingMaster(master);
417 }
418 void setReservoirCouplingSlave(ReservoirCouplingSlave<Scalar>* slave)
419 {
420 rescoup_.setSlave(slave);
421 this->guide_rate_handler_.setReservoirCouplingSlave(slave);
422 this->groupStateHelper().setReservoirCouplingSlave(slave);
423 }
424
426 void sendSlaveGroupDataToMaster();
427
429 void receiveSlaveGroupData();
430
431 void receiveGroupConstraintsFromMaster();
432 void sendMasterGroupConstraintsToSlaves();
433
442 std::optional<ReservoirCoupling::ScopedLoggerGuard>
443 setupRescoupScopedLogger(DeferredLogger& local_logger);
444#endif
445
446 bool updateWellControlsAndNetwork(const bool mandatory_network_balance,
447 const double dt,
448 DeferredLogger& local_deferredLogger);
449
450 // TODO: finding a better naming
451 void assembleWellEqWithoutIteration(const double dt);
452
453 const std::vector<Scalar>& B_avg() const
454 { return B_avg_; }
455
456 const ModelParameters& param() const
457 { return param_; }
458
459
460 template<class FluidState, class SingleWellState>
461 static Scalar computeTemperatureWeightFactor(const int perf_index, const int np, const FluidState& fs, const SingleWellState& ws)
462 {
463 const auto& perf_phase_rate = ws.perf_data.phase_rates;
464 // we only have one temperature pr cell any phaseIdx will do
465 Scalar cellTemperatures = fs.temperature(/*phaseIdx*/0).value();
466 Scalar weight_factor = 0.0;
467 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
468 if (!FluidSystem::phaseIsActive(phaseIdx)) {
469 continue;
470 }
471 Scalar cellInternalEnergy = fs.enthalpy(phaseIdx).value() -
472 fs.pressure(phaseIdx).value() / fs.density(phaseIdx).value();
473 Scalar cellBinv = fs.invB(phaseIdx).value();
474 Scalar cellDensity = fs.density(phaseIdx).value();
475 Scalar perfPhaseRate = perf_phase_rate[perf_index*np + phaseIdx];
476 weight_factor += cellDensity * (perfPhaseRate / cellBinv) * (cellInternalEnergy / cellTemperatures);
477 }
478 return (std::abs(weight_factor) + 1e-13);
479 }
480
481 protected:
483
484 // a vector of all the wells.
485 std::vector<WellInterfacePtr> well_container_{};
486
487 std::vector<bool> is_cell_perforated_{};
488
489 void initializeWellState(const int timeStepIdx);
490
491 // create the well container
492 void createWellContainer(const int report_step) override;
493
495 createWellPointer(const int wellID,
496 const int report_step) const;
497
498 template <typename WellType>
499 std::unique_ptr<WellType>
500 createTypedWellPointer(const int wellID,
501 const int time_step) const;
502
503 WellInterfacePtr createWellForWellTest(const std::string& well_name,
504 const int report_step,
505 DeferredLogger& deferred_logger) const;
506
508 std::size_t global_num_cells_{};
509 // the number of the cells in the local grid
510 std::size_t local_num_cells_{};
512 std::vector<Scalar> depth_{};
514 std::unique_ptr<RateConverterType> rateConverter_{};
515 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>> regionalAveragePressureCalculator_{};
516
520
521 // A flag to tell the convergence report whether we need to take another newton step
523
524 std::vector<Scalar> B_avg_{};
525
526 const EquilGrid& equilGrid() const
527 { return simulator_.vanguard().equilGrid(); }
528
529 const EclipseState& eclState() const
530 { return simulator_.vanguard().eclState(); }
531
532 // compute the well fluxes and assemble them in to the reservoir equations as source terms
533 // and in the well equations.
534 void assemble(const double dt);
535
536 // well controls and network pressures affect each other and are solved in an iterative manner.
537 // the function handles one iteration of updating well controls and network pressures.
538 // it is possible to decouple the update of well controls and network pressures further.
539 // the returned two booleans are {continue_due_to_network, well_group_control_changed}, respectively
540 std::tuple<bool, bool, Scalar> updateWellControlsAndNetworkIteration(const bool mandatory_network_balance,
541 const bool relax_network_tolerance,
542 const bool optimize_gas_lift,
543 const double dt,
544 DeferredLogger& local_deferredLogger);
545
554 void initializeLocalWellStructure(const int reportStepIdx,
555 const bool enableWellPIScaling);
556
560 void initializeGroupStructure(const int reportStepIdx);
561
562 // called at the end of a time step
563 void timeStepSucceeded(const double simulationTime, const double dt);
564
565 // called at the end of a report step
566 void endReportStep();
567
568 // setting the well_solutions_ based on well_state.
570
572
573 void computePotentials(const std::size_t widx,
574 const WellState<Scalar, IndexTraits>& well_state_copy,
575 std::string& exc_msg,
576 ExceptionType::ExcEnum& exc_type) override;
577
578 const std::vector<Scalar>& wellPerfEfficiencyFactors() const;
579
580 void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger& deferred_logger) override;
581 void calculateProductivityIndexValues(DeferredLogger& deferred_logger) override;
583 DeferredLogger& deferred_logger);
584
585 // The number of conservation quantities.
586 int numConservationQuantities() const;
587
588 int reportStepIndex() const;
589
590 void assembleWellEq(const double dt);
591
592 void prepareWellsBeforeAssembling(const double dt);
593
595
596 void extractLegacyDepth_();
597
599 void updateWellTestState(const double simulationTime, WellTestState& wellTestState);
600
601 void wellTesting(const int timeStepIdx, const double simulationTime, DeferredLogger& deferred_logger);
602
603 void calcResvCoeff(const int fipnum,
604 const int pvtreg,
605 const std::vector<Scalar>& production_rates,
606 std::vector<Scalar>& resv_coeff) const override;
607
608 void calcInjResvCoeff(const int fipnum,
609 const int pvtreg,
610 std::vector<Scalar>& resv_coeff) const override;
611
613
614 private:
617 BlackoilWellModelNldd<TypeTag>* nldd_ = nullptr;
618
619 // These members are used to avoid reallocation in specific functions
620 // instead of using local variables.
621 // Their state is not relevant between function calls, so they can
622 // (and must) be mutable, as the functions using them are const.
623 mutable BVector x_local_;
624
625 // Store cell rates after assembling to avoid iterating all wells and connections for every element
626 std::map<int, RateVector> cellRates_;
627
628 void assignWellTracerRates(data::Wells& wsrpt) const;
629 };
630
631} // namespace Opm
632
634
635#endif // OPM_BLACKOILWELLMODEL_HEADER_INCLUDED
Contains the high level supplements required to extend the black oil model by bioeffects.
Definition: blackoilbioeffectsmodules.hh:95
Contains the high level supplements required to extend the black oil model by polymer.
Definition: blackoilpolymermodules.hh:65
Class for handling the gaslift in the blackoil well model.
Definition: BlackoilWellModelGasLift.hpp:96
Class for handling the blackoil well model.
Definition: BlackoilWellModelGeneric.hpp:97
const GroupState< GetPropType< TypeTag, Properties::Scalar > > & groupState() const
Definition: BlackoilWellModelGeneric.hpp:142
BlackoilWellModelWBP< GetPropType< TypeTag, Properties::Scalar >, GetPropType< TypeTag, Properties::FluidSystem >::IndexTraitsType > wbp_
Definition: BlackoilWellModelGeneric.hpp:518
const WellState< GetPropType< TypeTag, Properties::Scalar >, GetPropType< TypeTag, Properties::FluidSystem >::IndexTraitsType > & wellState() const
Definition: BlackoilWellModelGeneric.hpp:152
void assignMassGasRate(data::Wells &wsrpt, const GetPropType< TypeTag, Properties::Scalar > gasDensity) const
Class for handling the guide rates in the blackoil well model.
Definition: BlackoilWellModelGuideRates.hpp:47
void assignWellGuideRates(data::Wells &wsrpt, const int reportStepIdx) const
Assign well guide rates.
Class for handling the blackoil well network model.
Definition: BlackoilWellModelNetwork.hpp:47
Class for handling the blackoil well model in a NLDD solver.
Definition: BlackoilWellModelNldd.hpp:80
Class for handling the blackoil well model.
Definition: BlackoilWellModel.hpp:98
ReservoirCoupling::Proxy< Scalar > & rescoup()
Get the reservoir coupling proxy.
Definition: BlackoilWellModel.hpp:375
void initializeGroupStructure(const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:294
bool network_needs_more_balancing_force_another_newton_iteration_
Definition: BlackoilWellModel.hpp:522
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:2155
ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster()
Get reference to reservoir coupling master.
Definition: BlackoilWellModel.hpp:401
void endTimeStep()
Definition: BlackoilWellModel.hpp:165
void prepareTimeStep(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1974
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:1288
WellInterfacePtr createWellPointer(const int wellID, const int report_step) const
Definition: BlackoilWellModel_impl.hpp:1084
GuideRateHandler< Scalar, IndexTraits > guide_rate_handler_
Definition: BlackoilWellModel.hpp:518
const EclipseState & eclState() const
Definition: BlackoilWellModel.hpp:529
void prepareWellsBeforeAssembling(const double dt)
Definition: BlackoilWellModel_impl.hpp:1364
const std::vector< WellInterfacePtr > & localNonshutWells() const
Get list of local nonshut wells.
Definition: BlackoilWellModel.hpp:306
void prepareDeserialize(const int report_step)
Definition: BlackoilWellModel.hpp:199
void init()
Definition: BlackoilWellModel_impl.hpp:163
void setNlddAdapter(BlackoilWellModelNldd< TypeTag > *mod)
Definition: BlackoilWellModel.hpp:369
const Simulator & simulator() const
Definition: BlackoilWellModel.hpp:366
std::vector< Scalar > depth_
Definition: BlackoilWellModel.hpp:512
static constexpr EnergyModules energyModuleType_
Definition: BlackoilWellModel.hpp:123
std::size_t global_num_cells_
Definition: BlackoilWellModel.hpp:508
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: BlackoilWellModel.hpp:107
const std::vector< Scalar > & B_avg() const
Definition: BlackoilWellModel.hpp:453
SimulatorReportSingle last_report_
Definition: BlackoilWellModel.hpp:517
static const int solventSaturationIdx
Definition: BlackoilWellModel.hpp:120
Scalar gravity_
Definition: BlackoilWellModel.hpp:511
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:2127
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:394
Dune::FieldVector< Scalar, numEq > VectorBlockType
Definition: BlackoilWellModel.hpp:129
void initFromRestartFile(const RestartValue &restartValues)
Definition: BlackoilWellModel.hpp:189
const ModelParameters param_
Definition: BlackoilWellModel.hpp:507
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:2099
const std::map< std::string, int > & well_domain() const
Definition: BlackoilWellModel.hpp:319
bool updateWellControls(DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:1651
void endIteration()
Definition: BlackoilWellModel.hpp:162
static constexpr bool has_polymer_
Definition: BlackoilWellModel.hpp:122
int reportStepIndex() const
Definition: BlackoilWellModel_impl.hpp:2143
void calculateProductivityIndexValues(DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1916
void extractLegacyDepth_()
Definition: BlackoilWellModel_impl.hpp:2115
void extractLegacyCellPvtRegionIndex_()
Definition: BlackoilWellModel_impl.hpp:2083
void recoverWellSolutionAndUpdateWellStateDomain(const BVector &x, const int domainIdx)
Definition: BlackoilWellModel_impl.hpp:1572
void updateAverageFormationFactor()
Definition: BlackoilWellModel_impl.hpp:2026
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: BlackoilWellModel.hpp:106
bool isReservoirCouplingSlave() const
Check if this process is a reservoir coupling slave.
Definition: BlackoilWellModel.hpp:397
void initializeWellState(const int timeStepIdx)
Definition: BlackoilWellModel_impl.hpp:832
const Grid & grid() const
Definition: BlackoilWellModel.hpp:363
void updatePrimaryVariables()
Definition: BlackoilWellModel_impl.hpp:2074
const ReservoirCoupling::Proxy< Scalar > & rescoup() const
Definition: BlackoilWellModel.hpp:376
void computeWellTemperature()
Definition: BlackoilWellModel_impl.hpp:2177
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:1473
const SimulatorReportSingle & lastReport() const
Definition: BlackoilWellModel_impl.hpp:711
bool updateWellControlsAndNetwork(const bool mandatory_network_balance, const double dt, DeferredLogger &local_deferredLogger)
Definition: BlackoilWellModel_impl.hpp:1229
void endEpisode()
Definition: BlackoilWellModel.hpp:171
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: BlackoilWellModel.hpp:110
void addWellContributions(SparseMatrixAdapter &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1463
const EquilGrid & equilGrid() const
Definition: BlackoilWellModel.hpp:526
void assembleWellEq(const double dt)
Definition: BlackoilWellModel_impl.hpp:1352
WellInterfacePtr createWellForWellTest(const std::string &well_name, const int report_step, DeferredLogger &deferred_logger) const
Definition: BlackoilWellModel_impl.hpp:1135
GetPropType< TypeTag, Properties::Indices > Indices
Definition: BlackoilWellModel.hpp:105
int compressedIndexForInteriorLGR(const std::string &lgr_tag, const Connection &conn) const override
Definition: BlackoilWellModel.hpp:343
void calculateExplicitQuantities() const
Definition: BlackoilWellModel_impl.hpp:1636
std::vector< WellInterfacePtr > well_container_
Definition: BlackoilWellModel.hpp:485
static constexpr std::size_t pressureVarIndex
Definition: BlackoilWellModel.hpp:117
void updateAndCommunicate(const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:1731
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: BlackoilWellModel.hpp:285
bool empty() const
Definition: BlackoilWellModel.hpp:330
void computeTotalRatesForDof(RateVector &rate, unsigned globalIdx) const
Definition: BlackoilWellModel_impl.hpp:795
bool addMatrixContributions() const
Definition: BlackoilWellModel.hpp:332
void beginTimeStep()
Definition: BlackoilWellModel_impl.hpp:326
std::vector< Scalar > B_avg_
Definition: BlackoilWellModel.hpp:524
const std::vector< Scalar > & wellPerfEfficiencyFactors() const
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: BlackoilWellModel.hpp:108
bool updateGroupControls(const Group &group, DeferredLogger &deferred_logger, const int reportStepIdx)
Definition: BlackoilWellModel_impl.hpp:1765
std::map< std::string, std::unique_ptr< AverageRegionalPressureType > > regionalAveragePressureCalculator_
Definition: BlackoilWellModel.hpp:515
void calcInjResvCoeff(const int fipnum, const int pvtreg, std::vector< Scalar > &resv_coeff) const override
Definition: BlackoilWellModel_impl.hpp:2166
static Scalar computeTemperatureWeightFactor(const int perf_index, const int np, const FluidState &fs, const SingleWellState &ws)
Definition: BlackoilWellModel.hpp:461
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:514
BlackoilWellModel(Simulator &simulator)
Definition: BlackoilWellModel_impl.hpp:76
const SparseTable< int > & well_local_cells() const
Definition: BlackoilWellModel.hpp:311
void wellTesting(const int timeStepIdx, const double simulationTime, DeferredLogger &deferred_logger)
Definition: BlackoilWellModel_impl.hpp:634
std::vector< bool > is_cell_perforated_
Definition: BlackoilWellModel.hpp:487
ConvergenceReport getWellConvergence(const std::vector< Scalar > &B_avg, const bool checkWellGroupControlsAndNetwork=false) const
Definition: BlackoilWellModel_impl.hpp:1585
int numStrictIterations() const
Definition: BlackoilWellModel.hpp:335
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:1410
data::WellBlockAveragePressures wellBlockAveragePressures() const
Definition: BlackoilWellModel.hpp:241
void assembleWellEqWithoutIteration(const double dt)
Definition: BlackoilWellModel_impl.hpp:1378
void updateCellRates()
Definition: BlackoilWellModel_impl.hpp:1398
void assemble(const double dt)
Definition: BlackoilWellModel_impl.hpp:1160
auto begin() const
Definition: BlackoilWellModel.hpp:328
std::size_t local_num_cells_
Definition: BlackoilWellModel.hpp:510
auto end() const
Definition: BlackoilWellModel.hpp:329
bool alternative_well_rate_init_
Definition: BlackoilWellModel.hpp:513
void timeStepSucceeded(const double simulationTime, const double dt)
Definition: BlackoilWellModel_impl.hpp:721
std::unique_ptr< WellType > createTypedWellPointer(const int wellID, const int time_step) const
Definition: BlackoilWellModel_impl.hpp:1104
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:1881
Simulator & simulator_
Definition: BlackoilWellModel.hpp:482
void updateGuideRates(const int report_step_idx, const double sim_time)
Update guide rates for all wells and groups.
Definition: BlackoilWellModel.hpp:385
ReservoirCoupling::Proxy< Scalar > rescoup_
Definition: BlackoilWellModel.hpp:519
void createWellContainer(const int report_step) override
Definition: BlackoilWellModel_impl.hpp:875
std::unique_ptr< WellInterface< TypeTag > > WellInterfacePtr
Definition: BlackoilWellModel.hpp:186
ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave()
Get reference to reservoir coupling slave.
Definition: BlackoilWellModel.hpp:407
data::Wells wellData() const
Definition: BlackoilWellModel.hpp:206
const ModelParameters & param() const
Definition: BlackoilWellModel.hpp:456
void updateWellTestState(const double simulationTime, WellTestState &wellTestState)
upate the wellTestState related to economic limits
Definition: BlackoilWellModel_impl.hpp:1810
void addWellPressureEquationsDomain(PressureMatrix &jacobian, const BVector &weights, const bool use_well_weights, const int domainIndex) const
Definition: BlackoilWellModel.hpp:291
void addReservoirSourceTerms(GlobalEqVector &residual, const std::vector< typename SparseMatrixAdapter::MatrixBlock * > &diagMatAddress) const
Definition: BlackoilWellModel_impl.hpp:1495
int compressedIndexForInterior(int cartesian_cell_idx) const override
get compressed index for interior cells (-1, otherwise
Definition: BlackoilWellModel.hpp:338
static constexpr bool has_energy_
Definition: BlackoilWellModel.hpp:124
void recoverWellSolutionAndUpdateWellState(const BVector &x)
Definition: BlackoilWellModel_impl.hpp:1548
void addWellPressureEquationsStruct(PressureMatrix &jacobian) const
Definition: BlackoilWellModel_impl.hpp:1525
void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger &deferred_logger) override
Definition: BlackoilWellModel_impl.hpp:1930
void endReportStep()
Definition: BlackoilWellModel_impl.hpp:694
data::WellBlockAveragePressures computeWellBlockAveragePressures(const Scalar gravity) const
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
Definition: GroupStateHelper.hpp:57
Handles computation and reporting of guide rates for wells and groups.
Definition: GuideRateHandler.hpp:46
Definition: RateConverter.hpp:70
Definition: RegionAverageCalculator.hpp:60
Thin proxy for reservoir coupling master/slave pointers.
Definition: RescoupProxy.hpp:54
Definition: ReservoirCouplingMaster.hpp:38
Definition: ReservoirCouplingSlave.hpp:40
Definition: SingleWellState.hpp:43
PerfData< Scalar > perf_data
Definition: SingleWellState.hpp:143
Definition: BlackoilWellModel.hpp:87
Definition: WellConnectionAuxiliaryModule.hpp:39
Definition: WellContributions.hpp:51
Definition: WellInterface.hpp:77
Definition: WellState.hpp:66
data::Wells report(const int *globalCellIdxMap, const std::function< bool(const int)> &wasDynamicallyClosed) const
ExcEnum
Definition: DeferredLogger.hpp:45
Definition: blackoilbioeffectsmodules.hh:45
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:194
bool matrix_add_well_contributions_
Whether to add influences of wells between cells to the matrix and preconditioner matrix.
Definition: BlackoilModelParameters.hpp:308
bool use_multisegment_well_
Definition: BlackoilModelParameters.hpp:302
int strict_outer_iter_wells_
Newton iteration where wells are stricly convergent.
Definition: BlackoilModelParameters.hpp:258
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34