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::ConstantTemperature,
133 (energyModuleType == EnergyModules::FullyImplicitThermal || energyModuleType == EnergyModules::SequentialImplicitThermal),
134 Indices::compositionSwitchIdx >= 0,
136 has_brine,
139 Indices::numPhases >;
141 WellInterface(const Well& well,
142 const ParallelWellInfo<Scalar>& pw_info,
143 const int time_step,
144 const ModelParameters& param,
145 const RateConverterType& rate_converter,
146 const int pvtRegionIdx,
147 const int num_conservation_quantities,
148 const int num_phases,
149 const int index_of_well,
150 const std::vector<PerforationData<Scalar>>& perf_data);
151
153 virtual ~WellInterface() = default;
154
155 virtual void init(const std::vector<Scalar>& depth_arg,
156 const Scalar gravity_arg,
157 const std::vector<Scalar>& B_avg,
158 const bool changed_to_open_this_step);
159
161 const std::vector<Scalar>& B_avg,
162 const bool relax_tolerance) const = 0;
163
164 virtual void solveEqAndUpdateWellState(const Simulator& simulator,
165 const GroupStateHelperType& groupStateHelper,
166 WellStateType& well_state) = 0;
167
168 void assembleWellEq(const Simulator& simulator,
169 const double dt,
170 const GroupStateHelperType& groupStateHelper,
171 WellStateType& well_state);
172
173 void assembleWellEqWithoutIteration(const Simulator& simulator,
174 const GroupStateHelperType& groupStateHelper,
175 const double dt,
176 WellStateType& well_state,
177 const bool solving_with_zero_rate);
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 const GroupStateHelperType& groupStateHelper,
183 WellStateType& well_state);
184
185 virtual void computeWellRatesWithBhp(const Simulator& ebosSimulator,
186 const Scalar& bhp,
187 std::vector<Scalar>& well_flux,
188 DeferredLogger& deferred_logger) const = 0;
189
190 virtual std::optional<Scalar>
192 const GroupStateHelperType& groupStateHelper,
193 const SummaryState& summary_state,
194 const Scalar alq_value,
195 bool iterate_if_no_solution) const = 0;
196
197 std::optional<Scalar>
199 const WellStateType& well_state,
200 Scalar bhp,
201 const SummaryState& summary_state,
202 const Scalar alq_value);
206 const BVector& x,
207 const GroupStateHelperType& groupStateHelper,
208 WellStateType& well_state) = 0;
209
211 virtual void apply(const BVector& x, BVector& Ax) const = 0;
212
214 virtual void apply(BVector& r) const = 0;
215
216 // TODO: before we decide to put more information under mutable, this function is not const
217 virtual void computeWellPotentials(const Simulator& simulator,
218 const WellStateType& well_state,
219 const GroupStateHelperType& groupStateHelper,
220 std::vector<Scalar>& well_potentials) = 0;
221
222 virtual void updateWellStateWithTarget(const Simulator& simulator,
223 const GroupStateHelperType& groupStateHelper,
224 WellStateType& well_state) const;
225
226 virtual void scaleSegmentRatesAndPressure(WellStateType& well_state) const;
227
228 virtual void computeWellRatesWithBhpIterations(const Simulator& simulator,
229 const Scalar& bhp,
230 const GroupStateHelperType& groupStateHelper,
231 std::vector<Scalar>& well_flux) const = 0;
232
233 bool wellUnderZeroRateTarget(const GroupStateHelperType& groupStateHelper) const;
234
235 bool stoppedOrZeroRateTarget(const GroupStateHelperType& groupStateHelper) const;
236
237 bool updateWellStateWithTHPTargetProd(const Simulator& simulator,
238 WellStateType& well_state,
239 const GroupStateHelperType& groupStateHelper) const;
240
241 bool wellUnderZeroGroupRateTarget(const GroupStateHelperType& groupStateHelper,
242 const std::optional<bool> group_control = std::nullopt) const;
243
244 enum class IndividualOrGroup { Individual, Group, Both };
245 bool updateWellControl(const Simulator& simulator,
246 const IndividualOrGroup iog,
247 const GroupStateHelperType& groupStateHelper,
248 WellStateType& well_state) /* const */;
249
251 const GroupStateHelperType& groupStateHelper,
252 const Well::InjectionControls& inj_controls,
253 const Well::ProductionControls& prod_controls,
254 const Scalar WQTotal,
255 WellStateType& well_state,
256 const bool fixed_control,
257 const bool fixed_status,
258 const bool solving_with_zero_rate);
259
260 virtual void updatePrimaryVariables(const GroupStateHelperType& groupStateHelper) = 0;
261
262 virtual void calculateExplicitQuantities(const Simulator& simulator,
263 const GroupStateHelperType& groupStateHelper) = 0; // should be const?
264
265 virtual void updateProductivityIndex(const Simulator& simulator,
266 const WellProdIndexCalculator<Scalar>& wellPICalc,
267 WellStateType& well_state,
268 DeferredLogger& deferred_logger) const = 0;
269
270 // Add well contributions to matrix
272
274 const BVector& x,
275 const int pressureVarIndex,
276 const bool use_well_weights,
277 const WellStateType& well_state) const = 0;
278
279 void addCellRates(std::map<int, RateVector>& cellRates_) const;
280
281 Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const;
282
283 // TODO: theoretically, it should be a const function
284 // Simulator is not const is because that assembleWellEq is non-const Simulator
285 void wellTesting(const Simulator& simulator,
286 const double simulation_time,
287 const GroupStateHelperType& groupStateHelper,
288 WellStateType& well_state,
289 WellTestState& welltest_state,
290 GLiftEclWells& ecl_well_map,
291 std::map<std::string, double>& open_times);
292
293 void checkWellOperability(const Simulator& simulator,
294 const WellStateType& well_state,
295 const GroupStateHelperType& groupStateHelper);
296
298 WellStateType& well_state,
299 const GroupState<Scalar>& group_state,
300 GLiftEclWells& ecl_well_map,
301 DeferredLogger& deferred_logger);
302
303 // check whether the well is operable under the current reservoir condition
304 // mostly related to BHP limit and THP limit
305 void updateWellOperability(const Simulator& simulator,
306 const WellStateType& well_state,
307 const GroupStateHelperType& groupStateHelper);
308
309 bool updateWellOperabilityFromWellEq(const Simulator& simulator,
310 const GroupStateHelperType& groupStateHelper);
311
312 // update perforation water throughput based on solved water rate
313 virtual void updateWaterThroughput(const double dt,
314 WellStateType& well_state) const = 0;
315
318 virtual std::vector<Scalar>
320 DeferredLogger& deferred_logger) const = 0;
321
325 void initializeProducerWellState(const Simulator& simulator,
326 WellStateType& well_state,
327 DeferredLogger& deferred_logger) const;
328
329 void solveWellEquation(const Simulator& simulator,
330 const GroupStateHelperType& groupStateHelper,
331 WellStateType& well_state);
332
333 const std::vector<RateVector>& connectionRates() const
334 {
335 return connectionRates_;
336 }
337
338 void updateConnectionDFactor(const Simulator& simulator,
339 SingleWellStateType& ws) const;
340
342 SingleWellStateType& ws) const;
343
344 virtual bool iterateWellEqWithSwitching(const Simulator& simulator,
345 const double dt,
346 const WellInjectionControls& inj_controls,
347 const WellProductionControls& prod_controls,
348 const GroupStateHelperType& groupStateHelper,
349 WellStateType& well_state,
350 const bool fixed_control,
351 const bool fixed_status,
352 const bool solving_with_zero_rate) = 0;
353
354 virtual Scalar maxPerfPress(const Simulator& simulator) const = 0;
355
356 virtual void updateIPRImplicit(const Simulator& simulator,
357 const GroupStateHelperType& groupStateHelper,
358 WellStateType& well_state) = 0;
359
360protected:
361 // simulation parameters
362 std::vector<RateVector> connectionRates_;
363 std::vector<Scalar> B_avg_;
367
368 Scalar wpolymer() const;
369 Scalar wfoam() const;
370 Scalar wsalt() const;
371 Scalar wmicrobes() const;
372 Scalar woxygen() const;
373 Scalar wurea() const;
374
375 virtual Scalar getRefDensity() const = 0;
376
377 std::vector<Scalar>
378 initialWellRateFractions(const Simulator& ebosSimulator,
379 const WellStateType& well_state) const;
380
381 // check whether the well is operable under BHP limit with current reservoir condition
382 virtual void checkOperabilityUnderBHPLimit(const WellStateType& well_state,
383 const Simulator& simulator,
384 DeferredLogger& deferred_logger) = 0;
385
386 // check whether the well is operable under THP limit with current reservoir condition
387 virtual void checkOperabilityUnderTHPLimit(const Simulator& simulator,
388 const WellStateType& well_state,
389 const GroupStateHelperType& groupStateHelper) = 0;
390
391 virtual void updateIPR(const Simulator& simulator,
392 DeferredLogger& deferred_logger) const = 0;
393
394 virtual void assembleWellEqWithoutIteration(const Simulator& simulator,
395 const GroupStateHelperType& groupStateHelper,
396 const double dt,
397 const WellInjectionControls& inj_controls,
398 const WellProductionControls& prod_controls,
399 WellStateType& well_state,
400 const bool solving_with_zero_rate) = 0;
401
402 // iterate well equations with the specified control until converged
403 virtual bool iterateWellEqWithControl(const Simulator& simulator,
404 const double dt,
405 const WellInjectionControls& inj_controls,
406 const WellProductionControls& prod_controls,
407 const GroupStateHelperType& groupStateHelper,
408 WellStateType& well_state) = 0;
409
410 bool iterateWellEquations(const Simulator& simulator,
411 const double dt,
412 const GroupStateHelperType& groupStateHelper,
413 WellStateType& well_state);
414
415 bool solveWellWithOperabilityCheck(const Simulator& simulator,
416 const double dt,
417 const Well::InjectionControls& inj_controls,
418 const Well::ProductionControls& prod_controls,
419 const GroupStateHelperType& groupStateHelper,
420 WellStateType& well_state);
421
422 std::optional<Scalar>
423 estimateOperableBhp(const Simulator& ebos_simulator,
424 const double dt,
425 const GroupStateHelperType& groupStateHelper,
426 const SummaryState& summary_state,
427 WellStateType& well_state);
428
429 bool solveWellWithBhp(const Simulator& simulator,
430 const double dt,
431 const Scalar bhp,
432 const GroupStateHelperType& groupStateHelper,
433 WellStateType& well_state);
434
435 bool solveWellWithZeroRate(const Simulator& simulator,
436 const double dt,
437 const GroupStateHelperType& groupStateHelper,
438 WellStateType& well_state);
439
440 bool solveWellForTesting(const Simulator& simulator,
441 const GroupStateHelperType& groupStateHelper,
442 WellStateType& well_state);
443
444
445 template<class GasLiftSingleWell>
446 std::unique_ptr<GasLiftSingleWell> initializeGliftWellTest_(const Simulator& simulator,
447 WellStateType& well_state,
448 const GroupState<Scalar>& group_state,
449 GLiftEclWells& ecl_well_map,
450 DeferredLogger& deferred_logger);
451
452 Eval getPerfCellPressure(const FluidState& fs) const;
453
454 // get the transmissibility multiplier for specific perforation
455 template<class Value, class Callback>
456 void getTransMult(Value& trans_mult,
457 const Simulator& simulator,
458 const int cell_idx,
459 Callback& extendEval) const;
460
461 // get the well transmissibility for specific perforation
462 template<class Value>
463 void getTw(std::vector<Value>& wi,
464 const int perf,
465 const IntensiveQuantities& intQuants,
466 const Value& trans_mult,
467 const SingleWellStateType& ws) 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: GroupStateHelper.hpp:54
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: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:1718
bool updateWellOperabilityFromWellEq(const Simulator &simulator, const GroupStateHelperType &groupStateHelper)
Definition: WellInterface_impl.hpp:1270
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:1142
void updateWellOperability(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper)
Definition: WellInterface_impl.hpp:1232
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:590
Scalar woxygen() const
Definition: WellInterface_impl.hpp:165
BlackOilFluidState< Eval, FluidSystem, energyModuleType==EnergyModules::ConstantTemperature,(energyModuleType==EnergyModules::FullyImplicitThermal||energyModuleType==EnergyModules::SequentialImplicitThermal), Indices::compositionSwitchIdx >=0, has_watVapor, has_brine, has_saltPrecip, has_disgas_in_water, Indices::numPhases > FluidState
Definition: WellInterface.hpp:139
virtual Scalar maxPerfPress(const Simulator &simulator) const =0
static constexpr bool has_brine
Definition: WellInterface.hpp:122
IndividualOrGroup
Definition: WellInterface.hpp:244
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:963
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:1965
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:333
Scalar wfoam() const
Definition: WellInterface_impl.hpp:127
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:277
void getTransMult(Value &trans_mult, const Simulator &simulator, const int cell_idx, Callback &extendEval) const
Definition: WellInterface_impl.hpp:2053
std::vector< RateVector > connectionRates_
Definition: WellInterface.hpp:362
bool solveWellForTesting(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:829
void computeConnLevelProdInd(const FluidState &fs, const std::function< Scalar(const Scalar)> &connPICalc, const std::vector< Scalar > &mobility, Scalar *connPI) const
Definition: WellInterface_impl.hpp:2211
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:1170
Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const
Definition: WellInterface_impl.hpp:1122
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
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:2176
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:1867
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:2066
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:1729
void solveWellEquation(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:885
static constexpr bool has_watVapor
Definition: WellInterface.hpp:123
void updateConnectionDFactor(const Simulator &simulator, SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:1945
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:2038
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:97
void initializeProducerWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1779
std::vector< Scalar > B_avg_
Definition: WellInterface.hpp:363
static constexpr bool has_bioeffects
Definition: WellInterface.hpp:126
bool changed_to_stopped_this_step_
Definition: WellInterface.hpp:364
virtual void updateWellStateWithTarget(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1300
static constexpr bool has_polymermw
Definition: WellInterface.hpp:120
void addCellRates(std::map< int, RateVector > &cellRates_) const
Definition: WellInterface_impl.hpp:1103
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:1686
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:2007
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:2245
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:2279
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:710
int number_of_well_reopenings_
Definition: WellInterface.hpp:366
Scalar wsalt() const
Definition: WellInterface_impl.hpp:141
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:797
bool wellUnderZeroGroupRateTarget(const GroupStateHelperType &groupStateHelper, const std::optional< bool > group_control=std::nullopt) const
Definition: WellInterface_impl.hpp:1703
bool solveWellWithBhp(const Simulator &simulator, const double dt, const Scalar bhp, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:747
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:983
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:395
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:58
WellState< Scalar, IndexTraits > WellStateType
Definition: WellInterface.hpp:100
bool updateWellStateWithTHPTargetProd(const Simulator &simulator, WellStateType &well_state, const GroupStateHelperType &groupStateHelper) const
Definition: WellInterface_impl.hpp:2144
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:521
static constexpr bool has_energy
Definition: WellInterface.hpp:117
Scalar wpolymer() const
Definition: WellInterface_impl.hpp:111
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:189
Scalar wurea() const
Definition: WellInterface_impl.hpp:177
bool thp_update_iterations
Definition: WellInterface.hpp:365
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:947
Scalar wmicrobes() const
Definition: WellInterface_impl.hpp:153
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:89
virtual void scaleSegmentRatesAndPressure(WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1292
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:148
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