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