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>
45#include <fmt/format.h>
50 template<
typename TypeTag>
57 const int pvtRegionIdx,
58 const int num_conservation_quantities,
60 const int index_of_well,
62 :
Base(well, pw_info, time_step, param, rate_converter, pvtRegionIdx, num_conservation_quantities, num_phases, index_of_well, perf_data)
73 template<
typename TypeTag>
76 init(
const std::vector<Scalar>& depth_arg,
78 const std::vector< Scalar >& B_avg,
79 const bool changed_to_open_this_step)
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);
89 template<
typename TypeTag>
94 const std::vector<Value>& mob,
96 const std::vector<Value>& Tw,
99 std::vector<Value>& cq_s,
103 auto obtain = [
this](
const Eval& value)
105 if constexpr (std::is_same_v<Value, Scalar>) {
106 static_cast<void>(
this);
107 return getValue(value);
109 return this->extendEval(value);
112 auto obtainN = [](
const auto& value)
114 if constexpr (std::is_same_v<Value, Scalar>) {
115 return getValue(value);
120 auto zeroElem = [
this]()
122 if constexpr (std::is_same_v<Value, Scalar>) {
123 static_cast<void>(
this);
126 return Value{this->primary_variables_.numWellEq() + Indices::numEq, 0.0};
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());
137 std::vector<Value> b_perfcells_dense(this->numConservationQuantities(), zeroElem());
138 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
139 if (!FluidSystem::phaseIsActive(phaseIdx)) {
142 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
143 b_perfcells_dense[compIdx] = obtain(fs.invB(phaseIdx));
145 if constexpr (has_solvent) {
146 b_perfcells_dense[Indices::contiSolventEqIdx] = obtain(intQuants.solventInverseFormationVolumeFactor());
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();
157 Value skin_pressure = zeroElem();
159 if (this->isInjector()) {
160 const int pskin_index = Bhp + 1 + this->numLocalPerfs() + perf;
161 skin_pressure = obtainN(this->primary_variables_.eval(pskin_index));
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));
191 template<
typename TypeTag>
192 template<
class Value>
196 const Value& pressure,
202 std::vector<Value>& b_perfcells_dense,
203 const std::vector<Value>& Tw,
206 const Value& skin_pressure,
207 const std::vector<Value>& cmix_s,
208 std::vector<Value>& cq_s,
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;
220 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)
221 ? FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx)
223 FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)
224 ? FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx)
226 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)
227 ? FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx)
235 if (!allow_cf && this->isInjector()) {
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;
245 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
246 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
248 ratioCalc.gasOilPerfRateProd(cq_s, perf_rates, rv, rs, rvw,
249 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx),
251 }
else if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) &&
252 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
254 ratioCalc.gasWaterPerfRateProd(cq_s, perf_rates, rvw, rsw, this->isProducer());
258 if (!allow_cf && this->isProducer()) {
263 Value total_mob_dense = mob[0];
264 for (
int componentIdx = 1; componentIdx < this->numConservationQuantities(); ++componentIdx) {
265 total_mob_dense += mob[componentIdx];
269 Value volumeRatio = bhp * 0.0;
271 if (FluidSystem::enableVaporizedWater() && FluidSystem::enableDissolvedGasInWater()) {
272 ratioCalc.disOilVapWatVolumeRatio(volumeRatio, rvw, rsw, pressure,
273 cmix_s, b_perfcells_dense, deferred_logger);
277 assert(FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx));
278 assert(FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx));
279 assert(!FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx));
282 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
283 const unsigned waterCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
284 volumeRatio += cmix_s[waterCompIdx] / b_perfcells_dense[waterCompIdx];
287 if constexpr (Indices::enableSolvent) {
288 volumeRatio += cmix_s[Indices::contiSolventEqIdx] / b_perfcells_dense[Indices::contiSolventEqIdx];
291 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
292 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
294 ratioCalc.gasOilVolumeRatio(volumeRatio, rv, rs, pressure,
295 cmix_s, b_perfcells_dense,
298 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
299 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
300 volumeRatio += cmix_s[oilCompIdx] / b_perfcells_dense[oilCompIdx];
302 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
303 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
304 volumeRatio += cmix_s[gasCompIdx] / b_perfcells_dense[gasCompIdx];
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;
317 if (this->isProducer()) {
318 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
319 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
321 ratioCalc.gasOilPerfRateInj(cq_s, perf_rates,
322 rv, rs, pressure, rvw,
323 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx),
326 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) &&
327 FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx))
330 ratioCalc.gasWaterPerfRateInj(cq_s, perf_rates, rvw, rsw,
331 pressure, deferred_logger);
338 template<
typename TypeTag>
344 const Well::InjectionControls& inj_controls,
345 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,
359 deferred_logger, solving_with_zero_rate);
365 template<
typename TypeTag>
371 const Well::InjectionControls& inj_controls,
372 const Well::ProductionControls& prod_controls,
375 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, deferred_logger);
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,
493 this->linSys_.invert();
495 OPM_DEFLOG_PROBLEM(NumericalProblem,
"Error when inverting local well equations for well " + name(), deferred_logger);
502 template<
typename TypeTag>
508 std::vector<RateVector>& connectionRates,
509 std::vector<EvalWell>& cq_s,
514 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
515 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
516 const int cell_idx = this->well_cells_[perf];
517 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
518 std::vector<EvalWell> mob(this->num_conservation_quantities_, {0.});
519 getMobility(simulator, perf, mob, deferred_logger);
523 getTransMult(trans_mult, simulator, cell_idx);
524 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
525 std::vector<EvalWell> Tw(this->num_conservation_quantities_, this->well_index_[perf] * trans_mult);
526 this->getTw(Tw, perf, intQuants, trans_mult, wellstate_nupcol);
527 computePerfRate(intQuants, mob, bhp, Tw, perf, allow_cf,
528 cq_s, perf_rates, deferred_logger);
530 auto& ws = well_state.
well(this->index_of_well_);
531 auto& perf_data = ws.perf_data;
532 if constexpr (has_polymer && Base::has_polymermw) {
533 if (this->isInjector()) {
536 const unsigned water_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
537 water_flux_s = cq_s[water_comp_idx];
540 handleInjectivityRate(simulator, perf, cq_s);
545 if (this->isProducer()) {
546 ws.phase_mixing_rates[ws.dissolved_gas] += perf_rates.
dis_gas;
547 ws.phase_mixing_rates[ws.dissolved_gas_in_water] += perf_rates.
dis_gas_in_water;
548 ws.phase_mixing_rates[ws.vaporized_oil] += perf_rates.
vap_oil;
549 ws.phase_mixing_rates[ws.vaporized_water] += perf_rates.
vap_wat;
550 perf_data.phase_mixing_rates[perf][ws.dissolved_gas] = perf_rates.
dis_gas;
551 perf_data.phase_mixing_rates[perf][ws.dissolved_gas_in_water] = perf_rates.
dis_gas_in_water;
552 perf_data.phase_mixing_rates[perf][ws.vaporized_oil] = perf_rates.
vap_oil;
553 perf_data.phase_mixing_rates[perf][ws.vaporized_water] = perf_rates.
vap_wat;
556 if constexpr (has_energy) {
557 connectionRates[perf][Indices::contiEnergyEqIdx] =
558 connectionRateEnergy(cq_s, intQuants, deferred_logger);
559 ws.energy_rate += getValue(connectionRates[perf][Indices::contiEnergyEqIdx]);
562 if constexpr (has_polymer) {
563 std::variant<Scalar,EvalWell> polymerConcentration;
564 if (this->isInjector()) {
565 polymerConcentration = this->wpolymer();
567 polymerConcentration = this->extendEval(intQuants.polymerConcentration() *
568 intQuants.polymerViscosityCorrection());
571 [[maybe_unused]]
EvalWell cq_s_poly;
572 std::tie(connectionRates[perf][Indices::contiPolymerEqIdx],
574 this->connections_.connectionRatePolymer(perf_data.polymer_rates[perf],
575 cq_s, polymerConcentration);
577 if constexpr (Base::has_polymermw) {
578 updateConnectionRatePolyMW(cq_s_poly, intQuants, well_state,
579 perf, connectionRates, deferred_logger);
583 if constexpr (has_foam) {
584 std::variant<Scalar,EvalWell> foamConcentration;
585 if (this->isInjector()) {
586 foamConcentration = this->wfoam();
588 foamConcentration = this->extendEval(intQuants.foamConcentration());
590 connectionRates[perf][Indices::contiFoamEqIdx] =
591 this->connections_.connectionRateFoam(cq_s, foamConcentration,
592 FoamModule::transportPhase(),
596 if constexpr (has_zFraction) {
597 std::variant<Scalar,std::array<EvalWell,2>> solventConcentration;
598 if (this->isInjector()) {
599 solventConcentration = this->wsolvent();
601 solventConcentration = std::array{this->extendEval(intQuants.xVolume()),
602 this->extendEval(intQuants.yVolume())};
604 std::tie(connectionRates[perf][Indices::contiZfracEqIdx],
605 cq_s_zfrac_effective) =
606 this->connections_.connectionRatezFraction(perf_data.solvent_rates[perf],
608 solventConcentration);
611 if constexpr (has_brine) {
612 std::variant<Scalar,EvalWell> saltConcentration;
613 if (this->isInjector()) {
614 saltConcentration = this->wsalt();
616 saltConcentration = this->extendEval(intQuants.fluidState().saltConcentration());
619 connectionRates[perf][Indices::contiBrineEqIdx] =
620 this->connections_.connectionRateBrine(perf_data.brine_rates[perf],
625 if constexpr (has_bioeffects) {
626 std::variant<Scalar,EvalWell> microbialConcentration;
627 if constexpr (has_micp) {
628 std::variant<Scalar,EvalWell> oxygenConcentration;
629 std::variant<Scalar,EvalWell> ureaConcentration;
630 if (this->isInjector()) {
631 microbialConcentration = this->wmicrobes();
632 oxygenConcentration = this->woxygen();
633 ureaConcentration = this->wurea();
635 microbialConcentration = this->extendEval(intQuants.microbialConcentration());
636 oxygenConcentration = this->extendEval(intQuants.oxygenConcentration());
637 ureaConcentration = this->extendEval(intQuants.ureaConcentration());
639 std::tie(connectionRates[perf][Indices::contiMicrobialEqIdx],
640 connectionRates[perf][Indices::contiOxygenEqIdx],
641 connectionRates[perf][Indices::contiUreaEqIdx]) =
642 this->connections_.connectionRatesMICP(perf_data.microbial_rates[perf],
643 perf_data.oxygen_rates[perf],
644 perf_data.urea_rates[perf],
646 microbialConcentration,
651 if (this->isProducer()) {
652 microbialConcentration = this->extendEval(intQuants.microbialConcentration());
653 connectionRates[perf][Indices::contiMicrobialEqIdx] =
654 this->connections_.connectionRateBioeffects(perf_data.microbial_rates[perf],
656 microbialConcentration);
662 perf_data.pressure[perf] = ws.bhp + this->connections_.pressure_diff(perf);
665 if (FluidSystem::phaseUsage().hasCO2orH2Store()) {
666 const unsigned gas_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
667 const Scalar rho = FluidSystem::referenceDensity( FluidSystem::gasPhaseIdx, Base::pvtRegionIdx() );
668 perf_data.gas_mass_rates[perf] = cq_s[gas_comp_idx].value() * rho;
672 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
673 const unsigned wat_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
674 const Scalar rho = FluidSystem::referenceDensity( FluidSystem::waterPhaseIdx, Base::pvtRegionIdx() );
675 perf_data.wat_mass_rates[perf] = cq_s[wat_comp_idx].value() * rho;
679 template<
typename TypeTag>
680 template<
class Value>
685 const int cell_idx)
const
687 auto obtain = [
this](
const Eval& value)
689 if constexpr (std::is_same_v<Value, Scalar>) {
690 static_cast<void>(
this);
691 return getValue(value);
693 return this->extendEval(value);
699 template<
typename TypeTag>
700 template<
class Value>
705 std::vector<Value>& mob,
708 auto obtain = [
this](
const Eval& value)
710 if constexpr (std::is_same_v<Value, Scalar>) {
711 static_cast<void>(
this);
712 return getValue(value);
714 return this->extendEval(value);
718 obtain, deferred_logger);
721 if constexpr (has_polymer) {
722 if (!FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
723 OPM_DEFLOG_THROW(std::runtime_error,
"Water is required when polymer is active", deferred_logger);
728 if constexpr (!Base::has_polymermw) {
729 if constexpr (std::is_same_v<Value, Scalar>) {
730 std::vector<EvalWell> mob_eval(this->num_conservation_quantities_, 0.);
731 for (std::size_t i = 0; i < mob.size(); ++i) {
732 mob_eval[i].setValue(mob[i]);
734 updateWaterMobilityWithPolymer(simulator, perf, mob_eval, deferred_logger);
735 for (std::size_t i = 0; i < mob.size(); ++i) {
736 mob[i] = getValue(mob_eval[i]);
739 updateWaterMobilityWithPolymer(simulator, perf, mob, deferred_logger);
746 const Scalar bhp = this->primary_variables_.value(Bhp);
747 const Scalar perf_press = bhp + this->connections_.pressure_diff(perf);
748 const Scalar multiplier = this->getInjMult(perf, bhp, perf_press, deferred_logger);
749 for (std::size_t i = 0; i < mob.size(); ++i) {
750 mob[i] *= multiplier;
756 template<
typename TypeTag>
765 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
767 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(groupStateHelper, deferred_logger);
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);
772 Base::calculateReservoirRates(simulator.vanguard().eclState().runspec().co2Storage(), well_state.
well(this->index_of_well_));
779 template<
typename TypeTag>
783 const bool stop_or_zero_rate_target,
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);
791 if constexpr (Base::has_polymermw) {
792 this->primary_variables_.updateNewtonPolyMW(dwells);
795 this->primary_variables_.checkFinite(deferred_logger,
"Newton update");
802 template<
typename TypeTag>
806 const SummaryState& summary_state,
809 this->primary_variables_.copyToWellState(well_state, deferred_logger);
812 updateThp(getRefDensity(),
813 [
this,&well_state]() {
return this->baseif_.getALQ(well_state); },
814 well_state, summary_state, deferred_logger);
817 if constexpr (Base::has_polymermw) {
818 this->primary_variables_.copyToWellStatePolyMW(well_state);
826 template<
typename TypeTag>
834 std::fill(this->ipr_a_.begin(), this->ipr_a_.end(), 0.);
835 std::fill(this->ipr_b_.begin(), this->ipr_b_.end(), 0.);
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);
841 const int cell_idx = this->well_cells_[perf];
842 const auto& int_quantities = simulator.model().intensiveQuantities(cell_idx, 0);
843 const auto& fs = int_quantities.fluidState();
845 Scalar p_r = this->getPerfCellPressure(fs).value();
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)) {
853 const unsigned comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phase));
854 b_perf[comp_idx] = fs.invB(phase).value();
856 if constexpr (has_solvent) {
857 b_perf[Indices::contiSolventEqIdx] = int_quantities.solventInverseFormationVolumeFactor().value();
861 const Scalar h_perf = this->connections_.pressure_diff(perf);
862 const Scalar pressure_diff = p_r - h_perf;
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");
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;
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();
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];
899 ipr_a_perf[gas_comp_idx] += dis_gas_a;
900 ipr_a_perf[oil_comp_idx] += vap_oil_a;
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];
905 ipr_b_perf[gas_comp_idx] += dis_gas_b;
906 ipr_b_perf[oil_comp_idx] += vap_oil_b;
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];
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());
918 template<
typename TypeTag>
933 auto rates = well_state.
well(this->index_of_well_).surface_rates;
935 for (std::size_t p = 0; p < rates.size(); ++p) {
936 zero_rates &= rates[p] == 0.0;
938 auto& ws = well_state.
well(this->index_of_well_);
940 const auto msg = fmt::format(
"updateIPRImplicit: Well {} has zero rate, IPRs might be problematic", this->name());
941 deferred_logger.
debug(msg);
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.);
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;
963 const auto cmode = ws.production_cmode;
964 ws.production_cmode = Well::ProducerCMode::BHP;
965 const double dt = simulator.timeStepSize();
966 assembleWellEqWithoutIteration(simulator, groupStateHelper, dt, inj_controls, prod_controls, well_state, deferred_logger,
969 const size_t nEq = this->primary_variables_.numWellEq();
973 for (
size_t i=0; i < nEq; ++i){
979 x_well[0].resize(nEq);
980 this->linSys_.solve(rhs, x_well);
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) {
987 ws.implicit_ipr_b[idx] -= x_well[0][pvIdx]*comp_rate.derivative(pvIdx+Indices::numEq);
989 ws.implicit_ipr_a[idx] = ws.implicit_ipr_b[idx]*ws.bhp - comp_rate.value();
992 ws.production_cmode = cmode;
995 template<
typename TypeTag>
1002 const auto& summaryState = simulator.vanguard().summaryState();
1006 const bool bhp_limit_not_defaulted = bhp_limit > 1.5 * unit::barsa;
1007 if ( bhp_limit_not_defaulted || !this->wellHasTHPConstraints(summaryState) ) {
1010 Scalar total_ipr_mass_rate = 0.0;
1011 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
1013 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1017 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
1018 const Scalar ipr_rate = this->ipr_a_[compIdx] - this->ipr_b_[compIdx] * bhp_limit;
1020 const Scalar rho = FluidSystem::referenceDensity( phaseIdx, Base::pvtRegionIdx() );
1021 total_ipr_mass_rate += ipr_rate * rho;
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;
1028 if (this->operability_status_.operable_under_only_bhp_limit && this->wellHasTHPConstraints(summaryState)) {
1032 std::vector<Scalar> well_rates_bhp_limit;
1033 computeWellRatesWithBhp(simulator, bhp_limit, well_rates_bhp_limit, deferred_logger);
1035 this->adaptRatesForVFP(well_rates_bhp_limit);
1036 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1039 this->getRefDensity(),
1040 this->getALQ(well_state),
1043 if ( (this->isProducer() && thp < thp_limit) || (this->isInjector() && thp > thp_limit) ) {
1044 this->operability_status_.obey_thp_limit_under_bhp_limit =
false;
1055 this->operability_status_.operable_under_only_bhp_limit =
true;
1056 this->operability_status_.obey_thp_limit_under_bhp_limit =
false;
1064 template<
typename TypeTag>
1072 const auto& summaryState = simulator.vanguard().summaryState();
1073 const auto obtain_bhp = this->isProducer() ? computeBhpAtThpLimitProd(well_state, simulator, groupStateHelper, summaryState, deferred_logger)
1074 : computeBhpAtThpLimitInj(simulator, groupStateHelper, summaryState, deferred_logger);
1077 this->operability_status_.can_obtain_bhp_with_thp_limit =
true;
1080 this->operability_status_.obey_bhp_limit_with_thp_limit = this->isProducer() ?
1081 *obtain_bhp >= bhp_limit : *obtain_bhp <= bhp_limit ;
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 "
1088 +
" bars as a producer for well " + name();
1089 deferred_logger.
debug(msg);
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 "
1095 +
" bars as a injector for well " + name();
1096 deferred_logger.
debug(msg);
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 "
1105 +
" bar for well " + name() +
", the well might need to be closed ");
1114 template<
typename TypeTag>
1119 bool all_drawdown_wrong_direction =
true;
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, 0);
1124 const auto& fs = intQuants.fluidState();
1126 const Scalar pressure = this->getPerfCellPressure(fs).value();
1127 const Scalar bhp = this->primary_variables_.eval(Bhp).value();
1130 const Scalar well_pressure = bhp + this->connections_.pressure_diff(perf);
1131 const Scalar drawdown = pressure - well_pressure;
1136 if ( (drawdown < 0. && this->isInjector()) ||
1137 (drawdown > 0. && this->isProducer()) ) {
1138 all_drawdown_wrong_direction =
false;
1143 const auto& comm = this->parallel_well_info_.communication();
1144 if (comm.size() > 1)
1146 all_drawdown_wrong_direction =
1147 (comm.min(all_drawdown_wrong_direction ? 1 : 0) == 1);
1150 return all_drawdown_wrong_direction;
1156 template<
typename TypeTag>
1161 return !this->getAllowCrossFlow() && allDrawDownWrongDirection(simulator);
1167 template<
typename TypeTag>
1173 auto prop_func =
typename StdWellEval::StdWellConnections::PressurePropertyFunctions {
1175 [&model = simulator.model()](
int cell_idx,
int phase_idx)
1177 return model.intensiveQuantities(cell_idx, 0)
1178 .fluidState().temperature(phase_idx).value();
1182 [&model = simulator.model()](
int cell_idx)
1184 return model.intensiveQuantities(cell_idx, 0)
1185 .fluidState().saltConcentration().value();
1189 [&model = simulator.model()](
int cell_idx)
1191 return model.intensiveQuantities(cell_idx, 0)
1192 .fluidState().pvtRegionIndex();
1196 if constexpr (Indices::enableSolvent) {
1197 prop_func.solventInverseFormationVolumeFactor =
1198 [&model = simulator.model()](
int cell_idx)
1200 return model.intensiveQuantities(cell_idx, 0)
1201 .solventInverseFormationVolumeFactor().value();
1204 prop_func.solventRefDensity = [&model = simulator.model()](
int cell_idx)
1206 return model.intensiveQuantities(cell_idx, 0)
1207 .solventRefDensity();
1211 return this->connections_.computePropertiesForPressures(well_state, prop_func);
1218 template<
typename TypeTag>
1222 const std::vector<Scalar>& B_avg,
1224 const bool relax_tolerance)
const
1228 assert((
int(B_avg.size()) == this->num_conservation_quantities_) || has_polymer || has_energy || has_foam || has_brine || has_zFraction || has_bioeffects);
1230 Scalar tol_wells = this->param_.tolerance_wells_;
1232 constexpr Scalar stopped_factor = 1.e-4;
1234 constexpr Scalar dynamic_thp_factor = 1.e-1;
1235 if (this->stoppedOrZeroRateTarget(groupStateHelper, deferred_logger)) {
1236 tol_wells = tol_wells*stopped_factor;
1237 }
else if (this->getDynamicThpLimit()) {
1238 tol_wells = tol_wells*dynamic_thp_factor;
1241 std::vector<Scalar> res;
1244 this->param_.max_residual_allowed_,
1246 this->param_.relaxed_tolerance_flow_well_,
1248 this->wellIsStopped(),
1252 checkConvergenceExtraEqs(res, report);
1261 template<
typename TypeTag>
1269 auto fluidState = [&simulator,
this](
const int perf)
1271 const auto cell_idx = this->well_cells_[perf];
1272 return simulator.model()
1273 .intensiveQuantities(cell_idx, 0).fluidState();
1276 const int np = this->number_of_phases_;
1277 auto setToZero = [np](
Scalar* x) ->
void
1279 std::fill_n(x, np, 0.0);
1282 auto addVector = [np](
const Scalar* src,
Scalar* dest) ->
void
1284 std::transform(src, src + np, dest, dest, std::plus<>{});
1287 auto& ws = well_state.
well(this->index_of_well_);
1288 auto& perf_data = ws.perf_data;
1289 auto* wellPI = ws.productivity_index.data();
1290 auto* connPI = perf_data.prod_index.data();
1294 const auto preferred_phase = this->well_ecl_.getPreferredPhase();
1295 auto subsetPerfID = 0;
1297 for (
const auto& perf : *this->perf_data_) {
1298 auto allPerfID = perf.ecl_index;
1300 auto connPICalc = [&wellPICalc, allPerfID](
const Scalar mobility) ->
Scalar
1305 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
1306 getMobility(simulator,
static_cast<int>(subsetPerfID), mob, deferred_logger);
1308 const auto& fs = fluidState(subsetPerfID);
1311 if (this->isInjector()) {
1312 this->computeConnLevelInjInd(fs, preferred_phase, connPICalc,
1313 mob, connPI, deferred_logger);
1316 this->computeConnLevelProdInd(fs, connPICalc, mob, connPI);
1319 addVector(connPI, wellPI);
1326 const auto& comm = this->parallel_well_info_.communication();
1327 if (comm.size() > 1) {
1328 comm.sum(wellPI, np);
1331 assert ((
static_cast<int>(subsetPerfID) == this->number_of_local_perforations_) &&
1332 "Internal logic error in processing connections for PI/II");
1337 template<
typename TypeTag>
1344 const auto& well_state = groupStateHelper.
wellState();
1348 const auto prop_func =
typename StdWellEval::StdWellConnections::DensityPropertyFunctions {
1353 [&model = simulator.model()](
const int cell,
1354 const std::vector<int>& phases,
1355 std::vector<Scalar>& mob)
1357 const auto& iq = model.intensiveQuantities(cell, 0);
1359 std::transform(phases.begin(), phases.end(), mob.begin(),
1360 [&iq](
const int phase) { return iq.mobility(phase).value(); });
1365 [&model = simulator.model()](
const int cell,
1366 const std::vector<int>& phases,
1367 std::vector<Scalar>& rho)
1369 const auto& fs = model.intensiveQuantities(cell, 0).fluidState();
1371 std::transform(phases.begin(), phases.end(), rho.begin(),
1372 [&fs](
const int phase) { return fs.density(phase).value(); });
1376 const auto stopped_or_zero_rate_target = this->
1377 stoppedOrZeroRateTarget(groupStateHelper, deferred_logger);
1380 .computeProperties(stopped_or_zero_rate_target, well_state,
1381 prop_func, props, deferred_logger);
1383 cachedRefDensity = this->connections_.rho(0);
1384 if (this->parallel_well_info_.communication().size() > 1) {
1385 cachedRefDensity = this->parallel_well_info_.broadcastFirstPerforationValue(cachedRefDensity);
1393 template<
typename TypeTag>
1400 const auto& well_state = groupStateHelper.
wellState();
1401 const auto props = computePropertiesForWellConnectionPressures
1402 (simulator, well_state);
1404 computeWellConnectionDensitesPressures(simulator, groupStateHelper,
1405 props, deferred_logger);
1412 template<
typename TypeTag>
1420 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1425 dx_well[0].resize(this->primary_variables_.numWellEq());
1426 this->linSys_.solve( dx_well);
1428 updateWellState(simulator, dx_well, groupStateHelper, well_state, deferred_logger);
1435 template<
typename TypeTag>
1442 updatePrimaryVariables(groupStateHelper, deferred_logger);
1443 computeWellConnectionPressures(simulator, groupStateHelper, deferred_logger);
1444 this->computeAccumWell();
1449 template<
typename TypeTag>
1454 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1456 if (this->param_.matrix_add_well_contributions_)
1462 this->linSys_.apply(x, Ax);
1468 template<
typename TypeTag>
1473 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1475 this->linSys_.apply(r);
1481 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, deferred_logger);
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,
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, deferred_logger
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, deferred_logger);
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 deferred_logger, potentials, this->getALQ(well_state)
1657 template<
typename TypeTag>
1662 std::vector<Scalar>& well_potentials,
1676 auto guard = groupStateHelper_copy.
pushWellState(well_state_copy);
1677 auto& ws = well_state_copy.
well(this->index_of_well_);
1680 const auto& summary_state = simulator.vanguard().summaryState();
1681 auto inj_controls = well_copy.
well_ecl_.isInjector()
1682 ? well_copy.
well_ecl_.injectionControls(summary_state)
1683 : Well::InjectionControls(0);
1684 auto prod_controls = well_copy.
well_ecl_.isProducer()
1685 ? well_copy.
well_ecl_.productionControls(summary_state) :
1686 Well::ProductionControls(0);
1692 const int num_perf = ws.perf_data.size();
1693 for (
int perf = 0; perf < num_perf; ++perf) {
1697 const int np = this->number_of_phases_;
1698 bool trivial =
true;
1699 for (
int phase = 0; phase < np; ++phase){
1700 trivial = trivial && (ws.well_potentials[phase] == 0.0) ;
1704 for (
int phase = 0; phase < np; ++phase) {
1705 ws.surface_rates[phase] = sign * ws.well_potentials[phase];
1710 const double dt = simulator.timeStepSize();
1712 bool converged =
false;
1713 if (this->well_ecl_.isProducer()) {
1715 simulator, dt, inj_controls, prod_controls, groupStateHelper_copy, well_state_copy, deferred_logger
1719 simulator, dt, inj_controls, prod_controls, groupStateHelper_copy, well_state_copy, deferred_logger,
1727 well_potentials.clear();
1728 well_potentials.resize(np, 0.0);
1729 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1730 if (has_solvent && comp_idx == Indices::contiSolventEqIdx)
continue;
1732 well_potentials[FluidSystem::activeCompToActivePhaseIdx(comp_idx)] = rate.value();
1736 if constexpr (has_solvent) {
1737 assert(FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx));
1739 const int gas_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
1741 well_potentials[gas_pos] += rate.value();
1747 template<
typename TypeTag>
1752 const SummaryState &summary_state,
1754 std::vector<Scalar>& potentials,
1758 auto bhp_at_thp_limit = computeBhpAtThpLimitProdWithAlq(
1759 simulator, groupStateHelper, summary_state, alq, deferred_logger,
true);
1760 if (bhp_at_thp_limit) {
1761 const auto& controls = this->well_ecl_.productionControls(summary_state);
1762 bhp = std::max(*bhp_at_thp_limit,
1763 static_cast<Scalar>(controls.bhp_limit));
1764 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1767 deferred_logger.
warning(
"FAILURE_GETTING_CONVERGED_POTENTIAL",
1768 "Failed in getting converged thp based potential calculation for well "
1769 + name() +
". Instead the bhp based value is used");
1770 const auto& controls = this->well_ecl_.productionControls(summary_state);
1771 bhp = controls.bhp_limit;
1772 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1777 template<
typename TypeTag>
1782 const SummaryState& summary_state,
1784 std::vector<Scalar>& potentials,
1788 computeWellRatesAndBhpWithThpAlqProd(simulator,
1796 template<
typename TypeTag>
1802 std::vector<Scalar>& well_potentials,
1805 const auto [compute_potential, bhp_controlled_well] =
1808 if (!compute_potential) {
1812 bool converged_implicit =
false;
1816 if (this->param_.local_well_solver_control_switching_ && !(this->changed_to_open_this_step_ && this->wellUnderZeroRateTarget(groupStateHelper, deferred_logger))) {
1817 converged_implicit = computeWellPotentialsImplicit(
1818 simulator, groupStateHelper, well_potentials, deferred_logger
1821 if (!converged_implicit) {
1823 const auto& summaryState = simulator.vanguard().summaryState();
1824 if (!Base::wellHasTHPConstraints(summaryState) || bhp_controlled_well) {
1834 const auto& ws = well_state.
well(this->index_of_well_);
1835 if (this->isInjector())
1836 bhp = std::max(ws.bhp, bhp);
1838 bhp = std::min(ws.bhp, bhp);
1840 assert(std::abs(bhp) != std::numeric_limits<Scalar>::max());
1841 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, well_potentials, deferred_logger);
1844 well_potentials = computeWellPotentialWithTHP(simulator, groupStateHelper, deferred_logger, well_state);
1848 this->checkNegativeWellPotentials(well_potentials,
1849 this->param_.check_well_operability_,
1859 template<
typename TypeTag>
1863 const int openConnIdx)
const
1865 return (openConnIdx < 0)
1867 : this->connections_.rho(openConnIdx);
1874 template<
typename TypeTag>
1880 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1882 const auto& well_state = groupStateHelper.
wellState();
1883 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(groupStateHelper, deferred_logger);
1884 this->primary_variables_.update(well_state, stop_or_zero_rate_target, deferred_logger);
1887 if constexpr (Base::has_polymermw) {
1888 this->primary_variables_.updatePolyMW(well_state);
1891 this->primary_variables_.checkFinite(deferred_logger,
"updating from well state");
1897 template<
typename TypeTag>
1902 return cachedRefDensity;
1908 template<
typename TypeTag>
1913 std::vector<EvalWell>& mob,
1916 const int cell_idx = this->well_cells_[perf];
1917 const auto& int_quant = simulator.model().intensiveQuantities(cell_idx, 0);
1918 const EvalWell polymer_concentration = this->extendEval(int_quant.polymerConcentration());
1922 if (this->isInjector()) {
1924 const auto& visc_mult_table = PolymerModule::plyviscViscosityMultiplierTable(int_quant.pvtRegionIndex());
1925 const unsigned waterCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
1926 mob[waterCompIdx] /= (this->extendEval(int_quant.waterViscosityCorrection()) * visc_mult_table.eval(polymer_concentration,
true) );
1929 if (PolymerModule::hasPlyshlog()) {
1932 if (this->isInjector() && this->wpolymer() == 0.) {
1937 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
1938 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
1940 std::vector<EvalWell> cq_s(this->num_conservation_quantities_, 0.);
1943 getTransMult(trans_mult, simulator, cell_idx);
1944 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1945 std::vector<EvalWell> Tw(this->num_conservation_quantities_, this->well_index_[perf] * trans_mult);
1946 this->getTw(Tw, perf, int_quant, trans_mult, wellstate_nupcol);
1947 computePerfRate(int_quant, mob, bhp, Tw, perf, allow_cf, cq_s,
1948 perf_rates, deferred_logger);
1950 const Scalar area = 2 * M_PI * this->perf_rep_radius_[perf] * this->perf_length_[perf];
1951 const auto& material_law_manager = simulator.problem().materialLawManager();
1952 const auto& scaled_drainage_info =
1953 material_law_manager->oilWaterScaledEpsInfoDrainage(cell_idx);
1954 const Scalar swcr = scaled_drainage_info.Swcr;
1955 const EvalWell poro = this->extendEval(int_quant.porosity());
1956 const EvalWell sw = this->extendEval(int_quant.fluidState().saturation(FluidSystem::waterPhaseIdx));
1958 const EvalWell denom = max( (area * poro * (sw - swcr)), 1e-12);
1959 const unsigned waterCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
1960 EvalWell water_velocity = cq_s[waterCompIdx] / denom * this->extendEval(int_quant.fluidState().invB(FluidSystem::waterPhaseIdx));
1962 if (PolymerModule::hasShrate()) {
1965 water_velocity *= PolymerModule::shrate( int_quant.pvtRegionIndex() ) / this->bore_diameters_[perf];
1967 const EvalWell shear_factor = PolymerModule::computeShearFactor(polymer_concentration,
1968 int_quant.pvtRegionIndex(),
1971 mob[waterCompIdx] /= shear_factor;
1975 template<
typename TypeTag>
1979 this->linSys_.extract(jacobian);
1983 template <
typename TypeTag>
1987 const int pressureVarIndex,
1988 const bool use_well_weights,
1991 this->linSys_.extractCPRPressureMatrix(jacobian,
2002 template<
typename TypeTag>
2009 if constexpr (Base::has_polymermw) {
2010 const int water_table_id = this->polymerWaterTable_();
2011 if (water_table_id <= 0) {
2013 fmt::format(
"Unused SKPRWAT table id used for well {}", name()),
2016 const auto& water_table_func = PolymerModule::getSkprwatTable(water_table_id);
2017 const EvalWell throughput_eval{throughput};
2019 EvalWell pskin_water = water_table_func.eval(throughput_eval, water_velocity);
2023 fmt::format(
"Polymermw is not activated, while injecting "
2024 "skin pressure is requested for well {}", name()),
2033 template<
typename TypeTag>
2041 if constexpr (Base::has_polymermw) {
2042 const Scalar sign = water_velocity >= 0. ? 1.0 : -1.0;
2043 const EvalWell water_velocity_abs = abs(water_velocity);
2044 if (poly_inj_conc == 0.) {
2045 return sign * pskinwater(throughput, water_velocity_abs, deferred_logger);
2047 const int polymer_table_id = this->polymerTable_();
2048 if (polymer_table_id <= 0) {
2050 fmt::format(
"Unavailable SKPRPOLY table id used for well {}", name()),
2053 const auto& skprpolytable = PolymerModule::getSkprpolyTable(polymer_table_id);
2054 const Scalar reference_concentration = skprpolytable.refConcentration;
2055 const EvalWell throughput_eval{throughput};
2057 const EvalWell pskin_poly = skprpolytable.table_func.eval(throughput_eval, water_velocity_abs);
2058 if (poly_inj_conc == reference_concentration) {
2059 return sign * pskin_poly;
2062 const EvalWell pskin_water = pskinwater(throughput, water_velocity_abs, deferred_logger);
2063 const EvalWell pskin = pskin_water + (pskin_poly - pskin_water) / reference_concentration * poly_inj_conc;
2064 return sign * pskin;
2067 fmt::format(
"Polymermw is not activated, while injecting "
2068 "skin pressure is requested for well {}", name()),
2077 template<
typename TypeTag>
2084 if constexpr (Base::has_polymermw) {
2085 const int table_id = this->polymerInjTable_();
2086 const auto& table_func = PolymerModule::getPlymwinjTable(table_id);
2087 const EvalWell throughput_eval{throughput};
2089 if (this->wpolymer() == 0.) {
2090 return molecular_weight;
2092 molecular_weight = table_func.eval(throughput_eval, abs(water_velocity));
2093 return molecular_weight;
2096 fmt::format(
"Polymermw is not activated, while injecting "
2097 "polymer molecular weight is requested for well {}", name()),
2106 template<
typename TypeTag>
2112 if constexpr (Base::has_polymermw) {
2113 if (!this->isInjector()) {
2117 auto& perf_water_throughput = well_state.
well(this->index_of_well_)
2118 .perf_data.water_throughput;
2120 for (
int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
2121 const Scalar perf_water_vel =
2122 this->primary_variables_.value(Bhp + 1 + perf);
2126 if (perf_water_vel >
Scalar{0}) {
2127 perf_water_throughput[perf] += perf_water_vel * dt;
2137 template<
typename TypeTag>
2142 std::vector<EvalWell>& cq_s)
const
2144 const int cell_idx = this->well_cells_[perf];
2145 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2146 const auto& fs = int_quants.fluidState();
2147 const EvalWell b_w = this->extendEval(fs.invB(FluidSystem::waterPhaseIdx));
2148 const Scalar area = M_PI * this->bore_diameters_[perf] * this->perf_length_[perf];
2149 const int wat_vel_index = Bhp + 1 + perf;
2150 const unsigned water_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
2154 cq_s[water_comp_idx] = area * this->primary_variables_.eval(wat_vel_index) * b_w;
2160 template<
typename TypeTag>
2169 const int cell_idx = this->well_cells_[perf];
2170 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2171 const auto& fs = int_quants.fluidState();
2172 const EvalWell b_w = this->extendEval(fs.invB(FluidSystem::waterPhaseIdx));
2173 const EvalWell water_flux_r = water_flux_s / b_w;
2174 const Scalar area = M_PI * this->bore_diameters_[perf] * this->perf_length_[perf];
2175 const EvalWell water_velocity = water_flux_r / area;
2176 const int wat_vel_index = Bhp + 1 + perf;
2179 const EvalWell eq_wat_vel = this->primary_variables_.eval(wat_vel_index) - water_velocity;
2181 const auto& ws = well_state.
well(this->index_of_well_);
2182 const auto& perf_data = ws.perf_data;
2183 const auto& perf_water_throughput = perf_data.water_throughput;
2184 const Scalar throughput = perf_water_throughput[perf];
2185 const int pskin_index = Bhp + 1 + this->number_of_local_perforations_ + perf;
2187 const EvalWell poly_conc(this->wpolymer());
2190 const EvalWell eq_pskin = this->primary_variables_.eval(pskin_index)
2191 - pskin(throughput, this->primary_variables_.eval(wat_vel_index), poly_conc, deferred_logger);
2194 assembleInjectivityEq(eq_pskin,
2199 this->primary_variables_.numWellEq(),
2207 template<
typename TypeTag>
2216 if constexpr (Base::has_polymermw) {
2218 checkConvergencePolyMW(res, Bhp, this->param_.max_residual_allowed_, report);
2226 template<
typename TypeTag>
2233 std::vector<RateVector>& connectionRates,
2238 if (this->isInjector()) {
2239 const int wat_vel_index = Bhp + 1 + perf;
2240 const EvalWell water_velocity = this->primary_variables_.eval(wat_vel_index);
2241 if (water_velocity > 0.) {
2242 const auto& ws = well_state.
well(this->index_of_well_);
2243 const auto& perf_water_throughput = ws.perf_data.water_throughput;
2244 const Scalar throughput = perf_water_throughput[perf];
2245 const EvalWell molecular_weight = wpolymermw(throughput, water_velocity, deferred_logger);
2246 cq_s_polymw *= molecular_weight;
2252 }
else if (this->isProducer()) {
2253 if (cq_s_polymw < 0.) {
2254 cq_s_polymw *= this->extendEval(int_quants.polymerMoleWeight() );
2261 connectionRates[perf][Indices::contiPolymerMWEqIdx] = Base::restrictEval(cq_s_polymw);
2268 template<
typename TypeTag>
2269 std::optional<typename StandardWell<TypeTag>::Scalar>
2274 const SummaryState& summary_state,
2277 return computeBhpAtThpLimitProdWithAlq(simulator,
2280 this->getALQ(well_state),
2285 template<
typename TypeTag>
2286 std::optional<typename StandardWell<TypeTag>::Scalar>
2290 const SummaryState& summary_state,
2293 bool iterate_if_no_solution)
const
2297 auto frates = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2303 std::vector<Scalar> rates(3);
2304 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2305 this->adaptRatesForVFP(rates);
2310 maxPerfPress(simulator),
2311 this->getRefDensity(),
2313 this->getTHPConstraint(summary_state),
2317 auto v = frates(*bhpAtLimit);
2318 if (std::all_of(v.cbegin(), v.cend(), [](
Scalar i){ return i <= 0; }) ) {
2323 if (!iterate_if_no_solution)
2324 return std::nullopt;
2326 auto fratesIter = [
this, &simulator, &groupStateHelper, &deferred_logger](
const Scalar bhp) {
2330 std::vector<Scalar> rates(3);
2331 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, rates, deferred_logger);
2332 this->adaptRatesForVFP(rates);
2338 maxPerfPress(simulator),
2339 this->getRefDensity(),
2341 this->getTHPConstraint(summary_state),
2347 auto v = frates(*bhpAtLimit);
2348 if (std::all_of(v.cbegin(), v.cend(), [](
Scalar i){ return i <= 0; }) ) {
2354 return std::nullopt;
2359 template<
typename TypeTag>
2360 std::optional<typename StandardWell<TypeTag>::Scalar>
2364 const SummaryState& summary_state,
2373 auto frates = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2379 std::vector<Scalar> rates(3);
2380 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2386 this->getRefDensity(),
2397 template<
typename TypeTag>
2402 const Well::InjectionControls& inj_controls,
2403 const Well::ProductionControls& prod_controls,
2408 updatePrimaryVariables(groupStateHelper, deferred_logger);
2410 const int max_iter = this->param_.max_inner_iter_wells_;
2413 bool relax_convergence =
false;
2414 this->regularize_ =
false;
2416 assembleWellEqWithoutIteration(simulator, groupStateHelper, dt, inj_controls, prod_controls, well_state, deferred_logger,
2419 if (it > this->param_.strict_inner_iter_wells_) {
2420 relax_convergence =
true;
2421 this->regularize_ =
true;
2424 auto report = getWellConvergence(groupStateHelper, Base::B_avg_, deferred_logger, relax_convergence);
2426 converged = report.converged();
2432 solveEqAndUpdateWellState(simulator, groupStateHelper, well_state, deferred_logger);
2439 }
while (it < max_iter);
2442 std::ostringstream sstr;
2443 sstr <<
" Well " << this->name() <<
" converged in " << it <<
" inner iterations.";
2444 if (relax_convergence)
2445 sstr <<
" (A relaxed tolerance was used after "<< this->param_.strict_inner_iter_wells_ <<
" iterations)";
2449 deferred_logger.
debug(sstr.str(), OpmLog::defaultDebugVerbosityLevel + (it == 0));
2451 std::ostringstream sstr;
2452 sstr <<
" Well " << this->name() <<
" did not converge in " << it <<
" inner iterations.";
2453 deferred_logger.
debug(sstr.str());
2460 template<
typename TypeTag>
2465 const Well::InjectionControls& inj_controls,
2466 const Well::ProductionControls& prod_controls,
2470 const bool fixed_control ,
2471 const bool fixed_status ,
2472 const bool solving_with_zero_rate )
2474 updatePrimaryVariables(groupStateHelper, deferred_logger);
2476 const int max_iter = this->param_.max_inner_iter_wells_;
2478 bool converged =
false;
2479 bool relax_convergence =
false;
2480 this->regularize_ =
false;
2481 const auto& summary_state = groupStateHelper.
summaryState();
2486 constexpr int min_its_after_switch = 4;
2488 const int max_status_switch = this->param_.max_well_status_switch_inner_iter_;
2489 int its_since_last_switch = min_its_after_switch;
2490 int switch_count= 0;
2492 const auto well_status_orig = this->wellStatus_;
2493 const auto operability_orig = this->operability_status_;
2494 auto well_status_cur = well_status_orig;
2495 int status_switch_count = 0;
2497 const bool allow_open = well_state.
well(this->index_of_well_).status == WellStatus::OPEN;
2499 const bool allow_switching =
2500 !this->wellUnderZeroRateTarget(groupStateHelper, deferred_logger) &&
2501 (!fixed_control || !fixed_status) && allow_open;
2503 bool changed =
false;
2504 bool final_check =
false;
2506 this->operability_status_.resetOperability();
2507 this->operability_status_.solvable =
true;
2509 its_since_last_switch++;
2510 if (allow_switching && its_since_last_switch >= min_its_after_switch && status_switch_count < max_status_switch){
2511 const Scalar wqTotal = this->primary_variables_.eval(WQTotal).value();
2512 changed = this->updateWellControlAndStatusLocalIteration(
2513 simulator, groupStateHelper, inj_controls, prod_controls, wqTotal,
2514 well_state, deferred_logger, fixed_control, fixed_status,
2515 solving_with_zero_rate
2518 its_since_last_switch = 0;
2520 if (well_status_cur != this->wellStatus_) {
2521 well_status_cur = this->wellStatus_;
2522 status_switch_count++;
2525 if (!changed && final_check) {
2528 final_check =
false;
2530 if (status_switch_count == max_status_switch) {
2531 this->wellStatus_ = well_status_orig;
2535 assembleWellEqWithoutIteration(simulator, groupStateHelper, dt, inj_controls, prod_controls, well_state, deferred_logger, solving_with_zero_rate);
2537 if (it > this->param_.strict_inner_iter_wells_) {
2538 relax_convergence =
true;
2539 this->regularize_ =
true;
2542 auto report = getWellConvergence(groupStateHelper, Base::B_avg_, deferred_logger, relax_convergence);
2544 converged = report.converged();
2548 if (switch_count > 0 && its_since_last_switch < min_its_after_switch) {
2550 its_since_last_switch = min_its_after_switch;
2557 solveEqAndUpdateWellState(simulator, groupStateHelper, well_state, deferred_logger);
2559 }
while (it < max_iter);
2562 if (allow_switching){
2564 const bool is_stopped = this->wellIsStopped();
2565 if (this->wellHasTHPConstraints(summary_state)){
2566 this->operability_status_.can_obtain_bhp_with_thp_limit = !is_stopped;
2567 this->operability_status_.obey_thp_limit_under_bhp_limit = !is_stopped;
2569 this->operability_status_.operable_under_only_bhp_limit = !is_stopped;
2572 std::string message = fmt::format(
" Well {} converged in {} inner iterations ("
2573 "{} control/status switches).", this->name(), it, switch_count);
2574 if (relax_convergence) {
2575 message.append(fmt::format(
" (A relaxed tolerance was used after {} iterations)",
2576 this->param_.strict_inner_iter_wells_));
2578 deferred_logger.
debug(message, OpmLog::defaultDebugVerbosityLevel + ((it == 0) && (switch_count == 0)));
2581 this->wellStatus_ = well_status_orig;
2582 this->operability_status_ = operability_orig;
2583 const std::string message = fmt::format(
" Well {} did not converge in {} inner iterations ("
2584 "{} switches, {} status changes).", this->name(), it, switch_count, status_switch_count);
2585 deferred_logger.
debug(message);
2591 template<
typename TypeTag>
2592 std::vector<typename StandardWell<TypeTag>::Scalar>
2598 std::vector<Scalar> well_q_s(this->num_conservation_quantities_, 0.);
2599 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
2600 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
2601 for (
int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
2602 const int cell_idx = this->well_cells_[perf];
2603 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
2604 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.);
2605 getMobility(simulator, perf, mob, deferred_logger);
2606 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.);
2608 getTransMult(trans_mult, simulator, cell_idx);
2609 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
2610 std::vector<Scalar> Tw(this->num_conservation_quantities_, this->well_index_[perf] * trans_mult);
2611 this->getTw(Tw, perf, intQuants, trans_mult, wellstate_nupcol);
2613 computePerfRate(intQuants, mob, bhp.value(), Tw, perf, allow_cf,
2614 cq_s, perf_rates, deferred_logger);
2615 for (
int comp = 0; comp < this->num_conservation_quantities_; ++comp) {
2616 well_q_s[comp] += cq_s[comp];
2619 const auto& comm = this->parallel_well_info_.communication();
2620 if (comm.size() > 1)
2622 comm.sum(well_q_s.data(), well_q_s.size());
2629 template <
typename TypeTag>
2630 std::vector<typename StandardWell<TypeTag>::Scalar>
2634 const int num_pri_vars = this->primary_variables_.numWellEq();
2635 std::vector<Scalar> retval(num_pri_vars);
2636 for (
int ii = 0; ii < num_pri_vars; ++ii) {
2637 retval[ii] = this->primary_variables_.value(ii);
2646 template <
typename TypeTag>
2651 const int num_pri_vars = this->primary_variables_.numWellEq();
2652 for (
int ii = 0; ii < num_pri_vars; ++ii) {
2653 this->primary_variables_.setValue(ii, it[ii]);
2655 return num_pri_vars;
2659 template <
typename TypeTag>
2663 const IntensiveQuantities& intQuants,
2666 auto fs = intQuants.fluidState();
2668 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
2669 if (!FluidSystem::phaseIsActive(phaseIdx)) {
2674 EvalWell cq_r_thermal{0.};
2675 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
2676 const bool both_oil_gas = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
2677 if (!both_oil_gas || FluidSystem::waterPhaseIdx == phaseIdx) {
2678 cq_r_thermal = cq_s[activeCompIdx] / this->extendEval(fs.invB(phaseIdx));
2681 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
2682 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
2687 const EvalWell d = this->extendEval(1.0 - fs.Rv() * fs.Rs());
2689 deferred_logger.
debug(
2690 fmt::format(
"Problematic d value {} obtained for well {}"
2691 " during calculateSinglePerf with rs {}"
2692 ", rv {}. Continue as if no dissolution (rs = 0) and"
2693 " vaporization (rv = 0) for this connection.",
2694 d, this->name(), fs.Rs(), fs.Rv()));
2695 cq_r_thermal = cq_s[activeCompIdx] / this->extendEval(fs.invB(phaseIdx));
2697 if (FluidSystem::gasPhaseIdx == phaseIdx) {
2698 cq_r_thermal = (cq_s[gasCompIdx] -
2699 this->extendEval(fs.Rs()) * cq_s[oilCompIdx]) /
2700 (d * this->extendEval(fs.invB(phaseIdx)) );
2701 }
else if (FluidSystem::oilPhaseIdx == phaseIdx) {
2703 cq_r_thermal = (cq_s[oilCompIdx] - this->extendEval(fs.Rv()) *
2705 (d * this->extendEval(fs.invB(phaseIdx)) );
2711 if (this->isInjector() && !this->wellIsStopped() && cq_r_thermal > 0.0){
2713 assert(this->well_ecl_.injectorType() != InjectorType::MULTI);
2714 fs.setTemperature(this->well_ecl_.inj_temperature());
2715 typedef typename std::decay<
decltype(fs)>::type::Scalar FsScalar;
2716 typename FluidSystem::template ParameterCache<FsScalar> paramCache;
2717 const unsigned pvtRegionIdx = intQuants.pvtRegionIndex();
2718 paramCache.setRegionIndex(pvtRegionIdx);
2719 paramCache.updatePhase(fs, phaseIdx);
2721 const auto& rho = FluidSystem::density(fs, paramCache, phaseIdx);
2722 fs.setDensity(phaseIdx, rho);
2723 const auto& h = FluidSystem::enthalpy(fs, paramCache, phaseIdx);
2724 fs.setEnthalpy(phaseIdx, h);
2725 cq_r_thermal *= this->extendEval(fs.enthalpy(phaseIdx)) * this->extendEval(fs.density(phaseIdx));
2726 result += getValue(cq_r_thermal);
2727 }
else if (cq_r_thermal > 0.0) {
2728 cq_r_thermal *= getValue(fs.enthalpy(phaseIdx)) * getValue(fs.density(phaseIdx));
2729 result += Base::restrictEval(cq_r_thermal);
2732 cq_r_thermal *= this->extendEval(fs.enthalpy(phaseIdx)) * this->extendEval(fs.density(phaseIdx));
2733 result += Base::restrictEval(cq_r_thermal);
2737 return result * this->well_efficiency_factor_;
2740 template <
typename TypeTag>
2744 Scalar max_pressure = 0.0;
2745 for (
int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
2746 const int cell_idx = this->well_cells_[perf];
2747 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2748 const auto& fs = int_quants.fluidState();
2749 Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2750 max_pressure = std::max(max_pressure, pressure_cell);
2752 const auto& comm = this->parallel_well_info_.communication();
2753 if (comm.size() > 1) {
2754 max_pressure = comm.max(max_pressure);
2756 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 warning(const std::string &tag, const std::string &message)
void debug(const std::string &tag, const std::string &message)
Definition: GroupStateHelper.hpp:53
GroupState< Scalar > & groupState() const
Definition: GroupStateHelper.hpp:180
const SummaryState & summaryState() const
Definition: GroupStateHelper.hpp:257
const WellState< Scalar, IndexTraits > & wellState() const
Definition: GroupStateHelper.hpp:306
WellStateGuard pushWellState(WellState< Scalar, IndexTraits > &well_state)
Definition: GroupStateHelper.hpp:195
GroupStateGuard pushGroupState(GroupState< Scalar > &group_state)
Definition: GroupStateHelper.hpp:190
Definition: GroupState.hpp:41
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:198
Definition: RatioCalculator.hpp:38
Class handling assemble of the equation system for StandardWell.
Definition: StandardWellAssemble.hpp:44
Scalar pressure_diff(const unsigned perf) const
Returns pressure drop for a given perforation.
Definition: StandardWellConnections.hpp:101
StdWellConnections connections_
Connection level values.
Definition: StandardWellEval.hpp:101
PrimaryVariables primary_variables_
Primary variables for well.
Definition: StandardWellEval.hpp:95
Definition: StandardWell.hpp:60
void updatePrimaryVariables(const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1877
EvalWell wpolymermw(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2080
void updateIPRImplicit(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:921
void assembleWellEqWithoutIteration(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, 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 solveEqAndUpdateWellState(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1415
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::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger, bool iterate_if_no_solution) const override
Definition: StandardWell_impl.hpp:2288
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2362
void addWellContributions(SparseMatrixAdapter &mat) const override
Definition: StandardWell_impl.hpp:1977
Scalar computeWellRatesAndBhpWithThpAlqProd(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1750
std::vector< Scalar > getPrimaryVars() const override
Definition: StandardWell_impl.hpp:2632
virtual ConvergenceReport getWellConvergence(const GroupStateHelperType &groupStateHelper, 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 computeWellRatesWithBhpIterations(const Simulator &ebosSimulator, const Scalar &bhp, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1553
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:1985
void updateWaterMobilityWithPolymer(const Simulator &simulator, const int perf, std::vector< EvalWell > &mob_water, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1911
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1484
std::vector< Scalar > computeCurrentWellRates(const Simulator &ebosSimulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:2594
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:505
void computeWellConnectionDensitesPressures(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const WellConnectionProps &props, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1339
std::vector< Scalar > computeWellPotentialWithTHP(const Simulator &ebosSimulator, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger, const WellStateType &well_state) const
Definition: StandardWell_impl.hpp:1624
bool computeWellPotentialsImplicit(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1660
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
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:52
typename StdWellEval::StdWellConnections::Properties WellConnectionProps
Definition: StandardWell.hpp:274
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:2229
void computeWellRatesWithBhp(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1505
void computeWellConnectionPressures(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1396
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:703
void getTransMult(Value &trans_mult, const Simulator &simulator, const int cell_indx) const
Definition: StandardWell_impl.hpp:683
void checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1067
void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:829
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:759
void handleInjectivityEquations(const Simulator &simulator, const WellStateType &well_state, const int perf, const EvalWell &water_flux_s, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:2163
void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) override
computing the well potentials for group control
Definition: StandardWell_impl.hpp:1799
virtual void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: StandardWell_impl.hpp:1452
void checkConvergenceExtraEqs(const std::vector< Scalar > &res, ConvergenceReport &report) const
Definition: StandardWell_impl.hpp:2210
void assembleWellEqWithoutIterationImpl(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, 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
typename StdWellEval::Eval Eval
Definition: StandardWell.hpp:120
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1159
Scalar maxPerfPress(const Simulator &simulator) const override
Definition: StandardWell_impl.hpp:2743
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:2036
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:2649
void computeWellRatesWithThpAlqProd(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1780
void updateWaterThroughput(const double dt, WellStateType &well_state) const override
Definition: StandardWell_impl.hpp:2109
void checkOperabilityUnderBHPLimit(const WellStateType &well_state, const Simulator &simulator, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:998
EvalWell pskinwater(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2005
bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:2400
void calculateExplicitQuantities(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1438
void handleInjectivityRate(const Simulator &simulator, const int perf, std::vector< EvalWell > &cq_s) const
Definition: StandardWell_impl.hpp:2140
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:1264
std::optional< Scalar > computeBhpAtThpLimitProd(const WellStateType &well_state, const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2271
Scalar getRefDensity() const override
Definition: StandardWell_impl.hpp:1900
bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, 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:2463
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: StandardWell_impl.hpp:1862
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: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
bool solveWellWithOperabilityCheck(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:589
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:2060
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:2073
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: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