22#ifndef OPM_STANDARDWELL_IMPL_HEADER_INCLUDED
23#define OPM_STANDARDWELL_IMPL_HEADER_INCLUDED
26#ifndef OPM_STANDARDWELL_HEADER_INCLUDED
31#include <opm/common/Exceptions.hpp>
33#include <opm/input/eclipse/Units/Units.hpp>
46#include <fmt/format.h>
51 template<
typename TypeTag>
58 const int pvtRegionIdx,
59 const int num_conservation_quantities,
61 const int index_of_well,
63 :
Base(well, pw_info, time_step, param, rate_converter, pvtRegionIdx, num_conservation_quantities, num_phases, index_of_well, perf_data)
74 template<
typename TypeTag>
77 init(
const std::vector<Scalar>& depth_arg,
79 const std::vector< Scalar >& B_avg,
80 const bool changed_to_open_this_step)
82 Base::init(depth_arg, gravity_arg, B_avg, changed_to_open_this_step);
83 this->StdWellEval::init(this->perf_depth_, depth_arg, Base::has_polymermw);
90 template<
typename TypeTag>
95 const std::vector<Value>& mob,
97 const std::vector<Value>& Tw,
100 std::vector<Value>& cq_s,
104 auto obtain = [
this](
const Eval& value)
106 if constexpr (std::is_same_v<Value, Scalar>) {
107 static_cast<void>(
this);
108 return getValue(value);
110 return this->extendEval(value);
113 auto obtainN = [](
const auto& value)
115 if constexpr (std::is_same_v<Value, Scalar>) {
116 return getValue(value);
121 auto zeroElem = [
this]()
123 if constexpr (std::is_same_v<Value, Scalar>) {
124 static_cast<void>(
this);
127 return Value{this->primary_variables_.numWellEq() + Indices::numEq, 0.0};
131 const auto& fs = intQuants.fluidState();
132 const Value pressure = obtain(this->getPerfCellPressure(fs));
133 const Value rs = obtain(fs.Rs());
134 const Value rv = obtain(fs.Rv());
135 const Value rvw = obtain(fs.Rvw());
136 const Value rsw = obtain(fs.Rsw());
138 std::vector<Value> b_perfcells_dense(this->numConservationQuantities(), zeroElem());
139 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
140 if (!FluidSystem::phaseIsActive(phaseIdx)) {
143 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
144 b_perfcells_dense[compIdx] = obtain(fs.invB(phaseIdx));
146 if constexpr (has_solvent) {
147 b_perfcells_dense[Indices::contiSolventEqIdx] = obtain(intQuants.solventInverseFormationVolumeFactor());
150 if constexpr (has_zFraction) {
151 if (this->isInjector()) {
152 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
153 b_perfcells_dense[gasCompIdx] *= (1.0 - this->wsolvent());
154 b_perfcells_dense[gasCompIdx] += this->wsolvent()*intQuants.zPureInvFormationVolumeFactor().value();
158 Value skin_pressure = zeroElem();
160 if (this->isInjector()) {
161 const int pskin_index = Bhp + 1 + this->numLocalPerfs() + perf;
162 skin_pressure = obtainN(this->primary_variables_.eval(pskin_index));
167 std::vector<Value> cmix_s(this->numConservationQuantities(), zeroElem());
168 for (
int componentIdx = 0; componentIdx < this->numConservationQuantities(); ++componentIdx) {
169 cmix_s[componentIdx] = obtainN(this->primary_variables_.surfaceVolumeFraction(componentIdx));
192 template<
typename TypeTag>
193 template<
class Value>
197 const Value& pressure,
203 std::vector<Value>& b_perfcells_dense,
204 const std::vector<Value>& Tw,
207 const Value& skin_pressure,
208 const std::vector<Value>& cmix_s,
209 std::vector<Value>& cq_s,
214 const Value well_pressure = bhp + this->connections_.pressure_diff(perf);
215 Value drawdown = pressure - well_pressure;
216 if (this->isInjector()) {
217 drawdown += skin_pressure;
221 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)
222 ? FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx)
224 FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)
225 ? FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx)
227 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
228 ? FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx)
236 if (!allow_cf && this->isInjector()) {
241 for (
int componentIdx = 0; componentIdx < this->numConservationQuantities(); ++componentIdx) {
242 const Value cq_p = - Tw[componentIdx] * (mob[componentIdx] * drawdown);
243 cq_s[componentIdx] = b_perfcells_dense[componentIdx] * cq_p;
246 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
247 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
249 ratioCalc.gasOilPerfRateProd(cq_s, perf_rates, rv, rs, rvw,
250 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx),
252 }
else if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) &&
253 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
255 ratioCalc.gasWaterPerfRateProd(cq_s, perf_rates, rvw, rsw, this->isProducer());
259 if (!allow_cf && this->isProducer()) {
264 Value total_mob_dense = mob[0];
265 for (
int componentIdx = 1; componentIdx < this->numConservationQuantities(); ++componentIdx) {
266 total_mob_dense += mob[componentIdx];
270 Value volumeRatio = bhp * 0.0;
272 if (FluidSystem::enableVaporizedWater() && FluidSystem::enableDissolvedGasInWater()) {
273 ratioCalc.disOilVapWatVolumeRatio(volumeRatio, rvw, rsw, pressure,
274 cmix_s, b_perfcells_dense, deferred_logger);
278 assert(FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx));
279 assert(FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx));
280 assert(!FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx));
283 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
284 const unsigned waterCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
285 volumeRatio += cmix_s[waterCompIdx] / b_perfcells_dense[waterCompIdx];
288 if constexpr (Indices::enableSolvent) {
289 volumeRatio += cmix_s[Indices::contiSolventEqIdx] / b_perfcells_dense[Indices::contiSolventEqIdx];
292 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
293 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
295 ratioCalc.gasOilVolumeRatio(volumeRatio, rv, rs, pressure,
296 cmix_s, b_perfcells_dense,
299 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
300 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
301 volumeRatio += cmix_s[oilCompIdx] / b_perfcells_dense[oilCompIdx];
303 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
304 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
305 volumeRatio += cmix_s[gasCompIdx] / b_perfcells_dense[gasCompIdx];
311 for (
int componentIdx = 0; componentIdx < this->numConservationQuantities(); ++componentIdx) {
312 const Value cqt_i = - Tw[componentIdx] * (total_mob_dense * drawdown);
313 Value cqt_is = cqt_i / volumeRatio;
314 cq_s[componentIdx] = cmix_s[componentIdx] * cqt_is;
318 if (this->isProducer()) {
319 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
320 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
322 ratioCalc.gasOilPerfRateInj(cq_s, perf_rates,
323 rv, rs, pressure, rvw,
324 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx),
327 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) &&
328 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx))
331 ratioCalc.gasWaterPerfRateInj(cq_s, perf_rates, rvw, rsw,
332 pressure, deferred_logger);
339 template<
typename TypeTag>
345 const Well::InjectionControls& inj_controls,
346 const Well::ProductionControls& prod_controls,
348 const bool solving_with_zero_rate)
352 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
355 this->linSys_.clear();
357 assembleWellEqWithoutIterationImpl(simulator, groupStateHelper, dt, inj_controls,
358 prod_controls, well_state, solving_with_zero_rate);
364 template<
typename TypeTag>
370 const Well::InjectionControls& inj_controls,
371 const Well::ProductionControls& prod_controls,
373 const bool solving_with_zero_rate)
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;
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;
388 const int np = this->number_of_phases_;
390 std::vector<RateVector> connectionRates = this->connectionRates_;
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) {
396 std::vector<EvalWell> cq_s(this->num_conservation_quantities_, 0.0);
399 calculateSinglePerf(simulator, perf, well_state, connectionRates,
400 cq_s, water_flux_s, cq_s_zfrac_effective, deferred_logger);
403 if constexpr (has_polymer && Base::has_polymermw) {
404 if (this->isInjector()) {
405 handleInjectivityEquations(simulator, well_state, perf,
406 water_flux_s, deferred_logger);
409 for (
int componentIdx = 0; componentIdx < this->num_conservation_quantities_; ++componentIdx) {
411 const EvalWell cq_s_effective = cq_s[componentIdx] * this->well_efficiency_factor_;
413 connectionRates[perf][componentIdx] = Base::restrictEval(cq_s_effective);
416 assemblePerforationEq(cq_s_effective,
419 this->primary_variables_.numWellEq(),
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();
427 perf_rates[perf*np + FluidSystem::activeCompToActivePhaseIdx(componentIdx)] = cq_s[componentIdx].value();
431 if constexpr (has_zFraction) {
433 assembleZFracEq(cq_s_zfrac_effective,
435 this->primary_variables_.numWellEq(),
440 this->connectionRates_ = connectionRates;
445 const auto& comm = this->parallel_well_info_.communication();
446 comm.sum(ws.phase_mixing_rates.data(), ws.phase_mixing_rates.size());
450 this->linSys_.sumDistributed(this->parallel_well_info_.communication());
453 for (
int componentIdx = 0; componentIdx < numWellConservationEq; ++componentIdx) {
457 if (FluidSystem::numActivePhases() > 1) {
459 resWell_loc += (this->primary_variables_.surfaceVolumeFraction(componentIdx) -
460 this->F0_[componentIdx]) * volume / dt;
462 resWell_loc -= this->primary_variables_.getQs(componentIdx) * this->well_efficiency_factor_;
464 assembleSourceEq(resWell_loc,
466 this->primary_variables_.numWellEq(),
470 const bool stopped_or_zero_target = this->stoppedOrZeroRateTarget(groupStateHelper);
476 auto& group_state = solving_with_zero_rate
480 auto group_guard = groupStateHelper_copy.
pushGroupState(group_state);
482 assembleControlEq(groupStateHelper_copy,
483 inj_controls, prod_controls,
484 this->primary_variables_,
485 this->getRefDensity(),
487 stopped_or_zero_target);
492 this->linSys_.invert();
494 OPM_DEFLOG_PROBLEM(NumericalProblem,
"Error when inverting local well equations for well " + name(), deferred_logger);
501 template<
typename TypeTag>
507 std::vector<RateVector>& connectionRates,
508 std::vector<EvalWell>& cq_s,
513 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
514 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
515 const int cell_idx = this->well_cells_[perf];
516 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
517 std::vector<EvalWell> mob(this->num_conservation_quantities_, {0.});
518 getMobility(simulator, perf, mob, deferred_logger);
522 getTransMult(trans_mult, simulator, cell_idx);
523 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
524 std::vector<EvalWell> Tw(this->num_conservation_quantities_, this->well_index_[perf] * trans_mult);
525 this->getTw(Tw, perf, intQuants, trans_mult, wellstate_nupcol);
526 computePerfRate(intQuants, mob, bhp, Tw, perf, allow_cf,
527 cq_s, perf_rates, deferred_logger);
529 auto& ws = well_state.
well(this->index_of_well_);
530 auto& perf_data = ws.perf_data;
531 if constexpr (has_polymer && Base::has_polymermw) {
532 if (this->isInjector()) {
535 const unsigned water_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
536 water_flux_s = cq_s[water_comp_idx];
539 handleInjectivityRate(simulator, perf, cq_s);
544 if (this->isProducer()) {
545 ws.phase_mixing_rates[ws.dissolved_gas] += perf_rates.
dis_gas;
546 ws.phase_mixing_rates[ws.dissolved_gas_in_water] += perf_rates.
dis_gas_in_water;
547 ws.phase_mixing_rates[ws.vaporized_oil] += perf_rates.
vap_oil;
548 ws.phase_mixing_rates[ws.vaporized_water] += perf_rates.
vap_wat;
549 perf_data.phase_mixing_rates[perf][ws.dissolved_gas] = perf_rates.
dis_gas;
550 perf_data.phase_mixing_rates[perf][ws.dissolved_gas_in_water] = perf_rates.
dis_gas_in_water;
551 perf_data.phase_mixing_rates[perf][ws.vaporized_oil] = perf_rates.
vap_oil;
552 perf_data.phase_mixing_rates[perf][ws.vaporized_water] = perf_rates.
vap_wat;
555 if constexpr (has_energy) {
556 connectionRates[perf][Indices::contiEnergyEqIdx] =
557 connectionRateEnergy(cq_s, intQuants, deferred_logger);
558 ws.energy_rate += getValue(connectionRates[perf][Indices::contiEnergyEqIdx]);
561 if constexpr (has_polymer) {
562 std::variant<Scalar,EvalWell> polymerConcentration;
563 if (this->isInjector()) {
564 polymerConcentration = this->wpolymer();
566 polymerConcentration = this->extendEval(intQuants.polymerConcentration() *
567 intQuants.polymerViscosityCorrection());
570 [[maybe_unused]]
EvalWell cq_s_poly;
571 std::tie(connectionRates[perf][Indices::contiPolymerEqIdx],
573 this->connections_.connectionRatePolymer(perf_data.polymer_rates[perf],
574 cq_s, polymerConcentration);
576 if constexpr (Base::has_polymermw) {
577 updateConnectionRatePolyMW(cq_s_poly, intQuants, well_state,
578 perf, connectionRates, deferred_logger);
582 if constexpr (has_foam) {
583 std::variant<Scalar,EvalWell> foamConcentration;
584 if (this->isInjector()) {
585 foamConcentration = this->wfoam();
587 foamConcentration = this->extendEval(intQuants.foamConcentration());
589 connectionRates[perf][Indices::contiFoamEqIdx] =
590 this->connections_.connectionRateFoam(cq_s, foamConcentration,
591 FoamModule::transportPhase(),
595 if constexpr (has_zFraction) {
596 std::variant<Scalar,std::array<EvalWell,2>> solventConcentration;
597 if (this->isInjector()) {
598 solventConcentration = this->wsolvent();
600 solventConcentration = std::array{this->extendEval(intQuants.xVolume()),
601 this->extendEval(intQuants.yVolume())};
603 std::tie(connectionRates[perf][Indices::contiZfracEqIdx],
604 cq_s_zfrac_effective) =
605 this->connections_.connectionRatezFraction(perf_data.solvent_rates[perf],
607 solventConcentration);
610 if constexpr (has_brine) {
611 std::variant<Scalar,EvalWell> saltConcentration;
612 if (this->isInjector()) {
613 saltConcentration = this->wsalt();
615 saltConcentration = this->extendEval(intQuants.fluidState().saltConcentration());
618 connectionRates[perf][Indices::contiBrineEqIdx] =
619 this->connections_.connectionRateBrine(perf_data.brine_rates[perf],
624 if constexpr (has_bioeffects) {
625 std::variant<Scalar,EvalWell> microbialConcentration;
626 if constexpr (has_micp) {
627 std::variant<Scalar,EvalWell> oxygenConcentration;
628 std::variant<Scalar,EvalWell> ureaConcentration;
629 if (this->isInjector()) {
630 microbialConcentration = this->wmicrobes();
631 oxygenConcentration = this->woxygen();
632 ureaConcentration = this->wurea();
634 microbialConcentration = this->extendEval(intQuants.microbialConcentration());
635 oxygenConcentration = this->extendEval(intQuants.oxygenConcentration());
636 ureaConcentration = this->extendEval(intQuants.ureaConcentration());
638 std::tie(connectionRates[perf][Indices::contiMicrobialEqIdx],
639 connectionRates[perf][Indices::contiOxygenEqIdx],
640 connectionRates[perf][Indices::contiUreaEqIdx]) =
641 this->connections_.connectionRatesMICP(perf_data.microbial_rates[perf],
642 perf_data.oxygen_rates[perf],
643 perf_data.urea_rates[perf],
645 microbialConcentration,
650 if (this->isProducer()) {
651 microbialConcentration = this->extendEval(intQuants.microbialConcentration());
652 connectionRates[perf][Indices::contiMicrobialEqIdx] =
653 this->connections_.connectionRateBioeffects(perf_data.microbial_rates[perf],
655 microbialConcentration);
661 perf_data.pressure[perf] = ws.bhp + this->connections_.pressure_diff(perf);
664 if (FluidSystem::phaseUsage().hasCO2orH2Store()) {
665 const unsigned gas_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
666 const Scalar rho = FluidSystem::referenceDensity( FluidSystem::gasPhaseIdx, Base::pvtRegionIdx() );
667 perf_data.gas_mass_rates[perf] = cq_s[gas_comp_idx].value() * rho;
671 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
672 const unsigned wat_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
673 const Scalar rho = FluidSystem::referenceDensity( FluidSystem::waterPhaseIdx, Base::pvtRegionIdx() );
674 perf_data.wat_mass_rates[perf] = cq_s[wat_comp_idx].value() * rho;
678 template<
typename TypeTag>
679 template<
class Value>
684 const int cell_idx)
const
686 auto obtain = [
this](
const Eval& value)
688 if constexpr (std::is_same_v<Value, Scalar>) {
689 static_cast<void>(
this);
690 return getValue(value);
692 return this->extendEval(value);
698 template<
typename TypeTag>
699 template<
class Value>
704 std::vector<Value>& mob,
707 auto obtain = [
this](
const Eval& value)
709 if constexpr (std::is_same_v<Value, Scalar>) {
710 static_cast<void>(
this);
711 return getValue(value);
713 return this->extendEval(value);
717 obtain, deferred_logger);
720 if constexpr (has_polymer) {
721 if (!FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
722 OPM_DEFLOG_THROW(std::runtime_error,
"Water is required when polymer is active", deferred_logger);
727 if constexpr (!Base::has_polymermw) {
728 if constexpr (std::is_same_v<Value, Scalar>) {
729 std::vector<EvalWell> mob_eval(this->num_conservation_quantities_, 0.);
730 for (std::size_t i = 0; i < mob.size(); ++i) {
731 mob_eval[i].setValue(mob[i]);
733 updateWaterMobilityWithPolymer(simulator, perf, mob_eval, deferred_logger);
734 for (std::size_t i = 0; i < mob.size(); ++i) {
735 mob[i] = getValue(mob_eval[i]);
738 updateWaterMobilityWithPolymer(simulator, perf, mob, deferred_logger);
745 const Scalar bhp = this->primary_variables_.value(Bhp);
746 const Scalar perf_press = bhp + this->connections_.pressure_diff(perf);
747 const Scalar multiplier = this->getInjMult(perf, bhp, perf_press, deferred_logger);
748 for (std::size_t i = 0; i < mob.size(); ++i) {
749 mob[i] *= multiplier;
755 template<
typename TypeTag>
763 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
767 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(groupStateHelper);
768 updatePrimaryVariablesNewton(dwells, stop_or_zero_rate_target, deferred_logger);
770 const auto& summary_state = simulator.vanguard().summaryState();
771 updateWellStateFromPrimaryVariables(well_state, summary_state, deferred_logger);
775 const bool isThermal = simulator.vanguard().eclState().getSimulationConfig().isThermal();
776 const bool co2store = simulator.vanguard().eclState().runspec().co2Storage();
777 Base::calculateReservoirRates( (isThermal || co2store), well_state.
well(this->index_of_well_));
784 template<
typename TypeTag>
788 const bool stop_or_zero_rate_target,
791 const Scalar dFLimit = this->param_.dwell_fraction_max_;
792 const Scalar dBHPLimit = this->param_.dbhp_max_rel_;
793 this->primary_variables_.updateNewton(dwells, stop_or_zero_rate_target, dFLimit, dBHPLimit, deferred_logger);
796 if constexpr (Base::has_polymermw) {
797 this->primary_variables_.updateNewtonPolyMW(dwells);
800 this->primary_variables_.checkFinite(deferred_logger,
"Newton update");
807 template<
typename TypeTag>
811 const SummaryState& summary_state,
814 this->primary_variables_.copyToWellState(well_state, deferred_logger);
817 updateThp(getRefDensity(),
818 [
this,&well_state]() {
return this->baseif_.getALQ(well_state); },
819 well_state, summary_state, deferred_logger);
822 if constexpr (Base::has_polymermw) {
823 this->primary_variables_.copyToWellStatePolyMW(well_state);
831 template<
typename TypeTag>
839 std::ranges::fill(this->ipr_a_, 0.0);
840 std::ranges::fill(this->ipr_b_, 0.0);
842 for (
int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
843 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
844 getMobility(simulator, perf, mob, deferred_logger);
846 const int cell_idx = this->well_cells_[perf];
847 const auto& int_quantities = simulator.model().intensiveQuantities(cell_idx, 0);
848 const auto& fs = int_quantities.fluidState();
850 Scalar p_r = this->getPerfCellPressure(fs).value();
853 std::vector<Scalar> b_perf(this->num_conservation_quantities_);
854 for (std::size_t phase = 0; phase < FluidSystem::numPhases; ++phase) {
855 if (!FluidSystem::phaseIsActive(phase)) {
858 const unsigned comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phase));
859 b_perf[comp_idx] = fs.invB(phase).value();
861 if constexpr (has_solvent) {
862 b_perf[Indices::contiSolventEqIdx] = int_quantities.solventInverseFormationVolumeFactor().value();
866 const Scalar h_perf = this->connections_.pressure_diff(perf);
867 const Scalar pressure_diff = p_r - h_perf;
872 if ( (this->isProducer() && pressure_diff < 0.) || (this->isInjector() && pressure_diff > 0.) ) {
873 deferred_logger.
debug(
"CROSSFLOW_IPR",
874 "cross flow found when updateIPR for well " + name()
875 +
" . The connection is ignored in IPR calculations");
882 getTransMult(trans_mult, simulator, cell_idx);
883 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
884 std::vector<Scalar> tw_perf(this->num_conservation_quantities_, this->well_index_[perf] * trans_mult);
885 this->getTw(tw_perf, perf, int_quantities, trans_mult, wellstate_nupcol);
886 std::vector<Scalar> ipr_a_perf(this->ipr_a_.size());
887 std::vector<Scalar> ipr_b_perf(this->ipr_b_.size());
888 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
889 const Scalar tw_mob = tw_perf[comp_idx] * mob[comp_idx] * b_perf[comp_idx];
890 ipr_a_perf[comp_idx] += tw_mob * pressure_diff;
891 ipr_b_perf[comp_idx] += tw_mob;
895 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
896 const unsigned oil_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
897 const unsigned gas_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
898 const Scalar rs = (fs.Rs()).value();
899 const Scalar rv = (fs.Rv()).value();
901 const Scalar dis_gas_a = rs * ipr_a_perf[oil_comp_idx];
902 const Scalar vap_oil_a = rv * ipr_a_perf[gas_comp_idx];
904 ipr_a_perf[gas_comp_idx] += dis_gas_a;
905 ipr_a_perf[oil_comp_idx] += vap_oil_a;
907 const Scalar dis_gas_b = rs * ipr_b_perf[oil_comp_idx];
908 const Scalar vap_oil_b = rv * ipr_b_perf[gas_comp_idx];
910 ipr_b_perf[gas_comp_idx] += dis_gas_b;
911 ipr_b_perf[oil_comp_idx] += vap_oil_b;
914 for (std::size_t comp_idx = 0; comp_idx < ipr_a_perf.size(); ++comp_idx) {
915 this->ipr_a_[comp_idx] += ipr_a_perf[comp_idx];
916 this->ipr_b_[comp_idx] += ipr_b_perf[comp_idx];
919 this->parallel_well_info_.communication().sum(this->ipr_a_.data(), this->ipr_a_.size());
920 this->parallel_well_info_.communication().sum(this->ipr_b_.data(), this->ipr_b_.size());
923 template<
typename TypeTag>
938 auto rates = well_state.
well(this->index_of_well_).surface_rates;
940 for (std::size_t p = 0; p < rates.size(); ++p) {
941 zero_rates &= rates[p] == 0.0;
943 auto& ws = well_state.
well(this->index_of_well_);
945 const auto msg = fmt::format(
"updateIPRImplicit: Well {} has zero rate, IPRs might be problematic", this->name());
946 deferred_logger.debug(msg);
959 std::ranges::fill(ws.implicit_ipr_a, 0.0);
960 std::ranges::fill(ws.implicit_ipr_b, 0.0);
962 auto inj_controls = Well::InjectionControls(0);
963 auto prod_controls = Well::ProductionControls(0);
964 prod_controls.addControl(Well::ProducerCMode::BHP);
965 prod_controls.bhp_limit = well_state.
well(this->index_of_well_).bhp;
968 const auto cmode = ws.production_cmode;
969 ws.production_cmode = Well::ProducerCMode::BHP;
970 const double dt = simulator.timeStepSize();
971 assembleWellEqWithoutIteration(simulator, groupStateHelper, dt, inj_controls, prod_controls, well_state,
974 const size_t nEq = this->primary_variables_.numWellEq();
978 for (
size_t i=0; i < nEq; ++i){
984 x_well[0].resize(nEq);
985 this->linSys_.solve(rhs, x_well);
987 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx){
988 EvalWell comp_rate = this->primary_variables_.getQs(comp_idx);
989 const int idx = FluidSystem::activeCompToActivePhaseIdx(comp_idx);
990 for (
size_t pvIdx = 0; pvIdx < nEq; ++pvIdx) {
992 ws.implicit_ipr_b[idx] -= x_well[0][pvIdx]*comp_rate.derivative(pvIdx+Indices::numEq);
994 ws.implicit_ipr_a[idx] = ws.implicit_ipr_b[idx]*ws.bhp - comp_rate.value();
997 ws.production_cmode = cmode;
1000 template<
typename TypeTag>
1007 const auto& summaryState = simulator.vanguard().summaryState();
1011 const bool bhp_limit_not_defaulted = bhp_limit > 1.5 * unit::barsa;
1012 if ( bhp_limit_not_defaulted || !this->wellHasTHPConstraints(summaryState) ) {
1015 Scalar total_ipr_mass_rate = 0.0;
1016 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
1018 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1022 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
1023 const Scalar ipr_rate = this->ipr_a_[compIdx] - this->ipr_b_[compIdx] * bhp_limit;
1025 const Scalar rho = FluidSystem::referenceDensity( phaseIdx, Base::pvtRegionIdx() );
1026 total_ipr_mass_rate += ipr_rate * rho;
1028 if ( (this->isProducer() && total_ipr_mass_rate < 0.) || (this->isInjector() && total_ipr_mass_rate > 0.) ) {
1029 this->operability_status_.operable_under_only_bhp_limit =
false;
1033 if (this->operability_status_.operable_under_only_bhp_limit && this->wellHasTHPConstraints(summaryState)) {
1037 std::vector<Scalar> well_rates_bhp_limit;
1038 computeWellRatesWithBhp(simulator, bhp_limit, well_rates_bhp_limit, deferred_logger);
1040 this->adaptRatesForVFP(well_rates_bhp_limit);
1041 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1044 this->getRefDensity(),
1045 this->getALQ(well_state),
1048 if ( (this->isProducer() && thp < thp_limit) || (this->isInjector() && thp > thp_limit) ) {
1049 this->operability_status_.obey_thp_limit_under_bhp_limit =
false;
1060 this->operability_status_.operable_under_only_bhp_limit =
true;
1061 this->operability_status_.obey_thp_limit_under_bhp_limit =
false;
1069 template<
typename TypeTag>
1077 const auto& summaryState = simulator.vanguard().summaryState();
1078 const auto obtain_bhp = this->isProducer() ? computeBhpAtThpLimitProd(well_state, simulator, groupStateHelper, summaryState)
1079 : computeBhpAtThpLimitInj(simulator, groupStateHelper, summaryState);
1082 this->operability_status_.can_obtain_bhp_with_thp_limit =
true;
1085 this->operability_status_.obey_bhp_limit_with_thp_limit = this->isProducer() ?
1086 *obtain_bhp >= bhp_limit : *obtain_bhp <= bhp_limit ;
1088 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1089 if (this->isProducer() && *obtain_bhp < thp_limit) {
1090 const std::string msg =
" obtained bhp " +
std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1091 +
" bars is SMALLER than thp limit "
1093 +
" bars as a producer for well " + name();
1094 deferred_logger.debug(msg);
1096 else if (this->isInjector() && *obtain_bhp > thp_limit) {
1097 const std::string msg =
" obtained bhp " +
std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1098 +
" bars is LARGER than thp limit "
1100 +
" bars as a injector for well " + name();
1101 deferred_logger.debug(msg);
1104 this->operability_status_.can_obtain_bhp_with_thp_limit =
false;
1105 this->operability_status_.obey_bhp_limit_with_thp_limit =
false;
1106 if (!this->wellIsStopped()) {
1107 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1108 deferred_logger.debug(
" could not find bhp value at thp limit "
1110 +
" bar for well " + name() +
", the well might need to be closed ");
1119 template<
typename TypeTag>
1124 bool all_drawdown_wrong_direction =
true;
1126 for (
int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
1127 const int cell_idx = this->well_cells_[perf];
1128 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
1129 const auto& fs = intQuants.fluidState();
1131 const Scalar pressure = this->getPerfCellPressure(fs).value();
1132 const Scalar bhp = this->primary_variables_.eval(Bhp).value();
1135 const Scalar well_pressure = bhp + this->connections_.pressure_diff(perf);
1136 const Scalar drawdown = pressure - well_pressure;
1141 if ( (drawdown < 0. && this->isInjector()) ||
1142 (drawdown > 0. && this->isProducer()) ) {
1143 all_drawdown_wrong_direction =
false;
1148 const auto& comm = this->parallel_well_info_.communication();
1149 if (comm.size() > 1)
1151 all_drawdown_wrong_direction =
1152 (comm.min(all_drawdown_wrong_direction ? 1 : 0) == 1);
1155 return all_drawdown_wrong_direction;
1161 template<
typename TypeTag>
1166 return !this->getAllowCrossFlow() && allDrawDownWrongDirection(simulator);
1172 template<
typename TypeTag>
1178 auto prop_func =
typename StdWellEval::StdWellConnections::PressurePropertyFunctions {
1180 [&model = simulator.model()](
int cell_idx,
int phase_idx)
1182 return model.intensiveQuantities(cell_idx, 0)
1183 .fluidState().temperature(phase_idx).value();
1187 [&model = simulator.model()](
int cell_idx)
1189 return model.intensiveQuantities(cell_idx, 0)
1190 .fluidState().saltConcentration().value();
1194 [&model = simulator.model()](
int cell_idx)
1196 return model.intensiveQuantities(cell_idx, 0)
1197 .fluidState().pvtRegionIndex();
1201 if constexpr (Indices::enableSolvent) {
1202 prop_func.solventInverseFormationVolumeFactor =
1203 [&model = simulator.model()](
int cell_idx)
1205 return model.intensiveQuantities(cell_idx, 0)
1206 .solventInverseFormationVolumeFactor().value();
1209 prop_func.solventRefDensity = [&model = simulator.model()](
int cell_idx)
1211 return model.intensiveQuantities(cell_idx, 0)
1212 .solventRefDensity();
1216 return this->connections_.computePropertiesForPressures(well_state, prop_func);
1223 template<
typename TypeTag>
1227 const std::vector<Scalar>& B_avg,
1228 const bool relax_tolerance)
const
1232 assert((
int(B_avg.size()) == this->num_conservation_quantities_) || has_polymer || has_energy || has_foam || has_brine || has_zFraction || has_bioeffects);
1235 Scalar tol_wells = this->param_.tolerance_wells_;
1237 constexpr Scalar stopped_factor = 1.e-4;
1239 constexpr Scalar dynamic_thp_factor = 1.e-1;
1240 if (this->stoppedOrZeroRateTarget(groupStateHelper)) {
1241 tol_wells = tol_wells*stopped_factor;
1242 }
else if (this->getDynamicThpLimit()) {
1243 tol_wells = tol_wells*dynamic_thp_factor;
1246 std::vector<Scalar> res;
1249 this->param_.max_residual_allowed_,
1251 this->param_.relaxed_tolerance_flow_well_,
1253 this->wellIsStopped(),
1257 checkConvergenceExtraEqs(res, report);
1266 template<
typename TypeTag>
1274 auto fluidState = [&simulator,
this](
const int perf)
1276 const auto cell_idx = this->well_cells_[perf];
1277 return simulator.model()
1278 .intensiveQuantities(cell_idx, 0).fluidState();
1281 const int np = this->number_of_phases_;
1282 auto setToZero = [np](
Scalar* x) ->
void
1284 std::fill_n(x, np, 0.0);
1287 auto addVector = [np](
const Scalar* src,
Scalar* dest) ->
void
1289 std::transform(src, src + np, dest, dest, std::plus<>{});
1292 auto& ws = well_state.
well(this->index_of_well_);
1293 auto& perf_data = ws.perf_data;
1294 auto* wellPI = ws.productivity_index.data();
1295 auto* connPI = perf_data.prod_index.data();
1299 const auto preferred_phase = this->well_ecl_.getPreferredPhase();
1300 auto subsetPerfID = 0;
1302 for (
const auto& perf : *this->perf_data_) {
1303 auto allPerfID = perf.ecl_index;
1305 auto connPICalc = [&wellPICalc, allPerfID](
const Scalar mobility) ->
Scalar
1310 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
1311 getMobility(simulator,
static_cast<int>(subsetPerfID), mob, deferred_logger);
1313 const auto& fs = fluidState(subsetPerfID);
1316 if (this->isInjector()) {
1317 this->computeConnLevelInjInd(fs, preferred_phase, connPICalc,
1318 mob, connPI, deferred_logger);
1321 this->computeConnLevelProdInd(fs, connPICalc, mob, connPI);
1324 addVector(connPI, wellPI);
1331 const auto& comm = this->parallel_well_info_.communication();
1332 if (comm.size() > 1) {
1333 comm.sum(wellPI, np);
1336 assert ((
static_cast<int>(subsetPerfID) == this->number_of_local_perforations_) &&
1337 "Internal logic error in processing connections for PI/II");
1342 template<
typename TypeTag>
1349 const auto& well_state = groupStateHelper.
wellState();
1353 const auto prop_func =
typename StdWellEval::StdWellConnections::DensityPropertyFunctions {
1358 [&model = simulator.model()](
const int cell,
1359 const std::vector<int>& phases,
1360 std::vector<Scalar>& mob)
1362 const auto& iq = model.intensiveQuantities(cell, 0);
1364 std::ranges::transform(phases, mob.begin(),
1365 [&iq](
const int phase) { return iq.mobility(phase).value(); });
1370 [&model = simulator.model()](
const int cell,
1371 const std::vector<int>& phases,
1372 std::vector<Scalar>& rho)
1374 const auto& fs = model.intensiveQuantities(cell, 0).fluidState();
1376 std::ranges::transform(phases, rho.begin(),
1377 [&fs](
const int phase) { return fs.density(phase).value(); });
1381 const auto stopped_or_zero_rate_target = this->
1382 stoppedOrZeroRateTarget(groupStateHelper);
1385 .computeProperties(stopped_or_zero_rate_target, well_state,
1386 prop_func, props, deferred_logger);
1388 cachedRefDensity = this->connections_.rho(0);
1389 if (this->parallel_well_info_.communication().size() > 1) {
1390 cachedRefDensity = this->parallel_well_info_.broadcastFirstPerforationValue(cachedRefDensity);
1398 template<
typename TypeTag>
1404 const auto& well_state = groupStateHelper.
wellState();
1405 const auto props = computePropertiesForWellConnectionPressures
1406 (simulator, well_state);
1408 computeWellConnectionDensitesPressures(simulator, groupStateHelper, props);
1415 template<
typename TypeTag>
1422 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1427 dx_well[0].resize(this->primary_variables_.numWellEq());
1428 this->linSys_.solve( dx_well);
1430 updateWellState(simulator, dx_well, groupStateHelper, well_state);
1437 template<
typename TypeTag>
1443 updatePrimaryVariables(groupStateHelper);
1444 computeWellConnectionPressures(simulator, groupStateHelper);
1445 this->computeAccumWell();
1450 template<
typename TypeTag>
1455 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1457 if (this->param_.matrix_add_well_contributions_)
1463 this->linSys_.apply(x, Ax);
1469 template<
typename TypeTag>
1474 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1476 this->linSys_.apply(r);
1482 template<
typename TypeTag>
1490 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1493 xw[0].resize(this->primary_variables_.numWellEq());
1495 this->linSys_.recoverSolutionWell(x, xw);
1496 updateWellState(simulator, xw, groupStateHelper, well_state);
1502 template<
typename TypeTag>
1507 std::vector<Scalar>& well_flux,
1511 const int np = this->number_of_phases_;
1512 well_flux.resize(np, 0.0);
1514 const bool allow_cf = this->getAllowCrossFlow();
1516 for (
int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
1517 const int cell_idx = this->well_cells_[perf];
1518 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
1520 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.);
1521 getMobility(simulator, perf, mob, deferred_logger);
1523 getTransMult(trans_mult, simulator, cell_idx);
1524 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1525 std::vector<Scalar> Tw(this->num_conservation_quantities_, this->well_index_[perf] * trans_mult);
1526 this->getTw(Tw, perf, intQuants, trans_mult, wellstate_nupcol);
1528 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.);
1530 computePerfRate(intQuants, mob, bhp, Tw, perf, allow_cf,
1531 cq_s, perf_rates, deferred_logger);
1533 for(
int p = 0; p < np; ++p) {
1534 well_flux[FluidSystem::activeCompToActivePhaseIdx(p)] += cq_s[p];
1538 if constexpr (has_solvent) {
1539 assert(FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx));
1541 const int gas_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
1542 well_flux[gas_pos] += cq_s[Indices::contiSolventEqIdx];
1545 this->parallel_well_info_.communication().sum(well_flux.data(), well_flux.size());
1550 template<
typename TypeTag>
1556 std::vector<Scalar>& well_flux)
const
1572 auto guard = groupStateHelper_copy.
pushWellState(well_state_copy);
1575 const auto& summary_state = simulator.vanguard().summaryState();
1576 auto inj_controls = well_copy.
well_ecl_.isInjector()
1577 ? well_copy.
well_ecl_.injectionControls(summary_state)
1578 : Well::InjectionControls(0);
1579 auto prod_controls = well_copy.
well_ecl_.isProducer()
1580 ? well_copy.
well_ecl_.productionControls(summary_state) :
1581 Well::ProductionControls(0);
1584 auto& ws = well_state_copy.
well(this->index_of_well_);
1586 inj_controls.bhp_limit = bhp;
1587 ws.injection_cmode = Well::InjectorCMode::BHP;
1589 prod_controls.bhp_limit = bhp;
1590 ws.production_cmode = Well::ProducerCMode::BHP;
1595 const int np = this->number_of_phases_;
1596 const Scalar sign = this->well_ecl_.isInjector() ? 1.0 : -1.0;
1597 for (
int phase = 0; phase < np; ++phase){
1598 well_state_copy.
wellRates(this->index_of_well_)[phase]
1599 = sign * ws.well_potentials[phase];
1604 const double dt = simulator.timeStepSize();
1606 simulator, dt, inj_controls, prod_controls, groupStateHelper_copy, well_state_copy
1609 const std::string msg =
" well " + name() +
" did not get converged during well potential calculations "
1610 " potentials are computed based on unconverged solution";
1611 deferred_logger.debug(msg);
1621 template<
typename TypeTag>
1622 std::vector<typename StandardWell<TypeTag>::Scalar>
1629 std::vector<Scalar> potentials(this->number_of_phases_, 0.0);
1630 const auto& summary_state = simulator.vanguard().summaryState();
1632 const auto& well = this->well_ecl_;
1633 if (well.isInjector()){
1634 const auto& controls = this->well_ecl_.injectionControls(summary_state);
1635 auto bhp_at_thp_limit = computeBhpAtThpLimitInj(simulator, groupStateHelper, summary_state);
1636 if (bhp_at_thp_limit) {
1637 const Scalar bhp = std::min(*bhp_at_thp_limit,
1638 static_cast<Scalar>(controls.bhp_limit));
1639 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1641 deferred_logger.warning(
"FAILURE_GETTING_CONVERGED_POTENTIAL",
1642 "Failed in getting converged thp based potential calculation for well "
1643 + name() +
". Instead the bhp based value is used");
1644 const Scalar bhp = controls.bhp_limit;
1645 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1648 computeWellRatesWithThpAlqProd(
1649 simulator, groupStateHelper, summary_state,
1650 potentials, this->getALQ(well_state)
1657 template<
typename TypeTag>
1662 std::vector<Scalar>& well_potentials)
const
1675 auto guard = groupStateHelper_copy.
pushWellState(well_state_copy);
1676 auto& ws = well_state_copy.
well(this->index_of_well_);
1679 const auto& summary_state = simulator.vanguard().summaryState();
1680 auto inj_controls = well_copy.
well_ecl_.isInjector()
1681 ? well_copy.
well_ecl_.injectionControls(summary_state)
1682 : Well::InjectionControls(0);
1683 auto prod_controls = well_copy.
well_ecl_.isProducer()
1684 ? well_copy.
well_ecl_.productionControls(summary_state) :
1685 Well::ProductionControls(0);
1691 const int num_perf = ws.perf_data.size();
1692 for (
int perf = 0; perf < num_perf; ++perf) {
1696 const int np = this->number_of_phases_;
1697 bool trivial =
true;
1698 for (
int phase = 0; phase < np; ++phase){
1699 trivial = trivial && (ws.well_potentials[phase] == 0.0) ;
1703 for (
int phase = 0; phase < np; ++phase) {
1704 ws.surface_rates[phase] = sign * ws.well_potentials[phase];
1709 const double dt = simulator.timeStepSize();
1711 bool converged =
false;
1712 if (this->well_ecl_.isProducer()) {
1714 simulator, dt, inj_controls, prod_controls, groupStateHelper_copy, well_state_copy
1718 simulator, dt, inj_controls, prod_controls, groupStateHelper_copy, well_state_copy,
1726 well_potentials.clear();
1727 well_potentials.resize(np, 0.0);
1728 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1729 if (has_solvent && comp_idx == Indices::contiSolventEqIdx)
continue;
1731 well_potentials[FluidSystem::activeCompToActivePhaseIdx(comp_idx)] = rate.value();
1735 if constexpr (has_solvent) {
1736 assert(FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx));
1738 const int gas_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
1740 well_potentials[gas_pos] += rate.value();
1746 template<
typename TypeTag>
1751 const SummaryState &summary_state,
1752 std::vector<Scalar>& potentials,
1757 auto bhp_at_thp_limit = computeBhpAtThpLimitProdWithAlq(
1758 simulator, groupStateHelper, summary_state, alq,
true);
1759 if (bhp_at_thp_limit) {
1760 const auto& controls = this->well_ecl_.productionControls(summary_state);
1761 bhp = std::max(*bhp_at_thp_limit,
1762 static_cast<Scalar>(controls.bhp_limit));
1763 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1766 deferred_logger.warning(
"FAILURE_GETTING_CONVERGED_POTENTIAL",
1767 "Failed in getting converged thp based potential calculation for well "
1768 + name() +
". Instead the bhp based value is used");
1769 const auto& controls = this->well_ecl_.productionControls(summary_state);
1770 bhp = controls.bhp_limit;
1771 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1776 template<
typename TypeTag>
1781 const SummaryState& summary_state,
1782 std::vector<Scalar>& potentials,
1786 computeWellRatesAndBhpWithThpAlqProd(simulator,
1793 template<
typename TypeTag>
1799 std::vector<Scalar>& well_potentials)
1802 const auto [compute_potential, bhp_controlled_well] =
1805 if (!compute_potential) {
1809 bool converged_implicit =
false;
1813 if (this->param_.local_well_solver_control_switching_ && !(this->changed_to_open_this_step_ && this->wellUnderZeroRateTarget(groupStateHelper))) {
1814 converged_implicit = computeWellPotentialsImplicit(
1815 simulator, groupStateHelper, well_potentials
1818 if (!converged_implicit) {
1820 const auto& summaryState = simulator.vanguard().summaryState();
1821 if (!Base::wellHasTHPConstraints(summaryState) || bhp_controlled_well) {
1831 const auto& ws = well_state.
well(this->index_of_well_);
1832 if (this->isInjector())
1833 bhp = std::max(ws.bhp, bhp);
1835 bhp = std::min(ws.bhp, bhp);
1837 assert(std::abs(bhp) != std::numeric_limits<Scalar>::max());
1838 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, well_potentials);
1841 well_potentials = computeWellPotentialWithTHP(simulator, groupStateHelper, well_state);
1845 this->checkNegativeWellPotentials(well_potentials,
1846 this->param_.check_well_operability_,
1856 template<
typename TypeTag>
1860 const int openConnIdx)
const
1862 return (openConnIdx < 0)
1864 : this->connections_.rho(openConnIdx);
1871 template<
typename TypeTag>
1876 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1879 const auto& well_state = groupStateHelper.
wellState();
1880 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(groupStateHelper);
1881 this->primary_variables_.update(well_state, stop_or_zero_rate_target, deferred_logger);
1884 if constexpr (Base::has_polymermw) {
1885 this->primary_variables_.updatePolyMW(well_state);
1888 this->primary_variables_.checkFinite(deferred_logger,
"updating from well state");
1894 template<
typename TypeTag>
1899 return cachedRefDensity;
1905 template<
typename TypeTag>
1910 std::vector<EvalWell>& mob,
1913 const int cell_idx = this->well_cells_[perf];
1914 const auto& int_quant = simulator.model().intensiveQuantities(cell_idx, 0);
1915 const EvalWell polymer_concentration = this->extendEval(int_quant.polymerConcentration());
1919 if (this->isInjector()) {
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,
true) );
1926 if (PolymerModule::hasPlyshlog()) {
1929 if (this->isInjector() && this->wpolymer() == 0.) {
1934 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
1935 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
1937 std::vector<EvalWell> cq_s(this->num_conservation_quantities_, 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);
1947 const Scalar area = 2 * std::numbers::pi_v<Scalar> * this->perf_rep_radius_[perf] * this->perf_length_[perf];
1948 const auto& material_law_manager = simulator.problem().materialLawManager();
1949 const auto& scaled_drainage_info =
1950 material_law_manager->oilWaterScaledEpsInfoDrainage(cell_idx);
1951 const Scalar swcr = scaled_drainage_info.Swcr;
1952 const EvalWell poro = this->extendEval(int_quant.porosity());
1953 const EvalWell sw = this->extendEval(int_quant.fluidState().saturation(FluidSystem::waterPhaseIdx));
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));
1959 if (PolymerModule::hasShrate()) {
1962 water_velocity *= PolymerModule::shrate( int_quant.pvtRegionIndex() ) / this->bore_diameters_[perf];
1964 const EvalWell shear_factor = PolymerModule::computeShearFactor(polymer_concentration,
1965 int_quant.pvtRegionIndex(),
1968 mob[waterCompIdx] /= shear_factor;
1972 template<
typename TypeTag>
1976 this->linSys_.extract(jacobian);
1980 template <
typename TypeTag>
1984 const int pressureVarIndex,
1985 const bool use_well_weights,
1988 this->linSys_.extractCPRPressureMatrix(jacobian,
1999 template<
typename TypeTag>
2006 if constexpr (Base::has_polymermw) {
2007 const int water_table_id = this->polymerWaterTable_();
2008 if (water_table_id <= 0) {
2010 fmt::format(
"Unused SKPRWAT table id used for well {}", name()),
2013 const auto& water_table_func = PolymerModule::getSkprwatTable(water_table_id);
2014 const EvalWell throughput_eval{throughput};
2016 EvalWell pskin_water = water_table_func.eval(throughput_eval, water_velocity);
2020 fmt::format(
"Polymermw is not activated, while injecting "
2021 "skin pressure is requested for well {}", name()),
2030 template<
typename TypeTag>
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);
2044 const int polymer_table_id = this->polymerTable_();
2045 if (polymer_table_id <= 0) {
2047 fmt::format(
"Unavailable SKPRPOLY table id used for well {}", name()),
2050 const auto& skprpolytable = PolymerModule::getSkprpolyTable(polymer_table_id);
2051 const Scalar reference_concentration = skprpolytable.refConcentration;
2052 const EvalWell throughput_eval{throughput};
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;
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;
2064 fmt::format(
"Polymermw is not activated, while injecting "
2065 "skin pressure is requested for well {}", name()),
2074 template<
typename TypeTag>
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};
2086 if (this->wpolymer() == 0.) {
2087 return molecular_weight;
2089 molecular_weight = table_func.eval(throughput_eval, abs(water_velocity));
2090 return molecular_weight;
2093 fmt::format(
"Polymermw is not activated, while injecting "
2094 "polymer molecular weight is requested for well {}", name()),
2103 template<
typename TypeTag>
2109 if constexpr (Base::has_polymermw) {
2110 if (!this->isInjector()) {
2114 auto& perf_water_throughput = well_state.
well(this->index_of_well_)
2115 .perf_data.water_throughput;
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);
2123 if (perf_water_vel >
Scalar{0}) {
2124 perf_water_throughput[perf] += perf_water_vel * dt;
2134 template<
typename TypeTag>
2139 std::vector<EvalWell>& cq_s)
const
2141 const int cell_idx = this->well_cells_[perf];
2142 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2143 const auto& fs = int_quants.fluidState();
2144 const EvalWell b_w = this->extendEval(fs.invB(FluidSystem::waterPhaseIdx));
2145 const Scalar area = std::numbers::pi_v<Scalar> * this->bore_diameters_[perf] * this->perf_length_[perf];
2146 const int wat_vel_index = Bhp + 1 + perf;
2147 const unsigned water_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
2151 cq_s[water_comp_idx] = area * this->primary_variables_.eval(wat_vel_index) * b_w;
2157 template<
typename TypeTag>
2166 const int cell_idx = this->well_cells_[perf];
2167 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2168 const auto& fs = int_quants.fluidState();
2169 const EvalWell b_w = this->extendEval(fs.invB(FluidSystem::waterPhaseIdx));
2170 const EvalWell water_flux_r = water_flux_s / b_w;
2171 const Scalar area = std::numbers::pi_v<Scalar> * this->bore_diameters_[perf] * this->perf_length_[perf];
2172 const EvalWell water_velocity = water_flux_r / area;
2173 const int wat_vel_index = Bhp + 1 + perf;
2176 const EvalWell eq_wat_vel = this->primary_variables_.eval(wat_vel_index) - water_velocity;
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;
2184 const EvalWell poly_conc(this->wpolymer());
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);
2191 assembleInjectivityEq(eq_pskin,
2196 this->primary_variables_.numWellEq(),
2204 template<
typename TypeTag>
2213 if constexpr (Base::has_polymermw) {
2215 checkConvergencePolyMW(res, Bhp, this->param_.max_residual_allowed_, report);
2223 template<
typename TypeTag>
2230 std::vector<RateVector>& connectionRates,
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.) {
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;
2249 }
else if (this->isProducer()) {
2250 if (cq_s_polymw < 0.) {
2251 cq_s_polymw *= this->extendEval(int_quants.polymerMoleWeight() );
2258 connectionRates[perf][Indices::contiPolymerMWEqIdx] = Base::restrictEval(cq_s_polymw);
2265 template<
typename TypeTag>
2266 std::optional<typename StandardWell<TypeTag>::Scalar>
2271 const SummaryState& summary_state)
const
2273 return computeBhpAtThpLimitProdWithAlq(simulator,
2276 this->getALQ(well_state),
2280 template<
typename TypeTag>
2281 std::optional<typename StandardWell<TypeTag>::Scalar>
2285 const SummaryState& summary_state,
2287 bool iterate_if_no_solution)
const
2292 auto frates = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2298 std::vector<Scalar> rates(3);
2299 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2300 this->adaptRatesForVFP(rates);
2305 maxPerfPress(simulator),
2306 this->getRefDensity(),
2308 this->getTHPConstraint(summary_state),
2312 auto v = frates(*bhpAtLimit);
2313 if (std::ranges::all_of(v, [](
Scalar i) {
return i <= 0; })) {
2318 if (!iterate_if_no_solution)
2319 return std::nullopt;
2321 auto fratesIter = [
this, &simulator, &groupStateHelper](
const Scalar bhp) {
2325 std::vector<Scalar> rates(3);
2326 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, rates);
2327 this->adaptRatesForVFP(rates);
2333 maxPerfPress(simulator),
2334 this->getRefDensity(),
2336 this->getTHPConstraint(summary_state),
2342 auto v = frates(*bhpAtLimit);
2343 if (std::ranges::all_of(v, [](
Scalar i) {
return i <= 0; })) {
2349 return std::nullopt;
2354 template<
typename TypeTag>
2355 std::optional<typename StandardWell<TypeTag>::Scalar>
2359 const SummaryState& summary_state)
const
2363 auto frates = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2369 std::vector<Scalar> rates(3);
2370 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2376 this->getRefDensity(),
2387 template<
typename TypeTag>
2392 const Well::InjectionControls& inj_controls,
2393 const Well::ProductionControls& prod_controls,
2399 updatePrimaryVariables(groupStateHelper);
2401 const int max_iter = this->param_.max_inner_iter_wells_;
2404 bool relax_convergence =
false;
2405 this->regularize_ =
false;
2407 assembleWellEqWithoutIteration(simulator, groupStateHelper, dt, inj_controls, prod_controls, well_state,
2410 if (it > this->param_.strict_inner_iter_wells_) {
2411 relax_convergence =
true;
2412 this->regularize_ =
true;
2415 auto report = getWellConvergence(groupStateHelper, Base::B_avg_, relax_convergence);
2417 converged = report.converged();
2423 solveEqAndUpdateWellState(simulator, groupStateHelper, well_state);
2430 }
while (it < max_iter);
2433 std::ostringstream sstr;
2434 sstr <<
" Well " << this->name() <<
" converged in " << it <<
" inner iterations.";
2435 if (relax_convergence)
2436 sstr <<
" (A relaxed tolerance was used after "<< this->param_.strict_inner_iter_wells_ <<
" iterations)";
2440 deferred_logger.debug(sstr.str(), OpmLog::defaultDebugVerbosityLevel + (it == 0));
2442 std::ostringstream sstr;
2443 sstr <<
" Well " << this->name() <<
" did not converge in " << it <<
" inner iterations.";
2444 deferred_logger.debug(sstr.str());
2451 template<
typename TypeTag>
2456 const Well::InjectionControls& inj_controls,
2457 const Well::ProductionControls& prod_controls,
2460 const bool fixed_control ,
2461 const bool fixed_status ,
2462 const bool solving_with_zero_rate )
2466 updatePrimaryVariables(groupStateHelper);
2468 const int max_iter = this->param_.max_inner_iter_wells_;
2470 bool converged =
false;
2471 bool relax_convergence =
false;
2472 this->regularize_ =
false;
2473 const auto& summary_state = groupStateHelper.
summaryState();
2478 constexpr int min_its_after_switch = 4;
2480 const int max_status_switch = this->param_.max_well_status_switch_inner_iter_;
2481 int its_since_last_switch = min_its_after_switch;
2482 int switch_count= 0;
2484 const auto well_status_orig = this->wellStatus_;
2485 const auto operability_orig = this->operability_status_;
2486 auto well_status_cur = well_status_orig;
2487 int status_switch_count = 0;
2489 const bool allow_open = well_state.
well(this->index_of_well_).status == WellStatus::OPEN;
2491 const bool allow_switching =
2492 !this->wellUnderZeroRateTarget(groupStateHelper) &&
2493 (!fixed_control || !fixed_status) && allow_open;
2495 bool changed =
false;
2496 bool final_check =
false;
2498 this->operability_status_.resetOperability();
2499 this->operability_status_.solvable =
true;
2501 its_since_last_switch++;
2502 if (allow_switching && its_since_last_switch >= min_its_after_switch && status_switch_count < max_status_switch){
2503 const Scalar wqTotal = this->primary_variables_.eval(WQTotal).value();
2504 changed = this->updateWellControlAndStatusLocalIteration(
2505 simulator, groupStateHelper, inj_controls, prod_controls, wqTotal,
2506 well_state, fixed_control, fixed_status,
2507 solving_with_zero_rate
2510 its_since_last_switch = 0;
2512 if (well_status_cur != this->wellStatus_) {
2513 well_status_cur = this->wellStatus_;
2514 status_switch_count++;
2517 if (!changed && final_check) {
2520 final_check =
false;
2522 if (status_switch_count == max_status_switch) {
2523 this->wellStatus_ = well_status_orig;
2527 assembleWellEqWithoutIteration(simulator, groupStateHelper, dt, inj_controls, prod_controls, well_state, solving_with_zero_rate);
2529 if (it > this->param_.strict_inner_iter_wells_) {
2530 relax_convergence =
true;
2531 this->regularize_ =
true;
2534 auto report = getWellConvergence(groupStateHelper, Base::B_avg_, relax_convergence);
2536 converged = report.converged();
2540 if (switch_count > 0 && its_since_last_switch < min_its_after_switch) {
2542 its_since_last_switch = min_its_after_switch;
2549 solveEqAndUpdateWellState(simulator, groupStateHelper, well_state);
2551 }
while (it < max_iter);
2554 if (allow_switching){
2556 const bool is_stopped = this->wellIsStopped();
2557 if (this->wellHasTHPConstraints(summary_state)){
2558 this->operability_status_.can_obtain_bhp_with_thp_limit = !is_stopped;
2559 this->operability_status_.obey_thp_limit_under_bhp_limit = !is_stopped;
2561 this->operability_status_.operable_under_only_bhp_limit = !is_stopped;
2564 std::string message = fmt::format(
" Well {} converged in {} inner iterations ("
2565 "{} control/status switches).", this->name(), it, switch_count);
2566 if (relax_convergence) {
2567 message.append(fmt::format(
" (A relaxed tolerance was used after {} iterations)",
2568 this->param_.strict_inner_iter_wells_));
2570 deferred_logger.debug(message, OpmLog::defaultDebugVerbosityLevel + ((it == 0) && (switch_count == 0)));
2573 this->wellStatus_ = well_status_orig;
2574 this->operability_status_ = operability_orig;
2575 const std::string message = fmt::format(
" Well {} did not converge in {} inner iterations ("
2576 "{} switches, {} status changes).", this->name(), it, switch_count, status_switch_count);
2577 deferred_logger.debug(message);
2583 template<
typename TypeTag>
2584 std::vector<typename StandardWell<TypeTag>::Scalar>
2590 std::vector<Scalar> well_q_s(this->num_conservation_quantities_, 0.);
2591 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
2592 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
2593 for (
int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
2594 const int cell_idx = this->well_cells_[perf];
2595 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
2596 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.);
2597 getMobility(simulator, perf, mob, deferred_logger);
2598 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.);
2600 getTransMult(trans_mult, simulator, cell_idx);
2601 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
2602 std::vector<Scalar> Tw(this->num_conservation_quantities_, this->well_index_[perf] * trans_mult);
2603 this->getTw(Tw, perf, intQuants, trans_mult, wellstate_nupcol);
2605 computePerfRate(intQuants, mob, bhp.value(), Tw, perf, allow_cf,
2606 cq_s, perf_rates, deferred_logger);
2607 for (
int comp = 0; comp < this->num_conservation_quantities_; ++comp) {
2608 well_q_s[comp] += cq_s[comp];
2611 const auto& comm = this->parallel_well_info_.communication();
2612 if (comm.size() > 1)
2614 comm.sum(well_q_s.data(), well_q_s.size());
2621 template <
typename TypeTag>
2622 std::vector<typename StandardWell<TypeTag>::Scalar>
2626 const int num_pri_vars = this->primary_variables_.numWellEq();
2627 std::vector<Scalar> retval(num_pri_vars);
2628 for (
int ii = 0; ii < num_pri_vars; ++ii) {
2629 retval[ii] = this->primary_variables_.value(ii);
2638 template <
typename TypeTag>
2643 const int num_pri_vars = this->primary_variables_.numWellEq();
2644 for (
int ii = 0; ii < num_pri_vars; ++ii) {
2645 this->primary_variables_.setValue(ii, it[ii]);
2647 return num_pri_vars;
2651 template <
typename TypeTag>
2655 const IntensiveQuantities& intQuants,
2658 auto fs = intQuants.fluidState();
2660 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
2661 if (!FluidSystem::phaseIsActive(phaseIdx)) {
2666 EvalWell cq_r_thermal{0.};
2667 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
2668 const bool both_oil_gas = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
2669 if (!both_oil_gas || FluidSystem::waterPhaseIdx == phaseIdx) {
2670 cq_r_thermal = cq_s[activeCompIdx] / this->extendEval(fs.invB(phaseIdx));
2673 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
2674 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
2679 const EvalWell d = this->extendEval(1.0 - fs.Rv() * fs.Rs());
2681 deferred_logger.
debug(
2682 fmt::format(
"Problematic d value {} obtained for well {}"
2683 " during calculateSinglePerf with rs {}"
2684 ", rv {}. Continue as if no dissolution (rs = 0) and"
2685 " vaporization (rv = 0) for this connection.",
2686 d, this->name(), fs.Rs(), fs.Rv()));
2687 cq_r_thermal = cq_s[activeCompIdx] / this->extendEval(fs.invB(phaseIdx));
2689 if (FluidSystem::gasPhaseIdx == phaseIdx) {
2690 cq_r_thermal = (cq_s[gasCompIdx] -
2691 this->extendEval(fs.Rs()) * cq_s[oilCompIdx]) /
2692 (d * this->extendEval(fs.invB(phaseIdx)) );
2693 }
else if (FluidSystem::oilPhaseIdx == phaseIdx) {
2695 cq_r_thermal = (cq_s[oilCompIdx] - this->extendEval(fs.Rv()) *
2697 (d * this->extendEval(fs.invB(phaseIdx)) );
2703 if (this->isInjector() && !this->wellIsStopped() && cq_r_thermal > 0.0){
2705 assert(this->well_ecl_.injectorType() != InjectorType::MULTI);
2706 fs.setTemperature(this->well_ecl_.inj_temperature());
2707 typedef typename std::decay<
decltype(fs)>::type::ValueType FsValueType;
2708 typename FluidSystem::template ParameterCache<FsValueType> paramCache;
2709 const unsigned pvtRegionIdx = intQuants.pvtRegionIndex();
2710 paramCache.setRegionIndex(pvtRegionIdx);
2711 paramCache.updatePhase(fs, phaseIdx);
2713 const auto& rho = FluidSystem::density(fs, paramCache, phaseIdx);
2714 fs.setDensity(phaseIdx, rho);
2715 const auto& h = FluidSystem::enthalpy(fs, paramCache, phaseIdx);
2716 fs.setEnthalpy(phaseIdx, h);
2717 cq_r_thermal *= this->extendEval(fs.enthalpy(phaseIdx)) * this->extendEval(fs.density(phaseIdx));
2718 result += getValue(cq_r_thermal);
2719 }
else if (cq_r_thermal > 0.0) {
2720 cq_r_thermal *= getValue(fs.enthalpy(phaseIdx)) * getValue(fs.density(phaseIdx));
2721 result += Base::restrictEval(cq_r_thermal);
2724 cq_r_thermal *= this->extendEval(fs.enthalpy(phaseIdx)) * this->extendEval(fs.density(phaseIdx));
2725 result += Base::restrictEval(cq_r_thermal);
2729 return result * this->well_efficiency_factor_;
2732 template <
typename TypeTag>
2736 Scalar max_pressure = 0.0;
2737 for (
int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
2738 const int cell_idx = this->well_cells_[perf];
2739 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2740 const auto& fs = int_quants.fluidState();
2741 Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2742 max_pressure = std::max(max_pressure, pressure_cell);
2744 const auto& comm = this->parallel_well_info_.communication();
2745 if (comm.size() > 1) {
2746 max_pressure = comm.max(max_pressure);
2748 return max_pressure;
#define OPM_DEFLOG_THROW(Exception, message, deferred_logger)
Definition: DeferredLoggingErrorHelpers.hpp:45
#define OPM_DEFLOG_PROBLEM(Exception, message, deferred_logger)
Definition: DeferredLoggingErrorHelpers.hpp:61
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
void debug(const std::string &tag, const std::string &message)
Definition: GroupStateHelper.hpp:57
GroupState< Scalar > & groupState() const
Definition: GroupStateHelper.hpp:298
const SummaryState & summaryState() const
Definition: GroupStateHelper.hpp:415
const WellState< Scalar, IndexTraits > & wellState() const
Definition: GroupStateHelper.hpp:494
DeferredLogger & deferredLogger() const
Get the deferred logger.
Definition: GroupStateHelper.hpp:234
WellStateGuard pushWellState(WellState< Scalar, IndexTraits > &well_state)
Definition: GroupStateHelper.hpp:354
GroupStateGuard pushGroupState(GroupState< Scalar > &group_state)
Definition: GroupStateHelper.hpp:331
Definition: GroupState.hpp:41
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:198
Definition: RatioCalculator.hpp:38
Class handling assemble of the equation system for StandardWell.
Definition: StandardWellAssemble.hpp:44
Scalar pressure_diff(const unsigned perf) const
Returns pressure drop for a given perforation.
Definition: StandardWellConnections.hpp:101
StdWellConnections connections_
Connection level values.
Definition: StandardWellEval.hpp:101
PrimaryVariables primary_variables_
Primary variables for well.
Definition: StandardWellEval.hpp:95
Definition: StandardWell.hpp:60
EvalWell wpolymermw(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2077
std::vector< Scalar > computeWellPotentialWithTHP(const Simulator &ebosSimulator, const GroupStateHelperType &groupStateHelper, const WellStateType &well_state) const
Definition: StandardWell_impl.hpp:1624
typename StdWellEval::EvalWell EvalWell
Definition: StandardWell.hpp:121
void updateWellStateFromPrimaryVariables(WellStateType &well_state, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:810
virtual ConvergenceReport getWellConvergence(const GroupStateHelperType &groupStateHelper, const std::vector< Scalar > &B_avg, const bool relax_tolerance) const override
check whether the well equations get converged for this well
Definition: StandardWell_impl.hpp:1226
WellConnectionProps computePropertiesForWellConnectionPressures(const Simulator &simulator, const WellStateType &well_state) const
Definition: StandardWell_impl.hpp:1175
std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, const Scalar alq_value, bool iterate_if_no_solution) const override
Definition: StandardWell_impl.hpp:2283
typename StdWellEval::BVectorWell BVectorWell
Definition: StandardWell.hpp:122
void addWellContributions(SparseMatrixAdapter &mat) const override
Definition: StandardWell_impl.hpp:1974
std::vector< Scalar > getPrimaryVars() const override
Definition: StandardWell_impl.hpp:2624
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: StandardWell_impl.hpp:758
void updatePrimaryVariables(const GroupStateHelperType &groupStateHelper) override
Definition: StandardWell_impl.hpp:1874
void solveEqAndUpdateWellState(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) override
Definition: StandardWell_impl.hpp:1418
void computeWellConnectionDensitesPressures(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const WellConnectionProps &props)
Definition: StandardWell_impl.hpp:1344
std::optional< Scalar > computeBhpAtThpLimitProd(const WellStateType &well_state, const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state) const
Definition: StandardWell_impl.hpp:2268
void addWellPressureEquations(PressureMatrix &mat, const BVector &x, const int pressureVarIndex, const bool use_well_weights, const WellStateType &well_state) const override
Definition: StandardWell_impl.hpp:1982
void updateWaterMobilityWithPolymer(const Simulator &simulator, const int perf, std::vector< EvalWell > &mob_water, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1908
bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) override
Definition: StandardWell_impl.hpp:2390
std::vector< Scalar > computeCurrentWellRates(const Simulator &ebosSimulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:2586
void calculateSinglePerf(const Simulator &simulator, const int perf, WellStateType &well_state, std::vector< RateVector > &connectionRates, std::vector< EvalWell > &cq_s, EvalWell &water_flux_s, EvalWell &cq_s_zfrac_effective, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:504
void computeWellConnectionPressures(const Simulator &simulator, const GroupStateHelperType &groupStateHelper)
Definition: StandardWell_impl.hpp:1401
void updatePrimaryVariablesNewton(const BVectorWell &dwells, const bool stop_or_zero_rate_target, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:787
void assembleWellEqWithoutIteration(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellStateType &well_state, const bool solving_with_zero_rate) override
Definition: StandardWell_impl.hpp:342
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_potentials) override
computing the well potentials for group control
Definition: StandardWell_impl.hpp:1796
StandardWell(const Well &well, const ParallelWellInfo< Scalar > &pw_info, const int time_step, const ModelParameters ¶m, const RateConverterType &rate_converter, const int pvtRegionIdx, const int num_conservation_quantities, const int num_phases, const int index_of_well, const std::vector< PerforationData< Scalar > > &perf_data)
Definition: StandardWell_impl.hpp:53
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state) const
Definition: StandardWell_impl.hpp:2357
typename StdWellEval::StdWellConnections::Properties WellConnectionProps
Definition: StandardWell.hpp:261
void computeWellRatesWithBhpIterations(const Simulator &ebosSimulator, const Scalar &bhp, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_flux) const override
Definition: StandardWell_impl.hpp:1553
void updateConnectionRatePolyMW(const EvalWell &cq_s_poly, const IntensiveQuantities &int_quants, const WellStateType &well_state, const int perf, std::vector< RateVector > &connectionRates, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2226
void computeWellRatesWithBhp(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1505
void updateIPRImplicit(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) override
Definition: StandardWell_impl.hpp:926
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:702
void getTransMult(Value &trans_mult, const Simulator &simulator, const int cell_indx) const
Definition: StandardWell_impl.hpp:682
void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:834
void handleInjectivityEquations(const Simulator &simulator, const WellStateType &well_state, const int perf, const EvalWell &water_flux_s, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:2160
virtual void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: StandardWell_impl.hpp:1453
void checkConvergenceExtraEqs(const std::vector< Scalar > &res, ConvergenceReport &report) const
Definition: StandardWell_impl.hpp:2207
void computeWellRatesWithThpAlqProd(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1779
typename StdWellEval::Eval Eval
Definition: StandardWell.hpp:120
Scalar computeWellRatesAndBhpWithThpAlqProd(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1749
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1164
bool computeWellPotentialsImplicit(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_potentials) const
Definition: StandardWell_impl.hpp:1660
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) override
Definition: StandardWell_impl.hpp:1485
Scalar maxPerfPress(const Simulator &simulator) const override
Definition: StandardWell_impl.hpp:2735
bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, const bool fixed_control, const bool fixed_status, const bool solving_with_zero_rate) override
Definition: StandardWell_impl.hpp:2454
void assembleWellEqWithoutIterationImpl(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellStateType &well_state, const bool solving_with_zero_rate)
Definition: StandardWell_impl.hpp:367
bool allDrawDownWrongDirection(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1122
EvalWell pskin(const Scalar throughput, const EvalWell &water_velocity, const EvalWell &poly_inj_conc, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2033
void computePerfRate(const IntensiveQuantities &intQuants, const std::vector< Value > &mob, const Value &bhp, const std::vector< Value > &Tw, const int perf, const bool allow_cf, std::vector< Value > &cq_s, PerforationRates< Scalar > &perf_rates, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:94
static constexpr int numWellConservationEq
Definition: StandardWell.hpp:97
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: StandardWell_impl.hpp:2641
void updateWaterThroughput(const double dt, WellStateType &well_state) const override
Definition: StandardWell_impl.hpp:2106
void checkOperabilityUnderBHPLimit(const WellStateType &well_state, const Simulator &simulator, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1003
EvalWell pskinwater(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2002
void handleInjectivityRate(const Simulator &simulator, const int perf, std::vector< EvalWell > &cq_s) const
Definition: StandardWell_impl.hpp:2137
virtual void init(const std::vector< Scalar > &depth_arg, const Scalar gravity_arg, const std::vector< Scalar > &B_avg, const bool changed_to_open_this_step) override
Definition: StandardWell_impl.hpp:77
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1269
void calculateExplicitQuantities(const Simulator &simulator, const GroupStateHelperType &groupStateHelper) override
Definition: StandardWell_impl.hpp:1440
Scalar getRefDensity() const override
Definition: StandardWell_impl.hpp:1897
void checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper) override
Definition: StandardWell_impl.hpp:1072
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: StandardWell_impl.hpp:1859
EvalWell getQs(const int compIdx) const
Returns scaled rate for a component.
Class for computing BHP limits.
Definition: WellBhpThpCalculator.hpp:41
Scalar calculateThpFromBhp(const std::vector< Scalar > &rates, const Scalar bhp, const Scalar rho, const std::optional< Scalar > &alq, const Scalar thp_limit, DeferredLogger &deferred_logger) const
Calculates THP from BHP.
std::optional< Scalar > computeBhpAtThpLimitProd(const std::function< std::vector< Scalar >(const Scalar)> &frates, const SummaryState &summary_state, const Scalar maxPerfPress, const Scalar rho, const Scalar alq_value, const Scalar thp_limit, DeferredLogger &deferred_logger) const
Compute BHP from THP limit for a producer.
Scalar mostStrictBhpFromBhpLimits(const SummaryState &summaryState) const
Obtain the most strict BHP from BHP limits.
std::optional< Scalar > computeBhpAtThpLimitInj(const std::function< std::vector< Scalar >(const Scalar)> &frates, const SummaryState &summary_state, const Scalar rho, const Scalar flo_rel_tol, const int max_iteration, const bool throwOnError, DeferredLogger &deferred_logger) const
Compute BHP from THP limit for an injector.
Definition: WellConvergence.hpp:38
const int num_conservation_quantities_
Definition: WellInterfaceGeneric.hpp:315
Well well_ecl_
Definition: WellInterfaceGeneric.hpp:305
void onlyKeepBHPandTHPcontrols(const SummaryState &summary_state, WellStateType &well_state, Well::InjectionControls &inj_controls, Well::ProductionControls &prod_controls) const
void resetDampening()
Definition: WellInterfaceGeneric.hpp:247
std::pair< bool, bool > computeWellPotentials(std::vector< Scalar > &well_potentials, const WellStateType &well_state)
Definition: WellInterfaceIndices.hpp:34
Definition: WellInterface.hpp:77
bool solveWellWithOperabilityCheck(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:601
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:82
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:105
void getTransMult(Value &trans_mult, const Simulator &simulator, const int cell_idx, Callback &extendEval) const
Definition: WellInterface_impl.hpp:2071
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:98
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:2084
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:87
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:97
typename Base::ModelParameters ModelParameters
Definition: WellInterface.hpp:111
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:84
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:86
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:89
Definition: WellProdIndexCalculator.hpp:37
Scalar connectionProdIndStandard(const std::size_t connIdx, const Scalar connMobility) const
Definition: WellState.hpp:66
const SingleWellState< Scalar, IndexTraits > & well(std::size_t well_index) const
Definition: WellState.hpp:290
std::vector< Scalar > & wellRates(std::size_t well_index)
One rate per well and phase.
Definition: WellState.hpp:255
@ NONE
Definition: DeferredLogger.hpp:46
Definition: blackoilbioeffectsmodules.hh:45
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Static data associated with a well perforation.
Definition: PerforationData.hpp:30
Definition: PerforationData.hpp:41
Scalar dis_gas
Definition: PerforationData.hpp:42
Scalar vap_wat
Definition: PerforationData.hpp:45
Scalar vap_oil
Definition: PerforationData.hpp:44
Scalar dis_gas_in_water
Definition: PerforationData.hpp:43