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
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:
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
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 using Base::name;
112 using Base::Water;
113 using Base::Oil;
114 using Base::Gas;
115
116 using typename Base::BVector;
117
118 using Eval = typename StdWellEval::Eval;
121
122 StandardWell(const Well& well,
123 const ParallelWellInfo<Scalar>& pw_info,
124 const int time_step,
125 const ModelParameters& param,
126 const RateConverterType& rate_converter,
127 const int pvtRegionIdx,
128 const int num_components,
129 const int num_phases,
130 const int index_of_well,
131 const std::vector<PerforationData<Scalar>>& perf_data);
132
133 virtual void init(const PhaseUsage* phase_usage_arg,
134 const std::vector<Scalar>& depth_arg,
135 const Scalar gravity_arg,
136 const int num_cells,
137 const std::vector<Scalar>& B_avg,
138 const bool changed_to_open_this_step) override;
139
140 void initPrimaryVariablesEvaluation() override;
141
143 virtual ConvergenceReport getWellConvergence(const Simulator& simulator,
144 const WellState<Scalar>& well_state,
145 const std::vector<Scalar>& B_avg,
146 DeferredLogger& deferred_logger,
147 const bool relax_tolerance) const override;
148
150 virtual void apply(const BVector& x, BVector& Ax) const override;
152 virtual void apply(BVector& r) const override;
153
157 const BVector& x,
158 WellState<Scalar>& well_state,
159 DeferredLogger& deferred_logger) override;
160
162 void computeWellPotentials(const Simulator& simulator,
163 const WellState<Scalar>& well_state,
164 std::vector<Scalar>& well_potentials,
165 DeferredLogger& deferred_logger) /* const */ override;
166
167 void updatePrimaryVariables(const Simulator& simulator,
168 const WellState<Scalar>& well_state,
169 DeferredLogger& deferred_logger) override;
170
171 void solveEqAndUpdateWellState(const Simulator& simulator,
172 WellState<Scalar>& well_state,
173 DeferredLogger& deferred_logger) override;
174
175 void calculateExplicitQuantities(const Simulator& simulator,
176 const WellState<Scalar>& well_state,
177 DeferredLogger& deferred_logger) override; // should be const?
178
179 void updateProductivityIndex(const Simulator& simulator,
180 const WellProdIndexCalculator<Scalar>& wellPICalc,
181 WellState<Scalar>& well_state,
182 DeferredLogger& deferred_logger) const override;
183
184 Scalar connectionDensity(const int globalConnIdx,
185 const int openConnIdx) const override;
186
187 void addWellContributions(SparseMatrixAdapter& mat) const override;
188
190 const BVector& x,
191 const int pressureVarIndex,
192 const bool use_well_weights,
193 const WellState<Scalar>& well_state) const override;
194
195 // iterate well equations with the specified control until converged
196 bool iterateWellEqWithControl(const Simulator& simulator,
197 const double dt,
198 const Well::InjectionControls& inj_controls,
199 const Well::ProductionControls& prod_controls,
200 WellState<Scalar>& well_state,
201 const GroupState<Scalar>& group_state,
202 DeferredLogger& deferred_logger) override;
203
204 // iterate well equations including control switching
205 bool iterateWellEqWithSwitching(const Simulator& simulator,
206 const double dt,
207 const Well::InjectionControls& inj_controls,
208 const Well::ProductionControls& prod_controls,
209 WellState<Scalar>& well_state,
210 const GroupState<Scalar>& group_state,
211 DeferredLogger& deferred_logger,
212 const bool fixed_control = false,
213 const bool fixed_status = false) override;
214
217 {
219 }
220
221 /* returns BHP */
223 const SummaryState &summary_state,
224 DeferredLogger& deferred_logger,
225 std::vector<Scalar>& potentials,
226 Scalar alq) const;
227
228 void computeWellRatesWithThpAlqProd(const Simulator& ebos_simulator,
229 const SummaryState& summary_state,
230 DeferredLogger& deferred_logger,
231 std::vector<Scalar>& potentials,
232 Scalar alq) const;
233
234 std::optional<Scalar>
235 computeBhpAtThpLimitProdWithAlq(const Simulator& ebos_simulator,
236 const SummaryState& summary_state,
237 const Scalar alq_value,
238 DeferredLogger& deferred_logger) const override;
239
240 void updateIPRImplicit(const Simulator& simulator,
241 WellState<Scalar>& well_state,
242 DeferredLogger& deferred_logger) override;
243
244 void computeWellRatesWithBhp(const Simulator& ebosSimulator,
245 const Scalar& bhp,
246 std::vector<Scalar>& well_flux,
247 DeferredLogger& deferred_logger) const override;
248
249 // NOTE: These cannot be protected since they are used by GasLiftRuntime
250 using Base::phaseUsage;
252
253 std::vector<Scalar>
254 computeCurrentWellRates(const Simulator& ebosSimulator,
255 DeferredLogger& deferred_logger) const override;
256
257 std::vector<Scalar> getPrimaryVars() const override;
258
259 int setPrimaryVars(typename std::vector<Scalar>::const_iterator it) override;
260
261 protected:
263
264 // updating the well_state based on well solution dwells
265 void updateWellState(const Simulator& simulator,
266 const BVectorWell& dwells,
267 WellState<Scalar>& well_state,
268 DeferredLogger& deferred_logger);
269
270 using WellConnectionProps = typename StdWellEval::StdWellConnections::Properties;
271
272 // Compute connection level PVT properties needed to calulate the
273 // pressure difference between well connections.
276 const WellState<Scalar>& well_state) 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<Scalar>& 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<Scalar>& perf_rates,
314 DeferredLogger& deferred_logger) const;
315
316 void computeWellRatesWithBhpIterations(const Simulator& ebosSimulator,
317 const Scalar& bhp,
318 std::vector<Scalar>& well_flux,
319 DeferredLogger& deferred_logger) const override;
320
321 std::vector<Scalar>
322 computeWellPotentialWithTHP(const Simulator& ebosSimulator,
323 DeferredLogger& deferred_logger,
324 const WellState<Scalar>& well_state) const;
325
326 bool computeWellPotentialsImplicit(const Simulator& ebos_simulator,
327 std::vector<Scalar>& well_potentials,
328 DeferredLogger& deferred_logger) const;
329
330 Scalar 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 Scalar throughput,
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 Scalar 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 Scalar 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<Scalar>& 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 std::optional<Scalar>
456 const Simulator& simulator,
457 const SummaryState& summary_state,
458 DeferredLogger& deferred_logger) const;
459
460 std::optional<Scalar>
461 computeBhpAtThpLimitInj(const Simulator& simulator,
462 const SummaryState& summary_state,
463 DeferredLogger& deferred_logger) const;
464
465 private:
466 Eval connectionRateEnergy(const Scalar maxOilSaturation,
467 const std::vector<EvalWell>& cq_s,
468 const IntensiveQuantities& intQuants,
469 DeferredLogger& deferred_logger) const;
470
471 template<class Value>
472 void gasOilPerfRateInj(const std::vector<Value>& cq_s,
473 PerforationRates<Scalar>& perf_rates,
474 const Value& rv,
475 const Value& rs,
476 const Value& pressure,
477 const Value& rvw,
478 DeferredLogger& deferred_logger) const;
479
480 template<class Value>
481 void gasOilPerfRateProd(std::vector<Value>& cq_s,
482 PerforationRates<Scalar>& perf_rates,
483 const Value& rv,
484 const Value& rs,
485 const Value& rvw) const;
486
487 template<class Value>
488 void gasWaterPerfRateProd(std::vector<Value>& cq_s,
489 PerforationRates<Scalar>& perf_rates,
490 const Value& rvw,
491 const Value& rsw) const;
492
493 template<class Value>
494 void gasWaterPerfRateInj(const std::vector<Value>& cq_s,
495 PerforationRates<Scalar>& perf_rates,
496 const Value& rvw,
497 const Value& rsw,
498 const Value& pressure,
499 DeferredLogger& deferred_logger) const;
500
501 template<class Value>
502 void disOilVapWatVolumeRatio(Value& volumeRatio,
503 const Value& rvw,
504 const Value& rsw,
505 const Value& pressure,
506 const std::vector<Value>& cmix_s,
507 const std::vector<Value>& b_perfcells_dense,
508 DeferredLogger& deferred_logger) const;
509
510 template<class Value>
511 void gasOilVolumeRatio(Value& volumeRatio,
512 const Value& rv,
513 const Value& rs,
514 const Value& pressure,
515 const std::vector<Value>& cmix_s,
516 const std::vector<Value>& b_perfcells_dense,
517 DeferredLogger& deferred_logger) const;
518 };
519
520}
521
522#ifndef OPM_STANDARDWELL_IMPL_HEADER_INCLUDED
523#include "StandardWell_impl.hpp"
524#endif
525
526#endif // OPM_STANDARDWELL_HEADER_INCLUDED
Contains the classes required to extend the black-oil model by brine.
Contains the classes required to extend the black-oil model by solvent component. For details,...
Contains the classes required to extend the black-oil model to include the effects of foam.
Contains the classes required to extend the black-oil model by MICP.
Contains the classes required to extend the black-oil model by polymer.
Contains the classes required to extend the black-oil model by solvents.
Contains the high level supplements required to extend the black oil model by brine.
Definition: blackoilbrinemodules.hh:50
Contains the high level supplements required to extend the black oil model to include the effects of ...
Definition: blackoilfoammodules.hh:60
Contains the high level supplements required to extend the black oil model by polymer.
Definition: blackoilpolymermodules.hh:54
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
Definition: GroupState.hpp:38
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:186
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:92
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:1365
EvalWell wpolymermw(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2017
void computeWellPotentials(const Simulator &simulator, const WellState< Scalar > &well_state, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) override
computing the well potentials for group control
Definition: StandardWell_impl.hpp:1742
typename StdWellEval::EvalWell EvalWell
Definition: StandardWell.hpp:119
static constexpr int numStaticWellEq
Definition: StandardWell.hpp:100
bool regularize_
Definition: StandardWell.hpp:262
std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &ebos_simulator, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:2224
typename StdWellEval::BVectorWell BVectorWell
Definition: StandardWell.hpp:120
void checkOperabilityUnderBHPLimit(const WellState< Scalar > &well_state, const Simulator &simulator, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:937
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:2334
void computeWellRatesWithThpAlqProd(const Simulator &ebos_simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1725
void addWellContributions(SparseMatrixAdapter &mat) const override
Definition: StandardWell_impl.hpp:1912
std::vector< Scalar > getPrimaryVars() const override
Definition: StandardWell_impl.hpp:2531
void initPrimaryVariablesEvaluation() override
Definition: StandardWell_impl.hpp:113
void updateIPRImplicit(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:861
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1239
void updateWaterMobilityWithPolymer(const Simulator &simulator, const int perf, std::vector< EvalWell > &mob_water, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1848
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:2380
StandardWell(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_components, const int num_phases, const int index_of_well, const std::vector< PerforationData< Scalar > > &perf_data)
Definition: StandardWell_impl.hpp:72
void updateWaterThroughput(const double dt, WellState< Scalar > &well_state) const override
Definition: StandardWell_impl.hpp:2046
std::vector< Scalar > computeCurrentWellRates(const Simulator &ebosSimulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:2495
virtual void init(const PhaseUsage *phase_usage_arg, const std::vector< Scalar > &depth_arg, const Scalar 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 updatePrimaryVariablesNewton(const BVectorWell &dwells, const bool stop_or_zero_rate_target, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:725
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:79
void calculateExplicitQuantities(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1405
void solveEqAndUpdateWellState(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1383
bool computeWellPotentialsImplicit(const Simulator &ebos_simulator, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1618
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1452
bool jacobianContainsWellContributions() const override
Wether the Jacobian will also have well contributions in it.
Definition: StandardWell.hpp:216
typename StdWellEval::StdWellConnections::Properties WellConnectionProps
Definition: StandardWell.hpp:270
void computeWellRatesWithBhp(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1472
WellConnectionProps computePropertiesForWellConnectionPressures(const Simulator &simulator, const WellState< Scalar > &well_state) const
Definition: StandardWell_impl.hpp:1144
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:647
void computeWellRatesWithBhpIterations(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1518
void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:768
std::vector< Scalar > computeWellPotentialWithTHP(const Simulator &ebosSimulator, DeferredLogger &deferred_logger, const WellState< Scalar > &well_state) const
Definition: StandardWell_impl.hpp:1583
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< Scalar > &perf_rates, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:126
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:703
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:344
void computeWellConnectionDensitesPressures(const Simulator &simulator, const WellState< Scalar > &well_state, const WellConnectionProps &props, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1314
virtual void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: StandardWell_impl.hpp:1420
std::optional< Scalar > computeBhpAtThpLimitProd(const WellState< Scalar > &well_state, const Simulator &simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2210
void checkConvergenceExtraEqs(const std::vector< Scalar > &res, ConvergenceReport &report) const
Definition: StandardWell_impl.hpp:2148
typename StdWellEval::Eval Eval
Definition: StandardWell.hpp:118
Scalar computeWellRatesAndBhpWithThpAlqProd(const Simulator &ebos_simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1696
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:2167
static constexpr int Bhp
Definition: StandardWell.hpp:105
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1133
void updatePrimaryVariables(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1814
virtual ConvergenceReport getWellConvergence(const Simulator &simulator, const WellState< Scalar > &well_state, const std::vector< Scalar > &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:1195
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2302
bool allDrawDownWrongDirection(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1055
EvalWell pskin(const Scalar throughput, const EvalWell &water_velocity, const EvalWell &poly_inj_conc, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1972
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:748
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:2100
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: StandardWell_impl.hpp:2548
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:370
void checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1006
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:1920
EvalWell pskinwater(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1940
bool canProduceInjectWithCurrentBhp(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1097
void handleInjectivityRate(const Simulator &simulator, const int perf, std::vector< EvalWell > &cq_s) const
Definition: StandardWell_impl.hpp:2077
static constexpr int numWellControlEq
Definition: StandardWell.hpp:97
Scalar getRefDensity() const override
Definition: StandardWell_impl.hpp:1837
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:497
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: StandardWell_impl.hpp:1799
static constexpr int Oil
Definition: WellInterfaceFluidSystem.hpp:64
static constexpr int Water
Definition: WellInterfaceFluidSystem.hpp:63
static constexpr int Gas
Definition: WellInterfaceFluidSystem.hpp:65
const std::string & name() const
Well name.
int pvtRegionIdx() const
Definition: WellInterfaceGeneric.hpp:119
const VFPProperties< FluidSystem::Scalar > * vfp_properties_
Definition: WellInterfaceGeneric.hpp:362
Definition: WellInterface.hpp:73
const ModelParameters & param_
Definition: WellInterface.hpp:388
static constexpr bool has_brine
Definition: WellInterface.hpp:116
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:78
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:100
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:132
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:97
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:82
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:79
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: WellInterface.hpp:83
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:96
static constexpr bool has_polymermw
Definition: WellInterface.hpp:114
static constexpr bool has_polymer
Definition: WellInterface.hpp:110
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:80
static constexpr bool has_solvent
Definition: WellInterface.hpp:108
static constexpr bool has_foam
Definition: WellInterface.hpp:115
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: WellInterface.hpp:85
static constexpr bool has_micp
Definition: WellInterface.hpp:120
BlackoilModelParameters< Scalar > ModelParameters
Definition: WellInterface.hpp:106
static constexpr bool has_energy
Definition: WellInterface.hpp:111
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:81
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:84
static constexpr bool has_zFraction
Definition: WellInterface.hpp:109
Definition: WellProdIndexCalculator.hpp:37
Definition: WellState.hpp:62
Definition: blackoilboundaryratevector.hh:37
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:235
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:147
bool matrix_add_well_contributions_
Whether to add influences of wells between cells to the matrix and preconditioner matrix.
Definition: BlackoilModelParameters.hpp:244
Static data associated with a well perforation.
Definition: PerforationData.hpp:30
Definition: PerforationData.hpp:41
Definition: BlackoilPhases.hpp:46