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_, 0.0);
395 EvalWell water_flux_s{0.0};
396 EvalWell cq_s_zfrac_effective{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(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_, {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_, 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 const WellGroupHelperType& wgHelper,
1040 DeferredLogger& deferred_logger)
1041 {
1042 const auto& summaryState = simulator.vanguard().summaryState();
1043 const auto obtain_bhp = this->isProducer() ? computeBhpAtThpLimitProd(well_state, simulator, wgHelper, summaryState, deferred_logger)
1044 : computeBhpAtThpLimitInj(simulator, wgHelper, summaryState, deferred_logger);
1045
1046 if (obtain_bhp) {
1047 this->operability_status_.can_obtain_bhp_with_thp_limit = true;
1048
1049 const Scalar bhp_limit = WellBhpThpCalculator(*this).mostStrictBhpFromBhpLimits(summaryState);
1050 this->operability_status_.obey_bhp_limit_with_thp_limit = this->isProducer() ?
1051 *obtain_bhp >= bhp_limit : *obtain_bhp <= bhp_limit ;
1052
1053 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1054 if (this->isProducer() && *obtain_bhp < thp_limit) {
1055 const std::string msg = " obtained bhp " + std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1056 + " bars is SMALLER than thp limit "
1057 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1058 + " bars as a producer for well " + name();
1059 deferred_logger.debug(msg);
1060 }
1061 else if (this->isInjector() && *obtain_bhp > thp_limit) {
1062 const std::string msg = " obtained bhp " + std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1063 + " bars is LARGER than thp limit "
1064 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1065 + " bars as a injector for well " + name();
1066 deferred_logger.debug(msg);
1067 }
1068 } else {
1069 this->operability_status_.can_obtain_bhp_with_thp_limit = false;
1070 this->operability_status_.obey_bhp_limit_with_thp_limit = false;
1071 if (!this->wellIsStopped()) {
1072 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1073 deferred_logger.debug(" could not find bhp value at thp limit "
1074 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1075 + " bar for well " + name() + ", the well might need to be closed ");
1076 }
1077 }
1078 }
1079
1080
1081
1082
1083
1084 template<typename TypeTag>
1085 bool
1087 allDrawDownWrongDirection(const Simulator& simulator) const
1088 {
1089 bool all_drawdown_wrong_direction = true;
1090
1091 for (int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
1092 const int cell_idx = this->well_cells_[perf];
1093 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/0);
1094 const auto& fs = intQuants.fluidState();
1095
1096 const Scalar pressure = this->getPerfCellPressure(fs).value();
1097 const Scalar bhp = this->primary_variables_.eval(Bhp).value();
1098
1099 // Pressure drawdown (also used to determine direction of flow)
1100 const Scalar well_pressure = bhp + this->connections_.pressure_diff(perf);
1101 const Scalar drawdown = pressure - well_pressure;
1102
1103 // for now, if there is one perforation can produce/inject in the correct
1104 // direction, we consider this well can still produce/inject.
1105 // TODO: it can be more complicated than this to cause wrong-signed rates
1106 if ( (drawdown < 0. && this->isInjector()) ||
1107 (drawdown > 0. && this->isProducer()) ) {
1108 all_drawdown_wrong_direction = false;
1109 break;
1110 }
1111 }
1112
1113 const auto& comm = this->parallel_well_info_.communication();
1114 if (comm.size() > 1)
1115 {
1116 all_drawdown_wrong_direction =
1117 (comm.min(all_drawdown_wrong_direction ? 1 : 0) == 1);
1118 }
1119
1120 return all_drawdown_wrong_direction;
1121 }
1122
1123
1124
1125
1126 template<typename TypeTag>
1127 bool
1129 openCrossFlowAvoidSingularity(const Simulator& simulator) const
1130 {
1131 return !this->getAllowCrossFlow() && allDrawDownWrongDirection(simulator);
1132 }
1133
1134
1135
1136
1137 template<typename TypeTag>
1141 const WellStateType& well_state) const
1142 {
1143 auto prop_func = typename StdWellEval::StdWellConnections::PressurePropertyFunctions {
1144 // getTemperature
1145 [&model = simulator.model()](int cell_idx, int phase_idx)
1146 {
1147 return model.intensiveQuantities(cell_idx, /* time_idx = */ 0)
1148 .fluidState().temperature(phase_idx).value();
1149 },
1150
1151 // getSaltConcentration
1152 [&model = simulator.model()](int cell_idx)
1153 {
1154 return model.intensiveQuantities(cell_idx, /* time_idx = */ 0)
1155 .fluidState().saltConcentration().value();
1156 },
1157
1158 // getPvtRegionIdx
1159 [&model = simulator.model()](int cell_idx)
1160 {
1161 return model.intensiveQuantities(cell_idx, /* time_idx = */ 0)
1162 .fluidState().pvtRegionIndex();
1163 }
1164 };
1165
1166 if constexpr (Indices::enableSolvent) {
1167 prop_func.solventInverseFormationVolumeFactor =
1168 [&model = simulator.model()](int cell_idx)
1169 {
1170 return model.intensiveQuantities(cell_idx, /* time_idx = */ 0)
1171 .solventInverseFormationVolumeFactor().value();
1172 };
1173
1174 prop_func.solventRefDensity = [&model = simulator.model()](int cell_idx)
1175 {
1176 return model.intensiveQuantities(cell_idx, /* time_idx = */ 0)
1177 .solventRefDensity();
1178 };
1179 }
1180
1181 return this->connections_.computePropertiesForPressures(well_state, prop_func);
1182 }
1183
1184
1185
1186
1187
1188 template<typename TypeTag>
1191 getWellConvergence(const Simulator& simulator,
1192 const WellStateType& well_state,
1193 const std::vector<Scalar>& B_avg,
1194 DeferredLogger& deferred_logger,
1195 const bool relax_tolerance) const
1196 {
1197 // the following implementation assume that the polymer is always after the w-o-g phases
1198 // For the polymer, energy and foam cases, there is one more mass balance equations of reservoir than wells
1199 assert((int(B_avg.size()) == this->num_conservation_quantities_) || has_polymer || has_energy || has_foam || has_brine || has_zFraction || has_bioeffects);
1200
1201 Scalar tol_wells = this->param_.tolerance_wells_;
1202 // use stricter tolerance for stopped wells and wells under zero rate target control.
1203 constexpr Scalar stopped_factor = 1.e-4;
1204 // use stricter tolerance for dynamic thp to ameliorate network convergence
1205 constexpr Scalar dynamic_thp_factor = 1.e-1;
1206 if (this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger)) {
1207 tol_wells = tol_wells*stopped_factor;
1208 } else if (this->getDynamicThpLimit()) {
1209 tol_wells = tol_wells*dynamic_thp_factor;
1210 }
1211
1212 std::vector<Scalar> res;
1213 ConvergenceReport report = this->StdWellEval::getWellConvergence(well_state,
1214 B_avg,
1215 this->param_.max_residual_allowed_,
1216 tol_wells,
1217 this->param_.relaxed_tolerance_flow_well_,
1218 relax_tolerance,
1219 this->wellIsStopped(),
1220 res,
1221 deferred_logger);
1222
1223 checkConvergenceExtraEqs(res, report);
1224
1225 return report;
1226 }
1227
1228
1229
1230
1231
1232 template<typename TypeTag>
1233 void
1235 updateProductivityIndex(const Simulator& simulator,
1236 const WellProdIndexCalculator<Scalar>& wellPICalc,
1237 WellStateType& well_state,
1238 DeferredLogger& deferred_logger) const
1239 {
1240 auto fluidState = [&simulator, this](const int perf)
1241 {
1242 const auto cell_idx = this->well_cells_[perf];
1243 return simulator.model()
1244 .intensiveQuantities(cell_idx, /*timeIdx=*/ 0).fluidState();
1245 };
1246
1247 const int np = this->number_of_phases_;
1248 auto setToZero = [np](Scalar* x) -> void
1249 {
1250 std::fill_n(x, np, 0.0);
1251 };
1252
1253 auto addVector = [np](const Scalar* src, Scalar* dest) -> void
1254 {
1255 std::transform(src, src + np, dest, dest, std::plus<>{});
1256 };
1257
1258 auto& ws = well_state.well(this->index_of_well_);
1259 auto& perf_data = ws.perf_data;
1260 auto* wellPI = ws.productivity_index.data();
1261 auto* connPI = perf_data.prod_index.data();
1262
1263 setToZero(wellPI);
1264
1265 const auto preferred_phase = this->well_ecl_.getPreferredPhase();
1266 auto subsetPerfID = 0;
1267
1268 for (const auto& perf : *this->perf_data_) {
1269 auto allPerfID = perf.ecl_index;
1270
1271 auto connPICalc = [&wellPICalc, allPerfID](const Scalar mobility) -> Scalar
1272 {
1273 return wellPICalc.connectionProdIndStandard(allPerfID, mobility);
1274 };
1275
1276 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
1277 getMobility(simulator, static_cast<int>(subsetPerfID), mob, deferred_logger);
1278
1279 const auto& fs = fluidState(subsetPerfID);
1280 setToZero(connPI);
1281
1282 if (this->isInjector()) {
1283 this->computeConnLevelInjInd(fs, preferred_phase, connPICalc,
1284 mob, connPI, deferred_logger);
1285 }
1286 else { // Production or zero flow rate
1287 this->computeConnLevelProdInd(fs, connPICalc, mob, connPI);
1288 }
1289
1290 addVector(connPI, wellPI);
1291
1292 ++subsetPerfID;
1293 connPI += np;
1294 }
1295
1296 // Sum with communication in case of distributed well.
1297 const auto& comm = this->parallel_well_info_.communication();
1298 if (comm.size() > 1) {
1299 comm.sum(wellPI, np);
1300 }
1301
1302 assert ((static_cast<int>(subsetPerfID) == this->number_of_local_perforations_) &&
1303 "Internal logic error in processing connections for PI/II");
1304 }
1305
1306
1307
1308 template<typename TypeTag>
1311 const WellStateType& well_state,
1312 const WellConnectionProps& props,
1313 DeferredLogger& deferred_logger)
1314 {
1315 // Cell level dynamic property call-back functions as fall-back
1316 // option for calculating connection level mixture densities in
1317 // stopped or zero-rate producer wells.
1318 const auto prop_func = typename StdWellEval::StdWellConnections::DensityPropertyFunctions {
1319 // This becomes slightly more palatable with C++20's designated
1320 // initialisers.
1321
1322 // mobility: Phase mobilities in specified cell.
1323 [&model = simulator.model()](const int cell,
1324 const std::vector<int>& phases,
1325 std::vector<Scalar>& mob)
1326 {
1327 const auto& iq = model.intensiveQuantities(cell, /* time_idx = */ 0);
1328
1329 std::transform(phases.begin(), phases.end(), mob.begin(),
1330 [&iq](const int phase) { return iq.mobility(phase).value(); });
1331 },
1332
1333 // densityInCell: Reservoir condition phase densities in
1334 // specified cell.
1335 [&model = simulator.model()](const int cell,
1336 const std::vector<int>& phases,
1337 std::vector<Scalar>& rho)
1338 {
1339 const auto& fs = model.intensiveQuantities(cell, /* time_idx = */ 0).fluidState();
1340
1341 std::transform(phases.begin(), phases.end(), rho.begin(),
1342 [&fs](const int phase) { return fs.density(phase).value(); });
1343 }
1344 };
1345
1346 const auto stopped_or_zero_rate_target = this->
1347 stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
1348
1349 this->connections_
1350 .computeProperties(stopped_or_zero_rate_target, well_state,
1351 prop_func, props, deferred_logger);
1352 // density was updated
1353 cachedRefDensity = this->connections_.rho(0);
1354 if (this->parallel_well_info_.communication().size() > 1) {
1355 cachedRefDensity = this->parallel_well_info_.broadcastFirstPerforationValue(cachedRefDensity);
1356 }
1357 }
1358
1359
1360
1361
1362
1363 template<typename TypeTag>
1364 void
1367 const WellStateType& well_state,
1368 DeferredLogger& deferred_logger)
1369 {
1370 const auto props = computePropertiesForWellConnectionPressures
1371 (simulator, well_state);
1372
1373 computeWellConnectionDensitesPressures(simulator, well_state,
1374 props, deferred_logger);
1375 }
1376
1377
1378
1379
1380
1381 template<typename TypeTag>
1382 void
1384 solveEqAndUpdateWellState(const Simulator& simulator,
1385 WellStateType& well_state,
1386 DeferredLogger& deferred_logger)
1387 {
1388 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1389
1390 // We assemble the well equations, then we check the convergence,
1391 // which is why we do not put the assembleWellEq here.
1392 BVectorWell dx_well(1);
1393 dx_well[0].resize(this->primary_variables_.numWellEq());
1394 this->linSys_.solve( dx_well);
1395
1396 updateWellState(simulator, dx_well, well_state, deferred_logger);
1397 }
1398
1399
1400
1401
1402
1403 template<typename TypeTag>
1404 void
1407 const WellStateType& well_state,
1408 DeferredLogger& deferred_logger)
1409 {
1410 updatePrimaryVariables(simulator, well_state, deferred_logger);
1411 computeWellConnectionPressures(simulator, well_state, deferred_logger);
1412 this->computeAccumWell();
1413 }
1414
1415
1416
1417 template<typename TypeTag>
1418 void
1420 apply(const BVector& x, BVector& Ax) const
1421 {
1422 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1423
1424 if (this->param_.matrix_add_well_contributions_)
1425 {
1426 // Contributions are already in the matrix itself
1427 return;
1428 }
1429
1430 this->linSys_.apply(x, Ax);
1431 }
1432
1433
1434
1435
1436 template<typename TypeTag>
1437 void
1439 apply(BVector& r) const
1440 {
1441 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1442
1443 this->linSys_.apply(r);
1444 }
1445
1446
1447
1448
1449 template<typename TypeTag>
1450 void
1453 const BVector& x,
1454 WellStateType& well_state,
1455 DeferredLogger& deferred_logger)
1456 {
1457 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1458
1459 BVectorWell xw(1);
1460 xw[0].resize(this->primary_variables_.numWellEq());
1461
1462 this->linSys_.recoverSolutionWell(x, xw);
1463 updateWellState(simulator, xw, well_state, deferred_logger);
1464 }
1465
1466
1467
1468
1469 template<typename TypeTag>
1470 void
1472 computeWellRatesWithBhp(const Simulator& simulator,
1473 const Scalar& bhp,
1474 std::vector<Scalar>& well_flux,
1475 DeferredLogger& deferred_logger) const
1476 {
1477 OPM_TIMEFUNCTION();
1478 const int np = this->number_of_phases_;
1479 well_flux.resize(np, 0.0);
1480
1481 const bool allow_cf = this->getAllowCrossFlow();
1482
1483 for (int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
1484 const int cell_idx = this->well_cells_[perf];
1485 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
1486 // flux for each perforation
1487 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.);
1488 getMobility(simulator, perf, mob, deferred_logger);
1489 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(intQuants, cell_idx);
1490 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1491 const std::vector<Scalar> Tw = this->wellIndex(perf, intQuants, trans_mult, wellstate_nupcol);
1492
1493 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.);
1494 PerforationRates<Scalar> perf_rates;
1495 computePerfRate(intQuants, mob, bhp, Tw, perf, allow_cf,
1496 cq_s, perf_rates, deferred_logger);
1497
1498 for(int p = 0; p < np; ++p) {
1499 well_flux[FluidSystem::activeCompToActivePhaseIdx(p)] += cq_s[p];
1500 }
1501
1502 // the solvent contribution is added to the gas potentials
1503 if constexpr (has_solvent) {
1504 assert(FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx));
1505 // TODO: should we use compIdx here?
1506 const int gas_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
1507 well_flux[gas_pos] += cq_s[Indices::contiSolventEqIdx];
1508 }
1509 }
1510 this->parallel_well_info_.communication().sum(well_flux.data(), well_flux.size());
1511 }
1512
1513
1514
1515 template<typename TypeTag>
1516 void
1519 const Scalar& bhp,
1520 const WellGroupHelperType& wgHelper,
1521 std::vector<Scalar>& well_flux,
1522 DeferredLogger& deferred_logger) const
1523 {
1524 // creating a copy of the well itself, to avoid messing up the explicit information
1525 // during this copy, the only information not copied properly is the well controls
1526 StandardWell<TypeTag> well_copy(*this);
1527 well_copy.resetDampening();
1528
1529 // iterate to get a more accurate well density
1530 // create a copy of the well_state to use. If the operability checking is sucessful, we use this one
1531 // to replace the original one
1532 WellGroupHelperType wgHelper_copy = wgHelper;
1533 WellStateType well_state_copy = wgHelper_copy.wellState();
1534 // Ensure that wgHelper_copy uses well_state_copy as WellState for the rest of this function,
1535 // and the guard ensures that the original well state is restored at scope exit, i.e. at
1536 // the end of this function.
1537 auto guard = wgHelper_copy.pushWellState(well_state_copy);
1538
1539 // Get the current controls.
1540 const auto& summary_state = simulator.vanguard().summaryState();
1541 auto inj_controls = well_copy.well_ecl_.isInjector()
1542 ? well_copy.well_ecl_.injectionControls(summary_state)
1543 : Well::InjectionControls(0);
1544 auto prod_controls = well_copy.well_ecl_.isProducer()
1545 ? well_copy.well_ecl_.productionControls(summary_state) :
1546 Well::ProductionControls(0);
1547
1548 // Set current control to bhp, and bhp value in state, modify bhp limit in control object.
1549 auto& ws = well_state_copy.well(this->index_of_well_);
1550 if (well_copy.well_ecl_.isInjector()) {
1551 inj_controls.bhp_limit = bhp;
1552 ws.injection_cmode = Well::InjectorCMode::BHP;
1553 } else {
1554 prod_controls.bhp_limit = bhp;
1555 ws.production_cmode = Well::ProducerCMode::BHP;
1556 }
1557 ws.bhp = bhp;
1558
1559 // initialized the well rates with the potentials i.e. the well rates based on bhp
1560 const int np = this->number_of_phases_;
1561 const Scalar sign = this->well_ecl_.isInjector() ? 1.0 : -1.0;
1562 for (int phase = 0; phase < np; ++phase){
1563 well_state_copy.wellRates(this->index_of_well_)[phase]
1564 = sign * ws.well_potentials[phase];
1565 }
1566 well_copy.updatePrimaryVariables(simulator, well_state_copy, deferred_logger);
1567 well_copy.computeAccumWell();
1568
1569 const double dt = simulator.timeStepSize();
1570 const bool converged = well_copy.iterateWellEqWithControl(
1571 simulator, dt, inj_controls, prod_controls, wgHelper_copy, well_state_copy, deferred_logger
1572 );
1573 if (!converged) {
1574 const std::string msg = " well " + name() + " did not get converged during well potential calculations "
1575 " potentials are computed based on unconverged solution";
1576 deferred_logger.debug(msg);
1577 }
1578 well_copy.updatePrimaryVariables(simulator, well_state_copy, deferred_logger);
1579 well_copy.computeWellConnectionPressures(simulator, well_state_copy, deferred_logger);
1580 well_copy.computeWellRatesWithBhp(simulator, bhp, well_flux, deferred_logger);
1581 }
1582
1583
1584
1585
1586 template<typename TypeTag>
1587 std::vector<typename StandardWell<TypeTag>::Scalar>
1590 const WellGroupHelperType& wgHelper,
1591 DeferredLogger& deferred_logger,
1592 const WellStateType& well_state) const
1593 {
1594 std::vector<Scalar> potentials(this->number_of_phases_, 0.0);
1595 const auto& summary_state = simulator.vanguard().summaryState();
1596
1597 const auto& well = this->well_ecl_;
1598 if (well.isInjector()){
1599 const auto& controls = this->well_ecl_.injectionControls(summary_state);
1600 auto bhp_at_thp_limit = computeBhpAtThpLimitInj(simulator, wgHelper, summary_state, deferred_logger);
1601 if (bhp_at_thp_limit) {
1602 const Scalar bhp = std::min(*bhp_at_thp_limit,
1603 static_cast<Scalar>(controls.bhp_limit));
1604 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1605 } else {
1606 deferred_logger.warning("FAILURE_GETTING_CONVERGED_POTENTIAL",
1607 "Failed in getting converged thp based potential calculation for well "
1608 + name() + ". Instead the bhp based value is used");
1609 const Scalar bhp = controls.bhp_limit;
1610 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1611 }
1612 } else {
1613 computeWellRatesWithThpAlqProd(
1614 simulator, wgHelper, summary_state,
1615 deferred_logger, potentials, this->getALQ(well_state)
1616 );
1617 }
1618
1619 return potentials;
1620 }
1621
1622 template<typename TypeTag>
1623 bool
1626 const WellGroupHelperType& wgHelper,
1627 std::vector<Scalar>& well_potentials,
1628 DeferredLogger& deferred_logger) const
1629 {
1630 // Create a copy of the well.
1631 // TODO: check if we can avoid taking multiple copies. Call from updateWellPotentials
1632 // is allready a copy, but not from other calls.
1633 StandardWell<TypeTag> well_copy(*this);
1634
1635 // store a copy of the well state, we don't want to update the real well state
1636 WellStateType well_state_copy = wgHelper.wellState();
1637 WellGroupHelperType wgHelper_copy = wgHelper;
1638 // Ensure that wgHelper_copy uses well_state_copy as WellState for the rest of this function,
1639 // and the guard ensures that the original well state is restored at scope exit, i.e. at
1640 // the end of this function.
1641 auto guard = wgHelper_copy.pushWellState(well_state_copy);
1642 auto& ws = well_state_copy.well(this->index_of_well_);
1643
1644 // get current controls
1645 const auto& summary_state = simulator.vanguard().summaryState();
1646 auto inj_controls = well_copy.well_ecl_.isInjector()
1647 ? well_copy.well_ecl_.injectionControls(summary_state)
1648 : Well::InjectionControls(0);
1649 auto prod_controls = well_copy.well_ecl_.isProducer()
1650 ? well_copy.well_ecl_.productionControls(summary_state) :
1651 Well::ProductionControls(0);
1652
1653 // prepare/modify well state and control
1654 well_copy.onlyKeepBHPandTHPcontrols(summary_state, well_state_copy, inj_controls, prod_controls);
1655
1656 // update connection pressures relative to updated bhp to get better estimate of connection dp
1657 const int num_perf = ws.perf_data.size();
1658 for (int perf = 0; perf < num_perf; ++perf) {
1659 ws.perf_data.pressure[perf] = ws.bhp + well_copy.connections_.pressure_diff(perf);
1660 }
1661 // initialize rates from previous potentials
1662 const int np = this->number_of_phases_;
1663 bool trivial = true;
1664 for (int phase = 0; phase < np; ++phase){
1665 trivial = trivial && (ws.well_potentials[phase] == 0.0) ;
1666 }
1667 if (!trivial) {
1668 const Scalar sign = well_copy.well_ecl_.isInjector() ? 1.0 : -1.0;
1669 for (int phase = 0; phase < np; ++phase) {
1670 ws.surface_rates[phase] = sign * ws.well_potentials[phase];
1671 }
1672 }
1673
1674 well_copy.calculateExplicitQuantities(simulator, well_state_copy, deferred_logger);
1675 const double dt = simulator.timeStepSize();
1676 // iterate to get a solution at the given bhp.
1677 bool converged = false;
1678 if (this->well_ecl_.isProducer()) {
1679 converged = well_copy.solveWellWithOperabilityCheck(
1680 simulator, dt, inj_controls, prod_controls, wgHelper_copy, well_state_copy, deferred_logger
1681 );
1682 } else {
1683 converged = well_copy.iterateWellEqWithSwitching(
1684 simulator, dt, inj_controls, prod_controls, wgHelper_copy, well_state_copy, deferred_logger
1685 );
1686 }
1687
1688 // fetch potentials (sign is updated on the outside).
1689 well_potentials.clear();
1690 well_potentials.resize(np, 0.0);
1691 for (int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1692 if (has_solvent && comp_idx == Indices::contiSolventEqIdx) continue; // we do not store the solvent in the well_potentials
1693 const EvalWell rate = well_copy.primary_variables_.getQs(comp_idx);
1694 well_potentials[FluidSystem::activeCompToActivePhaseIdx(comp_idx)] = rate.value();
1695 }
1696
1697 // the solvent contribution is added to the gas potentials
1698 if constexpr (has_solvent) {
1699 assert(FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx));
1700 // TODO: should we use compIdx here?
1701 const int gas_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
1702 const EvalWell rate = well_copy.primary_variables_.getQs(Indices::contiSolventEqIdx);
1703 well_potentials[gas_pos] += rate.value();
1704 }
1705 return converged;
1706 }
1707
1708
1709 template<typename TypeTag>
1713 const WellGroupHelperType& wgHelper,
1714 const SummaryState &summary_state,
1715 DeferredLogger& deferred_logger,
1716 std::vector<Scalar>& potentials,
1717 Scalar alq) const
1718 {
1719 Scalar bhp;
1720 auto bhp_at_thp_limit = computeBhpAtThpLimitProdWithAlq(
1721 simulator, wgHelper, summary_state, alq, deferred_logger, /*iterate_if_no_solution */ true);
1722 if (bhp_at_thp_limit) {
1723 const auto& controls = this->well_ecl_.productionControls(summary_state);
1724 bhp = std::max(*bhp_at_thp_limit,
1725 static_cast<Scalar>(controls.bhp_limit));
1726 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1727 }
1728 else {
1729 deferred_logger.warning("FAILURE_GETTING_CONVERGED_POTENTIAL",
1730 "Failed in getting converged thp based potential calculation for well "
1731 + name() + ". Instead the bhp based value is used");
1732 const auto& controls = this->well_ecl_.productionControls(summary_state);
1733 bhp = controls.bhp_limit;
1734 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1735 }
1736 return bhp;
1737 }
1738
1739 template<typename TypeTag>
1740 void
1743 const WellGroupHelperType& wgHelper,
1744 const SummaryState& summary_state,
1745 DeferredLogger& deferred_logger,
1746 std::vector<Scalar>& potentials,
1747 Scalar alq) const
1748 {
1749 /*double bhp =*/
1750 computeWellRatesAndBhpWithThpAlqProd(simulator,
1751 wgHelper,
1752 summary_state,
1753 deferred_logger,
1754 potentials,
1755 alq);
1756 }
1757
1758 template<typename TypeTag>
1759 void
1761 computeWellPotentials(const Simulator& simulator,
1762 const WellStateType& well_state,
1763 const WellGroupHelperType& wgHelper,
1764 std::vector<Scalar>& well_potentials,
1765 DeferredLogger& deferred_logger) // const
1766 {
1767 const auto [compute_potential, bhp_controlled_well] =
1769
1770 if (!compute_potential) {
1771 return;
1772 }
1773
1774 bool converged_implicit = false;
1775 // for newly opened wells we dont compute the potentials implicit
1776 // group controlled wells with defaulted guiderates will have zero targets as
1777 // the potentials are used to compute the well fractions.
1778 if (this->param_.local_well_solver_control_switching_ && !(this->changed_to_open_this_step_ && this->wellUnderZeroRateTarget(simulator, well_state, deferred_logger))) {
1779 converged_implicit = computeWellPotentialsImplicit(
1780 simulator, wgHelper, well_potentials, deferred_logger
1781 );
1782 }
1783 if (!converged_implicit) {
1784 // does the well have a THP related constraint?
1785 const auto& summaryState = simulator.vanguard().summaryState();
1786 if (!Base::wellHasTHPConstraints(summaryState) || bhp_controlled_well) {
1787 // get the bhp value based on the bhp constraints
1788 Scalar bhp = WellBhpThpCalculator(*this).mostStrictBhpFromBhpLimits(summaryState);
1789
1790 // In some very special cases the bhp pressure target are
1791 // temporary violated. This may lead to too small or negative potentials
1792 // that could lead to premature shutting of wells.
1793 // As a remedy the bhp that gives the largest potential is used.
1794 // For converged cases, ws.bhp <=bhp for injectors and ws.bhp >= bhp,
1795 // and the potentials will be computed using the limit as expected.
1796 const auto& ws = well_state.well(this->index_of_well_);
1797 if (this->isInjector())
1798 bhp = std::max(ws.bhp, bhp);
1799 else
1800 bhp = std::min(ws.bhp, bhp);
1801
1802 assert(std::abs(bhp) != std::numeric_limits<Scalar>::max());
1803 computeWellRatesWithBhpIterations(simulator, bhp, wgHelper, well_potentials, deferred_logger);
1804 } else {
1805 // the well has a THP related constraint
1806 well_potentials = computeWellPotentialWithTHP(simulator, wgHelper, deferred_logger, well_state);
1807 }
1808 }
1809
1810 this->checkNegativeWellPotentials(well_potentials,
1811 this->param_.check_well_operability_,
1812 deferred_logger);
1813 }
1814
1815
1816
1817
1818
1819
1820
1821 template<typename TypeTag>
1824 connectionDensity([[maybe_unused]] const int globalConnIdx,
1825 const int openConnIdx) const
1826 {
1827 return (openConnIdx < 0)
1828 ? 0.0
1829 : this->connections_.rho(openConnIdx);
1830 }
1831
1832
1833
1834
1835
1836 template<typename TypeTag>
1837 void
1839 updatePrimaryVariables(const Simulator& simulator,
1840 const WellStateType& well_state,
1841 DeferredLogger& deferred_logger)
1842 {
1843 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1844
1845 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
1846 this->primary_variables_.update(well_state, stop_or_zero_rate_target, deferred_logger);
1847
1848 // other primary variables related to polymer injection
1849 if constexpr (Base::has_polymermw) {
1850 this->primary_variables_.updatePolyMW(well_state);
1851 }
1852
1853 this->primary_variables_.checkFinite(deferred_logger);
1854 }
1855
1856
1857
1858
1859 template<typename TypeTag>
1862 getRefDensity() const
1863 {
1864 return cachedRefDensity;
1865 }
1866
1867
1868
1869
1870 template<typename TypeTag>
1871 void
1874 const int perf,
1875 std::vector<EvalWell>& mob,
1876 DeferredLogger& deferred_logger) const
1877 {
1878 const int cell_idx = this->well_cells_[perf];
1879 const auto& int_quant = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
1880 const EvalWell polymer_concentration = this->extendEval(int_quant.polymerConcentration());
1881
1882 // TODO: not sure should based on the well type or injecting/producing peforations
1883 // it can be different for crossflow
1884 if (this->isInjector()) {
1885 // assume fully mixing within injecting wellbore
1886 const auto& visc_mult_table = PolymerModule::plyviscViscosityMultiplierTable(int_quant.pvtRegionIndex());
1887 const unsigned waterCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
1888 mob[waterCompIdx] /= (this->extendEval(int_quant.waterViscosityCorrection()) * visc_mult_table.eval(polymer_concentration, /*extrapolate=*/true) );
1889 }
1890
1891 if (PolymerModule::hasPlyshlog()) {
1892 // we do not calculate the shear effects for injection wells when they do not
1893 // inject polymer.
1894 if (this->isInjector() && this->wpolymer() == 0.) {
1895 return;
1896 }
1897 // compute the well water velocity with out shear effects.
1898 // TODO: do we need to turn on crossflow here?
1899 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
1900 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
1901
1902 std::vector<EvalWell> cq_s(this->num_conservation_quantities_, 0.);
1903 PerforationRates<Scalar> perf_rates;
1904 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quant, cell_idx);
1905 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1906 const std::vector<Scalar> Tw = this->wellIndex(perf, int_quant, trans_mult, wellstate_nupcol);
1907 computePerfRate(int_quant, mob, bhp, Tw, perf, allow_cf, cq_s,
1908 perf_rates, deferred_logger);
1909 // TODO: make area a member
1910 const Scalar area = 2 * M_PI * this->perf_rep_radius_[perf] * this->perf_length_[perf];
1911 const auto& material_law_manager = simulator.problem().materialLawManager();
1912 const auto& scaled_drainage_info =
1913 material_law_manager->oilWaterScaledEpsInfoDrainage(cell_idx);
1914 const Scalar swcr = scaled_drainage_info.Swcr;
1915 const EvalWell poro = this->extendEval(int_quant.porosity());
1916 const EvalWell sw = this->extendEval(int_quant.fluidState().saturation(FluidSystem::waterPhaseIdx));
1917 // guard against zero porosity and no water
1918 const EvalWell denom = max( (area * poro * (sw - swcr)), 1e-12);
1919 const unsigned waterCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
1920 EvalWell water_velocity = cq_s[waterCompIdx] / denom * this->extendEval(int_quant.fluidState().invB(FluidSystem::waterPhaseIdx));
1921
1922 if (PolymerModule::hasShrate()) {
1923 // the equation for the water velocity conversion for the wells and reservoir are from different version
1924 // of implementation. It can be changed to be more consistent when possible.
1925 water_velocity *= PolymerModule::shrate( int_quant.pvtRegionIndex() ) / this->bore_diameters_[perf];
1926 }
1927 const EvalWell shear_factor = PolymerModule::computeShearFactor(polymer_concentration,
1928 int_quant.pvtRegionIndex(),
1929 water_velocity);
1930 // modify the mobility with the shear factor.
1931 mob[waterCompIdx] /= shear_factor;
1932 }
1933 }
1934
1935 template<typename TypeTag>
1936 void
1938 {
1939 this->linSys_.extract(jacobian);
1940 }
1941
1942
1943 template <typename TypeTag>
1944 void
1946 const BVector& weights,
1947 const int pressureVarIndex,
1948 const bool use_well_weights,
1949 const WellStateType& well_state) const
1950 {
1951 this->linSys_.extractCPRPressureMatrix(jacobian,
1952 weights,
1953 pressureVarIndex,
1954 use_well_weights,
1955 *this,
1956 Bhp,
1957 well_state);
1958 }
1959
1960
1961
1962 template<typename TypeTag>
1965 pskinwater(const Scalar throughput,
1966 const EvalWell& water_velocity,
1967 DeferredLogger& deferred_logger) const
1968 {
1969 if constexpr (Base::has_polymermw) {
1970 const int water_table_id = this->polymerWaterTable_();
1971 if (water_table_id <= 0) {
1972 OPM_DEFLOG_THROW(std::runtime_error,
1973 fmt::format("Unused SKPRWAT table id used for well {}", name()),
1974 deferred_logger);
1975 }
1976 const auto& water_table_func = PolymerModule::getSkprwatTable(water_table_id);
1977 const EvalWell throughput_eval{throughput};
1978 // the skin pressure when injecting water, which also means the polymer concentration is zero
1979 EvalWell pskin_water = water_table_func.eval(throughput_eval, water_velocity);
1980 return pskin_water;
1981 } else {
1982 OPM_DEFLOG_THROW(std::runtime_error,
1983 fmt::format("Polymermw is not activated, while injecting "
1984 "skin pressure is requested for well {}", name()),
1985 deferred_logger);
1986 }
1987 }
1988
1989
1990
1991
1992
1993 template<typename TypeTag>
1996 pskin(const Scalar throughput,
1997 const EvalWell& water_velocity,
1998 const EvalWell& poly_inj_conc,
1999 DeferredLogger& deferred_logger) const
2000 {
2001 if constexpr (Base::has_polymermw) {
2002 const Scalar sign = water_velocity >= 0. ? 1.0 : -1.0;
2003 const EvalWell water_velocity_abs = abs(water_velocity);
2004 if (poly_inj_conc == 0.) {
2005 return sign * pskinwater(throughput, water_velocity_abs, deferred_logger);
2006 }
2007 const int polymer_table_id = this->polymerTable_();
2008 if (polymer_table_id <= 0) {
2009 OPM_DEFLOG_THROW(std::runtime_error,
2010 fmt::format("Unavailable SKPRPOLY table id used for well {}", name()),
2011 deferred_logger);
2012 }
2013 const auto& skprpolytable = PolymerModule::getSkprpolyTable(polymer_table_id);
2014 const Scalar reference_concentration = skprpolytable.refConcentration;
2015 const EvalWell throughput_eval{throughput};
2016 // the skin pressure when injecting water, which also means the polymer concentration is zero
2017 const EvalWell pskin_poly = skprpolytable.table_func.eval(throughput_eval, water_velocity_abs);
2018 if (poly_inj_conc == reference_concentration) {
2019 return sign * pskin_poly;
2020 }
2021 // poly_inj_conc != reference concentration of the table, then some interpolation will be required
2022 const EvalWell pskin_water = pskinwater(throughput, water_velocity_abs, deferred_logger);
2023 const EvalWell pskin = pskin_water + (pskin_poly - pskin_water) / reference_concentration * poly_inj_conc;
2024 return sign * pskin;
2025 } else {
2026 OPM_DEFLOG_THROW(std::runtime_error,
2027 fmt::format("Polymermw is not activated, while injecting "
2028 "skin pressure is requested for well {}", name()),
2029 deferred_logger);
2030 }
2031 }
2032
2033
2034
2035
2036
2037 template<typename TypeTag>
2040 wpolymermw(const Scalar throughput,
2041 const EvalWell& water_velocity,
2042 DeferredLogger& deferred_logger) const
2043 {
2044 if constexpr (Base::has_polymermw) {
2045 const int table_id = this->polymerInjTable_();
2046 const auto& table_func = PolymerModule::getPlymwinjTable(table_id);
2047 const EvalWell throughput_eval{throughput};
2048 EvalWell molecular_weight{0.};
2049 if (this->wpolymer() == 0.) { // not injecting polymer
2050 return molecular_weight;
2051 }
2052 molecular_weight = table_func.eval(throughput_eval, abs(water_velocity));
2053 return molecular_weight;
2054 } else {
2055 OPM_DEFLOG_THROW(std::runtime_error,
2056 fmt::format("Polymermw is not activated, while injecting "
2057 "polymer molecular weight is requested for well {}", name()),
2058 deferred_logger);
2059 }
2060 }
2061
2062
2063
2064
2065
2066 template<typename TypeTag>
2067 void
2069 updateWaterThroughput([[maybe_unused]] const double dt,
2070 WellStateType& well_state) const
2071 {
2072 if constexpr (Base::has_polymermw) {
2073 if (!this->isInjector()) {
2074 return;
2075 }
2076
2077 auto& perf_water_throughput = well_state.well(this->index_of_well_)
2078 .perf_data.water_throughput;
2079
2080 for (int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
2081 const Scalar perf_water_vel =
2082 this->primary_variables_.value(Bhp + 1 + perf);
2083
2084 // we do not consider the formation damage due to water
2085 // flowing from reservoir into wellbore
2086 if (perf_water_vel > Scalar{0}) {
2087 perf_water_throughput[perf] += perf_water_vel * dt;
2088 }
2089 }
2090 }
2091 }
2092
2093
2094
2095
2096
2097 template<typename TypeTag>
2098 void
2100 handleInjectivityRate(const Simulator& simulator,
2101 const int perf,
2102 std::vector<EvalWell>& cq_s) const
2103 {
2104 const int cell_idx = this->well_cells_[perf];
2105 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2106 const auto& fs = int_quants.fluidState();
2107 const EvalWell b_w = this->extendEval(fs.invB(FluidSystem::waterPhaseIdx));
2108 const Scalar area = M_PI * this->bore_diameters_[perf] * this->perf_length_[perf];
2109 const int wat_vel_index = Bhp + 1 + perf;
2110 const unsigned water_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
2111
2112 // water rate is update to use the form from water velocity, since water velocity is
2113 // a primary variable now
2114 cq_s[water_comp_idx] = area * this->primary_variables_.eval(wat_vel_index) * b_w;
2115 }
2116
2117
2118
2119
2120 template<typename TypeTag>
2121 void
2123 handleInjectivityEquations(const Simulator& simulator,
2124 const WellStateType& well_state,
2125 const int perf,
2126 const EvalWell& water_flux_s,
2127 DeferredLogger& deferred_logger)
2128 {
2129 const int cell_idx = this->well_cells_[perf];
2130 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2131 const auto& fs = int_quants.fluidState();
2132 const EvalWell b_w = this->extendEval(fs.invB(FluidSystem::waterPhaseIdx));
2133 const EvalWell water_flux_r = water_flux_s / b_w;
2134 const Scalar area = M_PI * this->bore_diameters_[perf] * this->perf_length_[perf];
2135 const EvalWell water_velocity = water_flux_r / area;
2136 const int wat_vel_index = Bhp + 1 + perf;
2137
2138 // equation for the water velocity
2139 const EvalWell eq_wat_vel = this->primary_variables_.eval(wat_vel_index) - water_velocity;
2140
2141 const auto& ws = well_state.well(this->index_of_well_);
2142 const auto& perf_data = ws.perf_data;
2143 const auto& perf_water_throughput = perf_data.water_throughput;
2144 const Scalar throughput = perf_water_throughput[perf];
2145 const int pskin_index = Bhp + 1 + this->number_of_local_perforations_ + perf;
2146
2147 const EvalWell poly_conc(this->wpolymer());
2148
2149 // equation for the skin pressure
2150 const EvalWell eq_pskin = this->primary_variables_.eval(pskin_index)
2151 - pskin(throughput, this->primary_variables_.eval(wat_vel_index), poly_conc, deferred_logger);
2152
2154 assembleInjectivityEq(eq_pskin,
2155 eq_wat_vel,
2156 pskin_index,
2157 wat_vel_index,
2158 perf,
2159 this->primary_variables_.numWellEq(),
2160 this->linSys_);
2161 }
2162
2163
2164
2165
2166
2167 template<typename TypeTag>
2168 void
2170 checkConvergenceExtraEqs(const std::vector<Scalar>& res,
2171 ConvergenceReport& report) const
2172 {
2173 // if different types of extra equations are involved, this function needs to be refactored further
2174
2175 // checking the convergence of the extra equations related to polymer injectivity
2176 if constexpr (Base::has_polymermw) {
2177 WellConvergence(*this).
2178 checkConvergencePolyMW(res, Bhp, this->param_.max_residual_allowed_, report);
2179 }
2180 }
2181
2182
2183
2184
2185
2186 template<typename TypeTag>
2187 void
2189 updateConnectionRatePolyMW(const EvalWell& cq_s_poly,
2190 const IntensiveQuantities& int_quants,
2191 const WellStateType& well_state,
2192 const int perf,
2193 std::vector<RateVector>& connectionRates,
2194 DeferredLogger& deferred_logger) const
2195 {
2196 // the source term related to transport of molecular weight
2197 EvalWell cq_s_polymw = cq_s_poly;
2198 if (this->isInjector()) {
2199 const int wat_vel_index = Bhp + 1 + perf;
2200 const EvalWell water_velocity = this->primary_variables_.eval(wat_vel_index);
2201 if (water_velocity > 0.) { // injecting
2202 const auto& ws = well_state.well(this->index_of_well_);
2203 const auto& perf_water_throughput = ws.perf_data.water_throughput;
2204 const Scalar throughput = perf_water_throughput[perf];
2205 const EvalWell molecular_weight = wpolymermw(throughput, water_velocity, deferred_logger);
2206 cq_s_polymw *= molecular_weight;
2207 } else {
2208 // we do not consider the molecular weight from the polymer
2209 // going-back to the wellbore through injector
2210 cq_s_polymw *= 0.;
2211 }
2212 } else if (this->isProducer()) {
2213 if (cq_s_polymw < 0.) {
2214 cq_s_polymw *= this->extendEval(int_quants.polymerMoleWeight() );
2215 } else {
2216 // we do not consider the molecular weight from the polymer
2217 // re-injecting back through producer
2218 cq_s_polymw *= 0.;
2219 }
2220 }
2221 connectionRates[perf][Indices::contiPolymerMWEqIdx] = Base::restrictEval(cq_s_polymw);
2222 }
2223
2224
2225
2226
2227
2228 template<typename TypeTag>
2229 std::optional<typename StandardWell<TypeTag>::Scalar>
2232 const Simulator& simulator,
2233 const WellGroupHelperType& wgHelper,
2234 const SummaryState& summary_state,
2235 DeferredLogger& deferred_logger) const
2236 {
2237 return computeBhpAtThpLimitProdWithAlq(simulator,
2238 wgHelper,
2239 summary_state,
2240 this->getALQ(well_state),
2241 deferred_logger,
2242 /*iterate_if_no_solution */ true);
2243 }
2244
2245 template<typename TypeTag>
2246 std::optional<typename StandardWell<TypeTag>::Scalar>
2249 const WellGroupHelperType& wgHelper,
2250 const SummaryState& summary_state,
2251 const Scalar alq_value,
2252 DeferredLogger& deferred_logger,
2253 bool iterate_if_no_solution) const
2254 {
2255 OPM_TIMEFUNCTION();
2256 // Make the frates() function.
2257 auto frates = [this, &simulator, &deferred_logger](const Scalar bhp) {
2258 // Not solving the well equations here, which means we are
2259 // calculating at the current Fg/Fw values of the
2260 // well. This does not matter unless the well is
2261 // crossflowing, and then it is likely still a good
2262 // approximation.
2263 std::vector<Scalar> rates(3);
2264 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2265 this->adaptRatesForVFP(rates);
2266 return rates;
2267 };
2268
2269 Scalar max_pressure = 0.0;
2270 for (int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
2271 const int cell_idx = this->well_cells_[perf];
2272 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2273 const auto& fs = int_quants.fluidState();
2274 Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2275 max_pressure = std::max(max_pressure, pressure_cell);
2276 }
2277 const auto& comm = this->parallel_well_info_.communication();
2278 if (comm.size() > 1) {
2279 max_pressure = comm.max(max_pressure);
2280 }
2281 auto bhpAtLimit = WellBhpThpCalculator(*this).computeBhpAtThpLimitProd(frates,
2282 summary_state,
2283 max_pressure,
2284 this->getRefDensity(),
2285 alq_value,
2286 this->getTHPConstraint(summary_state),
2287 deferred_logger);
2288
2289 if (bhpAtLimit) {
2290 auto v = frates(*bhpAtLimit);
2291 if (std::all_of(v.cbegin(), v.cend(), [](Scalar i){ return i <= 0; }) ) {
2292 return bhpAtLimit;
2293 }
2294 }
2295
2296 if (!iterate_if_no_solution)
2297 return std::nullopt;
2298
2299 auto fratesIter = [this, &simulator, &wgHelper, &deferred_logger](const Scalar bhp) {
2300 // Solver the well iterations to see if we are
2301 // able to get a solution with an update
2302 // solution
2303 std::vector<Scalar> rates(3);
2304 computeWellRatesWithBhpIterations(simulator, bhp, wgHelper, rates, deferred_logger);
2305 this->adaptRatesForVFP(rates);
2306 return rates;
2307 };
2308
2309 bhpAtLimit = WellBhpThpCalculator(*this).computeBhpAtThpLimitProd(fratesIter,
2310 summary_state,
2311 max_pressure,
2312 this->getRefDensity(),
2313 alq_value,
2314 this->getTHPConstraint(summary_state),
2315 deferred_logger);
2316
2317
2318 if (bhpAtLimit) {
2319 // should we use fratesIter here since fratesIter is used in computeBhpAtThpLimitProd above?
2320 auto v = frates(*bhpAtLimit);
2321 if (std::all_of(v.cbegin(), v.cend(), [](Scalar i){ return i <= 0; }) ) {
2322 return bhpAtLimit;
2323 }
2324 }
2325
2326 // we still don't get a valied solution.
2327 return std::nullopt;
2328 }
2329
2330
2331
2332 template<typename TypeTag>
2333 std::optional<typename StandardWell<TypeTag>::Scalar>
2335 computeBhpAtThpLimitInj(const Simulator& simulator,
2336 [[maybe_unused]] const WellGroupHelperType& wgHelper,
2337 const SummaryState& summary_state,
2338 DeferredLogger& deferred_logger) const
2339 {
2340 // Note: wgHelper parameter is currently unused in StandardWell but kept for consistency
2341 // with MultisegmentWell::computeBhpAtThpLimitInj which uses it in fratesIter lambda.
2342 // This maintains parallel API structure between well types and allows for future
2343 // enhancements without breaking the interface.
2344
2345 // Make the frates() function.
2346 auto frates = [this, &simulator, &deferred_logger](const Scalar bhp) {
2347 // Not solving the well equations here, which means we are
2348 // calculating at the current Fg/Fw values of the
2349 // well. This does not matter unless the well is
2350 // crossflowing, and then it is likely still a good
2351 // approximation.
2352 std::vector<Scalar> rates(3);
2353 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2354 return rates;
2355 };
2356
2357 return WellBhpThpCalculator(*this).computeBhpAtThpLimitInj(frates,
2358 summary_state,
2359 this->getRefDensity(),
2360 1e-6,
2361 50,
2362 true,
2363 deferred_logger);
2364 }
2365
2366
2367
2368
2369
2370 template<typename TypeTag>
2371 bool
2373 iterateWellEqWithControl(const Simulator& simulator,
2374 const double dt,
2375 const Well::InjectionControls& inj_controls,
2376 const Well::ProductionControls& prod_controls,
2377 const WellGroupHelperType& wgHelper,
2378 WellStateType& well_state,
2379 DeferredLogger& deferred_logger)
2380 {
2381 const auto& group_state = wgHelper.groupState();
2382 updatePrimaryVariables(simulator, well_state, deferred_logger);
2383
2384 const int max_iter = this->param_.max_inner_iter_wells_;
2385 int it = 0;
2386 bool converged;
2387 bool relax_convergence = false;
2388 this->regularize_ = false;
2389 do {
2390 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
2391
2392 if (it > this->param_.strict_inner_iter_wells_) {
2393 relax_convergence = true;
2394 this->regularize_ = true;
2395 }
2396
2397 auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
2398
2399 converged = report.converged();
2400 if (converged) {
2401 break;
2402 }
2403
2404 ++it;
2405 solveEqAndUpdateWellState(simulator, well_state, deferred_logger);
2406
2407 // TODO: when this function is used for well testing purposes, will need to check the controls, so that we will obtain convergence
2408 // under the most restrictive control. Based on this converged results, we can check whether to re-open the well. Either we refactor
2409 // this function or we use different functions for the well testing purposes.
2410 // We don't allow for switching well controls while computing well potentials and testing wells
2411 // updateWellControl(simulator, well_state, deferred_logger);
2412 } while (it < max_iter);
2413
2414 if (converged) {
2415 std::ostringstream sstr;
2416 sstr << " Well " << this->name() << " converged in " << it << " inner iterations.";
2417 if (relax_convergence)
2418 sstr << " (A relaxed tolerance was used after "<< this->param_.strict_inner_iter_wells_ << " iterations)";
2419
2420 // Output "converged in 0 inner iterations" messages only at
2421 // elevated verbosity levels.
2422 deferred_logger.debug(sstr.str(), OpmLog::defaultDebugVerbosityLevel + (it == 0));
2423 } else {
2424 std::ostringstream sstr;
2425 sstr << " Well " << this->name() << " did not converge in " << it << " inner iterations.";
2426 deferred_logger.debug(sstr.str());
2427 }
2428
2429 return converged;
2430 }
2431
2432
2433 template<typename TypeTag>
2434 bool
2436 iterateWellEqWithSwitching(const Simulator& simulator,
2437 const double dt,
2438 const Well::InjectionControls& inj_controls,
2439 const Well::ProductionControls& prod_controls,
2440 const WellGroupHelperType& wgHelper,
2441 WellStateType& well_state,
2442 DeferredLogger& deferred_logger,
2443 const bool fixed_control /*false*/,
2444 const bool fixed_status /*false*/)
2445 {
2446 const auto& group_state = wgHelper.groupState();
2447 updatePrimaryVariables(simulator, well_state, deferred_logger);
2448
2449 const int max_iter = this->param_.max_inner_iter_wells_;
2450 int it = 0;
2451 bool converged = false;
2452 bool relax_convergence = false;
2453 this->regularize_ = false;
2454 const auto& summary_state = simulator.vanguard().summaryState();
2455
2456 // Always take a few (more than one) iterations after a switch before allowing a new switch
2457 // The optimal number here is subject to further investigation, but it has been observerved
2458 // that unless this number is >1, we may get stuck in a cycle
2459 constexpr int min_its_after_switch = 4;
2460 // We also want to restrict the number of status switches to avoid oscillation between STOP<->OPEN
2461 const int max_status_switch = this->param_.max_well_status_switch_inner_iter_;
2462 int its_since_last_switch = min_its_after_switch;
2463 int switch_count= 0;
2464 // if we fail to solve eqs, we reset status/operability before leaving
2465 const auto well_status_orig = this->wellStatus_;
2466 const auto operability_orig = this->operability_status_;
2467 auto well_status_cur = well_status_orig;
2468 int status_switch_count = 0;
2469 // don't allow opening wells that has a stopped well status
2470 const bool allow_open = well_state.well(this->index_of_well_).status == WellStatus::OPEN;
2471 // don't allow switcing for wells under zero rate target or requested fixed status and control
2472 const bool allow_switching =
2473 !this->wellUnderZeroRateTarget(simulator, well_state, deferred_logger) &&
2474 (!fixed_control || !fixed_status) && allow_open;
2475
2476 bool changed = false;
2477 bool final_check = false;
2478 // well needs to be set operable or else solving/updating of re-opened wells is skipped
2479 this->operability_status_.resetOperability();
2480 this->operability_status_.solvable = true;
2481 do {
2482 its_since_last_switch++;
2483 if (allow_switching && its_since_last_switch >= min_its_after_switch && status_switch_count < max_status_switch){
2484 const Scalar wqTotal = this->primary_variables_.eval(WQTotal).value();
2485 changed = this->updateWellControlAndStatusLocalIteration(
2486 simulator, wgHelper, inj_controls, prod_controls, wqTotal,
2487 well_state, deferred_logger, fixed_control, fixed_status
2488 );
2489 if (changed){
2490 its_since_last_switch = 0;
2491 switch_count++;
2492 if (well_status_cur != this->wellStatus_) {
2493 well_status_cur = this->wellStatus_;
2494 status_switch_count++;
2495 }
2496 }
2497 if (!changed && final_check) {
2498 break;
2499 } else {
2500 final_check = false;
2501 }
2502 if (status_switch_count == max_status_switch) {
2503 this->wellStatus_ = well_status_orig;
2504 }
2505 }
2506
2507 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
2508
2509 if (it > this->param_.strict_inner_iter_wells_) {
2510 relax_convergence = true;
2511 this->regularize_ = true;
2512 }
2513
2514 auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
2515
2516 converged = report.converged();
2517 if (converged) {
2518 // if equations are sufficiently linear they might converge in less than min_its_after_switch
2519 // in this case, make sure all constraints are satisfied before returning
2520 if (switch_count > 0 && its_since_last_switch < min_its_after_switch) {
2521 final_check = true;
2522 its_since_last_switch = min_its_after_switch;
2523 } else {
2524 break;
2525 }
2526 }
2527
2528 ++it;
2529 solveEqAndUpdateWellState(simulator, well_state, deferred_logger);
2530
2531 } while (it < max_iter);
2532
2533 if (converged) {
2534 if (allow_switching){
2535 // update operability if status change
2536 const bool is_stopped = this->wellIsStopped();
2537 if (this->wellHasTHPConstraints(summary_state)){
2538 this->operability_status_.can_obtain_bhp_with_thp_limit = !is_stopped;
2539 this->operability_status_.obey_thp_limit_under_bhp_limit = !is_stopped;
2540 } else {
2541 this->operability_status_.operable_under_only_bhp_limit = !is_stopped;
2542 }
2543 }
2544 std::string message = fmt::format(" Well {} converged in {} inner iterations ("
2545 "{} control/status switches).", this->name(), it, switch_count);
2546 if (relax_convergence) {
2547 message.append(fmt::format(" (A relaxed tolerance was used after {} iterations)",
2548 this->param_.strict_inner_iter_wells_));
2549 }
2550 deferred_logger.debug(message, OpmLog::defaultDebugVerbosityLevel + ((it == 0) && (switch_count == 0)));
2551
2552 } else {
2553 this->wellStatus_ = well_status_orig;
2554 this->operability_status_ = operability_orig;
2555 const std::string message = fmt::format(" Well {} did not converge in {} inner iterations ("
2556 "{} switches, {} status changes).", this->name(), it, switch_count, status_switch_count);
2557 deferred_logger.debug(message);
2558 // add operability here as well ?
2559 }
2560 return converged;
2561 }
2562
2563 template<typename TypeTag>
2564 std::vector<typename StandardWell<TypeTag>::Scalar>
2566 computeCurrentWellRates(const Simulator& simulator,
2567 DeferredLogger& deferred_logger) const
2568 {
2569 // Calculate the rates that follow from the current primary variables.
2570 std::vector<Scalar> well_q_s(this->num_conservation_quantities_, 0.);
2571 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
2572 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
2573 for (int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
2574 const int cell_idx = this->well_cells_[perf];
2575 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2576 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.);
2577 getMobility(simulator, perf, mob, deferred_logger);
2578 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.);
2579 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(intQuants, cell_idx);
2580 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
2581 const std::vector<Scalar> Tw = this->wellIndex(perf, intQuants, trans_mult, wellstate_nupcol);
2582 PerforationRates<Scalar> perf_rates;
2583 computePerfRate(intQuants, mob, bhp.value(), Tw, perf, allow_cf,
2584 cq_s, perf_rates, deferred_logger);
2585 for (int comp = 0; comp < this->num_conservation_quantities_; ++comp) {
2586 well_q_s[comp] += cq_s[comp];
2587 }
2588 }
2589 const auto& comm = this->parallel_well_info_.communication();
2590 if (comm.size() > 1)
2591 {
2592 comm.sum(well_q_s.data(), well_q_s.size());
2593 }
2594 return well_q_s;
2595 }
2596
2597
2598
2599 template <typename TypeTag>
2600 std::vector<typename StandardWell<TypeTag>::Scalar>
2602 getPrimaryVars() const
2603 {
2604 const int num_pri_vars = this->primary_variables_.numWellEq();
2605 std::vector<Scalar> retval(num_pri_vars);
2606 for (int ii = 0; ii < num_pri_vars; ++ii) {
2607 retval[ii] = this->primary_variables_.value(ii);
2608 }
2609 return retval;
2610 }
2611
2612
2613
2614
2615
2616 template <typename TypeTag>
2617 int
2619 setPrimaryVars(typename std::vector<Scalar>::const_iterator it)
2620 {
2621 const int num_pri_vars = this->primary_variables_.numWellEq();
2622 for (int ii = 0; ii < num_pri_vars; ++ii) {
2623 this->primary_variables_.setValue(ii, it[ii]);
2624 }
2625 return num_pri_vars;
2626 }
2627
2628
2629 template <typename TypeTag>
2632 connectionRateEnergy(const std::vector<EvalWell>& cq_s,
2633 const IntensiveQuantities& intQuants,
2634 DeferredLogger& deferred_logger) const
2635 {
2636 auto fs = intQuants.fluidState();
2637 Eval result = 0;
2638 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
2639 if (!FluidSystem::phaseIsActive(phaseIdx)) {
2640 continue;
2641 }
2642
2643 // convert to reservoir conditions
2644 EvalWell cq_r_thermal{0.};
2645 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
2646 const bool both_oil_gas = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
2647 if (!both_oil_gas || FluidSystem::waterPhaseIdx == phaseIdx) {
2648 cq_r_thermal = cq_s[activeCompIdx] / this->extendEval(fs.invB(phaseIdx));
2649 } else {
2650 // remove dissolved gas and vapporized oil
2651 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
2652 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
2653 // q_os = q_or * b_o + rv * q_gr * b_g
2654 // q_gs = q_gr * g_g + rs * q_or * b_o
2655 // q_gr = 1 / (b_g * d) * (q_gs - rs * q_os)
2656 // d = 1.0 - rs * rv
2657 const EvalWell d = this->extendEval(1.0 - fs.Rv() * fs.Rs());
2658 if (d <= 0.0) {
2659 deferred_logger.debug(
2660 fmt::format("Problematic d value {} obtained for well {}"
2661 " during calculateSinglePerf with rs {}"
2662 ", rv {}. Continue as if no dissolution (rs = 0) and"
2663 " vaporization (rv = 0) for this connection.",
2664 d, this->name(), fs.Rs(), fs.Rv()));
2665 cq_r_thermal = cq_s[activeCompIdx] / this->extendEval(fs.invB(phaseIdx));
2666 } else {
2667 if (FluidSystem::gasPhaseIdx == phaseIdx) {
2668 cq_r_thermal = (cq_s[gasCompIdx] -
2669 this->extendEval(fs.Rs()) * cq_s[oilCompIdx]) /
2670 (d * this->extendEval(fs.invB(phaseIdx)) );
2671 } else if (FluidSystem::oilPhaseIdx == phaseIdx) {
2672 // q_or = 1 / (b_o * d) * (q_os - rv * q_gs)
2673 cq_r_thermal = (cq_s[oilCompIdx] - this->extendEval(fs.Rv()) *
2674 cq_s[gasCompIdx]) /
2675 (d * this->extendEval(fs.invB(phaseIdx)) );
2676 }
2677 }
2678 }
2679
2680 // change temperature for injecting fluids
2681 if (this->isInjector() && !this->wellIsStopped() && cq_r_thermal > 0.0){
2682 // only handles single phase injection now
2683 assert(this->well_ecl_.injectorType() != InjectorType::MULTI);
2684 fs.setTemperature(this->well_ecl_.inj_temperature());
2685 typedef typename std::decay<decltype(fs)>::type::Scalar FsScalar;
2686 typename FluidSystem::template ParameterCache<FsScalar> paramCache;
2687 const unsigned pvtRegionIdx = intQuants.pvtRegionIndex();
2688 paramCache.setRegionIndex(pvtRegionIdx);
2689 paramCache.updatePhase(fs, phaseIdx);
2690
2691 const auto& rho = FluidSystem::density(fs, paramCache, phaseIdx);
2692 fs.setDensity(phaseIdx, rho);
2693 const auto& h = FluidSystem::enthalpy(fs, paramCache, phaseIdx);
2694 fs.setEnthalpy(phaseIdx, h);
2695 cq_r_thermal *= this->extendEval(fs.enthalpy(phaseIdx)) * this->extendEval(fs.density(phaseIdx));
2696 result += getValue(cq_r_thermal);
2697 } else if (cq_r_thermal > 0.0) {
2698 cq_r_thermal *= getValue(fs.enthalpy(phaseIdx)) * getValue(fs.density(phaseIdx));
2699 result += Base::restrictEval(cq_r_thermal);
2700 } else {
2701 // compute the thermal flux
2702 cq_r_thermal *= this->extendEval(fs.enthalpy(phaseIdx)) * this->extendEval(fs.density(phaseIdx));
2703 result += Base::restrictEval(cq_r_thermal);
2704 }
2705 }
2706
2707 return result * this->well_efficiency_factor_;
2708 }
2709} // namespace Opm
2710
2711#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:1406
EvalWell wpolymermw(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2040
typename StdWellEval::EvalWell EvalWell
Definition: StandardWell.hpp:121
void updateWellStateFromPrimaryVariables(WellStateType &well_state, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:775
WellConnectionProps computePropertiesForWellConnectionPressures(const Simulator &simulator, const WellStateType &well_state) const
Definition: StandardWell_impl.hpp:1140
typename StdWellEval::BVectorWell BVectorWell
Definition: StandardWell.hpp:122
std::vector< Scalar > computeWellPotentialWithTHP(const Simulator &ebosSimulator, const WellGroupHelperType &wgHelper, DeferredLogger &deferred_logger, const WellStateType &well_state) const
Definition: StandardWell_impl.hpp:1589
void addWellContributions(SparseMatrixAdapter &mat) const override
Definition: StandardWell_impl.hpp:1937
std::vector< Scalar > getPrimaryVars() const override
Definition: StandardWell_impl.hpp:2602
void addWellPressureEquations(PressureMatrix &mat, const BVector &x, const int pressureVarIndex, const bool use_well_weights, const WellStateType &well_state) const override
Definition: StandardWell_impl.hpp:1945
bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger, const bool fixed_control=false, const bool fixed_status=false) override
Definition: StandardWell_impl.hpp:2436
void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) override
computing the well potentials for group control
Definition: StandardWell_impl.hpp:1761
void updateWaterMobilityWithPolymer(const Simulator &simulator, const int perf, std::vector< EvalWell > &mob_water, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1873
void assembleWellEqWithoutIteration(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:341
std::vector< Scalar > computeCurrentWellRates(const Simulator &ebosSimulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:2566
void calculateSinglePerf(const Simulator &simulator, const int perf, WellStateType &well_state, std::vector< RateVector > &connectionRates, std::vector< EvalWell > &cq_s, EvalWell &water_flux_s, EvalWell &cq_s_zfrac_effective, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:496
bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:2373
void updatePrimaryVariablesNewton(const BVectorWell &dwells, const bool stop_or_zero_rate_target, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:752
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
bool computeWellPotentialsImplicit(const Simulator &ebos_simulator, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1625
void computeWellRatesWithBhpIterations(const Simulator &ebosSimulator, const Scalar &bhp, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1518
void computeWellConnectionPressures(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1366
StandardWell(const Well &well, const ParallelWellInfo< Scalar > &pw_info, const int time_step, const ModelParameters &param, const RateConverterType &rate_converter, const int pvtRegionIdx, const int num_conservation_quantities, const int num_phases, const int index_of_well, const std::vector< PerforationData< Scalar > > &perf_data)
Definition: StandardWell_impl.hpp:52
typename StdWellEval::StdWellConnections::Properties WellConnectionProps
Definition: StandardWell.hpp:269
void updateConnectionRatePolyMW(const EvalWell &cq_s_poly, const IntensiveQuantities &int_quants, const WellStateType &well_state, const int perf, std::vector< RateVector > &connectionRates, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2189
void computeWellRatesWithBhp(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1472
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:674
void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:799
void updateIPRImplicit(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:892
void computeWellConnectionDensitesPressures(const Simulator &simulator, const WellStateType &well_state, const WellConnectionProps &props, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1310
void computePerfRate(const IntensiveQuantities &intQuants, const std::vector< Value > &mob, const Value &bhp, const std::vector< Scalar > &Tw, const int perf, const bool allow_cf, std::vector< Value > &cq_s, PerforationRates< Scalar > &perf_rates, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:93
std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &ebos_simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger, bool iterate_if_no_solution) const override
Definition: StandardWell_impl.hpp:2248
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:730
void handleInjectivityEquations(const Simulator &simulator, const WellStateType &well_state, const int perf, const EvalWell &water_flux_s, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:2123
virtual void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: StandardWell_impl.hpp:1420
virtual ConvergenceReport getWellConvergence(const Simulator &simulator, const WellStateType &well_state, const std::vector< Scalar > &B_avg, DeferredLogger &deferred_logger, const bool relax_tolerance) const override
check whether the well equations get converged for this well
Definition: StandardWell_impl.hpp:1191
void checkConvergenceExtraEqs(const std::vector< Scalar > &res, ConvergenceReport &report) const
Definition: StandardWell_impl.hpp:2170
typename StdWellEval::Eval Eval
Definition: StandardWell.hpp:120
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1129
void computeWellRatesWithThpAlqProd(const Simulator &ebos_simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1742
Scalar computeWellRatesAndBhpWithThpAlqProd(const Simulator &ebos_simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1712
bool allDrawDownWrongDirection(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1087
EvalWell pskin(const Scalar throughput, const EvalWell &water_velocity, const EvalWell &poly_inj_conc, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1996
static constexpr int numWellConservationEq
Definition: StandardWell.hpp:97
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: StandardWell_impl.hpp:2619
void updateWaterThroughput(const double dt, WellStateType &well_state) const override
Definition: StandardWell_impl.hpp:2069
void checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellStateType &well_state, const WellGroupHelperType &wgHelper, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1037
void checkOperabilityUnderBHPLimit(const WellStateType &well_state, const Simulator &simulator, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:968
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2335
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1452
void assembleWellEqWithoutIterationImpl(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:367
EvalWell pskinwater(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1965
void solveEqAndUpdateWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1384
void handleInjectivityRate(const Simulator &simulator, const int perf, std::vector< EvalWell > &cq_s) const
Definition: StandardWell_impl.hpp:2100
virtual void init(const std::vector< Scalar > &depth_arg, const Scalar gravity_arg, const std::vector< Scalar > &B_avg, const bool changed_to_open_this_step) override
Definition: StandardWell_impl.hpp:76
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1235
void updatePrimaryVariables(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1839
std::optional< Scalar > computeBhpAtThpLimitProd(const WellStateType &well_state, const Simulator &simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2231
Scalar getRefDensity() const override
Definition: StandardWell_impl.hpp:1862
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: StandardWell_impl.hpp:1824
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
Definition: WellGroupHelper.hpp:51
const GroupState< Scalar > & groupState() const
Definition: WellGroupHelper.hpp:175
const WellState< Scalar, IndexTraits > & wellState() const
Definition: WellGroupHelper.hpp:256
WellStateGuard pushWellState(WellState< Scalar, IndexTraits > &well_state)
Definition: WellGroupHelper.hpp:190
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:77
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:82
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:105
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:98
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:2049
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:87
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:97
typename Base::ModelParameters ModelParameters
Definition: WellInterface.hpp:111
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:84
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:86
bool solveWellWithOperabilityCheck(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:580
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:89
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