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_conservation_quantities,
70 const int index_of_well,
72 :
Base(well, pw_info, time_step, param, rate_converter, pvtRegionIdx, num_conservation_quantities, num_phases, index_of_well, perf_data)
75 , segment_fluid_initial_(this->numberOfSegments(), std::vector<
Scalar>(this->num_conservation_quantities_, 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>
122 init(
const std::vector<Scalar>& depth_arg,
124 const std::vector< Scalar >& B_avg,
125 const bool changed_to_open_this_step)
127 Base::init(depth_arg, gravity_arg, B_avg, changed_to_open_this_step);
139 this->initMatrixAndVectors();
142 for (
int local_perf_index = 0; local_perf_index < this->number_of_local_perforations_; ++local_perf_index) {
144 const int cell_idx = this->well_cells_[local_perf_index];
146 this->cell_perforation_depth_diffs_[local_perf_index] = depth_arg[cell_idx] - this->perf_depth_[this->pw_info_.localPerfToActivePerf(local_perf_index)];
154 template <
typename TypeTag>
159 const auto& well_state = groupStateHelper.
wellState();
160 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(groupStateHelper);
161 this->primary_variables_.update(well_state, stop_or_zero_rate_target);
169 template <
typename TypeTag>
174 this->scaleSegmentRatesWithWellRates(this->segments_.inlets(),
175 this->segments_.perforations(),
177 this->scaleSegmentPressuresWithBhp(well_state);
180 template <
typename TypeTag>
187 Base::updateWellStateWithTarget(simulator, groupStateHelper, well_state);
190 this->scaleSegmentRatesWithWellRates(this->segments_.inlets(),
191 this->segments_.perforations(),
193 this->scaleSegmentPressuresWithBhp(well_state);
199 template <
typename TypeTag>
203 const std::vector<Scalar>& B_avg,
204 const bool relax_tolerance)
const
206 const auto& well_state = groupStateHelper.
wellState();
208 return this->MSWEval::getWellConvergence(well_state,
211 this->param_.max_residual_allowed_,
212 this->param_.tolerance_wells_,
213 this->param_.relaxed_tolerance_flow_well_,
214 this->param_.tolerance_pressure_ms_wells_,
215 this->param_.relaxed_tolerance_pressure_ms_well_,
217 this->wellIsStopped());
225 template <
typename TypeTag>
230 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) {
234 if (this->param_.matrix_add_well_contributions_) {
239 this->linSys_.apply(x, Ax);
246 template <
typename TypeTag>
251 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) {
255 this->linSys_.apply(r);
260 template <
typename TypeTag>
268 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) {
275 this->linSys_.recoverSolutionWell(x, xw);
277 updateWellState(simulator, xw, groupStateHelper, well_state);
279 catch (
const NumericalProblem& exp) {
283 deferred_logger.problem(
"In MultisegmentWell::recoverWellSolutionAndUpdateWellState for well "
284 + this->name() +
": "+exp.what());
293 template <
typename TypeTag>
299 std::vector<Scalar>& well_potentials)
302 const auto [compute_potential, bhp_controlled_well] =
305 if (!compute_potential) {
309 debug_cost_counter_ = 0;
310 bool converged_implicit =
false;
311 if (this->param_.local_well_solver_control_switching_) {
312 converged_implicit = computeWellPotentialsImplicit(simulator, groupStateHelper, well_potentials);
313 if (!converged_implicit) {
314 deferred_logger.debug(
"Implicit potential calculations failed for well "
315 + this->name() +
", reverting to original aproach.");
318 if (!converged_implicit) {
320 const auto& summaryState = simulator.vanguard().summaryState();
321 if (!Base::wellHasTHPConstraints(summaryState) || bhp_controlled_well) {
322 computeWellRatesAtBhpLimit(simulator, groupStateHelper, well_potentials);
324 well_potentials = computeWellPotentialWithTHP(
325 well_state, simulator, groupStateHelper);
328 deferred_logger.debug(
"Cost in iterations of finding well potential for well "
331 this->checkNegativeWellPotentials(well_potentials,
332 this->param_.check_well_operability_,
339 template<
typename TypeTag>
344 std::vector<Scalar>& well_flux)
const
346 if (this->well_ecl_.isInjector()) {
347 const auto controls = this->well_ecl_.injectionControls(simulator.vanguard().summaryState());
348 computeWellRatesWithBhpIterations(simulator, controls.bhp_limit, groupStateHelper, well_flux);
350 const auto controls = this->well_ecl_.productionControls(simulator.vanguard().summaryState());
351 computeWellRatesWithBhpIterations(simulator, controls.bhp_limit, groupStateHelper, well_flux);
355 template<
typename TypeTag>
360 std::vector<Scalar>& well_flux,
363 const int np = this->number_of_phases_;
365 well_flux.resize(np, 0.0);
366 const bool allow_cf = this->getAllowCrossFlow();
367 const int nseg = this->numberOfSegments();
368 const WellStateType& well_state = simulator.problem().wellModel().wellState();
369 const auto& ws = well_state.
well(this->indexOfWell());
370 auto segments_copy = ws.segments;
371 segments_copy.scale_pressure(bhp);
372 const auto& segment_pressure = segments_copy.pressure;
373 for (
int seg = 0; seg < nseg; ++seg) {
374 for (
const int perf : this->segments_.perforations()[seg]) {
375 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
376 if (local_perf_index < 0)
378 const int cell_idx = this->well_cells_[local_perf_index];
379 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
381 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.);
382 getMobility(simulator, local_perf_index, mob, deferred_logger);
384 getTransMult(trans_mult, simulator, cell_idx);
385 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
386 std::vector<Scalar> Tw(this->num_conservation_quantities_,
387 this->well_index_[local_perf_index] * trans_mult);
388 this->getTw(Tw, local_perf_index, intQuants, trans_mult, wellstate_nupcol);
389 const Scalar seg_pressure = segment_pressure[seg];
390 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.);
393 computePerfRate(intQuants, mob, Tw, seg, perf, seg_pressure,
394 allow_cf, cq_s, perf_press, perf_rates, deferred_logger);
396 for(
int p = 0; p < np; ++p) {
397 well_flux[FluidSystem::activeCompToActivePhaseIdx(p)] += cq_s[p];
401 this->parallel_well_info_.communication().sum(well_flux.data(), well_flux.size());
405 template<
typename TypeTag>
411 std::vector<Scalar>& well_flux)
const
424 auto guard = groupStateHelper_copy.
pushWellState(well_state_copy);
425 auto& ws = well_state_copy.
well(this->index_of_well_);
428 const auto& summary_state = simulator.vanguard().summaryState();
429 auto inj_controls = well_copy.
well_ecl_.isInjector()
430 ? well_copy.
well_ecl_.injectionControls(summary_state)
431 : Well::InjectionControls(0);
432 auto prod_controls = well_copy.
well_ecl_.isProducer()
433 ? well_copy.
well_ecl_.productionControls(summary_state) :
434 Well::ProductionControls(0);
438 inj_controls.bhp_limit = bhp;
439 ws.injection_cmode = Well::InjectorCMode::BHP;
441 prod_controls.bhp_limit = bhp;
442 ws.production_cmode = Well::ProducerCMode::BHP;
448 const int np = this->number_of_phases_;
450 for (
int phase = 0; phase < np; ++phase){
451 trivial = trivial && (ws.well_potentials[phase] == 0.0) ;
455 for (
int phase = 0; phase < np; ++phase) {
456 ws.surface_rates[phase] = sign * ws.well_potentials[phase];
460 this->segments_.perforations(),
464 const double dt = simulator.timeStepSize();
471 well_flux.resize(np, 0.0);
472 for (
int compIdx = 0; compIdx < this->num_conservation_quantities_; ++compIdx) {
474 well_flux[FluidSystem::activeCompToActivePhaseIdx(compIdx)] = rate.value();
481 template<
typename TypeTag>
482 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
489 std::vector<Scalar> potentials(this->number_of_phases_, 0.0);
490 const auto& summary_state = simulator.vanguard().summaryState();
492 const auto& well = this->well_ecl_;
493 if (well.isInjector()) {
494 auto bhp_at_thp_limit = computeBhpAtThpLimitInj(simulator, groupStateHelper, summary_state);
495 if (bhp_at_thp_limit) {
496 const auto& controls = well.injectionControls(summary_state);
497 const Scalar bhp = std::min(*bhp_at_thp_limit,
498 static_cast<Scalar>(controls.bhp_limit));
499 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, potentials);
500 deferred_logger.debug(
"Converged thp based potential calculation for well "
503 deferred_logger.warning(
"FAILURE_GETTING_CONVERGED_POTENTIAL",
504 "Failed in getting converged thp based potential calculation for well "
505 + this->name() +
". Instead the bhp based value is used");
506 const auto& controls = well.injectionControls(summary_state);
507 const Scalar bhp = controls.bhp_limit;
508 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, potentials);
511 auto bhp_at_thp_limit = computeBhpAtThpLimitProd(
512 well_state, simulator, groupStateHelper, summary_state);
513 if (bhp_at_thp_limit) {
514 const auto& controls = well.productionControls(summary_state);
515 const Scalar bhp = std::max(*bhp_at_thp_limit,
516 static_cast<Scalar>(controls.bhp_limit));
517 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, potentials);
518 deferred_logger.debug(
"Converged thp based potential calculation for well "
521 deferred_logger.warning(
"FAILURE_GETTING_CONVERGED_POTENTIAL",
522 "Failed in getting converged thp based potential calculation for well "
523 + this->name() +
". Instead the bhp based value is used");
524 const auto& controls = well.productionControls(summary_state);
525 const Scalar bhp = controls.bhp_limit;
526 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, potentials);
533 template<
typename TypeTag>
538 std::vector<Scalar>& well_potentials)
const
549 auto guard = groupStateHelper_copy.
pushWellState(well_state_copy);
550 auto& ws = well_state_copy.
well(this->index_of_well_);
553 const auto& summary_state = simulator.vanguard().summaryState();
554 auto inj_controls = well_copy.
well_ecl_.isInjector()
555 ? well_copy.
well_ecl_.injectionControls(summary_state)
556 : Well::InjectionControls(0);
557 auto prod_controls = well_copy.
well_ecl_.isProducer()
558 ? well_copy.
well_ecl_.productionControls(summary_state)
559 : Well::ProductionControls(0);
567 const int np = this->number_of_phases_;
569 for (
int phase = 0; phase < np; ++phase){
570 trivial = trivial && (ws.well_potentials[phase] == 0.0) ;
574 for (
int phase = 0; phase < np; ++phase) {
575 ws.surface_rates[phase] = sign * ws.well_potentials[phase];
579 this->segments_.perforations(),
583 const double dt = simulator.timeStepSize();
585 bool converged =
false;
586 if (this->well_ecl_.isProducer()) {
588 simulator, dt, inj_controls, prod_controls, groupStateHelper_copy, well_state_copy
592 simulator, dt, inj_controls, prod_controls, groupStateHelper_copy, well_state_copy,
600 well_potentials.clear();
601 well_potentials.resize(np, 0.0);
602 for (
int compIdx = 0; compIdx < this->num_conservation_quantities_; ++compIdx) {
604 well_potentials[FluidSystem::activeCompToActivePhaseIdx(compIdx)] = rate.value();
610 template <
typename TypeTag>
617 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
622 const BVectorWell dx_well = this->linSys_.solve();
623 updateWellState(simulator, dx_well, groupStateHelper, well_state);
625 catch(
const NumericalProblem& exp) {
630 deferred_logger.
problem(
"In MultisegmentWell::solveEqAndUpdateWellState for well "
631 + this->name() +
": "+exp.what());
640 template <
typename TypeTag>
647 for (
int local_perf_index = 0; local_perf_index < this->number_of_local_perforations_; ++local_perf_index) {
650 std::vector<Scalar> kr(this->number_of_phases_, 0.0);
651 std::vector<Scalar> density(this->number_of_phases_, 0.0);
653 const int cell_idx = this->well_cells_[local_perf_index];
654 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
655 const auto& fs = intQuants.fluidState();
659 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
660 const int water_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::waterPhaseIdx);
661 kr[water_pos] = intQuants.relativePermeability(FluidSystem::waterPhaseIdx).value();
662 sum_kr += kr[water_pos];
663 density[water_pos] = fs.density(FluidSystem::waterPhaseIdx).value();
666 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
667 const int oil_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::oilPhaseIdx);
668 kr[oil_pos] = intQuants.relativePermeability(FluidSystem::oilPhaseIdx).value();
669 sum_kr += kr[oil_pos];
670 density[oil_pos] = fs.density(FluidSystem::oilPhaseIdx).value();
673 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
674 const int gas_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
675 kr[gas_pos] = intQuants.relativePermeability(FluidSystem::gasPhaseIdx).value();
676 sum_kr += kr[gas_pos];
677 density[gas_pos] = fs.density(FluidSystem::gasPhaseIdx).value();
680 assert(sum_kr != 0.);
683 Scalar average_density = 0.;
684 for (
int p = 0; p < this->number_of_phases_; ++p) {
685 average_density += kr[p] * density[p];
687 average_density /= sum_kr;
689 this->cell_perforation_pressure_diffs_[local_perf_index] = this->gravity_ * average_density * this->cell_perforation_depth_diffs_[local_perf_index];
697 template <
typename TypeTag>
703 for (
int seg = 0; seg < this->numberOfSegments(); ++seg) {
705 const Scalar surface_volume = getSegmentSurfaceVolume(simulator, seg, deferred_logger).value();
706 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
707 segment_fluid_initial_[seg][comp_idx] = surface_volume * this->primary_variables_.surfaceVolumeFraction(seg, comp_idx).value();
716 template <
typename TypeTag>
720 const BVectorWell& dwells,
723 const Scalar relaxation_factor)
725 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
729 const Scalar dFLimit = this->param_.dwell_fraction_max_;
730 const Scalar max_pressure_change = this->param_.max_pressure_change_ms_wells_;
731 const bool stop_or_zero_rate_target =
732 this->stoppedOrZeroRateTarget(groupStateHelper);
733 this->primary_variables_.updateNewton(dwells,
736 stop_or_zero_rate_target,
737 max_pressure_change);
739 const auto& summary_state = simulator.vanguard().summaryState();
740 this->primary_variables_.copyToWellState(*
this, getRefDensity(),
746 auto& ws = well_state.
well(this->index_of_well_);
747 this->segments_.copyPhaseDensities(ws.segments);
750 Base::calculateReservoirRates(simulator.vanguard().eclState().runspec().co2Storage(), well_state.
well(this->index_of_well_));
757 template <
typename TypeTag>
764 updatePrimaryVariables(groupStateHelper);
765 computePerfCellPressDiffs(simulator);
766 computeInitialSegmentFluids(simulator, deferred_logger);
773 template<
typename TypeTag>
781 auto fluidState = [&simulator,
this](
const int local_perf_index)
783 const auto cell_idx = this->well_cells_[local_perf_index];
784 return simulator.model()
785 .intensiveQuantities(cell_idx, 0).fluidState();
788 const int np = this->number_of_phases_;
789 auto setToZero = [np](
Scalar* x) ->
void
791 std::fill_n(x, np, 0.0);
794 auto addVector = [np](
const Scalar* src,
Scalar* dest) ->
void
796 std::transform(src, src + np, dest, dest, std::plus<>{});
799 auto& ws = well_state.
well(this->index_of_well_);
800 auto& perf_data = ws.perf_data;
801 auto* connPI = perf_data.prod_index.data();
802 auto* wellPI = ws.productivity_index.data();
806 const auto preferred_phase = this->well_ecl_.getPreferredPhase();
807 auto subsetPerfID = 0;
809 for (
const auto& perf : *this->perf_data_){
810 auto allPerfID = perf.ecl_index;
812 auto connPICalc = [&wellPICalc, allPerfID](
const Scalar mobility) ->
Scalar
817 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
822 getMobility(simulator,
static_cast<int>(subsetPerfID), mob, deferred_logger);
824 const auto& fs = fluidState(subsetPerfID);
827 if (this->isInjector()) {
828 this->computeConnLevelInjInd(fs, preferred_phase, connPICalc,
829 mob, connPI, deferred_logger);
832 this->computeConnLevelProdInd(fs, connPICalc, mob, connPI);
835 addVector(connPI, wellPI);
842 const auto& comm = this->parallel_well_info_.communication();
843 if (comm.size() > 1) {
844 comm.sum(wellPI, np);
847 assert (
static_cast<int>(subsetPerfID) == this->number_of_local_perforations_ &&
848 "Internal logic error in processing connections for PI/II");
855 template<
typename TypeTag>
859 [[maybe_unused]]
const int openConnIdx)
const
864 const auto segNum = this->wellEcl()
865 .getConnections()[globalConnIdx].segment();
867 const auto segIdx = this->wellEcl()
868 .getSegments().segmentNumberToIndex(segNum);
870 return this->segments_.density(segIdx).value();
877 template<
typename TypeTag>
882 if (this->number_of_local_perforations_ == 0) {
886 this->linSys_.extract(jacobian);
890 template<
typename TypeTag>
895 const int pressureVarIndex,
896 const bool use_well_weights,
899 if (this->number_of_local_perforations_ == 0) {
904 this->linSys_.extractCPRPressureMatrix(jacobian,
914 template<
typename TypeTag>
915 template<
class Value>
921 const std::vector<Value>& b_perfcells,
922 const std::vector<Value>& mob_perfcells,
923 const std::vector<Value>& Tw,
925 const Value& segment_pressure,
926 const Value& segment_density,
927 const bool& allow_cf,
928 const std::vector<Value>& cmix_s,
929 std::vector<Value>& cq_s,
934 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
935 if (local_perf_index < 0)
939 const Value perf_seg_press_diff = this->gravity() * segment_density *
940 this->segments_.local_perforation_depth_diff(local_perf_index);
942 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[local_perf_index];
946 perf_press = segment_pressure + perf_seg_press_diff;
949 const Value cell_press_at_perf = pressure_cell - cell_perf_press_diff;
952 const Value drawdown = cell_press_at_perf - perf_press;
955 if (drawdown > 0.0) {
957 if (!allow_cf && this->isInjector()) {
962 for (
int comp_idx = 0; comp_idx < this->numConservationQuantities(); ++comp_idx) {
963 const Value cq_p = - Tw[comp_idx] * (mob_perfcells[comp_idx] * drawdown);
964 cq_s[comp_idx] = b_perfcells[comp_idx] * cq_p;
967 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
968 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
969 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
970 const Value cq_s_oil = cq_s[oilCompIdx];
971 const Value cq_s_gas = cq_s[gasCompIdx];
972 cq_s[gasCompIdx] += rs * cq_s_oil;
973 cq_s[oilCompIdx] += rv * cq_s_gas;
977 if (!allow_cf && this->isProducer()) {
982 Value total_mob = mob_perfcells[0];
983 for (
int comp_idx = 1; comp_idx < this->numConservationQuantities(); ++comp_idx) {
984 total_mob += mob_perfcells[comp_idx];
988 Value volume_ratio = 0.0;
989 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
990 const unsigned waterCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
991 volume_ratio += cmix_s[waterCompIdx] / b_perfcells[waterCompIdx];
994 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
995 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
996 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
1001 const Value d = 1.0 - rv * rs;
1003 if (getValue(d) == 0.0) {
1005 fmt::format(
"Zero d value obtained for well {} "
1006 "during flux calculation with rs {} and rv {}",
1007 this->name(), rs, rv),
1011 const Value tmp_oil = (cmix_s[oilCompIdx] - rv * cmix_s[gasCompIdx]) / d;
1012 volume_ratio += tmp_oil / b_perfcells[oilCompIdx];
1014 const Value tmp_gas = (cmix_s[gasCompIdx] - rs * cmix_s[oilCompIdx]) / d;
1015 volume_ratio += tmp_gas / b_perfcells[gasCompIdx];
1017 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
1018 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
1019 volume_ratio += cmix_s[oilCompIdx] / b_perfcells[oilCompIdx];
1021 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1022 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
1023 volume_ratio += cmix_s[gasCompIdx] / b_perfcells[gasCompIdx];
1027 for (
int componentIdx = 0; componentIdx < this->numConservationQuantities(); ++componentIdx) {
1028 const Value cqt_i = - Tw[componentIdx] * (total_mob * drawdown);
1029 Value cqt_is = cqt_i / volume_ratio;
1030 cq_s[componentIdx] = cmix_s[componentIdx] * cqt_is;
1035 if (this->isProducer()) {
1036 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1037 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
1038 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
1047 const Scalar d = 1.0 - getValue(rv) * getValue(rs);
1050 perf_rates.
vap_oil = getValue(rv) * (getValue(cq_s[gasCompIdx]) - getValue(rs) * getValue(cq_s[oilCompIdx])) / d;
1053 perf_rates.
dis_gas = getValue(rs) * (getValue(cq_s[oilCompIdx]) - getValue(rv) * getValue(cq_s[gasCompIdx])) / d;
1058 template <
typename TypeTag>
1059 template<
class Value>
1063 const std::vector<Value>& mob_perfcells,
1064 const std::vector<Value>& Tw,
1067 const Value& segment_pressure,
1068 const bool& allow_cf,
1069 std::vector<Value>& cq_s,
1075 auto obtain = [
this](
const Eval& value)
1077 if constexpr (std::is_same_v<Value, Scalar>) {
1078 static_cast<void>(
this);
1079 return getValue(value);
1081 return this->extendEval(value);
1084 auto obtainN = [](
const auto& value)
1086 if constexpr (std::is_same_v<Value, Scalar>) {
1087 return getValue(value);
1092 const auto& fs = int_quants.fluidState();
1094 const Value pressure_cell = obtain(this->getPerfCellPressure(fs));
1095 const Value rs = obtain(fs.Rs());
1096 const Value rv = obtain(fs.Rv());
1099 std::vector<Value> b_perfcells(this->num_conservation_quantities_, 0.0);
1101 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
1102 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1106 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
1107 b_perfcells[compIdx] = obtain(fs.invB(phaseIdx));
1110 std::vector<Value> cmix_s(this->numConservationQuantities(), 0.0);
1111 for (
int comp_idx = 0; comp_idx < this->numConservationQuantities(); ++comp_idx) {
1112 cmix_s[comp_idx] = obtainN(this->primary_variables_.surfaceVolumeFraction(seg, comp_idx));
1115 this->computePerfRate(pressure_cell,
1123 obtainN(this->segments_.density(seg)),
1132 template <
typename TypeTag>
1152 auto info = this->getFirstPerforationFluidStateInfo(simulator);
1153 temperature.setValue(std::get<0>(info));
1154 saltConcentration = this->extendEval(std::get<1>(info));
1156 this->segments_.computeFluidProperties(temperature,
1158 this->primary_variables_,
1162 template<
typename TypeTag>
1163 template<
class Value>
1168 const int cell_idx)
const
1170 auto obtain = [
this](
const Eval& value)
1172 if constexpr (std::is_same_v<Value, Scalar>) {
1173 static_cast<void>(
this);
1174 return getValue(value);
1176 return this->extendEval(value);
1182 template <
typename TypeTag>
1183 template<
class Value>
1187 const int local_perf_index,
1188 std::vector<Value>& mob,
1191 auto obtain = [
this](
const Eval& value)
1193 if constexpr (std::is_same_v<Value, Scalar>) {
1194 static_cast<void>(
this);
1195 return getValue(value);
1197 return this->extendEval(value);
1204 const auto perf_ecl_index = this->perforationData()[local_perf_index].ecl_index;
1205 const Connection& con = this->well_ecl_.getConnections()[perf_ecl_index];
1206 const int seg = this->segmentNumberToIndex(con.segment());
1210 const Scalar segment_pres = this->primary_variables_.getSegmentPressure(seg).value();
1211 const Scalar perf_seg_press_diff = this->gravity() * this->segments_.density(seg).value()
1212 * this->segments_.local_perforation_depth_diff(local_perf_index);
1213 const Scalar perf_press = segment_pres + perf_seg_press_diff;
1214 const Scalar multiplier = this->getInjMult(local_perf_index, segment_pres, perf_press, deferred_logger);
1215 for (std::size_t i = 0; i < mob.size(); ++i) {
1216 mob[i] *= multiplier;
1223 template<
typename TypeTag>
1228 return this->segments_.getRefDensity();
1231 template<
typename TypeTag>
1238 const auto& summaryState = simulator.vanguard().summaryState();
1242 const bool bhp_limit_not_defaulted = bhp_limit > 1.5 * unit::barsa;
1243 if ( bhp_limit_not_defaulted || !this->wellHasTHPConstraints(summaryState) ) {
1246 Scalar total_ipr_mass_rate = 0.0;
1247 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
1249 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1253 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
1254 const Scalar ipr_rate = this->ipr_a_[compIdx] - this->ipr_b_[compIdx] * bhp_limit;
1256 const Scalar rho = FluidSystem::referenceDensity( phaseIdx, Base::pvtRegionIdx() );
1257 total_ipr_mass_rate += ipr_rate * rho;
1259 if ( (this->isProducer() && total_ipr_mass_rate < 0.) || (this->isInjector() && total_ipr_mass_rate > 0.) ) {
1260 this->operability_status_.operable_under_only_bhp_limit =
false;
1264 if (this->operability_status_.operable_under_only_bhp_limit && this->wellHasTHPConstraints(summaryState)) {
1268 std::vector<Scalar> well_rates_bhp_limit;
1269 computeWellRatesWithBhp(simulator, bhp_limit, well_rates_bhp_limit, deferred_logger);
1271 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1274 this->getRefDensity(),
1275 this->wellEcl().alq_value(summaryState),
1278 if ( (this->isProducer() && thp < thp_limit) || (this->isInjector() && thp > thp_limit) ) {
1279 this->operability_status_.obey_thp_limit_under_bhp_limit =
false;
1290 this->operability_status_.operable_under_only_bhp_limit =
true;
1291 this->operability_status_.obey_thp_limit_under_bhp_limit =
false;
1297 template<
typename TypeTag>
1305 std::fill(this->ipr_a_.begin(), this->ipr_a_.end(), 0.);
1306 std::fill(this->ipr_b_.begin(), this->ipr_b_.end(), 0.);
1308 const int nseg = this->numberOfSegments();
1309 std::vector<Scalar> seg_dp(nseg, 0.0);
1310 for (
int seg = 0; seg < nseg; ++seg) {
1312 const Scalar dp = this->getSegmentDp(seg,
1313 this->segments_.density(seg).value(),
1316 for (
const int perf : this->segments_.perforations()[seg]) {
1317 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
1318 if (local_perf_index < 0)
1320 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
1323 getMobility(simulator, local_perf_index, mob, deferred_logger);
1325 const int cell_idx = this->well_cells_[local_perf_index];
1326 const auto& int_quantities = simulator.model().intensiveQuantities(cell_idx, 0);
1327 const auto& fs = int_quantities.fluidState();
1329 const Scalar perf_seg_press_diff = this->segments_.getPressureDiffSegLocalPerf(seg, local_perf_index);
1331 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[local_perf_index];
1332 const Scalar pressure_cell = this->getPerfCellPressure(fs).value();
1335 std::vector<Scalar> b_perf(this->num_conservation_quantities_);
1336 for (std::size_t phase = 0; phase < FluidSystem::numPhases; ++phase) {
1337 if (!FluidSystem::phaseIsActive(phase)) {
1340 const unsigned comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phase));
1341 b_perf[comp_idx] = fs.invB(phase).value();
1345 const Scalar h_perf = cell_perf_press_diff + perf_seg_press_diff + dp;
1346 const Scalar pressure_diff = pressure_cell - h_perf;
1349 if ( (this->isProducer() && pressure_diff < 0.) || (this->isInjector() && pressure_diff > 0.) ) {
1350 deferred_logger.
debug(
"CROSSFLOW_IPR",
1351 "cross flow found when updateIPR for well " + this->name());
1356 getTransMult(trans_mult, simulator, cell_idx);
1357 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1358 std::vector<Scalar> tw_perf(this->num_conservation_quantities_, this->well_index_[perf] * trans_mult);
1359 this->getTw(tw_perf, local_perf_index, int_quantities, trans_mult, wellstate_nupcol);
1360 std::vector<Scalar> ipr_a_perf(this->ipr_a_.size());
1361 std::vector<Scalar> ipr_b_perf(this->ipr_b_.size());
1362 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1363 const Scalar tw_mob = tw_perf[comp_idx] * mob[comp_idx] * b_perf[comp_idx];
1364 ipr_a_perf[comp_idx] += tw_mob * pressure_diff;
1365 ipr_b_perf[comp_idx] += tw_mob;
1369 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1370 const unsigned oil_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
1371 const unsigned gas_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
1372 const Scalar rs = (fs.Rs()).value();
1373 const Scalar rv = (fs.Rv()).value();
1375 const Scalar dis_gas_a = rs * ipr_a_perf[oil_comp_idx];
1376 const Scalar vap_oil_a = rv * ipr_a_perf[gas_comp_idx];
1378 ipr_a_perf[gas_comp_idx] += dis_gas_a;
1379 ipr_a_perf[oil_comp_idx] += vap_oil_a;
1381 const Scalar dis_gas_b = rs * ipr_b_perf[oil_comp_idx];
1382 const Scalar vap_oil_b = rv * ipr_b_perf[gas_comp_idx];
1384 ipr_b_perf[gas_comp_idx] += dis_gas_b;
1385 ipr_b_perf[oil_comp_idx] += vap_oil_b;
1388 for (std::size_t comp_idx = 0; comp_idx < ipr_a_perf.size(); ++comp_idx) {
1389 this->ipr_a_[comp_idx] += ipr_a_perf[comp_idx];
1390 this->ipr_b_[comp_idx] += ipr_b_perf[comp_idx];
1394 this->parallel_well_info_.communication().sum(this->ipr_a_.data(), this->ipr_a_.size());
1395 this->parallel_well_info_.communication().sum(this->ipr_b_.data(), this->ipr_b_.size());
1398 template<
typename TypeTag>
1413 auto rates = well_state.
well(this->index_of_well_).surface_rates;
1415 for (std::size_t p = 0; p < rates.size(); ++p) {
1416 zero_rates &= rates[p] == 0.0;
1418 auto& ws = well_state.
well(this->index_of_well_);
1420 const auto msg = fmt::format(
"updateIPRImplicit: Well {} has zero rate, IPRs might be problematic", this->name());
1421 deferred_logger.debug(msg);
1434 std::fill(ws.implicit_ipr_a.begin(), ws.implicit_ipr_a.end(), 0.);
1435 std::fill(ws.implicit_ipr_b.begin(), ws.implicit_ipr_b.end(), 0.);
1437 auto inj_controls = Well::InjectionControls(0);
1438 auto prod_controls = Well::ProductionControls(0);
1439 prod_controls.addControl(Well::ProducerCMode::BHP);
1440 prod_controls.bhp_limit = well_state.
well(this->index_of_well_).bhp;
1443 const auto cmode = ws.production_cmode;
1444 ws.production_cmode = Well::ProducerCMode::BHP;
1445 const double dt = simulator.timeStepSize();
1446 assembleWellEqWithoutIteration(simulator, groupStateHelper, dt, inj_controls, prod_controls, well_state,
1449 BVectorWell rhs(this->numberOfSegments());
1451 rhs[0][SPres] = -1.0;
1453 const BVectorWell x_well = this->linSys_.solve(rhs);
1454 constexpr int num_eq = MSWEval::numWellEq;
1455 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx){
1456 const EvalWell comp_rate = this->primary_variables_.getQs(comp_idx);
1457 const int idx = FluidSystem::activeCompToActivePhaseIdx(comp_idx);
1458 for (
size_t pvIdx = 0; pvIdx < num_eq; ++pvIdx) {
1460 ws.implicit_ipr_b[idx] -= x_well[0][pvIdx]*comp_rate.derivative(pvIdx+Indices::numEq);
1462 ws.implicit_ipr_a[idx] = ws.implicit_ipr_b[idx]*ws.bhp - comp_rate.value();
1465 ws.production_cmode = cmode;
1468 template<
typename TypeTag>
1476 const auto& summaryState = simulator.vanguard().summaryState();
1477 const auto obtain_bhp = this->isProducer()
1478 ? computeBhpAtThpLimitProd(
1479 well_state, simulator, groupStateHelper, summaryState)
1480 : computeBhpAtThpLimitInj(simulator, groupStateHelper, summaryState);
1483 this->operability_status_.can_obtain_bhp_with_thp_limit =
true;
1486 this->operability_status_.obey_bhp_limit_with_thp_limit = (*obtain_bhp >= bhp_limit);
1488 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1489 if (this->isProducer() && *obtain_bhp < thp_limit) {
1490 const std::string msg =
" obtained bhp " +
std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1491 +
" bars is SMALLER than thp limit "
1493 +
" bars as a producer for well " + this->name();
1494 deferred_logger.debug(msg);
1496 else if (this->isInjector() && *obtain_bhp > thp_limit) {
1497 const std::string msg =
" obtained bhp " +
std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1498 +
" bars is LARGER than thp limit "
1500 +
" bars as a injector for well " + this->name();
1501 deferred_logger.debug(msg);
1506 this->operability_status_.can_obtain_bhp_with_thp_limit =
false;
1507 this->operability_status_.obey_bhp_limit_with_thp_limit =
false;
1508 if (!this->wellIsStopped()) {
1509 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1510 deferred_logger.debug(
" could not find bhp value at thp limit "
1512 +
" bar for well " + this->name() +
", the well might need to be closed ");
1521 template<
typename TypeTag>
1526 const Well::InjectionControls& inj_controls,
1527 const Well::ProductionControls& prod_controls,
1531 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return true;
1535 const int max_iter_number = this->param_.max_inner_iter_ms_wells_;
1539 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1544 updatePrimaryVariables(groupStateHelper);
1546 std::vector<std::vector<Scalar> > residual_history;
1547 std::vector<Scalar> measure_history;
1550 Scalar relaxation_factor = 1.;
1551 bool converged =
false;
1552 bool relax_convergence =
false;
1553 this->regularize_ =
false;
1554 for (; it < max_iter_number; ++it, ++debug_cost_counter_) {
1556 if (it > this->param_.strict_inner_iter_wells_) {
1557 relax_convergence =
true;
1558 this->regularize_ =
true;
1561 assembleWellEqWithoutIteration(simulator, groupStateHelper, dt, inj_controls, prod_controls,
1565 const auto report = getWellConvergence(groupStateHelper, Base::B_avg_, relax_convergence);
1566 if (report.converged()) {
1573 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1577 residual_history.push_back(residuals);
1578 measure_history.push_back(this->getResidualMeasureValue(well_state,
1579 residual_history[it],
1580 this->param_.tolerance_wells_,
1581 this->param_.tolerance_pressure_ms_wells_,
1584 bool min_relaxation_reached = this->update_relaxation_factor(measure_history, relaxation_factor, this->regularize_, deferred_logger);
1585 if (min_relaxation_reached || this->repeatedStagnation(measure_history, this->regularize_, deferred_logger)) {
1587 const auto reportStag = getWellConvergence(groupStateHelper, Base::B_avg_,
true);
1588 if (reportStag.converged()) {
1590 std::string message = fmt::format(
"Well stagnates/oscillates but {} manages to get converged with relaxed tolerances in {} inner iterations."
1592 deferred_logger.debug(message);
1599 BVectorWell dx_well;
1601 dx_well = this->linSys_.solve();
1602 updateWellState(simulator, dx_well, groupStateHelper, well_state, relaxation_factor);
1604 catch(
const NumericalProblem& exp) {
1608 deferred_logger.problem(
"In MultisegmentWell::iterateWellEqWithControl for well "
1609 + this->name() +
": "+exp.what());
1616 std::ostringstream sstr;
1617 sstr <<
" Well " << this->name() <<
" converged in " << it <<
" inner iterations.";
1618 if (relax_convergence)
1619 sstr <<
" (A relaxed tolerance was used after "<< this->param_.strict_inner_iter_wells_ <<
" iterations)";
1623 deferred_logger.debug(sstr.str(), OpmLog::defaultDebugVerbosityLevel + (it == 0));
1625 std::ostringstream sstr;
1626 sstr <<
" Well " << this->name() <<
" did not converge in " << it <<
" inner iterations.";
1627#define EXTRA_DEBUG_MSW 0
1629 sstr <<
"***** Outputting the residual history for well " << this->name() <<
" during inner iterations:";
1630 for (
int i = 0; i < it; ++i) {
1631 const auto& residual = residual_history[i];
1632 sstr <<
" residual at " << i <<
"th iteration ";
1633 for (
const auto& res : residual) {
1636 sstr <<
" " << measure_history[i] <<
" \n";
1639#undef EXTRA_DEBUG_MSW
1640 deferred_logger.debug(sstr.str());
1647 template<
typename TypeTag>
1652 const Well::InjectionControls& inj_controls,
1653 const Well::ProductionControls& prod_controls,
1656 const bool fixed_control ,
1657 const bool fixed_status ,
1658 const bool solving_with_zero_rate )
1662 const int max_iter_number = this->param_.max_inner_iter_ms_wells_;
1666 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1671 updatePrimaryVariables(groupStateHelper);
1673 std::vector<std::vector<Scalar> > residual_history;
1674 std::vector<Scalar> measure_history;
1677 Scalar relaxation_factor = 1.;
1678 bool converged =
false;
1679 bool relax_convergence =
false;
1680 this->regularize_ =
false;
1681 const auto& summary_state = groupStateHelper.
summaryState();
1686 const int min_its_after_switch = 3;
1688 const int max_status_switch = this->param_.max_well_status_switch_inner_iter_;
1689 int its_since_last_switch = min_its_after_switch;
1690 int switch_count= 0;
1691 int status_switch_count = 0;
1693 const auto well_status_orig = this->wellStatus_;
1694 const auto operability_orig = this->operability_status_;
1695 auto well_status_cur = well_status_orig;
1697 const bool allow_open = well_state.
well(this->index_of_well_).status == WellStatus::OPEN;
1699 const bool allow_switching = !this->wellUnderZeroRateTarget(groupStateHelper) &&
1700 (!fixed_control || !fixed_status) && allow_open;
1701 bool final_check =
false;
1703 this->operability_status_.resetOperability();
1704 this->operability_status_.solvable =
true;
1706 for (; it < max_iter_number; ++it, ++debug_cost_counter_) {
1707 ++its_since_last_switch;
1708 if (allow_switching && its_since_last_switch >= min_its_after_switch && status_switch_count < max_status_switch){
1709 const Scalar wqTotal = this->primary_variables_.getWQTotal().value();
1710 bool changed = this->updateWellControlAndStatusLocalIteration(
1711 simulator, groupStateHelper, inj_controls, prod_controls, wqTotal,
1712 well_state, fixed_control, fixed_status,
1713 solving_with_zero_rate
1716 its_since_last_switch = 0;
1718 if (well_status_cur != this->wellStatus_) {
1719 well_status_cur = this->wellStatus_;
1720 status_switch_count++;
1723 if (!changed && final_check) {
1726 final_check =
false;
1728 if (status_switch_count == max_status_switch) {
1729 this->wellStatus_ = well_status_orig;
1733 if (it > this->param_.strict_inner_iter_wells_) {
1734 relax_convergence =
true;
1735 this->regularize_ =
true;
1738 assembleWellEqWithoutIteration(simulator, groupStateHelper, dt, inj_controls, prod_controls,
1739 well_state, solving_with_zero_rate);
1742 const auto report = getWellConvergence(groupStateHelper, Base::B_avg_, relax_convergence);
1743 converged = report.converged();
1744 if (this->parallel_well_info_.communication().size() > 1 &&
1745 this->parallel_well_info_.communication().max(converged) != this->parallel_well_info_.communication().min(converged)) {
1746 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()));
1751 if (switch_count > 0 && its_since_last_switch < min_its_after_switch) {
1753 its_since_last_switch = min_its_after_switch;
1761 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1767 residual_history.push_back(residuals);
1771 measure_history.push_back(this->getResidualMeasureValue(well_state,
1772 residual_history[it],
1773 this->param_.tolerance_wells_,
1774 this->param_.tolerance_pressure_ms_wells_,
1776 bool min_relaxation_reached = this->update_relaxation_factor(measure_history, relaxation_factor, this->regularize_, deferred_logger);
1777 if (min_relaxation_reached || this->repeatedStagnation(measure_history, this->regularize_, deferred_logger)) {
1783 const BVectorWell dx_well = this->linSys_.solve();
1784 updateWellState(simulator, dx_well, groupStateHelper, well_state, relaxation_factor);
1786 catch(
const NumericalProblem& exp) {
1790 deferred_logger.problem(
"In MultisegmentWell::iterateWellEqWithSwitching for well "
1791 + this->name() +
": "+exp.what());
1797 if (allow_switching){
1799 const bool is_stopped = this->wellIsStopped();
1800 if (this->wellHasTHPConstraints(summary_state)){
1801 this->operability_status_.can_obtain_bhp_with_thp_limit = !is_stopped;
1802 this->operability_status_.obey_thp_limit_under_bhp_limit = !is_stopped;
1804 this->operability_status_.operable_under_only_bhp_limit = !is_stopped;
1807 std::string message = fmt::format(
" Well {} converged in {} inner iterations ("
1808 "{} control/status switches).", this->name(), it, switch_count);
1809 if (relax_convergence) {
1810 message.append(fmt::format(
" (A relaxed tolerance was used after {} iterations)",
1811 this->param_.strict_inner_iter_wells_));
1813 deferred_logger.debug(message, OpmLog::defaultDebugVerbosityLevel + ((it == 0) && (switch_count == 0)));
1815 this->wellStatus_ = well_status_orig;
1816 this->operability_status_ = operability_orig;
1817 const std::string message = fmt::format(
" Well {} did not converge in {} inner iterations ("
1818 "{} switches, {} status changes).", this->name(), it, switch_count, status_switch_count);
1819 deferred_logger.debug(message);
1820 this->primary_variables_.outputLowLimitPressureSegments(deferred_logger);
1827 template<
typename TypeTag>
1833 const Well::InjectionControls& inj_controls,
1834 const Well::ProductionControls& prod_controls,
1836 const bool solving_with_zero_rate)
1838 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1843 this->segments_.updateUpwindingSegments(this->primary_variables_);
1846 computeSegmentFluidProperties(simulator, deferred_logger);
1849 this->linSys_.clear();
1851 auto& ws = well_state.
well(this->index_of_well_);
1852 ws.phase_mixing_rates.fill(0.0);
1859 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
1861 const int nseg = this->numberOfSegments();
1863 const Scalar rhow = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) ?
1864 FluidSystem::referenceDensity( FluidSystem::waterPhaseIdx, Base::pvtRegionIdx() ) : 0.0;
1865 const unsigned watCompIdx = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) ?
1866 FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx) : 0;
1868 for (
int seg = 0; seg < nseg; ++seg) {
1870 const EvalWell seg_pressure = this->primary_variables_.getSegmentPressure(seg);
1871 auto& perf_data = ws.perf_data;
1872 auto& perf_rates = perf_data.phase_rates;
1873 auto& perf_press_state = perf_data.pressure;
1874 for (
const int perf : this->segments_.perforations()[seg]) {
1875 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
1876 if (local_perf_index < 0)
1878 const int cell_idx = this->well_cells_[local_perf_index];
1879 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
1880 std::vector<EvalWell> mob(this->num_conservation_quantities_, 0.0);
1881 getMobility(simulator, local_perf_index, mob, deferred_logger);
1883 getTransMult(trans_mult, simulator, cell_idx);
1884 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1885 std::vector<EvalWell> Tw(this->num_conservation_quantities_, this->well_index_[local_perf_index] * trans_mult);
1886 this->getTw(Tw, local_perf_index, int_quants, trans_mult, wellstate_nupcol);
1887 std::vector<EvalWell> cq_s(this->num_conservation_quantities_, 0.0);
1890 computePerfRate(int_quants, mob, Tw, seg, perf, seg_pressure,
1891 allow_cf, cq_s, perf_press, perfRates, deferred_logger);
1894 if (this->isProducer()) {
1895 ws.phase_mixing_rates[ws.dissolved_gas] += perfRates.
dis_gas;
1896 ws.phase_mixing_rates[ws.vaporized_oil] += perfRates.
vap_oil;
1897 perf_data.phase_mixing_rates[local_perf_index][ws.dissolved_gas] = perfRates.
dis_gas;
1898 perf_data.phase_mixing_rates[local_perf_index][ws.vaporized_oil] = perfRates.
vap_oil;
1902 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1903 perf_rates[local_perf_index*this->number_of_phases_ + FluidSystem::activeCompToActivePhaseIdx(comp_idx)] = cq_s[comp_idx].value();
1905 perf_press_state[local_perf_index] = perf_press.value();
1908 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
1909 perf_data.wat_mass_rates[local_perf_index] = cq_s[watCompIdx].value() * rhow;
1912 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1914 const EvalWell cq_s_effective = cq_s[comp_idx] * this->well_efficiency_factor_;
1916 this->connectionRates_[local_perf_index][comp_idx] = Base::restrictEval(cq_s_effective);
1919 assemblePerforationEq(seg, local_perf_index, comp_idx, cq_s_effective, this->linSys_);
1925 const auto& comm = this->parallel_well_info_.communication();
1926 comm.sum(ws.phase_mixing_rates.data(), ws.phase_mixing_rates.size());
1929 if (this->parallel_well_info_.communication().size() > 1) {
1931 this->linSys_.sumDistributed(this->parallel_well_info_.communication());
1933 for (
int seg = 0; seg < nseg; ++seg) {
1937 const EvalWell segment_surface_volume = getSegmentSurfaceVolume(simulator, seg, deferred_logger);
1942 const Scalar regularization_factor = this->regularize_? this->param_.regularization_factor_wells_ : 1.0;
1944 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1945 const EvalWell accumulation_term = regularization_factor * (segment_surface_volume * this->primary_variables_.surfaceVolumeFraction(seg, comp_idx)
1946 - segment_fluid_initial_[seg][comp_idx]) / dt;
1948 assembleAccumulationTerm(seg, comp_idx, accumulation_term, this->linSys_);
1953 const int seg_upwind = this->segments_.upwinding_segment(seg);
1954 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1956 this->primary_variables_.getSegmentRateUpwinding(seg,
1959 this->well_efficiency_factor_;
1961 assembleOutflowTerm(seg, seg_upwind, comp_idx, segment_rate, this->linSys_);
1967 for (
const int inlet : this->segments_.inlets()[seg]) {
1968 const int inlet_upwind = this->segments_.upwinding_segment(inlet);
1969 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1971 this->primary_variables_.getSegmentRateUpwinding(inlet,
1974 this->well_efficiency_factor_;
1976 assembleInflowTerm(seg, inlet, inlet_upwind, comp_idx, inlet_rate, this->linSys_);
1983 const bool stopped_or_zero_target = this->stoppedOrZeroRateTarget(groupStateHelper);
1991 auto& group_state = solving_with_zero_rate
1995 auto group_guard = groupStateHelper_copy.
pushGroupState(group_state);
1997 assembleControlEq(groupStateHelper_copy,
2000 this->getRefDensity(),
2001 this->primary_variables_,
2003 stopped_or_zero_target);
2005 const UnitSystem& unit_system = simulator.vanguard().eclState().getDeckUnitSystem();
2006 const auto& summary_state = simulator.vanguard().summaryState();
2007 this->assemblePressureEq(seg, unit_system, well_state, summary_state, this->param_.use_average_density_ms_wells_, deferred_logger);
2011 this->parallel_well_info_.communication().sum(this->ipr_a_.data(), this->ipr_a_.size());
2012 this->linSys_.createSolver();
2018 template<
typename TypeTag>
2023 return !this->getAllowCrossFlow() && allDrawDownWrongDirection(simulator);
2027 template<
typename TypeTag>
2032 bool all_drawdown_wrong_direction =
true;
2033 const int nseg = this->numberOfSegments();
2035 for (
int seg = 0; seg < nseg; ++seg) {
2036 const EvalWell segment_pressure = this->primary_variables_.getSegmentPressure(seg);
2037 for (
const int perf : this->segments_.perforations()[seg]) {
2038 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
2039 if (local_perf_index < 0)
2042 const int cell_idx = this->well_cells_[local_perf_index];
2043 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
2044 const auto& fs = intQuants.fluidState();
2047 const EvalWell perf_seg_press_diff = this->segments_.getPressureDiffSegLocalPerf(seg, local_perf_index);
2049 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[local_perf_index];
2051 const Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2052 const Scalar perf_press = pressure_cell - cell_perf_press_diff;
2055 const EvalWell drawdown = perf_press - (segment_pressure + perf_seg_press_diff);
2060 if ( (drawdown < 0. && this->isInjector()) ||
2061 (drawdown > 0. && this->isProducer()) ) {
2062 all_drawdown_wrong_direction =
false;
2067 const auto& comm = this->parallel_well_info_.communication();
2068 if (comm.size() > 1)
2070 all_drawdown_wrong_direction =
2071 (comm.min(all_drawdown_wrong_direction ? 1 : 0) == 1);
2074 return all_drawdown_wrong_direction;
2080 template<
typename TypeTag>
2091 template<
typename TypeTag>
2101 auto info = this->getFirstPerforationFluidStateInfo(simulator);
2102 temperature.setValue(std::get<0>(info));
2103 saltConcentration = this->extendEval(std::get<1>(info));
2105 return this->segments_.getSurfaceVolume(temperature,
2107 this->primary_variables_,
2113 template<
typename TypeTag>
2114 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2119 const SummaryState& summary_state)
const
2125 this->getALQ(well_state),
2131 template<
typename TypeTag>
2132 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2136 const SummaryState& summary_state,
2138 bool iterate_if_no_solution)
const
2143 auto frates = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2149 std::vector<Scalar> rates(3);
2150 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2155 computeBhpAtThpLimitProd(frates,
2157 maxPerfPress(simulator),
2158 this->getRefDensity(),
2160 this->getTHPConstraint(summary_state),
2166 if (!iterate_if_no_solution)
2167 return std::nullopt;
2169 auto fratesIter = [
this, &simulator, &groupStateHelper, &deferred_logger](
const Scalar bhp) {
2173 std::vector<Scalar> rates(3);
2174 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, rates);
2179 computeBhpAtThpLimitProd(fratesIter,
2181 maxPerfPress(simulator),
2182 this->getRefDensity(),
2184 this->getTHPConstraint(summary_state),
2188 template<
typename TypeTag>
2189 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2193 const SummaryState& summary_state)
const
2197 auto frates = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2203 std::vector<Scalar> rates(3);
2204 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2209 computeBhpAtThpLimitInj(frates,
2211 this->getRefDensity(),
2220 auto fratesIter = [
this, &simulator, &groupStateHelper, &deferred_logger](
const Scalar bhp) {
2224 std::vector<Scalar> rates(3);
2225 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, rates);
2230 computeBhpAtThpLimitInj(fratesIter,
2232 this->getRefDensity(),
2243 template<
typename TypeTag>
2248 Scalar max_pressure = 0.0;
2249 const int nseg = this->numberOfSegments();
2250 for (
int seg = 0; seg < nseg; ++seg) {
2251 for (
const int perf : this->segments_.perforations()[seg]) {
2252 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
2253 if (local_perf_index < 0)
2256 const int cell_idx = this->well_cells_[local_perf_index];
2257 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2258 const auto& fs = int_quants.fluidState();
2259 Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2260 max_pressure = std::max(max_pressure, pressure_cell);
2263 max_pressure = this->parallel_well_info_.communication().max(max_pressure);
2264 return max_pressure;
2271 template<
typename TypeTag>
2272 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
2278 std::vector<Scalar> well_q_s(this->num_conservation_quantities_, 0.0);
2279 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
2280 const int nseg = this->numberOfSegments();
2281 for (
int seg = 0; seg < nseg; ++seg) {
2283 const Scalar seg_pressure = getValue(this->primary_variables_.getSegmentPressure(seg));
2284 for (
const int perf : this->segments_.perforations()[seg]) {
2285 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
2286 if (local_perf_index < 0)
2289 const int cell_idx = this->well_cells_[local_perf_index];
2290 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2291 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
2292 getMobility(simulator, local_perf_index, mob, deferred_logger);
2294 getTransMult(trans_mult, simulator, cell_idx);
2295 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
2296 std::vector<Scalar> Tw(this->num_conservation_quantities_, this->well_index_[local_perf_index] * trans_mult);
2297 this->getTw(Tw, local_perf_index, int_quants, trans_mult, wellstate_nupcol);
2298 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.0);
2301 computePerfRate(int_quants, mob, Tw, seg, perf, seg_pressure,
2302 allow_cf, cq_s, perf_press, perf_rates, deferred_logger);
2303 for (
int comp = 0; comp < this->num_conservation_quantities_; ++comp) {
2304 well_q_s[comp] += cq_s[comp];
2308 const auto& comm = this->parallel_well_info_.communication();
2309 if (comm.size() > 1)
2311 comm.sum(well_q_s.data(), well_q_s.size());
2317 template <
typename TypeTag>
2318 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
2322 const int num_seg = this->numberOfSegments();
2323 constexpr int num_eq = MSWEval::numWellEq;
2324 std::vector<Scalar> retval(num_seg * num_eq);
2325 for (
int ii = 0; ii < num_seg; ++ii) {
2326 const auto& pv = this->primary_variables_.value(ii);
2327 std::copy(pv.begin(), pv.end(), retval.begin() + ii * num_eq);
2335 template <
typename TypeTag>
2340 const int num_seg = this->numberOfSegments();
2341 constexpr int num_eq = MSWEval::numWellEq;
2342 std::array<Scalar, num_eq> tmp;
2343 for (
int ii = 0; ii < num_seg; ++ii) {
2344 const auto start = it + ii * num_eq;
2345 std::copy(start, start + num_eq, tmp.begin());
2346 this->primary_variables_.setValue(ii, tmp);
2348 return num_seg * num_eq;
2351 template <
typename TypeTag>
2356 Scalar fsTemperature = 0.0;
2357 using SaltConcType =
typename std::decay<
decltype(std::declval<
decltype(simulator.model().intensiveQuantities(0, 0).fluidState())>().saltConcentration())>::type;
2358 SaltConcType fsSaltConcentration{};
2361 if (this->well_cells_.size() > 0) {
2364 const int cell_idx = this->well_cells_[0];
2365 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
2366 const auto& fs = intQuants.fluidState();
2368 fsTemperature = fs.temperature(FluidSystem::oilPhaseIdx).value();
2369 fsSaltConcentration = fs.saltConcentration();
2372 auto info = std::make_tuple(fsTemperature, fsSaltConcentration);
2376 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 problem(const std::string &tag, const std::string &message)
void debug(const std::string &tag, const std::string &message)
Definition: GroupStateHelper.hpp:54
GroupState< Scalar > & groupState() const
Definition: GroupStateHelper.hpp:272
const SummaryState & summaryState() const
Definition: GroupStateHelper.hpp:387
const WellState< Scalar, IndexTraits > & wellState() const
Definition: GroupStateHelper.hpp:440
DeferredLogger & deferredLogger() const
Get the deferred logger.
Definition: GroupStateHelper.hpp:226
WellStateGuard pushWellState(WellState< Scalar, IndexTraits > &well_state)
Definition: GroupStateHelper.hpp:326
GroupStateGuard pushGroupState(GroupState< Scalar > &group_state)
Definition: GroupStateHelper.hpp:303
Definition: GroupState.hpp:41
Class handling assemble of the equation system for MultisegmentWell.
Definition: MultisegmentWellAssemble.hpp:45
typename PrimaryVariables::EvalWell EvalWell
Definition: MultisegmentWellEval.hpp:66
PrimaryVariables primary_variables_
The primary variables.
Definition: MultisegmentWellEval.hpp:140
void scaleSegmentRatesWithWellRates(const std::vector< std::vector< int > > &segment_inlets, const std::vector< std::vector< int > > &segment_perforations, WellState< Scalar, IndexTraits > &well_state) const
void scaleSegmentPressuresWithBhp(WellState< Scalar, IndexTraits > &well_state) const
Definition: MultisegmentWell.hpp:38
bool computeWellPotentialsImplicit(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_potentials) const
Definition: MultisegmentWell_impl.hpp:536
bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, const bool fixed_control, const bool fixed_status, const bool solving_with_zero_rate) override
Definition: MultisegmentWell_impl.hpp:1650
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, const Scalar relaxation_factor=1.0)
Definition: MultisegmentWell_impl.hpp:719
void updateWaterThroughput(const double dt, WellStateType &well_state) const override
Definition: MultisegmentWell_impl.hpp:2083
void addWellPressureEquations(PressureMatrix &mat, const BVector &x, const int pressureVarIndex, const bool use_well_weights, const WellStateType &well_state) const override
Definition: MultisegmentWell_impl.hpp:893
void assembleWellEqWithoutIteration(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellStateType &well_state, const bool solving_with_zero_rate) override
Definition: MultisegmentWell_impl.hpp:1830
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: MultisegmentWell_impl.hpp:858
void addWellContributions(SparseMatrixAdapter &jacobian) const override
Definition: MultisegmentWell_impl.hpp:880
void getTransMult(Value &trans_mult, const Simulator &simulator, const int cell_indx) const
Definition: MultisegmentWell_impl.hpp:1166
std::vector< Scalar > computeWellPotentialWithTHP(const WellStateType &well_state, const Simulator &simulator, const GroupStateHelperType &groupStateHelper) const
Definition: MultisegmentWell_impl.hpp:484
Scalar getRefDensity() const override
Definition: MultisegmentWell_impl.hpp:1226
void updateWellStateWithTarget(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) const override
updating the well state based the current control mode
Definition: MultisegmentWell_impl.hpp:183
std::vector< Scalar > getPrimaryVars() const override
Definition: MultisegmentWell_impl.hpp:2320
void computeWellRatesWithBhpIterations(const Simulator &simulator, const Scalar &bhp, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_flux) const override
Definition: MultisegmentWell_impl.hpp:408
void checkOperabilityUnderTHPLimit(const Simulator &ebos_simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper) override
Definition: MultisegmentWell_impl.hpp:1471
bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) override
Definition: MultisegmentWell_impl.hpp:1524
std::vector< Scalar > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:2274
void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: MultisegmentWell_impl.hpp:228
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_conservation_quantities, const int num_phases, const int index_of_well, const std::vector< PerforationData< Scalar > > &perf_data)
Definition: MultisegmentWell_impl.hpp:62
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
void checkOperabilityUnderBHPLimit(const WellStateType &well_state, const Simulator &ebos_simulator, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1234
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:1186
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) override
Definition: MultisegmentWell_impl.hpp:263
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2021
void computeSegmentFluidProperties(const Simulator &simulator, DeferredLogger &deferred_logger)
Definition: MultisegmentWell_impl.hpp:1135
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: MultisegmentWell_impl.hpp:2338
void computeWellRatesWithBhp(const Simulator &simulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:358
void scaleSegmentRatesAndPressure(WellStateType &well_state) const override
updating the segment pressure and rates based the current bhp and well rates
Definition: MultisegmentWell_impl.hpp:172
bool allDrawDownWrongDirection(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2030
int debug_cost_counter_
Definition: MultisegmentWell.hpp:174
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:776
static constexpr bool has_polymer
Definition: WellInterface.hpp:115
std::optional< Scalar > computeBhpAtThpLimitProd(const WellStateType &well_state, const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state) const
Definition: MultisegmentWell_impl.hpp:2116
Scalar maxPerfPress(const Simulator &simulator) const override
Definition: MultisegmentWell_impl.hpp:2246
void computePerfRate(const IntensiveQuantities &int_quants, const std::vector< Value > &mob_perfcells, const std::vector< Value > &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:1062
static constexpr bool has_solvent
Definition: WellInterface.hpp:113
void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_potentials) override
computing the well potentials for group control
Definition: MultisegmentWell_impl.hpp:296
void calculateExplicitQuantities(const Simulator &simulator, const GroupStateHelperType &groupStateHelper) override
Definition: MultisegmentWell_impl.hpp:760
void computePerfCellPressDiffs(const Simulator &simulator)
Definition: MultisegmentWell_impl.hpp:643
void solveEqAndUpdateWellState(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) override
Definition: MultisegmentWell_impl.hpp:613
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state) const
Definition: MultisegmentWell_impl.hpp:2191
void updateIPR(const Simulator &ebos_simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:1300
ConvergenceReport getWellConvergence(const GroupStateHelperType &groupStateHelper, const std::vector< Scalar > &B_avg, const bool relax_tolerance) const override
check whether the well equations get converged for this well
Definition: MultisegmentWell_impl.hpp:202
void updatePrimaryVariables(const GroupStateHelperType &groupStateHelper) override
Definition: MultisegmentWell_impl.hpp:157
void updateIPRImplicit(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) override
Definition: MultisegmentWell_impl.hpp:1401
void computeWellRatesAtBhpLimit(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_flux) const
Definition: MultisegmentWell_impl.hpp:342
FSInfo getFirstPerforationFluidStateInfo(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2354
std::tuple< Scalar, typename std::decay< decltype(std::declval< decltype(std::declval< const Simulator & >().model().intensiveQuantities(0, 0).fluidState())>().saltConcentration())>::type > FSInfo
Definition: MultisegmentWell.hpp:74
void init(const std::vector< Scalar > &depth_arg, const Scalar gravity_arg, const std::vector< Scalar > &B_avg, const bool changed_to_open_this_step) override
Definition: MultisegmentWell_impl.hpp:122
void computeInitialSegmentFluids(const Simulator &simulator, DeferredLogger &deferred_logger)
Definition: MultisegmentWell_impl.hpp:700
EvalWell getSegmentSurfaceVolume(const Simulator &simulator, const int seg_idx, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2094
std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, const Scalar alq_value, bool iterate_if_no_solution) const override
Definition: MultisegmentWell_impl.hpp:2134
EvalWell getQs(const int comp_idx) const
Returns scaled rate for a component.
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:198
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.
Well well_ecl_
Definition: WellInterfaceGeneric.hpp:304
void onlyKeepBHPandTHPcontrols(const SummaryState &summary_state, WellStateType &well_state, Well::InjectionControls &inj_controls, Well::ProductionControls &prod_controls) const
FluidSystem::Scalar rsRvInj() const
void resetDampening()
Definition: WellInterfaceGeneric.hpp:247
std::pair< bool, bool > computeWellPotentials(std::vector< Scalar > &well_potentials, const WellStateType &well_state)
Definition: WellInterfaceIndices.hpp:34
Definition: WellInterface.hpp:77
bool solveWellWithOperabilityCheck(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:590
static constexpr bool has_brine
Definition: WellInterface.hpp:122
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:82
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:105
void getTransMult(Value &trans_mult, const Simulator &simulator, const int cell_idx, Callback &extendEval) const
Definition: WellInterface_impl.hpp:2053
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:98
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:2066
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:87
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
static constexpr bool has_watVapor
Definition: WellInterface.hpp:123
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:97
typename Base::ModelParameters ModelParameters
Definition: WellInterface.hpp:111
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:84
static constexpr bool has_foam
Definition: WellInterface.hpp:121
static constexpr bool has_micp
Definition: WellInterface.hpp:127
typename Base::Eval Eval
Definition: WellInterface.hpp:96
static constexpr bool has_energy
Definition: WellInterface.hpp:117
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:86
bool thp_update_iterations
Definition: WellInterface.hpp:365
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:89
Definition: WellProdIndexCalculator.hpp:37
Scalar connectionProdIndStandard(const std::size_t connIdx, const Scalar connMobility) const
Definition: WellState.hpp:66
const SingleWellState< Scalar, IndexTraits > & well(std::size_t well_index) const
Definition: WellState.hpp:290
@ NONE
Definition: DeferredLogger.hpp:46
Definition: blackoilbioeffectsmodules.hh:43
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Static data associated with a well perforation.
Definition: PerforationData.hpp:30
Definition: PerforationData.hpp:41
Scalar dis_gas
Definition: PerforationData.hpp:42
Scalar vap_oil
Definition: PerforationData.hpp:44