opm-simulators
Loading...
Searching...
No Matches
BlackoilWellModel.hpp
1/*
2 Copyright 2016 SINTEF ICT, Applied Mathematics.
3 Copyright 2016 - 2017 Statoil ASA.
4 Copyright 2017 Dr. Blatt - HPC-Simulation-Software & Services
5 Copyright 2016 - 2018 IRIS AS
6
7 This file is part of the Open Porous Media project (OPM).
8
9 OPM is free software: you can redistribute it and/or modify
10 it under the terms of the GNU General Public License as published by
11 the Free Software Foundation, either version 3 of the License, or
12 (at your option) any later version.
13
14 OPM is distributed in the hope that it will be useful,
15 but WITHOUT ANY WARRANTY; without even the implied warranty of
16 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 GNU General Public License for more details.
18
19 You should have received a copy of the GNU General Public License
20 along with OPM. If not, see <http://www.gnu.org/licenses/>.
21*/
22
23#ifndef OPM_BLACKOILWELLMODEL_HEADER_INCLUDED
24#define OPM_BLACKOILWELLMODEL_HEADER_INCLUDED
25
26#include <dune/common/fmatrix.hh>
27
28#include <dune/istl/bcrsmatrix.hh>
29#include <dune/istl/matrixmatrix.hh>
30
31#include <opm/common/OpmLog/OpmLog.hpp>
32
33#include <opm/input/eclipse/Schedule/Group/Group.hpp>
34#include <opm/input/eclipse/Schedule/Group/GuideRate.hpp>
35#include <opm/input/eclipse/Schedule/Schedule.hpp>
36#include <opm/input/eclipse/Schedule/Well/WellTestState.hpp>
37
38#include <opm/material/densead/Math.hpp>
39
40#include <opm/simulators/flow/countGlobalCells.hpp>
42
43#include <opm/simulators/linalg/matrixblock.hh>
44
45#include <opm/simulators/timestepping/SimulatorReport.hpp>
46#include <opm/simulators/timestepping/gatherConvergenceReport.hpp>
47
48#include <opm/simulators/utils/DeferredLogger.hpp>
49
50#include <opm/simulators/wells/BlackoilWellModelGasLift.hpp>
51#include <opm/simulators/wells/BlackoilWellModelGeneric.hpp>
52#include <opm/simulators/wells/BlackoilWellModelGuideRates.hpp>
53#include <opm/simulators/wells/BlackoilWellModelNetwork.hpp>
54#include <opm/simulators/wells/GasLiftGroupInfo.hpp>
55#include <opm/simulators/wells/GasLiftSingleWell.hpp>
56#include <opm/simulators/wells/GasLiftSingleWellGeneric.hpp>
57#include <opm/simulators/wells/GasLiftWellState.hpp>
58#include <opm/simulators/wells/GroupStateHelper.hpp>
59#include <opm/simulators/wells/GuideRateHandler.hpp>
60#include <opm/simulators/wells/MultisegmentWell.hpp>
61#include <opm/simulators/wells/ParallelWBPCalculation.hpp>
62#include <opm/simulators/wells/ParallelWellInfo.hpp>
63#include <opm/simulators/wells/PerforationData.hpp>
66#include <opm/simulators/wells/StandardWell.hpp>
67#include <opm/simulators/wells/VFPInjProperties.hpp>
68#include <opm/simulators/wells/VFPProdProperties.hpp>
69#include <opm/simulators/wells/WGState.hpp>
70#include <opm/simulators/wells/WellConnectionAuxiliaryModule.hpp>
71#include <opm/simulators/wells/WellInterface.hpp>
72#include <opm/simulators/wells/WellProdIndexCalculator.hpp>
73#include <opm/simulators/wells/WellState.hpp>
74#include <opm/simulators/wells/rescoup/RescoupProxy.hpp>
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 ---------
111 using ModelParameters = BlackoilModelParameters<Scalar>;
112
113 using WellConnectionModule = WellConnectionAuxiliaryModule<TypeTag, BlackoilWellModel<TypeTag>>;
114 using IndexTraits = typename FluidSystem::IndexTraitsType;
115 using GroupStateHelperType = GroupStateHelper<Scalar, IndexTraits>;
116
117 constexpr static std::size_t pressureVarIndex = GetPropType<TypeTag, Properties::Indices>::pressureSwitchIdx;
118
119 static const int numEq = Indices::numEq;
120 static const int solventSaturationIdx = Indices::solventSaturationIdx;
121 static constexpr bool has_solvent_ = getPropValue<TypeTag, Properties::EnableSolvent>();
122 static constexpr bool has_polymer_ = getPropValue<TypeTag, Properties::EnablePolymer>();
123 static constexpr EnergyModules energyModuleType_ = getPropValue<TypeTag, Properties::EnergyModuleType>();
124 static constexpr bool has_energy_ = (energyModuleType_ == EnergyModules::FullyImplicitThermal);
125 static constexpr bool has_micp_ = Indices::enableMICP;
126 static constexpr bool has_geochem_ = getPropValue<TypeTag, Properties::EnableGeochemistry>();
127
128 // TODO: where we should put these types, WellInterface or Well Model?
129 // or there is some other strategy, like TypeTag
130 using VectorBlockType = Dune::FieldVector<Scalar, numEq>;
131 using BVector = Dune::BlockVector<VectorBlockType>;
132
133 using PolymerModule = BlackOilPolymerModule<TypeTag>;
134 using BioeffectsModule = BlackOilBioeffectsModule<TypeTag>;
135
136 // For the conversion between the surface volume rate and reservoir voidage rate
137 using RateConverterType = RateConverter::
138 SurfaceToReservoirVoidage<FluidSystem, std::vector<int> >;
139
140 // For computing average pressured used by gpmaint
141 using AverageRegionalPressureType = RegionAverageCalculator::
142 AverageRegionalPressure<FluidSystem, std::vector<int> >;
143
144 explicit BlackoilWellModel(Simulator& simulator);
145
146 void init();
147 void initWellContainer(const int reportStepIdx) override;
148
149 void beginEpisode()
150 {
151 OPM_TIMEBLOCK(beginEpsiode);
152 beginReportStep(simulator_.episodeIndex());
153 }
154
155 void beginTimeStep();
156
157 void beginIteration()
158 {
159 OPM_TIMEBLOCK(beginIteration);
160 assemble(simulator_.timeStepSize());
161 }
162
163 void endIteration()
164 { }
165
166 void endTimeStep()
167 {
168 OPM_TIMEBLOCK(endTimeStep);
169 timeStepSucceeded(simulator_.time(), simulator_.timeStepSize());
170 }
171
172 void endEpisode()
173 {
174 endReportStep();
175 }
176
177 void computeTotalRatesForDof(RateVector& rate,
178 unsigned globalIdx) const;
179
180 template <class Context>
181 void computeTotalRatesForDof(RateVector& rate,
182 const Context& context,
183 unsigned spaceIdx,
184 unsigned timeIdx) const;
185
186
187 using WellInterfacePtr = std::unique_ptr<WellInterface<TypeTag>>;
188
189 using BlackoilWellModelGeneric<Scalar, IndexTraits>::initFromRestartFile;
190 void initFromRestartFile(const RestartValue& restartValues)
191 {
192 initFromRestartFile(restartValues,
193 this->simulator_.vanguard().transferWTestState(),
194 grid().size(0),
195 param_.use_multisegment_well_,
196 this->simulator_.vanguard().enableDistributedWells());
197 }
198
199 using BlackoilWellModelGeneric<Scalar, IndexTraits>::prepareDeserialize;
200 void prepareDeserialize(const int report_step)
201 {
202 prepareDeserialize(report_step, grid().size(0),
203 param_.use_multisegment_well_,
204 this->simulator_.vanguard().enableDistributedWells());
205 }
206
207 data::Wells wellData() const
208 {
209 auto wsrpt = this->wellState()
210 .report(this->simulator_.vanguard().globalCell().data(),
211 [this](const int well_index)
212 { return this->wasDynamicallyShutThisTimeStep(well_index); },
213 this->rsConstInfo());
214
216 .assignWellGuideRates(wsrpt, this->reportStepIndex());
217
218 this->assignWellTracerRates(wsrpt);
219
220 if constexpr (has_geochem_) {
221 this->assignWellSpeciesRates(wsrpt);
222 }
223
224 if (const auto& rspec = eclState().runspec();
225 rspec.co2Storage() || rspec.h2Storage())
226 {
227 // The gas reference density (surface condition) is the
228 // same for all PVT regions in CO2STORE/H2STORE runs so,
229 // for simplicity, we use region zero (0) here.
230
231 this->assignMassGasRate(wsrpt, FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, 0));
232 }
233
234 this->assignWellTargets(wsrpt);
235
236 this->assignDynamicWellStatus(wsrpt);
237
238 // Assigning (a subset of the) property values in shut
239 // connections should be the last step of wellData().
240 this->assignShutConnections(wsrpt, this->reportStepIndex());
241
242 return wsrpt;
243 }
244
245 data::WellBlockAveragePressures wellBlockAveragePressures() const
246 {
247 return this->wbp_.computeWellBlockAveragePressures(this->gravity_);
248 }
249
250#if COMPILE_GPU_BRIDGE
251 // accumulate the contributions of all Wells in the WellContributions object
252 void getWellContributions(WellContributions<Scalar>& x) const;
253#endif
254
255 // Check if well equations is converged.
256 ConvergenceReport getWellConvergence(const std::vector<Scalar>& B_avg, const bool checkWellGroupControlsAndNetwork = false) const;
257
258 const SimulatorReportSingle& lastReport() const;
259
260 void addWellContributions(SparseMatrixAdapter& jacobian) const;
261
262 // add source from wells to the reservoir matrix
263 void addReservoirSourceTerms(GlobalEqVector& residual,
264 const std::vector<typename SparseMatrixAdapter::MatrixBlock*>& diagMatAddress) const;
265
266 // called at the beginning of a report step
267 void beginReportStep(const int time_step);
268
272 void calculateExplicitQuantities() const;
273
276 void prepareTimeStep(DeferredLogger& deferred_logger);
277
278 bool
279 updateWellControls(DeferredLogger& deferred_logger);
280
281 void updateAndCommunicate(const int reportStepIdx);
282
283 bool updateGroupControls(const Group& group,
284 DeferredLogger& deferred_logger,
285 const int reportStepIdx);
286
287 const WellInterface<TypeTag>& getWell(const std::string& well_name) const;
288
289 using PressureMatrix = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, 1, 1>>;
290
291 void addWellPressureEquations(PressureMatrix& jacobian,
292 const BVector& weights,
293 const bool use_well_weights) const;
294 void addWellPressureEquationsStruct(PressureMatrix& jacobian) const;
295 void addWellPressureEquationsDomain(PressureMatrix& jacobian,
296 const BVector& weights,
297 const bool use_well_weights,
298 const int domainIndex) const
299 {
300 if (!nldd_) {
301 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
302 }
303 return nldd_->addWellPressureEquations(jacobian,
304 weights,
305 use_well_weights,
306 domainIndex);
307 }
308
310 const std::vector<WellInterfacePtr>& localNonshutWells() const
311 {
312 return well_container_;
313 }
314
315 const SparseTable<int>& well_local_cells() const
316 {
317 if (!nldd_) {
318 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
319 }
320 return nldd_->well_local_cells();
321 }
322
323 const std::map<std::string, int>& well_domain() const
324 {
325 if (!nldd_) {
326 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
327 }
328
329 return nldd_->well_domain();
330 }
331
332 auto begin() const { return well_container_.begin(); }
333 auto end() const { return well_container_.end(); }
334 bool empty() const { return well_container_.empty(); }
335
336 bool addMatrixContributions() const
337 { return param_.matrix_add_well_contributions_; }
338
339 int numStrictIterations() const
340 { return param_.strict_outer_iter_wells_; }
341
342 int compressedIndexForInterior(int cartesian_cell_idx) const override
343 {
344 return simulator_.vanguard().compressedIndexForInterior(cartesian_cell_idx);
345 }
346
347 int compressedIndexForInteriorLGR(const std::string& lgr_tag, const Connection& conn) const override
348 {
349 return simulator_.vanguard().compressedIndexForInteriorLGR(lgr_tag, conn);
350 }
351
352 // using the solution x to recover the solution xw for wells and applying
353 // xw to update Well State
354 void recoverWellSolutionAndUpdateWellState(const BVector& x);
355
356 // using the solution x to recover the solution xw for wells and applying
357 // xw to update Well State
358 void recoverWellSolutionAndUpdateWellStateDomain(const BVector& x,
359 const int domainIdx);
360 // Update cellRates_ with contributions from all wells
361 void updateCellRates();
362
363 // Update cellRates_ with contributions from wells in a specific domain
364 void updateCellRatesForDomain(int domainIndex,
365 const std::map<std::string, int>& well_domain_map);
366
367 const Grid& grid() const
368 { return simulator_.vanguard().grid(); }
369
370 const Simulator& simulator() const
371 { return simulator_; }
372
373 void setNlddAdapter(BlackoilWellModelNldd<TypeTag>* mod)
374 { nldd_ = mod; }
375
376 // === Reservoir Coupling ===
377
380 const ReservoirCoupling::Proxy<Scalar>& rescoup() const { return rescoup_; }
381
389 void updateGuideRates(const int report_step_idx,
390 const double sim_time)
391 {
392 this->guide_rate_handler_.updateGuideRates(
393 report_step_idx, sim_time, this->wellState(), this->groupState()
394 );
395 }
396
398 bool isReservoirCouplingMaster() const { return rescoup_.isMaster(); }
399
401 bool isReservoirCouplingSlave() const { return rescoup_.isSlave(); }
402
406 return rescoup_.master();
407 }
408
412 return rescoup_.slave();
413 }
414
415#ifdef RESERVOIR_COUPLING_ENABLED
416 void setReservoirCouplingMaster(ReservoirCouplingMaster<Scalar>* master)
417 {
418 rescoup_.setMaster(master);
419 this->guide_rate_handler_.setReservoirCouplingMaster(master);
420 this->groupStateHelper().setReservoirCouplingMaster(master);
421 }
422 void setReservoirCouplingSlave(ReservoirCouplingSlave<Scalar>* slave)
423 {
424 rescoup_.setSlave(slave);
425 this->guide_rate_handler_.setReservoirCouplingSlave(slave);
426 this->groupStateHelper().setReservoirCouplingSlave(slave);
427 }
428
430 void sendSlaveGroupDataToMaster();
431
433 void receiveSlaveGroupData();
434
435 void receiveGroupConstraintsFromMaster();
436 void sendMasterGroupConstraintsToSlaves();
437 void rescoupSyncSummaryData();
438
447 std::optional<ReservoirCoupling::ScopedLoggerGuard>
448 setupRescoupScopedLogger(DeferredLogger& local_logger);
449#endif
450
451 bool updateWellControlsAndNetwork(const bool mandatory_network_balance,
452 const double dt,
453 DeferredLogger& local_deferredLogger);
454
455 // TODO: finding a better naming
456 void assembleWellEqWithoutIteration(const double dt);
457
458 const std::vector<Scalar>& B_avg() const
459 { return B_avg_; }
460
461 const ModelParameters& param() const
462 { return param_; }
463
464
465 template<class FluidState, class SingleWellState>
466 static Scalar computeTemperatureWeightFactor(const int perf_index, const int np, const FluidState& fs, const SingleWellState& ws)
467 {
468 const auto& perf_phase_rate = ws.perf_data.phase_rates;
469 // we only have one temperature pr cell any phaseIdx will do
470 Scalar cellTemperatures = fs.temperature(/*phaseIdx*/0).value();
471 Scalar weight_factor = 0.0;
472 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
473 if (!FluidSystem::phaseIsActive(phaseIdx)) {
474 continue;
475 }
476 Scalar cellInternalEnergy = fs.enthalpy(phaseIdx).value() -
477 fs.pressure(phaseIdx).value() / fs.density(phaseIdx).value();
478 Scalar cellBinv = fs.invB(phaseIdx).value();
479 Scalar cellDensity = fs.density(phaseIdx).value();
480 Scalar perfPhaseRate = perf_phase_rate[perf_index*np + phaseIdx];
481 weight_factor += cellDensity * (perfPhaseRate / cellBinv) * (cellInternalEnergy / cellTemperatures);
482 }
483 return (std::abs(weight_factor) + 1e-13);
484 }
485
486 protected:
487 Simulator& simulator_;
488
489 // a vector of all the wells.
490 std::vector<WellInterfacePtr> well_container_{};
491
492 std::vector<bool> is_cell_perforated_{};
493
494 void initializeWellState(const int timeStepIdx);
495
496 // create the well container
497 void createWellContainer(const int report_step) override;
498
499 WellInterfacePtr
500 createWellPointer(const int wellID,
501 const int report_step) const;
502
503 template <typename WellType>
504 std::unique_ptr<WellType>
505 createTypedWellPointer(const int wellID,
506 const int time_step) const;
507
508 WellInterfacePtr createWellForWellTest(const std::string& well_name,
509 const int report_step,
510 DeferredLogger& deferred_logger) const;
511
512 const ModelParameters param_;
513 std::size_t global_num_cells_{};
514 // the number of the cells in the local grid
515 std::size_t local_num_cells_{};
516 Scalar gravity_{};
517 std::vector<Scalar> depth_{};
518 bool alternative_well_rate_init_{};
519 std::unique_ptr<RateConverterType> rateConverter_{};
520 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>> regionalAveragePressureCalculator_{};
521
522 SimulatorReportSingle last_report_{};
523 GuideRateHandler<Scalar, IndexTraits> guide_rate_handler_{};
524 ReservoirCoupling::Proxy<Scalar> rescoup_{};
525
526 // A flag to tell the convergence report whether we need to take another newton step
527 bool network_needs_more_balancing_force_another_newton_iteration_{false};
528
529 std::vector<Scalar> B_avg_{};
530
531 const EquilGrid& equilGrid() const
532 { return simulator_.vanguard().equilGrid(); }
533
534 const EclipseState& eclState() const
535 { return simulator_.vanguard().eclState(); }
536
537 // compute the well fluxes and assemble them in to the reservoir equations as source terms
538 // and in the well equations.
539 void assemble(const double dt);
540
541 // well controls and network pressures affect each other and are solved in an iterative manner.
542 // the function handles one iteration of updating well controls and network pressures.
543 // it is possible to decouple the update of well controls and network pressures further.
544 // the returned two booleans are {continue_due_to_network, well_group_control_changed}, respectively
545 std::tuple<bool, bool, Scalar> updateWellControlsAndNetworkIteration(const bool mandatory_network_balance,
546 const bool relax_network_tolerance,
547 const bool optimize_gas_lift,
548 const double dt,
549 DeferredLogger& local_deferredLogger);
550
559 void initializeLocalWellStructure(const int reportStepIdx,
560 const bool enableWellPIScaling);
561
565 void initializeGroupStructure(const int reportStepIdx);
566
567 // called at the end of a time step
568 void timeStepSucceeded(const double simulationTime, const double dt);
569
570 // called at the end of a report step
571 void endReportStep();
572
573 // setting the well_solutions_ based on well_state.
574 void updatePrimaryVariables();
575
576 void updateAverageFormationFactor();
577
578 void computePotentials(const std::size_t widx,
579 const WellState<Scalar, IndexTraits>& well_state_copy,
580 std::string& exc_msg,
581 ExceptionType::ExcEnum& exc_type) override;
582
583 const std::vector<Scalar>& wellPerfEfficiencyFactors() const;
584
585 void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger& deferred_logger) override;
586 void calculateProductivityIndexValues(DeferredLogger& deferred_logger) override;
587 void calculateProductivityIndexValues(const WellInterface<TypeTag>* wellPtr,
588 DeferredLogger& deferred_logger);
589
590 // The number of conservation quantities.
591 int numConservationQuantities() const;
592
593 int reportStepIndex() const;
594
595 void assembleWellEq(const double dt);
596
597 void prepareWellsBeforeAssembling(const double dt);
598
599 void extractLegacyCellPvtRegionIndex_();
600
601 void extractLegacyDepth_();
602
604 void updateWellTestState(const double simulationTime, WellTestState& wellTestState);
605
606 void wellTesting(const int timeStepIdx, const double simulationTime, DeferredLogger& deferred_logger);
607
608 void calcResvCoeff(const int fipnum,
609 const int pvtreg,
610 const std::vector<Scalar>& production_rates,
611 std::vector<Scalar>& resv_coeff) const override;
612
613 void calcInjResvCoeff(const int fipnum,
614 const int pvtreg,
615 std::vector<Scalar>& resv_coeff) const override;
616
617 void computeWellTemperature();
618
619 private:
620 BlackoilWellModelGasLift<TypeTag> gaslift_;
621 BlackoilWellModelNetwork<TypeTag> network_;
622 BlackoilWellModelNldd<TypeTag>* nldd_ = nullptr;
623
624 // These members are used to avoid reallocation in specific functions
625 // instead of using local variables.
626 // Their state is not relevant between function calls, so they can
627 // (and must) be mutable, as the functions using them are const.
628 mutable BVector x_local_;
629
630 // Store cell rates after assembling to avoid iterating all wells and connections for every element
631 std::map<int, RateVector> cellRates_;
632
633 void assignWellTracerRates(data::Wells& wsrpt) const;
634 void assignWellSpeciesRates(data::Wells& wsrpt) const;
635
636 [[nodiscard]] auto rsConstInfo() const
637 -> typename WellState<Scalar,IndexTraits>::RsConstInfo;
638 };
639
640} // namespace Opm
641
642#include "BlackoilWellModel_impl.hpp"
643
644#endif // OPM_BLACKOILWELLMODEL_HEADER_INCLUDED
Helper class for grid instantiation of ECL file-format using problems.
Facility for converting component rates at surface conditions to phase (voidage) rates at reservoir c...
Facility for converting component rates at surface conditions to phase (voidage) rates at reservoir c...
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
void assignShutConnections(data::Wells &wsrpt, const int reportStepIndex) const
Definition BlackoilWellModelGeneric.cpp:1102
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.
Definition BlackoilWellModelGuideRates.cpp:478
Class for handling the blackoil well model in a NLDD solver.
Definition BlackoilWellModelNldd.hpp:80
ReservoirCoupling::Proxy< Scalar > & rescoup()
Get the reservoir coupling proxy.
Definition BlackoilWellModel.hpp:379
void initializeGroupStructure(const int reportStepIdx)
Initialize group control modes/constraints and group solution state.
Definition BlackoilWellModel_impl.hpp:297
ReservoirCouplingMaster< Scalar > & reservoirCouplingMaster()
Get reference to reservoir coupling master.
Definition BlackoilWellModel.hpp:405
void prepareTimeStep(DeferredLogger &deferred_logger)
One-time initialization at the start of each timestep.
Definition BlackoilWellModel_impl.hpp:2020
const std::vector< WellInterfacePtr > & localNonshutWells() const
Get list of local nonshut wells.
Definition BlackoilWellModel.hpp:310
bool isReservoirCouplingMaster() const
Check if this process is a reservoir coupling master.
Definition BlackoilWellModel.hpp:398
bool isReservoirCouplingSlave() const
Check if this process is a reservoir coupling slave.
Definition BlackoilWellModel.hpp:401
void calculateExplicitQuantities() const
Calculating the explicit quantities used in the well calculation.
Definition BlackoilWellModel_impl.hpp:1682
void initializeLocalWellStructure(const int reportStepIdx, const bool enableWellPIScaling)
Update rank's notion of intersecting wells and their associate solution variables.
Definition BlackoilWellModel_impl.hpp:251
void updateGuideRates(const int report_step_idx, const double sim_time)
Update guide rates for all wells and groups.
Definition BlackoilWellModel.hpp:389
ReservoirCouplingSlave< Scalar > & reservoirCouplingSlave()
Get reference to reservoir coupling slave.
Definition BlackoilWellModel.hpp:411
void updateWellTestState(const double simulationTime, WellTestState &wellTestState)
upate the wellTestState related to economic limits
Definition BlackoilWellModel_impl.hpp:1856
int compressedIndexForInterior(int cartesian_cell_idx) const override
get compressed index for interior cells (-1, otherwise
Definition BlackoilWellModel.hpp:342
Represents the convergence status of the whole simulator, to make it possible to query and store the ...
Definition ConvergenceReport.hpp:38
Definition DeferredLogger.hpp:57
Definition GroupStateHelper.hpp:57
Definition ReservoirCouplingMaster.hpp:38
Definition ReservoirCouplingSlave.hpp:40
Thin proxy for reservoir coupling master/slave pointers.
Definition RescoupProxy.hpp:54
Definition BlackoilWellModel.hpp:87
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:51
Definition WellInterface.hpp:77
This file contains a set of helper functions used by VFPProd / VFPInj.
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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
Solver parameters for the BlackoilModel.
Definition BlackoilModelParameters.hpp:194
A struct for returning timing data from a simulator to its caller.
Definition SimulatorReport.hpp:34