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 DeferredLogger& deferred_logger,
163 const bool relax_tolerance) const = 0;
164
165 virtual void solveEqAndUpdateWellState(const Simulator& simulator,
166 const GroupStateHelperType& groupStateHelper,
167 WellStateType& well_state,
168 DeferredLogger& deferred_logger) = 0;
169
170 void assembleWellEq(const Simulator& simulator,
171 const double dt,
172 const GroupStateHelperType& groupStateHelper,
173 WellStateType& well_state,
174 DeferredLogger& deferred_logger);
175
176 void assembleWellEqWithoutIteration(const Simulator& simulator,
177 const GroupStateHelperType& groupStateHelper,
178 const double dt,
179 WellStateType& well_state,
180 DeferredLogger& deferred_logger,
181 const bool solving_with_zero_rate);
182
183 // TODO: better name or further refactoring the function to make it more clear
184 void prepareWellBeforeAssembling(const Simulator& simulator,
185 const double dt,
186 const GroupStateHelperType& groupStateHelper,
187 WellStateType& well_state,
188 DeferredLogger& deferred_logger);
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 DeferredLogger& deferred_logger,
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,
209 DeferredLogger& deferred_logger);
213 const BVector& x,
214 const GroupStateHelperType& groupStateHelper,
215 WellStateType& well_state,
216 DeferredLogger& deferred_logger) = 0;
217
219 virtual void apply(const BVector& x, BVector& Ax) const = 0;
220
222 virtual void apply(BVector& r) const = 0;
223
224 // TODO: before we decide to put more information under mutable, this function is not const
225 virtual void computeWellPotentials(const Simulator& simulator,
226 const WellStateType& well_state,
227 const GroupStateHelperType& groupStateHelper,
228 std::vector<Scalar>& well_potentials,
229 DeferredLogger& deferred_logger) = 0;
230
231 virtual void updateWellStateWithTarget(const Simulator& simulator,
232 const GroupStateHelperType& groupStateHelper,
233 WellStateType& well_state,
234 DeferredLogger& deferred_logger) const;
235
236 virtual void scaleSegmentRatesAndPressure(WellStateType& well_state) const;
237
238 virtual void computeWellRatesWithBhpIterations(const Simulator& simulator,
239 const Scalar& bhp,
240 const GroupStateHelperType& groupStateHelper,
241 std::vector<Scalar>& well_flux,
242 DeferredLogger& deferred_logger) const = 0;
243
244 bool wellUnderZeroRateTarget(const GroupStateHelperType& groupStateHelper,
245 DeferredLogger& deferred_logger) const;
246
247 bool stoppedOrZeroRateTarget(const GroupStateHelperType& groupStateHelper,
248 DeferredLogger& deferred_logger) const;
249
250 bool updateWellStateWithTHPTargetProd(const Simulator& simulator,
251 WellStateType& well_state,
252 const GroupStateHelperType& groupStateHelper,
253 DeferredLogger& deferred_logger) const;
254
255 bool wellUnderZeroGroupRateTarget(const GroupStateHelperType& groupStateHelper,
256 DeferredLogger& deferred_logger,
257 const std::optional<bool> group_control = std::nullopt) const;
258
259 enum class IndividualOrGroup { Individual, Group, Both };
260 bool updateWellControl(const Simulator& simulator,
261 const IndividualOrGroup iog,
262 const GroupStateHelperType& groupStateHelper,
263 WellStateType& well_state,
264 DeferredLogger& deferred_logger) /* const */;
265
267 const GroupStateHelperType& groupStateHelper,
268 const Well::InjectionControls& inj_controls,
269 const Well::ProductionControls& prod_controls,
270 const Scalar WQTotal,
271 WellStateType& well_state,
272 DeferredLogger& deferred_logger,
273 const bool fixed_control,
274 const bool fixed_status,
275 const bool solving_with_zero_rate);
276
277 virtual void updatePrimaryVariables(const GroupStateHelperType& groupStateHelper,
278 DeferredLogger& deferred_logger) = 0;
279
280 virtual void calculateExplicitQuantities(const Simulator& simulator,
281 const GroupStateHelperType& groupStateHelper,
282 DeferredLogger& deferred_logger) = 0; // should be const?
283
284 virtual void updateProductivityIndex(const Simulator& simulator,
285 const WellProdIndexCalculator<Scalar>& wellPICalc,
286 WellStateType& well_state,
287 DeferredLogger& deferred_logger) const = 0;
288
289 // Add well contributions to matrix
291
293 const BVector& x,
294 const int pressureVarIndex,
295 const bool use_well_weights,
296 const WellStateType& well_state) const = 0;
297
298 void addCellRates(std::map<int, RateVector>& cellRates_) const;
299
300 Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const;
301
302 // TODO: theoretically, it should be a const function
303 // Simulator is not const is because that assembleWellEq is non-const Simulator
304 void wellTesting(const Simulator& simulator,
305 const double simulation_time,
306 const GroupStateHelperType& groupStateHelper,
307 WellStateType& well_state,
308 WellTestState& welltest_state,
309 GLiftEclWells& ecl_well_map,
310 std::map<std::string, double>& open_times,
311 DeferredLogger& deferred_logger);
312
313 void checkWellOperability(const Simulator& simulator,
314 const WellStateType& well_state,
315 const GroupStateHelperType& groupStateHelper,
316 DeferredLogger& deferred_logger);
317
319 WellStateType& well_state,
320 const GroupState<Scalar>& group_state,
321 GLiftEclWells& ecl_well_map,
322 DeferredLogger& deferred_logger);
323
324 // check whether the well is operable under the current reservoir condition
325 // mostly related to BHP limit and THP limit
326 void updateWellOperability(const Simulator& simulator,
327 const WellStateType& well_state,
328 const GroupStateHelperType& groupStateHelper,
329 DeferredLogger& deferred_logger);
330
331 bool updateWellOperabilityFromWellEq(const Simulator& simulator,
332 const GroupStateHelperType& groupStateHelper,
333 DeferredLogger& deferred_logger);
334
335 // update perforation water throughput based on solved water rate
336 virtual void updateWaterThroughput(const double dt,
337 WellStateType& well_state) const = 0;
338
341 virtual std::vector<Scalar>
343 DeferredLogger& deferred_logger) const = 0;
344
348 void initializeProducerWellState(const Simulator& simulator,
349 WellStateType& well_state,
350 DeferredLogger& deferred_logger) const;
351
352 void solveWellEquation(const Simulator& simulator,
353 const GroupStateHelperType& groupStateHelper,
354 WellStateType& well_state,
355 DeferredLogger& deferred_logger);
356
357 const std::vector<RateVector>& connectionRates() const
358 {
359 return connectionRates_;
360 }
361
362 void updateConnectionDFactor(const Simulator& simulator,
363 SingleWellStateType& ws) const;
364
366 SingleWellStateType& ws) const;
367
368 virtual bool iterateWellEqWithSwitching(const Simulator& simulator,
369 const double dt,
370 const WellInjectionControls& inj_controls,
371 const WellProductionControls& prod_controls,
372 const GroupStateHelperType& groupStateHelper,
373 WellStateType& well_state,
374 DeferredLogger& deferred_logger,
375 const bool fixed_control,
376 const bool fixed_status,
377 const bool solving_with_zero_rate) = 0;
378
379 virtual Scalar maxPerfPress(const Simulator& simulator) const = 0;
380
381 virtual void updateIPRImplicit(const Simulator& simulator,
382 const GroupStateHelperType& groupStateHelper,
383 WellStateType& well_state,
384 DeferredLogger& deferred_logger) = 0;
385
386protected:
387 // simulation parameters
388 std::vector<RateVector> connectionRates_;
389 std::vector<Scalar> B_avg_;
393
394 Scalar wpolymer() const;
395 Scalar wfoam() const;
396 Scalar wsalt() const;
397 Scalar wmicrobes() const;
398 Scalar woxygen() const;
399 Scalar wurea() const;
400
401 virtual Scalar getRefDensity() const = 0;
402
403 std::vector<Scalar>
404 initialWellRateFractions(const Simulator& ebosSimulator,
405 const WellStateType& well_state) const;
406
407 // check whether the well is operable under BHP limit with current reservoir condition
408 virtual void checkOperabilityUnderBHPLimit(const WellStateType& well_state,
409 const Simulator& simulator,
410 DeferredLogger& deferred_logger) = 0;
411
412 // check whether the well is operable under THP limit with current reservoir condition
413 virtual void checkOperabilityUnderTHPLimit(const Simulator& simulator,
414 const WellStateType& well_state,
415 const GroupStateHelperType& groupStateHelper,
416 DeferredLogger& deferred_logger) = 0;
417
418 virtual void updateIPR(const Simulator& simulator,
419 DeferredLogger& deferred_logger) const = 0;
420
421 virtual void assembleWellEqWithoutIteration(const Simulator& simulator,
422 const GroupStateHelperType& groupStateHelper,
423 const double dt,
424 const WellInjectionControls& inj_controls,
425 const WellProductionControls& prod_controls,
426 WellStateType& well_state,
427 DeferredLogger& deferred_logger,
428 const bool solving_with_zero_rate) = 0;
429
430 // iterate well equations with the specified control until converged
431 virtual bool iterateWellEqWithControl(const Simulator& simulator,
432 const double dt,
433 const WellInjectionControls& inj_controls,
434 const WellProductionControls& prod_controls,
435 const GroupStateHelperType& groupStateHelper,
436 WellStateType& well_state,
437 DeferredLogger& deferred_logger) = 0;
438
439 bool iterateWellEquations(const Simulator& simulator,
440 const double dt,
441 const GroupStateHelperType& groupStateHelper,
442 WellStateType& well_state,
443 DeferredLogger& deferred_logger);
444
445 bool solveWellWithOperabilityCheck(const Simulator& simulator,
446 const double dt,
447 const Well::InjectionControls& inj_controls,
448 const Well::ProductionControls& prod_controls,
449 const GroupStateHelperType& groupStateHelper,
450 WellStateType& well_state,
451 DeferredLogger& deferred_logger);
452
453 std::optional<Scalar>
454 estimateOperableBhp(const Simulator& ebos_simulator,
455 const double dt,
456 const GroupStateHelperType& groupStateHelper,
457 const SummaryState& summary_state,
458 WellStateType& well_state,
459 DeferredLogger& deferred_logger);
460
461 bool solveWellWithBhp(const Simulator& simulator,
462 const double dt,
463 const Scalar bhp,
464 const GroupStateHelperType& groupStateHelper,
465 WellStateType& well_state,
466 DeferredLogger& deferred_logger);
467
468 bool solveWellWithZeroRate(const Simulator& simulator,
469 const double dt,
470 const GroupStateHelperType& groupStateHelper,
471 WellStateType& well_state,
472 DeferredLogger& deferred_logger);
473
474 bool solveWellForTesting(const Simulator& simulator,
475 const GroupStateHelperType& groupStateHelper,
476 WellStateType& well_state,
477 DeferredLogger& deferred_logger);
478
479
480 template<class GasLiftSingleWell>
481 std::unique_ptr<GasLiftSingleWell> initializeGliftWellTest_(const Simulator& simulator,
482 WellStateType& well_state,
483 const GroupState<Scalar>& group_state,
484 GLiftEclWells& ecl_well_map,
485 DeferredLogger& deferred_logger);
486
487 Eval getPerfCellPressure(const FluidState& fs) const;
488
489 // get the transmissibility multiplier for specific perforation
490 template<class Value, class Callback>
491 void getTransMult(Value& trans_mult,
492 const Simulator& simulator,
493 const int cell_idx,
494 Callback& extendEval) const;
495
496 // get the well transmissibility for specific perforation
497 template<class Value>
498 void getTw(std::vector<Value>& wi,
499 const int perf,
500 const IntensiveQuantities& intQuants,
501 const Value& trans_mult,
502 const SingleWellStateType& ws) const;
503
504 // get the mobility for specific perforation
505 template<class Value, class Callback>
506 void getMobility(const Simulator& simulator,
507 const int local_perf_index,
508 std::vector<Value>& mob,
509 Callback& extendEval,
510 [[maybe_unused]] DeferredLogger& deferred_logger) const;
511
512 void computeConnLevelProdInd(const FluidState& fs,
513 const std::function<Scalar(const Scalar)>& connPICalc,
514 const std::vector<Scalar>& mobility,
515 Scalar* connPI) const;
516
517 void computeConnLevelInjInd(const FluidState& fs,
518 const Phase preferred_phase,
519 const std::function<Scalar(const Scalar)>& connIICalc,
520 const std::vector<Scalar>& mobility,
521 Scalar* connII,
522 DeferredLogger& deferred_logger) const;
523
524 Scalar computeConnectionDFactor(const int perf,
525 const IntensiveQuantities& intQuants,
526 const SingleWellStateType& ws) const;
527};
528
529} // namespace Opm
530
531#include "WellInterface_impl.hpp"
532
533#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:53
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 solveWellWithOperabilityCheck(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:589
virtual std::vector< Scalar > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const =0
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:259
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:82
virtual std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger, bool iterate_if_no_solution) const =0
virtual void solveEqAndUpdateWellState(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger)=0
void assembleWellEqWithoutIteration(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const double dt, WellStateType &well_state, DeferredLogger &deferred_logger, const bool solving_with_zero_rate)
Definition: WellInterface_impl.hpp:964
Scalar computeConnectionDFactor(const int perf, const IntensiveQuantities &intQuants, const SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:1972
virtual void updateIPRImplicit(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger)=0
void checkWellOperability(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1144
std::optional< Scalar > computeBhpAtThpLimitProdWithAlqUsingIPR(const Simulator &simulator, const WellStateType &well_state, Scalar bhp, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:2183
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
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, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:395
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:105
const std::vector< RateVector > & connectionRates() const
Definition: WellInterface.hpp:357
Scalar wfoam() const
Definition: WellInterface_impl.hpp:127
virtual void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const =0
virtual void updatePrimaryVariables(const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger)=0
void getTransMult(Value &trans_mult, const Simulator &simulator, const int cell_idx, Callback &extendEval) const
Definition: WellInterface_impl.hpp:2060
std::vector< RateVector > connectionRates_
Definition: WellInterface.hpp:388
void computeConnLevelProdInd(const FluidState &fs, const std::function< Scalar(const Scalar)> &connPICalc, const std::vector< Scalar > &mobility, Scalar *connPI) const
Definition: WellInterface_impl.hpp:2219
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:98
bool iterateWellEquations(const Simulator &simulator, const double dt, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:521
void gliftBeginTimeStepWellTestUpdateALQ(const Simulator &simulator, WellStateType &well_state, const GroupState< Scalar > &group_state, GLiftEclWells &ecl_well_map, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1172
Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const
Definition: WellInterface_impl.hpp:1124
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
virtual void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger)=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:1874
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:2073
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
bool wellUnderZeroGroupRateTarget(const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger, const std::optional< bool > group_control=std::nullopt) const
Definition: WellInterface_impl.hpp:1708
virtual void checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger)=0
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
bool updateWellStateWithTHPTargetProd(const Simulator &simulator, WellStateType &well_state, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:2151
std::vector< Scalar > initialWellRateFractions(const Simulator &ebosSimulator, const WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1736
virtual void computeWellRatesWithBhpIterations(const Simulator &simulator, const Scalar &bhp, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const =0
virtual void updateWellStateWithTarget(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1303
static constexpr bool has_watVapor
Definition: WellInterface.hpp:123
void updateConnectionDFactor(const Simulator &simulator, SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:1952
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:2045
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:97
void initializeProducerWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1786
void updateWellOperability(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1234
std::vector< Scalar > B_avg_
Definition: WellInterface.hpp:389
static constexpr bool has_bioeffects
Definition: WellInterface.hpp:126
bool changed_to_stopped_this_step_
Definition: WellInterface.hpp:390
bool updateWellControl(const Simulator &simulator, const IndividualOrGroup iog, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:189
static constexpr bool has_polymermw
Definition: WellInterface.hpp:120
void addCellRates(std::map< int, RateVector > &cellRates_) const
Definition: WellInterface_impl.hpp:1105
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 void assembleWellEqWithoutIteration(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const double dt, const WellInjectionControls &inj_controls, const WellProductionControls &prod_controls, WellStateType &well_state, DeferredLogger &deferred_logger, const bool solving_with_zero_rate)=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, DeferredLogger &deferred_logger, const bool fixed_control, const bool fixed_status, const bool solving_with_zero_rate)
Definition: WellInterface_impl.hpp:277
bool wellUnderZeroRateTarget(const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1691
virtual ~WellInterface()=default
Virtual destructor.
virtual void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger)=0
typename Base::ModelParameters ModelParameters
Definition: WellInterface.hpp:111
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
static constexpr bool has_saltPrecip
Definition: WellInterface.hpp:125
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:2014
virtual bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const WellInjectionControls &inj_controls, const WellProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger, const bool fixed_control, const bool fixed_status, const bool solving_with_zero_rate)=0
void computeConnLevelInjInd(const FluidState &fs, const Phase preferred_phase, const std::function< Scalar(const Scalar)> &connIICalc, const std::vector< Scalar > &mobility, Scalar *connII, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:2253
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:2287
virtual void addWellPressureEquations(PressureMatrix &mat, const BVector &x, const int pressureVarIndex, const bool use_well_weights, const WellStateType &well_state) const =0
int number_of_well_reopenings_
Definition: WellInterface.hpp:392
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, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:797
bool updateWellOperabilityFromWellEq(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1272
Dune::FieldMatrix< Scalar, Indices::numEq, Indices::numEq > MatrixBlockType
Definition: WellInterface.hpp:95
typename Base::Eval Eval
Definition: WellInterface.hpp:96
bool stoppedOrZeroRateTarget(const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1724
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
virtual ConvergenceReport getWellConvergence(const GroupStateHelperType &groupStateHelper, const std::vector< Scalar > &B_avg, DeferredLogger &deferred_logger, const bool relax_tolerance) const =0
void solveWellEquation(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:885
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 void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const =0
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:86
bool solveWellForTesting(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:830
Scalar wurea() const
Definition: WellInterface_impl.hpp:177
bool thp_update_iterations
Definition: WellInterface.hpp:391
Dune::FieldVector< Scalar, Indices::numEq > VectorBlockType
Definition: WellInterface.hpp:94
Scalar wmicrobes() const
Definition: WellInterface_impl.hpp:153
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:89
void prepareWellBeforeAssembling(const Simulator &simulator, const double dt, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:985
virtual void scaleSegmentRatesAndPressure(WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1295
virtual void calculateExplicitQuantities(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger)=0
std::optional< Scalar > estimateOperableBhp(const Simulator &ebos_simulator, const double dt, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:709
GroupStateHelper< Scalar, IndexTraits > GroupStateHelperType
Definition: WellInterface.hpp:102
bool solveWellWithBhp(const Simulator &simulator, const double dt, const Scalar bhp, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:747
virtual bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const WellInjectionControls &inj_controls, const WellProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger)=0
static constexpr bool has_zFraction
Definition: WellInterface.hpp:114
void assembleWellEq(const Simulator &simulator, const double dt, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:947
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:147
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