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::GroupStateHelperType;
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 GroupStateHelperType& groupStateHelper,
145 const std::vector<Scalar>& B_avg,
146 const bool relax_tolerance) const override;
147
149 virtual void apply(const BVector& x, BVector& Ax) const override;
151 virtual void apply(BVector& r) const override;
152
156 const BVector& x,
157 const GroupStateHelperType& groupStateHelper,
158 WellStateType& well_state) override;
159
161 void computeWellPotentials(const Simulator& simulator,
162 const WellStateType& well_state,
163 const GroupStateHelperType& groupStateHelper,
164 std::vector<Scalar>& well_potentials) /* const */ override;
165
166 void updatePrimaryVariables(const GroupStateHelperType& groupStateHelper) override;
167
168 void solveEqAndUpdateWellState(const Simulator& simulator,
169 const GroupStateHelperType& groupStateHelper,
170 WellStateType& well_state) override;
171
172 void calculateExplicitQuantities(const Simulator& simulator,
173 const GroupStateHelperType& groupStateHelper) override; // should be const?
174
175 void updateProductivityIndex(const Simulator& simulator,
176 const WellProdIndexCalculator<Scalar>& wellPICalc,
177 WellStateType& well_state,
178 DeferredLogger& deferred_logger) const override;
179
180 Scalar connectionDensity(const int globalConnIdx,
181 const int openConnIdx) const override;
182
183 void addWellContributions(SparseMatrixAdapter& mat) const override;
184
186 const BVector& x,
187 const int pressureVarIndex,
188 const bool use_well_weights,
189 const WellStateType& well_state) const override;
190
191 // iterate well equations with the specified control until converged
192 bool iterateWellEqWithControl(const Simulator& simulator,
193 const double dt,
194 const Well::InjectionControls& inj_controls,
195 const Well::ProductionControls& prod_controls,
196 const GroupStateHelperType& groupStateHelper,
197 WellStateType& well_state) override;
198
199 // iterate well equations including control switching
200 bool iterateWellEqWithSwitching(const Simulator& simulator,
201 const double dt,
202 const Well::InjectionControls& inj_controls,
203 const Well::ProductionControls& prod_controls,
204 const GroupStateHelperType& groupStateHelper,
205 WellStateType& well_state,
206 const bool fixed_control,
207 const bool fixed_status,
208 const bool solving_with_zero_rate) override;
209
210 /* returns BHP */
212 const GroupStateHelperType& groupStateHelper,
213 const SummaryState &summary_state,
214 std::vector<Scalar>& potentials,
215 Scalar alq) const;
216
217 void computeWellRatesWithThpAlqProd(const Simulator& ebos_simulator,
218 const GroupStateHelperType& groupStateHelper,
219 const SummaryState& summary_state,
220 std::vector<Scalar>& potentials,
221 Scalar alq) const;
222
223 std::optional<Scalar>
224 computeBhpAtThpLimitProdWithAlq(const Simulator& ebos_simulator,
225 const GroupStateHelperType& groupStateHelper,
226 const SummaryState& summary_state,
227 const Scalar alq_value,
228 bool iterate_if_no_solution) const override;
229
230 void updateIPRImplicit(const Simulator& simulator,
231 const GroupStateHelperType& groupStateHelper,
232 WellStateType& well_state) override;
233
234 void computeWellRatesWithBhp(const Simulator& ebosSimulator,
235 const Scalar& bhp,
236 std::vector<Scalar>& well_flux,
237 DeferredLogger& deferred_logger) const override;
238
239 Scalar maxPerfPress(const Simulator& simulator) const override;
240
241 // NOTE: These cannot be protected since they are used by GasLiftRuntime
243
244 std::vector<Scalar>
245 computeCurrentWellRates(const Simulator& ebosSimulator,
246 DeferredLogger& deferred_logger) const override;
247
248 std::vector<Scalar> getPrimaryVars() const override;
249
250 int setPrimaryVars(typename std::vector<Scalar>::const_iterator it) override;
251
252 protected:
254
255 // updating the well_state based on well solution dwells
256 void updateWellState(const Simulator& simulator,
257 const BVectorWell& dwells,
258 const GroupStateHelperType& groupStateHelper,
259 WellStateType& well_state);
260
261 using WellConnectionProps = typename StdWellEval::StdWellConnections::Properties;
262
263 // Compute connection level PVT properties needed to calulate the
264 // pressure difference between well connections.
267 const WellStateType& well_state) const;
268
270 const GroupStateHelperType& groupStateHelper,
271 const WellConnectionProps& props);
272
273 void computeWellConnectionPressures(const Simulator& simulator,
274 const GroupStateHelperType& groupStateHelper);
275
276 template<class Value>
277 void computePerfRate(const IntensiveQuantities& intQuants,
278 const std::vector<Value>& mob,
279 const Value& bhp,
280 const std::vector<Value>& Tw,
281 const int perf,
282 const bool allow_cf,
283 std::vector<Value>& cq_s,
284 PerforationRates<Scalar>& perf_rates,
285 DeferredLogger& deferred_logger) const;
286
287 template<class Value>
288 void computePerfRate(const std::vector<Value>& mob,
289 const Value& pressure,
290 const Value& bhp,
291 const Value& rs,
292 const Value& rv,
293 const Value& rvw,
294 const Value& rsw,
295 std::vector<Value>& b_perfcells_dense,
296 const std::vector<Value>& Tw,
297 const int perf,
298 const bool allow_cf,
299 const Value& skin_pressure,
300 const std::vector<Value>& cmix_s,
301 std::vector<Value>& cq_s,
302 PerforationRates<Scalar>& perf_rates,
303 DeferredLogger& deferred_logger) const;
304
305 void computeWellRatesWithBhpIterations(const Simulator& ebosSimulator,
306 const Scalar& bhp,
307 const GroupStateHelperType& groupStateHelper,
308 std::vector<Scalar>& well_flux) const override;
309
310 std::vector<Scalar>
311 computeWellPotentialWithTHP(const Simulator& ebosSimulator,
312 const GroupStateHelperType& groupStateHelper,
313 const WellStateType& well_state) const;
314
315 bool computeWellPotentialsImplicit(const Simulator& ebos_simulator,
316 const GroupStateHelperType& groupStateHelper,
317 std::vector<Scalar>& well_potentials) const;
318
319 // return the density at the perforation[0] of the rank owning this well,
320 // value is cached to minimize the number of broadcasts
321 Scalar getRefDensity() const override;
322
323 // get the transmissibility multiplier for specific perforation
324 template<class Value>
325 void getTransMult(Value& trans_mult,
326 const Simulator& simulator,
327 const int cell_indx) const;
328
329 // get the mobility for specific perforation
330 template<class Value>
331 void getMobility(const Simulator& simulator,
332 const int perf,
333 std::vector<Value>& mob,
334 DeferredLogger& deferred_logger) const;
335
336 void updateWaterMobilityWithPolymer(const Simulator& simulator,
337 const int perf,
338 std::vector<EvalWell>& mob_water,
339 DeferredLogger& deferred_logger) const;
340
341 void updatePrimaryVariablesNewton(const BVectorWell& dwells,
342 const bool stop_or_zero_rate_target,
343 DeferredLogger& deferred_logger);
344
346 const SummaryState& summary_state,
347 DeferredLogger& deferred_logger) const;
348
349 void assembleWellEqWithoutIteration(const Simulator& simulator,
350 const GroupStateHelperType& groupStateHelper,
351 const double dt,
352 const Well::InjectionControls& inj_controls,
353 const Well::ProductionControls& prod_controls,
354 WellStateType& well_state,
355 const bool solving_with_zero_rate) override;
356
357 void assembleWellEqWithoutIterationImpl(const Simulator& simulator,
358 const GroupStateHelperType& groupStateHelper,
359 const double dt,
360 const Well::InjectionControls& inj_controls,
361 const Well::ProductionControls& prod_controls,
362 WellStateType& well_state,
363 const bool solving_with_zero_rate);
364
365 void calculateSinglePerf(const Simulator& simulator,
366 const int perf,
367 WellStateType& well_state,
368 std::vector<RateVector>& connectionRates,
369 std::vector<EvalWell>& cq_s,
370 EvalWell& water_flux_s,
371 EvalWell& cq_s_zfrac_effective,
372 DeferredLogger& deferred_logger) const;
373
374 // check whether the well is operable under BHP limit with current reservoir condition
375 void checkOperabilityUnderBHPLimit(const WellStateType& well_state,
376 const Simulator& simulator,
377 DeferredLogger& deferred_logger) override;
378
379 // check whether the well is operable under THP limit with current reservoir condition
380 void checkOperabilityUnderTHPLimit(const Simulator& simulator,
381 const WellStateType& well_state,
382 const GroupStateHelperType& groupStateHelper) override;
383
384 // updating the inflow based on the current reservoir condition
385 void updateIPR(const Simulator& simulator,
386 DeferredLogger& deferred_logger) const override;
387
388 // for a well, when all drawdown are in the wrong direction, then this well will not
389 // be able to produce/inject .
390 bool allDrawDownWrongDirection(const Simulator& simulator) const;
391
392 // turn on crossflow to avoid singular well equations
393 // when the well is banned from cross-flow and the BHP is not properly initialized,
394 // we turn on crossflow to avoid singular well equations. It can result in wrong-signed
395 // well rates, it can cause problem for THP calculation
396 // TODO: looking for better alternative to avoid wrong-signed well rates
397 bool openCrossFlowAvoidSingularity(const Simulator& simulator) const;
398
399 // calculate the skin pressure based on water velocity, throughput and polymer concentration.
400 // throughput is used to describe the formation damage during water/polymer injection.
401 // calculated skin pressure will be applied to the drawdown during perforation rate calculation
402 // to handle the effect from formation damage.
403 EvalWell pskin(const Scalar throughput,
404 const EvalWell& water_velocity,
405 const EvalWell& poly_inj_conc,
406 DeferredLogger& deferred_logger) const;
407
408 // calculate the skin pressure based on water velocity, throughput during water injection.
409 EvalWell pskinwater(const Scalar throughput,
410 const EvalWell& water_velocity,
411 DeferredLogger& deferred_logger) const;
412
413 // calculate the injecting polymer molecular weight based on the througput and water velocity
414 EvalWell wpolymermw(const Scalar throughput,
415 const EvalWell& water_velocity,
416 DeferredLogger& deferred_logger) const;
417
418 // modify the water rate for polymer injectivity study
419 void handleInjectivityRate(const Simulator& simulator,
420 const int perf,
421 std::vector<EvalWell>& cq_s) const;
422
423 // handle the extra equations for polymer injectivity study
424 void handleInjectivityEquations(const Simulator& simulator,
425 const WellStateType& well_state,
426 const int perf,
427 const EvalWell& water_flux_s,
428 DeferredLogger& deferred_logger);
429
430 void updateWaterThroughput(const double dt,
431 WellStateType& well_state) const override;
432
433 // checking convergence of extra equations, if there are any
434 void checkConvergenceExtraEqs(const std::vector<Scalar>& res,
435 ConvergenceReport& report) const;
436
437 // updating the connectionRates_ related polymer molecular weight
438 void updateConnectionRatePolyMW(const EvalWell& cq_s_poly,
439 const IntensiveQuantities& int_quants,
440 const WellStateType& well_state,
441 const int perf,
442 std::vector<RateVector>& connectionRates,
443 DeferredLogger& deferred_logger) const;
444
445 std::optional<Scalar>
446 computeBhpAtThpLimitProd(const WellStateType& well_state,
447 const Simulator& simulator,
448 const GroupStateHelperType& groupStateHelper,
449 const SummaryState& summary_state) const;
450
451 std::optional<Scalar>
452 computeBhpAtThpLimitInj(const Simulator& simulator,
453 const GroupStateHelperType& groupStateHelper,
454 const SummaryState& summary_state) const;
455
456 private:
457 Eval connectionRateEnergy(const std::vector<EvalWell>& cq_s,
458 const IntensiveQuantities& intQuants,
459 DeferredLogger& deferred_logger) const;
460
461 // density of the first perforation, might not be from this rank
462 Scalar cachedRefDensity{0};
463 };
464
465}
466
467#include "StandardWell_impl.hpp"
468
469#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: GroupStateHelper.hpp:54
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:198
Definition: StandardWellEval.hpp:46
DenseAd::Evaluation< Scalar, Indices::numDerivatives > Eval
Definition: StandardWellEval.hpp:64
typename StandardWellEquations< Scalar, IndexTraits, Indices::numEq >::BVectorWell BVectorWell
Definition: StandardWellEval.hpp:65
Definition: StandardWell.hpp:60
EvalWell wpolymermw(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2071
std::vector< Scalar > computeWellPotentialWithTHP(const Simulator &ebosSimulator, const GroupStateHelperType &groupStateHelper, const WellStateType &well_state) const
Definition: StandardWell_impl.hpp:1618
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:804
static constexpr int numStaticWellEq
Definition: StandardWell.hpp:102
bool regularize_
Definition: StandardWell.hpp:253
virtual ConvergenceReport getWellConvergence(const GroupStateHelperType &groupStateHelper, const std::vector< Scalar > &B_avg, const bool relax_tolerance) const override
check whether the well equations get converged for this well
Definition: StandardWell_impl.hpp:1220
WellConnectionProps computePropertiesForWellConnectionPressures(const Simulator &simulator, const WellStateType &well_state) const
Definition: StandardWell_impl.hpp:1169
std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, const Scalar alq_value, bool iterate_if_no_solution) const override
Definition: StandardWell_impl.hpp:2277
typename StdWellEval::BVectorWell BVectorWell
Definition: StandardWell.hpp:122
void addWellContributions(SparseMatrixAdapter &mat) const override
Definition: StandardWell_impl.hpp:1968
std::vector< Scalar > getPrimaryVars() const override
Definition: StandardWell_impl.hpp:2618
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: StandardWell_impl.hpp:757
void updatePrimaryVariables(const GroupStateHelperType &groupStateHelper) override
Definition: StandardWell_impl.hpp:1868
void solveEqAndUpdateWellState(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) override
Definition: StandardWell_impl.hpp:1412
void computeWellConnectionDensitesPressures(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const WellConnectionProps &props)
Definition: StandardWell_impl.hpp:1338
std::optional< Scalar > computeBhpAtThpLimitProd(const WellStateType &well_state, const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state) const
Definition: StandardWell_impl.hpp:2262
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:1976
void updateWaterMobilityWithPolymer(const Simulator &simulator, const int perf, std::vector< EvalWell > &mob_water, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1902
bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) override
Definition: StandardWell_impl.hpp:2384
std::vector< Scalar > computeCurrentWellRates(const Simulator &ebosSimulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:2580
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:503
void computeWellConnectionPressures(const Simulator &simulator, const GroupStateHelperType &groupStateHelper)
Definition: StandardWell_impl.hpp:1395
void updatePrimaryVariablesNewton(const BVectorWell &dwells, const bool stop_or_zero_rate_target, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:781
void assembleWellEqWithoutIteration(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellStateType &well_state, const bool solving_with_zero_rate) override
Definition: StandardWell_impl.hpp:341
void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_potentials) override
computing the well potentials for group control
Definition: StandardWell_impl.hpp:1790
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
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state) const
Definition: StandardWell_impl.hpp:2351
typename StdWellEval::StdWellConnections::Properties WellConnectionProps
Definition: StandardWell.hpp:261
void computeWellRatesWithBhpIterations(const Simulator &ebosSimulator, const Scalar &bhp, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_flux) const override
Definition: StandardWell_impl.hpp:1547
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:2220
void computeWellRatesWithBhp(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1499
void updateIPRImplicit(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) override
Definition: StandardWell_impl.hpp:920
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:701
void getTransMult(Value &trans_mult, const Simulator &simulator, const int cell_indx) const
Definition: StandardWell_impl.hpp:681
void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:828
void handleInjectivityEquations(const Simulator &simulator, const WellStateType &well_state, const int perf, const EvalWell &water_flux_s, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:2154
virtual void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: StandardWell_impl.hpp:1447
void checkConvergenceExtraEqs(const std::vector< Scalar > &res, ConvergenceReport &report) const
Definition: StandardWell_impl.hpp:2201
void computeWellRatesWithThpAlqProd(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1773
typename StdWellEval::Eval Eval
Definition: StandardWell.hpp:120
Scalar computeWellRatesAndBhpWithThpAlqProd(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1743
static constexpr int Bhp
Definition: StandardWell.hpp:107
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1158
bool computeWellPotentialsImplicit(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_potentials) const
Definition: StandardWell_impl.hpp:1654
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) override
Definition: StandardWell_impl.hpp:1479
Scalar maxPerfPress(const Simulator &simulator) const override
Definition: StandardWell_impl.hpp:2729
bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, const bool fixed_control, const bool fixed_status, const bool solving_with_zero_rate) override
Definition: StandardWell_impl.hpp:2448
typename FluidSystem::IndexTraitsType IndexTraits
Definition: StandardWell.hpp:124
void assembleWellEqWithoutIterationImpl(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellStateType &well_state, const bool solving_with_zero_rate)
Definition: StandardWell_impl.hpp:366
bool allDrawDownWrongDirection(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1116
EvalWell pskin(const Scalar throughput, const EvalWell &water_velocity, const EvalWell &poly_inj_conc, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2027
void computePerfRate(const IntensiveQuantities &intQuants, const std::vector< Value > &mob, const Value &bhp, const std::vector< Value > &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
static constexpr int numWellConservationEq
Definition: StandardWell.hpp:97
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: StandardWell_impl.hpp:2635
void updateWaterThroughput(const double dt, WellStateType &well_state) const override
Definition: StandardWell_impl.hpp:2100
void checkOperabilityUnderBHPLimit(const WellStateType &well_state, const Simulator &simulator, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:997
EvalWell pskinwater(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1996
void handleInjectivityRate(const Simulator &simulator, const int perf, std::vector< EvalWell > &cq_s) const
Definition: StandardWell_impl.hpp:2131
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:1263
void calculateExplicitQuantities(const Simulator &simulator, const GroupStateHelperType &groupStateHelper) override
Definition: StandardWell_impl.hpp:1434
static constexpr int numWellControlEq
Definition: StandardWell.hpp:99
Scalar getRefDensity() const override
Definition: StandardWell_impl.hpp:1891
void checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper) override
Definition: StandardWell_impl.hpp:1066
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: StandardWell_impl.hpp:1853
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
BlackOilFluidState< Eval, FluidSystem, energyModuleType==EnergyModules::ConstantTemperature,(energyModuleType==EnergyModules::FullyImplicitThermal||energyModuleType==EnergyModules::SequentialImplicitThermal), Indices::compositionSwitchIdx >=0, has_watVapor, has_brine, has_saltPrecip, has_disgas_in_water, Indices::numPhases > FluidState
Definition: WellInterface.hpp:139
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:333
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
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:117
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:86
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:89
GroupStateHelper< Scalar, IndexTraits > GroupStateHelperType
Definition: WellInterface.hpp:102
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