StandardWell_impl.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#ifndef OPM_STANDARDWELL_IMPL_HEADER_INCLUDED
23#define OPM_STANDARDWELL_IMPL_HEADER_INCLUDED
24
25// Improve IDE experience
26#ifndef OPM_STANDARDWELL_HEADER_INCLUDED
27#include <config.h>
29#endif
30
31#include <opm/common/Exceptions.hpp>
32
33#include <opm/input/eclipse/Units/Units.hpp>
34
40
41#include <algorithm>
42#include <cstddef>
43#include <functional>
44
45#include <fmt/format.h>
46
47namespace Opm
48{
49
50 template<typename TypeTag>
52 StandardWell(const Well& well,
53 const ParallelWellInfo<Scalar>& pw_info,
54 const int time_step,
55 const ModelParameters& param,
56 const RateConverterType& rate_converter,
57 const int pvtRegionIdx,
58 const int num_conservation_quantities,
59 const int num_phases,
60 const int index_of_well,
61 const std::vector<PerforationData<Scalar>>& perf_data)
62 : Base(well, pw_info, time_step, param, rate_converter, pvtRegionIdx, num_conservation_quantities, num_phases, index_of_well, perf_data)
63 , StdWellEval(static_cast<const WellInterfaceIndices<FluidSystem,Indices>&>(*this))
64 , regularize_(false)
65 {
67 }
68
69
70
71
72
73 template<typename TypeTag>
74 void
76 init(const std::vector<Scalar>& depth_arg,
77 const Scalar gravity_arg,
78 const std::vector< Scalar >& B_avg,
79 const bool changed_to_open_this_step)
80 {
81 Base::init(depth_arg, gravity_arg, B_avg, changed_to_open_this_step);
82 this->StdWellEval::init(this->perf_depth_, depth_arg, Base::has_polymermw);
83 }
84
85
86
87
88
89 template<typename TypeTag>
90 template<class Value>
91 void
94 const std::vector<Value>& mob,
95 const Value& bhp,
96 const std::vector<Scalar>& Tw,
97 const int perf,
98 const bool allow_cf,
99 std::vector<Value>& cq_s,
100 PerforationRates<Scalar>& perf_rates,
101 DeferredLogger& deferred_logger) const
102 {
103 auto obtain = [this](const Eval& value)
104 {
105 if constexpr (std::is_same_v<Value, Scalar>) {
106 static_cast<void>(this); // suppress clang warning
107 return getValue(value);
108 } else {
109 return this->extendEval(value);
110 }
111 };
112 auto obtainN = [](const auto& value)
113 {
114 if constexpr (std::is_same_v<Value, Scalar>) {
115 return getValue(value);
116 } else {
117 return value;
118 }
119 };
120 auto zeroElem = [this]()
121 {
122 if constexpr (std::is_same_v<Value, Scalar>) {
123 static_cast<void>(this); // suppress clang warning
124 return 0.0;
125 } else {
126 return Value{this->primary_variables_.numWellEq() + Indices::numEq, 0.0};
127 }
128 };
129
130 const auto& fs = intQuants.fluidState();
131 const Value pressure = obtain(this->getPerfCellPressure(fs));
132 const Value rs = obtain(fs.Rs());
133 const Value rv = obtain(fs.Rv());
134 const Value rvw = obtain(fs.Rvw());
135 const Value rsw = obtain(fs.Rsw());
136
137 std::vector<Value> b_perfcells_dense(this->numConservationQuantities(), zeroElem());
138 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
139 if (!FluidSystem::phaseIsActive(phaseIdx)) {
140 continue;
141 }
142 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
143 b_perfcells_dense[compIdx] = obtain(fs.invB(phaseIdx));
144 }
145 if constexpr (has_solvent) {
146 b_perfcells_dense[Indices::contiSolventEqIdx] = obtain(intQuants.solventInverseFormationVolumeFactor());
147 }
148
149 if constexpr (has_zFraction) {
150 if (this->isInjector()) {
151 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
152 b_perfcells_dense[gasCompIdx] *= (1.0 - this->wsolvent());
153 b_perfcells_dense[gasCompIdx] += this->wsolvent()*intQuants.zPureInvFormationVolumeFactor().value();
154 }
155 }
156
157 Value skin_pressure = zeroElem();
158 if (has_polymermw) {
159 if (this->isInjector()) {
160 const int pskin_index = Bhp + 1 + this->numLocalPerfs() + perf;
161 skin_pressure = obtainN(this->primary_variables_.eval(pskin_index));
162 }
163 }
164
165 // surface volume fraction of fluids within wellbore
166 std::vector<Value> cmix_s(this->numConservationQuantities(), zeroElem());
167 for (int componentIdx = 0; componentIdx < this->numConservationQuantities(); ++componentIdx) {
168 cmix_s[componentIdx] = obtainN(this->primary_variables_.surfaceVolumeFraction(componentIdx));
169 }
170
171 computePerfRate(mob,
172 pressure,
173 bhp,
174 rs,
175 rv,
176 rvw,
177 rsw,
178 b_perfcells_dense,
179 Tw,
180 perf,
181 allow_cf,
182 skin_pressure,
183 cmix_s,
184 cq_s,
185 perf_rates,
186 deferred_logger);
187 }
188
189
190
191 template<typename TypeTag>
192 template<class Value>
193 void
195 computePerfRate(const std::vector<Value>& mob,
196 const Value& pressure,
197 const Value& bhp,
198 const Value& rs,
199 const Value& rv,
200 const Value& rvw,
201 const Value& rsw,
202 std::vector<Value>& b_perfcells_dense,
203 const std::vector<Scalar>& Tw,
204 const int perf,
205 const bool allow_cf,
206 const Value& skin_pressure,
207 const std::vector<Value>& cmix_s,
208 std::vector<Value>& cq_s,
209 PerforationRates<Scalar>& perf_rates,
210 DeferredLogger& deferred_logger) const
211 {
212 // Pressure drawdown (also used to determine direction of flow)
213 const Value well_pressure = bhp + this->connections_.pressure_diff(perf);
214 Value drawdown = pressure - well_pressure;
215 if (this->isInjector()) {
216 drawdown += skin_pressure;
217 }
218
219 RatioCalculator<Value> ratioCalc{
220 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)
221 ? FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx)
222 : -1,
223 FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)
224 ? FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx)
225 : -1,
226 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
227 ? FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx)
228 : -1,
229 this->name()
230 };
231
232 // producing perforations
233 if (drawdown > 0) {
234 // Do nothing if crossflow is not allowed
235 if (!allow_cf && this->isInjector()) {
236 return;
237 }
238
239 // compute component volumetric rates at standard conditions
240 for (int componentIdx = 0; componentIdx < this->numConservationQuantities(); ++componentIdx) {
241 const Value cq_p = - Tw[componentIdx] * (mob[componentIdx] * drawdown);
242 cq_s[componentIdx] = b_perfcells_dense[componentIdx] * cq_p;
243 }
244
245 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
246 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
247 {
248 ratioCalc.gasOilPerfRateProd(cq_s, perf_rates, rv, rs, rvw,
249 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx),
250 this->isProducer());
251 } else if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) &&
252 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
253 {
254 ratioCalc.gasWaterPerfRateProd(cq_s, perf_rates, rvw, rsw, this->isProducer());
255 }
256 } else {
257 // Do nothing if crossflow is not allowed
258 if (!allow_cf && this->isProducer()) {
259 return;
260 }
261
262 // Using total mobilities
263 Value total_mob_dense = mob[0];
264 for (int componentIdx = 1; componentIdx < this->numConservationQuantities(); ++componentIdx) {
265 total_mob_dense += mob[componentIdx];
266 }
267
268 // compute volume ratio between connection at standard conditions
269 Value volumeRatio = bhp * 0.0; // initialize it with the correct type
270
271 if (FluidSystem::enableVaporizedWater() && FluidSystem::enableDissolvedGasInWater()) {
272 ratioCalc.disOilVapWatVolumeRatio(volumeRatio, rvw, rsw, pressure,
273 cmix_s, b_perfcells_dense, deferred_logger);
274 // DISGASW only supported for gas-water CO2STORE/H2STORE case
275 // and the simulator will throw long before it reach to this point in the code
276 // For blackoil support of DISGASW we need to add the oil component here
277 assert(FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx));
278 assert(FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx));
279 assert(!FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx));
280 } else {
281
282 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
283 const unsigned waterCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
284 volumeRatio += cmix_s[waterCompIdx] / b_perfcells_dense[waterCompIdx];
285 }
286
287 if constexpr (Indices::enableSolvent) {
288 volumeRatio += cmix_s[Indices::contiSolventEqIdx] / b_perfcells_dense[Indices::contiSolventEqIdx];
289 }
290
291 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
292 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
293 {
294 ratioCalc.gasOilVolumeRatio(volumeRatio, rv, rs, pressure,
295 cmix_s, b_perfcells_dense,
296 deferred_logger);
297 } else {
298 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
299 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
300 volumeRatio += cmix_s[oilCompIdx] / b_perfcells_dense[oilCompIdx];
301 }
302 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
303 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
304 volumeRatio += cmix_s[gasCompIdx] / b_perfcells_dense[gasCompIdx];
305 }
306 }
307 }
308
309 // injecting connections total volumerates at standard conditions
310 for (int componentIdx = 0; componentIdx < this->numConservationQuantities(); ++componentIdx) {
311 const Value cqt_i = - Tw[componentIdx] * (total_mob_dense * drawdown);
312 Value cqt_is = cqt_i / volumeRatio;
313 cq_s[componentIdx] = cmix_s[componentIdx] * cqt_is;
314 }
315
316 // calculating the perforation solution gas rate and solution oil rates
317 if (this->isProducer()) {
318 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
319 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
320 {
321 ratioCalc.gasOilPerfRateInj(cq_s, perf_rates,
322 rv, rs, pressure, rvw,
323 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx),
324 deferred_logger);
325 }
326 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) &&
327 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx))
328 {
329 //no oil
330 ratioCalc.gasWaterPerfRateInj(cq_s, perf_rates, rvw, rsw,
331 pressure, deferred_logger);
332 }
333 }
334 }
335 }
336
337
338 template<typename TypeTag>
339 void
342 const double dt,
343 const Well::InjectionControls& inj_controls,
344 const Well::ProductionControls& prod_controls,
345 WellStateType& well_state,
346 const GroupState<Scalar>& group_state,
347 DeferredLogger& deferred_logger)
348 {
349 // TODO: only_wells should be put back to save some computation
350 // for example, the matrices B C does not need to update if only_wells
351 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
352
353 // clear all entries
354 this->linSys_.clear();
355
356 assembleWellEqWithoutIterationImpl(simulator, dt, inj_controls,
357 prod_controls, well_state,
358 group_state, deferred_logger);
359 }
360
361
362
363
364 template<typename TypeTag>
365 void
368 const double dt,
369 const Well::InjectionControls& inj_controls,
370 const Well::ProductionControls& prod_controls,
371 WellStateType& well_state,
372 const GroupState<Scalar>& group_state,
373 DeferredLogger& deferred_logger)
374 {
375 // try to regularize equation if the well does not converge
376 const Scalar regularization_factor = this->regularize_? this->param_.regularization_factor_wells_ : 1.0;
377 const Scalar volume = 0.1 * unit::cubic(unit::feet) * regularization_factor;
378
379 auto& ws = well_state.well(this->index_of_well_);
380 ws.phase_mixing_rates.fill(0.0);
381 if constexpr (has_energy) {
382 ws.energy_rate = 0.0;
383 }
384
385
386 const int np = this->number_of_phases_;
387
388 std::vector<RateVector> connectionRates = this->connectionRates_; // Copy to get right size.
389
390 auto& perf_data = ws.perf_data;
391 auto& perf_rates = perf_data.phase_rates;
392 for (int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
393 // Calculate perforation quantities.
394 std::vector<EvalWell> cq_s(this->num_conservation_quantities_, {this->primary_variables_.numWellEq() + Indices::numEq, 0.0});
395 EvalWell water_flux_s{this->primary_variables_.numWellEq() + Indices::numEq, 0.0};
396 EvalWell cq_s_zfrac_effective{this->primary_variables_.numWellEq() + Indices::numEq, 0.0};
397 calculateSinglePerf(simulator, perf, well_state, connectionRates,
398 cq_s, water_flux_s, cq_s_zfrac_effective, deferred_logger);
399
400 // Equation assembly for this perforation.
401 if constexpr (has_polymer && Base::has_polymermw) {
402 if (this->isInjector()) {
403 handleInjectivityEquations(simulator, well_state, perf,
404 water_flux_s, deferred_logger);
405 }
406 }
407 for (int componentIdx = 0; componentIdx < this->num_conservation_quantities_; ++componentIdx) {
408 // the cq_s entering mass balance equations need to consider the efficiency factors.
409 const EvalWell cq_s_effective = cq_s[componentIdx] * this->well_efficiency_factor_;
410
411 connectionRates[perf][componentIdx] = Base::restrictEval(cq_s_effective);
412
414 assemblePerforationEq(cq_s_effective,
415 componentIdx,
416 perf,
417 this->primary_variables_.numWellEq(),
418 this->linSys_);
419
420 // Store the perforation phase flux for later usage.
421 if (has_solvent && componentIdx == Indices::contiSolventEqIdx) {
422 auto& perf_rate_solvent = perf_data.solvent_rates;
423 perf_rate_solvent[perf] = cq_s[componentIdx].value();
424 } else {
425 perf_rates[perf*np + FluidSystem::activeCompToActivePhaseIdx(componentIdx)] = cq_s[componentIdx].value();
426 }
427 }
428
429 if constexpr (has_zFraction) {
431 assembleZFracEq(cq_s_zfrac_effective,
432 perf,
433 this->primary_variables_.numWellEq(),
434 this->linSys_);
435 }
436 }
437 // Update the connection
438 this->connectionRates_ = connectionRates;
439
440 // Accumulate dissolved gas and vaporized oil flow rates across all
441 // ranks sharing this well (this->index_of_well_).
442 {
443 const auto& comm = this->parallel_well_info_.communication();
444 comm.sum(ws.phase_mixing_rates.data(), ws.phase_mixing_rates.size());
445 }
446
447 // accumulate resWell_ and duneD_ in parallel to get effects of all perforations (might be distributed)
448 this->linSys_.sumDistributed(this->parallel_well_info_.communication());
449
450 // add vol * dF/dt + Q to the well equations;
451 for (int componentIdx = 0; componentIdx < numWellConservationEq; ++componentIdx) {
452 // TODO: following the development in MSW, we need to convert the volume of the wellbore to be surface volume
453 // since all the rates are under surface condition
454 EvalWell resWell_loc(this->primary_variables_.numWellEq() + Indices::numEq, 0.0);
455 if (FluidSystem::numActivePhases() > 1) {
456 assert(dt > 0);
457 resWell_loc += (this->primary_variables_.surfaceVolumeFraction(componentIdx) -
458 this->F0_[componentIdx]) * volume / dt;
459 }
460 resWell_loc -= this->primary_variables_.getQs(componentIdx) * this->well_efficiency_factor_;
462 assembleSourceEq(resWell_loc,
463 componentIdx,
464 this->primary_variables_.numWellEq(),
465 this->linSys_);
466 }
467
468 const auto& summaryState = simulator.vanguard().summaryState();
469 const Schedule& schedule = simulator.vanguard().schedule();
470 const bool stopped_or_zero_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
472 assembleControlEq(well_state, group_state,
473 schedule, summaryState,
474 inj_controls, prod_controls,
475 this->primary_variables_,
476 this->getRefDensity(),
477 this->linSys_,
478 stopped_or_zero_target,
479 deferred_logger);
480
481
482 // do the local inversion of D.
483 try {
484 this->linSys_.invert();
485 } catch( ... ) {
486 OPM_DEFLOG_PROBLEM(NumericalProblem, "Error when inverting local well equations for well " + name(), deferred_logger);
487 }
488 }
489
490
491
492
493 template<typename TypeTag>
494 void
496 calculateSinglePerf(const Simulator& simulator,
497 const int perf,
498 WellStateType& well_state,
499 std::vector<RateVector>& connectionRates,
500 std::vector<EvalWell>& cq_s,
501 EvalWell& water_flux_s,
502 EvalWell& cq_s_zfrac_effective,
503 DeferredLogger& deferred_logger) const
504 {
505 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
506 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
507 const int cell_idx = this->well_cells_[perf];
508 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
509 std::vector<EvalWell> mob(this->num_conservation_quantities_, {this->primary_variables_.numWellEq() + Indices::numEq, 0.});
510 getMobility(simulator, perf, mob, deferred_logger);
511
512 PerforationRates<Scalar> perf_rates;
513 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(intQuants, cell_idx);
514 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
515 const std::vector<Scalar> Tw = this->wellIndex(perf, intQuants, trans_mult, wellstate_nupcol);
516 computePerfRate(intQuants, mob, bhp, Tw, perf, allow_cf,
517 cq_s, perf_rates, deferred_logger);
518
519 auto& ws = well_state.well(this->index_of_well_);
520 auto& perf_data = ws.perf_data;
521 if constexpr (has_polymer && Base::has_polymermw) {
522 if (this->isInjector()) {
523 // Store the original water flux computed from the reservoir quantities.
524 // It will be required to assemble the injectivity equations.
525 const unsigned water_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
526 water_flux_s = cq_s[water_comp_idx];
527 // Modify the water flux for the rest of this function to depend directly on the
528 // local water velocity primary variable.
529 handleInjectivityRate(simulator, perf, cq_s);
530 }
531 }
532
533 // updating the solution gas rate and solution oil rate
534 if (this->isProducer()) {
535 ws.phase_mixing_rates[ws.dissolved_gas] += perf_rates.dis_gas;
536 ws.phase_mixing_rates[ws.dissolved_gas_in_water] += perf_rates.dis_gas_in_water;
537 ws.phase_mixing_rates[ws.vaporized_oil] += perf_rates.vap_oil;
538 ws.phase_mixing_rates[ws.vaporized_water] += perf_rates.vap_wat;
539 perf_data.phase_mixing_rates[perf][ws.dissolved_gas] = perf_rates.dis_gas;
540 perf_data.phase_mixing_rates[perf][ws.dissolved_gas_in_water] = perf_rates.dis_gas_in_water;
541 perf_data.phase_mixing_rates[perf][ws.vaporized_oil] = perf_rates.vap_oil;
542 perf_data.phase_mixing_rates[perf][ws.vaporized_water] = perf_rates.vap_wat;
543 }
544
545 if constexpr (has_energy) {
546 connectionRates[perf][Indices::contiEnergyEqIdx] =
547 connectionRateEnergy(cq_s, intQuants, deferred_logger);
548 ws.energy_rate += getValue(connectionRates[perf][Indices::contiEnergyEqIdx]);
549 }
550
551 if constexpr (has_polymer) {
552 std::variant<Scalar,EvalWell> polymerConcentration;
553 if (this->isInjector()) {
554 polymerConcentration = this->wpolymer();
555 } else {
556 polymerConcentration = this->extendEval(intQuants.polymerConcentration() *
557 intQuants.polymerViscosityCorrection());
558 }
559
560 [[maybe_unused]] EvalWell cq_s_poly;
561 std::tie(connectionRates[perf][Indices::contiPolymerEqIdx],
562 cq_s_poly) =
563 this->connections_.connectionRatePolymer(perf_data.polymer_rates[perf],
564 cq_s, polymerConcentration);
565
566 if constexpr (Base::has_polymermw) {
567 updateConnectionRatePolyMW(cq_s_poly, intQuants, well_state,
568 perf, connectionRates, deferred_logger);
569 }
570 }
571
572 if constexpr (has_foam) {
573 std::variant<Scalar,EvalWell> foamConcentration;
574 if (this->isInjector()) {
575 foamConcentration = this->wfoam();
576 } else {
577 foamConcentration = this->extendEval(intQuants.foamConcentration());
578 }
579 connectionRates[perf][Indices::contiFoamEqIdx] =
580 this->connections_.connectionRateFoam(cq_s, foamConcentration,
581 FoamModule::transportPhase(),
582 deferred_logger);
583 }
584
585 if constexpr (has_zFraction) {
586 std::variant<Scalar,std::array<EvalWell,2>> solventConcentration;
587 if (this->isInjector()) {
588 solventConcentration = this->wsolvent();
589 } else {
590 solventConcentration = std::array{this->extendEval(intQuants.xVolume()),
591 this->extendEval(intQuants.yVolume())};
592 }
593 std::tie(connectionRates[perf][Indices::contiZfracEqIdx],
594 cq_s_zfrac_effective) =
595 this->connections_.connectionRatezFraction(perf_data.solvent_rates[perf],
596 perf_rates.dis_gas, cq_s,
597 solventConcentration);
598 }
599
600 if constexpr (has_brine) {
601 std::variant<Scalar,EvalWell> saltConcentration;
602 if (this->isInjector()) {
603 saltConcentration = this->wsalt();
604 } else {
605 saltConcentration = this->extendEval(intQuants.fluidState().saltConcentration());
606 }
607
608 connectionRates[perf][Indices::contiBrineEqIdx] =
609 this->connections_.connectionRateBrine(perf_data.brine_rates[perf],
610 perf_rates.vap_wat, cq_s,
611 saltConcentration);
612 }
613
614 if constexpr (has_bioeffects) {
615 std::variant<Scalar,EvalWell> microbialConcentration;
616 if constexpr (has_micp) {
617 std::variant<Scalar,EvalWell> oxygenConcentration;
618 std::variant<Scalar,EvalWell> ureaConcentration;
619 if (this->isInjector()) {
620 microbialConcentration = this->wmicrobes();
621 oxygenConcentration = this->woxygen();
622 ureaConcentration = this->wurea();
623 } else {
624 microbialConcentration = this->extendEval(intQuants.microbialConcentration());
625 oxygenConcentration = this->extendEval(intQuants.oxygenConcentration());
626 ureaConcentration = this->extendEval(intQuants.ureaConcentration());
627 }
628 std::tie(connectionRates[perf][Indices::contiMicrobialEqIdx],
629 connectionRates[perf][Indices::contiOxygenEqIdx],
630 connectionRates[perf][Indices::contiUreaEqIdx]) =
631 this->connections_.connectionRatesMICP(perf_data.microbial_rates[perf],
632 perf_data.oxygen_rates[perf],
633 perf_data.urea_rates[perf],
634 cq_s,
635 microbialConcentration,
636 oxygenConcentration,
637 ureaConcentration);
638 }
639 else {
640 if (this->isProducer()) {
641 microbialConcentration = this->extendEval(intQuants.microbialConcentration());
642 connectionRates[perf][Indices::contiMicrobialEqIdx] =
643 this->connections_.connectionRateBioeffects(perf_data.microbial_rates[perf],
644 perf_rates.vap_wat, cq_s,
645 microbialConcentration);
646 }
647 }
648 }
649
650 // Store the perforation pressure for later usage.
651 perf_data.pressure[perf] = ws.bhp + this->connections_.pressure_diff(perf);
652
653 // Store the perforation gass mass rate.
654 if (FluidSystem::phaseUsage().hasCO2orH2Store()) {
655 const unsigned gas_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
656 const Scalar rho = FluidSystem::referenceDensity( FluidSystem::gasPhaseIdx, Base::pvtRegionIdx() );
657 perf_data.gas_mass_rates[perf] = cq_s[gas_comp_idx].value() * rho;
658 }
659
660 // Store the perforation water mass rate.
661 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
662 const unsigned wat_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
663 const Scalar rho = FluidSystem::referenceDensity( FluidSystem::waterPhaseIdx, Base::pvtRegionIdx() );
664 perf_data.wat_mass_rates[perf] = cq_s[wat_comp_idx].value() * rho;
665 }
666 }
667
668
669
670 template<typename TypeTag>
671 template<class Value>
672 void
674 getMobility(const Simulator& simulator,
675 const int perf,
676 std::vector<Value>& mob,
677 DeferredLogger& deferred_logger) const
678 {
679 auto obtain = [this](const Eval& value)
680 {
681 if constexpr (std::is_same_v<Value, Scalar>) {
682 static_cast<void>(this); // suppress clang warning
683 return getValue(value);
684 } else {
685 return this->extendEval(value);
686 }
687 };
688 WellInterface<TypeTag>::getMobility(simulator, perf, mob,
689 obtain, deferred_logger);
690
691 // modify the water mobility if polymer is present
692 if constexpr (has_polymer) {
693 if (!FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
694 OPM_DEFLOG_THROW(std::runtime_error, "Water is required when polymer is active", deferred_logger);
695 }
696
697 // for the cases related to polymer molecular weight, we assume fully mixing
698 // as a result, the polymer and water share the same viscosity
699 if constexpr (!Base::has_polymermw) {
700 if constexpr (std::is_same_v<Value, Scalar>) {
701 std::vector<EvalWell> mob_eval(this->num_conservation_quantities_, {this->primary_variables_.numWellEq() + Indices::numEq, 0.});
702 for (std::size_t i = 0; i < mob.size(); ++i) {
703 mob_eval[i].setValue(mob[i]);
704 }
705 updateWaterMobilityWithPolymer(simulator, perf, mob_eval, deferred_logger);
706 for (std::size_t i = 0; i < mob.size(); ++i) {
707 mob[i] = getValue(mob_eval[i]);
708 }
709 } else {
710 updateWaterMobilityWithPolymer(simulator, perf, mob, deferred_logger);
711 }
712 }
713 }
714
715 // if the injecting well has WINJMULT setup, we update the mobility accordingly
716 if (this->isInjector() && this->well_ecl_.getInjMultMode() != Well::InjMultMode::NONE) {
717 const Scalar bhp = this->primary_variables_.value(Bhp);
718 const Scalar perf_press = bhp + this->connections_.pressure_diff(perf);
719 const Scalar multiplier = this->getInjMult(perf, bhp, perf_press, deferred_logger);
720 for (std::size_t i = 0; i < mob.size(); ++i) {
721 mob[i] *= multiplier;
722 }
723 }
724 }
725
726
727 template<typename TypeTag>
728 void
730 updateWellState(const Simulator& simulator,
731 const BVectorWell& dwells,
732 WellStateType& well_state,
733 DeferredLogger& deferred_logger)
734 {
735 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
736
737 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
738 updatePrimaryVariablesNewton(dwells, stop_or_zero_rate_target, deferred_logger);
739
740 const auto& summary_state = simulator.vanguard().summaryState();
741 updateWellStateFromPrimaryVariables(well_state, summary_state, deferred_logger);
742 Base::calculateReservoirRates(simulator.vanguard().eclState().runspec().co2Storage(), well_state.well(this->index_of_well_));
743 }
744
745
746
747
748
749 template<typename TypeTag>
750 void
753 const bool stop_or_zero_rate_target,
754 DeferredLogger& deferred_logger)
755 {
756 const Scalar dFLimit = this->param_.dwell_fraction_max_;
757 const Scalar dBHPLimit = this->param_.dbhp_max_rel_;
758 this->primary_variables_.updateNewton(dwells, stop_or_zero_rate_target, dFLimit, dBHPLimit, deferred_logger);
759
760 // for the water velocity and skin pressure
761 if constexpr (Base::has_polymermw) {
762 this->primary_variables_.updateNewtonPolyMW(dwells);
763 }
764
765 this->primary_variables_.checkFinite(deferred_logger);
766 }
767
768
769
770
771
772 template<typename TypeTag>
773 void
776 const SummaryState& summary_state,
777 DeferredLogger& deferred_logger) const
778 {
779 this->primary_variables_.copyToWellState(well_state, deferred_logger);
780
781 WellBhpThpCalculator(this->baseif_).
782 updateThp(getRefDensity(),
783 [this,&well_state]() { return this->baseif_.getALQ(well_state); },
784 well_state, summary_state, deferred_logger);
785
786 // other primary variables related to polymer injectivity study
787 if constexpr (Base::has_polymermw) {
788 this->primary_variables_.copyToWellStatePolyMW(well_state);
789 }
790 }
791
792
793
794
795
796 template<typename TypeTag>
797 void
799 updateIPR(const Simulator& simulator, DeferredLogger& deferred_logger) const
800 {
801 // TODO: not handling solvent related here for now
802
803 // initialize all the values to be zero to begin with
804 std::fill(this->ipr_a_.begin(), this->ipr_a_.end(), 0.);
805 std::fill(this->ipr_b_.begin(), this->ipr_b_.end(), 0.);
806
807 for (int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
808 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
809 getMobility(simulator, perf, mob, deferred_logger);
810
811 const int cell_idx = this->well_cells_[perf];
812 const auto& int_quantities = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
813 const auto& fs = int_quantities.fluidState();
814 // the pressure of the reservoir grid block the well connection is in
815 Scalar p_r = this->getPerfCellPressure(fs).value();
816
817 // calculating the b for the connection
818 std::vector<Scalar> b_perf(this->num_conservation_quantities_);
819 for (std::size_t phase = 0; phase < FluidSystem::numPhases; ++phase) {
820 if (!FluidSystem::phaseIsActive(phase)) {
821 continue;
822 }
823 const unsigned comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phase));
824 b_perf[comp_idx] = fs.invB(phase).value();
825 }
826 if constexpr (has_solvent) {
827 b_perf[Indices::contiSolventEqIdx] = int_quantities.solventInverseFormationVolumeFactor().value();
828 }
829
830 // the pressure difference between the connection and BHP
831 const Scalar h_perf = this->connections_.pressure_diff(perf);
832 const Scalar pressure_diff = p_r - h_perf;
833
834 // Let us add a check, since the pressure is calculated based on zero value BHP
835 // it should not be negative anyway. If it is negative, we might need to re-formulate
836 // to taking into consideration the crossflow here.
837 if ( (this->isProducer() && pressure_diff < 0.) || (this->isInjector() && pressure_diff > 0.) ) {
838 deferred_logger.debug("CROSSFLOW_IPR",
839 "cross flow found when updateIPR for well " + name()
840 + " . The connection is ignored in IPR calculations");
841 // we ignore these connections for now
842 continue;
843 }
844
845 // the well index associated with the connection
846 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quantities, cell_idx);
847 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
848 const std::vector<Scalar> tw_perf = this->wellIndex(perf,
849 int_quantities,
850 trans_mult,
851 wellstate_nupcol);
852 std::vector<Scalar> ipr_a_perf(this->ipr_a_.size());
853 std::vector<Scalar> ipr_b_perf(this->ipr_b_.size());
854 for (int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
855 const Scalar tw_mob = tw_perf[comp_idx] * mob[comp_idx] * b_perf[comp_idx];
856 ipr_a_perf[comp_idx] += tw_mob * pressure_diff;
857 ipr_b_perf[comp_idx] += tw_mob;
858 }
859
860 // we need to handle the rs and rv when both oil and gas are present
861 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
862 const unsigned oil_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
863 const unsigned gas_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
864 const Scalar rs = (fs.Rs()).value();
865 const Scalar rv = (fs.Rv()).value();
866
867 const Scalar dis_gas_a = rs * ipr_a_perf[oil_comp_idx];
868 const Scalar vap_oil_a = rv * ipr_a_perf[gas_comp_idx];
869
870 ipr_a_perf[gas_comp_idx] += dis_gas_a;
871 ipr_a_perf[oil_comp_idx] += vap_oil_a;
872
873 const Scalar dis_gas_b = rs * ipr_b_perf[oil_comp_idx];
874 const Scalar vap_oil_b = rv * ipr_b_perf[gas_comp_idx];
875
876 ipr_b_perf[gas_comp_idx] += dis_gas_b;
877 ipr_b_perf[oil_comp_idx] += vap_oil_b;
878 }
879
880 for (std::size_t comp_idx = 0; comp_idx < ipr_a_perf.size(); ++comp_idx) {
881 this->ipr_a_[comp_idx] += ipr_a_perf[comp_idx];
882 this->ipr_b_[comp_idx] += ipr_b_perf[comp_idx];
883 }
884 }
885 this->parallel_well_info_.communication().sum(this->ipr_a_.data(), this->ipr_a_.size());
886 this->parallel_well_info_.communication().sum(this->ipr_b_.data(), this->ipr_b_.size());
887 }
888
889 template<typename TypeTag>
890 void
892 updateIPRImplicit(const Simulator& simulator,
893 WellStateType& well_state,
894 DeferredLogger& deferred_logger)
895 {
896 // Compute IPR based on *converged* well-equation:
897 // For a component rate r the derivative dr/dbhp is obtained by
898 // dr/dbhp = - (partial r/partial x) * inv(partial Eq/partial x) * (partial Eq/partial bhp_target)
899 // where Eq(x)=0 is the well equation setup with bhp control and primary variables x
900
901 // We shouldn't have zero rates at this stage, but check
902 bool zero_rates;
903 auto rates = well_state.well(this->index_of_well_).surface_rates;
904 zero_rates = true;
905 for (std::size_t p = 0; p < rates.size(); ++p) {
906 zero_rates &= rates[p] == 0.0;
907 }
908 auto& ws = well_state.well(this->index_of_well_);
909 if (zero_rates) {
910 const auto msg = fmt::format("updateIPRImplicit: Well {} has zero rate, IPRs might be problematic", this->name());
911 deferred_logger.debug(msg);
912 /*
913 // could revert to standard approach here:
914 updateIPR(simulator, deferred_logger);
915 for (int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx){
916 const int idx = this->activeCompToActivePhaseIdx(comp_idx);
917 ws.implicit_ipr_a[idx] = this->ipr_a_[comp_idx];
918 ws.implicit_ipr_b[idx] = this->ipr_b_[comp_idx];
919 }
920 return;
921 */
922 }
923 const auto& group_state = simulator.problem().wellModel().groupState();
924
925 std::fill(ws.implicit_ipr_a.begin(), ws.implicit_ipr_a.end(), 0.);
926 std::fill(ws.implicit_ipr_b.begin(), ws.implicit_ipr_b.end(), 0.);
927
928 auto inj_controls = Well::InjectionControls(0);
929 auto prod_controls = Well::ProductionControls(0);
930 prod_controls.addControl(Well::ProducerCMode::BHP);
931 prod_controls.bhp_limit = well_state.well(this->index_of_well_).bhp;
932
933 // Set current control to bhp, and bhp value in state, modify bhp limit in control object.
934 const auto cmode = ws.production_cmode;
935 ws.production_cmode = Well::ProducerCMode::BHP;
936 const double dt = simulator.timeStepSize();
937 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
938
939 const size_t nEq = this->primary_variables_.numWellEq();
940 BVectorWell rhs(1);
941 rhs[0].resize(nEq);
942 // rhs = 0 except -1 for control eq
943 for (size_t i=0; i < nEq; ++i){
944 rhs[0][i] = 0.0;
945 }
946 rhs[0][Bhp] = -1.0;
947
948 BVectorWell x_well(1);
949 x_well[0].resize(nEq);
950 this->linSys_.solve(rhs, x_well);
951
952 for (int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx){
953 EvalWell comp_rate = this->primary_variables_.getQs(comp_idx);
954 const int idx = FluidSystem::activeCompToActivePhaseIdx(comp_idx);
955 for (size_t pvIdx = 0; pvIdx < nEq; ++pvIdx) {
956 // well primary variable derivatives in EvalWell start at position Indices::numEq
957 ws.implicit_ipr_b[idx] -= x_well[0][pvIdx]*comp_rate.derivative(pvIdx+Indices::numEq);
958 }
959 ws.implicit_ipr_a[idx] = ws.implicit_ipr_b[idx]*ws.bhp - comp_rate.value();
960 }
961 // reset cmode
962 ws.production_cmode = cmode;
963 }
964
965 template<typename TypeTag>
966 void
969 const Simulator& simulator,
970 DeferredLogger& deferred_logger)
971 {
972 const auto& summaryState = simulator.vanguard().summaryState();
973 const Scalar bhp_limit = WellBhpThpCalculator(*this).mostStrictBhpFromBhpLimits(summaryState);
974 // Crude but works: default is one atmosphere.
975 // TODO: a better way to detect whether the BHP is defaulted or not
976 const bool bhp_limit_not_defaulted = bhp_limit > 1.5 * unit::barsa;
977 if ( bhp_limit_not_defaulted || !this->wellHasTHPConstraints(summaryState) ) {
978 // if the BHP limit is not defaulted or the well does not have a THP limit
979 // we need to check the BHP limit
980 Scalar total_ipr_mass_rate = 0.0;
981 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
982 {
983 if (!FluidSystem::phaseIsActive(phaseIdx)) {
984 continue;
985 }
986
987 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
988 const Scalar ipr_rate = this->ipr_a_[compIdx] - this->ipr_b_[compIdx] * bhp_limit;
989
990 const Scalar rho = FluidSystem::referenceDensity( phaseIdx, Base::pvtRegionIdx() );
991 total_ipr_mass_rate += ipr_rate * rho;
992 }
993 if ( (this->isProducer() && total_ipr_mass_rate < 0.) || (this->isInjector() && total_ipr_mass_rate > 0.) ) {
994 this->operability_status_.operable_under_only_bhp_limit = false;
995 }
996
997 // checking whether running under BHP limit will violate THP limit
998 if (this->operability_status_.operable_under_only_bhp_limit && this->wellHasTHPConstraints(summaryState)) {
999 // option 1: calculate well rates based on the BHP limit.
1000 // option 2: stick with the above IPR curve
1001 // we use IPR here
1002 std::vector<Scalar> well_rates_bhp_limit;
1003 computeWellRatesWithBhp(simulator, bhp_limit, well_rates_bhp_limit, deferred_logger);
1004
1005 this->adaptRatesForVFP(well_rates_bhp_limit);
1006 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1007 const Scalar thp = WellBhpThpCalculator(*this).calculateThpFromBhp(well_rates_bhp_limit,
1008 bhp_limit,
1009 this->getRefDensity(),
1010 this->getALQ(well_state),
1011 thp_limit,
1012 deferred_logger);
1013 if ( (this->isProducer() && thp < thp_limit) || (this->isInjector() && thp > thp_limit) ) {
1014 this->operability_status_.obey_thp_limit_under_bhp_limit = false;
1015 }
1016 }
1017 } else {
1018 // defaulted BHP and there is a THP constraint
1019 // default BHP limit is about 1 atm.
1020 // when applied the hydrostatic pressure correction dp,
1021 // most likely we get a negative value (bhp + dp)to search in the VFP table,
1022 // which is not desirable.
1023 // we assume we can operate under defaulted BHP limit and will violate the THP limit
1024 // when operating under defaulted BHP limit.
1025 this->operability_status_.operable_under_only_bhp_limit = true;
1026 this->operability_status_.obey_thp_limit_under_bhp_limit = false;
1027 }
1028 }
1029
1030
1031
1032
1033
1034 template<typename TypeTag>
1035 void
1038 const WellStateType& well_state,
1039 DeferredLogger& deferred_logger)
1040 {
1041 const auto& summaryState = simulator.vanguard().summaryState();
1042 const auto obtain_bhp = this->isProducer() ? computeBhpAtThpLimitProd(well_state, simulator, summaryState, deferred_logger)
1043 : computeBhpAtThpLimitInj(simulator, summaryState, deferred_logger);
1044
1045 if (obtain_bhp) {
1046 this->operability_status_.can_obtain_bhp_with_thp_limit = true;
1047
1048 const Scalar bhp_limit = WellBhpThpCalculator(*this).mostStrictBhpFromBhpLimits(summaryState);
1049 this->operability_status_.obey_bhp_limit_with_thp_limit = this->isProducer() ?
1050 *obtain_bhp >= bhp_limit : *obtain_bhp <= bhp_limit ;
1051
1052 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1053 if (this->isProducer() && *obtain_bhp < thp_limit) {
1054 const std::string msg = " obtained bhp " + std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1055 + " bars is SMALLER than thp limit "
1056 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1057 + " bars as a producer for well " + name();
1058 deferred_logger.debug(msg);
1059 }
1060 else if (this->isInjector() && *obtain_bhp > thp_limit) {
1061 const std::string msg = " obtained bhp " + std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1062 + " bars is LARGER than thp limit "
1063 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1064 + " bars as a injector for well " + name();
1065 deferred_logger.debug(msg);
1066 }
1067 } else {
1068 this->operability_status_.can_obtain_bhp_with_thp_limit = false;
1069 this->operability_status_.obey_bhp_limit_with_thp_limit = false;
1070 if (!this->wellIsStopped()) {
1071 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1072 deferred_logger.debug(" could not find bhp value at thp limit "
1073 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1074 + " bar for well " + name() + ", the well might need to be closed ");
1075 }
1076 }
1077 }
1078
1079
1080
1081
1082
1083 template<typename TypeTag>
1084 bool
1086 allDrawDownWrongDirection(const Simulator& simulator) const
1087 {
1088 bool all_drawdown_wrong_direction = true;
1089
1090 for (int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
1091 const int cell_idx = this->well_cells_[perf];
1092 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/0);
1093 const auto& fs = intQuants.fluidState();
1094
1095 const Scalar pressure = this->getPerfCellPressure(fs).value();
1096 const Scalar bhp = this->primary_variables_.eval(Bhp).value();
1097
1098 // Pressure drawdown (also used to determine direction of flow)
1099 const Scalar well_pressure = bhp + this->connections_.pressure_diff(perf);
1100 const Scalar drawdown = pressure - well_pressure;
1101
1102 // for now, if there is one perforation can produce/inject in the correct
1103 // direction, we consider this well can still produce/inject.
1104 // TODO: it can be more complicated than this to cause wrong-signed rates
1105 if ( (drawdown < 0. && this->isInjector()) ||
1106 (drawdown > 0. && this->isProducer()) ) {
1107 all_drawdown_wrong_direction = false;
1108 break;
1109 }
1110 }
1111
1112 const auto& comm = this->parallel_well_info_.communication();
1113 if (comm.size() > 1)
1114 {
1115 all_drawdown_wrong_direction =
1116 (comm.min(all_drawdown_wrong_direction ? 1 : 0) == 1);
1117 }
1118
1119 return all_drawdown_wrong_direction;
1120 }
1121
1122
1123
1124
1125 template<typename TypeTag>
1126 bool
1128 openCrossFlowAvoidSingularity(const Simulator& simulator) const
1129 {
1130 return !this->getAllowCrossFlow() && allDrawDownWrongDirection(simulator);
1131 }
1132
1133
1134
1135
1136 template<typename TypeTag>
1140 const WellStateType& well_state) const
1141 {
1142 auto prop_func = typename StdWellEval::StdWellConnections::PressurePropertyFunctions {
1143 // getTemperature
1144 [&model = simulator.model()](int cell_idx, int phase_idx)
1145 {
1146 return model.intensiveQuantities(cell_idx, /* time_idx = */ 0)
1147 .fluidState().temperature(phase_idx).value();
1148 },
1149
1150 // getSaltConcentration
1151 [&model = simulator.model()](int cell_idx)
1152 {
1153 return model.intensiveQuantities(cell_idx, /* time_idx = */ 0)
1154 .fluidState().saltConcentration().value();
1155 },
1156
1157 // getPvtRegionIdx
1158 [&model = simulator.model()](int cell_idx)
1159 {
1160 return model.intensiveQuantities(cell_idx, /* time_idx = */ 0)
1161 .fluidState().pvtRegionIndex();
1162 }
1163 };
1164
1165 if constexpr (Indices::enableSolvent) {
1166 prop_func.solventInverseFormationVolumeFactor =
1167 [&model = simulator.model()](int cell_idx)
1168 {
1169 return model.intensiveQuantities(cell_idx, /* time_idx = */ 0)
1170 .solventInverseFormationVolumeFactor().value();
1171 };
1172
1173 prop_func.solventRefDensity = [&model = simulator.model()](int cell_idx)
1174 {
1175 return model.intensiveQuantities(cell_idx, /* time_idx = */ 0)
1176 .solventRefDensity();
1177 };
1178 }
1179
1180 return this->connections_.computePropertiesForPressures(well_state, prop_func);
1181 }
1182
1183
1184
1185
1186
1187 template<typename TypeTag>
1190 getWellConvergence(const Simulator& simulator,
1191 const WellStateType& well_state,
1192 const std::vector<Scalar>& B_avg,
1193 DeferredLogger& deferred_logger,
1194 const bool relax_tolerance) const
1195 {
1196 // the following implementation assume that the polymer is always after the w-o-g phases
1197 // For the polymer, energy and foam cases, there is one more mass balance equations of reservoir than wells
1198 assert((int(B_avg.size()) == this->num_conservation_quantities_) || has_polymer || has_energy || has_foam || has_brine || has_zFraction || has_bioeffects);
1199
1200 Scalar tol_wells = this->param_.tolerance_wells_;
1201 // use stricter tolerance for stopped wells and wells under zero rate target control.
1202 constexpr Scalar stopped_factor = 1.e-4;
1203 // use stricter tolerance for dynamic thp to ameliorate network convergence
1204 constexpr Scalar dynamic_thp_factor = 1.e-1;
1205 if (this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger)) {
1206 tol_wells = tol_wells*stopped_factor;
1207 } else if (this->getDynamicThpLimit()) {
1208 tol_wells = tol_wells*dynamic_thp_factor;
1209 }
1210
1211 std::vector<Scalar> res;
1212 ConvergenceReport report = this->StdWellEval::getWellConvergence(well_state,
1213 B_avg,
1214 this->param_.max_residual_allowed_,
1215 tol_wells,
1216 this->param_.relaxed_tolerance_flow_well_,
1217 relax_tolerance,
1218 this->wellIsStopped(),
1219 res,
1220 deferred_logger);
1221
1222 checkConvergenceExtraEqs(res, report);
1223
1224 return report;
1225 }
1226
1227
1228
1229
1230
1231 template<typename TypeTag>
1232 void
1234 updateProductivityIndex(const Simulator& simulator,
1235 const WellProdIndexCalculator<Scalar>& wellPICalc,
1236 WellStateType& well_state,
1237 DeferredLogger& deferred_logger) const
1238 {
1239 auto fluidState = [&simulator, this](const int perf)
1240 {
1241 const auto cell_idx = this->well_cells_[perf];
1242 return simulator.model()
1243 .intensiveQuantities(cell_idx, /*timeIdx=*/ 0).fluidState();
1244 };
1245
1246 const int np = this->number_of_phases_;
1247 auto setToZero = [np](Scalar* x) -> void
1248 {
1249 std::fill_n(x, np, 0.0);
1250 };
1251
1252 auto addVector = [np](const Scalar* src, Scalar* dest) -> void
1253 {
1254 std::transform(src, src + np, dest, dest, std::plus<>{});
1255 };
1256
1257 auto& ws = well_state.well(this->index_of_well_);
1258 auto& perf_data = ws.perf_data;
1259 auto* wellPI = ws.productivity_index.data();
1260 auto* connPI = perf_data.prod_index.data();
1261
1262 setToZero(wellPI);
1263
1264 const auto preferred_phase = this->well_ecl_.getPreferredPhase();
1265 auto subsetPerfID = 0;
1266
1267 for (const auto& perf : *this->perf_data_) {
1268 auto allPerfID = perf.ecl_index;
1269
1270 auto connPICalc = [&wellPICalc, allPerfID](const Scalar mobility) -> Scalar
1271 {
1272 return wellPICalc.connectionProdIndStandard(allPerfID, mobility);
1273 };
1274
1275 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
1276 getMobility(simulator, static_cast<int>(subsetPerfID), mob, deferred_logger);
1277
1278 const auto& fs = fluidState(subsetPerfID);
1279 setToZero(connPI);
1280
1281 if (this->isInjector()) {
1282 this->computeConnLevelInjInd(fs, preferred_phase, connPICalc,
1283 mob, connPI, deferred_logger);
1284 }
1285 else { // Production or zero flow rate
1286 this->computeConnLevelProdInd(fs, connPICalc, mob, connPI);
1287 }
1288
1289 addVector(connPI, wellPI);
1290
1291 ++subsetPerfID;
1292 connPI += np;
1293 }
1294
1295 // Sum with communication in case of distributed well.
1296 const auto& comm = this->parallel_well_info_.communication();
1297 if (comm.size() > 1) {
1298 comm.sum(wellPI, np);
1299 }
1300
1301 assert ((static_cast<int>(subsetPerfID) == this->number_of_local_perforations_) &&
1302 "Internal logic error in processing connections for PI/II");
1303 }
1304
1305
1306
1307 template<typename TypeTag>
1310 const WellStateType& well_state,
1311 const WellConnectionProps& props,
1312 DeferredLogger& deferred_logger)
1313 {
1314 // Cell level dynamic property call-back functions as fall-back
1315 // option for calculating connection level mixture densities in
1316 // stopped or zero-rate producer wells.
1317 const auto prop_func = typename StdWellEval::StdWellConnections::DensityPropertyFunctions {
1318 // This becomes slightly more palatable with C++20's designated
1319 // initialisers.
1320
1321 // mobility: Phase mobilities in specified cell.
1322 [&model = simulator.model()](const int cell,
1323 const std::vector<int>& phases,
1324 std::vector<Scalar>& mob)
1325 {
1326 const auto& iq = model.intensiveQuantities(cell, /* time_idx = */ 0);
1327
1328 std::transform(phases.begin(), phases.end(), mob.begin(),
1329 [&iq](const int phase) { return iq.mobility(phase).value(); });
1330 },
1331
1332 // densityInCell: Reservoir condition phase densities in
1333 // specified cell.
1334 [&model = simulator.model()](const int cell,
1335 const std::vector<int>& phases,
1336 std::vector<Scalar>& rho)
1337 {
1338 const auto& fs = model.intensiveQuantities(cell, /* time_idx = */ 0).fluidState();
1339
1340 std::transform(phases.begin(), phases.end(), rho.begin(),
1341 [&fs](const int phase) { return fs.density(phase).value(); });
1342 }
1343 };
1344
1345 const auto stopped_or_zero_rate_target = this->
1346 stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
1347
1348 this->connections_
1349 .computeProperties(stopped_or_zero_rate_target, well_state,
1350 prop_func, props, deferred_logger);
1351 // density was updated
1352 cachedRefDensity = this->connections_.rho(0);
1353 if (this->parallel_well_info_.communication().size() > 1) {
1354 cachedRefDensity = this->parallel_well_info_.broadcastFirstPerforationValue(cachedRefDensity);
1355 }
1356 }
1357
1358
1359
1360
1361
1362 template<typename TypeTag>
1363 void
1366 const WellStateType& well_state,
1367 DeferredLogger& deferred_logger)
1368 {
1369 const auto props = computePropertiesForWellConnectionPressures
1370 (simulator, well_state);
1371
1372 computeWellConnectionDensitesPressures(simulator, well_state,
1373 props, deferred_logger);
1374 }
1375
1376
1377
1378
1379
1380 template<typename TypeTag>
1381 void
1383 solveEqAndUpdateWellState(const Simulator& simulator,
1384 WellStateType& well_state,
1385 DeferredLogger& deferred_logger)
1386 {
1387 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1388
1389 // We assemble the well equations, then we check the convergence,
1390 // which is why we do not put the assembleWellEq here.
1391 BVectorWell dx_well(1);
1392 dx_well[0].resize(this->primary_variables_.numWellEq());
1393 this->linSys_.solve( dx_well);
1394
1395 updateWellState(simulator, dx_well, well_state, deferred_logger);
1396 }
1397
1398
1399
1400
1401
1402 template<typename TypeTag>
1403 void
1406 const WellStateType& well_state,
1407 DeferredLogger& deferred_logger)
1408 {
1409 updatePrimaryVariables(simulator, well_state, deferred_logger);
1410 computeWellConnectionPressures(simulator, well_state, deferred_logger);
1411 this->computeAccumWell();
1412 }
1413
1414
1415
1416 template<typename TypeTag>
1417 void
1419 apply(const BVector& x, BVector& Ax) const
1420 {
1421 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1422
1423 if (this->param_.matrix_add_well_contributions_)
1424 {
1425 // Contributions are already in the matrix itself
1426 return;
1427 }
1428
1429 this->linSys_.apply(x, Ax);
1430 }
1431
1432
1433
1434
1435 template<typename TypeTag>
1436 void
1438 apply(BVector& r) const
1439 {
1440 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1441
1442 this->linSys_.apply(r);
1443 }
1444
1445
1446
1447
1448 template<typename TypeTag>
1449 void
1452 const BVector& x,
1453 WellStateType& well_state,
1454 DeferredLogger& deferred_logger)
1455 {
1456 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1457
1458 BVectorWell xw(1);
1459 xw[0].resize(this->primary_variables_.numWellEq());
1460
1461 this->linSys_.recoverSolutionWell(x, xw);
1462 updateWellState(simulator, xw, well_state, deferred_logger);
1463 }
1464
1465
1466
1467
1468 template<typename TypeTag>
1469 void
1471 computeWellRatesWithBhp(const Simulator& simulator,
1472 const Scalar& bhp,
1473 std::vector<Scalar>& well_flux,
1474 DeferredLogger& deferred_logger) const
1475 {
1476 OPM_TIMEFUNCTION();
1477 const int np = this->number_of_phases_;
1478 well_flux.resize(np, 0.0);
1479
1480 const bool allow_cf = this->getAllowCrossFlow();
1481
1482 for (int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
1483 const int cell_idx = this->well_cells_[perf];
1484 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
1485 // flux for each perforation
1486 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.);
1487 getMobility(simulator, perf, mob, deferred_logger);
1488 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(intQuants, cell_idx);
1489 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1490 const std::vector<Scalar> Tw = this->wellIndex(perf, intQuants, trans_mult, wellstate_nupcol);
1491
1492 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.);
1493 PerforationRates<Scalar> perf_rates;
1494 computePerfRate(intQuants, mob, bhp, Tw, perf, allow_cf,
1495 cq_s, perf_rates, deferred_logger);
1496
1497 for(int p = 0; p < np; ++p) {
1498 well_flux[FluidSystem::activeCompToActivePhaseIdx(p)] += cq_s[p];
1499 }
1500
1501 // the solvent contribution is added to the gas potentials
1502 if constexpr (has_solvent) {
1503 assert(FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx));
1504 // TODO: should we use compIdx here?
1505 const int gas_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
1506 well_flux[gas_pos] += cq_s[Indices::contiSolventEqIdx];
1507 }
1508 }
1509 this->parallel_well_info_.communication().sum(well_flux.data(), well_flux.size());
1510 }
1511
1512
1513
1514 template<typename TypeTag>
1515 void
1518 const Scalar& bhp,
1519 std::vector<Scalar>& well_flux,
1520 DeferredLogger& deferred_logger) const
1521 {
1522 // creating a copy of the well itself, to avoid messing up the explicit information
1523 // during this copy, the only information not copied properly is the well controls
1524 StandardWell<TypeTag> well_copy(*this);
1525 well_copy.resetDampening();
1526
1527 // iterate to get a more accurate well density
1528 // create a copy of the well_state to use. If the operability checking is sucessful, we use this one
1529 // to replace the original one
1530 WellStateType well_state_copy = simulator.problem().wellModel().wellState();
1531 const auto& group_state = simulator.problem().wellModel().groupState();
1532
1533 // Get the current controls.
1534 const auto& summary_state = simulator.vanguard().summaryState();
1535 auto inj_controls = well_copy.well_ecl_.isInjector()
1536 ? well_copy.well_ecl_.injectionControls(summary_state)
1537 : Well::InjectionControls(0);
1538 auto prod_controls = well_copy.well_ecl_.isProducer()
1539 ? well_copy.well_ecl_.productionControls(summary_state) :
1540 Well::ProductionControls(0);
1541
1542 // Set current control to bhp, and bhp value in state, modify bhp limit in control object.
1543 auto& ws = well_state_copy.well(this->index_of_well_);
1544 if (well_copy.well_ecl_.isInjector()) {
1545 inj_controls.bhp_limit = bhp;
1546 ws.injection_cmode = Well::InjectorCMode::BHP;
1547 } else {
1548 prod_controls.bhp_limit = bhp;
1549 ws.production_cmode = Well::ProducerCMode::BHP;
1550 }
1551 ws.bhp = bhp;
1552
1553 // initialized the well rates with the potentials i.e. the well rates based on bhp
1554 const int np = this->number_of_phases_;
1555 const Scalar sign = this->well_ecl_.isInjector() ? 1.0 : -1.0;
1556 for (int phase = 0; phase < np; ++phase){
1557 well_state_copy.wellRates(this->index_of_well_)[phase]
1558 = sign * ws.well_potentials[phase];
1559 }
1560 well_copy.updatePrimaryVariables(simulator, well_state_copy, deferred_logger);
1561 well_copy.computeAccumWell();
1562
1563 const double dt = simulator.timeStepSize();
1564 const bool converged = well_copy.iterateWellEqWithControl(simulator, dt, inj_controls, prod_controls, well_state_copy, group_state, deferred_logger);
1565 if (!converged) {
1566 const std::string msg = " well " + name() + " did not get converged during well potential calculations "
1567 " potentials are computed based on unconverged solution";
1568 deferred_logger.debug(msg);
1569 }
1570 well_copy.updatePrimaryVariables(simulator, well_state_copy, deferred_logger);
1571 well_copy.computeWellConnectionPressures(simulator, well_state_copy, deferred_logger);
1572 well_copy.computeWellRatesWithBhp(simulator, bhp, well_flux, deferred_logger);
1573 }
1574
1575
1576
1577
1578 template<typename TypeTag>
1579 std::vector<typename StandardWell<TypeTag>::Scalar>
1582 DeferredLogger& deferred_logger,
1583 const WellStateType& well_state) const
1584 {
1585 std::vector<Scalar> potentials(this->number_of_phases_, 0.0);
1586 const auto& summary_state = simulator.vanguard().summaryState();
1587
1588 const auto& well = this->well_ecl_;
1589 if (well.isInjector()){
1590 const auto& controls = this->well_ecl_.injectionControls(summary_state);
1591 auto bhp_at_thp_limit = computeBhpAtThpLimitInj(simulator, summary_state, deferred_logger);
1592 if (bhp_at_thp_limit) {
1593 const Scalar bhp = std::min(*bhp_at_thp_limit,
1594 static_cast<Scalar>(controls.bhp_limit));
1595 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1596 } else {
1597 deferred_logger.warning("FAILURE_GETTING_CONVERGED_POTENTIAL",
1598 "Failed in getting converged thp based potential calculation for well "
1599 + name() + ". Instead the bhp based value is used");
1600 const Scalar bhp = controls.bhp_limit;
1601 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1602 }
1603 } else {
1604 computeWellRatesWithThpAlqProd(
1605 simulator, summary_state,
1606 deferred_logger, potentials, this->getALQ(well_state)
1607 );
1608 }
1609
1610 return potentials;
1611 }
1612
1613 template<typename TypeTag>
1614 bool
1617 const WellStateType& well_state,
1618 std::vector<Scalar>& well_potentials,
1619 DeferredLogger& deferred_logger) const
1620 {
1621 // Create a copy of the well.
1622 // TODO: check if we can avoid taking multiple copies. Call from updateWellPotentials
1623 // is allready a copy, but not from other calls.
1624 StandardWell<TypeTag> well_copy(*this);
1625
1626 // store a copy of the well state, we don't want to update the real well state
1627 WellStateType well_state_copy = well_state;
1628 const auto& group_state = simulator.problem().wellModel().groupState();
1629 auto& ws = well_state_copy.well(this->index_of_well_);
1630
1631 // get current controls
1632 const auto& summary_state = simulator.vanguard().summaryState();
1633 auto inj_controls = well_copy.well_ecl_.isInjector()
1634 ? well_copy.well_ecl_.injectionControls(summary_state)
1635 : Well::InjectionControls(0);
1636 auto prod_controls = well_copy.well_ecl_.isProducer()
1637 ? well_copy.well_ecl_.productionControls(summary_state) :
1638 Well::ProductionControls(0);
1639
1640 // prepare/modify well state and control
1641 well_copy.onlyKeepBHPandTHPcontrols(summary_state, well_state_copy, inj_controls, prod_controls);
1642
1643 // update connection pressures relative to updated bhp to get better estimate of connection dp
1644 const int num_perf = ws.perf_data.size();
1645 for (int perf = 0; perf < num_perf; ++perf) {
1646 ws.perf_data.pressure[perf] = ws.bhp + well_copy.connections_.pressure_diff(perf);
1647 }
1648 // initialize rates from previous potentials
1649 const int np = this->number_of_phases_;
1650 bool trivial = true;
1651 for (int phase = 0; phase < np; ++phase){
1652 trivial = trivial && (ws.well_potentials[phase] == 0.0) ;
1653 }
1654 if (!trivial) {
1655 const Scalar sign = well_copy.well_ecl_.isInjector() ? 1.0 : -1.0;
1656 for (int phase = 0; phase < np; ++phase) {
1657 ws.surface_rates[phase] = sign * ws.well_potentials[phase];
1658 }
1659 }
1660
1661 well_copy.calculateExplicitQuantities(simulator, well_state_copy, deferred_logger);
1662 const double dt = simulator.timeStepSize();
1663 // iterate to get a solution at the given bhp.
1664 bool converged = false;
1665 if (this->well_ecl_.isProducer()) {
1666 converged = well_copy.solveWellWithOperabilityCheck(simulator, dt, inj_controls, prod_controls, well_state_copy, group_state, deferred_logger);
1667 } else {
1668 converged = well_copy.iterateWellEqWithSwitching(simulator, dt, inj_controls, prod_controls, well_state_copy, group_state, deferred_logger);
1669 }
1670
1671 // fetch potentials (sign is updated on the outside).
1672 well_potentials.clear();
1673 well_potentials.resize(np, 0.0);
1674 for (int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1675 if (has_solvent && comp_idx == Indices::contiSolventEqIdx) continue; // we do not store the solvent in the well_potentials
1676 const EvalWell rate = well_copy.primary_variables_.getQs(comp_idx);
1677 well_potentials[FluidSystem::activeCompToActivePhaseIdx(comp_idx)] = rate.value();
1678 }
1679
1680 // the solvent contribution is added to the gas potentials
1681 if constexpr (has_solvent) {
1682 assert(FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx));
1683 // TODO: should we use compIdx here?
1684 const int gas_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
1685 const EvalWell rate = well_copy.primary_variables_.getQs(Indices::contiSolventEqIdx);
1686 well_potentials[gas_pos] += rate.value();
1687 }
1688 return converged;
1689 }
1690
1691
1692 template<typename TypeTag>
1696 const SummaryState &summary_state,
1697 DeferredLogger& deferred_logger,
1698 std::vector<Scalar>& potentials,
1699 Scalar alq) const
1700 {
1701 Scalar bhp;
1702 auto bhp_at_thp_limit = computeBhpAtThpLimitProdWithAlq(
1703 simulator, summary_state, alq, deferred_logger, /*iterate_if_no_solution */ true);
1704 if (bhp_at_thp_limit) {
1705 const auto& controls = this->well_ecl_.productionControls(summary_state);
1706 bhp = std::max(*bhp_at_thp_limit,
1707 static_cast<Scalar>(controls.bhp_limit));
1708 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1709 }
1710 else {
1711 deferred_logger.warning("FAILURE_GETTING_CONVERGED_POTENTIAL",
1712 "Failed in getting converged thp based potential calculation for well "
1713 + name() + ". Instead the bhp based value is used");
1714 const auto& controls = this->well_ecl_.productionControls(summary_state);
1715 bhp = controls.bhp_limit;
1716 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1717 }
1718 return bhp;
1719 }
1720
1721 template<typename TypeTag>
1722 void
1725 const SummaryState& summary_state,
1726 DeferredLogger& deferred_logger,
1727 std::vector<Scalar>& potentials,
1728 Scalar alq) const
1729 {
1730 /*double bhp =*/
1731 computeWellRatesAndBhpWithThpAlqProd(simulator,
1732 summary_state,
1733 deferred_logger,
1734 potentials,
1735 alq);
1736 }
1737
1738 template<typename TypeTag>
1739 void
1741 computeWellPotentials(const Simulator& simulator,
1742 const WellStateType& well_state,
1743 std::vector<Scalar>& well_potentials,
1744 DeferredLogger& deferred_logger) // const
1745 {
1746 const auto [compute_potential, bhp_controlled_well] =
1748
1749 if (!compute_potential) {
1750 return;
1751 }
1752
1753 bool converged_implicit = false;
1754 // for newly opened wells we dont compute the potentials implicit
1755 // group controlled wells with defaulted guiderates will have zero targets as
1756 // the potentials are used to compute the well fractions.
1757 if (this->param_.local_well_solver_control_switching_ && !(this->changed_to_open_this_step_ && this->wellUnderZeroRateTarget(simulator, well_state, deferred_logger))) {
1758 converged_implicit = computeWellPotentialsImplicit(simulator, well_state, well_potentials, deferred_logger);
1759 }
1760 if (!converged_implicit) {
1761 // does the well have a THP related constraint?
1762 const auto& summaryState = simulator.vanguard().summaryState();
1763 if (!Base::wellHasTHPConstraints(summaryState) || bhp_controlled_well) {
1764 // get the bhp value based on the bhp constraints
1765 Scalar bhp = WellBhpThpCalculator(*this).mostStrictBhpFromBhpLimits(summaryState);
1766
1767 // In some very special cases the bhp pressure target are
1768 // temporary violated. This may lead to too small or negative potentials
1769 // that could lead to premature shutting of wells.
1770 // As a remedy the bhp that gives the largest potential is used.
1771 // For converged cases, ws.bhp <=bhp for injectors and ws.bhp >= bhp,
1772 // and the potentials will be computed using the limit as expected.
1773 const auto& ws = well_state.well(this->index_of_well_);
1774 if (this->isInjector())
1775 bhp = std::max(ws.bhp, bhp);
1776 else
1777 bhp = std::min(ws.bhp, bhp);
1778
1779 assert(std::abs(bhp) != std::numeric_limits<Scalar>::max());
1780 computeWellRatesWithBhpIterations(simulator, bhp, well_potentials, deferred_logger);
1781 } else {
1782 // the well has a THP related constraint
1783 well_potentials = computeWellPotentialWithTHP(simulator, deferred_logger, well_state);
1784 }
1785 }
1786
1787 this->checkNegativeWellPotentials(well_potentials,
1788 this->param_.check_well_operability_,
1789 deferred_logger);
1790 }
1791
1792
1793
1794
1795
1796
1797
1798 template<typename TypeTag>
1801 connectionDensity([[maybe_unused]] const int globalConnIdx,
1802 const int openConnIdx) const
1803 {
1804 return (openConnIdx < 0)
1805 ? 0.0
1806 : this->connections_.rho(openConnIdx);
1807 }
1808
1809
1810
1811
1812
1813 template<typename TypeTag>
1814 void
1816 updatePrimaryVariables(const Simulator& simulator,
1817 const WellStateType& well_state,
1818 DeferredLogger& deferred_logger)
1819 {
1820 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1821
1822 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
1823 this->primary_variables_.update(well_state, stop_or_zero_rate_target, deferred_logger);
1824
1825 // other primary variables related to polymer injection
1826 if constexpr (Base::has_polymermw) {
1827 this->primary_variables_.updatePolyMW(well_state);
1828 }
1829
1830 this->primary_variables_.checkFinite(deferred_logger);
1831 }
1832
1833
1834
1835
1836 template<typename TypeTag>
1839 getRefDensity() const
1840 {
1841 return cachedRefDensity;
1842 }
1843
1844
1845
1846
1847 template<typename TypeTag>
1848 void
1851 const int perf,
1852 std::vector<EvalWell>& mob,
1853 DeferredLogger& deferred_logger) const
1854 {
1855 const int cell_idx = this->well_cells_[perf];
1856 const auto& int_quant = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
1857 const EvalWell polymer_concentration = this->extendEval(int_quant.polymerConcentration());
1858
1859 // TODO: not sure should based on the well type or injecting/producing peforations
1860 // it can be different for crossflow
1861 if (this->isInjector()) {
1862 // assume fully mixing within injecting wellbore
1863 const auto& visc_mult_table = PolymerModule::plyviscViscosityMultiplierTable(int_quant.pvtRegionIndex());
1864 const unsigned waterCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
1865 mob[waterCompIdx] /= (this->extendEval(int_quant.waterViscosityCorrection()) * visc_mult_table.eval(polymer_concentration, /*extrapolate=*/true) );
1866 }
1867
1868 if (PolymerModule::hasPlyshlog()) {
1869 // we do not calculate the shear effects for injection wells when they do not
1870 // inject polymer.
1871 if (this->isInjector() && this->wpolymer() == 0.) {
1872 return;
1873 }
1874 // compute the well water velocity with out shear effects.
1875 // TODO: do we need to turn on crossflow here?
1876 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
1877 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
1878
1879 std::vector<EvalWell> cq_s(this->num_conservation_quantities_, {this->primary_variables_.numWellEq() + Indices::numEq, 0.});
1880 PerforationRates<Scalar> perf_rates;
1881 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quant, cell_idx);
1882 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1883 const std::vector<Scalar> Tw = this->wellIndex(perf, int_quant, trans_mult, wellstate_nupcol);
1884 computePerfRate(int_quant, mob, bhp, Tw, perf, allow_cf, cq_s,
1885 perf_rates, deferred_logger);
1886 // TODO: make area a member
1887 const Scalar area = 2 * M_PI * this->perf_rep_radius_[perf] * this->perf_length_[perf];
1888 const auto& material_law_manager = simulator.problem().materialLawManager();
1889 const auto& scaled_drainage_info =
1890 material_law_manager->oilWaterScaledEpsInfoDrainage(cell_idx);
1891 const Scalar swcr = scaled_drainage_info.Swcr;
1892 const EvalWell poro = this->extendEval(int_quant.porosity());
1893 const EvalWell sw = this->extendEval(int_quant.fluidState().saturation(FluidSystem::waterPhaseIdx));
1894 // guard against zero porosity and no water
1895 const EvalWell denom = max( (area * poro * (sw - swcr)), 1e-12);
1896 const unsigned waterCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
1897 EvalWell water_velocity = cq_s[waterCompIdx] / denom * this->extendEval(int_quant.fluidState().invB(FluidSystem::waterPhaseIdx));
1898
1899 if (PolymerModule::hasShrate()) {
1900 // the equation for the water velocity conversion for the wells and reservoir are from different version
1901 // of implementation. It can be changed to be more consistent when possible.
1902 water_velocity *= PolymerModule::shrate( int_quant.pvtRegionIndex() ) / this->bore_diameters_[perf];
1903 }
1904 const EvalWell shear_factor = PolymerModule::computeShearFactor(polymer_concentration,
1905 int_quant.pvtRegionIndex(),
1906 water_velocity);
1907 // modify the mobility with the shear factor.
1908 mob[waterCompIdx] /= shear_factor;
1909 }
1910 }
1911
1912 template<typename TypeTag>
1913 void
1915 {
1916 this->linSys_.extract(jacobian);
1917 }
1918
1919
1920 template <typename TypeTag>
1921 void
1923 const BVector& weights,
1924 const int pressureVarIndex,
1925 const bool use_well_weights,
1926 const WellStateType& well_state) const
1927 {
1928 this->linSys_.extractCPRPressureMatrix(jacobian,
1929 weights,
1930 pressureVarIndex,
1931 use_well_weights,
1932 *this,
1933 Bhp,
1934 well_state);
1935 }
1936
1937
1938
1939 template<typename TypeTag>
1942 pskinwater(const Scalar throughput,
1943 const EvalWell& water_velocity,
1944 DeferredLogger& deferred_logger) const
1945 {
1946 if constexpr (Base::has_polymermw) {
1947 const int water_table_id = this->polymerWaterTable_();
1948 if (water_table_id <= 0) {
1949 OPM_DEFLOG_THROW(std::runtime_error,
1950 fmt::format("Unused SKPRWAT table id used for well {}", name()),
1951 deferred_logger);
1952 }
1953 const auto& water_table_func = PolymerModule::getSkprwatTable(water_table_id);
1954 const EvalWell throughput_eval(this->primary_variables_.numWellEq() + Indices::numEq, throughput);
1955 // the skin pressure when injecting water, which also means the polymer concentration is zero
1956 EvalWell pskin_water(this->primary_variables_.numWellEq() + Indices::numEq, 0.0);
1957 pskin_water = water_table_func.eval(throughput_eval, water_velocity);
1958 return pskin_water;
1959 } else {
1960 OPM_DEFLOG_THROW(std::runtime_error,
1961 fmt::format("Polymermw is not activated, while injecting "
1962 "skin pressure is requested for well {}", name()),
1963 deferred_logger);
1964 }
1965 }
1966
1967
1968
1969
1970
1971 template<typename TypeTag>
1974 pskin(const Scalar throughput,
1975 const EvalWell& water_velocity,
1976 const EvalWell& poly_inj_conc,
1977 DeferredLogger& deferred_logger) const
1978 {
1979 if constexpr (Base::has_polymermw) {
1980 const Scalar sign = water_velocity >= 0. ? 1.0 : -1.0;
1981 const EvalWell water_velocity_abs = abs(water_velocity);
1982 if (poly_inj_conc == 0.) {
1983 return sign * pskinwater(throughput, water_velocity_abs, deferred_logger);
1984 }
1985 const int polymer_table_id = this->polymerTable_();
1986 if (polymer_table_id <= 0) {
1987 OPM_DEFLOG_THROW(std::runtime_error,
1988 fmt::format("Unavailable SKPRPOLY table id used for well {}", name()),
1989 deferred_logger);
1990 }
1991 const auto& skprpolytable = PolymerModule::getSkprpolyTable(polymer_table_id);
1992 const Scalar reference_concentration = skprpolytable.refConcentration;
1993 const EvalWell throughput_eval(this->primary_variables_.numWellEq() + Indices::numEq, throughput);
1994 // the skin pressure when injecting water, which also means the polymer concentration is zero
1995 EvalWell pskin_poly(this->primary_variables_.numWellEq() + Indices::numEq, 0.0);
1996 pskin_poly = skprpolytable.table_func.eval(throughput_eval, water_velocity_abs);
1997 if (poly_inj_conc == reference_concentration) {
1998 return sign * pskin_poly;
1999 }
2000 // poly_inj_conc != reference concentration of the table, then some interpolation will be required
2001 const EvalWell pskin_water = pskinwater(throughput, water_velocity_abs, deferred_logger);
2002 const EvalWell pskin = pskin_water + (pskin_poly - pskin_water) / reference_concentration * poly_inj_conc;
2003 return sign * pskin;
2004 } else {
2005 OPM_DEFLOG_THROW(std::runtime_error,
2006 fmt::format("Polymermw is not activated, while injecting "
2007 "skin pressure is requested for well {}", name()),
2008 deferred_logger);
2009 }
2010 }
2011
2012
2013
2014
2015
2016 template<typename TypeTag>
2019 wpolymermw(const Scalar throughput,
2020 const EvalWell& water_velocity,
2021 DeferredLogger& deferred_logger) const
2022 {
2023 if constexpr (Base::has_polymermw) {
2024 const int table_id = this->polymerInjTable_();
2025 const auto& table_func = PolymerModule::getPlymwinjTable(table_id);
2026 const EvalWell throughput_eval(this->primary_variables_.numWellEq() + Indices::numEq, throughput);
2027 EvalWell molecular_weight(this->primary_variables_.numWellEq() + Indices::numEq, 0.);
2028 if (this->wpolymer() == 0.) { // not injecting polymer
2029 return molecular_weight;
2030 }
2031 molecular_weight = table_func.eval(throughput_eval, abs(water_velocity));
2032 return molecular_weight;
2033 } else {
2034 OPM_DEFLOG_THROW(std::runtime_error,
2035 fmt::format("Polymermw is not activated, while injecting "
2036 "polymer molecular weight is requested for well {}", name()),
2037 deferred_logger);
2038 }
2039 }
2040
2041
2042
2043
2044
2045 template<typename TypeTag>
2046 void
2048 updateWaterThroughput([[maybe_unused]] const double dt,
2049 WellStateType& well_state) const
2050 {
2051 if constexpr (Base::has_polymermw) {
2052 if (!this->isInjector()) {
2053 return;
2054 }
2055
2056 auto& perf_water_throughput = well_state.well(this->index_of_well_)
2057 .perf_data.water_throughput;
2058
2059 for (int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
2060 const Scalar perf_water_vel =
2061 this->primary_variables_.value(Bhp + 1 + perf);
2062
2063 // we do not consider the formation damage due to water
2064 // flowing from reservoir into wellbore
2065 if (perf_water_vel > Scalar{0}) {
2066 perf_water_throughput[perf] += perf_water_vel * dt;
2067 }
2068 }
2069 }
2070 }
2071
2072
2073
2074
2075
2076 template<typename TypeTag>
2077 void
2079 handleInjectivityRate(const Simulator& simulator,
2080 const int perf,
2081 std::vector<EvalWell>& cq_s) const
2082 {
2083 const int cell_idx = this->well_cells_[perf];
2084 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2085 const auto& fs = int_quants.fluidState();
2086 const EvalWell b_w = this->extendEval(fs.invB(FluidSystem::waterPhaseIdx));
2087 const Scalar area = M_PI * this->bore_diameters_[perf] * this->perf_length_[perf];
2088 const int wat_vel_index = Bhp + 1 + perf;
2089 const unsigned water_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
2090
2091 // water rate is update to use the form from water velocity, since water velocity is
2092 // a primary variable now
2093 cq_s[water_comp_idx] = area * this->primary_variables_.eval(wat_vel_index) * b_w;
2094 }
2095
2096
2097
2098
2099 template<typename TypeTag>
2100 void
2102 handleInjectivityEquations(const Simulator& simulator,
2103 const WellStateType& well_state,
2104 const int perf,
2105 const EvalWell& water_flux_s,
2106 DeferredLogger& deferred_logger)
2107 {
2108 const int cell_idx = this->well_cells_[perf];
2109 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2110 const auto& fs = int_quants.fluidState();
2111 const EvalWell b_w = this->extendEval(fs.invB(FluidSystem::waterPhaseIdx));
2112 const EvalWell water_flux_r = water_flux_s / b_w;
2113 const Scalar area = M_PI * this->bore_diameters_[perf] * this->perf_length_[perf];
2114 const EvalWell water_velocity = water_flux_r / area;
2115 const int wat_vel_index = Bhp + 1 + perf;
2116
2117 // equation for the water velocity
2118 const EvalWell eq_wat_vel = this->primary_variables_.eval(wat_vel_index) - water_velocity;
2119
2120 const auto& ws = well_state.well(this->index_of_well_);
2121 const auto& perf_data = ws.perf_data;
2122 const auto& perf_water_throughput = perf_data.water_throughput;
2123 const Scalar throughput = perf_water_throughput[perf];
2124 const int pskin_index = Bhp + 1 + this->number_of_local_perforations_ + perf;
2125
2126 EvalWell poly_conc(this->primary_variables_.numWellEq() + Indices::numEq, 0.0);
2127 poly_conc.setValue(this->wpolymer());
2128
2129 // equation for the skin pressure
2130 const EvalWell eq_pskin = this->primary_variables_.eval(pskin_index)
2131 - pskin(throughput, this->primary_variables_.eval(wat_vel_index), poly_conc, deferred_logger);
2132
2134 assembleInjectivityEq(eq_pskin,
2135 eq_wat_vel,
2136 pskin_index,
2137 wat_vel_index,
2138 perf,
2139 this->primary_variables_.numWellEq(),
2140 this->linSys_);
2141 }
2142
2143
2144
2145
2146
2147 template<typename TypeTag>
2148 void
2150 checkConvergenceExtraEqs(const std::vector<Scalar>& res,
2151 ConvergenceReport& report) const
2152 {
2153 // if different types of extra equations are involved, this function needs to be refactored further
2154
2155 // checking the convergence of the extra equations related to polymer injectivity
2156 if constexpr (Base::has_polymermw) {
2157 WellConvergence(*this).
2158 checkConvergencePolyMW(res, Bhp, this->param_.max_residual_allowed_, report);
2159 }
2160 }
2161
2162
2163
2164
2165
2166 template<typename TypeTag>
2167 void
2169 updateConnectionRatePolyMW(const EvalWell& cq_s_poly,
2170 const IntensiveQuantities& int_quants,
2171 const WellStateType& well_state,
2172 const int perf,
2173 std::vector<RateVector>& connectionRates,
2174 DeferredLogger& deferred_logger) const
2175 {
2176 // the source term related to transport of molecular weight
2177 EvalWell cq_s_polymw = cq_s_poly;
2178 if (this->isInjector()) {
2179 const int wat_vel_index = Bhp + 1 + perf;
2180 const EvalWell water_velocity = this->primary_variables_.eval(wat_vel_index);
2181 if (water_velocity > 0.) { // injecting
2182 const auto& ws = well_state.well(this->index_of_well_);
2183 const auto& perf_water_throughput = ws.perf_data.water_throughput;
2184 const Scalar throughput = perf_water_throughput[perf];
2185 const EvalWell molecular_weight = wpolymermw(throughput, water_velocity, deferred_logger);
2186 cq_s_polymw *= molecular_weight;
2187 } else {
2188 // we do not consider the molecular weight from the polymer
2189 // going-back to the wellbore through injector
2190 cq_s_polymw *= 0.;
2191 }
2192 } else if (this->isProducer()) {
2193 if (cq_s_polymw < 0.) {
2194 cq_s_polymw *= this->extendEval(int_quants.polymerMoleWeight() );
2195 } else {
2196 // we do not consider the molecular weight from the polymer
2197 // re-injecting back through producer
2198 cq_s_polymw *= 0.;
2199 }
2200 }
2201 connectionRates[perf][Indices::contiPolymerMWEqIdx] = Base::restrictEval(cq_s_polymw);
2202 }
2203
2204
2205
2206
2207
2208 template<typename TypeTag>
2209 std::optional<typename StandardWell<TypeTag>::Scalar>
2212 const Simulator& simulator,
2213 const SummaryState& summary_state,
2214 DeferredLogger& deferred_logger) const
2215 {
2216 return computeBhpAtThpLimitProdWithAlq(simulator,
2217 summary_state,
2218 this->getALQ(well_state),
2219 deferred_logger,
2220 /*iterate_if_no_solution */ true);
2221 }
2222
2223 template<typename TypeTag>
2224 std::optional<typename StandardWell<TypeTag>::Scalar>
2227 const SummaryState& summary_state,
2228 const Scalar alq_value,
2229 DeferredLogger& deferred_logger,
2230 bool iterate_if_no_solution) const
2231 {
2232 OPM_TIMEFUNCTION();
2233 // Make the frates() function.
2234 auto frates = [this, &simulator, &deferred_logger](const Scalar bhp) {
2235 // Not solving the well equations here, which means we are
2236 // calculating at the current Fg/Fw values of the
2237 // well. This does not matter unless the well is
2238 // crossflowing, and then it is likely still a good
2239 // approximation.
2240 std::vector<Scalar> rates(3);
2241 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2242 this->adaptRatesForVFP(rates);
2243 return rates;
2244 };
2245
2246 Scalar max_pressure = 0.0;
2247 for (int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
2248 const int cell_idx = this->well_cells_[perf];
2249 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2250 const auto& fs = int_quants.fluidState();
2251 Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2252 max_pressure = std::max(max_pressure, pressure_cell);
2253 }
2254 const auto& comm = this->parallel_well_info_.communication();
2255 if (comm.size() > 1) {
2256 max_pressure = comm.max(max_pressure);
2257 }
2258 auto bhpAtLimit = WellBhpThpCalculator(*this).computeBhpAtThpLimitProd(frates,
2259 summary_state,
2260 max_pressure,
2261 this->getRefDensity(),
2262 alq_value,
2263 this->getTHPConstraint(summary_state),
2264 deferred_logger);
2265
2266 if (bhpAtLimit) {
2267 auto v = frates(*bhpAtLimit);
2268 if (std::all_of(v.cbegin(), v.cend(), [](Scalar i){ return i <= 0; }) ) {
2269 return bhpAtLimit;
2270 }
2271 }
2272
2273 if (!iterate_if_no_solution)
2274 return std::nullopt;
2275
2276 auto fratesIter = [this, &simulator, &deferred_logger](const Scalar bhp) {
2277 // Solver the well iterations to see if we are
2278 // able to get a solution with an update
2279 // solution
2280 std::vector<Scalar> rates(3);
2281 computeWellRatesWithBhpIterations(simulator, bhp, rates, deferred_logger);
2282 this->adaptRatesForVFP(rates);
2283 return rates;
2284 };
2285
2286 bhpAtLimit = WellBhpThpCalculator(*this).computeBhpAtThpLimitProd(fratesIter,
2287 summary_state,
2288 max_pressure,
2289 this->getRefDensity(),
2290 alq_value,
2291 this->getTHPConstraint(summary_state),
2292 deferred_logger);
2293
2294
2295 if (bhpAtLimit) {
2296 // should we use fratesIter here since fratesIter is used in computeBhpAtThpLimitProd above?
2297 auto v = frates(*bhpAtLimit);
2298 if (std::all_of(v.cbegin(), v.cend(), [](Scalar i){ return i <= 0; }) ) {
2299 return bhpAtLimit;
2300 }
2301 }
2302
2303 // we still don't get a valied solution.
2304 return std::nullopt;
2305 }
2306
2307
2308
2309 template<typename TypeTag>
2310 std::optional<typename StandardWell<TypeTag>::Scalar>
2312 computeBhpAtThpLimitInj(const Simulator& simulator,
2313 const SummaryState& summary_state,
2314 DeferredLogger& deferred_logger) const
2315 {
2316 // Make the frates() function.
2317 auto frates = [this, &simulator, &deferred_logger](const Scalar bhp) {
2318 // Not solving the well equations here, which means we are
2319 // calculating at the current Fg/Fw values of the
2320 // well. This does not matter unless the well is
2321 // crossflowing, and then it is likely still a good
2322 // approximation.
2323 std::vector<Scalar> rates(3);
2324 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2325 return rates;
2326 };
2327
2328 return WellBhpThpCalculator(*this).computeBhpAtThpLimitInj(frates,
2329 summary_state,
2330 this->getRefDensity(),
2331 1e-6,
2332 50,
2333 true,
2334 deferred_logger);
2335 }
2336
2337
2338
2339
2340
2341 template<typename TypeTag>
2342 bool
2344 iterateWellEqWithControl(const Simulator& simulator,
2345 const double dt,
2346 const Well::InjectionControls& inj_controls,
2347 const Well::ProductionControls& prod_controls,
2348 WellStateType& well_state,
2349 const GroupState<Scalar>& group_state,
2350 DeferredLogger& deferred_logger)
2351 {
2352 updatePrimaryVariables(simulator, well_state, deferred_logger);
2353
2354 const int max_iter = this->param_.max_inner_iter_wells_;
2355 int it = 0;
2356 bool converged;
2357 bool relax_convergence = false;
2358 this->regularize_ = false;
2359 do {
2360 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
2361
2362 if (it > this->param_.strict_inner_iter_wells_) {
2363 relax_convergence = true;
2364 this->regularize_ = true;
2365 }
2366
2367 auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
2368
2369 converged = report.converged();
2370 if (converged) {
2371 break;
2372 }
2373
2374 ++it;
2375 solveEqAndUpdateWellState(simulator, well_state, deferred_logger);
2376
2377 // TODO: when this function is used for well testing purposes, will need to check the controls, so that we will obtain convergence
2378 // under the most restrictive control. Based on this converged results, we can check whether to re-open the well. Either we refactor
2379 // this function or we use different functions for the well testing purposes.
2380 // We don't allow for switching well controls while computing well potentials and testing wells
2381 // updateWellControl(simulator, well_state, deferred_logger);
2382 } while (it < max_iter);
2383
2384 return converged;
2385 }
2386
2387
2388 template<typename TypeTag>
2389 bool
2391 iterateWellEqWithSwitching(const Simulator& simulator,
2392 const double dt,
2393 const Well::InjectionControls& inj_controls,
2394 const Well::ProductionControls& prod_controls,
2395 WellStateType& well_state,
2396 const GroupState<Scalar>& group_state,
2397 DeferredLogger& deferred_logger,
2398 const bool fixed_control /*false*/,
2399 const bool fixed_status /*false*/)
2400 {
2401 updatePrimaryVariables(simulator, well_state, deferred_logger);
2402
2403 const int max_iter = this->param_.max_inner_iter_wells_;
2404 int it = 0;
2405 bool converged = false;
2406 bool relax_convergence = false;
2407 this->regularize_ = false;
2408 const auto& summary_state = simulator.vanguard().summaryState();
2409
2410 // Always take a few (more than one) iterations after a switch before allowing a new switch
2411 // The optimal number here is subject to further investigation, but it has been observerved
2412 // that unless this number is >1, we may get stuck in a cycle
2413 constexpr int min_its_after_switch = 4;
2414 // We also want to restrict the number of status switches to avoid oscillation between STOP<->OPEN
2415 const int max_status_switch = this->param_.max_well_status_switch_inner_iter_;
2416 int its_since_last_switch = min_its_after_switch;
2417 int switch_count= 0;
2418 // if we fail to solve eqs, we reset status/operability before leaving
2419 const auto well_status_orig = this->wellStatus_;
2420 const auto operability_orig = this->operability_status_;
2421 auto well_status_cur = well_status_orig;
2422 int status_switch_count = 0;
2423 // don't allow opening wells that has a stopped well status
2424 const bool allow_open = well_state.well(this->index_of_well_).status == WellStatus::OPEN;
2425 // don't allow switcing for wells under zero rate target or requested fixed status and control
2426 const bool allow_switching =
2427 !this->wellUnderZeroRateTarget(simulator, well_state, deferred_logger) &&
2428 (!fixed_control || !fixed_status) && allow_open;
2429
2430 bool changed = false;
2431 bool final_check = false;
2432 // well needs to be set operable or else solving/updating of re-opened wells is skipped
2433 this->operability_status_.resetOperability();
2434 this->operability_status_.solvable = true;
2435 do {
2436 its_since_last_switch++;
2437 if (allow_switching && its_since_last_switch >= min_its_after_switch && status_switch_count < max_status_switch){
2438 const Scalar wqTotal = this->primary_variables_.eval(WQTotal).value();
2439 changed = this->updateWellControlAndStatusLocalIteration(simulator, well_state, group_state,
2440 inj_controls, prod_controls, wqTotal,
2441 deferred_logger, fixed_control, fixed_status);
2442 if (changed){
2443 its_since_last_switch = 0;
2444 switch_count++;
2445 if (well_status_cur != this->wellStatus_) {
2446 well_status_cur = this->wellStatus_;
2447 status_switch_count++;
2448 }
2449 }
2450 if (!changed && final_check) {
2451 break;
2452 } else {
2453 final_check = false;
2454 }
2455 if (status_switch_count == max_status_switch) {
2456 this->wellStatus_ = well_status_orig;
2457 }
2458 }
2459
2460 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
2461
2462 if (it > this->param_.strict_inner_iter_wells_) {
2463 relax_convergence = true;
2464 this->regularize_ = true;
2465 }
2466
2467 auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
2468
2469 converged = report.converged();
2470 if (converged) {
2471 // if equations are sufficiently linear they might converge in less than min_its_after_switch
2472 // in this case, make sure all constraints are satisfied before returning
2473 if (switch_count > 0 && its_since_last_switch < min_its_after_switch) {
2474 final_check = true;
2475 its_since_last_switch = min_its_after_switch;
2476 } else {
2477 break;
2478 }
2479 }
2480
2481 ++it;
2482 solveEqAndUpdateWellState(simulator, well_state, deferred_logger);
2483
2484 } while (it < max_iter);
2485
2486 if (converged) {
2487 if (allow_switching){
2488 // update operability if status change
2489 const bool is_stopped = this->wellIsStopped();
2490 if (this->wellHasTHPConstraints(summary_state)){
2491 this->operability_status_.can_obtain_bhp_with_thp_limit = !is_stopped;
2492 this->operability_status_.obey_thp_limit_under_bhp_limit = !is_stopped;
2493 } else {
2494 this->operability_status_.operable_under_only_bhp_limit = !is_stopped;
2495 }
2496 }
2497 } else {
2498 this->wellStatus_ = well_status_orig;
2499 this->operability_status_ = operability_orig;
2500 const std::string message = fmt::format(" Well {} did not converge in {} inner iterations ("
2501 "{} switches, {} status changes).", this->name(), it, switch_count, status_switch_count);
2502 deferred_logger.debug(message);
2503 // add operability here as well ?
2504 }
2505 return converged;
2506 }
2507
2508 template<typename TypeTag>
2509 std::vector<typename StandardWell<TypeTag>::Scalar>
2511 computeCurrentWellRates(const Simulator& simulator,
2512 DeferredLogger& deferred_logger) const
2513 {
2514 // Calculate the rates that follow from the current primary variables.
2515 std::vector<Scalar> well_q_s(this->num_conservation_quantities_, 0.);
2516 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
2517 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
2518 for (int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
2519 const int cell_idx = this->well_cells_[perf];
2520 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2521 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.);
2522 getMobility(simulator, perf, mob, deferred_logger);
2523 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.);
2524 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(intQuants, cell_idx);
2525 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
2526 const std::vector<Scalar> Tw = this->wellIndex(perf, intQuants, trans_mult, wellstate_nupcol);
2527 PerforationRates<Scalar> perf_rates;
2528 computePerfRate(intQuants, mob, bhp.value(), Tw, perf, allow_cf,
2529 cq_s, perf_rates, deferred_logger);
2530 for (int comp = 0; comp < this->num_conservation_quantities_; ++comp) {
2531 well_q_s[comp] += cq_s[comp];
2532 }
2533 }
2534 const auto& comm = this->parallel_well_info_.communication();
2535 if (comm.size() > 1)
2536 {
2537 comm.sum(well_q_s.data(), well_q_s.size());
2538 }
2539 return well_q_s;
2540 }
2541
2542
2543
2544 template <typename TypeTag>
2545 std::vector<typename StandardWell<TypeTag>::Scalar>
2547 getPrimaryVars() const
2548 {
2549 const int num_pri_vars = this->primary_variables_.numWellEq();
2550 std::vector<Scalar> retval(num_pri_vars);
2551 for (int ii = 0; ii < num_pri_vars; ++ii) {
2552 retval[ii] = this->primary_variables_.value(ii);
2553 }
2554 return retval;
2555 }
2556
2557
2558
2559
2560
2561 template <typename TypeTag>
2562 int
2564 setPrimaryVars(typename std::vector<Scalar>::const_iterator it)
2565 {
2566 const int num_pri_vars = this->primary_variables_.numWellEq();
2567 for (int ii = 0; ii < num_pri_vars; ++ii) {
2568 this->primary_variables_.setValue(ii, it[ii]);
2569 }
2570 return num_pri_vars;
2571 }
2572
2573
2574 template <typename TypeTag>
2577 connectionRateEnergy(const std::vector<EvalWell>& cq_s,
2578 const IntensiveQuantities& intQuants,
2579 DeferredLogger& deferred_logger) const
2580 {
2581 auto fs = intQuants.fluidState();
2582 Eval result = 0;
2583 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
2584 if (!FluidSystem::phaseIsActive(phaseIdx)) {
2585 continue;
2586 }
2587
2588 // convert to reservoir conditions
2589 EvalWell cq_r_thermal(this->primary_variables_.numWellEq() + Indices::numEq, 0.);
2590 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
2591 const bool both_oil_gas = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
2592 if (!both_oil_gas || FluidSystem::waterPhaseIdx == phaseIdx) {
2593 cq_r_thermal = cq_s[activeCompIdx] / this->extendEval(fs.invB(phaseIdx));
2594 } else {
2595 // remove dissolved gas and vapporized oil
2596 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
2597 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
2598 // q_os = q_or * b_o + rv * q_gr * b_g
2599 // q_gs = q_gr * g_g + rs * q_or * b_o
2600 // q_gr = 1 / (b_g * d) * (q_gs - rs * q_os)
2601 // d = 1.0 - rs * rv
2602 const EvalWell d = this->extendEval(1.0 - fs.Rv() * fs.Rs());
2603 if (d <= 0.0) {
2604 deferred_logger.debug(
2605 fmt::format("Problematic d value {} obtained for well {}"
2606 " during calculateSinglePerf with rs {}"
2607 ", rv {}. Continue as if no dissolution (rs = 0) and"
2608 " vaporization (rv = 0) for this connection.",
2609 d, this->name(), fs.Rs(), fs.Rv()));
2610 cq_r_thermal = cq_s[activeCompIdx] / this->extendEval(fs.invB(phaseIdx));
2611 } else {
2612 if (FluidSystem::gasPhaseIdx == phaseIdx) {
2613 cq_r_thermal = (cq_s[gasCompIdx] -
2614 this->extendEval(fs.Rs()) * cq_s[oilCompIdx]) /
2615 (d * this->extendEval(fs.invB(phaseIdx)) );
2616 } else if (FluidSystem::oilPhaseIdx == phaseIdx) {
2617 // q_or = 1 / (b_o * d) * (q_os - rv * q_gs)
2618 cq_r_thermal = (cq_s[oilCompIdx] - this->extendEval(fs.Rv()) *
2619 cq_s[gasCompIdx]) /
2620 (d * this->extendEval(fs.invB(phaseIdx)) );
2621 }
2622 }
2623 }
2624
2625 // change temperature for injecting fluids
2626 if (this->isInjector() && !this->wellIsStopped() && cq_r_thermal > 0.0){
2627 // only handles single phase injection now
2628 assert(this->well_ecl_.injectorType() != InjectorType::MULTI);
2629 fs.setTemperature(this->well_ecl_.inj_temperature());
2630 typedef typename std::decay<decltype(fs)>::type::Scalar FsScalar;
2631 typename FluidSystem::template ParameterCache<FsScalar> paramCache;
2632 const unsigned pvtRegionIdx = intQuants.pvtRegionIndex();
2633 paramCache.setRegionIndex(pvtRegionIdx);
2634 paramCache.updatePhase(fs, phaseIdx);
2635
2636 const auto& rho = FluidSystem::density(fs, paramCache, phaseIdx);
2637 fs.setDensity(phaseIdx, rho);
2638 const auto& h = FluidSystem::enthalpy(fs, paramCache, phaseIdx);
2639 fs.setEnthalpy(phaseIdx, h);
2640 cq_r_thermal *= this->extendEval(fs.enthalpy(phaseIdx)) * this->extendEval(fs.density(phaseIdx));
2641 result += getValue(cq_r_thermal);
2642 } else if (cq_r_thermal > 0.0) {
2643 cq_r_thermal *= getValue(fs.enthalpy(phaseIdx)) * getValue(fs.density(phaseIdx));
2644 result += Base::restrictEval(cq_r_thermal);
2645 } else {
2646 // compute the thermal flux
2647 cq_r_thermal *= this->extendEval(fs.enthalpy(phaseIdx)) * this->extendEval(fs.density(phaseIdx));
2648 result += Base::restrictEval(cq_r_thermal);
2649 }
2650 }
2651
2652 return result * this->well_efficiency_factor_;
2653 }
2654} // namespace Opm
2655
2656#endif
#define OPM_DEFLOG_THROW(Exception, message, deferred_logger)
Definition: DeferredLoggingErrorHelpers.hpp:45
#define OPM_DEFLOG_PROBLEM(Exception, message, deferred_logger)
Definition: DeferredLoggingErrorHelpers.hpp:61
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
void warning(const std::string &tag, const std::string &message)
void debug(const std::string &tag, const std::string &message)
Definition: GroupState.hpp:41
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:198
Definition: RatioCalculator.hpp:38
Class handling assemble of the equation system for StandardWell.
Definition: StandardWellAssemble.hpp:43
Scalar pressure_diff(const unsigned perf) const
Returns pressure drop for a given perforation.
Definition: StandardWellConnections.hpp:100
StdWellConnections connections_
Connection level values.
Definition: StandardWellEval.hpp:101
PrimaryVariables primary_variables_
Primary variables for well.
Definition: StandardWellEval.hpp:95
Definition: StandardWell.hpp:60
void calculateExplicitQuantities(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1405
EvalWell wpolymermw(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2019
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:775
WellConnectionProps computePropertiesForWellConnectionPressures(const Simulator &simulator, const WellStateType &well_state) const
Definition: StandardWell_impl.hpp:1139
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:2211
void computeWellRatesWithThpAlqProd(const Simulator &ebos_simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1724
void addWellContributions(SparseMatrixAdapter &mat) const override
Definition: StandardWell_impl.hpp:1914
std::vector< Scalar > getPrimaryVars() const override
Definition: StandardWell_impl.hpp:2547
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:1922
void updateWaterMobilityWithPolymer(const Simulator &simulator, const int perf, std::vector< EvalWell > &mob_water, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1850
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:2511
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
void updatePrimaryVariablesNewton(const BVectorWell &dwells, const bool stop_or_zero_rate_target, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:752
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:82
void computeWellConnectionPressures(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1365
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:2169
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:2391
void computeWellRatesWithBhp(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1471
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:674
void computeWellRatesWithBhpIterations(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1517
void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:799
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:2226
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:1309
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: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:2102
virtual void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: StandardWell_impl.hpp:1419
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:1190
void checkConvergenceExtraEqs(const std::vector< Scalar > &res, ConvergenceReport &report) const
Definition: StandardWell_impl.hpp:2150
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:1695
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1128
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2312
bool allDrawDownWrongDirection(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1086
EvalWell pskin(const Scalar throughput, const EvalWell &water_velocity, const EvalWell &poly_inj_conc, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1974
static constexpr int numWellConservationEq
Definition: StandardWell.hpp:96
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: StandardWell_impl.hpp:2564
void checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1037
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:1741
bool computeWellPotentialsImplicit(const Simulator &ebos_simulator, const WellStateType &well_state, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1616
void updateWaterThroughput(const double dt, WellStateType &well_state) const override
Definition: StandardWell_impl.hpp:2048
void checkOperabilityUnderBHPLimit(const WellStateType &well_state, const Simulator &simulator, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:968
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1451
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:2344
EvalWell pskinwater(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1942
void solveEqAndUpdateWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1383
void handleInjectivityRate(const Simulator &simulator, const int perf, std::vector< EvalWell > &cq_s) const
Definition: StandardWell_impl.hpp:2079
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:1581
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1234
void updatePrimaryVariables(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1816
Scalar getRefDensity() const override
Definition: StandardWell_impl.hpp:1839
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: StandardWell_impl.hpp:1801
EvalWell getQs(const int compIdx) const
Returns scaled rate for a component.
Class for computing BHP limits.
Definition: WellBhpThpCalculator.hpp:41
Scalar calculateThpFromBhp(const std::vector< Scalar > &rates, const Scalar bhp, const Scalar rho, const std::optional< Scalar > &alq, const Scalar thp_limit, DeferredLogger &deferred_logger) const
Calculates THP from BHP.
std::optional< Scalar > computeBhpAtThpLimitProd(const std::function< std::vector< Scalar >(const Scalar)> &frates, const SummaryState &summary_state, const Scalar maxPerfPress, const Scalar rho, const Scalar alq_value, const Scalar thp_limit, DeferredLogger &deferred_logger) const
Compute BHP from THP limit for a producer.
Scalar mostStrictBhpFromBhpLimits(const SummaryState &summaryState) const
Obtain the most strict BHP from BHP limits.
std::optional< Scalar > computeBhpAtThpLimitInj(const std::function< std::vector< Scalar >(const Scalar)> &frates, const SummaryState &summary_state, const Scalar rho, const Scalar flo_rel_tol, const int max_iteration, const bool throwOnError, DeferredLogger &deferred_logger) const
Compute BHP from THP limit for an injector.
Definition: WellConvergence.hpp:38
const int num_conservation_quantities_
Definition: WellInterfaceGeneric.hpp:314
Well well_ecl_
Definition: WellInterfaceGeneric.hpp:304
void onlyKeepBHPandTHPcontrols(const SummaryState &summary_state, WellStateType &well_state, Well::InjectionControls &inj_controls, Well::ProductionControls &prod_controls) const
void resetDampening()
Definition: WellInterfaceGeneric.hpp:247
std::pair< bool, bool > computeWellPotentials(std::vector< Scalar > &well_potentials, const WellStateType &well_state)
Definition: WellInterfaceIndices.hpp:34
Definition: WellInterface.hpp:76
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:81
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:103
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:97
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1960
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:86
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:82
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:96
typename Base::ModelParameters ModelParameters
Definition: WellInterface.hpp:109
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:83
bool solveWellWithOperabilityCheck(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: WellInterface_impl.hpp:566
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:85
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:88
Definition: WellProdIndexCalculator.hpp:37
Scalar connectionProdIndStandard(const std::size_t connIdx, const Scalar connMobility) const
Definition: WellState.hpp:66
const SingleWellState< Scalar, IndexTraits > & well(std::size_t well_index) const
Definition: WellState.hpp:290
std::vector< Scalar > & wellRates(std::size_t well_index)
One rate per well and phase.
Definition: WellState.hpp:255
@ NONE
Definition: DeferredLogger.hpp:46
Definition: blackoilbioeffectsmodules.hh:43
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Static data associated with a well perforation.
Definition: PerforationData.hpp:30
Definition: PerforationData.hpp:41
Scalar dis_gas
Definition: PerforationData.hpp:42
Scalar vap_wat
Definition: PerforationData.hpp:45
Scalar vap_oil
Definition: PerforationData.hpp:44
Scalar dis_gas_in_water
Definition: PerforationData.hpp:43