StandardWell.hpp
Go to the documentation of this file.
1/*
2 Copyright 2017 SINTEF Digital, Mathematics and Cybernetics.
3 Copyright 2017 Statoil ASA.
4 Copyright 2016 - 2017 IRIS AS.
5
6 This file is part of the Open Porous Media project (OPM).
7
8 OPM is free software: you can redistribute it and/or modify
9 it under the terms of the GNU General Public License as published by
10 the Free Software Foundation, either version 3 of the License, or
11 (at your option) any later version.
12
13 OPM is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License
19 along with OPM. If not, see <http://www.gnu.org/licenses/>.
20*/
21
22
23#ifndef OPM_STANDARDWELL_HEADER_INCLUDED
24#define OPM_STANDARDWELL_HEADER_INCLUDED
25
33
34#include <opm/models/blackoil/blackoilpolymermodules.hh>
35#include <opm/models/blackoil/blackoilsolventmodules.hh>
36#include <opm/models/blackoil/blackoilextbomodules.hh>
37#include <opm/models/blackoil/blackoilfoammodules.hh>
38#include <opm/models/blackoil/blackoilbrinemodules.hh>
39#include <opm/models/blackoil/blackoilmicpmodules.hh>
40
41#include <opm/material/densead/Evaluation.hpp>
42#include <opm/input/eclipse/Schedule/ScheduleTypes.hpp>
43
45
46#include <dune/common/dynvector.hh>
47#include <dune/common/dynmatrix.hh>
48
49#include <memory>
50#include <optional>
51
52namespace Opm
53{
54
55 template<typename TypeTag>
56 class StandardWell : public WellInterface<TypeTag>
57 , public StandardWellEval<GetPropType<TypeTag, Properties::FluidSystem>,
58 GetPropType<TypeTag, Properties::Indices>>
59 {
60
61 public:
64 GetPropType<TypeTag, Properties::Indices>>;
65
66 // TODO: some functions working with AD variables handles only with values (double) without
67 // dealing with derivatives. It can be beneficial to make functions can work with either AD or scalar value.
68 // And also, it can also be beneficial to make these functions hanle different types of AD variables.
69 using typename Base::Simulator;
70 using typename Base::IntensiveQuantities;
71 using typename Base::FluidSystem;
72 using typename Base::MaterialLaw;
73 using typename Base::ModelParameters;
74 using typename Base::Indices;
75 using typename Base::RateConverterType;
76 using typename Base::SparseMatrixAdapter;
77 using typename Base::FluidState;
78 using typename Base::RateVector;
79
84 using Base::has_foam;
85 using Base::has_brine;
86 using Base::has_energy;
87 using Base::has_micp;
88
89 using PolymerModule = BlackOilPolymerModule<TypeTag>;
90 using FoamModule = BlackOilFoamModule<TypeTag>;
91 using BrineModule = BlackOilBrineModule<TypeTag>;
92 using typename Base::PressureMatrix;
93
94 // number of the conservation equations
95 static constexpr int numWellConservationEq = Indices::numPhases + Indices::numSolvents;
96 // number of the well control equations
97 static constexpr int numWellControlEq = 1;
98 // number of the well equations that will always be used
99 // based on the solution strategy, there might be other well equations be introduced
101
102 // the index for Bhp in primary variables and also the index of well control equation
103 // they both will be the last one in their respective system.
104 // TODO: we should have indices for the well equations and well primary variables separately
105 static constexpr int Bhp = numStaticWellEq - numWellControlEq;
106
108
109 using typename Base::Scalar;
110
111
112 using Base::name;
113 using Base::Water;
114 using Base::Oil;
115 using Base::Gas;
116
117 using typename Base::BVector;
118
119 using Eval = typename StdWellEval::Eval;
122
123 StandardWell(const Well& well,
124 const ParallelWellInfo& pw_info,
125 const int time_step,
126 const ModelParameters& param,
127 const RateConverterType& rate_converter,
128 const int pvtRegionIdx,
129 const int num_components,
130 const int num_phases,
131 const int index_of_well,
132 const std::vector<PerforationData>& perf_data);
133
134 virtual void init(const PhaseUsage* phase_usage_arg,
135 const std::vector<double>& depth_arg,
136 const double gravity_arg,
137 const int num_cells,
138 const std::vector< Scalar >& B_avg,
139 const bool changed_to_open_this_step) override;
140
141
142 void initPrimaryVariablesEvaluation() override;
143
145 virtual ConvergenceReport getWellConvergence(const SummaryState& summary_state,
146 const WellState<Scalar>& well_state,
147 const std::vector<double>& B_avg,
148 DeferredLogger& deferred_logger,
149 const bool relax_tolerance) const override;
150
152 virtual void apply(const BVector& x, BVector& Ax) const override;
154 virtual void apply(BVector& r) const override;
155
158 void recoverWellSolutionAndUpdateWellState(const SummaryState& summary_state,
159 const BVector& x,
160 WellState<Scalar>& well_state,
161 DeferredLogger& deferred_logger) override;
162
164 void computeWellPotentials(const Simulator& simulator,
165 const WellState<Scalar>& well_state,
166 std::vector<double>& well_potentials,
167 DeferredLogger& deferred_logger) /* const */ override;
168
169 void updatePrimaryVariables(const SummaryState& summary_state,
170 const WellState<Scalar>& well_state,
171 DeferredLogger& deferred_logger) override;
172
173 void solveEqAndUpdateWellState(const SummaryState& summary_state,
174 WellState<Scalar>& well_state,
175 DeferredLogger& deferred_logger) override;
176
177 void calculateExplicitQuantities(const Simulator& simulator,
178 const WellState<Scalar>& well_state,
179 DeferredLogger& deferred_logger) override; // should be const?
180
181 void updateProductivityIndex(const Simulator& simulator,
182 const WellProdIndexCalculator& wellPICalc,
183 WellState<Scalar>& well_state,
184 DeferredLogger& deferred_logger) const override;
185
186 double connectionDensity(const int globalConnIdx,
187 const int openConnIdx) const override;
188
189 void addWellContributions(SparseMatrixAdapter& mat) const override;
190
192 const BVector& x,
193 const int pressureVarIndex,
194 const bool use_well_weights,
195 const WellState<Scalar>& well_state) const override;
196
197 // iterate well equations with the specified control until converged
198 bool iterateWellEqWithControl(const Simulator& simulator,
199 const double dt,
200 const Well::InjectionControls& inj_controls,
201 const Well::ProductionControls& prod_controls,
202 WellState<Scalar>& well_state,
203 const GroupState<Scalar>& group_state,
204 DeferredLogger& deferred_logger) override;
205
206 // iterate well equations including control switching
207 bool iterateWellEqWithSwitching(const Simulator& simulator,
208 const double dt,
209 const Well::InjectionControls& inj_controls,
210 const Well::ProductionControls& prod_controls,
211 WellState<Scalar>& well_state,
212 const GroupState<Scalar>& group_state,
213 DeferredLogger& deferred_logger,
214 const bool fixed_control = false,
215 const bool fixed_status = false) override;
216
219 {
221 }
222
223 /* returns BHP */
224 double computeWellRatesAndBhpWithThpAlqProd(const Simulator& simulator,
225 const SummaryState& summary_state,
226 DeferredLogger& deferred_logger,
227 std::vector<double>& potentials,
228 double alq) const;
229
230 void computeWellRatesWithThpAlqProd(const Simulator& simulator,
231 const SummaryState& summary_state,
232 DeferredLogger& deferred_logger,
233 std::vector<double>& potentials,
234 double alq) const;
235
236 std::optional<double>
238 const SummaryState& summary_state,
239 const double alq_value,
240 DeferredLogger& deferred_logger) const override;
241
242 void updateIPRImplicit(const Simulator& simulator,
243 WellState<Scalar>& well_state,
244 DeferredLogger& deferred_logger) override;
245
246 void computeWellRatesWithBhp(const Simulator& simulator,
247 const double& bhp,
248 std::vector<double>& well_flux,
249 DeferredLogger& deferred_logger) const override;
250
251 // NOTE: These cannot be protected since they are used by GasLiftRuntime
252 using Base::phaseUsage;
254
255 std::vector<double> computeCurrentWellRates(const Simulator& simulator,
256 DeferredLogger& deferred_logger) const override;
257
258 std::vector<double> getPrimaryVars() const override;
259
260 int setPrimaryVars(std::vector<double>::const_iterator it) override;
261
262 protected:
264
265 // updating the well_state based on well solution dwells
266 void updateWellState(const SummaryState& summary_state,
267 const BVectorWell& dwells,
268 WellState<Scalar>& well_state,
269 DeferredLogger& deferred_logger);
270
271 // calculate the properties for the well connections
272 // to calulate the pressure difference between well connections.
273 using WellConnectionProps = typename StdWellEval::StdWellConnections::Properties;
275 const WellState<Scalar>& well_state,
276 WellConnectionProps& props) const;
277
279 const WellState<Scalar>& well_state,
280 const WellConnectionProps& props,
281 DeferredLogger& deferred_logger);
282
283 void computeWellConnectionPressures(const Simulator& simulator,
284 const WellState<Scalar>& well_state,
285 DeferredLogger& deferred_logger);
286
287 template<class Value>
288 void computePerfRate(const IntensiveQuantities& intQuants,
289 const std::vector<Value>& mob,
290 const Value& bhp,
291 const std::vector<Scalar>& Tw,
292 const int perf,
293 const bool allow_cf,
294 std::vector<Value>& cq_s,
295 PerforationRates& perf_rates,
296 DeferredLogger& deferred_logger) const;
297
298 template<class Value>
299 void computePerfRate(const std::vector<Value>& mob,
300 const Value& pressure,
301 const Value& bhp,
302 const Value& rs,
303 const Value& rv,
304 const Value& rvw,
305 const Value& rsw,
306 std::vector<Value>& b_perfcells_dense,
307 const std::vector<Scalar>& Tw,
308 const int perf,
309 const bool allow_cf,
310 const Value& skin_pressure,
311 const std::vector<Value>& cmix_s,
312 std::vector<Value>& cq_s,
313 PerforationRates& perf_rates,
314 DeferredLogger& deferred_logger) const;
315
316 void computeWellRatesWithBhpIterations(const Simulator& simulator,
317 const double& bhp,
318 std::vector<double>& well_flux,
319 DeferredLogger& deferred_logger) const override;
320
321 std::vector<double> computeWellPotentialWithTHP(
322 const Simulator& simulator,
323 DeferredLogger& deferred_logger,
324 const WellState<Scalar>& well_state) const;
325
326 bool computeWellPotentialsImplicit(const Simulator& simulator,
327 std::vector<double>& well_potentials,
328 DeferredLogger& deferred_logger) const;
329
330 double getRefDensity() const override;
331
332 // get the mobility for specific perforation
333 template<class Value>
334 void getMobility(const Simulator& simulator,
335 const int perf,
336 std::vector<Value>& mob,
337 DeferredLogger& deferred_logger) const;
338
339 void updateWaterMobilityWithPolymer(const Simulator& simulator,
340 const int perf,
341 std::vector<EvalWell>& mob_water,
342 DeferredLogger& deferred_logger) const;
343
344 void updatePrimaryVariablesNewton(const BVectorWell& dwells,
345 const bool stop_or_zero_rate_target,
346 DeferredLogger& deferred_logger);
347
348 void updateWellStateFromPrimaryVariables(const bool stop_or_zero_rate_target,
349 WellState<Scalar>& well_state,
350 const SummaryState& summary_state,
351 DeferredLogger& deferred_logger) const;
352
353 void assembleWellEqWithoutIteration(const Simulator& simulator,
354 const double dt,
355 const Well::InjectionControls& inj_controls,
356 const Well::ProductionControls& prod_controls,
357 WellState<Scalar>& well_state,
358 const GroupState<Scalar>& group_state,
359 DeferredLogger& deferred_logger) override;
360
361 void assembleWellEqWithoutIterationImpl(const Simulator& simulator,
362 const double dt,
363 const Well::InjectionControls& inj_controls,
364 const Well::ProductionControls& prod_controls,
365 WellState<Scalar>& well_state,
366 const GroupState<Scalar>& group_state,
367 DeferredLogger& deferred_logger);
368
369 void calculateSinglePerf(const Simulator& simulator,
370 const int perf,
371 WellState<Scalar>& well_state,
372 std::vector<RateVector>& connectionRates,
373 std::vector<EvalWell>& cq_s,
374 EvalWell& water_flux_s,
375 EvalWell& cq_s_zfrac_effective,
376 DeferredLogger& deferred_logger) const;
377
378 // check whether the well is operable under BHP limit with current reservoir condition
380 const Simulator& simulator,
381 DeferredLogger& deferred_logger) override;
382
383 // check whether the well is operable under THP limit with current reservoir condition
384 void checkOperabilityUnderTHPLimit(const Simulator& simulator,
385 const WellState<Scalar>& well_state,
386 DeferredLogger& deferred_logger) override;
387
388 // updating the inflow based on the current reservoir condition
389 void updateIPR(const Simulator& simulator,
390 DeferredLogger& deferred_logger) const override;
391
392 // for a well, when all drawdown are in the wrong direction, then this well will not
393 // be able to produce/inject .
394 bool allDrawDownWrongDirection(const Simulator& simulator) const;
395
396 // whether the well can produce / inject based on the current well state (bhp)
397 bool canProduceInjectWithCurrentBhp(const Simulator& simulator,
398 const WellState<Scalar>& well_state,
399 DeferredLogger& deferred_logger);
400
401 // turn on crossflow to avoid singular well equations
402 // when the well is banned from cross-flow and the BHP is not properly initialized,
403 // we turn on crossflow to avoid singular well equations. It can result in wrong-signed
404 // well rates, it can cause problem for THP calculation
405 // TODO: looking for better alternative to avoid wrong-signed well rates
406 bool openCrossFlowAvoidSingularity(const Simulator& simulator) const;
407
408 // calculate the skin pressure based on water velocity, throughput and polymer concentration.
409 // throughput is used to describe the formation damage during water/polymer injection.
410 // calculated skin pressure will be applied to the drawdown during perforation rate calculation
411 // to handle the effect from formation damage.
412 EvalWell pskin(const double throuhgput,
413 const EvalWell& water_velocity,
414 const EvalWell& poly_inj_conc,
415 DeferredLogger& deferred_logger) const;
416
417 // calculate the skin pressure based on water velocity, throughput during water injection.
418 EvalWell pskinwater(const double throughput,
419 const EvalWell& water_velocity,
420 DeferredLogger& deferred_logger) const;
421
422 // calculate the injecting polymer molecular weight based on the througput and water velocity
423 EvalWell wpolymermw(const double throughput,
424 const EvalWell& water_velocity,
425 DeferredLogger& deferred_logger) const;
426
427 // modify the water rate for polymer injectivity study
428 void handleInjectivityRate(const Simulator& simulator,
429 const int perf,
430 std::vector<EvalWell>& cq_s) const;
431
432 // handle the extra equations for polymer injectivity study
433 void handleInjectivityEquations(const Simulator& simulator,
434 const WellState<Scalar>& well_state,
435 const int perf,
436 const EvalWell& water_flux_s,
437 DeferredLogger& deferred_logger);
438
439 void updateWaterThroughput(const double dt,
440 WellState<Scalar>& well_state) const override;
441
442 // checking convergence of extra equations, if there are any
443 void checkConvergenceExtraEqs(const std::vector<double>& res,
444 ConvergenceReport& report) const;
445
446 // updating the connectionRates_ related polymer molecular weight
447 void updateConnectionRatePolyMW(const EvalWell& cq_s_poly,
448 const IntensiveQuantities& int_quants,
449 const WellState<Scalar>& well_state,
450 const int perf,
451 std::vector<RateVector>& connectionRates,
452 DeferredLogger& deferred_logger) const;
453
454
455 std::optional<double>
457 const Simulator& simulator,
458 const SummaryState& summary_state,
459 DeferredLogger& deferred_logger) const;
460
461 std::optional<double>
462 computeBhpAtThpLimitInj(const Simulator& simulator,
463 const SummaryState& summary_state,
464 DeferredLogger& deferred_logger) const;
465
466 private:
467 Eval connectionRateEnergy(const double maxOilSaturation,
468 const std::vector<EvalWell>& cq_s,
469 const IntensiveQuantities& intQuants,
470 DeferredLogger& deferred_logger) const;
471
472 template<class Value>
473 void gasOilPerfRateInj(const std::vector<Value>& cq_s,
474 PerforationRates& perf_rates,
475 const Value& rv,
476 const Value& rs,
477 const Value& pressure,
478 const Value& rvw,
479 DeferredLogger& deferred_logger) const;
480
481 template<class Value>
482 void gasOilPerfRateProd(std::vector<Value>& cq_s,
483 PerforationRates& perf_rates,
484 const Value& rv,
485 const Value& rs,
486 const Value& rvw) const;
487
488 template<class Value>
489 void gasWaterPerfRateProd(std::vector<Value>& cq_s,
490 PerforationRates& perf_rates,
491 const Value& rvw,
492 const Value& rsw) const;
493
494 template<class Value>
495 void gasWaterPerfRateInj(const std::vector<Value>& cq_s,
496 PerforationRates& perf_rates,
497 const Value& rvw,
498 const Value& rsw,
499 const Value& pressure,
500 DeferredLogger& deferred_logger) const;
501
502 template<class Value>
503 void disOilVapWatVolumeRatio(Value& volumeRatio,
504 const Value& rvw,
505 const Value& rsw,
506 const Value& pressure,
507 const std::vector<Value>& cmix_s,
508 const std::vector<Value>& b_perfcells_dense,
509 DeferredLogger& deferred_logger) const;
510
511 template<class Value>
512 void gasOilVolumeRatio(Value& volumeRatio,
513 const Value& rv,
514 const Value& rs,
515 const Value& pressure,
516 const std::vector<Value>& cmix_s,
517 const std::vector<Value>& b_perfcells_dense,
518 DeferredLogger& deferred_logger) const;
519 };
520
521}
522
523#ifndef OPM_STANDARDWELL_IMPL_HEADER_INCLUDED
524#include "StandardWell_impl.hpp"
525#endif
526
527#endif // OPM_STANDARDWELL_HEADER_INCLUDED
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
Definition: GroupState.hpp:35
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:184
Definition: StandardWellEval.hpp:47
typename StandardWellEquations< Scalar, Indices::numEq >::BVectorWell BVectorWell
Definition: StandardWellEval.hpp:65
DenseAd::Evaluation< Scalar, Indices::numEq > Eval
Definition: StandardWellEval.hpp:64
Definition: StandardWell.hpp:59
void computeWellConnectionPressures(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1340
typename StdWellEval::EvalWell EvalWell
Definition: StandardWell.hpp:120
static constexpr int numStaticWellEq
Definition: StandardWell.hpp:100
bool regularize_
Definition: StandardWell.hpp:263
typename StdWellEval::BVectorWell BVectorWell
Definition: StandardWell.hpp:121
void checkOperabilityUnderBHPLimit(const WellState< Scalar > &well_state, const Simulator &simulator, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:925
bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:2288
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:82
std::vector< double > getPrimaryVars() const override
Definition: StandardWell_impl.hpp:2482
std::optional< double > computeBhpAtThpLimitProdWithAlq(const Simulator &simulator, const SummaryState &summary_state, const double alq_value, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:2178
double connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: StandardWell_impl.hpp:1760
void updatePrimaryVariables(const SummaryState &summary_state, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1775
void addWellContributions(SparseMatrixAdapter &mat) const override
Definition: StandardWell_impl.hpp:1873
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator &wellPICalc, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1222
EvalWell pskinwater(const double throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1901
void initPrimaryVariablesEvaluation() override
Definition: StandardWell_impl.hpp:113
void checkConvergenceExtraEqs(const std::vector< double > &res, ConvergenceReport &report) const
Definition: StandardWell_impl.hpp:2102
BlackOilFoamModule< TypeTag > FoamModule
Definition: StandardWell.hpp:90
void updateIPRImplicit(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:849
std::vector< double > computeWellPotentialWithTHP(const Simulator &simulator, DeferredLogger &deferred_logger, const WellState< Scalar > &well_state) const
Definition: StandardWell_impl.hpp:1561
void computePropertiesForWellConnectionPressures(const Simulator &simulator, const WellState< Scalar > &well_state, WellConnectionProps &props) const
Definition: StandardWell_impl.hpp:1132
EvalWell pskin(const double throuhgput, const EvalWell &water_velocity, const EvalWell &poly_inj_conc, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1933
void updateWaterMobilityWithPolymer(const Simulator &simulator, const int perf, std::vector< EvalWell > &mob_water, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1809
bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger, const bool fixed_control=false, const bool fixed_status=false) override
Definition: StandardWell_impl.hpp:2335
void updateWaterThroughput(const double dt, WellState< Scalar > &well_state) const override
Definition: StandardWell_impl.hpp:2007
int setPrimaryVars(std::vector< double >::const_iterator it) override
Definition: StandardWell_impl.hpp:2499
std::optional< double > computeBhpAtThpLimitInj(const Simulator &simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2256
void updatePrimaryVariablesNewton(const BVectorWell &dwells, const bool stop_or_zero_rate_target, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:716
BlackOilBrineModule< TypeTag > BrineModule
Definition: StandardWell.hpp:91
void calculateExplicitQuantities(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1382
void computeWellRatesWithThpAlqProd(const Simulator &simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger, std::vector< double > &potentials, double alq) const
Definition: StandardWell_impl.hpp:1686
std::vector< double > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:2446
StandardWell(const Well &well, const ParallelWellInfo &pw_info, const int time_step, const ModelParameters &param, const RateConverterType &rate_converter, const int pvtRegionIdx, const int num_components, const int num_phases, const int index_of_well, const std::vector< PerforationData > &perf_data)
Definition: StandardWell_impl.hpp:72
double computeWellRatesAndBhpWithThpAlqProd(const Simulator &simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger, std::vector< double > &potentials, double alq) const
Definition: StandardWell_impl.hpp:1658
bool jacobianContainsWellContributions() const override
Wether the Jacobian will also have well contributions in it.
Definition: StandardWell.hpp:218
typename StdWellEval::StdWellConnections::Properties WellConnectionProps
Definition: StandardWell.hpp:273
double getRefDensity() const override
Definition: StandardWell_impl.hpp:1798
void updateWellState(const SummaryState &summary_state, const BVectorWell &dwells, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:695
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:639
virtual void init(const PhaseUsage *phase_usage_arg, const std::vector< double > &depth_arg, const double gravity_arg, const int num_cells, const std::vector< Scalar > &B_avg, const bool changed_to_open_this_step) override
Definition: StandardWell_impl.hpp:96
void computePerfRate(const IntensiveQuantities &intQuants, const std::vector< Value > &mob, const Value &bhp, const std::vector< Scalar > &Tw, const int perf, const bool allow_cf, std::vector< Value > &cq_s, PerforationRates &perf_rates, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:126
void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:759
void assembleWellEqWithoutIteration(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:338
void computeWellConnectionDensitesPressures(const Simulator &simulator, const WellState< Scalar > &well_state, const WellConnectionProps &props, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1298
virtual void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: StandardWell_impl.hpp:1398
typename StdWellEval::Eval Eval
Definition: StandardWell.hpp:119
void computeWellPotentials(const Simulator &simulator, const WellState< Scalar > &well_state, std::vector< double > &well_potentials, DeferredLogger &deferred_logger) override
computing the well potentials for group control
Definition: StandardWell_impl.hpp:1703
void updateConnectionRatePolyMW(const EvalWell &cq_s_poly, const IntensiveQuantities &int_quants, const WellState< Scalar > &well_state, const int perf, std::vector< RateVector > &connectionRates, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2121
static constexpr int Bhp
Definition: StandardWell.hpp:105
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1121
std::optional< double > computeBhpAtThpLimitProd(const WellState< Scalar > &well_state, const Simulator &simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2164
bool allDrawDownWrongDirection(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1043
void recoverWellSolutionAndUpdateWellState(const SummaryState &summary_state, const BVector &x, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1430
void updateWellStateFromPrimaryVariables(const bool stop_or_zero_rate_target, WellState< Scalar > &well_state, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:739
static constexpr int numWellConservationEq
Definition: StandardWell.hpp:95
void handleInjectivityEquations(const Simulator &simulator, const WellState< Scalar > &well_state, const int perf, const EvalWell &water_flux_s, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:2054
virtual ConvergenceReport getWellConvergence(const SummaryState &summary_state, const WellState< Scalar > &well_state, const std::vector< double > &B_avg, DeferredLogger &deferred_logger, const bool relax_tolerance) const override
check whether the well equations get converged for this well
Definition: StandardWell_impl.hpp:1178
void computeWellRatesWithBhp(const Simulator &simulator, const double &bhp, std::vector< double > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1450
void assembleWellEqWithoutIterationImpl(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:364
void checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:994
void addWellPressureEquations(PressureMatrix &mat, const BVector &x, const int pressureVarIndex, const bool use_well_weights, const WellState< Scalar > &well_state) const override
Definition: StandardWell_impl.hpp:1881
BlackOilPolymerModule< TypeTag > PolymerModule
Definition: StandardWell.hpp:89
bool canProduceInjectWithCurrentBhp(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1085
void handleInjectivityRate(const Simulator &simulator, const int perf, std::vector< EvalWell > &cq_s) const
Definition: StandardWell_impl.hpp:2031
EvalWell wpolymermw(const double throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1978
bool computeWellPotentialsImplicit(const Simulator &simulator, std::vector< double > &well_potentials, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1595
void solveEqAndUpdateWellState(const SummaryState &summary_state, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1360
static constexpr int numWellControlEq
Definition: StandardWell.hpp:97
void computeWellRatesWithBhpIterations(const Simulator &simulator, const double &bhp, std::vector< double > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1496
void calculateSinglePerf(const Simulator &simulator, const int perf, WellState< Scalar > &well_state, std::vector< RateVector > &connectionRates, std::vector< EvalWell > &cq_s, EvalWell &water_flux_s, EvalWell &cq_s_zfrac_effective, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:489
static constexpr int Oil
Definition: WellInterfaceFluidSystem.hpp:58
static constexpr int Water
Definition: WellInterfaceFluidSystem.hpp:57
static constexpr int Gas
Definition: WellInterfaceFluidSystem.hpp:59
int pvtRegionIdx() const
Definition: WellInterfaceGeneric.hpp:126
const PhaseUsage & phaseUsage() const
const std::string & name() const
Well name.
const VFPProperties * vfp_properties_
Definition: WellInterfaceGeneric.hpp:391
Definition: WellInterface.hpp:75
Dune::BCRSMatrix< Opm::MatrixBlock< double, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:102
const ModelParameters & param_
Definition: WellInterface.hpp:373
static constexpr bool has_brine
Definition: WellInterface.hpp:119
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:82
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:105
const std::vector< RateVector > & connectionRates() const
Definition: WellInterface.hpp:344
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:135
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:85
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:96
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: WellInterface.hpp:86
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:101
static constexpr bool has_polymermw
Definition: WellInterface.hpp:117
static constexpr bool has_polymer
Definition: WellInterface.hpp:113
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:83
static constexpr bool has_solvent
Definition: WellInterface.hpp:111
BlackoilModelParameters< TypeTag > ModelParameters
Definition: WellInterface.hpp:79
static constexpr bool has_foam
Definition: WellInterface.hpp:118
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: WellInterface.hpp:88
static constexpr bool has_micp
Definition: WellInterface.hpp:123
static constexpr bool has_energy
Definition: WellInterface.hpp:114
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:84
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:87
static constexpr bool has_zFraction
Definition: WellInterface.hpp:112
Definition: WellProdIndexCalculator.hpp:36
Definition: WellState.hpp:62
VFPEvaluation bhp(const VFPProdTable &table, const double aqua, const double liquid, const double vapour, const double thp, const double alq, const double explicit_wfr, const double explicit_gfr, const bool use_vfpexplicit)
Definition: BlackoilPhases.hpp:27
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:484
bool matrix_add_well_contributions_
Whether to add influences of wells between cells to the matrix and preconditioner matrix.
Definition: BlackoilModelParameters.hpp:576
Definition: PerforationData.hpp:40
Definition: BlackoilPhases.hpp:46