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