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 WellGroupHelper;
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 bool has_energy = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal;
117 static const bool has_temperature = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::ConstantTemperature
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,
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 WellGroupHelperType& wgHelper,
173 WellStateType& well_state,
174 DeferredLogger& deferred_logger);
175
176 void assembleWellEqWithoutIteration(const Simulator& simulator,
177 const double dt,
178 WellStateType& well_state,
179 const GroupState<Scalar>& group_state,
180 DeferredLogger& deferred_logger);
181
182 // TODO: better name or further refactoring the function to make it more clear
183 void prepareWellBeforeAssembling(const Simulator& simulator,
184 const double dt,
185 const WellGroupHelperType& wgHelper,
186 WellStateType& well_state,
187 DeferredLogger& deferred_logger);
188
189 virtual void computeWellRatesWithBhp(const Simulator& ebosSimulator,
190 const Scalar& bhp,
191 std::vector<Scalar>& well_flux,
192 DeferredLogger& deferred_logger) const = 0;
193
194 virtual std::optional<Scalar>
196 const WellGroupHelperType& wgHelper,
197 const SummaryState& summary_state,
198 const Scalar alq_value,
199 DeferredLogger& deferred_logger,
200 bool iterate_if_no_solution) const = 0;
201
205 const BVector& x,
206 WellStateType& well_state,
207 DeferredLogger& deferred_logger) = 0;
208
210 virtual void apply(const BVector& x, BVector& Ax) const = 0;
211
213 virtual void apply(BVector& r) const = 0;
214
215 // TODO: before we decide to put more information under mutable, this function is not const
216 virtual void computeWellPotentials(const Simulator& simulator,
217 const WellStateType& well_state,
218 const WellGroupHelperType& wgHelper,
219 std::vector<Scalar>& well_potentials,
220 DeferredLogger& deferred_logger) = 0;
221
222 virtual void updateWellStateWithTarget(const Simulator& simulator,
223 const WellGroupHelperType& wgHelper,
224 WellStateType& well_state,
225 DeferredLogger& deferred_logger) const;
226
227 virtual void scaleSegmentRatesAndPressure(WellStateType& well_state) const;
228
229 virtual void computeWellRatesWithBhpIterations(const Simulator& simulator,
230 const Scalar& bhp,
231 const WellGroupHelperType& wgHelper,
232 std::vector<Scalar>& well_flux,
233 DeferredLogger& deferred_logger) const = 0;
234
235 bool wellUnderZeroRateTarget(const Simulator& simulator,
236 const WellStateType& well_state,
237 DeferredLogger& deferred_logger) const;
238
239 bool wellUnderZeroGroupRateTarget(const Simulator& simulator,
240 const WellStateType& well_state,
241 DeferredLogger& deferred_logger,
242 std::optional<bool> group_control = std::nullopt) const;
243
244 bool stoppedOrZeroRateTarget(const Simulator& simulator,
245 const WellStateType& well_state,
246 DeferredLogger& deferred_logger) const;
247
248 bool updateWellStateWithTHPTargetProd(const Simulator& simulator,
249 WellStateType& well_state,
250 const WellGroupHelperType& wgHelper,
251 DeferredLogger& deferred_logger) const;
252
253 enum class IndividualOrGroup { Individual, Group, Both };
254 bool updateWellControl(const Simulator& simulator,
255 const IndividualOrGroup iog,
256 const WellGroupHelperType& wgHelper,
257 WellStateType& well_state,
258 DeferredLogger& deferred_logger) /* const */;
259
261 const WellGroupHelperType& wgHelper,
262 const Well::InjectionControls& inj_controls,
263 const Well::ProductionControls& prod_controls,
264 const Scalar WQTotal,
265 WellStateType& well_state,
266 DeferredLogger& deferred_logger,
267 const bool fixed_control = false,
268 const bool fixed_status = false);
269
270 virtual void updatePrimaryVariables(const Simulator& simulator,
271 const WellStateType& well_state,
272 DeferredLogger& deferred_logger) = 0;
273
274 virtual void calculateExplicitQuantities(const Simulator& simulator,
275 const WellStateType& well_state,
276 DeferredLogger& deferred_logger) = 0; // should be const?
277
278 virtual void updateProductivityIndex(const Simulator& simulator,
279 const WellProdIndexCalculator<Scalar>& wellPICalc,
280 WellStateType& well_state,
281 DeferredLogger& deferred_logger) const = 0;
282
283 // Add well contributions to matrix
285
287 const BVector& x,
288 const int pressureVarIndex,
289 const bool use_well_weights,
290 const WellStateType& well_state) const = 0;
291
292 void addCellRates(std::map<int, RateVector>& cellRates_) const;
293
294 Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const;
295
296 // TODO: theoretically, it should be a const function
297 // Simulator is not const is because that assembleWellEq is non-const Simulator
298 void wellTesting(const Simulator& simulator,
299 const double simulation_time,
300 const WellGroupHelperType& wgHelper,
301 WellStateType& well_state,
302 WellTestState& welltest_state,
303 GLiftEclWells& ecl_well_map,
304 std::map<std::string, double>& open_times,
305 DeferredLogger& deferred_logger);
306
307 void checkWellOperability(const Simulator& simulator,
308 const WellStateType& well_state,
309 const WellGroupHelperType& wgHelper,
310 DeferredLogger& deferred_logger);
311
313 WellStateType& well_state,
314 const GroupState<Scalar>& group_state,
315 GLiftEclWells& ecl_well_map,
316 DeferredLogger& deferred_logger);
317
318 // check whether the well is operable under the current reservoir condition
319 // mostly related to BHP limit and THP limit
320 void updateWellOperability(const Simulator& simulator,
321 const WellStateType& well_state,
322 const WellGroupHelperType& wgHelper,
323 DeferredLogger& deferred_logger);
324
325 bool updateWellOperabilityFromWellEq(const Simulator& simulator,
326 const WellGroupHelperType& wgHelper,
327 DeferredLogger& deferred_logger);
328
329 // update perforation water throughput based on solved water rate
330 virtual void updateWaterThroughput(const double dt,
331 WellStateType& well_state) const = 0;
332
335 virtual std::vector<Scalar>
337 DeferredLogger& deferred_logger) const = 0;
338
342 void initializeProducerWellState(const Simulator& simulator,
343 WellStateType& well_state,
344 DeferredLogger& deferred_logger) const;
345
346 void solveWellEquation(const Simulator& simulator,
347 const WellGroupHelperType& wgHelper,
348 WellStateType& well_state,
349 DeferredLogger& deferred_logger);
350
351 const std::vector<RateVector>& connectionRates() const
352 {
353 return connectionRates_;
354 }
355
356 std::vector<Scalar> wellIndex(const int perf,
357 const IntensiveQuantities& intQuants,
358 const Scalar trans_mult,
359 const SingleWellStateType& ws) const;
360
361 void updateConnectionDFactor(const Simulator& simulator,
362 SingleWellStateType& ws) const;
363
365 SingleWellStateType& ws) const;
366
367 virtual bool iterateWellEqWithSwitching(const Simulator& simulator,
368 const double dt,
369 const WellInjectionControls& inj_controls,
370 const WellProductionControls& prod_controls,
371 const WellGroupHelperType& wgHelper,
372 WellStateType& well_state,
373 DeferredLogger& deferred_logger,
374 const bool fixed_control = false,
375 const bool fixed_status = false) = 0;
376protected:
377 // simulation parameters
378 std::vector<RateVector> connectionRates_;
379 std::vector<Scalar> B_avg_;
383
384 Scalar wpolymer() const;
385 Scalar wfoam() const;
386 Scalar wsalt() const;
387 Scalar wmicrobes() const;
388 Scalar woxygen() const;
389 Scalar wurea() const;
390
391 virtual Scalar getRefDensity() const = 0;
392
393 std::vector<Scalar>
394 initialWellRateFractions(const Simulator& ebosSimulator,
395 const WellStateType& well_state) const;
396
397 // check whether the well is operable under BHP limit with current reservoir condition
398 virtual void checkOperabilityUnderBHPLimit(const WellStateType& well_state,
399 const Simulator& simulator,
400 DeferredLogger& deferred_logger) = 0;
401
402 // check whether the well is operable under THP limit with current reservoir condition
403 virtual void checkOperabilityUnderTHPLimit(const Simulator& simulator,
404 const WellStateType& well_state,
405 const WellGroupHelperType& wgHelper,
406 DeferredLogger& deferred_logger) = 0;
407
408 virtual void updateIPR(const Simulator& simulator,
409 DeferredLogger& deferred_logger) const = 0;
410
411 virtual void assembleWellEqWithoutIteration(const Simulator& simulator,
412 const double dt,
413 const WellInjectionControls& inj_controls,
414 const WellProductionControls& prod_controls,
415 WellStateType& well_state,
416 const GroupState<Scalar>& group_state,
417 DeferredLogger& deferred_logger) = 0;
418
419 // iterate well equations with the specified control until converged
420 virtual bool iterateWellEqWithControl(const Simulator& simulator,
421 const double dt,
422 const WellInjectionControls& inj_controls,
423 const WellProductionControls& prod_controls,
424 const WellGroupHelperType& wgHelper,
425 WellStateType& well_state,
426 DeferredLogger& deferred_logger) = 0;
427
428 virtual void updateIPRImplicit(const Simulator& simulator,
429 WellStateType& well_state,
430 DeferredLogger& deferred_logger) = 0;
431
432 bool iterateWellEquations(const Simulator& simulator,
433 const double dt,
434 const WellGroupHelperType& wgHelper,
435 WellStateType& well_state,
436 DeferredLogger& deferred_logger);
437
438 bool solveWellWithOperabilityCheck(const Simulator& simulator,
439 const double dt,
440 const Well::InjectionControls& inj_controls,
441 const Well::ProductionControls& prod_controls,
442 const WellGroupHelperType& wgHelper,
443 WellStateType& well_state,
444 DeferredLogger& deferred_logger);
445
446 std::optional<Scalar>
447 estimateOperableBhp(const Simulator& ebos_simulator,
448 const double dt,
449 const WellGroupHelperType& wgHelper,
450 const SummaryState& summary_state,
451 WellStateType& well_state,
452 DeferredLogger& deferred_logger);
453
454 bool solveWellWithBhp(const Simulator& simulator,
455 const double dt,
456 const Scalar bhp,
457 const WellGroupHelperType& wgHelper,
458 WellStateType& well_state,
459 DeferredLogger& deferred_logger);
460
461 bool solveWellWithZeroRate(const Simulator& simulator,
462 const double dt,
463 const WellGroupHelperType& wgHelper,
464 WellStateType& well_state,
465 DeferredLogger& deferred_logger);
466
467 bool solveWellForTesting(const Simulator& simulator,
468 const WellGroupHelperType& wgHelper,
469 WellStateType& well_state,
470 DeferredLogger& deferred_logger);
471
472
473 template<class GasLiftSingleWell>
474 std::unique_ptr<GasLiftSingleWell> initializeGliftWellTest_(const Simulator& simulator,
475 WellStateType& well_state,
476 const GroupState<Scalar>& group_state,
477 GLiftEclWells& ecl_well_map,
478 DeferredLogger& deferred_logger);
479
480 Eval getPerfCellPressure(const FluidState& fs) const;
481
482 // get the mobility for specific perforation
483 template<class Value, class Callback>
484 void getMobility(const Simulator& simulator,
485 const int local_perf_index,
486 std::vector<Value>& mob,
487 Callback& extendEval,
488 [[maybe_unused]] DeferredLogger& deferred_logger) const;
489
490 void computeConnLevelProdInd(const FluidState& fs,
491 const std::function<Scalar(const Scalar)>& connPICalc,
492 const std::vector<Scalar>& mobility,
493 Scalar* connPI) const;
494
495 void computeConnLevelInjInd(const FluidState& fs,
496 const Phase preferred_phase,
497 const std::function<Scalar(const Scalar)>& connIICalc,
498 const std::vector<Scalar>& mobility,
499 Scalar* connII,
500 DeferredLogger& deferred_logger) const;
501
502 Scalar computeConnectionDFactor(const int perf,
503 const IntensiveQuantities& intQuants,
504 const SingleWellStateType& ws) const;
505};
506
507} // namespace Opm
508
509#include "WellInterface_impl.hpp"
510
511#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: 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: WellGroupHelper.hpp:51
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
const std::vector< FluidSystem::Scalar > & wellIndex() const
Definition: WellInterfaceGeneric.hpp:152
Definition: WellInterfaceIndices.hpp:34
typename WellInterfaceFluidSystem< GetPropType< TypeTag, Properties::FluidSystem > >::ModelParameters ModelParameters
Definition: WellInterfaceIndices.hpp:41
DenseAd::Evaluation< Scalar, Indices::numEq > Eval
Definition: WellInterfaceIndices.hpp:40
Definition: WellInterface.hpp:77
virtual std::vector< Scalar > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const =0
std::optional< Scalar > estimateOperableBhp(const Simulator &ebos_simulator, const double dt, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:695
Scalar woxygen() const
Definition: WellInterface_impl.hpp:165
static constexpr bool has_brine
Definition: WellInterface.hpp:122
IndividualOrGroup
Definition: WellInterface.hpp:253
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:82
bool updateWellOperabilityFromWellEq(const Simulator &simulator, const WellGroupHelperType &wgHelper, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1249
virtual void updateWellStateWithTarget(const Simulator &simulator, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1280
Scalar computeConnectionDFactor(const int perf, const IntensiveQuantities &intQuants, const SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:1963
virtual void computeWellRatesWithBhp(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const =0
virtual void addWellContributions(SparseMatrixAdapter &) const =0
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:105
const std::vector< RateVector > & connectionRates() const
Definition: WellInterface.hpp:351
BlackOilFluidState< Eval, FluidSystem, has_temperature, has_energy, Indices::compositionSwitchIdx >=0, has_watVapor, has_brine, has_saltPrecip, has_disgas_in_water, Indices::numPhases > FluidState
Definition: WellInterface.hpp:139
Scalar wfoam() const
Definition: WellInterface_impl.hpp:127
virtual void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const =0
std::vector< RateVector > connectionRates_
Definition: WellInterface.hpp:378
void computeConnLevelProdInd(const FluidState &fs, const std::function< Scalar(const Scalar)> &connPICalc, const std::vector< Scalar > &mobility, Scalar *connPI) const
Definition: WellInterface_impl.hpp:2159
bool solveWellForTesting(const Simulator &simulator, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:811
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:98
bool solveWellWithBhp(const Simulator &simulator, const double dt, const Scalar bhp, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:733
void gliftBeginTimeStepWellTestUpdateALQ(const Simulator &simulator, WellStateType &well_state, const GroupState< Scalar > &group_state, GLiftEclWells &ecl_well_map, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1149
Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const
Definition: WellInterface_impl.hpp:1101
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 checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellStateType &well_state, const WellGroupHelperType &wgHelper, DeferredLogger &deferred_logger)=0
bool iterateWellEquations(const Simulator &simulator, const double dt, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:513
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:2049
void solveWellEquation(const Simulator &simulator, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:863
void assembleWellEq(const Simulator &simulator, const double dt, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:925
void assembleWellEqWithoutIteration(const Simulator &simulator, const double dt, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:942
virtual void apply(const BVector &x, BVector &Ax) const =0
Ax = Ax - C D^-1 B x.
static const bool has_temperature
Definition: WellInterface.hpp:117
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:87
virtual void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger)=0
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
std::vector< Scalar > initialWellRateFractions(const Simulator &ebosSimulator, const WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1725
static constexpr bool has_watVapor
Definition: WellInterface.hpp:123
void updateConnectionDFactor(const Simulator &simulator, SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:1943
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:2034
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:97
void initializeProducerWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1775
std::vector< Scalar > B_avg_
Definition: WellInterface.hpp:379
static constexpr bool has_bioeffects
Definition: WellInterface.hpp:126
bool changed_to_stopped_this_step_
Definition: WellInterface.hpp:380
static constexpr bool has_polymermw
Definition: WellInterface.hpp:120
virtual void calculateExplicitQuantities(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger)=0
virtual bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const WellInjectionControls &inj_controls, const WellProductionControls &prod_controls, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger)=0
void addCellRates(std::map< int, RateVector > &cellRates_) const
Definition: WellInterface_impl.hpp:1082
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
bool wellUnderZeroGroupRateTarget(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger, std::optional< bool > group_control=std::nullopt) const
Definition: WellInterface_impl.hpp:1693
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
virtual void assembleWellEqWithoutIteration(const Simulator &simulator, const double dt, const WellInjectionControls &inj_controls, const WellProductionControls &prod_controls, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)=0
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
WellGroupHelper< Scalar, IndexTraits > WellGroupHelperType
Definition: WellInterface.hpp:102
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:2005
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:2193
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:2227
virtual void addWellPressureEquations(PressureMatrix &mat, const BVector &x, const int pressureVarIndex, const bool use_well_weights, const WellStateType &well_state) const =0
virtual void updateIPRImplicit(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger)=0
int number_of_well_reopenings_
Definition: WellInterface.hpp:382
virtual void solveEqAndUpdateWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger)=0
virtual void computeWellRatesWithBhpIterations(const Simulator &simulator, const Scalar &bhp, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const =0
bool updateWellStateWithTHPTargetProd(const Simulator &simulator, WellStateType &well_state, const WellGroupHelperType &wgHelper, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:2127
Scalar wsalt() const
Definition: WellInterface_impl.hpp:141
static constexpr bool has_micp
Definition: WellInterface.hpp:127
bool updateWellControlAndStatusLocalIteration(const Simulator &simulator, const WellGroupHelperType &wgHelper, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const Scalar WQTotal, WellStateType &well_state, DeferredLogger &deferred_logger, const bool fixed_control=false, const bool fixed_status=false)
Definition: WellInterface_impl.hpp:277
Dune::FieldMatrix< Scalar, Indices::numEq, Indices::numEq > MatrixBlockType
Definition: WellInterface.hpp:95
virtual std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &ebos_simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger, bool iterate_if_no_solution) const =0
typename Base::Eval Eval
Definition: WellInterface.hpp:96
bool solveWellWithZeroRate(const Simulator &simulator, const double dt, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:780
bool wellUnderZeroRateTarget(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1675
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
static constexpr bool has_energy
Definition: WellInterface.hpp:116
Scalar wpolymer() const
Definition: WellInterface_impl.hpp:111
SingleWellState< Scalar, IndexTraits > SingleWellStateType
Definition: WellInterface.hpp:101
bool updateWellControl(const Simulator &simulator, const IndividualOrGroup iog, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:189
virtual void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const =0
void wellTesting(const Simulator &simulator, const double simulation_time, const WellGroupHelperType &wgHelper, WellStateType &well_state, WellTestState &welltest_state, GLiftEclWells &ecl_well_map, std::map< std::string, double > &open_times, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:387
virtual bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const WellInjectionControls &inj_controls, const WellProductionControls &prod_controls, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger, const bool fixed_control=false, const bool fixed_status=false)=0
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:86
bool solveWellWithOperabilityCheck(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:580
void checkWellOperability(const Simulator &simulator, const WellStateType &well_state, const WellGroupHelperType &wgHelper, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1121
Scalar wurea() const
Definition: WellInterface_impl.hpp:177
void updateWellOperability(const Simulator &simulator, const WellStateType &well_state, const WellGroupHelperType &wgHelper, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1211
bool thp_update_iterations
Definition: WellInterface.hpp:381
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
virtual void scaleSegmentRatesAndPressure(WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1272
void prepareWellBeforeAssembling(const Simulator &simulator, const double dt, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:962
static constexpr bool has_zFraction
Definition: WellInterface.hpp:114
bool stoppedOrZeroRateTarget(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1712
Definition: WellProdIndexCalculator.hpp:37
Definition: WellState.hpp:66
Declares the properties required by the black oil model.
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