WellInterface.hpp
Go to the documentation of this file.
1/*
2 Copyright 2017 SINTEF Digital, Mathematics and Cybernetics.
3 Copyright 2017 Statoil ASA.
4 Copyright 2017 IRIS
5 Copyright 2019 Norce
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_WELLINTERFACE_HEADER_INCLUDED
24#define OPM_WELLINTERFACE_HEADER_INCLUDED
25
26// NOTE: GasLiftSingleWell.hpp includes StandardWell.hpp which includes ourself
27// (WellInterface.hpp), so we need to forward declare GasLiftSingleWell
28// for it to be defined in this file. Similar for BlackoilWellModel
29namespace Opm {
30 template<typename TypeTag> class GasLiftSingleWell;
31 template<typename TypeTag> class BlackoilWellModel;
32}
33
34#include <opm/common/OpmLog/OpmLog.hpp>
35#include <opm/common/ErrorMacros.hpp>
36#include <opm/common/Exceptions.hpp>
37
38#include <opm/input/eclipse/Schedule/Well/WellTestState.hpp>
39
41
43
52
54
56
57#include <dune/common/fmatrix.hh>
58#include <dune/istl/bcrsmatrix.hh>
59#include <dune/istl/matrixmatrix.hh>
60
61#include <opm/material/densead/Evaluation.hpp>
62
63#include <cassert>
64#include <vector>
65
66namespace Opm
67{
68
69class WellInjectionProperties;
70class WellProductionProperties;
71
72template<typename TypeTag>
73class WellInterface : public WellInterfaceIndices<GetPropType<TypeTag, Properties::FluidSystem>,
74 GetPropType<TypeTag, Properties::Indices>>
75{
77 GetPropType<TypeTag, Properties::Indices>>;
78public:
80
81 using Grid = GetPropType<TypeTag, Properties::Grid>;
82 using Simulator = GetPropType<TypeTag, Properties::Simulator>;
83 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
84 using Indices = GetPropType<TypeTag, Properties::Indices>;
85 using IntensiveQuantities = GetPropType<TypeTag, Properties::IntensiveQuantities>;
86 using MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>;
87 using SparseMatrixAdapter = GetPropType<TypeTag, Properties::SparseMatrixAdapter>;
88 using RateVector = GetPropType<TypeTag, Properties::RateVector>;
95
96 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
97
98 using VectorBlockType = Dune::FieldVector<Scalar, Indices::numEq>;
99 using MatrixBlockType = Dune::FieldMatrix<Scalar, Indices::numEq, Indices::numEq>;
100 using Eval = typename Base::Eval;
101 using BVector = Dune::BlockVector<VectorBlockType>;
102 using PressureMatrix = Dune::BCRSMatrix<Opm::MatrixBlock<double, 1, 1>>;
103
106
110
111 static constexpr bool has_solvent = getPropValue<TypeTag, Properties::EnableSolvent>();
112 static constexpr bool has_zFraction = getPropValue<TypeTag, Properties::EnableExtbo>();
113 static constexpr bool has_polymer = getPropValue<TypeTag, Properties::EnablePolymer>();
114 static constexpr bool has_energy = getPropValue<TypeTag, Properties::EnableEnergy>();
115 static const bool has_temperature = getPropValue<TypeTag, Properties::EnableTemperature>();
116 // flag for polymer molecular weight related
117 static constexpr bool has_polymermw = getPropValue<TypeTag, Properties::EnablePolymerMW>();
118 static constexpr bool has_foam = getPropValue<TypeTag, Properties::EnableFoam>();
119 static constexpr bool has_brine = getPropValue<TypeTag, Properties::EnableBrine>();
120 static constexpr bool has_watVapor = getPropValue<TypeTag, Properties::EnableVapwat>();
121 static constexpr bool has_disgas_in_water = getPropValue<TypeTag, Properties::EnableDisgasInWater>();
122 static constexpr bool has_saltPrecip = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
123 static constexpr bool has_micp = getPropValue<TypeTag, Properties::EnableMICP>();
124
125 // For the conversion between the surface volume rate and reservoir voidage rate
126 using FluidState = BlackOilFluidState<Eval,
130 Indices::compositionSwitchIdx >= 0,
132 has_brine,
135 Indices::numPhases >;
137 WellInterface(const Well& well,
138 const ParallelWellInfo& pw_info,
139 const int time_step,
140 const ModelParameters& param,
141 const RateConverterType& rate_converter,
142 const int pvtRegionIdx,
143 const int num_components,
144 const int num_phases,
145 const int index_of_well,
146 const std::vector<PerforationData>& perf_data);
147
149 virtual ~WellInterface() = default;
150
151 virtual void init(const PhaseUsage* phase_usage_arg,
152 const std::vector<double>& depth_arg,
153 const double gravity_arg,
154 const int num_cells,
155 const std::vector< Scalar >& B_avg,
156 const bool changed_to_open_this_step);
157
159
160 virtual ConvergenceReport getWellConvergence(const SummaryState& summary_state,
161 const WellState<Scalar>& well_state,
162 const std::vector<double>& B_avg,
163 DeferredLogger& deferred_logger,
164 const bool relax_tolerance) const = 0;
165
166 virtual void solveEqAndUpdateWellState(const SummaryState& summary_state,
167 WellState<Scalar>& well_state,
168 DeferredLogger& deferred_logger) = 0;
169
170 void assembleWellEq(const Simulator& simulator,
171 const double dt,
172 WellState<Scalar>& well_state,
173 const GroupState<Scalar>& group_state,
174 DeferredLogger& deferred_logger);
175
176 void assembleWellEqWithoutIteration(const Simulator& simulator,
177 const double dt,
178 WellState<Scalar>& well_state,
179 const GroupState<Scalar>& group_state,
180 DeferredLogger& deferred_logger);
181
182 // TODO: better name or further refactoring the function to make it more clear
183 void prepareWellBeforeAssembling(const Simulator& simulator,
184 const double dt,
185 WellState<Scalar>& well_state,
186 const GroupState<Scalar>& group_state,
187 DeferredLogger& deferred_logger);
188
189
191 const Simulator& simulator,
192 const double& bhp,
193 std::vector<double>& well_flux,
194 DeferredLogger& deferred_logger
195 ) const = 0;
196
197 virtual std::optional<double> computeBhpAtThpLimitProdWithAlq(
198 const Simulator& simulator,
199 const SummaryState& summary_state,
200 const double alq_value,
201 DeferredLogger& deferred_logger
202 ) const = 0;
203
206 virtual void recoverWellSolutionAndUpdateWellState(const SummaryState& summary_state,
207 const BVector& x,
208 WellState<Scalar>& well_state,
209 DeferredLogger& deferred_logger) = 0;
210
212 virtual void apply(const BVector& x, BVector& Ax) const = 0;
213
215 virtual void apply(BVector& r) const = 0;
216
217 // TODO: before we decide to put more information under mutable, this function is not const
218 virtual void computeWellPotentials(const Simulator& simulator,
219 const WellState<Scalar>& well_state,
220 std::vector<double>& well_potentials,
221 DeferredLogger& deferred_logger) = 0;
222
223 virtual void updateWellStateWithTarget(const Simulator& simulator,
224 const GroupState<Scalar>& group_state,
225 WellState<Scalar>& well_state,
226 DeferredLogger& deferred_logger) const;
227
228 virtual void computeWellRatesWithBhpIterations(const Simulator& simulator,
229 const Scalar& bhp,
230 std::vector<double>& well_flux,
231 DeferredLogger& deferred_logger) const = 0;
232
233 bool updateWellStateWithTHPTargetProd(const Simulator& simulator,
234 WellState<Scalar>& well_state,
235 DeferredLogger& deferred_logger) const;
236
237 enum class IndividualOrGroup { Individual, Group, Both };
238 bool updateWellControl(const Simulator& simulator,
239 const IndividualOrGroup iog,
240 WellState<Scalar>& well_state,
241 const GroupState<Scalar>& group_state,
242 DeferredLogger& deferred_logger) /* const */;
243
245 WellState<Scalar>& well_state,
246 const GroupState<Scalar>& group_state,
247 const Well::InjectionControls& inj_controls,
248 const Well::ProductionControls& prod_controls,
249 const double WQTotal,
250 DeferredLogger& deferred_logger,
251 const bool fixed_control = false,
252 const bool fixed_status = false);
253
254 virtual void updatePrimaryVariables(const SummaryState& summary_state,
255 const WellState<Scalar>& well_state,
256 DeferredLogger& deferred_logger) = 0;
257
258 virtual void calculateExplicitQuantities(const Simulator& simulator,
259 const WellState<Scalar>& well_state,
260 DeferredLogger& deferred_logger) = 0; // should be const?
261
262 virtual void updateProductivityIndex(const Simulator& simulator,
263 const WellProdIndexCalculator& wellPICalc,
264 WellState<Scalar>& well_state,
265 DeferredLogger& deferred_logger) const = 0;
266
267 virtual double connectionDensity(const int globalConnIdx,
268 const int openConnIdx) const = 0;
269
272 {
273 return false;
274 }
275
276 // Add well contributions to matrix
278
280 const BVector& x,
281 const int pressureVarIndex,
282 const bool use_well_weights,
283 const WellState<Scalar>& well_state) const = 0;
284
285 void addCellRates(RateVector& rates, int cellIdx) const;
286
287 Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const;
288
289 // TODO: theoretically, it should be a const function
290 // Simulator is not const is because that assembleWellEq is non-const Simulator
291 void wellTesting(const Simulator& simulator,
292 const double simulation_time,
293 /* const */ WellState<Scalar>& well_state,
294 const GroupState<Scalar>& group_state,
295 WellTestState& welltest_state,
296 DeferredLogger& deferred_logger);
297
298 void checkWellOperability(const Simulator& simulator,
299 const WellState<Scalar>& well_state,
300 DeferredLogger& deferred_logger);
301
303 const Simulator& simulator,
304 const double dt,
305 WellState<Scalar>& well_state,
306 const GroupState<Scalar>& group_state,
307 DeferredLogger& deferred_logger);
308
310 WellState<Scalar>& well_state,
311 DeferredLogger& deferred_logger);
312
313 // check whether the well is operable under the current reservoir condition
314 // mostly related to BHP limit and THP limit
315 void updateWellOperability(const Simulator& simulator,
316 const WellState<Scalar>& well_state,
317 DeferredLogger& deferred_logger);
318
319 bool updateWellOperabilityFromWellEq(const Simulator& simulator,
320 const WellState<Scalar>& well_state,
321 DeferredLogger& deferred_logger);
322
323 // update perforation water throughput based on solved water rate
324 virtual void updateWaterThroughput(const double dt,
325 WellState<Scalar>& well_state) const = 0;
326
329 virtual std::vector<double> computeCurrentWellRates(const Simulator& simulator,
330 DeferredLogger& deferred_logger) const = 0;
331
335 void updateWellStateRates(const Simulator& simulator,
336 WellState<Scalar>& well_state,
337 DeferredLogger& deferred_logger) const;
338
339 void solveWellEquation(const Simulator& simulator,
340 WellState<Scalar>& well_state,
341 const GroupState<Scalar>& group_state,
342 DeferredLogger& deferred_logger);
343
344 const std::vector<RateVector>& connectionRates() const
345 {
346 return connectionRates_;
347 }
348
349 virtual std::vector<double> getPrimaryVars() const
350 {
351 return {};
352 }
353
354 virtual int setPrimaryVars(std::vector<double>::const_iterator)
355 {
356 return 0;
357 }
358
359 std::vector<double> wellIndex(const int perf,
360 const IntensiveQuantities& intQuants,
361 const double trans_mult,
362 const SingleWellState<double>& ws) const;
363
364 void updateConnectionDFactor(const Simulator& simulator,
365 SingleWellState<double>& ws) const;
366
368 SingleWellState<double>& ws) const;
369
370
371protected:
372 // simulation parameters
374 std::vector<RateVector> connectionRates_;
375 std::vector< Scalar > B_avg_;
378
379 double wpolymer() const;
380
381 double wfoam() const;
382
383 double wsalt() const;
384
385 double wmicrobes() const;
386
387 double woxygen() const;
388
389 double wurea() const;
390
391 virtual double getRefDensity() const = 0;
392
393 // Component fractions for each phase for the well
394 const std::vector<double>& compFrac() const;
395
396 std::vector<double> initialWellRateFractions(const Simulator& simulator,
397 const WellState<Scalar>& well_state) const;
398
399 // check whether the well is operable under BHP limit with current reservoir condition
400 virtual void checkOperabilityUnderBHPLimit(const WellState<Scalar>& well_state,
401 const Simulator& simulator,
402 DeferredLogger& deferred_logger) = 0;
403
404 // check whether the well is operable under THP limit with current reservoir condition
405 virtual void checkOperabilityUnderTHPLimit(const Simulator& simulator,
406 const WellState<Scalar>& well_state,
407 DeferredLogger& deferred_logger) = 0;
408
409 virtual void updateIPR(const Simulator& simulator,
410 DeferredLogger& deferred_logger) const = 0;
411
412 virtual void assembleWellEqWithoutIteration(const Simulator& simulator,
413 const double dt,
414 const WellInjectionControls& inj_controls,
415 const WellProductionControls& prod_controls,
416 WellState<Scalar>& well_state,
417 const GroupState<Scalar>& group_state,
418 DeferredLogger& deferred_logger) = 0;
419
420 // iterate well equations with the specified control until converged
421 virtual bool iterateWellEqWithControl(const Simulator& simulator,
422 const double dt,
423 const WellInjectionControls& inj_controls,
424 const WellProductionControls& prod_controls,
425 WellState<Scalar>& well_state,
426 const GroupState<Scalar>& group_state,
427 DeferredLogger& deferred_logger) = 0;
428
429 virtual bool iterateWellEqWithSwitching(const Simulator& simulator,
430 const double dt,
431 const WellInjectionControls& inj_controls,
432 const WellProductionControls& prod_controls,
433 WellState<Scalar>& well_state,
434 const GroupState<Scalar>& group_state,
435 DeferredLogger& deferred_logger,
436 const bool fixed_control = false,
437 const bool fixed_status = false) = 0;
438
439 virtual void updateIPRImplicit(const Simulator& simulator,
440 WellState<Scalar>& well_state,
441 DeferredLogger& deferred_logger) = 0;
442
443 bool iterateWellEquations(const Simulator& simulator,
444 const double dt,
445 WellState<Scalar>& well_state,
446 const GroupState<Scalar>& group_state,
447 DeferredLogger& deferred_logger);
448
449 bool solveWellWithTHPConstraint(const Simulator& simulator,
450 const double dt,
451 const Well::InjectionControls& inj_controls,
452 const Well::ProductionControls& prod_controls,
453 WellState<Scalar>& well_state,
454 const GroupState<Scalar>& group_state,
455 DeferredLogger& deferred_logger);
456
457 std::optional<double> estimateOperableBhp(const Simulator& simulator,
458 const double dt,
459 WellState<Scalar>& well_state,
460 const SummaryState& summary_state,
461 DeferredLogger& deferred_logger);
462
463 bool solveWellWithBhp(const Simulator& simulator,
464 const double dt,
465 const double bhp,
466 WellState<Scalar>& well_state,
467 DeferredLogger& deferred_logger);
468
469 bool solveWellWithZeroRate(const Simulator& simulator,
470 const double dt,
471 WellState<Scalar>& well_state,
472 DeferredLogger& deferred_logger);
473
474 bool solveWellForTesting(const Simulator& simulator,
475 WellState<Scalar>& well_state,
476 const GroupState<Scalar>& group_state,
477 DeferredLogger& deferred_logger);
478
479 Eval getPerfCellPressure(const FluidState& fs) const;
480
481 // get the mobility for specific perforation
482 template<class Value, class Callback>
483 void getMobility(const Simulator& simulator,
484 const int perf,
485 std::vector<Value>& mob,
486 Callback& extendEval,
487 [[maybe_unused]] DeferredLogger& deferred_logger) const;
488
489 void computeConnLevelProdInd(const FluidState& fs,
490 const std::function<double(const double)>& connPICalc,
491 const std::vector<Scalar>& mobility,
492 double* connPI) const;
493
494 void computeConnLevelInjInd(const FluidState& fs,
495 const Phase preferred_phase,
496 const std::function<double(const double)>& connIICalc,
497 const std::vector<Scalar>& mobility,
498 double* connII,
499 DeferredLogger& deferred_logger) const;
500
501 double computeConnectionDFactor(const int perf,
502 const IntensiveQuantities& intQuants,
503 const SingleWellState<double>& ws) const;
504
505};
506
507} // namespace Opm
508
509#ifndef OPM_WELLINTERFACE_IMPL_HEADER_INCLUDED
510#include "WellInterface_impl.hpp"
511#endif
512
513#endif // OPM_WELLINTERFACE_HEADER_INCLUDED
typename BlackoilWellModelGeneric< Scalar >::GLiftProdWells GLiftProdWells
Definition: BlackoilWellModel.hpp:118
typename BlackoilWellModelGeneric< Scalar >::GLiftOptWells GLiftOptWells
Definition: BlackoilWellModel.hpp:117
typename BlackoilWellModelGeneric< Scalar >::GLiftWellStateMap GLiftWellStateMap
Definition: BlackoilWellModel.hpp:120
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
std::set< int > GLiftSyncGroups
Definition: GasLiftSingleWellGeneric.hpp:58
Definition: GasLiftSingleWell.hpp:38
Definition: GroupState.hpp:35
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:184
Definition: RateConverter.hpp:70
Definition: WellInterfaceFluidSystem.hpp:47
static constexpr int Oil
Definition: WellInterfaceFluidSystem.hpp:58
static constexpr int Water
Definition: WellInterfaceFluidSystem.hpp:57
static constexpr int Gas
Definition: WellInterfaceFluidSystem.hpp:59
int pvtRegionIdx() const
Definition: WellInterfaceGeneric.hpp:126
const std::vector< double > & wellIndex() const
Definition: WellInterfaceGeneric.hpp:170
Definition: WellInterfaceIndices.hpp:35
DenseAd::Evaluation< Scalar, Indices::numEq > Eval
Definition: WellInterfaceIndices.hpp:41
Definition: WellInterface.hpp:75
virtual void solveEqAndUpdateWellState(const SummaryState &summary_state, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
Dune::BCRSMatrix< Opm::MatrixBlock< double, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:102
virtual void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator &wellPICalc, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const =0
bool updateWellControl(const Simulator &simulator, const IndividualOrGroup iog, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:194
std::optional< double > estimateOperableBhp(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const SummaryState &summary_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:576
virtual std::vector< double > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const =0
virtual std::vector< double > getPrimaryVars() const
Definition: WellInterface.hpp:349
virtual double connectionDensity(const int globalConnIdx, const int openConnIdx) const =0
const ModelParameters & param_
Definition: WellInterface.hpp:373
bool solveWellWithBhp(const Simulator &simulator, const double dt, const double bhp, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:599
static constexpr bool has_brine
Definition: WellInterface.hpp:119
IndividualOrGroup
Definition: WellInterface.hpp:237
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:82
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1691
void assembleWellEq(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:747
virtual void addWellContributions(SparseMatrixAdapter &) const =0
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:105
const std::vector< RateVector > & connectionRates() const
Definition: WellInterface.hpp:344
BlackOilFluidState< Eval, FluidSystem, has_temperature, has_energy, Indices::compositionSwitchIdx >=0, has_watVapor, has_brine, has_saltPrecip, has_disgas_in_water, Indices::numPhases > FluidState
Definition: WellInterface.hpp:135
virtual void assembleWellEqWithoutIteration(const Simulator &simulator, const double dt, const WellInjectionControls &inj_controls, const WellProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)=0
virtual void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const =0
virtual void checkOperabilityUnderBHPLimit(const WellState< Scalar > &well_state, const Simulator &simulator, DeferredLogger &deferred_logger)=0
std::vector< RateVector > connectionRates_
Definition: WellInterface.hpp:374
bool updateWellStateWithTHPTargetProd(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1766
void addCellRates(RateVector &rates, int cellIdx) const
Definition: WellInterface_impl.hpp:839
virtual void calculateExplicitQuantities(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const
Definition: WellInterface_impl.hpp:855
void computeConnLevelInjInd(const FluidState &fs, const Phase preferred_phase, const std::function< double(const double)> &connIICalc, const std::vector< Scalar > &mobility, double *connII, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1831
virtual void addWellPressureEquations(PressureMatrix &mat, const BVector &x, const int pressureVarIndex, const bool use_well_weights, const WellState< Scalar > &well_state) const =0
virtual void apply(const BVector &x, BVector &Ax) const =0
Ax = Ax - C D^-1 B x.
static const bool has_temperature
Definition: WellInterface.hpp:115
bool iterateWellEquations(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:453
void computeConnLevelProdInd(const FluidState &fs, const std::function< double(const double)> &connPICalc, const std::vector< Scalar > &mobility, double *connPI) const
Definition: WellInterface_impl.hpp:1796
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:85
void updateConnectionDFactor(const Simulator &simulator, SingleWellState< double > &ws) const
Definition: WellInterface_impl.hpp:1585
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:96
void checkWellOperability(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:874
virtual double getRefDensity() const =0
void updateWellOperability(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:978
virtual void computeWellRatesWithBhpIterations(const Simulator &simulator, const Scalar &bhp, std::vector< double > &well_flux, DeferredLogger &deferred_logger) const =0
void updateWellStateRates(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1458
static constexpr bool has_watVapor
Definition: WellInterface.hpp:120
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: WellInterface.hpp:86
virtual void updateWellStateWithTarget(const Simulator &simulator, const GroupState< Scalar > &group_state, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1032
double wfoam() const
Definition: WellInterface_impl.hpp:126
double wsalt() const
Definition: WellInterface_impl.hpp:140
Eval getPerfCellPressure(const FluidState &fs) const
Definition: WellInterface_impl.hpp:1676
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:101
bool gliftBeginTimeStepWellTestIterateWellEquations(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:899
virtual ConvergenceReport getWellConvergence(const SummaryState &summary_state, const WellState< Scalar > &well_state, const std::vector< double > &B_avg, DeferredLogger &deferred_logger, const bool relax_tolerance) const =0
bool changed_to_stopped_this_step_
Definition: WellInterface.hpp:376
static constexpr bool has_polymermw
Definition: WellInterface.hpp:117
static constexpr bool has_polymer
Definition: WellInterface.hpp:113
virtual int setPrimaryVars(std::vector< double >::const_iterator)
Definition: WellInterface.hpp:354
virtual void apply(BVector &r) const =0
r = r - C D^-1 Rw
GetPropType< TypeTag, Properties::Grid > Grid
Definition: WellInterface.hpp:81
void updateConnectionTransmissibilityFactor(const Simulator &simulator, SingleWellState< double > &ws) const
Definition: WellInterface_impl.hpp:1647
virtual ~WellInterface()=default
Virtual destructor.
virtual void computeWellPotentials(const Simulator &simulator, const WellState< Scalar > &well_state, std::vector< double > &well_potentials, DeferredLogger &deferred_logger)=0
WellInterface(const Well &well, const ParallelWellInfo &pw_info, const int time_step, const ModelParameters &param, const RateConverterType &rate_converter, const int pvtRegionIdx, const int num_components, const int num_phases, const int index_of_well, const std::vector< PerforationData > &perf_data)
Constructor.
Definition: WellInterface_impl.hpp:54
typename BlackoilWellModel< TypeTag >::GLiftProdWells GLiftProdWells
Definition: WellInterface.hpp:91
virtual void updateIPRImplicit(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:83
static constexpr bool has_solvent
Definition: WellInterface.hpp:111
static constexpr bool has_disgas_in_water
Definition: WellInterface.hpp:121
void prepareWellBeforeAssembling(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:783
static constexpr bool has_saltPrecip
Definition: WellInterface.hpp:122
const std::vector< double > & compFrac() const
static constexpr bool has_foam
Definition: WellInterface.hpp:118
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: WellInterface.hpp:88
virtual void updateWaterThroughput(const double dt, WellState< Scalar > &well_state) const =0
void solveWellEquation(const Simulator &simulator, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:690
virtual bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const WellInjectionControls &inj_controls, const WellProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)=0
static constexpr bool has_micp
Definition: WellInterface.hpp:123
virtual bool jacobianContainsWellContributions() const
Wether the Jacobian will also have well contributions in it.
Definition: WellInterface.hpp:271
double woxygen() const
Definition: WellInterface_impl.hpp:164
bool updateWellOperabilityFromWellEq(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1014
double wurea() const
Definition: WellInterface_impl.hpp:182
virtual bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const WellInjectionControls &inj_controls, const WellProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger, const bool fixed_control=false, const bool fixed_status=false)=0
virtual std::optional< double > computeBhpAtThpLimitProdWithAlq(const Simulator &simulator, const SummaryState &summary_state, const double alq_value, DeferredLogger &deferred_logger) const =0
Dune::FieldMatrix< Scalar, Indices::numEq, Indices::numEq > MatrixBlockType
Definition: WellInterface.hpp:99
void gliftBeginTimeStepWellTestUpdateALQ(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:934
std::vector< double > initialWellRateFractions(const Simulator &simulator, const WellState< Scalar > &well_state) const
Definition: WellInterface_impl.hpp:1411
virtual void updatePrimaryVariables(const SummaryState &summary_state, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
typename Base::Eval Eval
Definition: WellInterface.hpp:100
double computeConnectionDFactor(const int perf, const IntensiveQuantities &intQuants, const SingleWellState< double > &ws) const
Definition: WellInterface_impl.hpp:1605
virtual void init(const PhaseUsage *phase_usage_arg, const std::vector< double > &depth_arg, const double gravity_arg, const int num_cells, const std::vector< Scalar > &B_avg, const bool changed_to_open_this_step)
Definition: WellInterface_impl.hpp:91
void wellTesting(const Simulator &simulator, const double simulation_time, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, WellTestState &welltest_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:355
static constexpr bool has_energy
Definition: WellInterface.hpp:114
virtual void computeWellRatesWithBhp(const Simulator &simulator, const double &bhp, std::vector< double > &well_flux, DeferredLogger &deferred_logger) const =0
double wpolymer() const
Definition: WellInterface_impl.hpp:110
bool solveWellWithZeroRate(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:635
bool solveWellForTesting(const Simulator &simulator, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:655
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:84
bool solveWellWithTHPConstraint(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:486
typename GasLiftSingleWellGeneric::GLiftSyncGroups GLiftSyncGroups
Definition: WellInterface.hpp:94
bool updateWellControlAndStatusLocalIteration(const Simulator &simulator, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const double WQTotal, DeferredLogger &deferred_logger, const bool fixed_control=false, const bool fixed_status=false)
Definition: WellInterface_impl.hpp:271
bool thp_update_iterations
Definition: WellInterface.hpp:377
Dune::FieldVector< Scalar, Indices::numEq > VectorBlockType
Definition: WellInterface.hpp:98
double wmicrobes() const
Definition: WellInterface_impl.hpp:152
void assembleWellEqWithoutIteration(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:764
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:87
typename BlackoilWellModel< TypeTag >::GLiftOptWells GLiftOptWells
Definition: WellInterface.hpp:90
virtual void recoverWellSolutionAndUpdateWellState(const SummaryState &summary_state, const BVector &x, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
virtual void checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
static constexpr bool has_zFraction
Definition: WellInterface.hpp:112
std::vector< Scalar > B_avg_
Definition: WellInterface.hpp:375
virtual void initPrimaryVariablesEvaluation()=0
typename BlackoilWellModel< TypeTag >::GLiftWellStateMap GLiftWellStateMap
Definition: WellInterface.hpp:93
Definition: WellProdIndexCalculator.hpp:36
Definition: WellState.hpp:62
VFPEvaluation bhp(const VFPProdTable &table, const double aqua, const double liquid, const double vapour, const double thp, const double alq, const double explicit_wfr, const double explicit_gfr, const bool use_vfpexplicit)
Definition: BlackoilPhases.hpp:27
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:484
Definition: BlackoilPhases.hpp:46