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 Scalar = GetPropType<TypeTag, Properties::Scalar>;
84 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
85 using Indices = GetPropType<TypeTag, Properties::Indices>;
86 using IntensiveQuantities = GetPropType<TypeTag, Properties::IntensiveQuantities>;
87 using MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>;
88 using SparseMatrixAdapter = GetPropType<TypeTag, Properties::SparseMatrixAdapter>;
89 using RateVector = GetPropType<TypeTag, Properties::RateVector>;
96
97 using VectorBlockType = Dune::FieldVector<Scalar, Indices::numEq>;
98 using MatrixBlockType = Dune::FieldMatrix<Scalar, Indices::numEq, Indices::numEq>;
99 using Eval = typename Base::Eval;
100 using BVector = Dune::BlockVector<VectorBlockType>;
101 using PressureMatrix = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, 1, 1>>;
102
105
109
110 static constexpr bool has_solvent = getPropValue<TypeTag, Properties::EnableSolvent>();
111 static constexpr bool has_zFraction = getPropValue<TypeTag, Properties::EnableExtbo>();
112 static constexpr bool has_polymer = getPropValue<TypeTag, Properties::EnablePolymer>();
113 static constexpr bool has_energy = getPropValue<TypeTag, Properties::EnableEnergy>();
114 static const bool has_temperature = getPropValue<TypeTag, Properties::EnableTemperature>();
115 // flag for polymer molecular weight related
116 static constexpr bool has_polymermw = getPropValue<TypeTag, Properties::EnablePolymerMW>();
117 static constexpr bool has_foam = getPropValue<TypeTag, Properties::EnableFoam>();
118 static constexpr bool has_brine = getPropValue<TypeTag, Properties::EnableBrine>();
119 static constexpr bool has_watVapor = getPropValue<TypeTag, Properties::EnableVapwat>();
120 static constexpr bool has_disgas_in_water = getPropValue<TypeTag, Properties::EnableDisgasInWater>();
121 static constexpr bool has_saltPrecip = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
122 static constexpr bool has_micp = getPropValue<TypeTag, Properties::EnableMICP>();
123
124 // For the conversion between the surface volume rate and reservoir voidage rate
125 using FluidState = BlackOilFluidState<Eval,
129 Indices::compositionSwitchIdx >= 0,
131 has_brine,
134 Indices::numPhases >;
136 WellInterface(const Well& well,
137 const ParallelWellInfo& pw_info,
138 const int time_step,
139 const ModelParameters& param,
140 const RateConverterType& rate_converter,
141 const int pvtRegionIdx,
142 const int num_components,
143 const int num_phases,
144 const int index_of_well,
145 const std::vector<PerforationData>& perf_data);
146
148 virtual ~WellInterface() = default;
149
150 virtual void init(const PhaseUsage* phase_usage_arg,
151 const std::vector<Scalar>& depth_arg,
152 const Scalar gravity_arg,
153 const int num_cells,
154 const std::vector<Scalar>& B_avg,
155 const bool changed_to_open_this_step);
156
158
159 virtual ConvergenceReport getWellConvergence(const SummaryState& summary_state,
160 const WellState<Scalar>& well_state,
161 const std::vector<Scalar>& B_avg,
162 DeferredLogger& deferred_logger,
163 const bool relax_tolerance) const = 0;
164
165 virtual void solveEqAndUpdateWellState(const SummaryState& summary_state,
166 WellState<Scalar>& well_state,
167 DeferredLogger& deferred_logger) = 0;
168
169 void assembleWellEq(const Simulator& simulator,
170 const double dt,
171 WellState<Scalar>& well_state,
172 const GroupState<Scalar>& group_state,
173 DeferredLogger& deferred_logger);
174
175 void assembleWellEqWithoutIteration(const Simulator& simulator,
176 const double dt,
177 WellState<Scalar>& well_state,
178 const GroupState<Scalar>& group_state,
179 DeferredLogger& deferred_logger);
180
181 // TODO: better name or further refactoring the function to make it more clear
182 void prepareWellBeforeAssembling(const Simulator& simulator,
183 const double dt,
184 WellState<Scalar>& well_state,
185 const GroupState<Scalar>& group_state,
186 DeferredLogger& deferred_logger);
187
188
189 virtual void computeWellRatesWithBhp(const Simulator& ebosSimulator,
190 const Scalar& bhp,
191 std::vector<Scalar>& well_flux,
192 DeferredLogger& deferred_logger) const = 0;
193
194 virtual std::optional<Scalar>
196 const SummaryState& summary_state,
197 const Scalar alq_value,
198 DeferredLogger& deferred_logger) const = 0;
199
202 virtual void recoverWellSolutionAndUpdateWellState(const SummaryState& summary_state,
203 const BVector& x,
204 WellState<Scalar>& well_state,
205 DeferredLogger& deferred_logger) = 0;
206
208 virtual void apply(const BVector& x, BVector& Ax) const = 0;
209
211 virtual void apply(BVector& r) const = 0;
212
213 // TODO: before we decide to put more information under mutable, this function is not const
214 virtual void computeWellPotentials(const Simulator& simulator,
215 const WellState<Scalar>& well_state,
216 std::vector<Scalar>& well_potentials,
217 DeferredLogger& deferred_logger) = 0;
218
219 virtual void updateWellStateWithTarget(const Simulator& simulator,
220 const GroupState<Scalar>& group_state,
221 WellState<Scalar>& well_state,
222 DeferredLogger& deferred_logger) const;
223
224 virtual void computeWellRatesWithBhpIterations(const Simulator& simulator,
225 const Scalar& bhp,
226 std::vector<Scalar>& well_flux,
227 DeferredLogger& deferred_logger) const = 0;
228
229 bool updateWellStateWithTHPTargetProd(const Simulator& simulator,
230 WellState<Scalar>& well_state,
231 DeferredLogger& deferred_logger) const;
232
233 enum class IndividualOrGroup { Individual, Group, Both };
234 bool updateWellControl(const Simulator& simulator,
235 const IndividualOrGroup iog,
236 WellState<Scalar>& well_state,
237 const GroupState<Scalar>& group_state,
238 DeferredLogger& deferred_logger) /* const */;
239
241 WellState<Scalar>& well_state,
242 const GroupState<Scalar>& group_state,
243 const Well::InjectionControls& inj_controls,
244 const Well::ProductionControls& prod_controls,
245 const Scalar WQTotal,
246 DeferredLogger& deferred_logger,
247 const bool fixed_control = false,
248 const bool fixed_status = false);
249
250 virtual void updatePrimaryVariables(const SummaryState& summary_state,
251 const WellState<Scalar>& well_state,
252 DeferredLogger& deferred_logger) = 0;
253
254 virtual void calculateExplicitQuantities(const Simulator& simulator,
255 const WellState<Scalar>& well_state,
256 DeferredLogger& deferred_logger) = 0; // should be const?
257
258 virtual void updateProductivityIndex(const Simulator& simulator,
259 const WellProdIndexCalculator& wellPICalc,
260 WellState<Scalar>& well_state,
261 DeferredLogger& deferred_logger) const = 0;
262
263 virtual Scalar connectionDensity(const int globalConnIdx,
264 const int openConnIdx) const = 0;
265
268 {
269 return false;
270 }
271
272 // Add well contributions to matrix
274
276 const BVector& x,
277 const int pressureVarIndex,
278 const bool use_well_weights,
279 const WellState<Scalar>& well_state) const = 0;
280
281 void addCellRates(RateVector& rates, int cellIdx) const;
282
283 Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const;
284
285 // TODO: theoretically, it should be a const function
286 // Simulator is not const is because that assembleWellEq is non-const Simulator
287 void wellTesting(const Simulator& simulator,
288 const double simulation_time,
289 /* const */ WellState<Scalar>& well_state,
290 const GroupState<Scalar>& group_state,
291 WellTestState& welltest_state,
292 DeferredLogger& deferred_logger);
293
294 void checkWellOperability(const Simulator& simulator,
295 const WellState<Scalar>& well_state,
296 DeferredLogger& deferred_logger);
297
299 const double dt,
300 WellState<Scalar>& well_state,
301 const GroupState<Scalar>& group_state,
302 DeferredLogger& deferred_logger);
303
305 WellState<Scalar>& well_state,
306 DeferredLogger& deferred_logger);
307
308 // check whether the well is operable under the current reservoir condition
309 // mostly related to BHP limit and THP limit
310 void updateWellOperability(const Simulator& simulator,
311 const WellState<Scalar>& well_state,
312 DeferredLogger& deferred_logger);
313
314 bool updateWellOperabilityFromWellEq(const Simulator& simulator,
315 const WellState<Scalar>& well_state,
316 DeferredLogger& deferred_logger);
317
318 // update perforation water throughput based on solved water rate
319 virtual void updateWaterThroughput(const double dt,
320 WellState<Scalar>& well_state) const = 0;
321
324 virtual std::vector<Scalar>
326 DeferredLogger& deferred_logger) const = 0;
327
331 void updateWellStateRates(const Simulator& simulator,
332 WellState<Scalar>& well_state,
333 DeferredLogger& deferred_logger) const;
334
335 void solveWellEquation(const Simulator& simulator,
336 WellState<Scalar>& well_state,
337 const GroupState<Scalar>& group_state,
338 DeferredLogger& deferred_logger);
339
340 const std::vector<RateVector>& connectionRates() const
341 {
342 return connectionRates_;
343 }
344
345 virtual std::vector<Scalar> getPrimaryVars() const
346 {
347 return {};
348 }
349
350 virtual int setPrimaryVars(typename std::vector<Scalar>::const_iterator)
351 {
352 return 0;
353 }
354
355 std::vector<Scalar> wellIndex(const int perf,
356 const IntensiveQuantities& intQuants,
357 const Scalar trans_mult,
358 const SingleWellState<Scalar>& ws) const;
359
360 void updateConnectionDFactor(const Simulator& simulator,
361 SingleWellState<Scalar>& ws) const;
362
364 SingleWellState<Scalar>& ws) const;
365
366protected:
367 // simulation parameters
369 std::vector<RateVector> connectionRates_;
370 std::vector<Scalar> B_avg_;
373
374 Scalar wpolymer() const;
375 Scalar wfoam() const;
376 Scalar wsalt() const;
377 Scalar wmicrobes() const;
378 Scalar woxygen() const;
379 Scalar wurea() const;
380
381 virtual Scalar getRefDensity() const = 0;
382
383 // Component fractions for each phase for the well
384 const std::vector<Scalar>& compFrac() const;
385
386 std::vector<Scalar>
387 initialWellRateFractions(const Simulator& ebosSimulator,
388 const WellState<Scalar>& well_state) const;
389
390 // check whether the well is operable under BHP limit with current reservoir condition
391 virtual void checkOperabilityUnderBHPLimit(const WellState<Scalar>& well_state,
392 const Simulator& simulator,
393 DeferredLogger& deferred_logger) = 0;
394
395 // check whether the well is operable under THP limit with current reservoir condition
396 virtual void checkOperabilityUnderTHPLimit(const Simulator& simulator,
397 const WellState<Scalar>& well_state,
398 DeferredLogger& deferred_logger) = 0;
399
400 virtual void updateIPR(const Simulator& simulator,
401 DeferredLogger& deferred_logger) const = 0;
402
403 virtual void assembleWellEqWithoutIteration(const Simulator& simulator,
404 const double dt,
405 const WellInjectionControls& inj_controls,
406 const WellProductionControls& prod_controls,
407 WellState<Scalar>& well_state,
408 const GroupState<Scalar>& group_state,
409 DeferredLogger& deferred_logger) = 0;
410
411 // iterate well equations with the specified control until converged
412 virtual bool iterateWellEqWithControl(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 virtual bool iterateWellEqWithSwitching(const Simulator& simulator,
421 const double dt,
422 const WellInjectionControls& inj_controls,
423 const WellProductionControls& prod_controls,
424 WellState<Scalar>& well_state,
425 const GroupState<Scalar>& group_state,
426 DeferredLogger& deferred_logger,
427 const bool fixed_control = false,
428 const bool fixed_status = false) = 0;
429
430 virtual void updateIPRImplicit(const Simulator& simulator,
431 WellState<Scalar>& well_state,
432 DeferredLogger& deferred_logger) = 0;
433
434 bool iterateWellEquations(const Simulator& simulator,
435 const double dt,
436 WellState<Scalar>& well_state,
437 const GroupState<Scalar>& group_state,
438 DeferredLogger& deferred_logger);
439
440 bool solveWellWithTHPConstraint(const Simulator& simulator,
441 const double dt,
442 const Well::InjectionControls& inj_controls,
443 const Well::ProductionControls& prod_controls,
444 WellState<Scalar>& well_state,
445 const GroupState<Scalar>& group_state,
446 DeferredLogger& deferred_logger);
447
448 std::optional<Scalar>
449 estimateOperableBhp(const Simulator& ebos_simulator,
450 const double dt,
451 WellState<Scalar>& well_state,
452 const SummaryState& summary_state,
453 DeferredLogger& deferred_logger);
454
455 bool solveWellWithBhp(const Simulator& simulator,
456 const double dt,
457 const Scalar bhp,
458 WellState<Scalar>& well_state,
459 DeferredLogger& deferred_logger);
460
461 bool solveWellWithZeroRate(const Simulator& simulator,
462 const double dt,
463 WellState<Scalar>& well_state,
464 DeferredLogger& deferred_logger);
465
466 bool solveWellForTesting(const Simulator& simulator,
467 WellState<Scalar>& well_state,
468 const GroupState<Scalar>& group_state,
469 DeferredLogger& deferred_logger);
470
471 Eval getPerfCellPressure(const FluidState& fs) const;
472
473 // get the mobility for specific perforation
474 template<class Value, class Callback>
475 void getMobility(const Simulator& simulator,
476 const int perf,
477 std::vector<Value>& mob,
478 Callback& extendEval,
479 [[maybe_unused]] DeferredLogger& deferred_logger) const;
480
481 void computeConnLevelProdInd(const FluidState& fs,
482 const std::function<Scalar(const Scalar)>& connPICalc,
483 const std::vector<Scalar>& mobility,
484 Scalar* connPI) const;
485
486 void computeConnLevelInjInd(const FluidState& fs,
487 const Phase preferred_phase,
488 const std::function<Scalar(const Scalar)>& connIICalc,
489 const std::vector<Scalar>& mobility,
490 Scalar* connII,
491 DeferredLogger& deferred_logger) const;
492
493 Scalar computeConnectionDFactor(const int perf,
494 const IntensiveQuantities& intQuants,
495 const SingleWellState<Scalar>& ws) const;
496};
497
498} // namespace Opm
499
500#ifndef OPM_WELLINTERFACE_IMPL_HEADER_INCLUDED
501#include "WellInterface_impl.hpp"
502#endif
503
504#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:36
Definition: GroupState.hpp:35
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:184
Definition: RateConverter.hpp:70
Definition: SingleWellState.hpp:41
Definition: WellInterfaceFluidSystem.hpp:48
static constexpr int Oil
Definition: WellInterfaceFluidSystem.hpp:61
static constexpr int Water
Definition: WellInterfaceFluidSystem.hpp:60
static constexpr int Gas
Definition: WellInterfaceFluidSystem.hpp:62
int pvtRegionIdx() const
Definition: WellInterfaceGeneric.hpp:118
const std::vector< FluidSystem::Scalar > & wellIndex() const
Definition: WellInterfaceGeneric.hpp:140
Definition: WellInterfaceIndices.hpp:34
DenseAd::Evaluation< Scalar, Indices::numEq > Eval
Definition: WellInterfaceIndices.hpp:40
Definition: WellInterface.hpp:75
const std::vector< Scalar > & compFrac() const
virtual void solveEqAndUpdateWellState(const SummaryState &summary_state, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
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
virtual std::vector< Scalar > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const =0
const ModelParameters & param_
Definition: WellInterface.hpp:368
Scalar woxygen() const
Definition: WellInterface_impl.hpp:164
static constexpr bool has_brine
Definition: WellInterface.hpp:118
IndividualOrGroup
Definition: WellInterface.hpp:233
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:1712
virtual void computeWellPotentials(const Simulator &simulator, const WellState< Scalar > &well_state, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger)=0
void assembleWellEq(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:770
typename GasLiftSingleWellGeneric< Scalar >::GLiftSyncGroups GLiftSyncGroups
Definition: WellInterface.hpp:95
virtual void computeWellRatesWithBhp(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const =0
virtual void addWellContributions(SparseMatrixAdapter &) const =0
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:104
const std::vector< RateVector > & connectionRates() const
Definition: WellInterface.hpp:340
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:134
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
Scalar wfoam() const
Definition: WellInterface_impl.hpp:126
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:369
void computeConnLevelProdInd(const FluidState &fs, const std::function< Scalar(const Scalar)> &connPICalc, const std::vector< Scalar > &mobility, Scalar *connPI) const
Definition: WellInterface_impl.hpp:1817
bool updateWellStateWithTHPTargetProd(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1787
void addCellRates(RateVector &rates, int cellIdx) const
Definition: WellInterface_impl.hpp:862
virtual void calculateExplicitQuantities(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:101
Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const
Definition: WellInterface_impl.hpp:878
virtual std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &ebos_simulator, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger) const =0
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 init(const PhaseUsage *phase_usage_arg, const std::vector< Scalar > &depth_arg, const Scalar gravity_arg, const int num_cells, const std::vector< Scalar > &B_avg, const bool changed_to_open_this_step)
Definition: WellInterface_impl.hpp:91
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:114
bool iterateWellEquations(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:469
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:86
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
void checkWellOperability(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:898
std::optional< Scalar > estimateOperableBhp(const Simulator &ebos_simulator, const double dt, WellState< Scalar > &well_state, const SummaryState &summary_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:599
void updateWellOperability(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1001
void updateWellStateRates(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1479
static constexpr bool has_watVapor
Definition: WellInterface.hpp:119
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: WellInterface.hpp:87
virtual void updateWellStateWithTarget(const Simulator &simulator, const GroupState< Scalar > &group_state, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1055
virtual Scalar getRefDensity() const =0
Eval getPerfCellPressure(const FluidState &fs) const
Definition: WellInterface_impl.hpp:1697
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:100
std::vector< Scalar > B_avg_
Definition: WellInterface.hpp:370
bool changed_to_stopped_this_step_
Definition: WellInterface.hpp:371
static constexpr bool has_polymermw
Definition: WellInterface.hpp:116
static constexpr bool has_polymer
Definition: WellInterface.hpp:112
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< Scalar > &ws) const
Definition: WellInterface_impl.hpp:1668
virtual ConvergenceReport getWellConvergence(const SummaryState &summary_state, const WellState< Scalar > &well_state, const std::vector< Scalar > &B_avg, DeferredLogger &deferred_logger, const bool relax_tolerance) const =0
virtual ~WellInterface()=default
Virtual destructor.
virtual void computeWellRatesWithBhpIterations(const Simulator &simulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const =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:92
virtual void updateIPRImplicit(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:84
static constexpr bool has_solvent
Definition: WellInterface.hpp:110
static constexpr bool has_disgas_in_water
Definition: WellInterface.hpp:120
void prepareWellBeforeAssembling(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:804
static constexpr bool has_saltPrecip
Definition: WellInterface.hpp:121
bool gliftBeginTimeStepWellTestIterateWellEquations(const Simulator &ebos_simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:923
virtual int setPrimaryVars(typename std::vector< Scalar >::const_iterator)
Definition: WellInterface.hpp:350
static constexpr bool has_foam
Definition: WellInterface.hpp:117
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: WellInterface.hpp:89
virtual void updateWaterThroughput(const double dt, WellState< Scalar > &well_state) const =0
void computeConnLevelInjInd(const FluidState &fs, const Phase preferred_phase, const std::function< Scalar(const Scalar)> &connIICalc, const std::vector< Scalar > &mobility, Scalar *connII, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1852
void solveWellEquation(const Simulator &simulator, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:713
void updateConnectionDFactor(const Simulator &simulator, SingleWellState< Scalar > &ws) const
Definition: WellInterface_impl.hpp:1606
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
std::vector< Scalar > initialWellRateFractions(const Simulator &ebosSimulator, const WellState< Scalar > &well_state) const
Definition: WellInterface_impl.hpp:1432
Scalar wsalt() const
Definition: WellInterface_impl.hpp:140
static constexpr bool has_micp
Definition: WellInterface.hpp:122
virtual Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const =0
virtual bool jacobianContainsWellContributions() const
Wether the Jacobian will also have well contributions in it.
Definition: WellInterface.hpp:267
bool updateWellOperabilityFromWellEq(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1037
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
Dune::FieldMatrix< Scalar, Indices::numEq, Indices::numEq > MatrixBlockType
Definition: WellInterface.hpp:98
virtual std::vector< Scalar > getPrimaryVars() const
Definition: WellInterface.hpp:345
void gliftBeginTimeStepWellTestUpdateALQ(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:957
virtual void updatePrimaryVariables(const SummaryState &summary_state, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
typename Base::Eval Eval
Definition: WellInterface.hpp:99
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:365
static constexpr bool has_energy
Definition: WellInterface.hpp:113
Scalar wpolymer() const
Definition: WellInterface_impl.hpp:110
Scalar computeConnectionDFactor(const int perf, const IntensiveQuantities &intQuants, const SingleWellState< Scalar > &ws) const
Definition: WellInterface_impl.hpp:1626
bool solveWellWithBhp(const Simulator &simulator, const double dt, const Scalar bhp, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:622
bool solveWellWithZeroRate(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:658
bool solveWellForTesting(const Simulator &simulator, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:678
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:85
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:502
Scalar wurea() const
Definition: WellInterface_impl.hpp:182
bool thp_update_iterations
Definition: WellInterface.hpp:372
Dune::FieldVector< Scalar, Indices::numEq > VectorBlockType
Definition: WellInterface.hpp:97
Scalar 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:785
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:88
typename BlackoilWellModel< TypeTag >::GLiftOptWells GLiftOptWells
Definition: WellInterface.hpp:91
virtual void recoverWellSolutionAndUpdateWellState(const SummaryState &summary_state, const BVector &x, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
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 Scalar WQTotal, DeferredLogger &deferred_logger, const bool fixed_control=false, const bool fixed_status=false)
Definition: WellInterface_impl.hpp:271
virtual void checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
static constexpr bool has_zFraction
Definition: WellInterface.hpp:111
virtual void initPrimaryVariablesEvaluation()=0
typename BlackoilWellModel< TypeTag >::GLiftWellStateMap GLiftWellStateMap
Definition: WellInterface.hpp:94
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:485
Definition: BlackoilPhases.hpp:46