22#ifndef OPM_STANDARDWELL_IMPL_HEADER_INCLUDED
23#define OPM_STANDARDWELL_IMPL_HEADER_INCLUDED
26#ifndef OPM_STANDARDWELL_HEADER_INCLUDED
28#include <opm/simulators/wells/StandardWell.hpp>
31#include <opm/common/Exceptions.hpp>
33#include <opm/input/eclipse/Units/Units.hpp>
35#include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>
36#include <opm/simulators/wells/StandardWellAssemble.hpp>
37#include <opm/simulators/wells/VFPHelpers.hpp>
38#include <opm/simulators/wells/WellBhpThpCalculator.hpp>
39#include <opm/simulators/wells/WellConvergence.hpp>
45#include <fmt/format.h>
50 template<
typename TypeTag>
55 const ModelParameters& param,
56 const RateConverterType& rate_converter,
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)
66 assert(this->num_conservation_quantities_ == numWellConservationEq);
73 template<
typename TypeTag>
75 StandardWell<TypeTag>::
76 init(
const std::vector<Scalar>& depth_arg,
77 const Scalar gravity_arg,
78 const std::vector< Scalar >& B_avg,
79 const bool changed_to_open_this_step)
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>
92 StandardWell<TypeTag>::
93 computePerfRate(
const IntensiveQuantities& intQuants,
94 const std::vector<Value>& mob,
96 const std::vector<Scalar>& Tw,
99 std::vector<Value>& cq_s,
100 PerforationRates<Scalar>& perf_rates,
101 DeferredLogger& deferred_logger)
const
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>
194 StandardWell<TypeTag>::
195 computePerfRate(
const std::vector<Value>& mob,
196 const Value& pressure,
202 std::vector<Value>& b_perfcells_dense,
203 const std::vector<Scalar>& Tw,
206 const Value& skin_pressure,
207 const std::vector<Value>& cmix_s,
208 std::vector<Value>& cq_s,
209 PerforationRates<Scalar>& perf_rates,
210 DeferredLogger& deferred_logger)
const
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;
219 RatioCalculator<Value> ratioCalc{
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>
340 StandardWell<TypeTag>::
341 assembleWellEqWithoutIteration(
const Simulator& simulator,
343 const Well::InjectionControls& inj_controls,
344 const Well::ProductionControls& prod_controls,
345 WellStateType& well_state,
346 const GroupState<Scalar>& group_state,
347 DeferredLogger& deferred_logger)
351 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
354 this->linSys_.clear();
356 assembleWellEqWithoutIterationImpl(simulator, dt, inj_controls,
357 prod_controls, well_state,
358 group_state, deferred_logger);
364 template<
typename TypeTag>
366 StandardWell<TypeTag>::
367 assembleWellEqWithoutIterationImpl(
const Simulator& simulator,
369 const Well::InjectionControls& inj_controls,
370 const Well::ProductionControls& prod_controls,
371 WellStateType& well_state,
372 const GroupState<Scalar>& group_state,
373 DeferredLogger& deferred_logger)
376 const Scalar regularization_factor = this->regularize_? this->param_.regularization_factor_wells_ : 1.0;
377 const Scalar volume = 0.1 * unit::cubic(unit::feet) * regularization_factor;
379 auto& ws = well_state.well(this->index_of_well_);
380 ws.phase_mixing_rates.fill(0.0);
383 const int np = this->number_of_phases_;
385 std::vector<RateVector> connectionRates = this->connectionRates_;
387 auto& perf_data = ws.perf_data;
388 auto& perf_rates = perf_data.phase_rates;
389 for (
int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
391 std::vector<EvalWell> cq_s(this->num_conservation_quantities_, {this->primary_variables_.numWellEq() + Indices::numEq, 0.0});
392 EvalWell water_flux_s{this->primary_variables_.numWellEq() + Indices::numEq, 0.0};
393 EvalWell cq_s_zfrac_effective{this->primary_variables_.numWellEq() + Indices::numEq, 0.0};
394 calculateSinglePerf(simulator, perf, well_state, connectionRates,
395 cq_s, water_flux_s, cq_s_zfrac_effective, deferred_logger);
398 if constexpr (has_polymer && Base::has_polymermw) {
399 if (this->isInjector()) {
400 handleInjectivityEquations(simulator, well_state, perf,
401 water_flux_s, deferred_logger);
404 for (
int componentIdx = 0; componentIdx < this->num_conservation_quantities_; ++componentIdx) {
406 const EvalWell cq_s_effective = cq_s[componentIdx] * this->well_efficiency_factor_;
408 connectionRates[perf][componentIdx] = Base::restrictEval(cq_s_effective);
410 StandardWellAssemble<FluidSystem,Indices>(*this).
411 assemblePerforationEq(cq_s_effective,
414 this->primary_variables_.numWellEq(),
418 if (has_solvent && componentIdx == Indices::contiSolventEqIdx) {
419 auto& perf_rate_solvent = perf_data.solvent_rates;
420 perf_rate_solvent[perf] = cq_s[componentIdx].value();
422 perf_rates[perf*np + FluidSystem::activeCompToActivePhaseIdx(componentIdx)] = cq_s[componentIdx].value();
426 if constexpr (has_zFraction) {
427 StandardWellAssemble<FluidSystem,Indices>(*this).
428 assembleZFracEq(cq_s_zfrac_effective,
430 this->primary_variables_.numWellEq(),
435 this->connectionRates_ = connectionRates;
440 const auto& comm = this->parallel_well_info_.communication();
441 comm.sum(ws.phase_mixing_rates.data(), ws.phase_mixing_rates.size());
445 this->linSys_.sumDistributed(this->parallel_well_info_.communication());
448 for (
int componentIdx = 0; componentIdx < numWellConservationEq; ++componentIdx) {
451 EvalWell resWell_loc(this->primary_variables_.numWellEq() + Indices::numEq, 0.0);
452 if (FluidSystem::numActivePhases() > 1) {
454 resWell_loc += (this->primary_variables_.surfaceVolumeFraction(componentIdx) -
455 this->F0_[componentIdx]) * volume / dt;
457 resWell_loc -= this->primary_variables_.getQs(componentIdx) * this->well_efficiency_factor_;
458 StandardWellAssemble<FluidSystem,Indices>(*this).
459 assembleSourceEq(resWell_loc,
461 this->primary_variables_.numWellEq(),
465 const auto& summaryState = simulator.vanguard().summaryState();
466 const Schedule& schedule = simulator.vanguard().schedule();
467 const bool stopped_or_zero_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
468 StandardWellAssemble<FluidSystem,Indices>(*this).
469 assembleControlEq(well_state, group_state,
470 schedule, summaryState,
471 inj_controls, prod_controls,
472 this->primary_variables_,
473 this->getRefDensity(),
475 stopped_or_zero_target,
481 this->linSys_.invert();
483 OPM_DEFLOG_PROBLEM(NumericalProblem,
"Error when inverting local well equations for well " + name(), deferred_logger);
490 template<
typename TypeTag>
492 StandardWell<TypeTag>::
493 calculateSinglePerf(
const Simulator& simulator,
495 WellStateType& well_state,
496 std::vector<RateVector>& connectionRates,
497 std::vector<EvalWell>& cq_s,
498 EvalWell& water_flux_s,
499 EvalWell& cq_s_zfrac_effective,
500 DeferredLogger& deferred_logger)
const
502 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
503 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
504 const int cell_idx = this->well_cells_[perf];
505 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
506 std::vector<EvalWell> mob(this->num_conservation_quantities_, {this->primary_variables_.numWellEq() + Indices::numEq, 0.});
507 getMobility(simulator, perf, mob, deferred_logger);
509 PerforationRates<Scalar> perf_rates;
510 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(intQuants, cell_idx);
511 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
512 const std::vector<Scalar> Tw = this->wellIndex(perf, intQuants, trans_mult, wellstate_nupcol);
513 computePerfRate(intQuants, mob, bhp, Tw, perf, allow_cf,
514 cq_s, perf_rates, deferred_logger);
516 auto& ws = well_state.well(this->index_of_well_);
517 auto& perf_data = ws.perf_data;
518 if constexpr (has_polymer && Base::has_polymermw) {
519 if (this->isInjector()) {
522 const unsigned water_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
523 water_flux_s = cq_s[water_comp_idx];
526 handleInjectivityRate(simulator, perf, cq_s);
531 if (this->isProducer()) {
532 ws.phase_mixing_rates[ws.dissolved_gas] += perf_rates.dis_gas;
533 ws.phase_mixing_rates[ws.dissolved_gas_in_water] += perf_rates.dis_gas_in_water;
534 ws.phase_mixing_rates[ws.vaporized_oil] += perf_rates.vap_oil;
535 ws.phase_mixing_rates[ws.vaporized_water] += perf_rates.vap_wat;
536 perf_data.phase_mixing_rates[perf][ws.dissolved_gas] = perf_rates.dis_gas;
537 perf_data.phase_mixing_rates[perf][ws.dissolved_gas_in_water] = perf_rates.dis_gas_in_water;
538 perf_data.phase_mixing_rates[perf][ws.vaporized_oil] = perf_rates.vap_oil;
539 perf_data.phase_mixing_rates[perf][ws.vaporized_water] = perf_rates.vap_wat;
542 if constexpr (has_energy) {
543 connectionRates[perf][Indices::contiEnergyEqIdx] =
544 connectionRateEnergy(cq_s, intQuants, deferred_logger);
547 if constexpr (has_polymer) {
548 std::variant<Scalar,EvalWell> polymerConcentration;
549 if (this->isInjector()) {
550 polymerConcentration = this->wpolymer();
552 polymerConcentration = this->extendEval(intQuants.polymerConcentration() *
553 intQuants.polymerViscosityCorrection());
556 [[maybe_unused]] EvalWell cq_s_poly;
557 std::tie(connectionRates[perf][Indices::contiPolymerEqIdx],
559 this->connections_.connectionRatePolymer(perf_data.polymer_rates[perf],
560 cq_s, polymerConcentration);
562 if constexpr (Base::has_polymermw) {
563 updateConnectionRatePolyMW(cq_s_poly, intQuants, well_state,
564 perf, connectionRates, deferred_logger);
568 if constexpr (has_foam) {
569 std::variant<Scalar,EvalWell> foamConcentration;
570 if (this->isInjector()) {
571 foamConcentration = this->wfoam();
573 foamConcentration = this->extendEval(intQuants.foamConcentration());
575 connectionRates[perf][Indices::contiFoamEqIdx] =
576 this->connections_.connectionRateFoam(cq_s, foamConcentration,
577 FoamModule::transportPhase(),
581 if constexpr (has_zFraction) {
582 std::variant<Scalar,std::array<EvalWell,2>> solventConcentration;
583 if (this->isInjector()) {
584 solventConcentration = this->wsolvent();
586 solventConcentration = std::array{this->extendEval(intQuants.xVolume()),
587 this->extendEval(intQuants.yVolume())};
589 std::tie(connectionRates[perf][Indices::contiZfracEqIdx],
590 cq_s_zfrac_effective) =
591 this->connections_.connectionRatezFraction(perf_data.solvent_rates[perf],
592 perf_rates.dis_gas, cq_s,
593 solventConcentration);
596 if constexpr (has_brine) {
597 std::variant<Scalar,EvalWell> saltConcentration;
598 if (this->isInjector()) {
599 saltConcentration = this->wsalt();
601 saltConcentration = this->extendEval(intQuants.fluidState().saltConcentration());
604 connectionRates[perf][Indices::contiBrineEqIdx] =
605 this->connections_.connectionRateBrine(perf_data.brine_rates[perf],
606 perf_rates.vap_wat, cq_s,
610 if constexpr (has_bioeffects) {
611 std::variant<Scalar,EvalWell> microbialConcentration;
612 if constexpr (has_micp) {
613 std::variant<Scalar,EvalWell> oxygenConcentration;
614 std::variant<Scalar,EvalWell> ureaConcentration;
615 if (this->isInjector()) {
616 microbialConcentration = this->wmicrobes();
617 oxygenConcentration = this->woxygen();
618 ureaConcentration = this->wurea();
620 microbialConcentration = this->extendEval(intQuants.microbialConcentration());
621 oxygenConcentration = this->extendEval(intQuants.oxygenConcentration());
622 ureaConcentration = this->extendEval(intQuants.ureaConcentration());
624 std::tie(connectionRates[perf][Indices::contiMicrobialEqIdx],
625 connectionRates[perf][Indices::contiOxygenEqIdx],
626 connectionRates[perf][Indices::contiUreaEqIdx]) =
627 this->connections_.connectionRatesMICP(perf_data.microbial_rates[perf],
628 perf_data.oxygen_rates[perf],
629 perf_data.urea_rates[perf],
631 microbialConcentration,
636 if (this->isProducer()) {
637 microbialConcentration = this->extendEval(intQuants.microbialConcentration());
638 connectionRates[perf][Indices::contiMicrobialEqIdx] =
639 this->connections_.connectionRateBioeffects(perf_data.microbial_rates[perf],
640 perf_rates.vap_wat, cq_s,
641 microbialConcentration);
647 perf_data.pressure[perf] = ws.bhp + this->connections_.pressure_diff(perf);
650 if (FluidSystem::phaseUsage().hasCO2orH2Store()) {
651 const unsigned gas_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
652 const Scalar rho = FluidSystem::referenceDensity( FluidSystem::gasPhaseIdx, Base::pvtRegionIdx() );
653 perf_data.gas_mass_rates[perf] = cq_s[gas_comp_idx].value() * rho;
657 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
658 const unsigned wat_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
659 const Scalar rho = FluidSystem::referenceDensity( FluidSystem::waterPhaseIdx, Base::pvtRegionIdx() );
660 perf_data.wat_mass_rates[perf] = cq_s[wat_comp_idx].value() * rho;
666 template<
typename TypeTag>
667 template<
class Value>
669 StandardWell<TypeTag>::
670 getMobility(
const Simulator& simulator,
672 std::vector<Value>& mob,
673 DeferredLogger& deferred_logger)
const
675 auto obtain = [
this](
const Eval& value)
677 if constexpr (std::is_same_v<Value, Scalar>) {
678 static_cast<void>(
this);
679 return getValue(value);
681 return this->extendEval(value);
684 WellInterface<TypeTag>::getMobility(simulator, perf, mob,
685 obtain, deferred_logger);
688 if constexpr (has_polymer) {
689 if (!FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
690 OPM_DEFLOG_THROW(std::runtime_error,
"Water is required when polymer is active", deferred_logger);
695 if constexpr (!Base::has_polymermw) {
696 if constexpr (std::is_same_v<Value, Scalar>) {
697 std::vector<EvalWell> mob_eval(this->num_conservation_quantities_, {this->primary_variables_.numWellEq() + Indices::numEq, 0.});
698 for (std::size_t i = 0; i < mob.size(); ++i) {
699 mob_eval[i].setValue(mob[i]);
701 updateWaterMobilityWithPolymer(simulator, perf, mob_eval, deferred_logger);
702 for (std::size_t i = 0; i < mob.size(); ++i) {
703 mob[i] = getValue(mob_eval[i]);
706 updateWaterMobilityWithPolymer(simulator, perf, mob, deferred_logger);
712 if (this->isInjector() && this->well_ecl_.getInjMultMode() != Well::InjMultMode::NONE) {
713 const Scalar bhp = this->primary_variables_.value(Bhp);
714 const Scalar perf_press = bhp + this->connections_.pressure_diff(perf);
715 const Scalar multiplier = this->getInjMult(perf, bhp, perf_press, deferred_logger);
716 for (std::size_t i = 0; i < mob.size(); ++i) {
717 mob[i] *= multiplier;
723 template<
typename TypeTag>
725 StandardWell<TypeTag>::
726 updateWellState(
const Simulator& simulator,
727 const BVectorWell& dwells,
728 WellStateType& well_state,
729 DeferredLogger& deferred_logger)
731 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
733 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
734 updatePrimaryVariablesNewton(dwells, stop_or_zero_rate_target, deferred_logger);
736 const auto& summary_state = simulator.vanguard().summaryState();
737 updateWellStateFromPrimaryVariables(well_state, summary_state, deferred_logger);
738 Base::calculateReservoirRates(simulator.vanguard().eclState().runspec().co2Storage(), well_state.well(this->index_of_well_));
745 template<
typename TypeTag>
747 StandardWell<TypeTag>::
748 updatePrimaryVariablesNewton(
const BVectorWell& dwells,
749 const bool stop_or_zero_rate_target,
750 DeferredLogger& deferred_logger)
752 const Scalar dFLimit = this->param_.dwell_fraction_max_;
753 const Scalar dBHPLimit = this->param_.dbhp_max_rel_;
754 this->primary_variables_.updateNewton(dwells, stop_or_zero_rate_target, dFLimit, dBHPLimit, deferred_logger);
757 if constexpr (Base::has_polymermw) {
758 this->primary_variables_.updateNewtonPolyMW(dwells);
761 this->primary_variables_.checkFinite(deferred_logger);
768 template<
typename TypeTag>
770 StandardWell<TypeTag>::
771 updateWellStateFromPrimaryVariables(WellStateType& well_state,
772 const SummaryState& summary_state,
773 DeferredLogger& deferred_logger)
const
775 this->primary_variables_.copyToWellState(well_state, deferred_logger);
777 WellBhpThpCalculator(this->baseif_).
778 updateThp(getRefDensity(),
779 [
this,&well_state]() {
return this->baseif_.getALQ(well_state); },
780 well_state, summary_state, deferred_logger);
783 if constexpr (Base::has_polymermw) {
784 this->primary_variables_.copyToWellStatePolyMW(well_state);
792 template<
typename TypeTag>
794 StandardWell<TypeTag>::
795 updateIPR(
const Simulator& simulator, DeferredLogger& deferred_logger)
const
800 std::fill(this->ipr_a_.begin(), this->ipr_a_.end(), 0.);
801 std::fill(this->ipr_b_.begin(), this->ipr_b_.end(), 0.);
803 for (
int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
804 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
805 getMobility(simulator, perf, mob, deferred_logger);
807 const int cell_idx = this->well_cells_[perf];
808 const auto& int_quantities = simulator.model().intensiveQuantities(cell_idx, 0);
809 const auto& fs = int_quantities.fluidState();
811 Scalar p_r = this->getPerfCellPressure(fs).value();
814 std::vector<Scalar> b_perf(this->num_conservation_quantities_);
815 for (std::size_t phase = 0; phase < FluidSystem::numPhases; ++phase) {
816 if (!FluidSystem::phaseIsActive(phase)) {
819 const unsigned comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phase));
820 b_perf[comp_idx] = fs.invB(phase).value();
822 if constexpr (has_solvent) {
823 b_perf[Indices::contiSolventEqIdx] = int_quantities.solventInverseFormationVolumeFactor().value();
827 const Scalar h_perf = this->connections_.pressure_diff(perf);
828 const Scalar pressure_diff = p_r - h_perf;
833 if ( (this->isProducer() && pressure_diff < 0.) || (this->isInjector() && pressure_diff > 0.) ) {
834 deferred_logger.debug(
"CROSSFLOW_IPR",
835 "cross flow found when updateIPR for well " + name()
836 +
" . The connection is ignored in IPR calculations");
842 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quantities, cell_idx);
843 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
844 const std::vector<Scalar> tw_perf = this->wellIndex(perf,
848 std::vector<Scalar> ipr_a_perf(this->ipr_a_.size());
849 std::vector<Scalar> ipr_b_perf(this->ipr_b_.size());
850 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
851 const Scalar tw_mob = tw_perf[comp_idx] * mob[comp_idx] * b_perf[comp_idx];
852 ipr_a_perf[comp_idx] += tw_mob * pressure_diff;
853 ipr_b_perf[comp_idx] += tw_mob;
857 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
858 const unsigned oil_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
859 const unsigned gas_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
860 const Scalar rs = (fs.Rs()).value();
861 const Scalar rv = (fs.Rv()).value();
863 const Scalar dis_gas_a = rs * ipr_a_perf[oil_comp_idx];
864 const Scalar vap_oil_a = rv * ipr_a_perf[gas_comp_idx];
866 ipr_a_perf[gas_comp_idx] += dis_gas_a;
867 ipr_a_perf[oil_comp_idx] += vap_oil_a;
869 const Scalar dis_gas_b = rs * ipr_b_perf[oil_comp_idx];
870 const Scalar vap_oil_b = rv * ipr_b_perf[gas_comp_idx];
872 ipr_b_perf[gas_comp_idx] += dis_gas_b;
873 ipr_b_perf[oil_comp_idx] += vap_oil_b;
876 for (std::size_t comp_idx = 0; comp_idx < ipr_a_perf.size(); ++comp_idx) {
877 this->ipr_a_[comp_idx] += ipr_a_perf[comp_idx];
878 this->ipr_b_[comp_idx] += ipr_b_perf[comp_idx];
881 this->parallel_well_info_.communication().sum(this->ipr_a_.data(), this->ipr_a_.size());
882 this->parallel_well_info_.communication().sum(this->ipr_b_.data(), this->ipr_b_.size());
885 template<
typename TypeTag>
887 StandardWell<TypeTag>::
888 updateIPRImplicit(
const Simulator& simulator,
889 WellStateType& well_state,
890 DeferredLogger& deferred_logger)
899 auto rates = well_state.well(this->index_of_well_).surface_rates;
901 for (std::size_t p = 0; p < rates.size(); ++p) {
902 zero_rates &= rates[p] == 0.0;
904 auto& ws = well_state.well(this->index_of_well_);
906 const auto msg = fmt::format(
"updateIPRImplicit: Well {} has zero rate, IPRs might be problematic", this->name());
907 deferred_logger.debug(msg);
919 const auto& group_state = simulator.problem().wellModel().groupState();
921 std::fill(ws.implicit_ipr_a.begin(), ws.implicit_ipr_a.end(), 0.);
922 std::fill(ws.implicit_ipr_b.begin(), ws.implicit_ipr_b.end(), 0.);
924 auto inj_controls = Well::InjectionControls(0);
925 auto prod_controls = Well::ProductionControls(0);
926 prod_controls.addControl(Well::ProducerCMode::BHP);
927 prod_controls.bhp_limit = well_state.well(this->index_of_well_).bhp;
930 const auto cmode = ws.production_cmode;
931 ws.production_cmode = Well::ProducerCMode::BHP;
932 const double dt = simulator.timeStepSize();
933 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
935 const size_t nEq = this->primary_variables_.numWellEq();
939 for (
size_t i=0; i < nEq; ++i){
944 BVectorWell x_well(1);
945 x_well[0].resize(nEq);
946 this->linSys_.solve(rhs, x_well);
948 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx){
949 EvalWell comp_rate = this->primary_variables_.getQs(comp_idx);
950 const int idx = FluidSystem::activeCompToActivePhaseIdx(comp_idx);
951 for (
size_t pvIdx = 0; pvIdx < nEq; ++pvIdx) {
953 ws.implicit_ipr_b[idx] -= x_well[0][pvIdx]*comp_rate.derivative(pvIdx+Indices::numEq);
955 ws.implicit_ipr_a[idx] = ws.implicit_ipr_b[idx]*ws.bhp - comp_rate.value();
958 ws.production_cmode = cmode;
961 template<
typename TypeTag>
963 StandardWell<TypeTag>::
964 checkOperabilityUnderBHPLimit(
const WellStateType& well_state,
965 const Simulator& simulator,
966 DeferredLogger& deferred_logger)
968 const auto& summaryState = simulator.vanguard().summaryState();
969 const Scalar bhp_limit = WellBhpThpCalculator(*this).mostStrictBhpFromBhpLimits(summaryState);
972 const bool bhp_limit_not_defaulted = bhp_limit > 1.5 * unit::barsa;
973 if ( bhp_limit_not_defaulted || !this->wellHasTHPConstraints(summaryState) ) {
976 Scalar total_ipr_mass_rate = 0.0;
977 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
979 if (!FluidSystem::phaseIsActive(phaseIdx)) {
983 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
984 const Scalar ipr_rate = this->ipr_a_[compIdx] - this->ipr_b_[compIdx] * bhp_limit;
986 const Scalar rho = FluidSystem::referenceDensity( phaseIdx, Base::pvtRegionIdx() );
987 total_ipr_mass_rate += ipr_rate * rho;
989 if ( (this->isProducer() && total_ipr_mass_rate < 0.) || (this->isInjector() && total_ipr_mass_rate > 0.) ) {
990 this->operability_status_.operable_under_only_bhp_limit =
false;
994 if (this->operability_status_.operable_under_only_bhp_limit && this->wellHasTHPConstraints(summaryState)) {
998 std::vector<Scalar> well_rates_bhp_limit;
999 computeWellRatesWithBhp(simulator, bhp_limit, well_rates_bhp_limit, deferred_logger);
1001 this->adaptRatesForVFP(well_rates_bhp_limit);
1002 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1003 const Scalar thp = WellBhpThpCalculator(*this).calculateThpFromBhp(well_rates_bhp_limit,
1005 this->getRefDensity(),
1006 this->getALQ(well_state),
1009 if ( (this->isProducer() && thp < thp_limit) || (this->isInjector() && thp > thp_limit) ) {
1010 this->operability_status_.obey_thp_limit_under_bhp_limit =
false;
1021 this->operability_status_.operable_under_only_bhp_limit =
true;
1022 this->operability_status_.obey_thp_limit_under_bhp_limit =
false;
1030 template<
typename TypeTag>
1032 StandardWell<TypeTag>::
1033 checkOperabilityUnderTHPLimit(
const Simulator& simulator,
1034 const WellStateType& well_state,
1035 DeferredLogger& deferred_logger)
1037 const auto& summaryState = simulator.vanguard().summaryState();
1038 const auto obtain_bhp = this->isProducer() ? computeBhpAtThpLimitProd(well_state, simulator, summaryState, deferred_logger)
1039 : computeBhpAtThpLimitInj(simulator, summaryState, deferred_logger);
1042 this->operability_status_.can_obtain_bhp_with_thp_limit =
true;
1044 const Scalar bhp_limit = WellBhpThpCalculator(*this).mostStrictBhpFromBhpLimits(summaryState);
1045 this->operability_status_.obey_bhp_limit_with_thp_limit = this->isProducer() ?
1046 *obtain_bhp >= bhp_limit : *obtain_bhp <= bhp_limit ;
1048 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1049 if (this->isProducer() && *obtain_bhp < thp_limit) {
1050 const std::string msg =
" obtained bhp " + std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1051 +
" bars is SMALLER than thp limit "
1052 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1053 +
" bars as a producer for well " + name();
1054 deferred_logger.debug(msg);
1056 else if (this->isInjector() && *obtain_bhp > thp_limit) {
1057 const std::string msg =
" obtained bhp " + std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1058 +
" bars is LARGER than thp limit "
1059 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1060 +
" bars as a injector for well " + name();
1061 deferred_logger.debug(msg);
1064 this->operability_status_.can_obtain_bhp_with_thp_limit =
false;
1065 this->operability_status_.obey_bhp_limit_with_thp_limit =
false;
1066 if (!this->wellIsStopped()) {
1067 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1068 deferred_logger.debug(
" could not find bhp value at thp limit "
1069 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1070 +
" bar for well " + name() +
", the well might need to be closed ");
1079 template<
typename TypeTag>
1081 StandardWell<TypeTag>::
1082 allDrawDownWrongDirection(
const Simulator& simulator)
const
1084 bool all_drawdown_wrong_direction =
true;
1086 for (
int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
1087 const int cell_idx = this->well_cells_[perf];
1088 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
1089 const auto& fs = intQuants.fluidState();
1091 const Scalar pressure = this->getPerfCellPressure(fs).value();
1092 const Scalar bhp = this->primary_variables_.eval(Bhp).value();
1095 const Scalar well_pressure = bhp + this->connections_.pressure_diff(perf);
1096 const Scalar drawdown = pressure - well_pressure;
1101 if ( (drawdown < 0. && this->isInjector()) ||
1102 (drawdown > 0. && this->isProducer()) ) {
1103 all_drawdown_wrong_direction =
false;
1108 const auto& comm = this->parallel_well_info_.communication();
1109 if (comm.size() > 1)
1111 all_drawdown_wrong_direction =
1112 (comm.min(all_drawdown_wrong_direction ? 1 : 0) == 1);
1115 return all_drawdown_wrong_direction;
1121 template<
typename TypeTag>
1123 StandardWell<TypeTag>::
1124 openCrossFlowAvoidSingularity(
const Simulator& simulator)
const
1126 return !this->getAllowCrossFlow() && allDrawDownWrongDirection(simulator);
1132 template<
typename TypeTag>
1133 typename StandardWell<TypeTag>::WellConnectionProps
1134 StandardWell<TypeTag>::
1135 computePropertiesForWellConnectionPressures(
const Simulator& simulator,
1136 const WellStateType& well_state)
const
1138 auto prop_func =
typename StdWellEval::StdWellConnections::PressurePropertyFunctions {
1140 [&model = simulator.model()](
int cell_idx,
int phase_idx)
1142 return model.intensiveQuantities(cell_idx, 0)
1143 .fluidState().temperature(phase_idx).value();
1147 [&model = simulator.model()](
int cell_idx)
1149 return model.intensiveQuantities(cell_idx, 0)
1150 .fluidState().saltConcentration().value();
1154 [&model = simulator.model()](
int cell_idx)
1156 return model.intensiveQuantities(cell_idx, 0)
1157 .fluidState().pvtRegionIndex();
1161 if constexpr (Indices::enableSolvent) {
1162 prop_func.solventInverseFormationVolumeFactor =
1163 [&model = simulator.model()](
int cell_idx)
1165 return model.intensiveQuantities(cell_idx, 0)
1166 .solventInverseFormationVolumeFactor().value();
1169 prop_func.solventRefDensity = [&model = simulator.model()](
int cell_idx)
1171 return model.intensiveQuantities(cell_idx, 0)
1172 .solventRefDensity();
1176 return this->connections_.computePropertiesForPressures(well_state, prop_func);
1183 template<
typename TypeTag>
1185 StandardWell<TypeTag>::
1186 getWellConvergence(
const Simulator& simulator,
1187 const WellStateType& well_state,
1188 const std::vector<Scalar>& B_avg,
1190 const bool relax_tolerance)
const
1194 assert((
int(B_avg.size()) == this->num_conservation_quantities_) || has_polymer || has_energy || has_foam || has_brine || has_zFraction || has_bioeffects);
1196 Scalar tol_wells = this->param_.tolerance_wells_;
1198 constexpr Scalar stopped_factor = 1.e-4;
1200 constexpr Scalar dynamic_thp_factor = 1.e-1;
1201 if (this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger)) {
1202 tol_wells = tol_wells*stopped_factor;
1203 }
else if (this->getDynamicThpLimit()) {
1204 tol_wells = tol_wells*dynamic_thp_factor;
1207 std::vector<Scalar> res;
1210 this->param_.max_residual_allowed_,
1212 this->param_.relaxed_tolerance_flow_well_,
1214 this->wellIsStopped(),
1218 checkConvergenceExtraEqs(res, report);
1227 template<
typename TypeTag>
1232 WellStateType& well_state,
1235 auto fluidState = [&simulator,
this](
const int perf)
1237 const auto cell_idx = this->well_cells_[perf];
1238 return simulator.
model()
1239 .intensiveQuantities(cell_idx, 0).fluidState();
1242 const int np = this->number_of_phases_;
1243 auto setToZero = [np](Scalar* x) ->
void
1245 std::fill_n(x, np, 0.0);
1248 auto addVector = [np](
const Scalar* src, Scalar* dest) ->
void
1250 std::transform(src, src + np, dest, dest, std::plus<>{});
1253 auto& ws = well_state.well(this->index_of_well_);
1254 auto& perf_data = ws.perf_data;
1255 auto* wellPI = ws.productivity_index.data();
1256 auto* connPI = perf_data.prod_index.data();
1260 const auto preferred_phase = this->well_ecl_.getPreferredPhase();
1261 auto subsetPerfID = 0;
1263 for (
const auto& perf : *this->perf_data_) {
1264 auto allPerfID = perf.ecl_index;
1266 auto connPICalc = [&wellPICalc, allPerfID](
const Scalar mobility) -> Scalar
1271 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
1272 getMobility(simulator,
static_cast<int>(subsetPerfID), mob, deferred_logger);
1274 const auto& fs = fluidState(subsetPerfID);
1277 if (this->isInjector()) {
1278 this->computeConnLevelInjInd(fs, preferred_phase, connPICalc,
1279 mob, connPI, deferred_logger);
1282 this->computeConnLevelProdInd(fs, connPICalc, mob, connPI);
1285 addVector(connPI, wellPI);
1292 const auto& comm = this->parallel_well_info_.communication();
1293 if (comm.size() > 1) {
1294 comm.sum(wellPI, np);
1297 assert ((
static_cast<int>(subsetPerfID) == this->number_of_local_perforations_) &&
1298 "Internal logic error in processing connections for PI/II");
1303 template<
typename TypeTag>
1304 void StandardWell<TypeTag>::
1305 computeWellConnectionDensitesPressures(
const Simulator& simulator,
1306 const WellStateType& well_state,
1307 const WellConnectionProps& props,
1308 DeferredLogger& deferred_logger)
1313 const auto prop_func =
typename StdWellEval::StdWellConnections::DensityPropertyFunctions {
1318 [&model = simulator.model()](
const int cell,
1319 const std::vector<int>& phases,
1320 std::vector<Scalar>& mob)
1322 const auto& iq = model.intensiveQuantities(cell, 0);
1324 std::transform(phases.begin(), phases.end(), mob.begin(),
1325 [&iq](
const int phase) { return iq.mobility(phase).value(); });
1330 [&model = simulator.model()](
const int cell,
1331 const std::vector<int>& phases,
1332 std::vector<Scalar>& rho)
1334 const auto& fs = model.intensiveQuantities(cell, 0).fluidState();
1336 std::transform(phases.begin(), phases.end(), rho.begin(),
1337 [&fs](
const int phase) { return fs.density(phase).value(); });
1341 const auto stopped_or_zero_rate_target = this->
1342 stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
1345 .computeProperties(stopped_or_zero_rate_target, well_state,
1346 prop_func, props, deferred_logger);
1348 cachedRefDensity = this->connections_.rho(0);
1349 if (this->parallel_well_info_.communication().size() > 1) {
1350 cachedRefDensity = this->parallel_well_info_.broadcastFirstPerforationValue(cachedRefDensity);
1358 template<
typename TypeTag>
1360 StandardWell<TypeTag>::
1361 computeWellConnectionPressures(
const Simulator& simulator,
1362 const WellStateType& well_state,
1363 DeferredLogger& deferred_logger)
1365 const auto props = computePropertiesForWellConnectionPressures
1366 (simulator, well_state);
1368 computeWellConnectionDensitesPressures(simulator, well_state,
1369 props, deferred_logger);
1376 template<
typename TypeTag>
1378 StandardWell<TypeTag>::
1379 solveEqAndUpdateWellState(
const Simulator& simulator,
1380 WellStateType& well_state,
1381 DeferredLogger& deferred_logger)
1383 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1387 BVectorWell dx_well(1);
1388 dx_well[0].resize(this->primary_variables_.numWellEq());
1389 this->linSys_.solve( dx_well);
1391 updateWellState(simulator, dx_well, well_state, deferred_logger);
1398 template<
typename TypeTag>
1400 StandardWell<TypeTag>::
1401 calculateExplicitQuantities(
const Simulator& simulator,
1402 const WellStateType& well_state,
1403 DeferredLogger& deferred_logger)
1405 updatePrimaryVariables(simulator, well_state, deferred_logger);
1406 computeWellConnectionPressures(simulator, well_state, deferred_logger);
1407 this->computeAccumWell();
1412 template<
typename TypeTag>
1414 StandardWell<TypeTag>::
1415 apply(
const BVector& x, BVector& Ax)
const
1417 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1419 if (this->param_.matrix_add_well_contributions_)
1431 template<
typename TypeTag>
1433 StandardWell<TypeTag>::
1434 apply(BVector& r)
const
1436 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1444 template<
typename TypeTag>
1446 StandardWell<TypeTag>::
1447 recoverWellSolutionAndUpdateWellState(
const Simulator& simulator,
1449 WellStateType& well_state,
1452 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1457 this->
linSys_.recoverSolutionWell(x, xw);
1458 updateWellState(simulator, xw, well_state, deferred_logger);
1464 template<
typename TypeTag>
1469 std::vector<Scalar>& well_flux,
1473 const int np = this->number_of_phases_;
1474 well_flux.resize(np, 0.0);
1476 const bool allow_cf = this->getAllowCrossFlow();
1478 for (
int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
1479 const int cell_idx = this->well_cells_[perf];
1480 const auto& intQuants = simulator.
model().intensiveQuantities(cell_idx, 0);
1482 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.);
1483 getMobility(simulator, perf, mob, deferred_logger);
1484 Scalar trans_mult = simulator.
problem().template wellTransMultiplier<Scalar>(intQuants, cell_idx);
1485 const auto& wellstate_nupcol = simulator.
problem().wellModel().nupcolWellState().well(this->index_of_well_);
1486 const std::vector<Scalar> Tw = this->wellIndex(perf, intQuants, trans_mult, wellstate_nupcol);
1488 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.);
1490 computePerfRate(intQuants, mob, bhp, Tw, perf, allow_cf,
1491 cq_s, perf_rates, deferred_logger);
1493 for(
int p = 0; p < np; ++p) {
1494 well_flux[FluidSystem::activeCompToActivePhaseIdx(p)] += cq_s[p];
1498 if constexpr (has_solvent) {
1499 assert(FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx));
1501 const int gas_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
1502 well_flux[gas_pos] += cq_s[Indices::contiSolventEqIdx];
1505 this->parallel_well_info_.communication().sum(well_flux.data(), well_flux.size());
1510 template<
typename TypeTag>
1512 StandardWell<TypeTag>::
1513 computeWellRatesWithBhpIterations(
const Simulator& simulator,
1515 std::vector<Scalar>& well_flux,
1516 DeferredLogger& deferred_logger)
const
1520 StandardWell<TypeTag> well_copy(*
this);
1521 well_copy.resetDampening();
1526 WellStateType well_state_copy = simulator.problem().wellModel().wellState();
1527 const auto& group_state = simulator.problem().wellModel().groupState();
1530 const auto& summary_state = simulator.vanguard().summaryState();
1531 auto inj_controls = well_copy.well_ecl_.isInjector()
1532 ? well_copy.well_ecl_.injectionControls(summary_state)
1533 : Well::InjectionControls(0);
1534 auto prod_controls = well_copy.well_ecl_.isProducer()
1535 ? well_copy.well_ecl_.productionControls(summary_state) :
1536 Well::ProductionControls(0);
1539 auto& ws = well_state_copy.well(this->index_of_well_);
1540 if (well_copy.well_ecl_.isInjector()) {
1541 inj_controls.bhp_limit = bhp;
1542 ws.injection_cmode = Well::InjectorCMode::BHP;
1544 prod_controls.bhp_limit = bhp;
1545 ws.production_cmode = Well::ProducerCMode::BHP;
1550 const int np = this->number_of_phases_;
1551 const Scalar sign = this->well_ecl_.isInjector() ? 1.0 : -1.0;
1552 for (
int phase = 0; phase < np; ++phase){
1553 well_state_copy.wellRates(this->index_of_well_)[phase]
1554 = sign * ws.well_potentials[phase];
1556 well_copy.updatePrimaryVariables(simulator, well_state_copy, deferred_logger);
1557 well_copy.computeAccumWell();
1559 const double dt = simulator.timeStepSize();
1560 const bool converged = well_copy.iterateWellEqWithControl(simulator, dt, inj_controls, prod_controls, well_state_copy, group_state, deferred_logger);
1562 const std::string msg =
" well " + name() +
" did not get converged during well potential calculations "
1563 " potentials are computed based on unconverged solution";
1564 deferred_logger.debug(msg);
1566 well_copy.updatePrimaryVariables(simulator, well_state_copy, deferred_logger);
1567 well_copy.computeWellConnectionPressures(simulator, well_state_copy, deferred_logger);
1568 well_copy.computeWellRatesWithBhp(simulator, bhp, well_flux, deferred_logger);
1574 template<
typename TypeTag>
1575 std::vector<typename StandardWell<TypeTag>::Scalar>
1576 StandardWell<TypeTag>::
1577 computeWellPotentialWithTHP(
const Simulator& simulator,
1578 DeferredLogger& deferred_logger,
1579 const WellStateType& well_state)
const
1581 std::vector<Scalar> potentials(this->number_of_phases_, 0.0);
1582 const auto& summary_state = simulator.vanguard().summaryState();
1584 const auto& well = this->well_ecl_;
1585 if (well.isInjector()){
1586 const auto& controls = this->well_ecl_.injectionControls(summary_state);
1587 auto bhp_at_thp_limit = computeBhpAtThpLimitInj(simulator, summary_state, deferred_logger);
1588 if (bhp_at_thp_limit) {
1589 const Scalar bhp = std::min(*bhp_at_thp_limit,
1590 static_cast<Scalar
>(controls.bhp_limit));
1591 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1593 deferred_logger.warning(
"FAILURE_GETTING_CONVERGED_POTENTIAL",
1594 "Failed in getting converged thp based potential calculation for well "
1595 + name() +
". Instead the bhp based value is used");
1596 const Scalar bhp = controls.bhp_limit;
1597 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1600 computeWellRatesWithThpAlqProd(
1601 simulator, summary_state,
1602 deferred_logger, potentials, this->getALQ(well_state)
1609 template<
typename TypeTag>
1611 StandardWell<TypeTag>::
1612 computeWellPotentialsImplicit(
const Simulator& simulator,
1613 const WellStateType& well_state,
1614 std::vector<Scalar>& well_potentials,
1615 DeferredLogger& deferred_logger)
const
1620 StandardWell<TypeTag> well_copy(*
this);
1623 WellStateType well_state_copy = well_state;
1624 const auto& group_state = simulator.problem().wellModel().groupState();
1625 auto& ws = well_state_copy.well(this->index_of_well_);
1628 const auto& summary_state = simulator.vanguard().summaryState();
1629 auto inj_controls = well_copy.well_ecl_.isInjector()
1630 ? well_copy.well_ecl_.injectionControls(summary_state)
1631 : Well::InjectionControls(0);
1632 auto prod_controls = well_copy.well_ecl_.isProducer()
1633 ? well_copy.well_ecl_.productionControls(summary_state) :
1634 Well::ProductionControls(0);
1637 well_copy.onlyKeepBHPandTHPcontrols(summary_state, well_state_copy, inj_controls, prod_controls);
1640 const int num_perf = ws.perf_data.size();
1641 for (
int perf = 0; perf < num_perf; ++perf) {
1642 ws.perf_data.pressure[perf] = ws.bhp + well_copy.connections_.pressure_diff(perf);
1645 const int np = this->number_of_phases_;
1646 bool trivial =
true;
1647 for (
int phase = 0; phase < np; ++phase){
1648 trivial = trivial && (ws.well_potentials[phase] == 0.0) ;
1651 const Scalar sign = well_copy.well_ecl_.isInjector() ? 1.0 : -1.0;
1652 for (
int phase = 0; phase < np; ++phase) {
1653 ws.surface_rates[phase] = sign * ws.well_potentials[phase];
1657 well_copy.calculateExplicitQuantities(simulator, well_state_copy, deferred_logger);
1658 const double dt = simulator.timeStepSize();
1660 bool converged =
false;
1661 if (this->well_ecl_.isProducer()) {
1662 converged = well_copy.solveWellWithOperabilityCheck(simulator, dt, inj_controls, prod_controls, well_state_copy, group_state, deferred_logger);
1664 converged = well_copy.iterateWellEqWithSwitching(simulator, dt, inj_controls, prod_controls, well_state_copy, group_state, deferred_logger);
1668 well_potentials.clear();
1669 well_potentials.resize(np, 0.0);
1670 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1671 if (has_solvent && comp_idx == Indices::contiSolventEqIdx)
continue;
1672 const EvalWell rate = well_copy.primary_variables_.getQs(comp_idx);
1673 well_potentials[FluidSystem::activeCompToActivePhaseIdx(comp_idx)] = rate.value();
1677 if constexpr (has_solvent) {
1678 assert(FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx));
1680 const int gas_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
1681 const EvalWell rate = well_copy.primary_variables_.getQs(Indices::contiSolventEqIdx);
1682 well_potentials[gas_pos] += rate.value();
1688 template<
typename TypeTag>
1689 typename StandardWell<TypeTag>::Scalar
1690 StandardWell<TypeTag>::
1691 computeWellRatesAndBhpWithThpAlqProd(
const Simulator &simulator,
1692 const SummaryState &summary_state,
1693 DeferredLogger& deferred_logger,
1694 std::vector<Scalar>& potentials,
1698 auto bhp_at_thp_limit = computeBhpAtThpLimitProdWithAlq(
1699 simulator, summary_state, alq, deferred_logger,
true);
1700 if (bhp_at_thp_limit) {
1701 const auto& controls = this->well_ecl_.productionControls(summary_state);
1702 bhp = std::max(*bhp_at_thp_limit,
1703 static_cast<Scalar
>(controls.bhp_limit));
1704 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1707 deferred_logger.warning(
"FAILURE_GETTING_CONVERGED_POTENTIAL",
1708 "Failed in getting converged thp based potential calculation for well "
1709 + name() +
". Instead the bhp based value is used");
1710 const auto& controls = this->well_ecl_.productionControls(summary_state);
1711 bhp = controls.bhp_limit;
1712 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1717 template<
typename TypeTag>
1719 StandardWell<TypeTag>::
1720 computeWellRatesWithThpAlqProd(
const Simulator& simulator,
1721 const SummaryState& summary_state,
1722 DeferredLogger& deferred_logger,
1723 std::vector<Scalar>& potentials,
1727 computeWellRatesAndBhpWithThpAlqProd(simulator,
1734 template<
typename TypeTag>
1736 StandardWell<TypeTag>::
1737 computeWellPotentials(
const Simulator& simulator,
1738 const WellStateType& well_state,
1739 std::vector<Scalar>& well_potentials,
1742 const auto [compute_potential, bhp_controlled_well] =
1743 this->WellInterfaceGeneric<Scalar, IndexTraits>::computeWellPotentials(well_potentials, well_state);
1745 if (!compute_potential) {
1749 bool converged_implicit =
false;
1753 if (this->param_.local_well_solver_control_switching_ && !(this->changed_to_open_this_step_ && this->wellUnderZeroRateTarget(simulator, well_state, deferred_logger))) {
1754 converged_implicit = computeWellPotentialsImplicit(simulator, well_state, well_potentials, deferred_logger);
1756 if (!converged_implicit) {
1758 const auto& summaryState = simulator.vanguard().summaryState();
1769 const auto& ws = well_state.well(this->index_of_well_);
1771 bhp = std::max(ws.bhp, bhp);
1773 bhp = std::min(ws.bhp, bhp);
1775 assert(std::abs(bhp) != std::numeric_limits<Scalar>::max());
1776 computeWellRatesWithBhpIterations(simulator, bhp, well_potentials, deferred_logger);
1779 well_potentials = computeWellPotentialWithTHP(simulator, deferred_logger, well_state);
1783 this->checkNegativeWellPotentials(well_potentials,
1784 this->param_.check_well_operability_,
1794 template<
typename TypeTag>
1795 typename StandardWell<TypeTag>::Scalar
1798 const int openConnIdx)
const
1800 return (openConnIdx < 0)
1802 : this->connections_.rho(openConnIdx);
1809 template<
typename TypeTag>
1811 StandardWell<TypeTag>::
1812 updatePrimaryVariables(
const Simulator& simulator,
1813 const WellStateType& well_state,
1814 DeferredLogger& deferred_logger)
1816 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1818 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
1819 this->primary_variables_.update(well_state, stop_or_zero_rate_target, deferred_logger);
1822 if constexpr (Base::has_polymermw) {
1823 this->primary_variables_.updatePolyMW(well_state);
1826 this->primary_variables_.checkFinite(deferred_logger);
1832 template<
typename TypeTag>
1833 typename StandardWell<TypeTag>::Scalar
1834 StandardWell<TypeTag>::
1835 getRefDensity()
const
1837 return cachedRefDensity;
1843 template<
typename TypeTag>
1845 StandardWell<TypeTag>::
1846 updateWaterMobilityWithPolymer(
const Simulator& simulator,
1848 std::vector<EvalWell>& mob,
1849 DeferredLogger& deferred_logger)
const
1851 const int cell_idx = this->well_cells_[perf];
1852 const auto& int_quant = simulator.model().intensiveQuantities(cell_idx, 0);
1853 const EvalWell polymer_concentration = this->extendEval(int_quant.polymerConcentration());
1857 if (this->isInjector()) {
1859 const auto& visc_mult_table = PolymerModule::plyviscViscosityMultiplierTable(int_quant.pvtRegionIndex());
1860 const unsigned waterCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
1861 mob[waterCompIdx] /= (this->extendEval(int_quant.waterViscosityCorrection()) * visc_mult_table.eval(polymer_concentration,
true) );
1864 if (PolymerModule::hasPlyshlog()) {
1867 if (this->isInjector() && this->wpolymer() == 0.) {
1872 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
1873 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
1875 std::vector<EvalWell> cq_s(this->num_conservation_quantities_, {this->primary_variables_.numWellEq() + Indices::numEq, 0.});
1876 PerforationRates<Scalar> perf_rates;
1877 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quant, cell_idx);
1878 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1879 const std::vector<Scalar> Tw = this->wellIndex(perf, int_quant, trans_mult, wellstate_nupcol);
1880 computePerfRate(int_quant, mob, bhp, Tw, perf, allow_cf, cq_s,
1881 perf_rates, deferred_logger);
1883 const Scalar area = 2 * M_PI * this->perf_rep_radius_[perf] * this->perf_length_[perf];
1884 const auto& material_law_manager = simulator.problem().materialLawManager();
1885 const auto& scaled_drainage_info =
1886 material_law_manager->oilWaterScaledEpsInfoDrainage(cell_idx);
1887 const Scalar swcr = scaled_drainage_info.Swcr;
1888 const EvalWell poro = this->extendEval(int_quant.porosity());
1889 const EvalWell sw = this->extendEval(int_quant.fluidState().saturation(FluidSystem::waterPhaseIdx));
1891 const EvalWell denom = max( (area * poro * (sw - swcr)), 1e-12);
1892 const unsigned waterCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
1893 EvalWell water_velocity = cq_s[waterCompIdx] / denom * this->extendEval(int_quant.fluidState().invB(FluidSystem::waterPhaseIdx));
1895 if (PolymerModule::hasShrate()) {
1898 water_velocity *= PolymerModule::shrate( int_quant.pvtRegionIndex() ) / this->bore_diameters_[perf];
1900 const EvalWell shear_factor = PolymerModule::computeShearFactor(polymer_concentration,
1901 int_quant.pvtRegionIndex(),
1904 mob[waterCompIdx] /= shear_factor;
1908 template<
typename TypeTag>
1910 StandardWell<TypeTag>::addWellContributions(SparseMatrixAdapter& jacobian)
const
1912 this->linSys_.extract(jacobian);
1916 template <
typename TypeTag>
1918 StandardWell<TypeTag>::addWellPressureEquations(PressureMatrix& jacobian,
1919 const BVector& weights,
1920 const int pressureVarIndex,
1921 const bool use_well_weights,
1922 const WellStateType& well_state)
const
1924 this->linSys_.extractCPRPressureMatrix(jacobian,
1935 template<
typename TypeTag>
1936 typename StandardWell<TypeTag>::EvalWell
1937 StandardWell<TypeTag>::
1938 pskinwater(
const Scalar throughput,
1939 const EvalWell& water_velocity,
1940 DeferredLogger& deferred_logger)
const
1942 if constexpr (Base::has_polymermw) {
1943 const int water_table_id = this->polymerWaterTable_();
1944 if (water_table_id <= 0) {
1945 OPM_DEFLOG_THROW(std::runtime_error,
1946 fmt::format(
"Unused SKPRWAT table id used for well {}", name()),
1949 const auto& water_table_func = PolymerModule::getSkprwatTable(water_table_id);
1950 const EvalWell throughput_eval(this->primary_variables_.numWellEq() + Indices::numEq, throughput);
1952 EvalWell pskin_water(this->primary_variables_.numWellEq() + Indices::numEq, 0.0);
1953 pskin_water = water_table_func.eval(throughput_eval, water_velocity);
1956 OPM_DEFLOG_THROW(std::runtime_error,
1957 fmt::format(
"Polymermw is not activated, while injecting "
1958 "skin pressure is requested for well {}", name()),
1967 template<
typename TypeTag>
1968 typename StandardWell<TypeTag>::EvalWell
1969 StandardWell<TypeTag>::
1970 pskin(
const Scalar throughput,
1971 const EvalWell& water_velocity,
1972 const EvalWell& poly_inj_conc,
1973 DeferredLogger& deferred_logger)
const
1975 if constexpr (Base::has_polymermw) {
1976 const Scalar sign = water_velocity >= 0. ? 1.0 : -1.0;
1977 const EvalWell water_velocity_abs = abs(water_velocity);
1978 if (poly_inj_conc == 0.) {
1979 return sign * pskinwater(throughput, water_velocity_abs, deferred_logger);
1981 const int polymer_table_id = this->polymerTable_();
1982 if (polymer_table_id <= 0) {
1983 OPM_DEFLOG_THROW(std::runtime_error,
1984 fmt::format(
"Unavailable SKPRPOLY table id used for well {}", name()),
1987 const auto& skprpolytable = PolymerModule::getSkprpolyTable(polymer_table_id);
1988 const Scalar reference_concentration = skprpolytable.refConcentration;
1989 const EvalWell throughput_eval(this->primary_variables_.numWellEq() + Indices::numEq, throughput);
1991 EvalWell pskin_poly(this->primary_variables_.numWellEq() + Indices::numEq, 0.0);
1992 pskin_poly = skprpolytable.table_func.eval(throughput_eval, water_velocity_abs);
1993 if (poly_inj_conc == reference_concentration) {
1994 return sign * pskin_poly;
1997 const EvalWell pskin_water = pskinwater(throughput, water_velocity_abs, deferred_logger);
1998 const EvalWell pskin = pskin_water + (pskin_poly - pskin_water) / reference_concentration * poly_inj_conc;
1999 return sign * pskin;
2001 OPM_DEFLOG_THROW(std::runtime_error,
2002 fmt::format(
"Polymermw is not activated, while injecting "
2003 "skin pressure is requested for well {}", name()),
2012 template<
typename TypeTag>
2013 typename StandardWell<TypeTag>::EvalWell
2014 StandardWell<TypeTag>::
2015 wpolymermw(
const Scalar throughput,
2016 const EvalWell& water_velocity,
2017 DeferredLogger& deferred_logger)
const
2019 if constexpr (Base::has_polymermw) {
2020 const int table_id = this->polymerInjTable_();
2021 const auto& table_func = PolymerModule::getPlymwinjTable(table_id);
2022 const EvalWell throughput_eval(this->primary_variables_.numWellEq() + Indices::numEq, throughput);
2023 EvalWell molecular_weight(this->primary_variables_.numWellEq() + Indices::numEq, 0.);
2024 if (this->wpolymer() == 0.) {
2025 return molecular_weight;
2027 molecular_weight = table_func.eval(throughput_eval, abs(water_velocity));
2028 return molecular_weight;
2030 OPM_DEFLOG_THROW(std::runtime_error,
2031 fmt::format(
"Polymermw is not activated, while injecting "
2032 "polymer molecular weight is requested for well {}", name()),
2041 template<
typename TypeTag>
2043 StandardWell<TypeTag>::
2044 updateWaterThroughput([[maybe_unused]]
const double dt,
2045 WellStateType& well_state)
const
2047 if constexpr (Base::has_polymermw) {
2048 if (!this->isInjector()) {
2052 auto& perf_water_throughput = well_state.well(this->index_of_well_)
2053 .perf_data.water_throughput;
2055 for (
int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
2056 const Scalar perf_water_vel =
2057 this->primary_variables_.value(Bhp + 1 + perf);
2061 if (perf_water_vel > Scalar{0}) {
2062 perf_water_throughput[perf] += perf_water_vel * dt;
2072 template<
typename TypeTag>
2074 StandardWell<TypeTag>::
2075 handleInjectivityRate(
const Simulator& simulator,
2077 std::vector<EvalWell>& cq_s)
const
2079 const int cell_idx = this->well_cells_[perf];
2080 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2081 const auto& fs = int_quants.fluidState();
2082 const EvalWell b_w = this->extendEval(fs.invB(FluidSystem::waterPhaseIdx));
2083 const Scalar area = M_PI * this->bore_diameters_[perf] * this->perf_length_[perf];
2084 const int wat_vel_index = Bhp + 1 + perf;
2085 const unsigned water_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
2089 cq_s[water_comp_idx] = area * this->primary_variables_.eval(wat_vel_index) * b_w;
2095 template<
typename TypeTag>
2097 StandardWell<TypeTag>::
2098 handleInjectivityEquations(
const Simulator& simulator,
2099 const WellStateType& well_state,
2101 const EvalWell& water_flux_s,
2102 DeferredLogger& deferred_logger)
2104 const int cell_idx = this->well_cells_[perf];
2105 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2106 const auto& fs = int_quants.fluidState();
2107 const EvalWell b_w = this->extendEval(fs.invB(FluidSystem::waterPhaseIdx));
2108 const EvalWell water_flux_r = water_flux_s / b_w;
2109 const Scalar area = M_PI * this->bore_diameters_[perf] * this->perf_length_[perf];
2110 const EvalWell water_velocity = water_flux_r / area;
2111 const int wat_vel_index = Bhp + 1 + perf;
2114 const EvalWell eq_wat_vel = this->primary_variables_.eval(wat_vel_index) - water_velocity;
2116 const auto& ws = well_state.well(this->index_of_well_);
2117 const auto& perf_data = ws.perf_data;
2118 const auto& perf_water_throughput = perf_data.water_throughput;
2119 const Scalar throughput = perf_water_throughput[perf];
2120 const int pskin_index = Bhp + 1 + this->number_of_local_perforations_ + perf;
2122 EvalWell poly_conc(this->primary_variables_.numWellEq() + Indices::numEq, 0.0);
2123 poly_conc.setValue(this->wpolymer());
2126 const EvalWell eq_pskin = this->primary_variables_.eval(pskin_index)
2127 - pskin(throughput, this->primary_variables_.eval(wat_vel_index), poly_conc, deferred_logger);
2129 StandardWellAssemble<FluidSystem,Indices>(*this).
2130 assembleInjectivityEq(eq_pskin,
2135 this->primary_variables_.numWellEq(),
2143 template<
typename TypeTag>
2145 StandardWell<TypeTag>::
2146 checkConvergenceExtraEqs(
const std::vector<Scalar>& res,
2147 ConvergenceReport& report)
const
2152 if constexpr (Base::has_polymermw) {
2153 WellConvergence(*this).
2154 checkConvergencePolyMW(res, Bhp, this->param_.max_residual_allowed_, report);
2162 template<
typename TypeTag>
2164 StandardWell<TypeTag>::
2165 updateConnectionRatePolyMW(
const EvalWell& cq_s_poly,
2166 const IntensiveQuantities& int_quants,
2167 const WellStateType& well_state,
2169 std::vector<RateVector>& connectionRates,
2170 DeferredLogger& deferred_logger)
const
2173 EvalWell cq_s_polymw = cq_s_poly;
2174 if (this->isInjector()) {
2175 const int wat_vel_index = Bhp + 1 + perf;
2176 const EvalWell water_velocity = this->primary_variables_.eval(wat_vel_index);
2177 if (water_velocity > 0.) {
2178 const auto& ws = well_state.well(this->index_of_well_);
2179 const auto& perf_water_throughput = ws.perf_data.water_throughput;
2180 const Scalar throughput = perf_water_throughput[perf];
2181 const EvalWell molecular_weight = wpolymermw(throughput, water_velocity, deferred_logger);
2182 cq_s_polymw *= molecular_weight;
2188 }
else if (this->isProducer()) {
2189 if (cq_s_polymw < 0.) {
2190 cq_s_polymw *= this->extendEval(int_quants.polymerMoleWeight() );
2197 connectionRates[perf][Indices::contiPolymerMWEqIdx] = Base::restrictEval(cq_s_polymw);
2204 template<
typename TypeTag>
2205 std::optional<typename StandardWell<TypeTag>::Scalar>
2206 StandardWell<TypeTag>::
2207 computeBhpAtThpLimitProd(
const WellStateType& well_state,
2208 const Simulator& simulator,
2209 const SummaryState& summary_state,
2210 DeferredLogger& deferred_logger)
const
2212 return computeBhpAtThpLimitProdWithAlq(simulator,
2214 this->getALQ(well_state),
2219 template<
typename TypeTag>
2220 std::optional<typename StandardWell<TypeTag>::Scalar>
2221 StandardWell<TypeTag>::
2222 computeBhpAtThpLimitProdWithAlq(
const Simulator& simulator,
2223 const SummaryState& summary_state,
2224 const Scalar alq_value,
2225 DeferredLogger& deferred_logger,
2226 bool iterate_if_no_solution)
const
2230 auto frates = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2236 std::vector<Scalar> rates(3);
2237 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2238 this->adaptRatesForVFP(rates);
2242 Scalar max_pressure = 0.0;
2243 for (
int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
2244 const int cell_idx = this->well_cells_[perf];
2245 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2246 const auto& fs = int_quants.fluidState();
2247 Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2248 max_pressure = std::max(max_pressure, pressure_cell);
2250 const auto& comm = this->parallel_well_info_.communication();
2251 if (comm.size() > 1) {
2252 max_pressure = comm.max(max_pressure);
2254 auto bhpAtLimit = WellBhpThpCalculator(*this).computeBhpAtThpLimitProd(frates,
2257 this->getRefDensity(),
2259 this->getTHPConstraint(summary_state),
2263 auto v = frates(*bhpAtLimit);
2264 if (std::all_of(v.cbegin(), v.cend(), [](Scalar i){ return i <= 0; }) ) {
2269 if (!iterate_if_no_solution)
2270 return std::nullopt;
2272 auto fratesIter = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2276 std::vector<Scalar> rates(3);
2277 computeWellRatesWithBhpIterations(simulator, bhp, rates, deferred_logger);
2278 this->adaptRatesForVFP(rates);
2282 bhpAtLimit = WellBhpThpCalculator(*this).computeBhpAtThpLimitProd(fratesIter,
2285 this->getRefDensity(),
2287 this->getTHPConstraint(summary_state),
2293 auto v = frates(*bhpAtLimit);
2294 if (std::all_of(v.cbegin(), v.cend(), [](Scalar i){ return i <= 0; }) ) {
2300 return std::nullopt;
2305 template<
typename TypeTag>
2306 std::optional<typename StandardWell<TypeTag>::Scalar>
2307 StandardWell<TypeTag>::
2308 computeBhpAtThpLimitInj(
const Simulator& simulator,
2309 const SummaryState& summary_state,
2310 DeferredLogger& deferred_logger)
const
2313 auto frates = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2319 std::vector<Scalar> rates(3);
2320 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2324 return WellBhpThpCalculator(*this).computeBhpAtThpLimitInj(frates,
2326 this->getRefDensity(),
2337 template<
typename TypeTag>
2339 StandardWell<TypeTag>::
2340 iterateWellEqWithControl(
const Simulator& simulator,
2342 const Well::InjectionControls& inj_controls,
2343 const Well::ProductionControls& prod_controls,
2344 WellStateType& well_state,
2345 const GroupState<Scalar>& group_state,
2346 DeferredLogger& deferred_logger)
2348 updatePrimaryVariables(simulator, well_state, deferred_logger);
2350 const int max_iter = this->param_.max_inner_iter_wells_;
2353 bool relax_convergence =
false;
2354 this->regularize_ =
false;
2356 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
2358 if (it > this->param_.strict_inner_iter_wells_) {
2359 relax_convergence =
true;
2360 this->regularize_ =
true;
2363 auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
2365 converged = report.converged();
2371 solveEqAndUpdateWellState(simulator, well_state, deferred_logger);
2378 }
while (it < max_iter);
2384 template<
typename TypeTag>
2386 StandardWell<TypeTag>::
2387 iterateWellEqWithSwitching(
const Simulator& simulator,
2389 const Well::InjectionControls& inj_controls,
2390 const Well::ProductionControls& prod_controls,
2391 WellStateType& well_state,
2392 const GroupState<Scalar>& group_state,
2393 DeferredLogger& deferred_logger,
2394 const bool fixed_control ,
2395 const bool fixed_status )
2397 updatePrimaryVariables(simulator, well_state, deferred_logger);
2399 const int max_iter = this->param_.max_inner_iter_wells_;
2401 bool converged =
false;
2402 bool relax_convergence =
false;
2403 this->regularize_ =
false;
2404 const auto& summary_state = simulator.vanguard().summaryState();
2409 constexpr int min_its_after_switch = 4;
2411 const int max_status_switch = this->param_.max_well_status_switch_inner_iter_;
2412 int its_since_last_switch = min_its_after_switch;
2413 int switch_count= 0;
2415 const auto well_status_orig = this->wellStatus_;
2416 const auto operability_orig = this->operability_status_;
2417 auto well_status_cur = well_status_orig;
2418 int status_switch_count = 0;
2420 const bool allow_open = well_state.well(this->index_of_well_).status == WellStatus::OPEN;
2422 const bool allow_switching =
2423 !this->wellUnderZeroRateTarget(simulator, well_state, deferred_logger) &&
2424 (!fixed_control || !fixed_status) && allow_open;
2426 bool changed =
false;
2427 bool final_check =
false;
2429 this->operability_status_.resetOperability();
2430 this->operability_status_.solvable =
true;
2432 its_since_last_switch++;
2433 if (allow_switching && its_since_last_switch >= min_its_after_switch && status_switch_count < max_status_switch){
2434 const Scalar wqTotal = this->primary_variables_.eval(WQTotal).value();
2435 changed = this->updateWellControlAndStatusLocalIteration(simulator, well_state, group_state,
2436 inj_controls, prod_controls, wqTotal,
2437 deferred_logger, fixed_control, fixed_status);
2439 its_since_last_switch = 0;
2441 if (well_status_cur != this->wellStatus_) {
2442 well_status_cur = this->wellStatus_;
2443 status_switch_count++;
2446 if (!changed && final_check) {
2449 final_check =
false;
2451 if (status_switch_count == max_status_switch) {
2452 this->wellStatus_ = well_status_orig;
2456 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
2458 if (it > this->param_.strict_inner_iter_wells_) {
2459 relax_convergence =
true;
2460 this->regularize_ =
true;
2463 auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
2465 converged = report.converged();
2469 if (switch_count > 0 && its_since_last_switch < min_its_after_switch) {
2471 its_since_last_switch = min_its_after_switch;
2478 solveEqAndUpdateWellState(simulator, well_state, deferred_logger);
2480 }
while (it < max_iter);
2483 if (allow_switching){
2485 const bool is_stopped = this->wellIsStopped();
2486 if (this->wellHasTHPConstraints(summary_state)){
2487 this->operability_status_.can_obtain_bhp_with_thp_limit = !is_stopped;
2488 this->operability_status_.obey_thp_limit_under_bhp_limit = !is_stopped;
2490 this->operability_status_.operable_under_only_bhp_limit = !is_stopped;
2494 this->wellStatus_ = well_status_orig;
2495 this->operability_status_ = operability_orig;
2496 const std::string message = fmt::format(
" Well {} did not converge in {} inner iterations ("
2497 "{} switches, {} status changes).", this->name(), it, switch_count, status_switch_count);
2498 deferred_logger.debug(message);
2504 template<
typename TypeTag>
2505 std::vector<typename StandardWell<TypeTag>::Scalar>
2506 StandardWell<TypeTag>::
2507 computeCurrentWellRates(
const Simulator& simulator,
2511 std::vector<Scalar> well_q_s(this->num_conservation_quantities_, 0.);
2513 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
2514 for (
int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
2515 const int cell_idx = this->well_cells_[perf];
2516 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
2517 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.);
2518 getMobility(simulator, perf, mob, deferred_logger);
2519 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.);
2520 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(intQuants, cell_idx);
2521 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
2522 const std::vector<Scalar> Tw = this->wellIndex(perf, intQuants, trans_mult, wellstate_nupcol);
2524 computePerfRate(intQuants, mob, bhp.value(), Tw, perf, allow_cf,
2525 cq_s, perf_rates, deferred_logger);
2526 for (
int comp = 0; comp < this->num_conservation_quantities_; ++comp) {
2527 well_q_s[comp] += cq_s[comp];
2530 const auto& comm = this->parallel_well_info_.communication();
2531 if (comm.size() > 1)
2533 comm.sum(well_q_s.data(), well_q_s.size());
2540 template <
typename TypeTag>
2541 std::vector<typename StandardWell<TypeTag>::Scalar>
2545 const int num_pri_vars = this->primary_variables_.numWellEq();
2546 std::vector<Scalar> retval(num_pri_vars);
2547 for (
int ii = 0; ii < num_pri_vars; ++ii) {
2548 retval[ii] = this->primary_variables_.value(ii);
2557 template <
typename TypeTag>
2559 StandardWell<TypeTag>::
2560 setPrimaryVars(
typename std::vector<Scalar>::const_iterator it)
2562 const int num_pri_vars = this->primary_variables_.numWellEq();
2563 for (
int ii = 0; ii < num_pri_vars; ++ii) {
2564 this->primary_variables_.setValue(ii, it[ii]);
2566 return num_pri_vars;
2570 template <
typename TypeTag>
2571 typename StandardWell<TypeTag>::Eval
2572 StandardWell<TypeTag>::
2573 connectionRateEnergy(
const std::vector<EvalWell>& cq_s,
2574 const IntensiveQuantities& intQuants,
2575 DeferredLogger& deferred_logger)
const
2577 auto fs = intQuants.fluidState();
2579 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
2580 if (!FluidSystem::phaseIsActive(phaseIdx)) {
2585 EvalWell cq_r_thermal(this->primary_variables_.numWellEq() + Indices::numEq, 0.);
2586 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
2587 const bool both_oil_gas = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
2588 if (!both_oil_gas || FluidSystem::waterPhaseIdx == phaseIdx) {
2589 cq_r_thermal = cq_s[activeCompIdx] / this->extendEval(fs.invB(phaseIdx));
2592 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
2593 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
2598 const EvalWell d = this->extendEval(1.0 - fs.Rv() * fs.Rs());
2600 deferred_logger.debug(
2601 fmt::format(
"Problematic d value {} obtained for well {}"
2602 " during calculateSinglePerf with rs {}"
2603 ", rv {}. Continue as if no dissolution (rs = 0) and"
2604 " vaporization (rv = 0) for this connection.",
2605 d, this->name(), fs.Rs(), fs.Rv()));
2606 cq_r_thermal = cq_s[activeCompIdx] / this->extendEval(fs.invB(phaseIdx));
2608 if (FluidSystem::gasPhaseIdx == phaseIdx) {
2609 cq_r_thermal = (cq_s[gasCompIdx] -
2610 this->extendEval(fs.Rs()) * cq_s[oilCompIdx]) /
2611 (d * this->extendEval(fs.invB(phaseIdx)) );
2612 }
else if (FluidSystem::oilPhaseIdx == phaseIdx) {
2614 cq_r_thermal = (cq_s[oilCompIdx] - this->extendEval(fs.Rv()) *
2616 (d * this->extendEval(fs.invB(phaseIdx)) );
2622 if (this->isInjector() && !this->wellIsStopped() && cq_r_thermal > 0.0){
2624 assert(this->well_ecl_.injectorType() != InjectorType::MULTI);
2625 fs.setTemperature(this->well_ecl_.inj_temperature());
2626 typedef typename std::decay<
decltype(fs)>::type::Scalar FsScalar;
2627 typename FluidSystem::template ParameterCache<FsScalar> paramCache;
2628 const unsigned pvtRegionIdx = intQuants.pvtRegionIndex();
2629 paramCache.setRegionIndex(pvtRegionIdx);
2630 paramCache.updatePhase(fs, phaseIdx);
2632 const auto& rho = FluidSystem::density(fs, paramCache, phaseIdx);
2633 fs.setDensity(phaseIdx, rho);
2634 const auto& h = FluidSystem::enthalpy(fs, paramCache, phaseIdx);
2635 fs.setEnthalpy(phaseIdx, h);
2636 cq_r_thermal *= this->extendEval(fs.enthalpy(phaseIdx)) * this->extendEval(fs.density(phaseIdx));
2637 result += getValue(cq_r_thermal);
2638 }
else if (cq_r_thermal > 0.0) {
2639 cq_r_thermal *= getValue(fs.enthalpy(phaseIdx)) * getValue(fs.density(phaseIdx));
2640 result += Base::restrictEval(cq_r_thermal);
2643 cq_r_thermal *= this->extendEval(fs.enthalpy(phaseIdx)) * this->extendEval(fs.density(phaseIdx));
2644 result += Base::restrictEval(cq_r_thermal);
2648 return result * this->well_efficiency_factor_;
Represents the convergence status of the whole simulator, to make it possible to query and store the ...
Definition ConvergenceReport.hpp:38
Definition DeferredLogger.hpp:57
Class encapsulating some information about parallel wells.
Definition ParallelWellInfo.hpp:198
Manages the initializing and running of time dependent problems.
Definition simulator.hh:84
Problem & problem()
Return the object which specifies the pysical setup of the simulation.
Definition simulator.hh:265
Model & model()
Return the physical model used in the simulation.
Definition simulator.hh:252
StandardWellEquations< Scalar, IndexTraits, Indices::numEq > linSys_
Definition StandardWellEval.hpp:100
PrimaryVariables primary_variables_
Definition StandardWellEval.hpp:95
Definition StandardWell.hpp:60
Class for computing BHP limits.
Definition WellBhpThpCalculator.hpp:41
Scalar mostStrictBhpFromBhpLimits(const SummaryState &summaryState) const
Obtain the most strict BHP from BHP limits.
Definition WellBhpThpCalculator.cpp:93
bool isInjector() const
Definition WellInterfaceGeneric.cpp:174
bool wellHasTHPConstraints(const SummaryState &summaryState) const
Definition WellInterfaceGeneric.cpp:325
Definition WellInterfaceIndices.hpp:34
Collect per-connection static information to enable calculating connection-level or well-level produc...
Definition WellProdIndexCalculator.hpp:37
Scalar connectionProdIndStandard(const std::size_t connIdx, const Scalar connMobility) const
Compute connection-level steady-state productivity index value using dynamic phase mobility.
Definition WellProdIndexCalculator.cpp:121
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
Static data associated with a well perforation.
Definition PerforationData.hpp:30
Definition PerforationData.hpp:41