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;
72template<typename Scalar, typename IndexTraits> class GroupStateHelper;
73
74template<typename TypeTag>
75class WellInterface : public WellInterfaceIndices<GetPropType<TypeTag, Properties::FluidSystem>,
76 GetPropType<TypeTag, Properties::Indices>>
77{
80public:
85 using IndexTraits = typename FluidSystem::IndexTraitsType;
93
94 using VectorBlockType = Dune::FieldVector<Scalar, Indices::numEq>;
95 using MatrixBlockType = Dune::FieldMatrix<Scalar, Indices::numEq, Indices::numEq>;
96 using Eval = typename Base::Eval;
97 using BVector = Dune::BlockVector<VectorBlockType>;
98 using PressureMatrix = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, 1, 1>>;
99
103
106
110
112
113 static constexpr bool has_solvent = getPropValue<TypeTag, Properties::EnableSolvent>();
114 static constexpr bool has_zFraction = getPropValue<TypeTag, Properties::EnableExtbo>();
115 static constexpr bool has_polymer = getPropValue<TypeTag, Properties::EnablePolymer>();
116 static constexpr EnergyModules energyModuleType = getPropValue<TypeTag, Properties::EnergyModuleType>();
117 static constexpr bool has_energy = energyModuleType == EnergyModules::FullyImplicitThermal;
118
119 // flag for polymer molecular weight related
120 static constexpr bool has_polymermw = getPropValue<TypeTag, Properties::EnablePolymerMW>();
121 static constexpr bool has_foam = getPropValue<TypeTag, Properties::EnableFoam>();
122 static constexpr bool has_brine = getPropValue<TypeTag, Properties::EnableBrine>();
123 static constexpr bool has_watVapor = getPropValue<TypeTag, Properties::EnableVapwat>();
124 static constexpr bool has_disgas_in_water = getPropValue<TypeTag, Properties::EnableDisgasInWater>();
125 static constexpr bool has_saltPrecip = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
126 static constexpr bool has_bioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
127 static constexpr bool has_micp = Indices::enableMICP;
128
129 // For the conversion between the surface volume rate and reservoir voidage rate
130 using FluidState = BlackOilFluidState<Eval,
132 energyModuleType != EnergyModules::NoTemperature,
133 energyModuleType == EnergyModules::FullyImplicitThermal,
134 Indices::compositionSwitchIdx >= 0,
136 has_brine,
140 Indices::numPhases >;
142 WellInterface(const Well& well,
143 const ParallelWellInfo<Scalar>& pw_info,
144 const int time_step,
145 const ModelParameters& param,
146 const RateConverterType& rate_converter,
147 const int pvtRegionIdx,
148 const int num_conservation_quantities,
149 const int num_phases,
150 const int index_of_well,
151 const std::vector<PerforationData<Scalar>>& perf_data);
152
154 virtual ~WellInterface() = default;
155
156 virtual void init(const std::vector<Scalar>& depth_arg,
157 const Scalar gravity_arg,
158 const std::vector<Scalar>& B_avg,
159 const bool changed_to_open_this_step);
160
162 const std::vector<Scalar>& B_avg,
163 const bool relax_tolerance) const = 0;
164
165 virtual void solveEqAndUpdateWellState(const Simulator& simulator,
166 const GroupStateHelperType& groupStateHelper,
167 WellStateType& well_state) = 0;
168
169 void assembleWellEq(const Simulator& simulator,
170 const double dt,
171 const GroupStateHelperType& groupStateHelper,
172 WellStateType& well_state);
173
174 void assembleWellEqWithoutIteration(const Simulator& simulator,
175 const GroupStateHelperType& groupStateHelper,
176 const double dt,
177 WellStateType& well_state,
178 const bool solving_with_zero_rate);
179
180 // Convenience overload that gets scaled fractions internally
182 DeferredLogger& deferred_logger) const;
183
184 // TODO: better name or further refactoring the function to make it more clear
185 void prepareWellBeforeAssembling(const Simulator& simulator,
186 const double dt,
187 const GroupStateHelperType& groupStateHelper,
188 WellStateType& well_state);
189
190 virtual void computeWellRatesWithBhp(const Simulator& ebosSimulator,
191 const Scalar& bhp,
192 std::vector<Scalar>& well_flux,
193 DeferredLogger& deferred_logger) const = 0;
194
195 virtual std::optional<Scalar>
197 const GroupStateHelperType& groupStateHelper,
198 const SummaryState& summary_state,
199 const Scalar alq_value,
200 bool iterate_if_no_solution) const = 0;
201
202 std::optional<Scalar>
204 const WellStateType& well_state,
205 Scalar bhp,
206 const SummaryState& summary_state,
207 const Scalar alq_value);
211 const BVector& x,
212 const GroupStateHelperType& groupStateHelper,
213 WellStateType& well_state) = 0;
214
216 virtual void apply(const BVector& x, BVector& Ax) const = 0;
217
219 virtual void apply(BVector& r) const = 0;
220
221 // TODO: before we decide to put more information under mutable, this function is not const
222 virtual void computeWellPotentials(const Simulator& simulator,
223 const WellStateType& well_state,
224 const GroupStateHelperType& groupStateHelper,
225 std::vector<Scalar>& well_potentials) = 0;
226
227 virtual void updateWellStateWithTarget(const Simulator& simulator,
228 const GroupStateHelperType& groupStateHelper,
229 WellStateType& well_state) const;
230
231 virtual void scaleSegmentRatesAndPressure(WellStateType& well_state) const;
232
233 virtual void computeWellRatesWithBhpIterations(const Simulator& simulator,
234 const Scalar& bhp,
235 const GroupStateHelperType& groupStateHelper,
236 std::vector<Scalar>& well_flux) const = 0;
237
238 bool wellUnderZeroRateTarget(const GroupStateHelperType& groupStateHelper) const;
239
240 bool stoppedOrZeroRateTarget(const GroupStateHelperType& groupStateHelper) const;
241
242 bool updateWellStateWithTHPTargetProd(const Simulator& simulator,
243 WellStateType& well_state,
244 const GroupStateHelperType& groupStateHelper) const;
245
246 bool wellUnderZeroGroupRateTarget(const GroupStateHelperType& groupStateHelper,
247 const std::optional<bool> group_control = std::nullopt) const;
248
249 enum class IndividualOrGroup { Individual, Group, Both };
250 bool updateWellControl(const Simulator& simulator,
251 const IndividualOrGroup iog,
252 const GroupStateHelperType& groupStateHelper,
253 WellStateType& well_state) /* const */;
254
256 const GroupStateHelperType& groupStateHelper,
257 const Well::InjectionControls& inj_controls,
258 const Well::ProductionControls& prod_controls,
259 const Scalar WQTotal,
260 WellStateType& well_state,
261 const bool fixed_control,
262 const bool fixed_status,
263 const bool solving_with_zero_rate);
264
265 virtual void updatePrimaryVariables(const GroupStateHelperType& groupStateHelper) = 0;
266
267 virtual void calculateExplicitQuantities(const Simulator& simulator,
268 const GroupStateHelperType& groupStateHelper) = 0; // should be const?
269
270 virtual void updateProductivityIndex(const Simulator& simulator,
271 const WellProdIndexCalculator<Scalar>& wellPICalc,
272 WellStateType& well_state,
273 DeferredLogger& deferred_logger) const = 0;
274
275 // Add well contributions to matrix
277
279 const BVector& x,
280 const int pressureVarIndex,
281 const bool use_well_weights,
282 const WellStateType& well_state) const = 0;
283
284 void addCellRates(std::map<int, RateVector>& cellRates_) const;
285
286 Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const;
287
288 // TODO: theoretically, it should be a const function
289 // Simulator is not const is because that assembleWellEq is non-const Simulator
290 void wellTesting(const Simulator& simulator,
291 const double simulation_time,
292 const GroupStateHelperType& groupStateHelper,
293 WellStateType& well_state,
294 WellTestState& welltest_state,
295 GLiftEclWells& ecl_well_map,
296 std::map<std::string, double>& open_times);
297
298 void checkWellOperability(const Simulator& simulator,
299 const WellStateType& well_state,
300 const GroupStateHelperType& groupStateHelper);
301
303 WellStateType& well_state,
304 const GroupState<Scalar>& group_state,
305 GLiftEclWells& ecl_well_map,
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 WellStateType& well_state,
312 const GroupStateHelperType& groupStateHelper);
313
314 bool updateWellOperabilityFromWellEq(const Simulator& simulator,
315 const GroupStateHelperType& groupStateHelper);
316
317 // update perforation water throughput based on solved water rate
318 virtual void updateWaterThroughput(const double dt,
319 WellStateType& well_state) const = 0;
320
323 virtual std::vector<Scalar>
325 DeferredLogger& deferred_logger) const = 0;
326
329 virtual void getScaledWellFractions(std::vector<Scalar>& scaled_fractions,
330 DeferredLogger& deferred_logger) const = 0;
331
335 void initializeProducerWellState(const Simulator& simulator,
336 WellStateType& well_state,
337 DeferredLogger& deferred_logger) const;
338
339 void solveWellEquation(const Simulator& simulator,
340 const GroupStateHelperType& groupStateHelper,
341 WellStateType& well_state);
342
343 const std::vector<RateVector>& connectionRates() const
344 {
345 return connectionRates_;
346 }
347
348 void updateConnectionDFactor(const Simulator& simulator,
349 SingleWellStateType& ws) const;
350
352 SingleWellStateType& ws) const;
353
354 virtual bool iterateWellEqWithSwitching(const Simulator& simulator,
355 const double dt,
356 const WellInjectionControls& inj_controls,
357 const WellProductionControls& prod_controls,
358 const GroupStateHelperType& groupStateHelper,
359 WellStateType& well_state,
360 const bool fixed_control,
361 const bool fixed_status,
362 const bool solving_with_zero_rate) = 0;
363
364 virtual Scalar maxPerfPress(const Simulator& simulator) const = 0;
365
366 virtual void updateIPRImplicit(const Simulator& simulator,
367 const GroupStateHelperType& groupStateHelper,
368 WellStateType& well_state) = 0;
369
370protected:
371 // simulation parameters
372 std::vector<RateVector> connectionRates_;
373 std::vector<Scalar> B_avg_;
377
378 Scalar wpolymer() const;
379 Scalar wfoam() const;
380 Scalar wsalt() const;
381 Scalar wmicrobes() const;
382 Scalar woxygen() const;
383 Scalar wurea() const;
384
385 virtual Scalar getRefDensity() const = 0;
386
387 std::vector<Scalar>
388 initialWellRateFractions(const Simulator& ebosSimulator,
389 const WellStateType& well_state) const;
390
391 // check whether the well is operable under BHP limit with current reservoir condition
392 virtual void checkOperabilityUnderBHPLimit(const WellStateType& well_state,
393 const Simulator& simulator,
394 DeferredLogger& deferred_logger) = 0;
395
396 // check whether the well is operable under THP limit with current reservoir condition
397 virtual void checkOperabilityUnderTHPLimit(const Simulator& simulator,
398 const WellStateType& well_state,
399 const GroupStateHelperType& groupStateHelper) = 0;
400
401 virtual void updateIPR(const Simulator& simulator,
402 DeferredLogger& deferred_logger) const = 0;
403
404 virtual void assembleWellEqWithoutIteration(const Simulator& simulator,
405 const GroupStateHelperType& groupStateHelper,
406 const double dt,
407 const WellInjectionControls& inj_controls,
408 const WellProductionControls& prod_controls,
409 WellStateType& well_state,
410 const bool solving_with_zero_rate) = 0;
411
412 // iterate well equations with the specified control until converged
413 virtual bool iterateWellEqWithControl(const Simulator& simulator,
414 const double dt,
415 const WellInjectionControls& inj_controls,
416 const WellProductionControls& prod_controls,
417 const GroupStateHelperType& groupStateHelper,
418 WellStateType& well_state) = 0;
419
420 bool iterateWellEquations(const Simulator& simulator,
421 const double dt,
422 const GroupStateHelperType& groupStateHelper,
423 WellStateType& well_state);
424
425 bool solveWellWithOperabilityCheck(const Simulator& simulator,
426 const double dt,
427 const Well::InjectionControls& inj_controls,
428 const Well::ProductionControls& prod_controls,
429 const GroupStateHelperType& groupStateHelper,
430 WellStateType& well_state);
431
432 std::optional<Scalar>
433 estimateOperableBhp(const Simulator& ebos_simulator,
434 const double dt,
435 const GroupStateHelperType& groupStateHelper,
436 const SummaryState& summary_state,
437 WellStateType& well_state);
438
439 bool solveWellWithBhp(const Simulator& simulator,
440 const double dt,
441 const Scalar bhp,
442 const GroupStateHelperType& groupStateHelper,
443 WellStateType& well_state);
444
445 bool solveWellWithZeroRate(const Simulator& simulator,
446 const double dt,
447 const GroupStateHelperType& groupStateHelper,
448 WellStateType& well_state);
449
450 bool solveWellForTesting(const Simulator& simulator,
451 const GroupStateHelperType& groupStateHelper,
452 WellStateType& well_state);
453
454
455 template<class GasLiftSingleWell>
456 std::unique_ptr<GasLiftSingleWell> initializeGliftWellTest_(const Simulator& simulator,
457 WellStateType& well_state,
458 const GroupState<Scalar>& group_state,
459 GLiftEclWells& ecl_well_map,
460 DeferredLogger& deferred_logger);
461
462 Eval getPerfCellPressure(const FluidState& fs) const;
463
464 // get the transmissibility multiplier for specific perforation
465 template<class Value, class Callback>
466 void getTransMult(Value& trans_mult,
467 const Simulator& simulator,
468 const int cell_idx,
469 Callback& extendEval) const;
470
471 // get the well transmissibility for specific perforation
472 template<class Value>
473 void getTw(std::vector<Value>& wi,
474 const int perf,
475 const IntensiveQuantities& intQuants,
476 const Value& trans_mult,
477 const SingleWellStateType& ws) const;
478
479 // get the mobility for specific perforation
480 template<class Value, class Callback>
481 void getMobility(const Simulator& simulator,
482 const int local_perf_index,
483 std::vector<Value>& mob,
484 Callback& extendEval,
485 [[maybe_unused]] DeferredLogger& deferred_logger) const;
486
487 void computeConnLevelProdInd(const FluidState& fs,
488 const std::function<Scalar(const Scalar)>& connPICalc,
489 const std::vector<Scalar>& mobility,
490 Scalar* connPI) const;
491
492 void computeConnLevelInjInd(const FluidState& fs,
493 const Phase preferred_phase,
494 const std::function<Scalar(const Scalar)>& connIICalc,
495 const std::vector<Scalar>& mobility,
496 Scalar* connII,
497 DeferredLogger& deferred_logger) const;
498
499 Scalar computeConnectionDFactor(const int perf,
500 const IntensiveQuantities& intQuants,
501 const SingleWellStateType& ws) const;
502};
503
504} // namespace Opm
505
506#include "WellInterface_impl.hpp"
507
508#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:65
Definition: GasLiftSingleWell.hpp:39
Definition: GroupStateHelper.hpp:56
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:44
Definition: WellInterfaceFluidSystem.hpp:50
static constexpr int Oil
Definition: WellInterfaceFluidSystem.hpp:65
static constexpr int Water
Definition: WellInterfaceFluidSystem.hpp:64
static constexpr int Gas
Definition: WellInterfaceFluidSystem.hpp:66
int pvtRegionIdx() const
Definition: WellInterfaceGeneric.hpp:130
Definition: WellInterfaceIndices.hpp:34
DenseAd::Evaluation< Scalar, Indices::numDerivatives > Eval
Definition: WellInterfaceIndices.hpp:40
typename WellInterfaceFluidSystem< GetPropType< TypeTag, Properties::FluidSystem > >::ModelParameters ModelParameters
Definition: WellInterfaceIndices.hpp:41
Definition: WellInterface.hpp:77
bool stoppedOrZeroRateTarget(const GroupStateHelperType &groupStateHelper) const
Definition: WellInterface_impl.hpp:1761
bool updateWellOperabilityFromWellEq(const Simulator &simulator, const GroupStateHelperType &groupStateHelper)
Definition: WellInterface_impl.hpp:1312
virtual void checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper)=0
virtual std::vector< Scalar > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const =0
void checkWellOperability(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper)
Definition: WellInterface_impl.hpp:1184
void updateWellOperability(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper)
Definition: WellInterface_impl.hpp:1274
bool solveWellWithOperabilityCheck(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:610
Scalar woxygen() const
Definition: WellInterface_impl.hpp:166
virtual Scalar maxPerfPress(const Simulator &simulator) const =0
static constexpr bool has_brine
Definition: WellInterface.hpp:122
IndividualOrGroup
Definition: WellInterface.hpp:249
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:82
void assembleWellEqWithoutIteration(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const double dt, WellStateType &well_state, const bool solving_with_zero_rate)
Definition: WellInterface_impl.hpp:983
virtual void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)=0
Scalar computeConnectionDFactor(const int perf, const IntensiveQuantities &intQuants, const SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:2008
virtual std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, const Scalar alq_value, bool iterate_if_no_solution) const =0
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:105
const std::vector< RateVector > & connectionRates() const
Definition: WellInterface.hpp:343
Scalar wfoam() const
Definition: WellInterface_impl.hpp:128
virtual void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const =0
bool updateWellControlAndStatusLocalIteration(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const Scalar WQTotal, WellStateType &well_state, const bool fixed_control, const bool fixed_status, const bool solving_with_zero_rate)
Definition: WellInterface_impl.hpp:287
void getTransMult(Value &trans_mult, const Simulator &simulator, const int cell_idx, Callback &extendEval) const
Definition: WellInterface_impl.hpp:2096
std::vector< RateVector > connectionRates_
Definition: WellInterface.hpp:372
bool solveWellForTesting(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:849
void computeConnLevelProdInd(const FluidState &fs, const std::function< Scalar(const Scalar)> &connPICalc, const std::vector< Scalar > &mobility, Scalar *connPI) const
Definition: WellInterface_impl.hpp:2262
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:98
void gliftBeginTimeStepWellTestUpdateALQ(const Simulator &simulator, WellStateType &well_state, const GroupState< Scalar > &group_state, GLiftEclWells &ecl_well_map, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1212
Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const
Definition: WellInterface_impl.hpp:1164
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:96
std::optional< Scalar > computeBhpAtThpLimitProdWithAlqUsingIPR(const Simulator &simulator, const WellStateType &well_state, Scalar bhp, const SummaryState &summary_state, const Scalar alq_value)
Definition: WellInterface_impl.hpp:2227
virtual void updatePrimaryVariables(const GroupStateHelperType &groupStateHelper)=0
void getTw(std::vector< Value > &wi, const int perf, const IntensiveQuantities &intQuants, const Value &trans_mult, const SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:1910
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:2109
virtual void apply(const BVector &x, BVector &Ax) const =0
Ax = Ax - C D^-1 B x.
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:87
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
virtual void updateIPRImplicit(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)=0
std::vector< Scalar > initialWellRateFractions(const Simulator &ebosSimulator, const WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1772
void solveWellEquation(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:905
static constexpr bool has_watVapor
Definition: WellInterface.hpp:123
void updateConnectionDFactor(const Simulator &simulator, SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:1988
virtual void updateWaterThroughput(const double dt, WellStateType &well_state) const =0
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: WellInterface.hpp:88
virtual Scalar getRefDensity() const =0
typename FluidSystem::IndexTraitsType IndexTraits
Definition: WellInterface.hpp:85
Eval getPerfCellPressure(const FluidState &fs) const
Definition: WellInterface_impl.hpp:2081
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:97
void initializeProducerWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1822
std::vector< Scalar > B_avg_
Definition: WellInterface.hpp:373
static constexpr bool has_bioeffects
Definition: WellInterface.hpp:126
bool changed_to_stopped_this_step_
Definition: WellInterface.hpp:374
virtual void updateWellStateWithTarget(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1342
static constexpr bool has_polymermw
Definition: WellInterface.hpp:120
void addCellRates(std::map< int, RateVector > &cellRates_) const
Definition: WellInterface_impl.hpp:1145
static constexpr bool has_polymer
Definition: WellInterface.hpp:115
virtual void apply(BVector &r) const =0
r = r - C D^-1 Rw
GetPropType< TypeTag, Properties::Grid > Grid
Definition: WellInterface.hpp:81
virtual ~WellInterface()=default
Virtual destructor.
typename Base::ModelParameters ModelParameters
Definition: WellInterface.hpp:111
virtual void solveEqAndUpdateWellState(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)=0
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:84
static constexpr bool has_solvent
Definition: WellInterface.hpp:113
static constexpr bool has_disgas_in_water
Definition: WellInterface.hpp:124
virtual bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const WellInjectionControls &inj_controls, const WellProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, const bool fixed_control, const bool fixed_status, const bool solving_with_zero_rate)=0
static constexpr bool has_saltPrecip
Definition: WellInterface.hpp:125
bool wellUnderZeroRateTarget(const GroupStateHelperType &groupStateHelper) const
Definition: WellInterface_impl.hpp:1729
virtual void assembleWellEqWithoutIteration(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const double dt, const WellInjectionControls &inj_controls, const WellProductionControls &prod_controls, WellStateType &well_state, const bool solving_with_zero_rate)=0
virtual void checkOperabilityUnderBHPLimit(const WellStateType &well_state, const Simulator &simulator, DeferredLogger &deferred_logger)=0
static constexpr bool has_foam
Definition: WellInterface.hpp:121
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: WellInterface.hpp:90
void updateConnectionTransmissibilityFactor(const Simulator &simulator, SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:2050
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:2296
typename GasLiftGroupInfo< Scalar, IndexTraits >::GLiftEclWells GLiftEclWells
Definition: WellInterface.hpp:92
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:2330
virtual void addWellPressureEquations(PressureMatrix &mat, const BVector &x, const int pressureVarIndex, const bool use_well_weights, const WellStateType &well_state) const =0
std::optional< Scalar > estimateOperableBhp(const Simulator &ebos_simulator, const double dt, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, WellStateType &well_state)
Definition: WellInterface_impl.hpp:730
int number_of_well_reopenings_
Definition: WellInterface.hpp:376
Scalar wsalt() const
Definition: WellInterface_impl.hpp:142
static constexpr bool has_micp
Definition: WellInterface.hpp:127
bool solveWellWithZeroRate(const Simulator &simulator, const double dt, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:817
bool wellUnderZeroGroupRateTarget(const GroupStateHelperType &groupStateHelper, const std::optional< bool > group_control=std::nullopt) const
Definition: WellInterface_impl.hpp:1746
bool solveWellWithBhp(const Simulator &simulator, const double dt, const Scalar bhp, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:767
Dune::FieldMatrix< Scalar, Indices::numEq, Indices::numEq > MatrixBlockType
Definition: WellInterface.hpp:95
void prepareWellBeforeAssembling(const Simulator &simulator, const double dt, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:1018
void wellTesting(const Simulator &simulator, const double simulation_time, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, WellTestState &welltest_state, GLiftEclWells &ecl_well_map, std::map< std::string, double > &open_times)
Definition: WellInterface_impl.hpp:406
virtual ConvergenceReport getWellConvergence(const GroupStateHelperType &groupStateHelper, const std::vector< Scalar > &B_avg, const bool relax_tolerance) const =0
typename Base::Eval Eval
Definition: WellInterface.hpp:96
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:59
WellState< Scalar, IndexTraits > WellStateType
Definition: WellInterface.hpp:100
bool updateWellStateWithTHPTargetProd(const Simulator &simulator, WellStateType &well_state, const GroupStateHelperType &groupStateHelper) const
Definition: WellInterface_impl.hpp:2195
virtual void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_potentials)=0
bool iterateWellEquations(const Simulator &simulator, const double dt, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:532
static constexpr bool has_energy
Definition: WellInterface.hpp:117
Scalar wpolymer() const
Definition: WellInterface_impl.hpp:112
SingleWellState< Scalar, IndexTraits > SingleWellStateType
Definition: WellInterface.hpp:101
virtual bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const WellInjectionControls &inj_controls, const WellProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)=0
virtual void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const =0
virtual void calculateExplicitQuantities(const Simulator &simulator, const GroupStateHelperType &groupStateHelper)=0
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:86
bool updateWellControl(const Simulator &simulator, const IndividualOrGroup iog, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:190
Scalar wurea() const
Definition: WellInterface_impl.hpp:178
void updateGroupTargetFallbackFlag(WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1003
virtual void getScaledWellFractions(std::vector< Scalar > &scaled_fractions, DeferredLogger &deferred_logger) const =0
bool thp_update_iterations
Definition: WellInterface.hpp:375
Dune::FieldVector< Scalar, Indices::numEq > VectorBlockType
Definition: WellInterface.hpp:94
void assembleWellEq(const Simulator &simulator, const double dt, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:967
BlackOilFluidState< Eval, FluidSystem, energyModuleType !=EnergyModules::NoTemperature, energyModuleType==EnergyModules::FullyImplicitThermal, Indices::compositionSwitchIdx >=0, has_watVapor, has_brine, has_saltPrecip, has_disgas_in_water, has_solvent, Indices::numPhases > FluidState
Definition: WellInterface.hpp:140
Scalar wmicrobes() const
Definition: WellInterface_impl.hpp:154
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:89
virtual void scaleSegmentRatesAndPressure(WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1334
GroupStateHelper< Scalar, IndexTraits > GroupStateHelperType
Definition: WellInterface.hpp:102
virtual void computeWellRatesWithBhpIterations(const Simulator &simulator, const Scalar &bhp, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_flux) const =0
static constexpr bool has_zFraction
Definition: WellInterface.hpp:114
static constexpr EnergyModules energyModuleType
Definition: WellInterface.hpp:116
Definition: WellProdIndexCalculator.hpp:37
Definition: WellState.hpp:66
Declares the properties required by the black oil model.
Phase
Phase indices for reservoir coupling, we currently only support black-oil phases (oil,...
Definition: ReservoirCoupling.hpp:156
Definition: blackoilbioeffectsmodules.hh:45
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
Static data associated with a well perforation.
Definition: PerforationData.hpp:30