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 <limits>
66#include <vector>
67
68namespace Opm
69{
70
71class WellInjectionProperties;
72class WellProductionProperties;
73template<typename Scalar, typename IndexTraits> class GroupStateHelper;
74
75template<typename TypeTag>
76class WellInterface : public WellInterfaceIndices<GetPropType<TypeTag, Properties::FluidSystem>,
77 GetPropType<TypeTag, Properties::Indices>>
78{
81public:
86 using IndexTraits = typename FluidSystem::IndexTraitsType;
94
95 using VectorBlockType = Dune::FieldVector<Scalar, Indices::numEq>;
96 using MatrixBlockType = Dune::FieldMatrix<Scalar, Indices::numEq, Indices::numEq>;
97 using Eval = typename Base::Eval;
98 using BVector = Dune::BlockVector<VectorBlockType>;
99 using PressureMatrix = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, 1, 1>>;
100
104
107
111
113
114 static constexpr bool has_solvent = getPropValue<TypeTag, Properties::EnableSolvent>();
115 static constexpr bool has_zFraction = getPropValue<TypeTag, Properties::EnableExtbo>();
116 static constexpr bool has_polymer = getPropValue<TypeTag, Properties::EnablePolymer>();
117 static constexpr EnergyModules energyModuleType = getPropValue<TypeTag, Properties::EnergyModuleType>();
118 static constexpr bool has_energy = energyModuleType == EnergyModules::FullyImplicitThermal;
119
120 // flag for polymer molecular weight related
121 static constexpr bool has_polymermw = getPropValue<TypeTag, Properties::EnablePolymerMW>();
122 static constexpr bool has_foam = getPropValue<TypeTag, Properties::EnableFoam>();
123 static constexpr bool has_brine = getPropValue<TypeTag, Properties::EnableBrine>();
124 static constexpr bool has_watVapor = getPropValue<TypeTag, Properties::EnableVapwat>();
125 static constexpr bool has_disgas_in_water = getPropValue<TypeTag, Properties::EnableDisgasInWater>();
126 static constexpr bool has_saltPrecip = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
127 static constexpr bool has_bioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
128 static constexpr bool has_micp = Indices::enableMICP;
129
130 // For the conversion between the surface volume rate and reservoir voidage rate
131 using FluidState = BlackOilFluidState<Eval,
133 energyModuleType != EnergyModules::NoTemperature,
134 energyModuleType == EnergyModules::FullyImplicitThermal,
135 Indices::compositionSwitchIdx != std::numeric_limits<unsigned>::max(),
137 has_brine,
141 Indices::numPhases >;
143 WellInterface(const Well& well,
144 const ParallelWellInfo<Scalar>& pw_info,
145 const int time_step,
146 const ModelParameters& param,
147 const RateConverterType& rate_converter,
148 const int pvtRegionIdx,
149 const int num_conservation_quantities,
150 const int num_phases,
151 const int index_of_well,
152 const std::vector<PerforationData<Scalar>>& perf_data);
153
155 virtual ~WellInterface() = default;
156
157 virtual void init(const std::vector<Scalar>& depth_arg,
158 const Scalar gravity_arg,
159 const std::vector<Scalar>& B_avg,
160 const bool changed_to_open_this_step);
161
163 const std::vector<Scalar>& B_avg,
164 const bool relax_tolerance) const = 0;
165
166 virtual void solveEqAndUpdateWellState(const Simulator& simulator,
167 const GroupStateHelperType& groupStateHelper,
168 WellStateType& well_state) = 0;
169
170 void assembleWellEq(const Simulator& simulator,
171 const double dt,
172 const GroupStateHelperType& groupStateHelper,
173 WellStateType& well_state);
174
175 void assembleWellEqWithoutIteration(const Simulator& simulator,
176 const GroupStateHelperType& groupStateHelper,
177 const double dt,
178 WellStateType& well_state,
179 const bool solving_with_zero_rate);
180
181 // Convenience overload that gets scaled fractions internally
183 DeferredLogger& deferred_logger) const;
184
185 // TODO: better name or further refactoring the function to make it more clear
186 void prepareWellBeforeAssembling(const Simulator& simulator,
187 const double dt,
188 const GroupStateHelperType& groupStateHelper,
189 WellStateType& well_state);
190
191 virtual void computeWellRatesWithBhp(const Simulator& ebosSimulator,
192 const Scalar& bhp,
193 std::vector<Scalar>& well_flux,
194 DeferredLogger& deferred_logger) const = 0;
195
196 virtual std::optional<Scalar>
198 const GroupStateHelperType& groupStateHelper,
199 const SummaryState& summary_state,
200 const Scalar alq_value,
201 bool iterate_if_no_solution) const = 0;
202
203 std::optional<Scalar>
205 const WellStateType& well_state,
206 Scalar bhp,
207 const SummaryState& summary_state,
208 const Scalar alq_value);
212 const BVector& x,
213 const GroupStateHelperType& groupStateHelper,
214 WellStateType& well_state) = 0;
215
217 virtual void apply(const BVector& x, BVector& Ax) const = 0;
218
220 virtual void apply(BVector& r) const = 0;
221
222 // TODO: before we decide to put more information under mutable, this function is not const
223 virtual void computeWellPotentials(const Simulator& simulator,
224 const WellStateType& well_state,
225 const GroupStateHelperType& groupStateHelper,
226 std::vector<Scalar>& well_potentials) = 0;
227
228 virtual void updateWellStateWithTarget(const Simulator& simulator,
229 const GroupStateHelperType& groupStateHelper,
230 WellStateType& well_state) const;
231
232 virtual void scaleSegmentRatesAndPressure(WellStateType& well_state) const;
233
234 virtual void computeWellRatesWithBhpIterations(const Simulator& simulator,
235 const Scalar& bhp,
236 const GroupStateHelperType& groupStateHelper,
237 std::vector<Scalar>& well_flux) const = 0;
238
239 bool wellUnderZeroRateTarget(const GroupStateHelperType& groupStateHelper) const;
240
241 bool stoppedOrZeroRateTarget(const GroupStateHelperType& groupStateHelper) const;
242
243 bool updateWellStateWithTHPTargetProd(const Simulator& simulator,
244 WellStateType& well_state,
245 const GroupStateHelperType& groupStateHelper) const;
246
247 bool wellUnderZeroGroupRateTarget(const GroupStateHelperType& groupStateHelper,
248 const std::optional<bool> group_control = std::nullopt) const;
249
250 enum class IndividualOrGroup { Individual, Group, Both };
251 bool updateWellControl(const Simulator& simulator,
252 const IndividualOrGroup iog,
253 const GroupStateHelperType& groupStateHelper,
254 WellStateType& well_state) /* const */;
255
257 const GroupStateHelperType& groupStateHelper,
258 const Well::InjectionControls& inj_controls,
259 const Well::ProductionControls& prod_controls,
260 const Scalar WQTotal,
261 WellStateType& well_state,
262 const bool fixed_control,
263 const bool fixed_status,
264 const bool solving_with_zero_rate);
265
266 virtual void updatePrimaryVariables(const GroupStateHelperType& groupStateHelper) = 0;
267
268 virtual void calculateExplicitQuantities(const Simulator& simulator,
269 const GroupStateHelperType& groupStateHelper) = 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(std::map<int, RateVector>& cellRates_) 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 GroupStateHelperType& groupStateHelper,
294 WellStateType& well_state,
295 WellTestState& welltest_state,
296 GLiftEclWells& ecl_well_map,
297 std::map<std::string, double>& open_times);
298
299 void checkWellOperability(const Simulator& simulator,
300 const WellStateType& well_state,
301 const GroupStateHelperType& groupStateHelper);
302
304 WellStateType& well_state,
305 const GroupState<Scalar>& group_state,
306 GLiftEclWells& ecl_well_map,
307 DeferredLogger& deferred_logger);
308
309 // check whether the well is operable under the current reservoir condition
310 // mostly related to BHP limit and THP limit
311 void updateWellOperability(const Simulator& simulator,
312 const WellStateType& well_state,
313 const GroupStateHelperType& groupStateHelper);
314
315 bool updateWellOperabilityFromWellEq(const Simulator& simulator,
316 const GroupStateHelperType& groupStateHelper);
317
318 // update perforation water throughput based on solved water rate
319 virtual void updateWaterThroughput(const double dt,
320 WellStateType& well_state) const = 0;
321
324 virtual std::vector<Scalar>
326 DeferredLogger& deferred_logger) const = 0;
327
330 virtual void getScaledWellFractions(std::vector<Scalar>& scaled_fractions,
331 DeferredLogger& deferred_logger) const = 0;
332
336 void initializeProducerWellState(const Simulator& simulator,
337 WellStateType& well_state,
338 DeferredLogger& deferred_logger) const;
339
340 void solveWellEquation(const Simulator& simulator,
341 const GroupStateHelperType& groupStateHelper,
342 WellStateType& well_state);
343
344 const std::vector<RateVector>& connectionRates() const
345 {
346 return connectionRates_;
347 }
348
349 void updateConnectionDFactor(const Simulator& simulator,
350 SingleWellStateType& ws) const;
351
353 SingleWellStateType& ws) const;
354
355 virtual bool iterateWellEqWithSwitching(const Simulator& simulator,
356 const double dt,
357 const WellInjectionControls& inj_controls,
358 const WellProductionControls& prod_controls,
359 const GroupStateHelperType& groupStateHelper,
360 WellStateType& well_state,
361 const bool fixed_control,
362 const bool fixed_status,
363 const bool solving_with_zero_rate) = 0;
364
365 virtual Scalar maxPerfPress(const Simulator& simulator) const = 0;
366
367 virtual void updateIPRImplicit(const Simulator& simulator,
368 const GroupStateHelperType& groupStateHelper,
369 WellStateType& well_state) = 0;
370
371protected:
372 // simulation parameters
373 std::vector<RateVector> connectionRates_;
374 std::vector<Scalar> B_avg_;
378
379 Scalar wpolymer() const;
380 Scalar wfoam() const;
381 Scalar wsalt() const;
382 Scalar wmicrobes() const;
383 Scalar woxygen() const;
384 Scalar wurea() const;
385
386 virtual Scalar getRefDensity() const = 0;
387
388 std::vector<Scalar>
389 initialWellRateFractions(const Simulator& ebosSimulator,
390 const WellStateType& well_state) const;
391
392 // check whether the well is operable under BHP limit with current reservoir condition
393 virtual void checkOperabilityUnderBHPLimit(const WellStateType& well_state,
394 const Simulator& simulator,
395 DeferredLogger& deferred_logger) = 0;
396
397 // check whether the well is operable under THP limit with current reservoir condition
398 virtual void checkOperabilityUnderTHPLimit(const Simulator& simulator,
399 const WellStateType& well_state,
400 const GroupStateHelperType& groupStateHelper) = 0;
401
402 virtual void updateIPR(const Simulator& simulator,
403 DeferredLogger& deferred_logger) const = 0;
404
405 virtual void assembleWellEqWithoutIteration(const Simulator& simulator,
406 const GroupStateHelperType& groupStateHelper,
407 const double dt,
408 const WellInjectionControls& inj_controls,
409 const WellProductionControls& prod_controls,
410 WellStateType& well_state,
411 const bool solving_with_zero_rate) = 0;
412
413 // iterate well equations with the specified control until converged
414 virtual bool iterateWellEqWithControl(const Simulator& simulator,
415 const double dt,
416 const WellInjectionControls& inj_controls,
417 const WellProductionControls& prod_controls,
418 const GroupStateHelperType& groupStateHelper,
419 WellStateType& well_state) = 0;
420
421 bool iterateWellEquations(const Simulator& simulator,
422 const double dt,
423 const GroupStateHelperType& groupStateHelper,
424 WellStateType& well_state);
425
426 bool solveWellWithOperabilityCheck(const Simulator& simulator,
427 const double dt,
428 const Well::InjectionControls& inj_controls,
429 const Well::ProductionControls& prod_controls,
430 const GroupStateHelperType& groupStateHelper,
431 WellStateType& well_state);
432
433 std::optional<Scalar>
434 estimateOperableBhp(const Simulator& ebos_simulator,
435 const double dt,
436 const GroupStateHelperType& groupStateHelper,
437 const SummaryState& summary_state,
438 WellStateType& well_state);
439
440 bool solveWellWithBhp(const Simulator& simulator,
441 const double dt,
442 const Scalar bhp,
443 const GroupStateHelperType& groupStateHelper,
444 WellStateType& well_state);
445
446 bool solveWellWithZeroRate(const Simulator& simulator,
447 const double dt,
448 const GroupStateHelperType& groupStateHelper,
449 WellStateType& well_state);
450
451 bool solveWellForTesting(const Simulator& simulator,
452 const GroupStateHelperType& groupStateHelper,
453 WellStateType& well_state);
454
455
456 template<class GasLiftSingleWell>
457 std::unique_ptr<GasLiftSingleWell> initializeGliftWellTest_(const Simulator& simulator,
458 WellStateType& well_state,
459 const GroupState<Scalar>& group_state,
460 GLiftEclWells& ecl_well_map,
461 DeferredLogger& deferred_logger);
462
463 Eval getPerfCellPressure(const FluidState& fs) const;
464
465 // get the transmissibility multiplier for specific perforation
466 template<class Value, class Callback>
467 void getTransMult(Value& trans_mult,
468 const Simulator& simulator,
469 const int cell_idx,
470 Callback& extendEval) const;
471
472 // get the well transmissibility for specific perforation
473 template<class Value>
474 void getTw(std::vector<Value>& wi,
475 const int perf,
476 const IntensiveQuantities& intQuants,
477 const Value& trans_mult,
478 const SingleWellStateType& ws) const;
479
480 // get the mobility for specific perforation
481 template<class Value, class Callback>
482 void getMobility(const Simulator& simulator,
483 const int local_perf_index,
484 std::vector<Value>& mob,
485 Callback& extendEval,
486 [[maybe_unused]] DeferredLogger& deferred_logger) const;
487
488 void computeConnLevelProdInd(const FluidState& fs,
489 const std::function<Scalar(const Scalar)>& connPICalc,
490 const std::vector<Scalar>& mobility,
491 Scalar* connPI) const;
492
493 void computeConnLevelInjInd(const FluidState& fs,
494 const Phase preferred_phase,
495 const std::function<Scalar(const Scalar)>& connIICalc,
496 const std::vector<Scalar>& mobility,
497 Scalar* connII,
498 DeferredLogger& deferred_logger) const;
499
500 Scalar computeConnectionDFactor(const int perf,
501 const IntensiveQuantities& intQuants,
502 const SingleWellStateType& ws) const;
503};
504
505} // namespace Opm
506
507#include "WellInterface_impl.hpp"
508
509#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:78
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
BlackOilFluidState< Eval, FluidSystem, energyModuleType !=EnergyModules::NoTemperature, energyModuleType==EnergyModules::FullyImplicitThermal, Indices::compositionSwitchIdx !=std::numeric_limits< unsigned >::max(), has_watVapor, has_brine, has_saltPrecip, has_disgas_in_water, has_solvent, Indices::numPhases > FluidState
Definition: WellInterface.hpp:141
Scalar woxygen() const
Definition: WellInterface_impl.hpp:166
virtual Scalar maxPerfPress(const Simulator &simulator) const =0
static constexpr bool has_brine
Definition: WellInterface.hpp:123
IndividualOrGroup
Definition: WellInterface.hpp:250
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:83
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:106
const std::vector< RateVector > & connectionRates() const
Definition: WellInterface.hpp:344
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:373
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:99
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:88
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:84
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:124
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:89
virtual Scalar getRefDensity() const =0
typename FluidSystem::IndexTraitsType IndexTraits
Definition: WellInterface.hpp:86
Eval getPerfCellPressure(const FluidState &fs) const
Definition: WellInterface_impl.hpp:2081
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:98
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:374
static constexpr bool has_bioeffects
Definition: WellInterface.hpp:127
bool changed_to_stopped_this_step_
Definition: WellInterface.hpp:375
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:121
void addCellRates(std::map< int, RateVector > &cellRates_) const
Definition: WellInterface_impl.hpp:1145
static constexpr bool has_polymer
Definition: WellInterface.hpp:116
virtual void apply(BVector &r) const =0
r = r - C D^-1 Rw
GetPropType< TypeTag, Properties::Grid > Grid
Definition: WellInterface.hpp:82
virtual ~WellInterface()=default
Virtual destructor.
typename Base::ModelParameters ModelParameters
Definition: WellInterface.hpp:112
virtual void solveEqAndUpdateWellState(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)=0
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:85
static constexpr bool has_solvent
Definition: WellInterface.hpp:114
static constexpr bool has_disgas_in_water
Definition: WellInterface.hpp:125
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:126
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:122
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: WellInterface.hpp:91
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:93
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:377
Scalar wsalt() const
Definition: WellInterface_impl.hpp:142
static constexpr bool has_micp
Definition: WellInterface.hpp:128
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:96
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:97
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:101
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:118
Scalar wpolymer() const
Definition: WellInterface_impl.hpp:112
SingleWellState< Scalar, IndexTraits > SingleWellStateType
Definition: WellInterface.hpp:102
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:87
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:376
Dune::FieldVector< Scalar, Indices::numEq > VectorBlockType
Definition: WellInterface.hpp:95
void assembleWellEq(const Simulator &simulator, const double dt, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:967
Scalar wmicrobes() const
Definition: WellInterface_impl.hpp:154
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:90
virtual void scaleSegmentRatesAndPressure(WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1334
GroupStateHelper< Scalar, IndexTraits > GroupStateHelperType
Definition: WellInterface.hpp:103
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:115
static constexpr EnergyModules energyModuleType
Definition: WellInterface.hpp:117
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:165
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