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
50
52
55
56#include <dune/common/fmatrix.hh>
57#include <dune/istl/bcrsmatrix.hh>
58#include <dune/istl/matrixmatrix.hh>
59
60#include <opm/material/densead/Evaluation.hpp>
61
62#include <vector>
63
64namespace Opm
65{
66
67class WellInjectionProperties;
68class WellProductionProperties;
69
70template<typename TypeTag>
71class WellInterface : public WellInterfaceIndices<GetPropType<TypeTag, Properties::FluidSystem>,
72 GetPropType<TypeTag, Properties::Indices>>
73{
76public:
92
93 using VectorBlockType = Dune::FieldVector<Scalar, Indices::numEq>;
94 using MatrixBlockType = Dune::FieldMatrix<Scalar, Indices::numEq, Indices::numEq>;
95 using Eval = typename Base::Eval;
96 using BVector = Dune::BlockVector<VectorBlockType>;
97 using PressureMatrix = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, 1, 1>>;
98
101
105
107
108 static constexpr bool has_solvent = getPropValue<TypeTag, Properties::EnableSolvent>();
109 static constexpr bool has_zFraction = getPropValue<TypeTag, Properties::EnableExtbo>();
110 static constexpr bool has_polymer = getPropValue<TypeTag, Properties::EnablePolymer>();
111 static constexpr bool has_energy = getPropValue<TypeTag, Properties::EnableEnergy>();
112 static const bool has_temperature = getPropValue<TypeTag, Properties::EnableTemperature>();
113 // flag for polymer molecular weight related
114 static constexpr bool has_polymermw = getPropValue<TypeTag, Properties::EnablePolymerMW>();
115 static constexpr bool has_foam = getPropValue<TypeTag, Properties::EnableFoam>();
116 static constexpr bool has_brine = getPropValue<TypeTag, Properties::EnableBrine>();
117 static constexpr bool has_watVapor = getPropValue<TypeTag, Properties::EnableVapwat>();
118 static constexpr bool has_disgas_in_water = getPropValue<TypeTag, Properties::EnableDisgasInWater>();
119 static constexpr bool has_saltPrecip = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
120 static constexpr bool has_micp = getPropValue<TypeTag, Properties::EnableMICP>();
121
122 // For the conversion between the surface volume rate and reservoir voidage rate
123 using FluidState = BlackOilFluidState<Eval,
127 Indices::compositionSwitchIdx >= 0,
129 has_brine,
132 Indices::numPhases >;
134 WellInterface(const Well& well,
135 const ParallelWellInfo<Scalar>& pw_info,
136 const int time_step,
137 const ModelParameters& param,
138 const RateConverterType& rate_converter,
139 const int pvtRegionIdx,
140 const int num_components,
141 const int num_phases,
142 const int index_of_well,
143 const std::vector<PerforationData<Scalar>>& perf_data);
144
146 virtual ~WellInterface() = default;
147
148 virtual void init(const PhaseUsage* phase_usage_arg,
149 const std::vector<Scalar>& depth_arg,
150 const Scalar gravity_arg,
151 const int num_cells,
152 const std::vector<Scalar>& B_avg,
153 const bool changed_to_open_this_step);
154
156
158 const WellState<Scalar>& well_state,
159 const std::vector<Scalar>& B_avg,
160 DeferredLogger& deferred_logger,
161 const bool relax_tolerance) const = 0;
162
163 virtual void solveEqAndUpdateWellState(const Simulator& simulator,
164 WellState<Scalar>& well_state,
165 DeferredLogger& deferred_logger) = 0;
166
167 void assembleWellEq(const Simulator& simulator,
168 const double dt,
169 WellState<Scalar>& well_state,
170 const GroupState<Scalar>& group_state,
171 DeferredLogger& deferred_logger);
172
173 void assembleWellEqWithoutIteration(const Simulator& simulator,
174 const double dt,
175 WellState<Scalar>& well_state,
176 const GroupState<Scalar>& group_state,
177 DeferredLogger& deferred_logger);
178
179 // TODO: better name or further refactoring the function to make it more clear
180 void prepareWellBeforeAssembling(const Simulator& simulator,
181 const double dt,
182 WellState<Scalar>& well_state,
183 const GroupState<Scalar>& group_state,
184 DeferredLogger& deferred_logger);
185
186
187 virtual void computeWellRatesWithBhp(const Simulator& ebosSimulator,
188 const Scalar& bhp,
189 std::vector<Scalar>& well_flux,
190 DeferredLogger& deferred_logger) const = 0;
191
192 virtual std::optional<Scalar>
194 const SummaryState& summary_state,
195 const Scalar alq_value,
196 DeferredLogger& deferred_logger) const = 0;
197
201 const BVector& x,
202 WellState<Scalar>& well_state,
203 DeferredLogger& deferred_logger) = 0;
204
206 virtual void apply(const BVector& x, BVector& Ax) const = 0;
207
209 virtual void apply(BVector& r) const = 0;
210
211 // TODO: before we decide to put more information under mutable, this function is not const
212 virtual void computeWellPotentials(const Simulator& simulator,
213 const WellState<Scalar>& well_state,
214 std::vector<Scalar>& well_potentials,
215 DeferredLogger& deferred_logger) = 0;
216
217 virtual void updateWellStateWithTarget(const Simulator& simulator,
218 const GroupState<Scalar>& group_state,
219 WellState<Scalar>& well_state,
220 DeferredLogger& deferred_logger) const;
221
222 virtual void computeWellRatesWithBhpIterations(const Simulator& simulator,
223 const Scalar& bhp,
224 std::vector<Scalar>& well_flux,
225 DeferredLogger& deferred_logger) const = 0;
226
227 bool wellUnderZeroRateTarget(const Simulator& simulator,
228 const WellState<Scalar>& well_state,
229 DeferredLogger& deferred_logger) const;
230
231 bool wellUnderZeroGroupRateTarget(const Simulator& simulator,
232 const WellState<Scalar>& well_state,
233 DeferredLogger& deferred_logger,
234 std::optional<bool> group_control = std::nullopt) const;
235
236 bool stoppedOrZeroRateTarget(const Simulator& simulator,
237 const WellState<Scalar>& well_state,
238 DeferredLogger& deferred_logger) const;
239
240 bool updateWellStateWithTHPTargetProd(const Simulator& simulator,
241 WellState<Scalar>& well_state,
242 DeferredLogger& deferred_logger) const;
243
244 enum class IndividualOrGroup { Individual, Group, Both };
245 bool updateWellControl(const Simulator& simulator,
246 const IndividualOrGroup iog,
247 WellState<Scalar>& well_state,
248 const GroupState<Scalar>& group_state,
249 DeferredLogger& deferred_logger) /* const */;
250
252 WellState<Scalar>& well_state,
253 const GroupState<Scalar>& group_state,
254 const Well::InjectionControls& inj_controls,
255 const Well::ProductionControls& prod_controls,
256 const Scalar WQTotal,
257 DeferredLogger& deferred_logger,
258 const bool fixed_control = false,
259 const bool fixed_status = false);
260
261 virtual void updatePrimaryVariables(const Simulator& simulator,
262 const WellState<Scalar>& well_state,
263 DeferredLogger& deferred_logger) = 0;
264
265 virtual void calculateExplicitQuantities(const Simulator& simulator,
266 const WellState<Scalar>& well_state,
267 DeferredLogger& deferred_logger) = 0; // should be const?
268
269 virtual void updateProductivityIndex(const Simulator& simulator,
270 const WellProdIndexCalculator<Scalar>& wellPICalc,
271 WellState<Scalar>& well_state,
272 DeferredLogger& deferred_logger) const = 0;
273
274 virtual Scalar connectionDensity(const int globalConnIdx,
275 const int openConnIdx) const = 0;
276
279 {
280 return false;
281 }
282
283 // Add well contributions to matrix
285
287 const BVector& x,
288 const int pressureVarIndex,
289 const bool use_well_weights,
290 const WellState<Scalar>& well_state) const = 0;
291
292 void addCellRates(RateVector& rates, int cellIdx) const;
293
294 Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const;
295
296 // TODO: theoretically, it should be a const function
297 // Simulator is not const is because that assembleWellEq is non-const Simulator
298 void wellTesting(const Simulator& simulator,
299 const double simulation_time,
300 /* const */ WellState<Scalar>& well_state,
301 const GroupState<Scalar>& group_state,
302 WellTestState& welltest_state,
303 DeferredLogger& deferred_logger);
304
305 void checkWellOperability(const Simulator& simulator,
306 const WellState<Scalar>& well_state,
307 DeferredLogger& deferred_logger);
308
310 const double dt,
311 WellState<Scalar>& well_state,
312 const GroupState<Scalar>& group_state,
313 DeferredLogger& deferred_logger);
314
316 WellState<Scalar>& well_state,
317 DeferredLogger& deferred_logger);
318
319 // check whether the well is operable under the current reservoir condition
320 // mostly related to BHP limit and THP limit
321 void updateWellOperability(const Simulator& simulator,
322 const WellState<Scalar>& well_state,
323 DeferredLogger& deferred_logger);
324
325 bool updateWellOperabilityFromWellEq(const Simulator& simulator,
326 const WellState<Scalar>& well_state,
327 DeferredLogger& deferred_logger);
328
329 // update perforation water throughput based on solved water rate
330 virtual void updateWaterThroughput(const double dt,
331 WellState<Scalar>& well_state) const = 0;
332
335 virtual std::vector<Scalar>
337 DeferredLogger& deferred_logger) const = 0;
338
342 void updateWellStateRates(const Simulator& simulator,
343 WellState<Scalar>& well_state,
344 DeferredLogger& deferred_logger) const;
345
346 void solveWellEquation(const Simulator& simulator,
347 WellState<Scalar>& well_state,
348 const GroupState<Scalar>& group_state,
349 DeferredLogger& deferred_logger);
350
351 const std::vector<RateVector>& connectionRates() const
352 {
353 return connectionRates_;
354 }
355
356 virtual std::vector<Scalar> getPrimaryVars() const
357 {
358 return {};
359 }
360
361 virtual int setPrimaryVars(typename std::vector<Scalar>::const_iterator)
362 {
363 return 0;
364 }
365
366 std::vector<Scalar> wellIndex(const int perf,
367 const IntensiveQuantities& intQuants,
368 const Scalar trans_mult,
369 const SingleWellState<Scalar>& ws) const;
370
371 void updateConnectionDFactor(const Simulator& simulator,
372 SingleWellState<Scalar>& ws) const;
373
375 SingleWellState<Scalar>& ws) const;
376
377 virtual bool iterateWellEqWithSwitching(const Simulator& simulator,
378 const double dt,
379 const WellInjectionControls& inj_controls,
380 const WellProductionControls& prod_controls,
381 WellState<Scalar>& well_state,
382 const GroupState<Scalar>& group_state,
383 DeferredLogger& deferred_logger,
384 const bool fixed_control = false,
385 const bool fixed_status = false) = 0;
386protected:
387 // simulation parameters
389 std::vector<RateVector> connectionRates_;
390 std::vector<Scalar> B_avg_;
393
394 Scalar wpolymer() const;
395 Scalar wfoam() const;
396 Scalar wsalt() const;
397 Scalar wmicrobes() const;
398 Scalar woxygen() const;
399 Scalar wurea() const;
400
401 virtual Scalar getRefDensity() const = 0;
402
403 // Component fractions for each phase for the well
404 const std::vector<Scalar>& compFrac() const;
405
406 std::vector<Scalar>
407 initialWellRateFractions(const Simulator& ebosSimulator,
408 const WellState<Scalar>& well_state) const;
409
410 // check whether the well is operable under BHP limit with current reservoir condition
411 virtual void checkOperabilityUnderBHPLimit(const WellState<Scalar>& well_state,
412 const Simulator& simulator,
413 DeferredLogger& deferred_logger) = 0;
414
415 // check whether the well is operable under THP limit with current reservoir condition
416 virtual void checkOperabilityUnderTHPLimit(const Simulator& simulator,
417 const WellState<Scalar>& well_state,
418 DeferredLogger& deferred_logger) = 0;
419
420 virtual void updateIPR(const Simulator& simulator,
421 DeferredLogger& deferred_logger) const = 0;
422
423 virtual void assembleWellEqWithoutIteration(const Simulator& simulator,
424 const double dt,
425 const WellInjectionControls& inj_controls,
426 const WellProductionControls& prod_controls,
427 WellState<Scalar>& well_state,
428 const GroupState<Scalar>& group_state,
429 DeferredLogger& deferred_logger) = 0;
430
431 // iterate well equations with the specified control until converged
432 virtual bool iterateWellEqWithControl(const Simulator& simulator,
433 const double dt,
434 const WellInjectionControls& inj_controls,
435 const WellProductionControls& prod_controls,
436 WellState<Scalar>& well_state,
437 const GroupState<Scalar>& group_state,
438 DeferredLogger& deferred_logger) = 0;
439
440 virtual void updateIPRImplicit(const Simulator& simulator,
441 WellState<Scalar>& well_state,
442 DeferredLogger& deferred_logger) = 0;
443
444 bool iterateWellEquations(const Simulator& simulator,
445 const double dt,
446 WellState<Scalar>& well_state,
447 const GroupState<Scalar>& group_state,
448 DeferredLogger& deferred_logger);
449
450 bool solveWellWithTHPConstraint(const Simulator& simulator,
451 const double dt,
452 const Well::InjectionControls& inj_controls,
453 const Well::ProductionControls& prod_controls,
454 WellState<Scalar>& well_state,
455 const GroupState<Scalar>& group_state,
456 DeferredLogger& deferred_logger);
457
458 std::optional<Scalar>
459 estimateOperableBhp(const Simulator& ebos_simulator,
460 const double dt,
461 WellState<Scalar>& well_state,
462 const SummaryState& summary_state,
463 DeferredLogger& deferred_logger);
464
465 bool solveWellWithBhp(const Simulator& simulator,
466 const double dt,
467 const Scalar bhp,
468 WellState<Scalar>& well_state,
469 DeferredLogger& deferred_logger);
470
471 bool solveWellWithZeroRate(const Simulator& simulator,
472 const double dt,
473 WellState<Scalar>& well_state,
474 DeferredLogger& deferred_logger);
475
476 bool solveWellForTesting(const Simulator& simulator,
477 WellState<Scalar>& well_state,
478 const GroupState<Scalar>& group_state,
479 DeferredLogger& deferred_logger);
480
481 Eval getPerfCellPressure(const FluidState& fs) const;
482
483 // get the mobility for specific perforation
484 template<class Value, class Callback>
485 void getMobility(const Simulator& simulator,
486 const int perf,
487 std::vector<Value>& mob,
488 Callback& extendEval,
489 [[maybe_unused]] DeferredLogger& deferred_logger) const;
490
491 void computeConnLevelProdInd(const FluidState& fs,
492 const std::function<Scalar(const Scalar)>& connPICalc,
493 const std::vector<Scalar>& mobility,
494 Scalar* connPI) const;
495
496 void computeConnLevelInjInd(const FluidState& fs,
497 const Phase preferred_phase,
498 const std::function<Scalar(const Scalar)>& connIICalc,
499 const std::vector<Scalar>& mobility,
500 Scalar* connII,
501 DeferredLogger& deferred_logger) const;
502
503 Scalar computeConnectionDFactor(const int perf,
504 const IntensiveQuantities& intQuants,
505 const SingleWellState<Scalar>& ws) const;
506};
507
508} // namespace Opm
509
510#ifndef OPM_WELLINTERFACE_IMPL_HEADER_INCLUDED
511#include "WellInterface_impl.hpp"
512#endif
513
514#endif // OPM_WELLINTERFACE_HEADER_INCLUDED
typename BlackoilWellModelGeneric< Scalar >::GLiftProdWells GLiftProdWells
Definition: BlackoilWellModel.hpp:113
typename BlackoilWellModelGeneric< Scalar >::GLiftOptWells GLiftOptWells
Definition: BlackoilWellModel.hpp:112
typename BlackoilWellModelGeneric< Scalar >::GLiftWellStateMap GLiftWellStateMap
Definition: BlackoilWellModel.hpp:115
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
std::set< int > GLiftSyncGroups
Definition: GasLiftSingleWellGeneric.hpp:59
Definition: GasLiftSingleWell.hpp:36
Definition: GroupState.hpp:38
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:186
Definition: RateConverter.hpp:71
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:92
Definition: SingleWellState.hpp:42
Definition: WellInterfaceFluidSystem.hpp:51
static constexpr int Oil
Definition: WellInterfaceFluidSystem.hpp:64
static constexpr int Water
Definition: WellInterfaceFluidSystem.hpp:63
static constexpr int Gas
Definition: WellInterfaceFluidSystem.hpp:65
int pvtRegionIdx() const
Definition: WellInterfaceGeneric.hpp:119
const std::vector< FluidSystem::Scalar > & wellIndex() const
Definition: WellInterfaceGeneric.hpp:141
Definition: WellInterfaceIndices.hpp:34
DenseAd::Evaluation< Scalar, Indices::numEq > Eval
Definition: WellInterfaceIndices.hpp:40
Definition: WellInterface.hpp:73
const std::vector< Scalar > & compFrac() const
bool updateWellControl(const Simulator &simulator, const IndividualOrGroup iog, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:195
virtual std::vector< Scalar > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const =0
const ModelParameters & param_
Definition: WellInterface.hpp:388
Scalar woxygen() const
Definition: WellInterface_impl.hpp:165
static constexpr bool has_brine
Definition: WellInterface.hpp:116
IndividualOrGroup
Definition: WellInterface.hpp:244
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:78
virtual ConvergenceReport getWellConvergence(const Simulator &simulator, const WellState< Scalar > &well_state, const std::vector< Scalar > &B_avg, DeferredLogger &deferred_logger, const bool relax_tolerance) const =0
WellInterface(const Well &well, const ParallelWellInfo< Scalar > &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< Scalar > > &perf_data)
Constructor.
Definition: WellInterface_impl.hpp:55
virtual void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1773
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:782
typename GasLiftSingleWellGeneric< Scalar >::GLiftSyncGroups GLiftSyncGroups
Definition: WellInterface.hpp:91
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:100
const std::vector< RateVector > & connectionRates() const
Definition: WellInterface.hpp:351
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:132
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:127
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:389
bool wellUnderZeroRateTarget(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1444
void computeConnLevelProdInd(const FluidState &fs, const std::function< Scalar(const Scalar)> &connPICalc, const std::vector< Scalar > &mobility, Scalar *connPI) const
Definition: WellInterface_impl.hpp:1878
bool updateWellStateWithTHPTargetProd(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1848
void addCellRates(RateVector &rates, int cellIdx) const
Definition: WellInterface_impl.hpp:874
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:97
Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const
Definition: WellInterface_impl.hpp:890
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:92
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:112
virtual void updatePrimaryVariables(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
bool iterateWellEquations(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:481
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:82
virtual void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const =0
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:79
void checkWellOperability(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:910
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:611
void updateWellOperability(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1013
void updateWellStateRates(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1540
static constexpr bool has_watVapor
Definition: WellInterface.hpp:117
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: WellInterface.hpp:83
virtual void updateWellStateWithTarget(const Simulator &simulator, const GroupState< Scalar > &group_state, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1067
virtual Scalar getRefDensity() const =0
Eval getPerfCellPressure(const FluidState &fs) const
Definition: WellInterface_impl.hpp:1758
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:96
std::vector< Scalar > B_avg_
Definition: WellInterface.hpp:390
bool changed_to_stopped_this_step_
Definition: WellInterface.hpp:391
virtual void solveEqAndUpdateWellState(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
static constexpr bool has_polymermw
Definition: WellInterface.hpp:114
static constexpr bool has_polymer
Definition: WellInterface.hpp:110
virtual void apply(BVector &r) const =0
r = r - C D^-1 Rw
GetPropType< TypeTag, Properties::Grid > Grid
Definition: WellInterface.hpp:77
void updateConnectionTransmissibilityFactor(const Simulator &simulator, SingleWellState< Scalar > &ws) const
Definition: WellInterface_impl.hpp:1729
virtual ~WellInterface()=default
Virtual destructor.
bool stoppedOrZeroRateTarget(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1480
virtual void computeWellRatesWithBhpIterations(const Simulator &simulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const =0
typename BlackoilWellModel< TypeTag >::GLiftProdWells GLiftProdWells
Definition: WellInterface.hpp:88
virtual void updateIPRImplicit(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)=0
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:80
static constexpr bool has_solvent
Definition: WellInterface.hpp:108
static constexpr bool has_disgas_in_water
Definition: WellInterface.hpp:118
void prepareWellBeforeAssembling(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:816
static constexpr bool has_saltPrecip
Definition: WellInterface.hpp:119
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:935
virtual int setPrimaryVars(typename std::vector< Scalar >::const_iterator)
Definition: WellInterface.hpp:361
static constexpr bool has_foam
Definition: WellInterface.hpp:115
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: WellInterface.hpp:85
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:1913
void solveWellEquation(const Simulator &simulator, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:725
void updateConnectionDFactor(const Simulator &simulator, SingleWellState< Scalar > &ws) const
Definition: WellInterface_impl.hpp:1667
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:1493
Scalar wsalt() const
Definition: WellInterface_impl.hpp:141
static constexpr bool has_micp
Definition: WellInterface.hpp:120
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:278
bool updateWellOperabilityFromWellEq(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1049
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:94
virtual std::vector< Scalar > getPrimaryVars() const
Definition: WellInterface.hpp:356
void gliftBeginTimeStepWellTestUpdateALQ(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:969
typename Base::Eval Eval
Definition: WellInterface.hpp:95
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:376
static constexpr bool has_energy
Definition: WellInterface.hpp:111
Scalar wpolymer() const
Definition: WellInterface_impl.hpp:111
Scalar computeConnectionDFactor(const int perf, const IntensiveQuantities &intQuants, const SingleWellState< Scalar > &ws) const
Definition: WellInterface_impl.hpp:1687
bool solveWellWithBhp(const Simulator &simulator, const double dt, const Scalar bhp, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:634
bool solveWellWithZeroRate(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:670
bool solveWellForTesting(const Simulator &simulator, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:690
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:81
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:514
Scalar wurea() const
Definition: WellInterface_impl.hpp:183
bool wellUnderZeroGroupRateTarget(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger, std::optional< bool > group_control=std::nullopt) const
Definition: WellInterface_impl.hpp:1461
bool thp_update_iterations
Definition: WellInterface.hpp:392
Dune::FieldVector< Scalar, Indices::numEq > VectorBlockType
Definition: WellInterface.hpp:93
Scalar wmicrobes() const
Definition: WellInterface_impl.hpp:153
void assembleWellEqWithoutIteration(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:797
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:84
typename BlackoilWellModel< TypeTag >::GLiftOptWells GLiftOptWells
Definition: WellInterface.hpp:87
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:109
virtual void initPrimaryVariablesEvaluation()=0
typename BlackoilWellModel< TypeTag >::GLiftWellStateMap GLiftWellStateMap
Definition: WellInterface.hpp:90
Definition: WellProdIndexCalculator.hpp:37
Definition: WellState.hpp:62
Definition: blackoilboundaryratevector.hh:37
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:235
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:147
Static data associated with a well perforation.
Definition: PerforationData.hpp:30
Definition: BlackoilPhases.hpp:46