22#ifndef OPM_MULTISEGMENTWELL_IMPL_HEADER_INCLUDED
23#define OPM_MULTISEGMENTWELL_IMPL_HEADER_INCLUDED
25#ifndef OPM_MULTISEGMENTWELL_HEADER_INCLUDED
30#include <opm/common/Exceptions.hpp>
31#include <opm/common/OpmLog/OpmLog.hpp>
33#include <opm/input/eclipse/Schedule/MSW/Segment.hpp>
34#include <opm/input/eclipse/Schedule/MSW/Valve.hpp>
35#include <opm/input/eclipse/Schedule/MSW/WellSegments.hpp>
36#include <opm/input/eclipse/Schedule/Well/Connection.hpp>
37#include <opm/input/eclipse/Schedule/Well/WellConnections.hpp>
39#include <opm/input/eclipse/Units/Units.hpp>
41#include <opm/material/densead/EvaluationFormat.hpp>
52#if COMPILE_GPU_BRIDGE && (HAVE_CUDA || HAVE_OPENCL)
60 template <
typename TypeTag>
67 const int pvtRegionIdx,
68 const int num_components,
70 const int index_of_well,
72 :
Base(well, pw_info, time_step, param, rate_converter, pvtRegionIdx, num_components, num_phases, index_of_well, perf_data)
75 , segment_fluid_initial_(this->numberOfSegments(), std::vector<
Scalar>(this->num_components_, 0.0))
79 OPM_THROW(std::runtime_error,
"solvent is not supported by multisegment well yet");
83 OPM_THROW(std::runtime_error,
"polymer is not supported by multisegment well yet");
87 OPM_THROW(std::runtime_error,
"energy is not supported by multisegment well yet");
91 OPM_THROW(std::runtime_error,
"foam is not supported by multisegment well yet");
95 OPM_THROW(std::runtime_error,
"brine is not supported by multisegment well yet");
99 OPM_THROW(std::runtime_error,
"water evaporation is not supported by multisegment well yet");
103 OPM_THROW(std::runtime_error,
"MICP is not supported by multisegment well yet");
107 OPM_THROW(std::runtime_error,
108 "dissolved gas/ vapporized oil in injected oil/gas not supported by multisegment well yet."
109 " \n See (WCONINJE item 10 / WCONHIST item 8)");
119 template <
typename TypeTag>
123 const std::vector<Scalar>& depth_arg,
125 const std::vector< Scalar >& B_avg,
126 const bool changed_to_open_this_step)
128 Base::init(phase_usage_arg, depth_arg, gravity_arg, B_avg, changed_to_open_this_step);
140 this->initMatrixAndVectors();
143 for (
int local_perf_index = 0; local_perf_index < this->number_of_local_perforations_; ++local_perf_index) {
145 const int cell_idx = this->well_cells_[local_perf_index];
147 this->cell_perforation_depth_diffs_[local_perf_index] = depth_arg[cell_idx] - this->perf_depth_[this->pw_info_.localPerfToActivePerf(local_perf_index)];
155 template <
typename TypeTag>
162 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
163 this->primary_variables_.update(well_state, stop_or_zero_rate_target);
171 template <
typename TypeTag>
179 Base::updateWellStateWithTarget(simulator, group_state, well_state, deferred_logger);
182 this->scaleSegmentRatesWithWellRates(this->segments_.inlets(),
183 this->segments_.perforations(),
185 this->scaleSegmentPressuresWithBhp(well_state);
192 template <
typename TypeTag>
197 const std::vector<Scalar>& B_avg,
199 const bool relax_tolerance)
const
201 return this->MSWEval::getWellConvergence(well_state,
204 this->param_.max_residual_allowed_,
205 this->param_.tolerance_wells_,
206 this->param_.relaxed_tolerance_flow_well_,
207 this->param_.tolerance_pressure_ms_wells_,
208 this->param_.relaxed_tolerance_pressure_ms_well_,
210 this->wellIsStopped());
218 template <
typename TypeTag>
223 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) {
227 if (this->param_.matrix_add_well_contributions_) {
232 this->linSys_.apply(x, Ax);
239 template <
typename TypeTag>
244 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) {
248 this->linSys_.apply(r);
253 template <
typename TypeTag>
261 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) {
267 this->linSys_.recoverSolutionWell(x, xw);
269 updateWellState(simulator, xw, well_state, deferred_logger);
271 catch (
const NumericalProblem& exp) {
275 deferred_logger.
problem(
"In MultisegmentWell::recoverWellSolutionAndUpdateWellState for well "
276 + this->name() +
": "+exp.what());
285 template <
typename TypeTag>
290 std::vector<Scalar>& well_potentials,
293 const auto [compute_potential, bhp_controlled_well] =
296 if (!compute_potential) {
300 debug_cost_counter_ = 0;
301 bool converged_implicit =
false;
302 if (this->param_.local_well_solver_control_switching_) {
303 converged_implicit = computeWellPotentialsImplicit(simulator, well_state, well_potentials, deferred_logger);
304 if (!converged_implicit) {
305 deferred_logger.
debug(
"Implicit potential calculations failed for well "
306 + this->name() +
", reverting to original aproach.");
309 if (!converged_implicit) {
311 const auto& summaryState = simulator.vanguard().summaryState();
312 if (!Base::wellHasTHPConstraints(summaryState) || bhp_controlled_well) {
313 computeWellRatesAtBhpLimit(simulator, well_potentials, deferred_logger);
315 well_potentials = computeWellPotentialWithTHP(
316 well_state, simulator, deferred_logger);
319 deferred_logger.
debug(
"Cost in iterations of finding well potential for well "
322 this->checkNegativeWellPotentials(well_potentials,
323 this->param_.check_well_operability_,
330 template<
typename TypeTag>
334 std::vector<Scalar>& well_flux,
337 if (this->well_ecl_.isInjector()) {
338 const auto controls = this->well_ecl_.injectionControls(simulator.vanguard().summaryState());
339 computeWellRatesWithBhpIterations(simulator, controls.bhp_limit, well_flux, deferred_logger);
341 const auto controls = this->well_ecl_.productionControls(simulator.vanguard().summaryState());
342 computeWellRatesWithBhpIterations(simulator, controls.bhp_limit, well_flux, deferred_logger);
346 template<
typename TypeTag>
351 std::vector<Scalar>& well_flux,
354 const int np = this->number_of_phases_;
356 well_flux.resize(np, 0.0);
357 const bool allow_cf = this->getAllowCrossFlow();
358 const int nseg = this->numberOfSegments();
359 const WellState<Scalar>& well_state = simulator.problem().wellModel().wellState();
360 const auto& ws = well_state.
well(this->indexOfWell());
361 auto segments_copy = ws.segments;
362 segments_copy.scale_pressure(bhp);
363 const auto& segment_pressure = segments_copy.pressure;
364 for (
int seg = 0; seg < nseg; ++seg) {
365 for (
const int perf : this->segments_.perforations()[seg]) {
366 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
367 if (local_perf_index < 0)
369 const int cell_idx = this->well_cells_[local_perf_index];
370 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
372 std::vector<Scalar> mob(this->num_components_, 0.);
373 getMobility(simulator, local_perf_index, mob, deferred_logger);
374 const Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(intQuants, cell_idx);
375 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
376 const std::vector<Scalar> Tw = this->wellIndex(local_perf_index, intQuants, trans_mult, wellstate_nupcol);
377 const Scalar seg_pressure = segment_pressure[seg];
378 std::vector<Scalar> cq_s(this->num_components_, 0.);
381 computePerfRate(intQuants, mob, Tw, seg, perf, seg_pressure,
382 allow_cf, cq_s, perf_press, perf_rates, deferred_logger);
384 for(
int p = 0; p < np; ++p) {
385 well_flux[this->modelCompIdxToFlowCompIdx(p)] += cq_s[p];
389 this->parallel_well_info_.communication().sum(well_flux.data(), well_flux.size());
393 template<
typename TypeTag>
398 std::vector<Scalar>& well_flux,
411 const auto& group_state = simulator.problem().wellModel().groupState();
412 auto& ws = well_state_copy.
well(this->index_of_well_);
415 const auto& summary_state = simulator.vanguard().summaryState();
416 auto inj_controls = well_copy.
well_ecl_.isInjector()
417 ? well_copy.
well_ecl_.injectionControls(summary_state)
418 : Well::InjectionControls(0);
419 auto prod_controls = well_copy.
well_ecl_.isProducer()
420 ? well_copy.
well_ecl_.productionControls(summary_state) :
421 Well::ProductionControls(0);
425 inj_controls.bhp_limit = bhp;
426 ws.injection_cmode = Well::InjectorCMode::BHP;
428 prod_controls.bhp_limit = bhp;
429 ws.production_cmode = Well::ProducerCMode::BHP;
435 const int np = this->number_of_phases_;
437 for (
int phase = 0; phase < np; ++phase){
438 trivial = trivial && (ws.well_potentials[phase] == 0.0) ;
442 for (
int phase = 0; phase < np; ++phase) {
443 ws.surface_rates[phase] = sign * ws.well_potentials[phase];
447 this->segments_.perforations(),
451 const double dt = simulator.timeStepSize();
458 well_flux.resize(np, 0.0);
459 for (
int compIdx = 0; compIdx < this->num_components_; ++compIdx) {
461 well_flux[this->modelCompIdxToFlowCompIdx(compIdx)] = rate.value();
468 template<
typename TypeTag>
469 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
475 std::vector<Scalar> potentials(this->number_of_phases_, 0.0);
476 const auto& summary_state = simulator.vanguard().summaryState();
478 const auto& well = this->well_ecl_;
479 if (well.isInjector()) {
480 auto bhp_at_thp_limit = computeBhpAtThpLimitInj(simulator, summary_state, deferred_logger);
481 if (bhp_at_thp_limit) {
482 const auto& controls = well.injectionControls(summary_state);
483 const Scalar bhp = std::min(*bhp_at_thp_limit,
484 static_cast<Scalar>(controls.bhp_limit));
485 computeWellRatesWithBhpIterations(simulator, bhp, potentials, deferred_logger);
486 deferred_logger.
debug(
"Converged thp based potential calculation for well "
489 deferred_logger.
warning(
"FAILURE_GETTING_CONVERGED_POTENTIAL",
490 "Failed in getting converged thp based potential calculation for well "
491 + this->name() +
". Instead the bhp based value is used");
492 const auto& controls = well.injectionControls(summary_state);
493 const Scalar bhp = controls.bhp_limit;
494 computeWellRatesWithBhpIterations(simulator, bhp, potentials, deferred_logger);
497 auto bhp_at_thp_limit = computeBhpAtThpLimitProd(
498 well_state, simulator, summary_state, deferred_logger);
499 if (bhp_at_thp_limit) {
500 const auto& controls = well.productionControls(summary_state);
501 const Scalar bhp = std::max(*bhp_at_thp_limit,
502 static_cast<Scalar>(controls.bhp_limit));
503 computeWellRatesWithBhpIterations(simulator, bhp, potentials, deferred_logger);
504 deferred_logger.
debug(
"Converged thp based potential calculation for well "
507 deferred_logger.
warning(
"FAILURE_GETTING_CONVERGED_POTENTIAL",
508 "Failed in getting converged thp based potential calculation for well "
509 + this->name() +
". Instead the bhp based value is used");
510 const auto& controls = well.productionControls(summary_state);
511 const Scalar bhp = controls.bhp_limit;
512 computeWellRatesWithBhpIterations(simulator, bhp, potentials, deferred_logger);
519 template<
typename TypeTag>
524 std::vector<Scalar>& well_potentials,
535 const auto& group_state = simulator.problem().wellModel().groupState();
536 auto& ws = well_state_copy.
well(this->index_of_well_);
539 const auto& summary_state = simulator.vanguard().summaryState();
540 auto inj_controls = well_copy.
well_ecl_.isInjector()
541 ? well_copy.
well_ecl_.injectionControls(summary_state)
542 : Well::InjectionControls(0);
543 auto prod_controls = well_copy.
well_ecl_.isProducer()
544 ? well_copy.
well_ecl_.productionControls(summary_state)
545 : Well::ProductionControls(0);
553 const int np = this->number_of_phases_;
555 for (
int phase = 0; phase < np; ++phase){
556 trivial = trivial && (ws.well_potentials[phase] == 0.0) ;
560 for (
int phase = 0; phase < np; ++phase) {
561 ws.surface_rates[phase] = sign * ws.well_potentials[phase];
565 this->segments_.perforations(),
569 const double dt = simulator.timeStepSize();
571 bool converged =
false;
572 if (this->well_ecl_.isProducer()) {
573 converged = well_copy.
solveWellWithOperabilityCheck(simulator, dt, inj_controls, prod_controls, well_state_copy, group_state, deferred_logger);
575 converged = well_copy.
iterateWellEqWithSwitching(simulator, dt, inj_controls, prod_controls, well_state_copy, group_state, deferred_logger);
579 well_potentials.clear();
580 well_potentials.resize(np, 0.0);
581 for (
int compIdx = 0; compIdx < this->num_components_; ++compIdx) {
583 well_potentials[this->modelCompIdxToFlowCompIdx(compIdx)] = rate.value();
589 template <
typename TypeTag>
596 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
601 const BVectorWell dx_well = this->linSys_.solve();
602 updateWellState(simulator, dx_well, well_state, deferred_logger);
604 catch(
const NumericalProblem& exp) {
608 deferred_logger.
problem(
"In MultisegmentWell::solveEqAndUpdateWellState for well "
609 + this->name() +
": "+exp.what());
618 template <
typename TypeTag>
625 for (
int local_perf_index = 0; local_perf_index < this->number_of_local_perforations_; ++local_perf_index) {
628 std::vector<Scalar> kr(this->number_of_phases_, 0.0);
629 std::vector<Scalar> density(this->number_of_phases_, 0.0);
631 const int cell_idx = this->well_cells_[local_perf_index];
632 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
633 const auto& fs = intQuants.fluidState();
639 const int water_pos = pu.
phase_pos[Water];
640 kr[water_pos] = intQuants.relativePermeability(FluidSystem::waterPhaseIdx).value();
641 sum_kr += kr[water_pos];
642 density[water_pos] = fs.density(FluidSystem::waterPhaseIdx).value();
647 kr[oil_pos] = intQuants.relativePermeability(FluidSystem::oilPhaseIdx).value();
648 sum_kr += kr[oil_pos];
649 density[oil_pos] = fs.density(FluidSystem::oilPhaseIdx).value();
654 kr[gas_pos] = intQuants.relativePermeability(FluidSystem::gasPhaseIdx).value();
655 sum_kr += kr[gas_pos];
656 density[gas_pos] = fs.density(FluidSystem::gasPhaseIdx).value();
659 assert(sum_kr != 0.);
662 Scalar average_density = 0.;
663 for (
int p = 0; p < this->number_of_phases_; ++p) {
664 average_density += kr[p] * density[p];
666 average_density /= sum_kr;
668 this->cell_perforation_pressure_diffs_[local_perf_index] = this->gravity_ * average_density * this->cell_perforation_depth_diffs_[local_perf_index];
676 template <
typename TypeTag>
681 for (
int seg = 0; seg < this->numberOfSegments(); ++seg) {
683 const Scalar surface_volume = getSegmentSurfaceVolume(simulator, seg).value();
684 for (
int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
685 segment_fluid_initial_[seg][comp_idx] = surface_volume * this->primary_variables_.surfaceVolumeFraction(seg, comp_idx).value();
694 template <
typename TypeTag>
698 const BVectorWell& dwells,
701 const Scalar relaxation_factor)
703 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
705 const Scalar dFLimit = this->param_.dwell_fraction_max_;
706 const Scalar max_pressure_change = this->param_.max_pressure_change_ms_wells_;
707 const bool stop_or_zero_rate_target =
708 this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
709 this->primary_variables_.updateNewton(dwells,
712 stop_or_zero_rate_target,
713 max_pressure_change);
715 const auto& summary_state = simulator.vanguard().summaryState();
716 this->primary_variables_.copyToWellState(*
this, getRefDensity(),
722 auto& ws = well_state.
well(this->index_of_well_);
723 this->segments_.copyPhaseDensities(ws.pu, ws.segments);
726 Base::calculateReservoirRates(simulator.vanguard().eclState().runspec().co2Storage(), well_state.
well(this->index_of_well_));
733 template <
typename TypeTag>
740 updatePrimaryVariables(simulator, well_state, deferred_logger);
741 computePerfCellPressDiffs(simulator);
742 computeInitialSegmentFluids(simulator);
749 template<
typename TypeTag>
757 auto fluidState = [&simulator,
this](
const int local_perf_index)
759 const auto cell_idx = this->well_cells_[local_perf_index];
760 return simulator.model()
761 .intensiveQuantities(cell_idx, 0).fluidState();
764 const int np = this->number_of_phases_;
765 auto setToZero = [np](
Scalar* x) ->
void
767 std::fill_n(x, np, 0.0);
770 auto addVector = [np](
const Scalar* src,
Scalar* dest) ->
void
772 std::transform(src, src + np, dest, dest, std::plus<>{});
775 auto& ws = well_state.
well(this->index_of_well_);
776 auto& perf_data = ws.perf_data;
777 auto* connPI = perf_data.prod_index.data();
778 auto* wellPI = ws.productivity_index.data();
782 const auto preferred_phase = this->well_ecl_.getPreferredPhase();
783 auto subsetPerfID = 0;
785 for (
const auto& perf : *this->perf_data_){
786 auto allPerfID = perf.ecl_index;
788 auto connPICalc = [&wellPICalc, allPerfID](
const Scalar mobility) ->
Scalar
793 std::vector<Scalar> mob(this->num_components_, 0.0);
798 getMobility(simulator,
static_cast<int>(subsetPerfID), mob, deferred_logger);
800 const auto& fs = fluidState(subsetPerfID);
803 if (this->isInjector()) {
804 this->computeConnLevelInjInd(fs, preferred_phase, connPICalc,
805 mob, connPI, deferred_logger);
808 this->computeConnLevelProdInd(fs, connPICalc, mob, connPI);
811 addVector(connPI, wellPI);
818 const auto& comm = this->parallel_well_info_.communication();
819 if (comm.size() > 1) {
820 comm.sum(wellPI, np);
823 assert (
static_cast<int>(subsetPerfID) == this->number_of_local_perforations_ &&
824 "Internal logic error in processing connections for PI/II");
831 template<
typename TypeTag>
835 [[maybe_unused]]
const int openConnIdx)
const
840 const auto segNum = this->wellEcl()
841 .getConnections()[globalConnIdx].segment();
843 const auto segIdx = this->wellEcl()
844 .getSegments().segmentNumberToIndex(segNum);
846 return this->segments_.density(segIdx).value();
853 template<
typename TypeTag>
858 if (this->number_of_local_perforations_ == 0) {
862 this->linSys_.extract(jacobian);
866 template<
typename TypeTag>
871 const int pressureVarIndex,
872 const bool use_well_weights,
875 if (this->number_of_local_perforations_ == 0) {
880 this->linSys_.extractCPRPressureMatrix(jacobian,
890 template<
typename TypeTag>
891 template<
class Value>
897 const std::vector<Value>& b_perfcells,
898 const std::vector<Value>& mob_perfcells,
899 const std::vector<Scalar>& Tw,
901 const Value& segment_pressure,
902 const Value& segment_density,
903 const bool& allow_cf,
904 const std::vector<Value>& cmix_s,
905 std::vector<Value>& cq_s,
910 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
911 if (local_perf_index < 0)
915 const Value perf_seg_press_diff = this->gravity() * segment_density *
916 this->segments_.local_perforation_depth_diff(local_perf_index);
918 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[local_perf_index];
922 perf_press = segment_pressure + perf_seg_press_diff;
925 const Value cell_press_at_perf = pressure_cell - cell_perf_press_diff;
928 const Value drawdown = cell_press_at_perf - perf_press;
931 if (drawdown > 0.0) {
933 if (!allow_cf && this->isInjector()) {
938 for (
int comp_idx = 0; comp_idx < this->numComponents(); ++comp_idx) {
939 const Value cq_p = - Tw[comp_idx] * (mob_perfcells[comp_idx] * drawdown);
940 cq_s[comp_idx] = b_perfcells[comp_idx] * cq_p;
943 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
944 const unsigned oilCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
945 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
946 const Value cq_s_oil = cq_s[oilCompIdx];
947 const Value cq_s_gas = cq_s[gasCompIdx];
948 cq_s[gasCompIdx] += rs * cq_s_oil;
949 cq_s[oilCompIdx] += rv * cq_s_gas;
953 if (!allow_cf && this->isProducer()) {
958 Value total_mob = mob_perfcells[0];
959 for (
int comp_idx = 1; comp_idx < this->numComponents(); ++comp_idx) {
960 total_mob += mob_perfcells[comp_idx];
964 Value volume_ratio = 0.0;
965 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
966 const unsigned waterCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
967 volume_ratio += cmix_s[waterCompIdx] / b_perfcells[waterCompIdx];
970 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
971 const unsigned oilCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
972 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
977 const Value d = 1.0 - rv * rs;
979 if (getValue(d) == 0.0) {
981 fmt::format(
"Zero d value obtained for well {} "
982 "during flux calculation with rs {} and rv {}",
983 this->name(), rs, rv),
987 const Value tmp_oil = (cmix_s[oilCompIdx] - rv * cmix_s[gasCompIdx]) / d;
988 volume_ratio += tmp_oil / b_perfcells[oilCompIdx];
990 const Value tmp_gas = (cmix_s[gasCompIdx] - rs * cmix_s[oilCompIdx]) / d;
991 volume_ratio += tmp_gas / b_perfcells[gasCompIdx];
993 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
994 const unsigned oilCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
995 volume_ratio += cmix_s[oilCompIdx] / b_perfcells[oilCompIdx];
997 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
998 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
999 volume_ratio += cmix_s[gasCompIdx] / b_perfcells[gasCompIdx];
1003 for (
int componentIdx = 0; componentIdx < this->numComponents(); ++componentIdx) {
1004 const Value cqt_i = - Tw[componentIdx] * (total_mob * drawdown);
1005 Value cqt_is = cqt_i / volume_ratio;
1006 cq_s[componentIdx] = cmix_s[componentIdx] * cqt_is;
1011 if (this->isProducer()) {
1012 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1013 const unsigned oilCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
1014 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
1023 const Scalar d = 1.0 - getValue(rv) * getValue(rs);
1026 perf_rates.
vap_oil = getValue(rv) * (getValue(cq_s[gasCompIdx]) - getValue(rs) * getValue(cq_s[oilCompIdx])) / d;
1029 perf_rates.
dis_gas = getValue(rs) * (getValue(cq_s[oilCompIdx]) - getValue(rv) * getValue(cq_s[gasCompIdx])) / d;
1034 template <
typename TypeTag>
1035 template<
class Value>
1039 const std::vector<Value>& mob_perfcells,
1040 const std::vector<Scalar>& Tw,
1043 const Value& segment_pressure,
1044 const bool& allow_cf,
1045 std::vector<Value>& cq_s,
1051 auto obtain = [
this](
const Eval& value)
1053 if constexpr (std::is_same_v<Value, Scalar>) {
1054 static_cast<void>(
this);
1055 return getValue(value);
1057 return this->extendEval(value);
1060 auto obtainN = [](
const auto& value)
1062 if constexpr (std::is_same_v<Value, Scalar>) {
1063 return getValue(value);
1068 const auto& fs = int_quants.fluidState();
1070 const Value pressure_cell = obtain(this->getPerfCellPressure(fs));
1071 const Value rs = obtain(fs.Rs());
1072 const Value rv = obtain(fs.Rv());
1075 std::vector<Value> b_perfcells(this->num_components_, 0.0);
1077 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
1078 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1082 const unsigned compIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
1083 b_perfcells[compIdx] = obtain(fs.invB(phaseIdx));
1086 std::vector<Value> cmix_s(this->numComponents(), 0.0);
1087 for (
int comp_idx = 0; comp_idx < this->numComponents(); ++comp_idx) {
1088 cmix_s[comp_idx] = obtainN(this->primary_variables_.surfaceVolumeFraction(seg, comp_idx));
1091 this->computePerfRate(pressure_cell,
1099 obtainN(this->segments_.density(seg)),
1108 template <
typename TypeTag>
1128 auto info = this->getFirstPerforationFluidStateInfo(simulator);
1129 temperature.setValue(std::get<0>(info));
1130 saltConcentration = this->extendEval(std::get<1>(info));
1132 this->segments_.computeFluidProperties(temperature,
1134 this->primary_variables_,
1139 template <
typename TypeTag>
1140 template<
class Value>
1144 const int local_perf_index,
1145 std::vector<Value>& mob,
1148 auto obtain = [
this](
const Eval& value)
1150 if constexpr (std::is_same_v<Value, Scalar>) {
1151 static_cast<void>(
this);
1152 return getValue(value);
1154 return this->extendEval(value);
1161 const auto perf_ecl_index = this->perforationData()[local_perf_index].ecl_index;
1162 const Connection& con = this->well_ecl_.getConnections()[perf_ecl_index];
1163 const int seg = this->segmentNumberToIndex(con.segment());
1167 const Scalar segment_pres = this->primary_variables_.getSegmentPressure(seg).value();
1168 const Scalar perf_seg_press_diff = this->gravity() * this->segments_.density(seg).value()
1169 * this->segments_.local_perforation_depth_diff(local_perf_index);
1170 const Scalar perf_press = segment_pres + perf_seg_press_diff;
1171 const Scalar multiplier = this->getInjMult(local_perf_index, segment_pres, perf_press, deferred_logger);
1172 for (std::size_t i = 0; i < mob.size(); ++i) {
1173 mob[i] *= multiplier;
1180 template<
typename TypeTag>
1185 return this->segments_.getRefDensity();
1188 template<
typename TypeTag>
1195 const auto& summaryState = simulator.vanguard().summaryState();
1199 const bool bhp_limit_not_defaulted = bhp_limit > 1.5 * unit::barsa;
1200 if ( bhp_limit_not_defaulted || !this->wellHasTHPConstraints(summaryState) ) {
1203 Scalar total_ipr_mass_rate = 0.0;
1204 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
1206 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1210 const unsigned compIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
1211 const Scalar ipr_rate = this->ipr_a_[compIdx] - this->ipr_b_[compIdx] * bhp_limit;
1213 const Scalar rho = FluidSystem::referenceDensity( phaseIdx, Base::pvtRegionIdx() );
1214 total_ipr_mass_rate += ipr_rate * rho;
1216 if ( (this->isProducer() && total_ipr_mass_rate < 0.) || (this->isInjector() && total_ipr_mass_rate > 0.) ) {
1217 this->operability_status_.operable_under_only_bhp_limit =
false;
1221 if (this->operability_status_.operable_under_only_bhp_limit && this->wellHasTHPConstraints(summaryState)) {
1225 std::vector<Scalar> well_rates_bhp_limit;
1226 computeWellRatesWithBhp(simulator, bhp_limit, well_rates_bhp_limit, deferred_logger);
1228 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1231 this->getRefDensity(),
1232 this->wellEcl().alq_value(summaryState),
1235 if ( (this->isProducer() && thp < thp_limit) || (this->isInjector() && thp > thp_limit) ) {
1236 this->operability_status_.obey_thp_limit_under_bhp_limit =
false;
1247 this->operability_status_.operable_under_only_bhp_limit =
true;
1248 this->operability_status_.obey_thp_limit_under_bhp_limit =
false;
1254 template<
typename TypeTag>
1262 std::fill(this->ipr_a_.begin(), this->ipr_a_.end(), 0.);
1263 std::fill(this->ipr_b_.begin(), this->ipr_b_.end(), 0.);
1265 const int nseg = this->numberOfSegments();
1266 std::vector<Scalar> seg_dp(nseg, 0.0);
1267 for (
int seg = 0; seg < nseg; ++seg) {
1269 const Scalar dp = this->getSegmentDp(seg,
1270 this->segments_.density(seg).value(),
1273 for (
const int perf : this->segments_.perforations()[seg]) {
1274 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
1275 if (local_perf_index < 0)
1277 std::vector<Scalar> mob(this->num_components_, 0.0);
1280 getMobility(simulator, local_perf_index, mob, deferred_logger);
1282 const int cell_idx = this->well_cells_[local_perf_index];
1283 const auto& int_quantities = simulator.model().intensiveQuantities(cell_idx, 0);
1284 const auto& fs = int_quantities.fluidState();
1286 const Scalar perf_seg_press_diff = this->segments_.getPressureDiffSegLocalPerf(seg, local_perf_index);
1288 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[local_perf_index];
1289 const Scalar pressure_cell = this->getPerfCellPressure(fs).value();
1292 std::vector<Scalar> b_perf(this->num_components_);
1293 for (std::size_t phase = 0; phase < FluidSystem::numPhases; ++phase) {
1294 if (!FluidSystem::phaseIsActive(phase)) {
1297 const unsigned comp_idx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phase));
1298 b_perf[comp_idx] = fs.invB(phase).value();
1302 const Scalar h_perf = cell_perf_press_diff + perf_seg_press_diff + dp;
1303 const Scalar pressure_diff = pressure_cell - h_perf;
1306 if ( (this->isProducer() && pressure_diff < 0.) || (this->isInjector() && pressure_diff > 0.) ) {
1307 deferred_logger.
debug(
"CROSSFLOW_IPR",
1308 "cross flow found when updateIPR for well " + this->name());
1312 const Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quantities, cell_idx);
1313 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1314 const std::vector<Scalar> tw_perf = this->wellIndex(local_perf_index, int_quantities, trans_mult, wellstate_nupcol);
1315 std::vector<Scalar> ipr_a_perf(this->ipr_a_.size());
1316 std::vector<Scalar> ipr_b_perf(this->ipr_b_.size());
1317 for (
int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
1318 const Scalar tw_mob = tw_perf[comp_idx] * mob[comp_idx] * b_perf[comp_idx];
1319 ipr_a_perf[comp_idx] += tw_mob * pressure_diff;
1320 ipr_b_perf[comp_idx] += tw_mob;
1324 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1325 const unsigned oil_comp_idx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
1326 const unsigned gas_comp_idx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
1327 const Scalar rs = (fs.Rs()).value();
1328 const Scalar rv = (fs.Rv()).value();
1330 const Scalar dis_gas_a = rs * ipr_a_perf[oil_comp_idx];
1331 const Scalar vap_oil_a = rv * ipr_a_perf[gas_comp_idx];
1333 ipr_a_perf[gas_comp_idx] += dis_gas_a;
1334 ipr_a_perf[oil_comp_idx] += vap_oil_a;
1336 const Scalar dis_gas_b = rs * ipr_b_perf[oil_comp_idx];
1337 const Scalar vap_oil_b = rv * ipr_b_perf[gas_comp_idx];
1339 ipr_b_perf[gas_comp_idx] += dis_gas_b;
1340 ipr_b_perf[oil_comp_idx] += vap_oil_b;
1343 for (std::size_t comp_idx = 0; comp_idx < ipr_a_perf.size(); ++comp_idx) {
1344 this->ipr_a_[comp_idx] += ipr_a_perf[comp_idx];
1345 this->ipr_b_[comp_idx] += ipr_b_perf[comp_idx];
1349 this->parallel_well_info_.communication().sum(this->ipr_a_.data(), this->ipr_a_.size());
1350 this->parallel_well_info_.communication().sum(this->ipr_b_.data(), this->ipr_b_.size());
1353 template<
typename TypeTag>
1367 auto rates = well_state.
well(this->index_of_well_).surface_rates;
1369 for (std::size_t p = 0; p < rates.size(); ++p) {
1370 zero_rates &= rates[p] == 0.0;
1372 auto& ws = well_state.
well(this->index_of_well_);
1374 const auto msg = fmt::format(
"updateIPRImplicit: Well {} has zero rate, IPRs might be problematic", this->name());
1375 deferred_logger.
debug(msg);
1387 const auto& group_state = simulator.problem().wellModel().groupState();
1389 std::fill(ws.implicit_ipr_a.begin(), ws.implicit_ipr_a.end(), 0.);
1390 std::fill(ws.implicit_ipr_b.begin(), ws.implicit_ipr_b.end(), 0.);
1392 auto inj_controls = Well::InjectionControls(0);
1393 auto prod_controls = Well::ProductionControls(0);
1394 prod_controls.addControl(Well::ProducerCMode::BHP);
1395 prod_controls.bhp_limit = well_state.
well(this->index_of_well_).bhp;
1398 const auto cmode = ws.production_cmode;
1399 ws.production_cmode = Well::ProducerCMode::BHP;
1400 const double dt = simulator.timeStepSize();
1401 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
1403 BVectorWell rhs(this->numberOfSegments());
1405 rhs[0][SPres] = -1.0;
1407 const BVectorWell x_well = this->linSys_.solve(rhs);
1408 constexpr int num_eq = MSWEval::numWellEq;
1409 for (
int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx){
1410 const EvalWell comp_rate = this->primary_variables_.getQs(comp_idx);
1411 const int idx = this->modelCompIdxToFlowCompIdx(comp_idx);
1412 for (
size_t pvIdx = 0; pvIdx < num_eq; ++pvIdx) {
1414 ws.implicit_ipr_b[idx] -= x_well[0][pvIdx]*comp_rate.derivative(pvIdx+Indices::numEq);
1416 ws.implicit_ipr_a[idx] = ws.implicit_ipr_b[idx]*ws.bhp - comp_rate.value();
1419 ws.production_cmode = cmode;
1422 template<
typename TypeTag>
1429 const auto& summaryState = simulator.vanguard().summaryState();
1430 const auto obtain_bhp = this->isProducer()
1431 ? computeBhpAtThpLimitProd(
1432 well_state, simulator, summaryState, deferred_logger)
1433 : computeBhpAtThpLimitInj(simulator, summaryState, deferred_logger);
1436 this->operability_status_.can_obtain_bhp_with_thp_limit =
true;
1439 this->operability_status_.obey_bhp_limit_with_thp_limit = (*obtain_bhp >= bhp_limit);
1441 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1442 if (this->isProducer() && *obtain_bhp < thp_limit) {
1443 const std::string msg =
" obtained bhp " +
std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1444 +
" bars is SMALLER than thp limit "
1446 +
" bars as a producer for well " + this->name();
1447 deferred_logger.
debug(msg);
1449 else if (this->isInjector() && *obtain_bhp > thp_limit) {
1450 const std::string msg =
" obtained bhp " +
std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1451 +
" bars is LARGER than thp limit "
1453 +
" bars as a injector for well " + this->name();
1454 deferred_logger.
debug(msg);
1459 this->operability_status_.can_obtain_bhp_with_thp_limit =
false;
1460 this->operability_status_.obey_bhp_limit_with_thp_limit =
false;
1461 if (!this->wellIsStopped()) {
1462 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1463 deferred_logger.
debug(
" could not find bhp value at thp limit "
1465 +
" bar for well " + this->name() +
", the well might need to be closed ");
1474 template<
typename TypeTag>
1479 const Well::InjectionControls& inj_controls,
1480 const Well::ProductionControls& prod_controls,
1485 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return true;
1487 const int max_iter_number = this->param_.max_inner_iter_ms_wells_;
1491 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1496 updatePrimaryVariables(simulator, well_state, deferred_logger);
1498 std::vector<std::vector<Scalar> > residual_history;
1499 std::vector<Scalar> measure_history;
1502 Scalar relaxation_factor = 1.;
1503 bool converged =
false;
1504 bool relax_convergence =
false;
1505 this->regularize_ =
false;
1506 for (; it < max_iter_number; ++it, ++debug_cost_counter_) {
1508 if (it > this->param_.strict_inner_iter_wells_) {
1509 relax_convergence =
true;
1510 this->regularize_ =
true;
1513 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls,
1514 well_state, group_state, deferred_logger);
1516 const auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
1517 if (report.converged()) {
1524 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1528 residual_history.push_back(residuals);
1529 measure_history.push_back(this->getResidualMeasureValue(well_state,
1530 residual_history[it],
1531 this->param_.tolerance_wells_,
1532 this->param_.tolerance_pressure_ms_wells_,
1535 bool min_relaxation_reached = this->update_relaxation_factor(measure_history, relaxation_factor, this->regularize_, deferred_logger);
1536 if (min_relaxation_reached || this->repeatedStagnation(measure_history, this->regularize_, deferred_logger)) {
1538 const auto reportStag = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger,
true);
1539 if (reportStag.converged()) {
1541 std::string message = fmt::format(
"Well stagnates/oscillates but {} manages to get converged with relaxed tolerances in {} inner iterations."
1543 deferred_logger.
debug(message);
1550 BVectorWell dx_well;
1552 dx_well = this->linSys_.solve();
1553 updateWellState(simulator, dx_well, well_state, deferred_logger, relaxation_factor);
1555 catch(
const NumericalProblem& exp) {
1559 deferred_logger.
problem(
"In MultisegmentWell::iterateWellEqWithControl for well "
1560 + this->name() +
": "+exp.what());
1567 std::ostringstream sstr;
1568 sstr <<
" Well " << this->name() <<
" converged in " << it <<
" inner iterations.";
1569 if (relax_convergence)
1570 sstr <<
" (A relaxed tolerance was used after "<< this->param_.strict_inner_iter_wells_ <<
" iterations)";
1574 deferred_logger.
debug(sstr.str(), OpmLog::defaultDebugVerbosityLevel + (it == 0));
1576 std::ostringstream sstr;
1577 sstr <<
" Well " << this->name() <<
" did not converge in " << it <<
" inner iterations.";
1578#define EXTRA_DEBUG_MSW 0
1580 sstr <<
"***** Outputting the residual history for well " << this->name() <<
" during inner iterations:";
1581 for (
int i = 0; i < it; ++i) {
1582 const auto& residual = residual_history[i];
1583 sstr <<
" residual at " << i <<
"th iteration ";
1584 for (
const auto& res : residual) {
1587 sstr <<
" " << measure_history[i] <<
" \n";
1590#undef EXTRA_DEBUG_MSW
1591 deferred_logger.
debug(sstr.str());
1598 template<
typename TypeTag>
1603 const Well::InjectionControls& inj_controls,
1604 const Well::ProductionControls& prod_controls,
1608 const bool fixed_control ,
1609 const bool fixed_status )
1611 const int max_iter_number = this->param_.max_inner_iter_ms_wells_;
1615 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1620 updatePrimaryVariables(simulator, well_state, deferred_logger);
1622 std::vector<std::vector<Scalar> > residual_history;
1623 std::vector<Scalar> measure_history;
1626 Scalar relaxation_factor = 1.;
1627 bool converged =
false;
1628 bool relax_convergence =
false;
1629 this->regularize_ =
false;
1630 const auto& summary_state = simulator.vanguard().summaryState();
1635 const int min_its_after_switch = 3;
1636 int its_since_last_switch = min_its_after_switch;
1637 int switch_count= 0;
1639 const auto well_status_orig = this->wellStatus_;
1640 const auto operability_orig = this->operability_status_;
1642 const bool allow_open = this->well_ecl_.getStatus() == WellStatus::OPEN &&
1643 well_state.
well(this->index_of_well_).status == WellStatus::OPEN;
1645 const bool allow_switching = !this->wellUnderZeroRateTarget(simulator, well_state, deferred_logger) &&
1646 (!fixed_control || !fixed_status) && allow_open;
1647 bool final_check =
false;
1649 this->operability_status_.resetOperability();
1650 this->operability_status_.solvable =
true;
1652 for (; it < max_iter_number; ++it, ++debug_cost_counter_) {
1653 ++its_since_last_switch;
1654 if (allow_switching && its_since_last_switch >= min_its_after_switch){
1655 const Scalar wqTotal = this->primary_variables_.getWQTotal().value();
1656 bool changed = this->updateWellControlAndStatusLocalIteration(simulator, well_state, group_state,
1657 inj_controls, prod_controls, wqTotal,
1658 deferred_logger, fixed_control,
1661 its_since_last_switch = 0;
1664 if (!changed && final_check) {
1667 final_check =
false;
1671 if (it > this->param_.strict_inner_iter_wells_) {
1672 relax_convergence =
true;
1673 this->regularize_ =
true;
1676 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls,
1677 well_state, group_state, deferred_logger);
1680 const auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
1681 converged = report.converged();
1682 if (this->parallel_well_info_.communication().size() > 1 &&
1683 this->parallel_well_info_.communication().max(converged) != this->parallel_well_info_.communication().min(converged)) {
1684 OPM_THROW(std::runtime_error, fmt::format(
"Misalignment of the parallel simulation run in iterateWellEqWithSwitching - the well calculation for well {} succeeded some ranks but failed on other ranks.", this->name()));
1689 if (switch_count > 0 && its_since_last_switch < min_its_after_switch) {
1691 its_since_last_switch = min_its_after_switch;
1699 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1703 residual_history.push_back(residuals);
1707 measure_history.push_back(this->getResidualMeasureValue(well_state,
1708 residual_history[it],
1709 this->param_.tolerance_wells_,
1710 this->param_.tolerance_pressure_ms_wells_,
1712 bool min_relaxation_reached = this->update_relaxation_factor(measure_history, relaxation_factor, this->regularize_, deferred_logger);
1713 if (min_relaxation_reached || this->repeatedStagnation(measure_history, this->regularize_, deferred_logger)) {
1719 const BVectorWell dx_well = this->linSys_.solve();
1720 updateWellState(simulator, dx_well, well_state, deferred_logger, relaxation_factor);
1722 catch(
const NumericalProblem& exp) {
1726 deferred_logger.
problem(
"In MultisegmentWell::iterateWellEqWithSwitching for well "
1727 + this->name() +
": "+exp.what());
1733 if (allow_switching){
1735 const bool is_stopped = this->wellIsStopped();
1736 if (this->wellHasTHPConstraints(summary_state)){
1737 this->operability_status_.can_obtain_bhp_with_thp_limit = !is_stopped;
1738 this->operability_status_.obey_thp_limit_under_bhp_limit = !is_stopped;
1740 this->operability_status_.operable_under_only_bhp_limit = !is_stopped;
1743 std::string message = fmt::format(
" Well {} converged in {} inner iterations ("
1744 "{} control/status switches).", this->name(), it, switch_count);
1745 if (relax_convergence) {
1746 message.append(fmt::format(
" (A relaxed tolerance was used after {} iterations)",
1747 this->param_.strict_inner_iter_wells_));
1749 deferred_logger.
debug(message, OpmLog::defaultDebugVerbosityLevel + ((it == 0) && (switch_count == 0)));
1751 this->wellStatus_ = well_status_orig;
1752 this->operability_status_ = operability_orig;
1753 const std::string message = fmt::format(
" Well {} did not converge in {} inner iterations ("
1754 "{} control/status switches).", this->name(), it, switch_count);
1755 deferred_logger.
debug(message);
1756 this->primary_variables_.outputLowLimitPressureSegments(deferred_logger);
1763 template<
typename TypeTag>
1768 const Well::InjectionControls& inj_controls,
1769 const Well::ProductionControls& prod_controls,
1774 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1777 this->segments_.updateUpwindingSegments(this->primary_variables_);
1780 computeSegmentFluidProperties(simulator, deferred_logger);
1783 this->linSys_.clear();
1785 auto& ws = well_state.
well(this->index_of_well_);
1786 ws.phase_mixing_rates.fill(0.0);
1793 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
1795 const int nseg = this->numberOfSegments();
1797 const Scalar rhow = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) ?
1798 FluidSystem::referenceDensity( FluidSystem::waterPhaseIdx, Base::pvtRegionIdx() ) : 0.0;
1799 const unsigned watCompIdx = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) ?
1800 Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx) : 0;
1802 for (
int seg = 0; seg < nseg; ++seg) {
1804 const EvalWell seg_pressure = this->primary_variables_.getSegmentPressure(seg);
1805 auto& perf_data = ws.perf_data;
1806 auto& perf_rates = perf_data.phase_rates;
1807 auto& perf_press_state = perf_data.pressure;
1808 for (
const int perf : this->segments_.perforations()[seg]) {
1809 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
1810 if (local_perf_index < 0)
1812 const int cell_idx = this->well_cells_[local_perf_index];
1813 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
1814 std::vector<EvalWell> mob(this->num_components_, 0.0);
1815 getMobility(simulator, local_perf_index, mob, deferred_logger);
1816 const Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quants, cell_idx);
1817 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1818 const std::vector<Scalar> Tw = this->wellIndex(local_perf_index, int_quants, trans_mult, wellstate_nupcol);
1819 std::vector<EvalWell> cq_s(this->num_components_, 0.0);
1822 computePerfRate(int_quants, mob, Tw, seg, perf, seg_pressure,
1823 allow_cf, cq_s, perf_press, perfRates, deferred_logger);
1826 if (this->isProducer()) {
1827 ws.phase_mixing_rates[ws.dissolved_gas] += perfRates.
dis_gas;
1828 ws.phase_mixing_rates[ws.vaporized_oil] += perfRates.
vap_oil;
1829 perf_data.phase_mixing_rates[local_perf_index][ws.dissolved_gas] = perfRates.
dis_gas;
1830 perf_data.phase_mixing_rates[local_perf_index][ws.vaporized_oil] = perfRates.
vap_oil;
1834 for (
int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
1835 perf_rates[local_perf_index*this->number_of_phases_ + this->modelCompIdxToFlowCompIdx(comp_idx)] = cq_s[comp_idx].value();
1837 perf_press_state[local_perf_index] = perf_press.value();
1840 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
1841 perf_data.wat_mass_rates[local_perf_index] = cq_s[watCompIdx].value() * rhow;
1844 for (
int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
1846 const EvalWell cq_s_effective = cq_s[comp_idx] * this->well_efficiency_factor_;
1848 this->connectionRates_[local_perf_index][comp_idx] = Base::restrictEval(cq_s_effective);
1851 assemblePerforationEq(seg, local_perf_index, comp_idx, cq_s_effective, this->linSys_);
1857 const auto& comm = this->parallel_well_info_.communication();
1858 comm.sum(ws.phase_mixing_rates.data(), ws.phase_mixing_rates.size());
1861 if (this->parallel_well_info_.communication().size() > 1) {
1863 this->linSys_.sumDistributed(this->parallel_well_info_.communication());
1865 for (
int seg = 0; seg < nseg; ++seg) {
1869 const EvalWell segment_surface_volume = getSegmentSurfaceVolume(simulator, seg);
1874 const Scalar regularization_factor = this->regularize_? this->param_.regularization_factor_wells_ : 1.0;
1876 for (
int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
1877 const EvalWell accumulation_term = regularization_factor * (segment_surface_volume * this->primary_variables_.surfaceVolumeFraction(seg, comp_idx)
1878 - segment_fluid_initial_[seg][comp_idx]) / dt;
1880 assembleAccumulationTerm(seg, comp_idx, accumulation_term, this->linSys_);
1885 const int seg_upwind = this->segments_.upwinding_segment(seg);
1886 for (
int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
1888 this->primary_variables_.getSegmentRateUpwinding(seg,
1891 this->well_efficiency_factor_;
1893 assembleOutflowTerm(seg, seg_upwind, comp_idx, segment_rate, this->linSys_);
1899 for (
const int inlet : this->segments_.inlets()[seg]) {
1900 const int inlet_upwind = this->segments_.upwinding_segment(inlet);
1901 for (
int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
1903 this->primary_variables_.getSegmentRateUpwinding(inlet,
1906 this->well_efficiency_factor_;
1908 assembleInflowTerm(seg, inlet, inlet_upwind, comp_idx, inlet_rate, this->linSys_);
1915 const auto& summaryState = simulator.vanguard().summaryState();
1916 const Schedule& schedule = simulator.vanguard().schedule();
1917 const bool stopped_or_zero_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
1919 assembleControlEq(well_state,
1926 this->primary_variables_,
1928 stopped_or_zero_target,
1931 const UnitSystem& unit_system = simulator.vanguard().eclState().getDeckUnitSystem();
1932 const auto& summary_state = simulator.vanguard().summaryState();
1933 this->assemblePressureEq(seg, unit_system, well_state, summary_state, this->param_.use_average_density_ms_wells_, deferred_logger);
1937 this->parallel_well_info_.communication().sum(this->ipr_a_.data(), this->ipr_a_.size());
1938 this->linSys_.createSolver();
1944 template<
typename TypeTag>
1949 return !this->getAllowCrossFlow() && allDrawDownWrongDirection(simulator);
1953 template<
typename TypeTag>
1958 bool all_drawdown_wrong_direction =
true;
1959 const int nseg = this->numberOfSegments();
1961 for (
int seg = 0; seg < nseg; ++seg) {
1962 const EvalWell segment_pressure = this->primary_variables_.getSegmentPressure(seg);
1963 for (
const int perf : this->segments_.perforations()[seg]) {
1964 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
1965 if (local_perf_index < 0)
1968 const int cell_idx = this->well_cells_[local_perf_index];
1969 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
1970 const auto& fs = intQuants.fluidState();
1973 const EvalWell perf_seg_press_diff = this->segments_.getPressureDiffSegLocalPerf(seg, local_perf_index);
1975 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[local_perf_index];
1977 const Scalar pressure_cell = this->getPerfCellPressure(fs).value();
1978 const Scalar perf_press = pressure_cell - cell_perf_press_diff;
1981 const EvalWell drawdown = perf_press - (segment_pressure + perf_seg_press_diff);
1986 if ( (drawdown < 0. && this->isInjector()) ||
1987 (drawdown > 0. && this->isProducer()) ) {
1988 all_drawdown_wrong_direction =
false;
1993 const auto& comm = this->parallel_well_info_.communication();
1994 if (comm.size() > 1)
1996 all_drawdown_wrong_direction =
1997 (comm.min(all_drawdown_wrong_direction ? 1 : 0) == 1);
2000 return all_drawdown_wrong_direction;
2006 template<
typename TypeTag>
2017 template<
typename TypeTag>
2025 auto info = this->getFirstPerforationFluidStateInfo(simulator);
2026 temperature.setValue(std::get<0>(info));
2027 saltConcentration = this->extendEval(std::get<1>(info));
2029 return this->segments_.getSurfaceVolume(temperature,
2031 this->primary_variables_,
2037 template<
typename TypeTag>
2038 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2042 const SummaryState& summary_state,
2048 this->getALQ(well_state),
2055 template<
typename TypeTag>
2056 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2059 const SummaryState& summary_state,
2062 bool iterate_if_no_solution)
const
2066 auto frates = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2072 std::vector<Scalar> rates(3);
2073 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2078 computeBhpAtThpLimitProd(frates,
2080 this->maxPerfPress(simulator),
2081 this->getRefDensity(),
2083 this->getTHPConstraint(summary_state),
2089 if (!iterate_if_no_solution)
2090 return std::nullopt;
2092 auto fratesIter = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2096 std::vector<Scalar> rates(3);
2097 computeWellRatesWithBhpIterations(simulator, bhp, rates, deferred_logger);
2102 computeBhpAtThpLimitProd(fratesIter,
2104 this->maxPerfPress(simulator),
2105 this->getRefDensity(),
2107 this->getTHPConstraint(summary_state),
2111 template<
typename TypeTag>
2112 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2115 const SummaryState& summary_state,
2119 auto frates = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2125 std::vector<Scalar> rates(3);
2126 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2131 computeBhpAtThpLimitInj(frates,
2133 this->getRefDensity(),
2142 auto fratesIter = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2146 std::vector<Scalar> rates(3);
2147 computeWellRatesWithBhpIterations(simulator, bhp, rates, deferred_logger);
2152 computeBhpAtThpLimitInj(fratesIter,
2154 this->getRefDensity(),
2165 template<
typename TypeTag>
2170 Scalar max_pressure = 0.0;
2171 const int nseg = this->numberOfSegments();
2172 for (
int seg = 0; seg < nseg; ++seg) {
2173 for (
const int perf : this->segments_.perforations()[seg]) {
2174 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
2175 if (local_perf_index < 0)
2178 const int cell_idx = this->well_cells_[local_perf_index];
2179 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2180 const auto& fs = int_quants.fluidState();
2181 Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2182 max_pressure = std::max(max_pressure, pressure_cell);
2185 max_pressure = this->parallel_well_info_.communication().max(max_pressure);
2186 return max_pressure;
2193 template<
typename TypeTag>
2194 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
2200 std::vector<Scalar> well_q_s(this->num_components_, 0.0);
2201 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
2202 const int nseg = this->numberOfSegments();
2203 for (
int seg = 0; seg < nseg; ++seg) {
2205 const Scalar seg_pressure = getValue(this->primary_variables_.getSegmentPressure(seg));
2206 for (
const int perf : this->segments_.perforations()[seg]) {
2207 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
2208 if (local_perf_index < 0)
2211 const int cell_idx = this->well_cells_[local_perf_index];
2212 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2213 std::vector<Scalar> mob(this->num_components_, 0.0);
2214 getMobility(simulator, local_perf_index, mob, deferred_logger);
2215 const Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quants, cell_idx);
2216 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
2217 const std::vector<Scalar> Tw = this->wellIndex(local_perf_index, int_quants, trans_mult, wellstate_nupcol);
2218 std::vector<Scalar> cq_s(this->num_components_, 0.0);
2221 computePerfRate(int_quants, mob, Tw, seg, perf, seg_pressure,
2222 allow_cf, cq_s, perf_press, perf_rates, deferred_logger);
2223 for (
int comp = 0; comp < this->num_components_; ++comp) {
2224 well_q_s[comp] += cq_s[comp];
2228 const auto& comm = this->parallel_well_info_.communication();
2229 if (comm.size() > 1)
2231 comm.sum(well_q_s.data(), well_q_s.size());
2237 template <
typename TypeTag>
2238 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
2242 const int num_seg = this->numberOfSegments();
2243 constexpr int num_eq = MSWEval::numWellEq;
2244 std::vector<Scalar> retval(num_seg * num_eq);
2245 for (
int ii = 0; ii < num_seg; ++ii) {
2246 const auto& pv = this->primary_variables_.value(ii);
2247 std::copy(pv.begin(), pv.end(), retval.begin() + ii * num_eq);
2255 template <
typename TypeTag>
2260 const int num_seg = this->numberOfSegments();
2261 constexpr int num_eq = MSWEval::numWellEq;
2262 std::array<Scalar, num_eq> tmp;
2263 for (
int ii = 0; ii < num_seg; ++ii) {
2264 const auto start = it + ii * num_eq;
2265 std::copy(start, start + num_eq, tmp.begin());
2266 this->primary_variables_.setValue(ii, tmp);
2268 return num_seg * num_eq;
2271 template <
typename TypeTag>
2275 Scalar fsTemperature = 0.0;
2276 using SaltConcType =
typename std::decay<
decltype(std::declval<
decltype(simulator.model().intensiveQuantities(0, 0).fluidState())>().saltConcentration())>::type;
2277 SaltConcType fsSaltConcentration{};
2278 int pvt_region_index = 0;
2281 if (this->well_cells_.size() > 0) {
2284 const int cell_idx = this->well_cells_[0];
2285 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
2286 const auto& fs = intQuants.fluidState();
2288 fsTemperature = fs.temperature(FluidSystem::oilPhaseIdx).value();
2289 fsSaltConcentration = fs.saltConcentration();
2290 pvt_region_index = fs.pvtRegionIndex();
2293 auto info = std::make_tuple(fsTemperature, fsSaltConcentration, pvt_region_index);
2297 return this->parallel_well_info_.communication().size() == 1 ? info : this->pw_info_.broadcastFirstPerforationValue(info);
#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 problem(const std::string &tag, const std::string &message)
void debug(const std::string &tag, const std::string &message)
Definition: GroupState.hpp:43
Class handling assemble of the equation system for MultisegmentWell.
Definition: MultisegmentWellAssemble.hpp:44
typename PrimaryVariables::EvalWell EvalWell
Definition: MultisegmentWellEval.hpp:64
PrimaryVariables primary_variables_
The primary variables.
Definition: MultisegmentWellEval.hpp:141
void scaleSegmentPressuresWithBhp(WellState< Scalar > &well_state) const
void scaleSegmentRatesWithWellRates(const std::vector< std::vector< int > > &segment_inlets, const std::vector< std::vector< int > > &segment_perforations, WellState< Scalar > &well_state) const
Definition: MultisegmentWell.hpp:38
std::vector< Scalar > computeWellPotentialWithTHP(const WellState< Scalar > &well_state, const Simulator &simulator, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:471
void computeWellPotentials(const Simulator &simulator, const WellState< Scalar > &well_state, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) override
computing the well potentials for group control
Definition: MultisegmentWell_impl.hpp:288
void calculateExplicitQuantities(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:736
std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &simulator, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger, bool iterate_if_no_solution) const override
Definition: MultisegmentWell_impl.hpp:2058
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: MultisegmentWell_impl.hpp:834
void addWellContributions(SparseMatrixAdapter &jacobian) const override
Definition: MultisegmentWell_impl.hpp:856
std::optional< Scalar > computeBhpAtThpLimitProd(const WellState< Scalar > &well_state, const Simulator &ebos_simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2040
Scalar getRefDensity() const override
Definition: MultisegmentWell_impl.hpp:1183
bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger, const bool fixed_control=false, const bool fixed_status=false) override
Definition: MultisegmentWell_impl.hpp:1601
std::vector< Scalar > getPrimaryVars() const override
Definition: MultisegmentWell_impl.hpp:2240
void solveEqAndUpdateWellState(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:592
ConvergenceReport getWellConvergence(const Simulator &simulator, const WellState< Scalar > &well_state, const std::vector< Scalar > &B_avg, DeferredLogger &deferred_logger, const bool relax_tolerance) const override
check whether the well equations get converged for this well
Definition: MultisegmentWell_impl.hpp:195
void updatePrimaryVariables(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:158
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:752
void computeWellRatesWithBhpIterations(const Simulator &simulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:396
std::vector< Scalar > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:2196
void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: MultisegmentWell_impl.hpp:221
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, WellState< Scalar > &well_state, DeferredLogger &deferred_logger, const Scalar relaxation_factor=1.0)
Definition: MultisegmentWell_impl.hpp:697
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:1143
void computePerfRate(const IntensiveQuantities &int_quants, const std::vector< Value > &mob_perfcells, const std::vector< Scalar > &Tw, const int seg, const int perf, const Value &segment_pressure, const bool &allow_cf, std::vector< Value > &cq_s, Value &perf_press, PerforationRates< Scalar > &perf_rates, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:1038
bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1477
void checkOperabilityUnderBHPLimit(const WellState< Scalar > &well_state, const Simulator &ebos_simulator, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1191
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:1947
void computeSegmentFluidProperties(const Simulator &simulator, DeferredLogger &deferred_logger)
Definition: MultisegmentWell_impl.hpp:1111
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:256
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: MultisegmentWell_impl.hpp:2258
void computeWellRatesWithBhp(const Simulator &simulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:349
EvalWell getSegmentSurfaceVolume(const Simulator &simulator, const int seg_idx) const
Definition: MultisegmentWell_impl.hpp:2020
bool allDrawDownWrongDirection(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:1956
bool computeWellPotentialsImplicit(const Simulator &simulator, const WellState< Scalar > &well_state, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:522
void init(const PhaseUsage *phase_usage_arg, 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: MultisegmentWell_impl.hpp:122
void updateWellStateWithTarget(const Simulator &simulator, const GroupState< Scalar > &group_state, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const override
updating the well state based the current control mode
Definition: MultisegmentWell_impl.hpp:174
int debug_cost_counter_
Definition: MultisegmentWell.hpp:175
void addWellPressureEquations(PressureMatrix &mat, const BVector &x, const int pressureVarIndex, const bool use_well_weights, const WellState< Scalar > &well_state) const override
Definition: MultisegmentWell_impl.hpp:869
static constexpr bool has_polymer
Definition: WellInterface.hpp:110
void assembleWellEqWithoutIteration(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1766
static constexpr bool has_solvent
Definition: WellInterface.hpp:108
void checkOperabilityUnderTHPLimit(const Simulator &ebos_simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1425
void computePerfCellPressDiffs(const Simulator &simulator)
Definition: MultisegmentWell_impl.hpp:621
void updateIPR(const Simulator &ebos_simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:1257
void updateWaterThroughput(const double dt, WellState< Scalar > &well_state) const override
Definition: MultisegmentWell_impl.hpp:2009
std::tuple< Scalar, typename std::decay< decltype(std::declval< decltype(std::declval< const Simulator & >().model().intensiveQuantities(0, 0).fluidState())>().saltConcentration())>::type, int > FSInfo
Definition: MultisegmentWell.hpp:71
FSInfo getFirstPerforationFluidStateInfo(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2273
void computeWellRatesAtBhpLimit(const Simulator &simulator, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:333
MultisegmentWell(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_components, const int num_phases, const int index_of_well, const std::vector< PerforationData< Scalar > > &perf_data)
Definition: MultisegmentWell_impl.hpp:62
Scalar maxPerfPress(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2168
void updateIPRImplicit(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1356
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &ebos_simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2114
void computeInitialSegmentFluids(const Simulator &simulator)
Definition: MultisegmentWell_impl.hpp:679
EvalWell getQs(const int comp_idx) const
Returns scaled rate for a component.
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:195
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.
Scalar mostStrictBhpFromBhpLimits(const SummaryState &summaryState) const
Obtain the most strict BHP from BHP limits.
void prepareForPotentialCalculations(const SummaryState &summary_state, WellState< Scalar > &well_state, Well::InjectionControls &inj_controls, Well::ProductionControls &prod_controls) const
std::pair< bool, bool > computeWellPotentials(std::vector< Scalar > &well_potentials, const WellState< Scalar > &well_state)
void resetDampening()
Definition: WellInterfaceGeneric.hpp:240
FluidSystem::Scalar rsRvInj() const
Well well_ecl_
Definition: WellInterfaceGeneric.hpp:297
Definition: WellInterfaceIndices.hpp:34
Definition: WellInterface.hpp:77
static constexpr bool has_brine
Definition: WellInterface.hpp:116
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:82
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:100
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:97
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1895
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:86
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
static constexpr bool has_watVapor
Definition: WellInterface.hpp:117
bool solveWellWithOperabilityCheck(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:566
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:96
typename Base::ModelParameters ModelParameters
Definition: WellInterface.hpp:106
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:84
static constexpr bool has_foam
Definition: WellInterface.hpp:115
static constexpr bool has_micp
Definition: WellInterface.hpp:120
typename Base::Eval Eval
Definition: WellInterface.hpp:95
static constexpr bool has_energy
Definition: WellInterface.hpp:111
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:85
bool thp_update_iterations
Definition: WellInterface.hpp:370
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:88
Definition: WellProdIndexCalculator.hpp:37
Scalar connectionProdIndStandard(const std::size_t connIdx, const Scalar connMobility) const
Definition: WellState.hpp:65
const SingleWellState< Scalar > & well(std::size_t well_index) const
Definition: WellState.hpp:285
@ NONE
Definition: DeferredLogger.hpp:46
Definition: blackoilboundaryratevector.hh:39
PhaseUsage phaseUsage(const Phases &phases)
Determine the active phases.
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_oil
Definition: PerforationData.hpp:44
Definition: BlackoilPhases.hpp:46
std::array< int, MaxNumPhases+NumCryptoPhases > phase_used
Definition: BlackoilPhases.hpp:51
std::array< int, MaxNumPhases+NumCryptoPhases > phase_pos
Definition: BlackoilPhases.hpp:52