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