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
40#include <opm/material/fluidstates/BlackOilFluidState.hpp>
41
43
45
54
56
58
59#include <dune/common/fmatrix.hh>
60#include <dune/istl/bcrsmatrix.hh>
61#include <dune/istl/matrixmatrix.hh>
62
63#include <opm/material/densead/Evaluation.hpp>
64
65#include <vector>
66
67namespace Opm
68{
69
70class WellInjectionProperties;
71class WellProductionProperties;
72
73template<typename TypeTag>
74class WellInterface : public WellInterfaceIndices<GetPropType<TypeTag, Properties::FluidSystem>,
75 GetPropType<TypeTag, Properties::Indices>>
76{
79public:
84 using IndexTraits = typename FluidSystem::IndexTraitsType;
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
104
108
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_bioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
124 static constexpr bool has_micp = Indices::enableMICP;
125
126 // For the conversion between the surface volume rate and reservoir voidage rate
127 using FluidState = BlackOilFluidState<Eval,
131 Indices::compositionSwitchIdx >= 0,
133 has_brine,
136 Indices::numPhases >;
138 WellInterface(const Well& well,
139 const ParallelWellInfo<Scalar>& pw_info,
140 const int time_step,
141 const ModelParameters& param,
142 const RateConverterType& rate_converter,
143 const int pvtRegionIdx,
144 const int num_conservation_quantities,
145 const int num_phases,
146 const int index_of_well,
147 const std::vector<PerforationData<Scalar>>& perf_data);
148
150 virtual ~WellInterface() = default;
151
152 virtual void init(const std::vector<Scalar>& depth_arg,
153 const Scalar gravity_arg,
154 const std::vector<Scalar>& B_avg,
155 const bool changed_to_open_this_step);
156
158 const WellStateType& 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 WellStateType& well_state,
165 DeferredLogger& deferred_logger) = 0;
166
167 void assembleWellEq(const Simulator& simulator,
168 const double dt,
169 WellStateType& well_state,
170 const GroupState<Scalar>& group_state,
171 DeferredLogger& deferred_logger);
172
173 void assembleWellEqWithoutIteration(const Simulator& simulator,
174 const double dt,
175 WellStateType& 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 WellStateType& well_state,
183 const GroupState<Scalar>& group_state,
184 DeferredLogger& deferred_logger);
185
186 virtual void computeWellRatesWithBhp(const Simulator& ebosSimulator,
187 const Scalar& bhp,
188 std::vector<Scalar>& well_flux,
189 DeferredLogger& deferred_logger) const = 0;
190
191 virtual std::optional<Scalar>
193 const SummaryState& summary_state,
194 const Scalar alq_value,
195 DeferredLogger& deferred_logger,
196 bool iterate_if_no_solution) const = 0;
197
201 const BVector& x,
202 WellStateType& 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 WellStateType& 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 WellStateType& well_state,
220 DeferredLogger& deferred_logger) const;
221
222 virtual void scaleSegmentRatesAndPressure(WellStateType& well_state) 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 wellUnderZeroRateTarget(const Simulator& simulator,
230 const WellStateType& well_state,
231 DeferredLogger& deferred_logger) const;
232
233 bool wellUnderZeroGroupRateTarget(const Simulator& simulator,
234 const WellStateType& well_state,
235 DeferredLogger& deferred_logger,
236 std::optional<bool> group_control = std::nullopt) const;
237
238 bool stoppedOrZeroRateTarget(const Simulator& simulator,
239 const WellStateType& well_state,
240 DeferredLogger& deferred_logger) const;
241
242 bool updateWellStateWithTHPTargetProd(const Simulator& simulator,
243 WellStateType& well_state,
244 DeferredLogger& deferred_logger) const;
245
246 enum class IndividualOrGroup { Individual, Group, Both };
247 bool updateWellControl(const Simulator& simulator,
248 const IndividualOrGroup iog,
249 WellStateType& well_state,
250 const GroupState<Scalar>& group_state,
251 DeferredLogger& deferred_logger) /* const */;
252
254 WellStateType& well_state,
255 const GroupState<Scalar>& group_state,
256 const Well::InjectionControls& inj_controls,
257 const Well::ProductionControls& prod_controls,
258 const Scalar WQTotal,
259 DeferredLogger& deferred_logger,
260 const bool fixed_control = false,
261 const bool fixed_status = false);
262
263 virtual void updatePrimaryVariables(const Simulator& simulator,
264 const WellStateType& well_state,
265 DeferredLogger& deferred_logger) = 0;
266
267 virtual void calculateExplicitQuantities(const Simulator& simulator,
268 const WellStateType& well_state,
269 DeferredLogger& deferred_logger) = 0; // should be const?
270
271 virtual void updateProductivityIndex(const Simulator& simulator,
272 const WellProdIndexCalculator<Scalar>& wellPICalc,
273 WellStateType& well_state,
274 DeferredLogger& deferred_logger) const = 0;
275
276 // Add well contributions to matrix
278
280 const BVector& x,
281 const int pressureVarIndex,
282 const bool use_well_weights,
283 const WellStateType& 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 */ WellStateType& well_state,
294 const GroupState<Scalar>& group_state,
295 WellTestState& welltest_state,
296 GLiftEclWells& ecl_well_map,
297 std::map<std::string, double>& open_times,
298 DeferredLogger& deferred_logger);
299
300 void checkWellOperability(const Simulator& simulator,
301 const WellStateType& well_state,
302 DeferredLogger& deferred_logger);
303
305 WellStateType& well_state,
306 const GroupState<Scalar>& group_state,
307 GLiftEclWells& ecl_well_map,
308 DeferredLogger& deferred_logger);
309
310 // check whether the well is operable under the current reservoir condition
311 // mostly related to BHP limit and THP limit
312 void updateWellOperability(const Simulator& simulator,
313 const WellStateType& well_state,
314 DeferredLogger& deferred_logger);
315
316 bool updateWellOperabilityFromWellEq(const Simulator& simulator,
317 const WellStateType& well_state,
318 DeferredLogger& deferred_logger);
319
320 // update perforation water throughput based on solved water rate
321 virtual void updateWaterThroughput(const double dt,
322 WellStateType& well_state) const = 0;
323
326 virtual std::vector<Scalar>
328 DeferredLogger& deferred_logger) const = 0;
329
333 void initializeProducerWellState(const Simulator& simulator,
334 WellStateType& well_state,
335 DeferredLogger& deferred_logger) const;
336
337 void solveWellEquation(const Simulator& simulator,
338 WellStateType& well_state,
339 const GroupState<Scalar>& group_state,
340 DeferredLogger& deferred_logger);
341
342 const std::vector<RateVector>& connectionRates() const
343 {
344 return connectionRates_;
345 }
346
347 std::vector<Scalar> wellIndex(const int perf,
348 const IntensiveQuantities& intQuants,
349 const Scalar trans_mult,
350 const SingleWellStateType& ws) const;
351
352 void updateConnectionDFactor(const Simulator& simulator,
353 SingleWellStateType& ws) const;
354
356 SingleWellStateType& ws) const;
357
358 virtual bool iterateWellEqWithSwitching(const Simulator& simulator,
359 const double dt,
360 const WellInjectionControls& inj_controls,
361 const WellProductionControls& prod_controls,
362 WellStateType& well_state,
363 const GroupState<Scalar>& group_state,
364 DeferredLogger& deferred_logger,
365 const bool fixed_control = false,
366 const bool fixed_status = false) = 0;
367protected:
368 // simulation parameters
369 std::vector<RateVector> connectionRates_;
370 std::vector<Scalar> B_avg_;
374
375 Scalar wpolymer() const;
376 Scalar wfoam() const;
377 Scalar wsalt() const;
378 Scalar wmicrobes() const;
379 Scalar woxygen() const;
380 Scalar wurea() const;
381
382 virtual Scalar getRefDensity() const = 0;
383
384 std::vector<Scalar>
385 initialWellRateFractions(const Simulator& ebosSimulator,
386 const WellStateType& well_state) const;
387
388 // check whether the well is operable under BHP limit with current reservoir condition
389 virtual void checkOperabilityUnderBHPLimit(const WellStateType& well_state,
390 const Simulator& simulator,
391 DeferredLogger& deferred_logger) = 0;
392
393 // check whether the well is operable under THP limit with current reservoir condition
394 virtual void checkOperabilityUnderTHPLimit(const Simulator& simulator,
395 const WellStateType& well_state,
396 DeferredLogger& deferred_logger) = 0;
397
398 virtual void updateIPR(const Simulator& simulator,
399 DeferredLogger& deferred_logger) const = 0;
400
401 virtual void assembleWellEqWithoutIteration(const Simulator& simulator,
402 const double dt,
403 const WellInjectionControls& inj_controls,
404 const WellProductionControls& prod_controls,
405 WellStateType& well_state,
406 const GroupState<Scalar>& group_state,
407 DeferredLogger& deferred_logger) = 0;
408
409 // iterate well equations with the specified control until converged
410 virtual bool iterateWellEqWithControl(const Simulator& simulator,
411 const double dt,
412 const WellInjectionControls& inj_controls,
413 const WellProductionControls& prod_controls,
414 WellStateType& well_state,
415 const GroupState<Scalar>& group_state,
416 DeferredLogger& deferred_logger) = 0;
417
418 virtual void updateIPRImplicit(const Simulator& simulator,
419 WellStateType& well_state,
420 DeferredLogger& deferred_logger) = 0;
421
422 bool iterateWellEquations(const Simulator& simulator,
423 const double dt,
424 WellStateType& well_state,
425 const GroupState<Scalar>& group_state,
426 DeferredLogger& deferred_logger);
427
428 bool solveWellWithOperabilityCheck(const Simulator& simulator,
429 const double dt,
430 const Well::InjectionControls& inj_controls,
431 const Well::ProductionControls& prod_controls,
432 WellStateType& well_state,
433 const GroupState<Scalar>& group_state,
434 DeferredLogger& deferred_logger);
435
436 std::optional<Scalar>
437 estimateOperableBhp(const Simulator& ebos_simulator,
438 const double dt,
439 WellStateType& well_state,
440 const SummaryState& summary_state,
441 DeferredLogger& deferred_logger);
442
443 bool solveWellWithBhp(const Simulator& simulator,
444 const double dt,
445 const Scalar bhp,
446 WellStateType& well_state,
447 DeferredLogger& deferred_logger);
448
449 bool solveWellWithZeroRate(const Simulator& simulator,
450 const double dt,
451 WellStateType& well_state,
452 DeferredLogger& deferred_logger);
453
454 bool solveWellForTesting(const Simulator& simulator,
455 WellStateType& well_state,
456 const GroupState<Scalar>& group_state,
457 DeferredLogger& deferred_logger);
458
459
460 template<class GasLiftSingleWell>
461 std::unique_ptr<GasLiftSingleWell> initializeGliftWellTest_(const Simulator& simulator,
462 WellStateType& well_state,
463 const GroupState<Scalar>& group_state,
464 GLiftEclWells& ecl_well_map,
465 DeferredLogger& deferred_logger);
466
467 Eval getPerfCellPressure(const FluidState& fs) const;
468
469 // get the mobility for specific perforation
470 template<class Value, class Callback>
471 void getMobility(const Simulator& simulator,
472 const int local_perf_index,
473 std::vector<Value>& mob,
474 Callback& extendEval,
475 [[maybe_unused]] DeferredLogger& deferred_logger) const;
476
477 void computeConnLevelProdInd(const FluidState& fs,
478 const std::function<Scalar(const Scalar)>& connPICalc,
479 const std::vector<Scalar>& mobility,
480 Scalar* connPI) const;
481
482 void computeConnLevelInjInd(const FluidState& fs,
483 const Phase preferred_phase,
484 const std::function<Scalar(const Scalar)>& connIICalc,
485 const std::vector<Scalar>& mobility,
486 Scalar* connII,
487 DeferredLogger& deferred_logger) const;
488
489 Scalar computeConnectionDFactor(const int perf,
490 const IntensiveQuantities& intQuants,
491 const SingleWellStateType& ws) const;
492};
493
494} // namespace Opm
495
496#include "WellInterface_impl.hpp"
497
498#endif // OPM_WELLINTERFACE_HEADER_INCLUDED
Declares the properties required by the black oil model.
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
std::map< std::string, std::pair< const Well *, int > > GLiftEclWells
Definition: GasLiftGroupInfo.hpp:64
Definition: GasLiftSingleWell.hpp:39
Definition: GroupState.hpp:41
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:198
Definition: RateConverter.hpp:70
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:84
Definition: SingleWellState.hpp:43
Definition: WellInterfaceFluidSystem.hpp:49
static constexpr int Oil
Definition: WellInterfaceFluidSystem.hpp:62
static constexpr int Water
Definition: WellInterfaceFluidSystem.hpp:61
static constexpr int Gas
Definition: WellInterfaceFluidSystem.hpp:63
int pvtRegionIdx() const
Definition: WellInterfaceGeneric.hpp:128
const std::vector< FluidSystem::Scalar > & wellIndex() const
Definition: WellInterfaceGeneric.hpp:150
Definition: WellInterfaceIndices.hpp:34
typename WellInterfaceFluidSystem< GetPropType< TypeTag, Properties::FluidSystem > >::ModelParameters ModelParameters
Definition: WellInterfaceIndices.hpp:41
DenseAd::Evaluation< Scalar, Indices::numEq > Eval
Definition: WellInterfaceIndices.hpp:40
Definition: WellInterface.hpp:76
bool solveWellWithBhp(const Simulator &simulator, const double dt, const Scalar bhp, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:701
virtual std::vector< Scalar > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const =0
void checkWellOperability(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1039
Scalar woxygen() const
Definition: WellInterface_impl.hpp:165
static constexpr bool has_brine
Definition: WellInterface.hpp:119
IndividualOrGroup
Definition: WellInterface.hpp:246
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:81
Scalar computeConnectionDFactor(const int perf, const IntensiveQuantities &intQuants, const SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:1874
void solveWellEquation(const Simulator &simulator, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:805
virtual void computeWellRatesWithBhp(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const =0
bool updateWellControlAndStatusLocalIteration(const Simulator &simulator, WellStateType &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:275
virtual void addWellContributions(SparseMatrixAdapter &) const =0
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:103
bool iterateWellEquations(const Simulator &simulator, const double dt, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:503
const std::vector< RateVector > & connectionRates() const
Definition: WellInterface.hpp:342
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:136
Scalar wfoam() const
Definition: WellInterface_impl.hpp:127
virtual void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const =0
std::optional< Scalar > estimateOperableBhp(const Simulator &ebos_simulator, const double dt, WellStateType &well_state, const SummaryState &summary_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:668
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:2069
void addCellRates(RateVector &rates, int cellIdx) const
Definition: WellInterface_impl.hpp:1003
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:97
void gliftBeginTimeStepWellTestUpdateALQ(const Simulator &simulator, WellStateType &well_state, const GroupState< Scalar > &group_state, GLiftEclWells &ecl_well_map, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1066
void assembleWellEq(const Simulator &simulator, const double dt, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:867
Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const
Definition: WellInterface_impl.hpp:1019
virtual void init(const std::vector< Scalar > &depth_arg, const Scalar gravity_arg, const std::vector< Scalar > &B_avg, const bool changed_to_open_this_step)
Definition: WellInterface_impl.hpp:95
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1960
void assembleWellEqWithoutIteration(const Simulator &simulator, const double dt, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:883
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
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:86
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:82
virtual void checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger)=0
std::vector< Scalar > initialWellRateFractions(const Simulator &ebosSimulator, const WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1636
bool updateWellStateWithTHPTargetProd(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:2038
virtual bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const WellInjectionControls &inj_controls, const WellProductionControls &prod_controls, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)=0
static constexpr bool has_watVapor
Definition: WellInterface.hpp:120
void updateConnectionDFactor(const Simulator &simulator, SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:1854
virtual void updateWaterThroughput(const double dt, WellStateType &well_state) const =0
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: WellInterface.hpp:87
virtual Scalar getRefDensity() const =0
typename FluidSystem::IndexTraitsType IndexTraits
Definition: WellInterface.hpp:84
virtual void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger)=0
Eval getPerfCellPressure(const FluidState &fs) const
Definition: WellInterface_impl.hpp:1945
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:96
void initializeProducerWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1686
std::vector< Scalar > B_avg_
Definition: WellInterface.hpp:370
static constexpr bool has_bioeffects
Definition: WellInterface.hpp:123
bool changed_to_stopped_this_step_
Definition: WellInterface.hpp:371
static constexpr bool has_polymermw
Definition: WellInterface.hpp:117
virtual void calculateExplicitQuantities(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger)=0
static constexpr bool has_polymer
Definition: WellInterface.hpp:113
virtual void apply(BVector &r) const =0
r = r - C D^-1 Rw
GetPropType< TypeTag, Properties::Grid > Grid
Definition: WellInterface.hpp:80
bool wellUnderZeroGroupRateTarget(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger, std::optional< bool > group_control=std::nullopt) const
Definition: WellInterface_impl.hpp:1604
virtual ~WellInterface()=default
Virtual destructor.
virtual std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &ebos_simulator, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger, bool iterate_if_no_solution) const =0
typename Base::ModelParameters ModelParameters
Definition: WellInterface.hpp:109
virtual void updatePrimaryVariables(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger)=0
bool updateWellControl(const Simulator &simulator, const IndividualOrGroup iog, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:189
virtual void computeWellRatesWithBhpIterations(const Simulator &simulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const =0
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:83
static constexpr bool has_solvent
Definition: WellInterface.hpp:111
virtual void assembleWellEqWithoutIteration(const Simulator &simulator, const double dt, const WellInjectionControls &inj_controls, const WellProductionControls &prod_controls, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)=0
static constexpr bool has_disgas_in_water
Definition: WellInterface.hpp:121
static constexpr bool has_saltPrecip
Definition: WellInterface.hpp:122
virtual void checkOperabilityUnderBHPLimit(const WellStateType &well_state, const Simulator &simulator, DeferredLogger &deferred_logger)=0
static constexpr bool has_foam
Definition: WellInterface.hpp:118
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: WellInterface.hpp:89
virtual void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, WellStateType &well_state, DeferredLogger &deferred_logger)=0
void updateConnectionTransmissibilityFactor(const Simulator &simulator, SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:1916
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:2103
typename GasLiftGroupInfo< Scalar, IndexTraits >::GLiftEclWells GLiftEclWells
Definition: WellInterface.hpp:91
std::unique_ptr< GasLiftSingleWell > initializeGliftWellTest_(const Simulator &simulator, WellStateType &well_state, const GroupState< Scalar > &group_state, GLiftEclWells &ecl_well_map, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:2137
virtual void addWellPressureEquations(PressureMatrix &mat, const BVector &x, const int pressureVarIndex, const bool use_well_weights, const WellStateType &well_state) const =0
virtual void updateIPRImplicit(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger)=0
int number_of_well_reopenings_
Definition: WellInterface.hpp:373
virtual void solveEqAndUpdateWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger)=0
void wellTesting(const Simulator &simulator, const double simulation_time, WellStateType &well_state, const GroupState< Scalar > &group_state, WellTestState &welltest_state, GLiftEclWells &ecl_well_map, std::map< std::string, double > &open_times, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:383
Scalar wsalt() const
Definition: WellInterface_impl.hpp:141
static constexpr bool has_micp
Definition: WellInterface.hpp:124
Dune::FieldMatrix< Scalar, Indices::numEq, Indices::numEq > MatrixBlockType
Definition: WellInterface.hpp:94
void prepareWellBeforeAssembling(const Simulator &simulator, const double dt, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:903
typename Base::Eval Eval
Definition: WellInterface.hpp:95
bool wellUnderZeroRateTarget(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1586
virtual ConvergenceReport getWellConvergence(const Simulator &simulator, const WellStateType &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_conservation_quantities, const int num_phases, const int index_of_well, const std::vector< PerforationData< Scalar > > &perf_data)
Constructor.
Definition: WellInterface_impl.hpp:58
WellState< Scalar, IndexTraits > WellStateType
Definition: WellInterface.hpp:99
bool solveWellForTesting(const Simulator &simulator, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:759
void updateWellOperability(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1128
static constexpr bool has_energy
Definition: WellInterface.hpp:114
bool solveWellWithOperabilityCheck(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:566
Scalar wpolymer() const
Definition: WellInterface_impl.hpp:111
SingleWellState< Scalar, IndexTraits > SingleWellStateType
Definition: WellInterface.hpp:100
virtual void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const =0
bool updateWellOperabilityFromWellEq(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1165
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:85
Scalar wurea() const
Definition: WellInterface_impl.hpp:177
virtual bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const WellInjectionControls &inj_controls, const WellProductionControls &prod_controls, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger, const bool fixed_control=false, const bool fixed_status=false)=0
bool solveWellWithZeroRate(const Simulator &simulator, const double dt, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:738
bool thp_update_iterations
Definition: WellInterface.hpp:372
Dune::FieldVector< Scalar, Indices::numEq > VectorBlockType
Definition: WellInterface.hpp:93
Scalar wmicrobes() const
Definition: WellInterface_impl.hpp:153
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:88
virtual void scaleSegmentRatesAndPressure(WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1184
virtual void updateWellStateWithTarget(const Simulator &simulator, const GroupState< Scalar > &group_state, WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1192
static constexpr bool has_zFraction
Definition: WellInterface.hpp:112
bool stoppedOrZeroRateTarget(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1623
Definition: WellProdIndexCalculator.hpp:37
Definition: WellState.hpp:66
Declares the properties required by the black oil model.
Definition: blackoilbioeffectsmodules.hh:43
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
Static data associated with a well perforation.
Definition: PerforationData.hpp:30