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
34
41
42#include <opm/material/densead/Evaluation.hpp>
43#include <opm/input/eclipse/Schedule/ScheduleTypes.hpp>
44
46
47#include <dune/common/dynvector.hh>
48#include <dune/common/dynmatrix.hh>
49
50#include <memory>
51#include <optional>
52
53namespace Opm
54{
55
56 template<typename TypeTag>
57 class StandardWell : public WellInterface<TypeTag>
58 , public StandardWellEval<GetPropType<TypeTag, Properties::FluidSystem>,
59 GetPropType<TypeTag, Properties::Indices>>
60 {
61
62 public:
66
67 // TODO: some functions working with AD variables handles only with values (double) without
68 // dealing with derivatives. It can be beneficial to make functions can work with either AD or scalar value.
69 // And also, it can also be beneficial to make these functions hanle different types of AD variables.
70 using typename Base::Simulator;
71 using typename Base::IntensiveQuantities;
72 using typename Base::FluidSystem;
73 using typename Base::MaterialLaw;
74 using typename Base::ModelParameters;
75 using typename Base::Indices;
76 using typename Base::RateConverterType;
77 using typename Base::SparseMatrixAdapter;
78 using typename Base::FluidState;
79 using typename Base::RateVector;
80 using typename Base::WellGroupHelperType;
81
86 using Base::has_foam;
87 using Base::has_brine;
88 using Base::has_energy;
90 using Base::has_micp;
91
94 using typename Base::PressureMatrix;
95
96 // number of the conservation equations
97 static constexpr int numWellConservationEq = Indices::numPhases + Indices::numSolvents;
98 // number of the well control equations
99 static constexpr int numWellControlEq = 1;
100 // number of the well equations that will always be used
101 // based on the solution strategy, there might be other well equations be introduced
103
104 // the index for Bhp in primary variables and also the index of well control equation
105 // they both will be the last one in their respective system.
106 // TODO: we should have indices for the well equations and well primary variables separately
107 static constexpr int Bhp = numStaticWellEq - numWellControlEq;
108
110
111 using typename Base::Scalar;
112
113 using Base::name;
114 using Base::Water;
115 using Base::Oil;
116 using Base::Gas;
117
118 using typename Base::BVector;
119
120 using Eval = typename StdWellEval::Eval;
123
124 using IndexTraits = typename FluidSystem::IndexTraitsType;
126
127 StandardWell(const Well& well,
128 const ParallelWellInfo<Scalar>& pw_info,
129 const int time_step,
130 const ModelParameters& param,
131 const RateConverterType& rate_converter,
132 const int pvtRegionIdx,
133 const int num_conservation_quantities,
134 const int num_phases,
135 const int index_of_well,
136 const std::vector<PerforationData<Scalar>>& perf_data);
137
138 virtual void init(const std::vector<Scalar>& depth_arg,
139 const Scalar gravity_arg,
140 const std::vector<Scalar>& B_avg,
141 const bool changed_to_open_this_step) override;
142
144 virtual ConvergenceReport getWellConvergence(const Simulator& simulator,
145 const WellStateType& well_state,
146 const std::vector<Scalar>& B_avg,
147 DeferredLogger& deferred_logger,
148 const bool relax_tolerance) const override;
149
151 virtual void apply(const BVector& x, BVector& Ax) const override;
153 virtual void apply(BVector& r) const override;
154
158 const BVector& x,
159 WellStateType& well_state,
160 DeferredLogger& deferred_logger) override;
161
163 void computeWellPotentials(const Simulator& simulator,
164 const WellStateType& well_state,
165 const WellGroupHelperType& wgHelper,
166 std::vector<Scalar>& well_potentials,
167 DeferredLogger& deferred_logger) /* const */ override;
168
169 void updatePrimaryVariables(const Simulator& simulator,
170 const WellStateType& well_state,
171 DeferredLogger& deferred_logger) override;
172
173 void solveEqAndUpdateWellState(const Simulator& simulator,
174 WellStateType& well_state,
175 DeferredLogger& deferred_logger) override;
176
177 void calculateExplicitQuantities(const Simulator& simulator,
178 const WellStateType& well_state,
179 DeferredLogger& deferred_logger) override; // should be const?
180
181 void updateProductivityIndex(const Simulator& simulator,
182 const WellProdIndexCalculator<Scalar>& wellPICalc,
183 WellStateType& well_state,
184 DeferredLogger& deferred_logger) const override;
185
186 Scalar 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 WellStateType& 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 const WellGroupHelperType& wgHelper,
203 WellStateType& well_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 const WellGroupHelperType& wgHelper,
212 WellStateType& well_state,
213 DeferredLogger& deferred_logger,
214 const bool fixed_control = false,
215 const bool fixed_status = false) override;
216
217 /* returns BHP */
219 const WellGroupHelperType& wgHelper,
220 const SummaryState &summary_state,
221 DeferredLogger& deferred_logger,
222 std::vector<Scalar>& potentials,
223 Scalar alq) const;
224
225 void computeWellRatesWithThpAlqProd(const Simulator& ebos_simulator,
226 const WellGroupHelperType& wgHelper,
227 const SummaryState& summary_state,
228 DeferredLogger& deferred_logger,
229 std::vector<Scalar>& potentials,
230 Scalar alq) const;
231
232 std::optional<Scalar>
233 computeBhpAtThpLimitProdWithAlq(const Simulator& ebos_simulator,
234 const WellGroupHelperType& wgHelper,
235 const SummaryState& summary_state,
236 const Scalar alq_value,
237 DeferredLogger& deferred_logger,
238 bool iterate_if_no_solution) const override;
239
240 void updateIPRImplicit(const Simulator& simulator,
241 WellStateType& 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
251
252 std::vector<Scalar>
253 computeCurrentWellRates(const Simulator& ebosSimulator,
254 DeferredLogger& deferred_logger) const override;
255
256 std::vector<Scalar> getPrimaryVars() const override;
257
258 int setPrimaryVars(typename std::vector<Scalar>::const_iterator it) override;
259
260 protected:
262
263 // updating the well_state based on well solution dwells
264 void updateWellState(const Simulator& simulator,
265 const BVectorWell& dwells,
266 WellStateType& well_state,
267 DeferredLogger& deferred_logger);
268
269 using WellConnectionProps = typename StdWellEval::StdWellConnections::Properties;
270
271 // Compute connection level PVT properties needed to calulate the
272 // pressure difference between well connections.
275 const WellStateType& well_state) const;
276
278 const WellStateType& well_state,
279 const WellConnectionProps& props,
280 DeferredLogger& deferred_logger);
281
282 void computeWellConnectionPressures(const Simulator& simulator,
283 const WellStateType& well_state,
284 DeferredLogger& deferred_logger);
285
286 template<class Value>
287 void computePerfRate(const IntensiveQuantities& intQuants,
288 const std::vector<Value>& mob,
289 const Value& bhp,
290 const std::vector<Scalar>& Tw,
291 const int perf,
292 const bool allow_cf,
293 std::vector<Value>& cq_s,
294 PerforationRates<Scalar>& perf_rates,
295 DeferredLogger& deferred_logger) const;
296
297 template<class Value>
298 void computePerfRate(const std::vector<Value>& mob,
299 const Value& pressure,
300 const Value& bhp,
301 const Value& rs,
302 const Value& rv,
303 const Value& rvw,
304 const Value& rsw,
305 std::vector<Value>& b_perfcells_dense,
306 const std::vector<Scalar>& Tw,
307 const int perf,
308 const bool allow_cf,
309 const Value& skin_pressure,
310 const std::vector<Value>& cmix_s,
311 std::vector<Value>& cq_s,
312 PerforationRates<Scalar>& perf_rates,
313 DeferredLogger& deferred_logger) const;
314
315 void computeWellRatesWithBhpIterations(const Simulator& ebosSimulator,
316 const Scalar& bhp,
317 const WellGroupHelperType& wgHelper,
318 std::vector<Scalar>& well_flux,
319 DeferredLogger& deferred_logger) const override;
320
321 std::vector<Scalar>
322 computeWellPotentialWithTHP(const Simulator& ebosSimulator,
323 const WellGroupHelperType& wgHelper,
324 DeferredLogger& deferred_logger,
325 const WellStateType& well_state) const;
326
327 bool computeWellPotentialsImplicit(const Simulator& ebos_simulator,
328 const WellGroupHelperType& wgHelper,
329 std::vector<Scalar>& well_potentials,
330 DeferredLogger& deferred_logger) const;
331
332 // return the density at the perforation[0] of the rank owning this well,
333 // value is cached to minimize the number of broadcasts
334 Scalar getRefDensity() const override;
335
336 // get the mobility for specific perforation
337 template<class Value>
338 void getMobility(const Simulator& simulator,
339 const int perf,
340 std::vector<Value>& mob,
341 DeferredLogger& deferred_logger) const;
342
343 void updateWaterMobilityWithPolymer(const Simulator& simulator,
344 const int perf,
345 std::vector<EvalWell>& mob_water,
346 DeferredLogger& deferred_logger) const;
347
348 void updatePrimaryVariablesNewton(const BVectorWell& dwells,
349 const bool stop_or_zero_rate_target,
350 DeferredLogger& deferred_logger);
351
353 const SummaryState& summary_state,
354 DeferredLogger& deferred_logger) const;
355
356 void assembleWellEqWithoutIteration(const Simulator& simulator,
357 const double dt,
358 const Well::InjectionControls& inj_controls,
359 const Well::ProductionControls& prod_controls,
360 WellStateType& well_state,
361 const GroupState<Scalar>& group_state,
362 DeferredLogger& deferred_logger) override;
363
364 void assembleWellEqWithoutIterationImpl(const Simulator& simulator,
365 const double dt,
366 const Well::InjectionControls& inj_controls,
367 const Well::ProductionControls& prod_controls,
368 WellStateType& well_state,
369 const GroupState<Scalar>& group_state,
370 DeferredLogger& deferred_logger);
371
372 void calculateSinglePerf(const Simulator& simulator,
373 const int perf,
374 WellStateType& well_state,
375 std::vector<RateVector>& connectionRates,
376 std::vector<EvalWell>& cq_s,
377 EvalWell& water_flux_s,
378 EvalWell& cq_s_zfrac_effective,
379 DeferredLogger& deferred_logger) const;
380
381 // check whether the well is operable under BHP limit with current reservoir condition
382 void checkOperabilityUnderBHPLimit(const WellStateType& well_state,
383 const Simulator& simulator,
384 DeferredLogger& deferred_logger) override;
385
386 // check whether the well is operable under THP limit with current reservoir condition
387 void checkOperabilityUnderTHPLimit(const Simulator& simulator,
388 const WellStateType& well_state,
389 const WellGroupHelperType& wgHelper,
390 DeferredLogger& deferred_logger) override;
391
392 // updating the inflow based on the current reservoir condition
393 void updateIPR(const Simulator& simulator,
394 DeferredLogger& deferred_logger) const override;
395
396 // for a well, when all drawdown are in the wrong direction, then this well will not
397 // be able to produce/inject .
398 bool allDrawDownWrongDirection(const Simulator& simulator) const;
399
400 // turn on crossflow to avoid singular well equations
401 // when the well is banned from cross-flow and the BHP is not properly initialized,
402 // we turn on crossflow to avoid singular well equations. It can result in wrong-signed
403 // well rates, it can cause problem for THP calculation
404 // TODO: looking for better alternative to avoid wrong-signed well rates
405 bool openCrossFlowAvoidSingularity(const Simulator& simulator) const;
406
407 // calculate the skin pressure based on water velocity, throughput and polymer concentration.
408 // throughput is used to describe the formation damage during water/polymer injection.
409 // calculated skin pressure will be applied to the drawdown during perforation rate calculation
410 // to handle the effect from formation damage.
411 EvalWell pskin(const Scalar throughput,
412 const EvalWell& water_velocity,
413 const EvalWell& poly_inj_conc,
414 DeferredLogger& deferred_logger) const;
415
416 // calculate the skin pressure based on water velocity, throughput during water injection.
417 EvalWell pskinwater(const Scalar throughput,
418 const EvalWell& water_velocity,
419 DeferredLogger& deferred_logger) const;
420
421 // calculate the injecting polymer molecular weight based on the througput and water velocity
422 EvalWell wpolymermw(const Scalar throughput,
423 const EvalWell& water_velocity,
424 DeferredLogger& deferred_logger) const;
425
426 // modify the water rate for polymer injectivity study
427 void handleInjectivityRate(const Simulator& simulator,
428 const int perf,
429 std::vector<EvalWell>& cq_s) const;
430
431 // handle the extra equations for polymer injectivity study
432 void handleInjectivityEquations(const Simulator& simulator,
433 const WellStateType& well_state,
434 const int perf,
435 const EvalWell& water_flux_s,
436 DeferredLogger& deferred_logger);
437
438 void updateWaterThroughput(const double dt,
439 WellStateType& well_state) const override;
440
441 // checking convergence of extra equations, if there are any
442 void checkConvergenceExtraEqs(const std::vector<Scalar>& res,
443 ConvergenceReport& report) const;
444
445 // updating the connectionRates_ related polymer molecular weight
446 void updateConnectionRatePolyMW(const EvalWell& cq_s_poly,
447 const IntensiveQuantities& int_quants,
448 const WellStateType& well_state,
449 const int perf,
450 std::vector<RateVector>& connectionRates,
451 DeferredLogger& deferred_logger) const;
452
453 std::optional<Scalar>
454 computeBhpAtThpLimitProd(const WellStateType& well_state,
455 const Simulator& simulator,
456 const WellGroupHelperType& wgHelper,
457 const SummaryState& summary_state,
458 DeferredLogger& deferred_logger) const;
459
460 std::optional<Scalar>
461 computeBhpAtThpLimitInj(const Simulator& simulator,
462 const WellGroupHelperType& wgHelper,
463 const SummaryState& summary_state,
464 DeferredLogger& deferred_logger) const;
465
466 private:
467 Eval connectionRateEnergy(const std::vector<EvalWell>& cq_s,
468 const IntensiveQuantities& intQuants,
469 DeferredLogger& deferred_logger) const;
470
471 // density of the first perforation, might not be from this rank
472 Scalar cachedRefDensity{0};
473 };
474
475}
476
477#include "StandardWell_impl.hpp"
478
479#endif // OPM_STANDARDWELL_HEADER_INCLUDED
Contains the classes required to extend the black-oil model by bioeffects.
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 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 to include the effects of ...
Definition: blackoilfoammodules.hh:58
Contains the high level supplements required to extend the black oil model by polymer.
Definition: blackoilpolymermodules.hh:64
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
Definition: GroupState.hpp:41
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:198
Definition: StandardWellEval.hpp:46
typename StandardWellEquations< Scalar, IndexTraits, Indices::numEq >::BVectorWell BVectorWell
Definition: StandardWellEval.hpp:65
DenseAd::Evaluation< Scalar, Indices::numEq > Eval
Definition: StandardWellEval.hpp:64
Definition: StandardWell.hpp:60
void calculateExplicitQuantities(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1406
EvalWell wpolymermw(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2040
typename StdWellEval::EvalWell EvalWell
Definition: StandardWell.hpp:121
void updateWellStateFromPrimaryVariables(WellStateType &well_state, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:775
static constexpr int numStaticWellEq
Definition: StandardWell.hpp:102
bool regularize_
Definition: StandardWell.hpp:261
WellConnectionProps computePropertiesForWellConnectionPressures(const Simulator &simulator, const WellStateType &well_state) const
Definition: StandardWell_impl.hpp:1140
typename StdWellEval::BVectorWell BVectorWell
Definition: StandardWell.hpp:122
std::vector< Scalar > computeWellPotentialWithTHP(const Simulator &ebosSimulator, const WellGroupHelperType &wgHelper, DeferredLogger &deferred_logger, const WellStateType &well_state) const
Definition: StandardWell_impl.hpp:1589
void addWellContributions(SparseMatrixAdapter &mat) const override
Definition: StandardWell_impl.hpp:1937
std::vector< Scalar > getPrimaryVars() const override
Definition: StandardWell_impl.hpp:2602
void addWellPressureEquations(PressureMatrix &mat, const BVector &x, const int pressureVarIndex, const bool use_well_weights, const WellStateType &well_state) const override
Definition: StandardWell_impl.hpp:1945
bool iterateWellEqWithSwitching(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, const bool fixed_control=false, const bool fixed_status=false) override
Definition: StandardWell_impl.hpp:2436
void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) override
computing the well potentials for group control
Definition: StandardWell_impl.hpp:1761
void updateWaterMobilityWithPolymer(const Simulator &simulator, const int perf, std::vector< EvalWell > &mob_water, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1873
void assembleWellEqWithoutIteration(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:341
std::vector< Scalar > computeCurrentWellRates(const Simulator &ebosSimulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:2566
void calculateSinglePerf(const Simulator &simulator, const int perf, WellStateType &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:496
bool iterateWellEqWithControl(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) override
Definition: StandardWell_impl.hpp:2373
void updatePrimaryVariablesNewton(const BVectorWell &dwells, const bool stop_or_zero_rate_target, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:752
bool computeWellPotentialsImplicit(const Simulator &ebos_simulator, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1625
void computeWellRatesWithBhpIterations(const Simulator &ebosSimulator, const Scalar &bhp, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1518
void computeWellConnectionPressures(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1366
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_conservation_quantities, const int num_phases, const int index_of_well, const std::vector< PerforationData< Scalar > > &perf_data)
Definition: StandardWell_impl.hpp:52
typename StdWellEval::StdWellConnections::Properties WellConnectionProps
Definition: StandardWell.hpp:269
void updateConnectionRatePolyMW(const EvalWell &cq_s_poly, const IntensiveQuantities &int_quants, const WellStateType &well_state, const int perf, std::vector< RateVector > &connectionRates, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2189
void computeWellRatesWithBhp(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1472
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:674
void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:799
void updateIPRImplicit(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:892
void computeWellConnectionDensitesPressures(const Simulator &simulator, const WellStateType &well_state, const WellConnectionProps &props, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1310
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:93
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 override
Definition: StandardWell_impl.hpp:2248
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:730
void handleInjectivityEquations(const Simulator &simulator, const WellStateType &well_state, const int perf, const EvalWell &water_flux_s, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:2123
virtual void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: StandardWell_impl.hpp:1420
virtual ConvergenceReport getWellConvergence(const Simulator &simulator, const WellStateType &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:1191
void checkConvergenceExtraEqs(const std::vector< Scalar > &res, ConvergenceReport &report) const
Definition: StandardWell_impl.hpp:2170
typename StdWellEval::Eval Eval
Definition: StandardWell.hpp:120
static constexpr int Bhp
Definition: StandardWell.hpp:107
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1129
void computeWellRatesWithThpAlqProd(const Simulator &ebos_simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1742
typename FluidSystem::IndexTraitsType IndexTraits
Definition: StandardWell.hpp:124
Scalar computeWellRatesAndBhpWithThpAlqProd(const Simulator &ebos_simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1712
bool allDrawDownWrongDirection(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1087
EvalWell pskin(const Scalar throughput, const EvalWell &water_velocity, const EvalWell &poly_inj_conc, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1996
static constexpr int numWellConservationEq
Definition: StandardWell.hpp:97
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: StandardWell_impl.hpp:2619
void updateWaterThroughput(const double dt, WellStateType &well_state) const override
Definition: StandardWell_impl.hpp:2069
void checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellStateType &well_state, const WellGroupHelperType &wgHelper, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1037
void checkOperabilityUnderBHPLimit(const WellStateType &well_state, const Simulator &simulator, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:968
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2335
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1452
void assembleWellEqWithoutIterationImpl(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:367
EvalWell pskinwater(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1965
void solveEqAndUpdateWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1384
void handleInjectivityRate(const Simulator &simulator, const int perf, std::vector< EvalWell > &cq_s) const
Definition: StandardWell_impl.hpp:2100
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) override
Definition: StandardWell_impl.hpp:76
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1235
void updatePrimaryVariables(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1839
std::optional< Scalar > computeBhpAtThpLimitProd(const WellStateType &well_state, const Simulator &simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2231
static constexpr int numWellControlEq
Definition: StandardWell.hpp:99
Scalar getRefDensity() const override
Definition: StandardWell_impl.hpp:1862
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: StandardWell_impl.hpp:1824
Definition: WellGroupHelper.hpp:51
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 VFPProperties< FluidSystem::Scalar, FluidSystem::IndexTraitsType > * vfp_properties_
Definition: WellInterfaceGeneric.hpp:399
Definition: WellInterface.hpp:77
static constexpr bool has_brine
Definition: WellInterface.hpp:122
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: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
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:98
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:87
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: WellInterface.hpp:88
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:97
static constexpr bool has_bioeffects
Definition: WellInterface.hpp:126
static constexpr bool has_polymermw
Definition: WellInterface.hpp:120
static constexpr bool has_polymer
Definition: WellInterface.hpp:115
typename Base::ModelParameters ModelParameters
Definition: WellInterface.hpp:111
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:84
static constexpr bool has_solvent
Definition: WellInterface.hpp:113
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
static constexpr bool has_micp
Definition: WellInterface.hpp:127
static constexpr bool has_energy
Definition: WellInterface.hpp:116
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:86
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:89
static constexpr bool has_zFraction
Definition: WellInterface.hpp:114
Definition: WellProdIndexCalculator.hpp:37
Definition: WellState.hpp:66
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
Definition: PerforationData.hpp:41