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>
161 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
162 this->primary_variables_.update(well_state, stop_or_zero_rate_target);
170 template <
typename TypeTag>
175 this->scaleSegmentRatesWithWellRates(this->segments_.inlets(),
176 this->segments_.perforations(),
178 this->scaleSegmentPressuresWithBhp(well_state);
181 template <
typename TypeTag>
189 Base::updateWellStateWithTarget(simulator, wgHelper, well_state, deferred_logger);
192 this->scaleSegmentRatesWithWellRates(this->segments_.inlets(),
193 this->segments_.perforations(),
195 this->scaleSegmentPressuresWithBhp(well_state);
201 template <
typename TypeTag>
206 const std::vector<Scalar>& B_avg,
208 const bool relax_tolerance)
const
210 return this->MSWEval::getWellConvergence(well_state,
213 this->param_.max_residual_allowed_,
214 this->param_.tolerance_wells_,
215 this->param_.relaxed_tolerance_flow_well_,
216 this->param_.tolerance_pressure_ms_wells_,
217 this->param_.relaxed_tolerance_pressure_ms_well_,
219 this->wellIsStopped());
227 template <
typename TypeTag>
232 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) {
236 if (this->param_.matrix_add_well_contributions_) {
241 this->linSys_.apply(x, Ax);
248 template <
typename TypeTag>
253 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) {
257 this->linSys_.apply(r);
262 template <
typename TypeTag>
270 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) {
276 this->linSys_.recoverSolutionWell(x, xw);
278 updateWellState(simulator, xw, well_state, deferred_logger);
280 catch (
const NumericalProblem& exp) {
284 deferred_logger.
problem(
"In MultisegmentWell::recoverWellSolutionAndUpdateWellState for well "
285 + this->name() +
": "+exp.what());
294 template <
typename TypeTag>
300 std::vector<Scalar>& well_potentials,
303 const auto [compute_potential, bhp_controlled_well] =
306 if (!compute_potential) {
310 debug_cost_counter_ = 0;
311 bool converged_implicit =
false;
312 if (this->param_.local_well_solver_control_switching_) {
313 converged_implicit = computeWellPotentialsImplicit(simulator, wgHelper, well_potentials, deferred_logger);
314 if (!converged_implicit) {
315 deferred_logger.
debug(
"Implicit potential calculations failed for well "
316 + this->name() +
", reverting to original aproach.");
319 if (!converged_implicit) {
321 const auto& summaryState = simulator.vanguard().summaryState();
322 if (!Base::wellHasTHPConstraints(summaryState) || bhp_controlled_well) {
323 computeWellRatesAtBhpLimit(simulator, wgHelper, well_potentials, deferred_logger);
325 well_potentials = computeWellPotentialWithTHP(
326 well_state, simulator, wgHelper, deferred_logger);
329 deferred_logger.
debug(
"Cost in iterations of finding well potential for well "
332 this->checkNegativeWellPotentials(well_potentials,
333 this->param_.check_well_operability_,
340 template<
typename TypeTag>
345 std::vector<Scalar>& well_flux,
348 if (this->well_ecl_.isInjector()) {
349 const auto controls = this->well_ecl_.injectionControls(simulator.vanguard().summaryState());
350 computeWellRatesWithBhpIterations(simulator, controls.bhp_limit, wgHelper, well_flux, deferred_logger);
352 const auto controls = this->well_ecl_.productionControls(simulator.vanguard().summaryState());
353 computeWellRatesWithBhpIterations(simulator, controls.bhp_limit, wgHelper, well_flux, deferred_logger);
357 template<
typename TypeTag>
362 std::vector<Scalar>& well_flux,
365 const int np = this->number_of_phases_;
367 well_flux.resize(np, 0.0);
368 const bool allow_cf = this->getAllowCrossFlow();
369 const int nseg = this->numberOfSegments();
370 const WellStateType& well_state = simulator.problem().wellModel().wellState();
371 const auto& ws = well_state.
well(this->indexOfWell());
372 auto segments_copy = ws.segments;
373 segments_copy.scale_pressure(bhp);
374 const auto& segment_pressure = segments_copy.pressure;
375 for (
int seg = 0; seg < nseg; ++seg) {
376 for (
const int perf : this->segments_.perforations()[seg]) {
377 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
378 if (local_perf_index < 0)
380 const int cell_idx = this->well_cells_[local_perf_index];
381 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
383 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.);
384 getMobility(simulator, local_perf_index, mob, deferred_logger);
385 const Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(intQuants, cell_idx);
386 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
387 const std::vector<Scalar> Tw = this->wellIndex(local_perf_index, intQuants, trans_mult, wellstate_nupcol);
388 const Scalar seg_pressure = segment_pressure[seg];
389 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.);
392 computePerfRate(intQuants, mob, Tw, seg, perf, seg_pressure,
393 allow_cf, cq_s, perf_press, perf_rates, deferred_logger);
395 for(
int p = 0; p < np; ++p) {
396 well_flux[FluidSystem::activeCompToActivePhaseIdx(p)] += cq_s[p];
400 this->parallel_well_info_.communication().sum(well_flux.data(), well_flux.size());
404 template<
typename TypeTag>
410 std::vector<Scalar>& well_flux,
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();
467 well_state_copy, deferred_logger);
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, wgHelper, summary_state, deferred_logger);
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, wgHelper, potentials, deferred_logger);
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, wgHelper, potentials, deferred_logger);
511 auto bhp_at_thp_limit = computeBhpAtThpLimitProd(
512 well_state, simulator, wgHelper, summary_state, deferred_logger);
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, wgHelper, potentials, deferred_logger);
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, wgHelper, potentials, deferred_logger);
533 template<
typename TypeTag>
538 std::vector<Scalar>& well_potentials,
551 auto& ws = well_state_copy.
well(this->index_of_well_);
554 const auto& summary_state = simulator.vanguard().summaryState();
555 auto inj_controls = well_copy.
well_ecl_.isInjector()
556 ? well_copy.
well_ecl_.injectionControls(summary_state)
557 : Well::InjectionControls(0);
558 auto prod_controls = well_copy.
well_ecl_.isProducer()
559 ? well_copy.
well_ecl_.productionControls(summary_state)
560 : Well::ProductionControls(0);
568 const int np = this->number_of_phases_;
570 for (
int phase = 0; phase < np; ++phase){
571 trivial = trivial && (ws.well_potentials[phase] == 0.0) ;
575 for (
int phase = 0; phase < np; ++phase) {
576 ws.surface_rates[phase] = sign * ws.well_potentials[phase];
580 this->segments_.perforations(),
584 const double dt = simulator.timeStepSize();
586 bool converged =
false;
587 if (this->well_ecl_.isProducer()) {
589 simulator, dt, inj_controls, prod_controls, wgHelper_copy, well_state_copy, deferred_logger
593 simulator, dt, inj_controls, prod_controls, wgHelper_copy, well_state_copy, deferred_logger
598 well_potentials.clear();
599 well_potentials.resize(np, 0.0);
600 for (
int compIdx = 0; compIdx < this->num_conservation_quantities_; ++compIdx) {
602 well_potentials[FluidSystem::activeCompToActivePhaseIdx(compIdx)] = rate.value();
608 template <
typename TypeTag>
615 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
620 const BVectorWell dx_well = this->linSys_.solve();
621 updateWellState(simulator, dx_well, well_state, deferred_logger);
623 catch(
const NumericalProblem& exp) {
627 deferred_logger.
problem(
"In MultisegmentWell::solveEqAndUpdateWellState for well "
628 + this->name() +
": "+exp.what());
637 template <
typename TypeTag>
644 for (
int local_perf_index = 0; local_perf_index < this->number_of_local_perforations_; ++local_perf_index) {
647 std::vector<Scalar> kr(this->number_of_phases_, 0.0);
648 std::vector<Scalar> density(this->number_of_phases_, 0.0);
650 const int cell_idx = this->well_cells_[local_perf_index];
651 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
652 const auto& fs = intQuants.fluidState();
656 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
657 const int water_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::waterPhaseIdx);
658 kr[water_pos] = intQuants.relativePermeability(FluidSystem::waterPhaseIdx).value();
659 sum_kr += kr[water_pos];
660 density[water_pos] = fs.density(FluidSystem::waterPhaseIdx).value();
663 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
664 const int oil_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::oilPhaseIdx);
665 kr[oil_pos] = intQuants.relativePermeability(FluidSystem::oilPhaseIdx).value();
666 sum_kr += kr[oil_pos];
667 density[oil_pos] = fs.density(FluidSystem::oilPhaseIdx).value();
670 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
671 const int gas_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
672 kr[gas_pos] = intQuants.relativePermeability(FluidSystem::gasPhaseIdx).value();
673 sum_kr += kr[gas_pos];
674 density[gas_pos] = fs.density(FluidSystem::gasPhaseIdx).value();
677 assert(sum_kr != 0.);
680 Scalar average_density = 0.;
681 for (
int p = 0; p < this->number_of_phases_; ++p) {
682 average_density += kr[p] * density[p];
684 average_density /= sum_kr;
686 this->cell_perforation_pressure_diffs_[local_perf_index] = this->gravity_ * average_density * this->cell_perforation_depth_diffs_[local_perf_index];
694 template <
typename TypeTag>
700 for (
int seg = 0; seg < this->numberOfSegments(); ++seg) {
702 const Scalar surface_volume = getSegmentSurfaceVolume(simulator, seg, deferred_logger).value();
703 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
704 segment_fluid_initial_[seg][comp_idx] = surface_volume * this->primary_variables_.surfaceVolumeFraction(seg, comp_idx).value();
713 template <
typename TypeTag>
717 const BVectorWell& dwells,
720 const Scalar relaxation_factor)
722 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
724 const Scalar dFLimit = this->param_.dwell_fraction_max_;
725 const Scalar max_pressure_change = this->param_.max_pressure_change_ms_wells_;
726 const bool stop_or_zero_rate_target =
727 this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
728 this->primary_variables_.updateNewton(dwells,
731 stop_or_zero_rate_target,
732 max_pressure_change);
734 const auto& summary_state = simulator.vanguard().summaryState();
735 this->primary_variables_.copyToWellState(*
this, getRefDensity(),
741 auto& ws = well_state.
well(this->index_of_well_);
742 this->segments_.copyPhaseDensities(ws.segments);
745 Base::calculateReservoirRates(simulator.vanguard().eclState().runspec().co2Storage(), well_state.
well(this->index_of_well_));
752 template <
typename TypeTag>
759 updatePrimaryVariables(simulator, well_state, deferred_logger);
760 computePerfCellPressDiffs(simulator);
761 computeInitialSegmentFluids(simulator, deferred_logger);
768 template<
typename TypeTag>
776 auto fluidState = [&simulator,
this](
const int local_perf_index)
778 const auto cell_idx = this->well_cells_[local_perf_index];
779 return simulator.model()
780 .intensiveQuantities(cell_idx, 0).fluidState();
783 const int np = this->number_of_phases_;
784 auto setToZero = [np](
Scalar* x) ->
void
786 std::fill_n(x, np, 0.0);
789 auto addVector = [np](
const Scalar* src,
Scalar* dest) ->
void
791 std::transform(src, src + np, dest, dest, std::plus<>{});
794 auto& ws = well_state.
well(this->index_of_well_);
795 auto& perf_data = ws.perf_data;
796 auto* connPI = perf_data.prod_index.data();
797 auto* wellPI = ws.productivity_index.data();
801 const auto preferred_phase = this->well_ecl_.getPreferredPhase();
802 auto subsetPerfID = 0;
804 for (
const auto& perf : *this->perf_data_){
805 auto allPerfID = perf.ecl_index;
807 auto connPICalc = [&wellPICalc, allPerfID](
const Scalar mobility) ->
Scalar
812 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
817 getMobility(simulator,
static_cast<int>(subsetPerfID), mob, deferred_logger);
819 const auto& fs = fluidState(subsetPerfID);
822 if (this->isInjector()) {
823 this->computeConnLevelInjInd(fs, preferred_phase, connPICalc,
824 mob, connPI, deferred_logger);
827 this->computeConnLevelProdInd(fs, connPICalc, mob, connPI);
830 addVector(connPI, wellPI);
837 const auto& comm = this->parallel_well_info_.communication();
838 if (comm.size() > 1) {
839 comm.sum(wellPI, np);
842 assert (
static_cast<int>(subsetPerfID) == this->number_of_local_perforations_ &&
843 "Internal logic error in processing connections for PI/II");
850 template<
typename TypeTag>
854 [[maybe_unused]]
const int openConnIdx)
const
859 const auto segNum = this->wellEcl()
860 .getConnections()[globalConnIdx].segment();
862 const auto segIdx = this->wellEcl()
863 .getSegments().segmentNumberToIndex(segNum);
865 return this->segments_.density(segIdx).value();
872 template<
typename TypeTag>
877 if (this->number_of_local_perforations_ == 0) {
881 this->linSys_.extract(jacobian);
885 template<
typename TypeTag>
890 const int pressureVarIndex,
891 const bool use_well_weights,
894 if (this->number_of_local_perforations_ == 0) {
899 this->linSys_.extractCPRPressureMatrix(jacobian,
909 template<
typename TypeTag>
910 template<
class Value>
916 const std::vector<Value>& b_perfcells,
917 const std::vector<Value>& mob_perfcells,
918 const std::vector<Scalar>& Tw,
920 const Value& segment_pressure,
921 const Value& segment_density,
922 const bool& allow_cf,
923 const std::vector<Value>& cmix_s,
924 std::vector<Value>& cq_s,
929 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
930 if (local_perf_index < 0)
934 const Value perf_seg_press_diff = this->gravity() * segment_density *
935 this->segments_.local_perforation_depth_diff(local_perf_index);
937 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[local_perf_index];
941 perf_press = segment_pressure + perf_seg_press_diff;
944 const Value cell_press_at_perf = pressure_cell - cell_perf_press_diff;
947 const Value drawdown = cell_press_at_perf - perf_press;
950 if (drawdown > 0.0) {
952 if (!allow_cf && this->isInjector()) {
957 for (
int comp_idx = 0; comp_idx < this->numConservationQuantities(); ++comp_idx) {
958 const Value cq_p = - Tw[comp_idx] * (mob_perfcells[comp_idx] * drawdown);
959 cq_s[comp_idx] = b_perfcells[comp_idx] * cq_p;
962 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
963 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
964 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
965 const Value cq_s_oil = cq_s[oilCompIdx];
966 const Value cq_s_gas = cq_s[gasCompIdx];
967 cq_s[gasCompIdx] += rs * cq_s_oil;
968 cq_s[oilCompIdx] += rv * cq_s_gas;
972 if (!allow_cf && this->isProducer()) {
977 Value total_mob = mob_perfcells[0];
978 for (
int comp_idx = 1; comp_idx < this->numConservationQuantities(); ++comp_idx) {
979 total_mob += mob_perfcells[comp_idx];
983 Value volume_ratio = 0.0;
984 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
985 const unsigned waterCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
986 volume_ratio += cmix_s[waterCompIdx] / b_perfcells[waterCompIdx];
989 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
990 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
991 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
996 const Value d = 1.0 - rv * rs;
998 if (getValue(d) == 0.0) {
1000 fmt::format(
"Zero d value obtained for well {} "
1001 "during flux calculation with rs {} and rv {}",
1002 this->name(), rs, rv),
1006 const Value tmp_oil = (cmix_s[oilCompIdx] - rv * cmix_s[gasCompIdx]) / d;
1007 volume_ratio += tmp_oil / b_perfcells[oilCompIdx];
1009 const Value tmp_gas = (cmix_s[gasCompIdx] - rs * cmix_s[oilCompIdx]) / d;
1010 volume_ratio += tmp_gas / b_perfcells[gasCompIdx];
1012 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
1013 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
1014 volume_ratio += cmix_s[oilCompIdx] / b_perfcells[oilCompIdx];
1016 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1017 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
1018 volume_ratio += cmix_s[gasCompIdx] / b_perfcells[gasCompIdx];
1022 for (
int componentIdx = 0; componentIdx < this->numConservationQuantities(); ++componentIdx) {
1023 const Value cqt_i = - Tw[componentIdx] * (total_mob * drawdown);
1024 Value cqt_is = cqt_i / volume_ratio;
1025 cq_s[componentIdx] = cmix_s[componentIdx] * cqt_is;
1030 if (this->isProducer()) {
1031 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1032 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
1033 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
1042 const Scalar d = 1.0 - getValue(rv) * getValue(rs);
1045 perf_rates.
vap_oil = getValue(rv) * (getValue(cq_s[gasCompIdx]) - getValue(rs) * getValue(cq_s[oilCompIdx])) / d;
1048 perf_rates.
dis_gas = getValue(rs) * (getValue(cq_s[oilCompIdx]) - getValue(rv) * getValue(cq_s[gasCompIdx])) / d;
1053 template <
typename TypeTag>
1054 template<
class Value>
1058 const std::vector<Value>& mob_perfcells,
1059 const std::vector<Scalar>& Tw,
1062 const Value& segment_pressure,
1063 const bool& allow_cf,
1064 std::vector<Value>& cq_s,
1070 auto obtain = [
this](
const Eval& value)
1072 if constexpr (std::is_same_v<Value, Scalar>) {
1073 static_cast<void>(
this);
1074 return getValue(value);
1076 return this->extendEval(value);
1079 auto obtainN = [](
const auto& value)
1081 if constexpr (std::is_same_v<Value, Scalar>) {
1082 return getValue(value);
1087 const auto& fs = int_quants.fluidState();
1089 const Value pressure_cell = obtain(this->getPerfCellPressure(fs));
1090 const Value rs = obtain(fs.Rs());
1091 const Value rv = obtain(fs.Rv());
1094 std::vector<Value> b_perfcells(this->num_conservation_quantities_, 0.0);
1096 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
1097 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1101 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
1102 b_perfcells[compIdx] = obtain(fs.invB(phaseIdx));
1105 std::vector<Value> cmix_s(this->numConservationQuantities(), 0.0);
1106 for (
int comp_idx = 0; comp_idx < this->numConservationQuantities(); ++comp_idx) {
1107 cmix_s[comp_idx] = obtainN(this->primary_variables_.surfaceVolumeFraction(seg, comp_idx));
1110 this->computePerfRate(pressure_cell,
1118 obtainN(this->segments_.density(seg)),
1127 template <
typename TypeTag>
1147 auto info = this->getFirstPerforationFluidStateInfo(simulator);
1148 temperature.setValue(std::get<0>(info));
1149 saltConcentration = this->extendEval(std::get<1>(info));
1151 this->segments_.computeFluidProperties(temperature,
1153 this->primary_variables_,
1157 template <
typename TypeTag>
1158 template<
class Value>
1162 const int local_perf_index,
1163 std::vector<Value>& mob,
1166 auto obtain = [
this](
const Eval& value)
1168 if constexpr (std::is_same_v<Value, Scalar>) {
1169 static_cast<void>(
this);
1170 return getValue(value);
1172 return this->extendEval(value);
1179 const auto perf_ecl_index = this->perforationData()[local_perf_index].ecl_index;
1180 const Connection& con = this->well_ecl_.getConnections()[perf_ecl_index];
1181 const int seg = this->segmentNumberToIndex(con.segment());
1185 const Scalar segment_pres = this->primary_variables_.getSegmentPressure(seg).value();
1186 const Scalar perf_seg_press_diff = this->gravity() * this->segments_.density(seg).value()
1187 * this->segments_.local_perforation_depth_diff(local_perf_index);
1188 const Scalar perf_press = segment_pres + perf_seg_press_diff;
1189 const Scalar multiplier = this->getInjMult(local_perf_index, segment_pres, perf_press, deferred_logger);
1190 for (std::size_t i = 0; i < mob.size(); ++i) {
1191 mob[i] *= multiplier;
1198 template<
typename TypeTag>
1203 return this->segments_.getRefDensity();
1206 template<
typename TypeTag>
1213 const auto& summaryState = simulator.vanguard().summaryState();
1217 const bool bhp_limit_not_defaulted = bhp_limit > 1.5 * unit::barsa;
1218 if ( bhp_limit_not_defaulted || !this->wellHasTHPConstraints(summaryState) ) {
1221 Scalar total_ipr_mass_rate = 0.0;
1222 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
1224 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1228 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
1229 const Scalar ipr_rate = this->ipr_a_[compIdx] - this->ipr_b_[compIdx] * bhp_limit;
1231 const Scalar rho = FluidSystem::referenceDensity( phaseIdx, Base::pvtRegionIdx() );
1232 total_ipr_mass_rate += ipr_rate * rho;
1234 if ( (this->isProducer() && total_ipr_mass_rate < 0.) || (this->isInjector() && total_ipr_mass_rate > 0.) ) {
1235 this->operability_status_.operable_under_only_bhp_limit =
false;
1239 if (this->operability_status_.operable_under_only_bhp_limit && this->wellHasTHPConstraints(summaryState)) {
1243 std::vector<Scalar> well_rates_bhp_limit;
1244 computeWellRatesWithBhp(simulator, bhp_limit, well_rates_bhp_limit, deferred_logger);
1246 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1249 this->getRefDensity(),
1250 this->wellEcl().alq_value(summaryState),
1253 if ( (this->isProducer() && thp < thp_limit) || (this->isInjector() && thp > thp_limit) ) {
1254 this->operability_status_.obey_thp_limit_under_bhp_limit =
false;
1265 this->operability_status_.operable_under_only_bhp_limit =
true;
1266 this->operability_status_.obey_thp_limit_under_bhp_limit =
false;
1272 template<
typename TypeTag>
1280 std::fill(this->ipr_a_.begin(), this->ipr_a_.end(), 0.);
1281 std::fill(this->ipr_b_.begin(), this->ipr_b_.end(), 0.);
1283 const int nseg = this->numberOfSegments();
1284 std::vector<Scalar> seg_dp(nseg, 0.0);
1285 for (
int seg = 0; seg < nseg; ++seg) {
1287 const Scalar dp = this->getSegmentDp(seg,
1288 this->segments_.density(seg).value(),
1291 for (
const int perf : this->segments_.perforations()[seg]) {
1292 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
1293 if (local_perf_index < 0)
1295 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
1298 getMobility(simulator, local_perf_index, mob, deferred_logger);
1300 const int cell_idx = this->well_cells_[local_perf_index];
1301 const auto& int_quantities = simulator.model().intensiveQuantities(cell_idx, 0);
1302 const auto& fs = int_quantities.fluidState();
1304 const Scalar perf_seg_press_diff = this->segments_.getPressureDiffSegLocalPerf(seg, local_perf_index);
1306 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[local_perf_index];
1307 const Scalar pressure_cell = this->getPerfCellPressure(fs).value();
1310 std::vector<Scalar> b_perf(this->num_conservation_quantities_);
1311 for (std::size_t phase = 0; phase < FluidSystem::numPhases; ++phase) {
1312 if (!FluidSystem::phaseIsActive(phase)) {
1315 const unsigned comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phase));
1316 b_perf[comp_idx] = fs.invB(phase).value();
1320 const Scalar h_perf = cell_perf_press_diff + perf_seg_press_diff + dp;
1321 const Scalar pressure_diff = pressure_cell - h_perf;
1324 if ( (this->isProducer() && pressure_diff < 0.) || (this->isInjector() && pressure_diff > 0.) ) {
1325 deferred_logger.
debug(
"CROSSFLOW_IPR",
1326 "cross flow found when updateIPR for well " + this->name());
1330 const Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quantities, cell_idx);
1331 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1332 const std::vector<Scalar> tw_perf = this->wellIndex(local_perf_index, int_quantities, trans_mult, wellstate_nupcol);
1333 std::vector<Scalar> ipr_a_perf(this->ipr_a_.size());
1334 std::vector<Scalar> ipr_b_perf(this->ipr_b_.size());
1335 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1336 const Scalar tw_mob = tw_perf[comp_idx] * mob[comp_idx] * b_perf[comp_idx];
1337 ipr_a_perf[comp_idx] += tw_mob * pressure_diff;
1338 ipr_b_perf[comp_idx] += tw_mob;
1342 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1343 const unsigned oil_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
1344 const unsigned gas_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
1345 const Scalar rs = (fs.Rs()).value();
1346 const Scalar rv = (fs.Rv()).value();
1348 const Scalar dis_gas_a = rs * ipr_a_perf[oil_comp_idx];
1349 const Scalar vap_oil_a = rv * ipr_a_perf[gas_comp_idx];
1351 ipr_a_perf[gas_comp_idx] += dis_gas_a;
1352 ipr_a_perf[oil_comp_idx] += vap_oil_a;
1354 const Scalar dis_gas_b = rs * ipr_b_perf[oil_comp_idx];
1355 const Scalar vap_oil_b = rv * ipr_b_perf[gas_comp_idx];
1357 ipr_b_perf[gas_comp_idx] += dis_gas_b;
1358 ipr_b_perf[oil_comp_idx] += vap_oil_b;
1361 for (std::size_t comp_idx = 0; comp_idx < ipr_a_perf.size(); ++comp_idx) {
1362 this->ipr_a_[comp_idx] += ipr_a_perf[comp_idx];
1363 this->ipr_b_[comp_idx] += ipr_b_perf[comp_idx];
1367 this->parallel_well_info_.communication().sum(this->ipr_a_.data(), this->ipr_a_.size());
1368 this->parallel_well_info_.communication().sum(this->ipr_b_.data(), this->ipr_b_.size());
1371 template<
typename TypeTag>
1385 auto rates = well_state.
well(this->index_of_well_).surface_rates;
1387 for (std::size_t p = 0; p < rates.size(); ++p) {
1388 zero_rates &= rates[p] == 0.0;
1390 auto& ws = well_state.
well(this->index_of_well_);
1392 const auto msg = fmt::format(
"updateIPRImplicit: Well {} has zero rate, IPRs might be problematic", this->name());
1393 deferred_logger.
debug(msg);
1405 const auto& group_state = simulator.problem().wellModel().groupState();
1407 std::fill(ws.implicit_ipr_a.begin(), ws.implicit_ipr_a.end(), 0.);
1408 std::fill(ws.implicit_ipr_b.begin(), ws.implicit_ipr_b.end(), 0.);
1410 auto inj_controls = Well::InjectionControls(0);
1411 auto prod_controls = Well::ProductionControls(0);
1412 prod_controls.addControl(Well::ProducerCMode::BHP);
1413 prod_controls.bhp_limit = well_state.
well(this->index_of_well_).bhp;
1416 const auto cmode = ws.production_cmode;
1417 ws.production_cmode = Well::ProducerCMode::BHP;
1418 const double dt = simulator.timeStepSize();
1419 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
1421 BVectorWell rhs(this->numberOfSegments());
1423 rhs[0][SPres] = -1.0;
1425 const BVectorWell x_well = this->linSys_.solve(rhs);
1426 constexpr int num_eq = MSWEval::numWellEq;
1427 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx){
1428 const EvalWell comp_rate = this->primary_variables_.getQs(comp_idx);
1429 const int idx = FluidSystem::activeCompToActivePhaseIdx(comp_idx);
1430 for (
size_t pvIdx = 0; pvIdx < num_eq; ++pvIdx) {
1432 ws.implicit_ipr_b[idx] -= x_well[0][pvIdx]*comp_rate.derivative(pvIdx+Indices::numEq);
1434 ws.implicit_ipr_a[idx] = ws.implicit_ipr_b[idx]*ws.bhp - comp_rate.value();
1437 ws.production_cmode = cmode;
1440 template<
typename TypeTag>
1448 const auto& summaryState = simulator.vanguard().summaryState();
1449 const auto obtain_bhp = this->isProducer()
1450 ? computeBhpAtThpLimitProd(
1451 well_state, simulator, wgHelper, summaryState, deferred_logger)
1452 : computeBhpAtThpLimitInj(simulator, wgHelper, summaryState, deferred_logger);
1455 this->operability_status_.can_obtain_bhp_with_thp_limit =
true;
1458 this->operability_status_.obey_bhp_limit_with_thp_limit = (*obtain_bhp >= bhp_limit);
1460 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1461 if (this->isProducer() && *obtain_bhp < thp_limit) {
1462 const std::string msg =
" obtained bhp " +
std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1463 +
" bars is SMALLER than thp limit "
1465 +
" bars as a producer for well " + this->name();
1466 deferred_logger.
debug(msg);
1468 else if (this->isInjector() && *obtain_bhp > thp_limit) {
1469 const std::string msg =
" obtained bhp " +
std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1470 +
" bars is LARGER than thp limit "
1472 +
" bars as a injector for well " + this->name();
1473 deferred_logger.
debug(msg);
1478 this->operability_status_.can_obtain_bhp_with_thp_limit =
false;
1479 this->operability_status_.obey_bhp_limit_with_thp_limit =
false;
1480 if (!this->wellIsStopped()) {
1481 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1482 deferred_logger.
debug(
" could not find bhp value at thp limit "
1484 +
" bar for well " + this->name() +
", the well might need to be closed ");
1493 template<
typename TypeTag>
1498 const Well::InjectionControls& inj_controls,
1499 const Well::ProductionControls& prod_controls,
1504 const auto& group_state = wgHelper.
groupState();
1505 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return true;
1507 const int max_iter_number = this->param_.max_inner_iter_ms_wells_;
1511 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1516 updatePrimaryVariables(simulator, well_state, deferred_logger);
1518 std::vector<std::vector<Scalar> > residual_history;
1519 std::vector<Scalar> measure_history;
1522 Scalar relaxation_factor = 1.;
1523 bool converged =
false;
1524 bool relax_convergence =
false;
1525 this->regularize_ =
false;
1526 for (; it < max_iter_number; ++it, ++debug_cost_counter_) {
1528 if (it > this->param_.strict_inner_iter_wells_) {
1529 relax_convergence =
true;
1530 this->regularize_ =
true;
1533 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls,
1534 well_state, group_state, deferred_logger);
1536 const auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
1537 if (report.converged()) {
1544 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1548 residual_history.push_back(residuals);
1549 measure_history.push_back(this->getResidualMeasureValue(well_state,
1550 residual_history[it],
1551 this->param_.tolerance_wells_,
1552 this->param_.tolerance_pressure_ms_wells_,
1555 bool min_relaxation_reached = this->update_relaxation_factor(measure_history, relaxation_factor, this->regularize_, deferred_logger);
1556 if (min_relaxation_reached || this->repeatedStagnation(measure_history, this->regularize_, deferred_logger)) {
1558 const auto reportStag = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger,
true);
1559 if (reportStag.converged()) {
1561 std::string message = fmt::format(
"Well stagnates/oscillates but {} manages to get converged with relaxed tolerances in {} inner iterations."
1563 deferred_logger.
debug(message);
1570 BVectorWell dx_well;
1572 dx_well = this->linSys_.solve();
1573 updateWellState(simulator, dx_well, well_state, deferred_logger, relaxation_factor);
1575 catch(
const NumericalProblem& exp) {
1579 deferred_logger.
problem(
"In MultisegmentWell::iterateWellEqWithControl for well "
1580 + this->name() +
": "+exp.what());
1587 std::ostringstream sstr;
1588 sstr <<
" Well " << this->name() <<
" converged in " << it <<
" inner iterations.";
1589 if (relax_convergence)
1590 sstr <<
" (A relaxed tolerance was used after "<< this->param_.strict_inner_iter_wells_ <<
" iterations)";
1594 deferred_logger.
debug(sstr.str(), OpmLog::defaultDebugVerbosityLevel + (it == 0));
1596 std::ostringstream sstr;
1597 sstr <<
" Well " << this->name() <<
" did not converge in " << it <<
" inner iterations.";
1598#define EXTRA_DEBUG_MSW 0
1600 sstr <<
"***** Outputting the residual history for well " << this->name() <<
" during inner iterations:";
1601 for (
int i = 0; i < it; ++i) {
1602 const auto& residual = residual_history[i];
1603 sstr <<
" residual at " << i <<
"th iteration ";
1604 for (
const auto& res : residual) {
1607 sstr <<
" " << measure_history[i] <<
" \n";
1610#undef EXTRA_DEBUG_MSW
1611 deferred_logger.
debug(sstr.str());
1618 template<
typename TypeTag>
1623 const Well::InjectionControls& inj_controls,
1624 const Well::ProductionControls& prod_controls,
1628 const bool fixed_control ,
1629 const bool fixed_status )
1631 const auto& group_state = wgHelper.
groupState();
1632 const int max_iter_number = this->param_.max_inner_iter_ms_wells_;
1636 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1641 updatePrimaryVariables(simulator, well_state, deferred_logger);
1643 std::vector<std::vector<Scalar> > residual_history;
1644 std::vector<Scalar> measure_history;
1647 Scalar relaxation_factor = 1.;
1648 bool converged =
false;
1649 bool relax_convergence =
false;
1650 this->regularize_ =
false;
1651 const auto& summary_state = simulator.vanguard().summaryState();
1656 const int min_its_after_switch = 3;
1658 const int max_status_switch = this->param_.max_well_status_switch_inner_iter_;
1659 int its_since_last_switch = min_its_after_switch;
1660 int switch_count= 0;
1661 int status_switch_count = 0;
1663 const auto well_status_orig = this->wellStatus_;
1664 const auto operability_orig = this->operability_status_;
1665 auto well_status_cur = well_status_orig;
1667 const bool allow_open = well_state.
well(this->index_of_well_).status == WellStatus::OPEN;
1669 const bool allow_switching = !this->wellUnderZeroRateTarget(simulator, well_state, deferred_logger) &&
1670 (!fixed_control || !fixed_status) && allow_open;
1671 bool final_check =
false;
1673 this->operability_status_.resetOperability();
1674 this->operability_status_.solvable =
true;
1676 for (; it < max_iter_number; ++it, ++debug_cost_counter_) {
1677 ++its_since_last_switch;
1678 if (allow_switching && its_since_last_switch >= min_its_after_switch && status_switch_count < max_status_switch){
1679 const Scalar wqTotal = this->primary_variables_.getWQTotal().value();
1680 bool changed = this->updateWellControlAndStatusLocalIteration(
1681 simulator, wgHelper, inj_controls, prod_controls, wqTotal,
1682 well_state, deferred_logger, fixed_control, fixed_status
1685 its_since_last_switch = 0;
1687 if (well_status_cur != this->wellStatus_) {
1688 well_status_cur = this->wellStatus_;
1689 status_switch_count++;
1692 if (!changed && final_check) {
1695 final_check =
false;
1697 if (status_switch_count == max_status_switch) {
1698 this->wellStatus_ = well_status_orig;
1702 if (it > this->param_.strict_inner_iter_wells_) {
1703 relax_convergence =
true;
1704 this->regularize_ =
true;
1707 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls,
1708 well_state, group_state, deferred_logger);
1711 const auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
1712 converged = report.converged();
1713 if (this->parallel_well_info_.communication().size() > 1 &&
1714 this->parallel_well_info_.communication().max(converged) != this->parallel_well_info_.communication().min(converged)) {
1715 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()));
1720 if (switch_count > 0 && its_since_last_switch < min_its_after_switch) {
1722 its_since_last_switch = min_its_after_switch;
1730 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1736 residual_history.push_back(residuals);
1740 measure_history.push_back(this->getResidualMeasureValue(well_state,
1741 residual_history[it],
1742 this->param_.tolerance_wells_,
1743 this->param_.tolerance_pressure_ms_wells_,
1745 bool min_relaxation_reached = this->update_relaxation_factor(measure_history, relaxation_factor, this->regularize_, deferred_logger);
1746 if (min_relaxation_reached || this->repeatedStagnation(measure_history, this->regularize_, deferred_logger)) {
1752 const BVectorWell dx_well = this->linSys_.solve();
1753 updateWellState(simulator, dx_well, well_state, deferred_logger, relaxation_factor);
1755 catch(
const NumericalProblem& exp) {
1759 deferred_logger.
problem(
"In MultisegmentWell::iterateWellEqWithSwitching for well "
1760 + this->name() +
": "+exp.what());
1766 if (allow_switching){
1768 const bool is_stopped = this->wellIsStopped();
1769 if (this->wellHasTHPConstraints(summary_state)){
1770 this->operability_status_.can_obtain_bhp_with_thp_limit = !is_stopped;
1771 this->operability_status_.obey_thp_limit_under_bhp_limit = !is_stopped;
1773 this->operability_status_.operable_under_only_bhp_limit = !is_stopped;
1776 std::string message = fmt::format(
" Well {} converged in {} inner iterations ("
1777 "{} control/status switches).", this->name(), it, switch_count);
1778 if (relax_convergence) {
1779 message.append(fmt::format(
" (A relaxed tolerance was used after {} iterations)",
1780 this->param_.strict_inner_iter_wells_));
1782 deferred_logger.
debug(message, OpmLog::defaultDebugVerbosityLevel + ((it == 0) && (switch_count == 0)));
1784 this->wellStatus_ = well_status_orig;
1785 this->operability_status_ = operability_orig;
1786 const std::string message = fmt::format(
" Well {} did not converge in {} inner iterations ("
1787 "{} switches, {} status changes).", this->name(), it, switch_count, status_switch_count);
1788 deferred_logger.
debug(message);
1789 this->primary_variables_.outputLowLimitPressureSegments(deferred_logger);
1796 template<
typename TypeTag>
1801 const Well::InjectionControls& inj_controls,
1802 const Well::ProductionControls& prod_controls,
1807 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1810 this->segments_.updateUpwindingSegments(this->primary_variables_);
1813 computeSegmentFluidProperties(simulator, deferred_logger);
1816 this->linSys_.clear();
1818 auto& ws = well_state.
well(this->index_of_well_);
1819 ws.phase_mixing_rates.fill(0.0);
1826 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
1828 const int nseg = this->numberOfSegments();
1830 const Scalar rhow = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) ?
1831 FluidSystem::referenceDensity( FluidSystem::waterPhaseIdx, Base::pvtRegionIdx() ) : 0.0;
1832 const unsigned watCompIdx = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) ?
1833 FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx) : 0;
1835 for (
int seg = 0; seg < nseg; ++seg) {
1837 const EvalWell seg_pressure = this->primary_variables_.getSegmentPressure(seg);
1838 auto& perf_data = ws.perf_data;
1839 auto& perf_rates = perf_data.phase_rates;
1840 auto& perf_press_state = perf_data.pressure;
1841 for (
const int perf : this->segments_.perforations()[seg]) {
1842 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
1843 if (local_perf_index < 0)
1845 const int cell_idx = this->well_cells_[local_perf_index];
1846 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
1847 std::vector<EvalWell> mob(this->num_conservation_quantities_, 0.0);
1848 getMobility(simulator, local_perf_index, mob, deferred_logger);
1849 const Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quants, cell_idx);
1850 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1851 const std::vector<Scalar> Tw = this->wellIndex(local_perf_index, int_quants, trans_mult, wellstate_nupcol);
1852 std::vector<EvalWell> cq_s(this->num_conservation_quantities_, 0.0);
1855 computePerfRate(int_quants, mob, Tw, seg, perf, seg_pressure,
1856 allow_cf, cq_s, perf_press, perfRates, deferred_logger);
1859 if (this->isProducer()) {
1860 ws.phase_mixing_rates[ws.dissolved_gas] += perfRates.
dis_gas;
1861 ws.phase_mixing_rates[ws.vaporized_oil] += perfRates.
vap_oil;
1862 perf_data.phase_mixing_rates[local_perf_index][ws.dissolved_gas] = perfRates.
dis_gas;
1863 perf_data.phase_mixing_rates[local_perf_index][ws.vaporized_oil] = perfRates.
vap_oil;
1867 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1868 perf_rates[local_perf_index*this->number_of_phases_ + FluidSystem::activeCompToActivePhaseIdx(comp_idx)] = cq_s[comp_idx].value();
1870 perf_press_state[local_perf_index] = perf_press.value();
1873 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
1874 perf_data.wat_mass_rates[local_perf_index] = cq_s[watCompIdx].value() * rhow;
1877 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1879 const EvalWell cq_s_effective = cq_s[comp_idx] * this->well_efficiency_factor_;
1881 this->connectionRates_[local_perf_index][comp_idx] = Base::restrictEval(cq_s_effective);
1884 assemblePerforationEq(seg, local_perf_index, comp_idx, cq_s_effective, this->linSys_);
1890 const auto& comm = this->parallel_well_info_.communication();
1891 comm.sum(ws.phase_mixing_rates.data(), ws.phase_mixing_rates.size());
1894 if (this->parallel_well_info_.communication().size() > 1) {
1896 this->linSys_.sumDistributed(this->parallel_well_info_.communication());
1898 for (
int seg = 0; seg < nseg; ++seg) {
1902 const EvalWell segment_surface_volume = getSegmentSurfaceVolume(simulator, seg, deferred_logger);
1907 const Scalar regularization_factor = this->regularize_? this->param_.regularization_factor_wells_ : 1.0;
1909 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1910 const EvalWell accumulation_term = regularization_factor * (segment_surface_volume * this->primary_variables_.surfaceVolumeFraction(seg, comp_idx)
1911 - segment_fluid_initial_[seg][comp_idx]) / dt;
1913 assembleAccumulationTerm(seg, comp_idx, accumulation_term, this->linSys_);
1918 const int seg_upwind = this->segments_.upwinding_segment(seg);
1919 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1921 this->primary_variables_.getSegmentRateUpwinding(seg,
1924 this->well_efficiency_factor_;
1926 assembleOutflowTerm(seg, seg_upwind, comp_idx, segment_rate, this->linSys_);
1932 for (
const int inlet : this->segments_.inlets()[seg]) {
1933 const int inlet_upwind = this->segments_.upwinding_segment(inlet);
1934 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1936 this->primary_variables_.getSegmentRateUpwinding(inlet,
1939 this->well_efficiency_factor_;
1941 assembleInflowTerm(seg, inlet, inlet_upwind, comp_idx, inlet_rate, this->linSys_);
1948 const auto& summaryState = simulator.vanguard().summaryState();
1949 const Schedule& schedule = simulator.vanguard().schedule();
1950 const bool stopped_or_zero_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
1952 assembleControlEq(well_state,
1959 this->primary_variables_,
1961 stopped_or_zero_target,
1964 const UnitSystem& unit_system = simulator.vanguard().eclState().getDeckUnitSystem();
1965 const auto& summary_state = simulator.vanguard().summaryState();
1966 this->assemblePressureEq(seg, unit_system, well_state, summary_state, this->param_.use_average_density_ms_wells_, deferred_logger);
1970 this->parallel_well_info_.communication().sum(this->ipr_a_.data(), this->ipr_a_.size());
1971 this->linSys_.createSolver();
1977 template<
typename TypeTag>
1982 return !this->getAllowCrossFlow() && allDrawDownWrongDirection(simulator);
1986 template<
typename TypeTag>
1991 bool all_drawdown_wrong_direction =
true;
1992 const int nseg = this->numberOfSegments();
1994 for (
int seg = 0; seg < nseg; ++seg) {
1995 const EvalWell segment_pressure = this->primary_variables_.getSegmentPressure(seg);
1996 for (
const int perf : this->segments_.perforations()[seg]) {
1997 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
1998 if (local_perf_index < 0)
2001 const int cell_idx = this->well_cells_[local_perf_index];
2002 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
2003 const auto& fs = intQuants.fluidState();
2006 const EvalWell perf_seg_press_diff = this->segments_.getPressureDiffSegLocalPerf(seg, local_perf_index);
2008 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[local_perf_index];
2010 const Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2011 const Scalar perf_press = pressure_cell - cell_perf_press_diff;
2014 const EvalWell drawdown = perf_press - (segment_pressure + perf_seg_press_diff);
2019 if ( (drawdown < 0. && this->isInjector()) ||
2020 (drawdown > 0. && this->isProducer()) ) {
2021 all_drawdown_wrong_direction =
false;
2026 const auto& comm = this->parallel_well_info_.communication();
2027 if (comm.size() > 1)
2029 all_drawdown_wrong_direction =
2030 (comm.min(all_drawdown_wrong_direction ? 1 : 0) == 1);
2033 return all_drawdown_wrong_direction;
2039 template<
typename TypeTag>
2050 template<
typename TypeTag>
2060 auto info = this->getFirstPerforationFluidStateInfo(simulator);
2061 temperature.setValue(std::get<0>(info));
2062 saltConcentration = this->extendEval(std::get<1>(info));
2064 return this->segments_.getSurfaceVolume(temperature,
2066 this->primary_variables_,
2072 template<
typename TypeTag>
2073 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2078 const SummaryState& summary_state,
2085 this->getALQ(well_state),
2092 template<
typename TypeTag>
2093 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2097 const SummaryState& summary_state,
2100 bool iterate_if_no_solution)
const
2104 auto frates = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2110 std::vector<Scalar> rates(3);
2111 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2116 computeBhpAtThpLimitProd(frates,
2118 this->maxPerfPress(simulator),
2119 this->getRefDensity(),
2121 this->getTHPConstraint(summary_state),
2127 if (!iterate_if_no_solution)
2128 return std::nullopt;
2130 auto fratesIter = [
this, &simulator, &wgHelper, &deferred_logger](
const Scalar bhp) {
2134 std::vector<Scalar> rates(3);
2135 computeWellRatesWithBhpIterations(simulator, bhp, wgHelper, rates, deferred_logger);
2140 computeBhpAtThpLimitProd(fratesIter,
2142 this->maxPerfPress(simulator),
2143 this->getRefDensity(),
2145 this->getTHPConstraint(summary_state),
2149 template<
typename TypeTag>
2150 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2154 const SummaryState& summary_state,
2158 auto frates = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2164 std::vector<Scalar> rates(3);
2165 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2170 computeBhpAtThpLimitInj(frates,
2172 this->getRefDensity(),
2181 auto fratesIter = [
this, &simulator, &wgHelper, &deferred_logger](
const Scalar bhp) {
2185 std::vector<Scalar> rates(3);
2186 computeWellRatesWithBhpIterations(simulator, bhp, wgHelper, rates, deferred_logger);
2191 computeBhpAtThpLimitInj(fratesIter,
2193 this->getRefDensity(),
2204 template<
typename TypeTag>
2209 Scalar max_pressure = 0.0;
2210 const int nseg = this->numberOfSegments();
2211 for (
int seg = 0; seg < nseg; ++seg) {
2212 for (
const int perf : this->segments_.perforations()[seg]) {
2213 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
2214 if (local_perf_index < 0)
2217 const int cell_idx = this->well_cells_[local_perf_index];
2218 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2219 const auto& fs = int_quants.fluidState();
2220 Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2221 max_pressure = std::max(max_pressure, pressure_cell);
2224 max_pressure = this->parallel_well_info_.communication().max(max_pressure);
2225 return max_pressure;
2232 template<
typename TypeTag>
2233 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
2239 std::vector<Scalar> well_q_s(this->num_conservation_quantities_, 0.0);
2240 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
2241 const int nseg = this->numberOfSegments();
2242 for (
int seg = 0; seg < nseg; ++seg) {
2244 const Scalar seg_pressure = getValue(this->primary_variables_.getSegmentPressure(seg));
2245 for (
const int perf : this->segments_.perforations()[seg]) {
2246 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
2247 if (local_perf_index < 0)
2250 const int cell_idx = this->well_cells_[local_perf_index];
2251 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2252 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
2253 getMobility(simulator, local_perf_index, mob, deferred_logger);
2254 const Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quants, cell_idx);
2255 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
2256 const std::vector<Scalar> Tw = this->wellIndex(local_perf_index, int_quants, trans_mult, wellstate_nupcol);
2257 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.0);
2260 computePerfRate(int_quants, mob, Tw, seg, perf, seg_pressure,
2261 allow_cf, cq_s, perf_press, perf_rates, deferred_logger);
2262 for (
int comp = 0; comp < this->num_conservation_quantities_; ++comp) {
2263 well_q_s[comp] += cq_s[comp];
2267 const auto& comm = this->parallel_well_info_.communication();
2268 if (comm.size() > 1)
2270 comm.sum(well_q_s.data(), well_q_s.size());
2276 template <
typename TypeTag>
2277 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
2281 const int num_seg = this->numberOfSegments();
2282 constexpr int num_eq = MSWEval::numWellEq;
2283 std::vector<Scalar> retval(num_seg * num_eq);
2284 for (
int ii = 0; ii < num_seg; ++ii) {
2285 const auto& pv = this->primary_variables_.value(ii);
2286 std::copy(pv.begin(), pv.end(), retval.begin() + ii * num_eq);
2294 template <
typename TypeTag>
2299 const int num_seg = this->numberOfSegments();
2300 constexpr int num_eq = MSWEval::numWellEq;
2301 std::array<Scalar, num_eq> tmp;
2302 for (
int ii = 0; ii < num_seg; ++ii) {
2303 const auto start = it + ii * num_eq;
2304 std::copy(start, start + num_eq, tmp.begin());
2305 this->primary_variables_.setValue(ii, tmp);
2307 return num_seg * num_eq;
2310 template <
typename TypeTag>
2315 Scalar fsTemperature = 0.0;
2316 using SaltConcType =
typename std::decay<
decltype(std::declval<
decltype(simulator.model().intensiveQuantities(0, 0).fluidState())>().saltConcentration())>::type;
2317 SaltConcType fsSaltConcentration{};
2320 if (this->well_cells_.size() > 0) {
2323 const int cell_idx = this->well_cells_[0];
2324 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
2325 const auto& fs = intQuants.fluidState();
2327 fsTemperature = fs.temperature(FluidSystem::oilPhaseIdx).value();
2328 fsSaltConcentration = fs.saltConcentration();
2331 auto info = std::make_tuple(fsTemperature, fsSaltConcentration);
2335 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:41
Class handling assemble of the equation system for MultisegmentWell.
Definition: MultisegmentWellAssemble.hpp:44
typename PrimaryVariables::EvalWell EvalWell
Definition: MultisegmentWellEval.hpp:66
PrimaryVariables primary_variables_
The primary variables.
Definition: MultisegmentWellEval.hpp:143
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
std::optional< Scalar > computeBhpAtThpLimitProd(const WellStateType &well_state, const Simulator &ebos_simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2075
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, WellStateType &well_state, DeferredLogger &deferred_logger, const Scalar relaxation_factor=1.0)
Definition: MultisegmentWell_impl.hpp:716
void updateWaterThroughput(const double dt, WellStateType &well_state) const override
Definition: MultisegmentWell_impl.hpp:2042
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:888
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: MultisegmentWell_impl.hpp:853
void addWellContributions(SparseMatrixAdapter &jacobian) const override
Definition: MultisegmentWell_impl.hpp:875
bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1496
void solveEqAndUpdateWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:611
Scalar getRefDensity() const override
Definition: MultisegmentWell_impl.hpp:1201
std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger, bool iterate_if_no_solution) const override
Definition: MultisegmentWell_impl.hpp:2095
std::vector< Scalar > getPrimaryVars() const override
Definition: MultisegmentWell_impl.hpp:2279
bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger, const bool fixed_control=false, const bool fixed_status=false) override
Definition: MultisegmentWell_impl.hpp:1621
void computeWellRatesWithBhpIterations(const Simulator &simulator, const Scalar &bhp, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:407
void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) override
computing the well potentials for group control
Definition: MultisegmentWell_impl.hpp:297
std::vector< Scalar > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:2235
void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: MultisegmentWell_impl.hpp:230
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:1209
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:1161
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:1057
void updateIPRImplicit(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1374
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:1980
void computeSegmentFluidProperties(const Simulator &simulator, DeferredLogger &deferred_logger)
Definition: MultisegmentWell_impl.hpp:1130
bool computeWellPotentialsImplicit(const Simulator &simulator, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:536
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: MultisegmentWell_impl.hpp:2297
void computeWellRatesWithBhp(const Simulator &simulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:360
void scaleSegmentRatesAndPressure(WellStateType &well_state) const override
updating the segment pressure and rates based the current bhp and well rates
Definition: MultisegmentWell_impl.hpp:173
bool allDrawDownWrongDirection(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:1989
int debug_cost_counter_
Definition: MultisegmentWell.hpp:182
ConvergenceReport getWellConvergence(const Simulator &simulator, const WellStateType &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:204
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:771
static constexpr bool has_polymer
Definition: WellInterface.hpp:115
static constexpr bool has_solvent
Definition: WellInterface.hpp:113
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &ebos_simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2152
void computePerfCellPressDiffs(const Simulator &simulator)
Definition: MultisegmentWell_impl.hpp:640
void computeWellRatesAtBhpLimit(const Simulator &simulator, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:343
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:265
std::vector< Scalar > computeWellPotentialWithTHP(const WellStateType &well_state, const Simulator &simulator, const WellGroupHelperType &wgHelper, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:484
void updateIPR(const Simulator &ebos_simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:1275
FSInfo getFirstPerforationFluidStateInfo(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2313
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 assembleWellEqWithoutIteration(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1799
Scalar maxPerfPress(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2207
void computeInitialSegmentFluids(const Simulator &simulator, DeferredLogger &deferred_logger)
Definition: MultisegmentWell_impl.hpp:697
void checkOperabilityUnderTHPLimit(const Simulator &ebos_simulator, const WellStateType &well_state, const WellGroupHelperType &wgHelper, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1443
void updatePrimaryVariables(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:157
EvalWell getSegmentSurfaceVolume(const Simulator &simulator, const int seg_idx, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2053
void calculateExplicitQuantities(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:755
void updateWellStateWithTarget(const Simulator &simulator, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger) const override
updating the well state based the current control mode
Definition: MultisegmentWell_impl.hpp:184
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.
Definition: WellGroupHelper.hpp:51
const GroupState< Scalar > & groupState() const
Definition: WellGroupHelper.hpp:175
const WellState< Scalar, IndexTraits > & wellState() const
Definition: WellGroupHelper.hpp:256
WellStateGuard pushWellState(WellState< Scalar, IndexTraits > &well_state)
Definition: WellGroupHelper.hpp:190
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
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
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:2049
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:116
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:86
bool solveWellWithOperabilityCheck(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:580
bool thp_update_iterations
Definition: WellInterface.hpp:381
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