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 WellStateType& well_state,
162 const std::vector<Scalar>& B_avg,
163 DeferredLogger& deferred_logger,
164 const bool relax_tolerance) const = 0;
165
166 virtual void solveEqAndUpdateWellState(const Simulator& simulator,
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 WellStateType& well_state,
215 DeferredLogger& deferred_logger) = 0;
216
218 virtual void apply(const BVector& x, BVector& Ax) const = 0;
219
221 virtual void apply(BVector& r) const = 0;
222
223 // TODO: before we decide to put more information under mutable, this function is not const
224 virtual void computeWellPotentials(const Simulator& simulator,
225 const WellStateType& well_state,
226 const GroupStateHelperType& groupStateHelper,
227 std::vector<Scalar>& well_potentials,
228 DeferredLogger& deferred_logger) = 0;
229
230 virtual void updateWellStateWithTarget(const Simulator& simulator,
231 const GroupStateHelperType& groupStateHelper,
232 WellStateType& well_state,
233 DeferredLogger& deferred_logger) const;
234
235 virtual void scaleSegmentRatesAndPressure(WellStateType& well_state) const;
236
237 virtual void computeWellRatesWithBhpIterations(const Simulator& simulator,
238 const Scalar& bhp,
239 const GroupStateHelperType& groupStateHelper,
240 std::vector<Scalar>& well_flux,
241 DeferredLogger& deferred_logger) const = 0;
242
243 bool wellUnderZeroRateTarget(const Simulator& simulator,
244 const WellStateType& well_state,
245 DeferredLogger& deferred_logger) const;
246
247 bool wellUnderZeroGroupRateTarget(const Simulator& simulator,
248 const WellStateType& well_state,
249 DeferredLogger& deferred_logger,
250 std::optional<bool> group_control = std::nullopt) const;
251
252 bool stoppedOrZeroRateTarget(const Simulator& simulator,
253 const WellStateType& well_state,
254 DeferredLogger& deferred_logger) const;
255
256 bool updateWellStateWithTHPTargetProd(const Simulator& simulator,
257 WellStateType& well_state,
258 const GroupStateHelperType& groupStateHelper,
259 DeferredLogger& deferred_logger) const;
260
261 enum class IndividualOrGroup { Individual, Group, Both };
262 bool updateWellControl(const Simulator& simulator,
263 const IndividualOrGroup iog,
264 const GroupStateHelperType& groupStateHelper,
265 WellStateType& well_state,
266 DeferredLogger& deferred_logger) /* const */;
267
269 const GroupStateHelperType& groupStateHelper,
270 const Well::InjectionControls& inj_controls,
271 const Well::ProductionControls& prod_controls,
272 const Scalar WQTotal,
273 WellStateType& well_state,
274 DeferredLogger& deferred_logger,
275 const bool fixed_control,
276 const bool fixed_status,
277 const bool solving_with_zero_rate);
278
279 virtual void updatePrimaryVariables(const Simulator& simulator,
280 const WellStateType& well_state,
281 DeferredLogger& deferred_logger) = 0;
282
283 virtual void calculateExplicitQuantities(const Simulator& simulator,
284 const WellStateType& well_state,
285 DeferredLogger& deferred_logger) = 0; // should be const?
286
287 virtual void updateProductivityIndex(const Simulator& simulator,
288 const WellProdIndexCalculator<Scalar>& wellPICalc,
289 WellStateType& well_state,
290 DeferredLogger& deferred_logger) const = 0;
291
292 // Add well contributions to matrix
294
296 const BVector& x,
297 const int pressureVarIndex,
298 const bool use_well_weights,
299 const WellStateType& well_state) const = 0;
300
301 void addCellRates(std::map<int, RateVector>& cellRates_) const;
302
303 Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const;
304
305 // TODO: theoretically, it should be a const function
306 // Simulator is not const is because that assembleWellEq is non-const Simulator
307 void wellTesting(const Simulator& simulator,
308 const double simulation_time,
309 const GroupStateHelperType& groupStateHelper,
310 WellStateType& well_state,
311 WellTestState& welltest_state,
312 GLiftEclWells& ecl_well_map,
313 std::map<std::string, double>& open_times,
314 DeferredLogger& deferred_logger);
315
316 void checkWellOperability(const Simulator& simulator,
317 const WellStateType& well_state,
318 const GroupStateHelperType& groupStateHelper,
319 DeferredLogger& deferred_logger);
320
322 WellStateType& well_state,
323 const GroupState<Scalar>& group_state,
324 GLiftEclWells& ecl_well_map,
325 DeferredLogger& deferred_logger);
326
327 // check whether the well is operable under the current reservoir condition
328 // mostly related to BHP limit and THP limit
329 void updateWellOperability(const Simulator& simulator,
330 const WellStateType& well_state,
331 const GroupStateHelperType& groupStateHelper,
332 DeferredLogger& deferred_logger);
333
334 bool updateWellOperabilityFromWellEq(const Simulator& simulator,
335 const GroupStateHelperType& groupStateHelper,
336 DeferredLogger& deferred_logger);
337
338 // update perforation water throughput based on solved water rate
339 virtual void updateWaterThroughput(const double dt,
340 WellStateType& well_state) const = 0;
341
344 virtual std::vector<Scalar>
346 DeferredLogger& deferred_logger) const = 0;
347
351 void initializeProducerWellState(const Simulator& simulator,
352 WellStateType& well_state,
353 DeferredLogger& deferred_logger) const;
354
355 void solveWellEquation(const Simulator& simulator,
356 const GroupStateHelperType& groupStateHelper,
357 WellStateType& well_state,
358 DeferredLogger& deferred_logger);
359
360 const std::vector<RateVector>& connectionRates() const
361 {
362 return connectionRates_;
363 }
364
365 void updateConnectionDFactor(const Simulator& simulator,
366 SingleWellStateType& ws) const;
367
369 SingleWellStateType& ws) const;
370
371 virtual bool iterateWellEqWithSwitching(const Simulator& simulator,
372 const double dt,
373 const WellInjectionControls& inj_controls,
374 const WellProductionControls& prod_controls,
375 const GroupStateHelperType& groupStateHelper,
376 WellStateType& well_state,
377 DeferredLogger& deferred_logger,
378 const bool fixed_control,
379 const bool fixed_status,
380 const bool solving_with_zero_rate) = 0;
381
382 virtual Scalar maxPerfPress(const Simulator& simulator) const = 0;
383
384 virtual void updateIPRImplicit(const Simulator& simulator,
385 const GroupStateHelperType& groupStateHelper,
386 WellStateType& well_state,
387 DeferredLogger& deferred_logger) = 0;
388
389protected:
390 // simulation parameters
391 std::vector<RateVector> connectionRates_;
392 std::vector<Scalar> B_avg_;
396
397 Scalar wpolymer() const;
398 Scalar wfoam() const;
399 Scalar wsalt() const;
400 Scalar wmicrobes() const;
401 Scalar woxygen() const;
402 Scalar wurea() const;
403
404 virtual Scalar getRefDensity() const = 0;
405
406 std::vector<Scalar>
407 initialWellRateFractions(const Simulator& ebosSimulator,
408 const WellStateType& well_state) const;
409
410 // check whether the well is operable under BHP limit with current reservoir condition
411 virtual void checkOperabilityUnderBHPLimit(const WellStateType& well_state,
412 const Simulator& simulator,
413 DeferredLogger& deferred_logger) = 0;
414
415 // check whether the well is operable under THP limit with current reservoir condition
416 virtual void checkOperabilityUnderTHPLimit(const Simulator& simulator,
417 const WellStateType& well_state,
418 const GroupStateHelperType& groupStateHelper,
419 DeferredLogger& deferred_logger) = 0;
420
421 virtual void updateIPR(const Simulator& simulator,
422 DeferredLogger& deferred_logger) const = 0;
423
424 virtual void assembleWellEqWithoutIteration(const Simulator& simulator,
425 const GroupStateHelperType& groupStateHelper,
426 const double dt,
427 const WellInjectionControls& inj_controls,
428 const WellProductionControls& prod_controls,
429 WellStateType& well_state,
430 DeferredLogger& deferred_logger,
431 const bool solving_with_zero_rate) = 0;
432
433 // iterate well equations with the specified control until converged
434 virtual bool iterateWellEqWithControl(const Simulator& simulator,
435 const double dt,
436 const WellInjectionControls& inj_controls,
437 const WellProductionControls& prod_controls,
438 const GroupStateHelperType& groupStateHelper,
439 WellStateType& well_state,
440 DeferredLogger& deferred_logger) = 0;
441
442 bool iterateWellEquations(const Simulator& simulator,
443 const double dt,
444 const GroupStateHelperType& groupStateHelper,
445 WellStateType& well_state,
446 DeferredLogger& deferred_logger);
447
448 bool solveWellWithOperabilityCheck(const Simulator& simulator,
449 const double dt,
450 const Well::InjectionControls& inj_controls,
451 const Well::ProductionControls& prod_controls,
452 const GroupStateHelperType& groupStateHelper,
453 WellStateType& well_state,
454 DeferredLogger& deferred_logger);
455
456 std::optional<Scalar>
457 estimateOperableBhp(const Simulator& ebos_simulator,
458 const double dt,
459 const GroupStateHelperType& groupStateHelper,
460 const SummaryState& summary_state,
461 WellStateType& well_state,
462 DeferredLogger& deferred_logger);
463
464 bool solveWellWithBhp(const Simulator& simulator,
465 const double dt,
466 const Scalar bhp,
467 const GroupStateHelperType& groupStateHelper,
468 WellStateType& well_state,
469 DeferredLogger& deferred_logger);
470
471 bool solveWellWithZeroRate(const Simulator& simulator,
472 const double dt,
473 const GroupStateHelperType& groupStateHelper,
474 WellStateType& well_state,
475 DeferredLogger& deferred_logger);
476
477 bool solveWellForTesting(const Simulator& simulator,
478 const GroupStateHelperType& groupStateHelper,
479 WellStateType& well_state,
480 DeferredLogger& deferred_logger);
481
482
483 template<class GasLiftSingleWell>
484 std::unique_ptr<GasLiftSingleWell> initializeGliftWellTest_(const Simulator& simulator,
485 WellStateType& well_state,
486 const GroupState<Scalar>& group_state,
487 GLiftEclWells& ecl_well_map,
488 DeferredLogger& deferred_logger);
489
490 Eval getPerfCellPressure(const FluidState& fs) const;
491
492 // get the transmissibility multiplier for specific perforation
493 template<class Value, class Callback>
494 void getTransMult(Value& trans_mult,
495 const Simulator& simulator,
496 const int cell_idx,
497 Callback& extendEval) const;
498
499 // get the well transmissibility for specific perforation
500 template<class Value>
501 void getTw(std::vector<Value>& wi,
502 const int perf,
503 const IntensiveQuantities& intQuants,
504 const Value& trans_mult,
505 const SingleWellStateType& ws) const;
506
507 // get the mobility for specific perforation
508 template<class Value, class Callback>
509 void getMobility(const Simulator& simulator,
510 const int local_perf_index,
511 std::vector<Value>& mob,
512 Callback& extendEval,
513 [[maybe_unused]] DeferredLogger& deferred_logger) const;
514
515 void computeConnLevelProdInd(const FluidState& fs,
516 const std::function<Scalar(const Scalar)>& connPICalc,
517 const std::vector<Scalar>& mobility,
518 Scalar* connPI) const;
519
520 void computeConnLevelInjInd(const FluidState& fs,
521 const Phase preferred_phase,
522 const std::function<Scalar(const Scalar)>& connIICalc,
523 const std::vector<Scalar>& mobility,
524 Scalar* connII,
525 DeferredLogger& deferred_logger) const;
526
527 Scalar computeConnectionDFactor(const int perf,
528 const IntensiveQuantities& intQuants,
529 const SingleWellStateType& ws) const;
530};
531
532} // namespace Opm
533
534#include "WellInterface_impl.hpp"
535
536#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:261
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
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:1984
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:2195
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:360
Scalar wfoam() const
Definition: WellInterface_impl.hpp:127
virtual void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const =0
void getTransMult(Value &trans_mult, const Simulator &simulator, const int cell_idx, Callback &extendEval) const
Definition: WellInterface_impl.hpp:2072
std::vector< RateVector > connectionRates_
Definition: WellInterface.hpp:391
void computeConnLevelProdInd(const FluidState &fs, const std::function< Scalar(const Scalar)> &connPICalc, const std::vector< Scalar > &mobility, Scalar *connPI) const
Definition: WellInterface_impl.hpp:2231
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:1886
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:2085
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
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:2163
std::vector< Scalar > initialWellRateFractions(const Simulator &ebosSimulator, const WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1748
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:1964
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:2057
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:97
void initializeProducerWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1798
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:392
static constexpr bool has_bioeffects
Definition: WellInterface.hpp:126
bool changed_to_stopped_this_step_
Definition: WellInterface.hpp:393
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
virtual void calculateExplicitQuantities(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger)=0
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 wellUnderZeroGroupRateTarget(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger, std::optional< bool > group_control=std::nullopt) const
Definition: WellInterface_impl.hpp:1716
virtual ~WellInterface()=default
Virtual destructor.
typename Base::ModelParameters ModelParameters
Definition: WellInterface.hpp:111
virtual void updatePrimaryVariables(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger)=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
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
virtual void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, WellStateType &well_state, DeferredLogger &deferred_logger)=0
void updateConnectionTransmissibilityFactor(const Simulator &simulator, SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:2026
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:2265
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:2299
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:395
virtual void solveEqAndUpdateWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger)=0
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 wellUnderZeroRateTarget(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1698
virtual ConvergenceReport getWellConvergence(const Simulator &simulator, const WellStateType &well_state, const std::vector< Scalar > &B_avg, DeferredLogger &deferred_logger, const bool relax_tolerance) const =0
WellInterface(const Well &well, const ParallelWellInfo< Scalar > &pw_info, const int time_step, const ModelParameters &param, const RateConverterType &rate_converter, const int pvtRegionIdx, const int num_conservation_quantities, const int num_phases, const int index_of_well, const std::vector< PerforationData< Scalar > > &perf_data)
Constructor.
Definition: WellInterface_impl.hpp:58
WellState< Scalar, IndexTraits > WellStateType
Definition: WellInterface.hpp:100
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:394
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
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
bool stoppedOrZeroRateTarget(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1735
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:141
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