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
85 using Base::has_foam;
86 using Base::has_brine;
87 using Base::has_energy;
89 using Base::has_micp;
90
93 using typename Base::PressureMatrix;
94
95 // number of the conservation equations
96 static constexpr int numWellConservationEq = Indices::numPhases + Indices::numSolvents;
97 // number of the well control equations
98 static constexpr int numWellControlEq = 1;
99 // number of the well equations that will always be used
100 // based on the solution strategy, there might be other well equations be introduced
102
103 // the index for Bhp in primary variables and also the index of well control equation
104 // they both will be the last one in their respective system.
105 // TODO: we should have indices for the well equations and well primary variables separately
106 static constexpr int Bhp = numStaticWellEq - numWellControlEq;
107
109
110 using typename Base::Scalar;
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 using IndexTraits = typename FluidSystem::IndexTraitsType;
125
126 StandardWell(const Well& well,
127 const ParallelWellInfo<Scalar>& pw_info,
128 const int time_step,
129 const ModelParameters& param,
130 const RateConverterType& rate_converter,
131 const int pvtRegionIdx,
132 const int num_conservation_quantities,
133 const int num_phases,
134 const int index_of_well,
135 const std::vector<PerforationData<Scalar>>& perf_data);
136
137 virtual void init(const std::vector<Scalar>& depth_arg,
138 const Scalar gravity_arg,
139 const std::vector<Scalar>& B_avg,
140 const bool changed_to_open_this_step) override;
141
143 virtual ConvergenceReport getWellConvergence(const Simulator& simulator,
144 const WellStateType& 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 WellStateType& well_state,
159 DeferredLogger& deferred_logger) override;
160
162 void computeWellPotentials(const Simulator& simulator,
163 const WellStateType& well_state,
164 std::vector<Scalar>& well_potentials,
165 DeferredLogger& deferred_logger) /* const */ override;
166
167 void updatePrimaryVariables(const Simulator& simulator,
168 const WellStateType& well_state,
169 DeferredLogger& deferred_logger) override;
170
171 void solveEqAndUpdateWellState(const Simulator& simulator,
172 WellStateType& well_state,
173 DeferredLogger& deferred_logger) override;
174
175 void calculateExplicitQuantities(const Simulator& simulator,
176 const WellStateType& well_state,
177 DeferredLogger& deferred_logger) override; // should be const?
178
179 void updateProductivityIndex(const Simulator& simulator,
180 const WellProdIndexCalculator<Scalar>& wellPICalc,
181 WellStateType& 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 WellStateType& 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 WellStateType& 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 WellStateType& 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
215 /* returns BHP */
217 const SummaryState &summary_state,
218 DeferredLogger& deferred_logger,
219 std::vector<Scalar>& potentials,
220 Scalar alq) const;
221
222 void computeWellRatesWithThpAlqProd(const Simulator& ebos_simulator,
223 const SummaryState& summary_state,
224 DeferredLogger& deferred_logger,
225 std::vector<Scalar>& potentials,
226 Scalar alq) const;
227
228 std::optional<Scalar>
229 computeBhpAtThpLimitProdWithAlq(const Simulator& ebos_simulator,
230 const SummaryState& summary_state,
231 const Scalar alq_value,
232 DeferredLogger& deferred_logger,
233 bool iterate_if_no_solution) const override;
234
235 void updateIPRImplicit(const Simulator& simulator,
236 WellStateType& well_state,
237 DeferredLogger& deferred_logger) override;
238
239 void computeWellRatesWithBhp(const Simulator& ebosSimulator,
240 const Scalar& bhp,
241 std::vector<Scalar>& well_flux,
242 DeferredLogger& deferred_logger) const override;
243
244 // NOTE: These cannot be protected since they are used by GasLiftRuntime
246
247 std::vector<Scalar>
248 computeCurrentWellRates(const Simulator& ebosSimulator,
249 DeferredLogger& deferred_logger) const override;
250
251 std::vector<Scalar> getPrimaryVars() const override;
252
253 int setPrimaryVars(typename std::vector<Scalar>::const_iterator it) override;
254
255 protected:
257
258 // updating the well_state based on well solution dwells
259 void updateWellState(const Simulator& simulator,
260 const BVectorWell& dwells,
261 WellStateType& well_state,
262 DeferredLogger& deferred_logger);
263
264 using WellConnectionProps = typename StdWellEval::StdWellConnections::Properties;
265
266 // Compute connection level PVT properties needed to calulate the
267 // pressure difference between well connections.
270 const WellStateType& well_state) const;
271
273 const WellStateType& well_state,
274 const WellConnectionProps& props,
275 DeferredLogger& deferred_logger);
276
277 void computeWellConnectionPressures(const Simulator& simulator,
278 const WellStateType& well_state,
279 DeferredLogger& deferred_logger);
280
281 template<class Value>
282 void computePerfRate(const IntensiveQuantities& intQuants,
283 const std::vector<Value>& mob,
284 const Value& bhp,
285 const std::vector<Scalar>& Tw,
286 const int perf,
287 const bool allow_cf,
288 std::vector<Value>& cq_s,
289 PerforationRates<Scalar>& perf_rates,
290 DeferredLogger& deferred_logger) const;
291
292 template<class Value>
293 void computePerfRate(const std::vector<Value>& mob,
294 const Value& pressure,
295 const Value& bhp,
296 const Value& rs,
297 const Value& rv,
298 const Value& rvw,
299 const Value& rsw,
300 std::vector<Value>& b_perfcells_dense,
301 const std::vector<Scalar>& Tw,
302 const int perf,
303 const bool allow_cf,
304 const Value& skin_pressure,
305 const std::vector<Value>& cmix_s,
306 std::vector<Value>& cq_s,
307 PerforationRates<Scalar>& perf_rates,
308 DeferredLogger& deferred_logger) const;
309
310 void computeWellRatesWithBhpIterations(const Simulator& ebosSimulator,
311 const Scalar& bhp,
312 std::vector<Scalar>& well_flux,
313 DeferredLogger& deferred_logger) const override;
314
315 std::vector<Scalar>
316 computeWellPotentialWithTHP(const Simulator& ebosSimulator,
317 DeferredLogger& deferred_logger,
318 const WellStateType& well_state) const;
319
320 bool computeWellPotentialsImplicit(const Simulator& ebos_simulator,
321 const WellStateType& well_state,
322 std::vector<Scalar>& well_potentials,
323 DeferredLogger& deferred_logger) const;
324
325 Scalar getRefDensity() const override;
326
327 // get the mobility for specific perforation
328 template<class Value>
329 void getMobility(const Simulator& simulator,
330 const int perf,
331 std::vector<Value>& mob,
332 DeferredLogger& deferred_logger) const;
333
334 void updateWaterMobilityWithPolymer(const Simulator& simulator,
335 const int perf,
336 std::vector<EvalWell>& mob_water,
337 DeferredLogger& deferred_logger) const;
338
339 void updatePrimaryVariablesNewton(const BVectorWell& dwells,
340 const bool stop_or_zero_rate_target,
341 DeferredLogger& deferred_logger);
342
344 const SummaryState& summary_state,
345 DeferredLogger& deferred_logger) const;
346
347 void assembleWellEqWithoutIteration(const Simulator& simulator,
348 const double dt,
349 const Well::InjectionControls& inj_controls,
350 const Well::ProductionControls& prod_controls,
351 WellStateType& well_state,
352 const GroupState<Scalar>& group_state,
353 DeferredLogger& deferred_logger) override;
354
355 void assembleWellEqWithoutIterationImpl(const Simulator& simulator,
356 const double dt,
357 const Well::InjectionControls& inj_controls,
358 const Well::ProductionControls& prod_controls,
359 WellStateType& well_state,
360 const GroupState<Scalar>& group_state,
361 DeferredLogger& deferred_logger);
362
363 void calculateSinglePerf(const Simulator& simulator,
364 const int perf,
365 WellStateType& well_state,
366 std::vector<RateVector>& connectionRates,
367 std::vector<EvalWell>& cq_s,
368 EvalWell& water_flux_s,
369 EvalWell& cq_s_zfrac_effective,
370 DeferredLogger& deferred_logger) const;
371
372 // check whether the well is operable under BHP limit with current reservoir condition
373 void checkOperabilityUnderBHPLimit(const WellStateType& well_state,
374 const Simulator& simulator,
375 DeferredLogger& deferred_logger) override;
376
377 // check whether the well is operable under THP limit with current reservoir condition
378 void checkOperabilityUnderTHPLimit(const Simulator& simulator,
379 const WellStateType& well_state,
380 DeferredLogger& deferred_logger) override;
381
382 // updating the inflow based on the current reservoir condition
383 void updateIPR(const Simulator& simulator,
384 DeferredLogger& deferred_logger) const override;
385
386 // for a well, when all drawdown are in the wrong direction, then this well will not
387 // be able to produce/inject .
388 bool allDrawDownWrongDirection(const Simulator& simulator) const;
389
390 // turn on crossflow to avoid singular well equations
391 // when the well is banned from cross-flow and the BHP is not properly initialized,
392 // we turn on crossflow to avoid singular well equations. It can result in wrong-signed
393 // well rates, it can cause problem for THP calculation
394 // TODO: looking for better alternative to avoid wrong-signed well rates
395 bool openCrossFlowAvoidSingularity(const Simulator& simulator) const;
396
397 // calculate the skin pressure based on water velocity, throughput and polymer concentration.
398 // throughput is used to describe the formation damage during water/polymer injection.
399 // calculated skin pressure will be applied to the drawdown during perforation rate calculation
400 // to handle the effect from formation damage.
401 EvalWell pskin(const Scalar throughput,
402 const EvalWell& water_velocity,
403 const EvalWell& poly_inj_conc,
404 DeferredLogger& deferred_logger) const;
405
406 // calculate the skin pressure based on water velocity, throughput during water injection.
407 EvalWell pskinwater(const Scalar throughput,
408 const EvalWell& water_velocity,
409 DeferredLogger& deferred_logger) const;
410
411 // calculate the injecting polymer molecular weight based on the througput and water velocity
412 EvalWell wpolymermw(const Scalar throughput,
413 const EvalWell& water_velocity,
414 DeferredLogger& deferred_logger) const;
415
416 // modify the water rate for polymer injectivity study
417 void handleInjectivityRate(const Simulator& simulator,
418 const int perf,
419 std::vector<EvalWell>& cq_s) const;
420
421 // handle the extra equations for polymer injectivity study
422 void handleInjectivityEquations(const Simulator& simulator,
423 const WellStateType& well_state,
424 const int perf,
425 const EvalWell& water_flux_s,
426 DeferredLogger& deferred_logger);
427
428 void updateWaterThroughput(const double dt,
429 WellStateType& well_state) const override;
430
431 // checking convergence of extra equations, if there are any
432 void checkConvergenceExtraEqs(const std::vector<Scalar>& res,
433 ConvergenceReport& report) const;
434
435 // updating the connectionRates_ related polymer molecular weight
436 void updateConnectionRatePolyMW(const EvalWell& cq_s_poly,
437 const IntensiveQuantities& int_quants,
438 const WellStateType& well_state,
439 const int perf,
440 std::vector<RateVector>& connectionRates,
441 DeferredLogger& deferred_logger) const;
442
443 std::optional<Scalar>
444 computeBhpAtThpLimitProd(const WellStateType& well_state,
445 const Simulator& simulator,
446 const SummaryState& summary_state,
447 DeferredLogger& deferred_logger) const;
448
449 std::optional<Scalar>
450 computeBhpAtThpLimitInj(const Simulator& simulator,
451 const SummaryState& summary_state,
452 DeferredLogger& deferred_logger) const;
453
454 private:
455 Eval connectionRateEnergy(const std::vector<EvalWell>& cq_s,
456 const IntensiveQuantities& intQuants,
457 DeferredLogger& deferred_logger) const;
458 };
459
460}
461
462#include "StandardWell_impl.hpp"
463
464#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:1391
EvalWell wpolymermw(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2005
typename StdWellEval::EvalWell EvalWell
Definition: StandardWell.hpp:120
void updateWellStateFromPrimaryVariables(WellStateType &well_state, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:771
static constexpr int numStaticWellEq
Definition: StandardWell.hpp:101
bool regularize_
Definition: StandardWell.hpp:256
WellConnectionProps computePropertiesForWellConnectionPressures(const Simulator &simulator, const WellStateType &well_state) const
Definition: StandardWell_impl.hpp:1130
typename StdWellEval::BVectorWell BVectorWell
Definition: StandardWell.hpp:121
std::optional< Scalar > computeBhpAtThpLimitProd(const WellStateType &well_state, const Simulator &simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2197
void computeWellRatesWithThpAlqProd(const Simulator &ebos_simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1710
void addWellContributions(SparseMatrixAdapter &mat) const override
Definition: StandardWell_impl.hpp:1900
std::vector< Scalar > getPrimaryVars() const override
Definition: StandardWell_impl.hpp:2529
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:1908
void updateWaterMobilityWithPolymer(const Simulator &simulator, const int perf, std::vector< EvalWell > &mob_water, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1836
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:2493
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:493
void updatePrimaryVariablesNewton(const BVectorWell &dwells, const bool stop_or_zero_rate_target, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:748
void computeWellConnectionPressures(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1351
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:264
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:2155
bool iterateWellEqWithSwitching(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, const bool fixed_control=false, const bool fixed_status=false) override
Definition: StandardWell_impl.hpp:2373
void computeWellRatesWithBhp(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1457
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:670
void computeWellRatesWithBhpIterations(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1503
void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:790
std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &ebos_simulator, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger, bool iterate_if_no_solution) const override
Definition: StandardWell_impl.hpp:2212
void updateIPRImplicit(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:883
void computeWellConnectionDensitesPressures(const Simulator &simulator, const WellStateType &well_state, const WellConnectionProps &props, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1300
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
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:726
void handleInjectivityEquations(const Simulator &simulator, const WellStateType &well_state, const int perf, const EvalWell &water_flux_s, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:2088
virtual void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: StandardWell_impl.hpp:1405
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:1181
void checkConvergenceExtraEqs(const std::vector< Scalar > &res, ConvergenceReport &report) const
Definition: StandardWell_impl.hpp:2136
typename StdWellEval::Eval Eval
Definition: StandardWell.hpp:119
Scalar computeWellRatesAndBhpWithThpAlqProd(const Simulator &ebos_simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1681
static constexpr int Bhp
Definition: StandardWell.hpp:106
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1119
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2294
typename FluidSystem::IndexTraitsType IndexTraits
Definition: StandardWell.hpp:123
bool allDrawDownWrongDirection(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1077
EvalWell pskin(const Scalar throughput, const EvalWell &water_velocity, const EvalWell &poly_inj_conc, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1960
static constexpr int numWellConservationEq
Definition: StandardWell.hpp:96
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: StandardWell_impl.hpp:2546
void checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1028
void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) override
computing the well potentials for group control
Definition: StandardWell_impl.hpp:1727
bool computeWellPotentialsImplicit(const Simulator &ebos_simulator, const WellStateType &well_state, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1602
void updateWaterThroughput(const double dt, WellStateType &well_state) const override
Definition: StandardWell_impl.hpp:2034
void checkOperabilityUnderBHPLimit(const WellStateType &well_state, const Simulator &simulator, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:959
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1437
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
bool iterateWellEqWithControl(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:2326
EvalWell pskinwater(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1928
void solveEqAndUpdateWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1369
void handleInjectivityRate(const Simulator &simulator, const int perf, std::vector< EvalWell > &cq_s) const
Definition: StandardWell_impl.hpp:2065
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
std::vector< Scalar > computeWellPotentialWithTHP(const Simulator &ebosSimulator, DeferredLogger &deferred_logger, const WellStateType &well_state) const
Definition: StandardWell_impl.hpp:1567
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1225
void updatePrimaryVariables(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1802
static constexpr int numWellControlEq
Definition: StandardWell.hpp:98
Scalar getRefDensity() const override
Definition: StandardWell_impl.hpp:1825
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: StandardWell_impl.hpp:1787
static constexpr int Oil
Definition: WellInterfaceFluidSystem.hpp:62
static constexpr int Water
Definition: WellInterfaceFluidSystem.hpp:61
static constexpr int Gas
Definition: WellInterfaceFluidSystem.hpp:63
int pvtRegionIdx() const
Definition: WellInterfaceGeneric.hpp:128
const VFPProperties< FluidSystem::Scalar, FluidSystem::IndexTraitsType > * vfp_properties_
Definition: WellInterfaceGeneric.hpp:397
Definition: WellInterface.hpp:76
static constexpr bool has_brine
Definition: WellInterface.hpp:119
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:81
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:103
const std::vector< RateVector > & connectionRates() const
Definition: WellInterface.hpp:342
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:136
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:97
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:86
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:82
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: WellInterface.hpp:87
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:96
static constexpr bool has_bioeffects
Definition: WellInterface.hpp:123
static constexpr bool has_polymermw
Definition: WellInterface.hpp:117
static constexpr bool has_polymer
Definition: WellInterface.hpp:113
typename Base::ModelParameters ModelParameters
Definition: WellInterface.hpp:109
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:83
static constexpr bool has_solvent
Definition: WellInterface.hpp:111
static constexpr bool has_foam
Definition: WellInterface.hpp:118
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: WellInterface.hpp:89
static constexpr bool has_micp
Definition: WellInterface.hpp:124
static constexpr bool has_energy
Definition: WellInterface.hpp:114
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:85
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:88
static constexpr bool has_zFraction
Definition: WellInterface.hpp:112
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