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