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>
160 const auto& well_state = groupStateHelper.
wellState();
161 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(groupStateHelper, 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, groupStateHelper, well_state, deferred_logger);
192 this->scaleSegmentRatesWithWellRates(this->segments_.inlets(),
193 this->segments_.perforations(),
195 this->scaleSegmentPressuresWithBhp(well_state);
201 template <
typename TypeTag>
205 const std::vector<Scalar>& B_avg,
207 const bool relax_tolerance)
const
209 const auto& well_state = groupStateHelper.
wellState();
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>
271 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) {
277 this->linSys_.recoverSolutionWell(x, xw);
279 updateWellState(simulator, xw, groupStateHelper, well_state, deferred_logger);
281 catch (
const NumericalProblem& exp) {
285 deferred_logger.
problem(
"In MultisegmentWell::recoverWellSolutionAndUpdateWellState for well "
286 + this->name() +
": "+exp.what());
295 template <
typename TypeTag>
301 std::vector<Scalar>& well_potentials,
304 const auto [compute_potential, bhp_controlled_well] =
307 if (!compute_potential) {
311 debug_cost_counter_ = 0;
312 bool converged_implicit =
false;
313 if (this->param_.local_well_solver_control_switching_) {
314 converged_implicit = computeWellPotentialsImplicit(simulator, groupStateHelper, well_potentials, deferred_logger);
315 if (!converged_implicit) {
316 deferred_logger.
debug(
"Implicit potential calculations failed for well "
317 + this->name() +
", reverting to original aproach.");
320 if (!converged_implicit) {
322 const auto& summaryState = simulator.vanguard().summaryState();
323 if (!Base::wellHasTHPConstraints(summaryState) || bhp_controlled_well) {
324 computeWellRatesAtBhpLimit(simulator, groupStateHelper, well_potentials, deferred_logger);
326 well_potentials = computeWellPotentialWithTHP(
327 well_state, simulator, groupStateHelper, deferred_logger);
330 deferred_logger.
debug(
"Cost in iterations of finding well potential for well "
333 this->checkNegativeWellPotentials(well_potentials,
334 this->param_.check_well_operability_,
341 template<
typename TypeTag>
346 std::vector<Scalar>& well_flux,
349 if (this->well_ecl_.isInjector()) {
350 const auto controls = this->well_ecl_.injectionControls(simulator.vanguard().summaryState());
351 computeWellRatesWithBhpIterations(simulator, controls.bhp_limit, groupStateHelper, well_flux, deferred_logger);
353 const auto controls = this->well_ecl_.productionControls(simulator.vanguard().summaryState());
354 computeWellRatesWithBhpIterations(simulator, controls.bhp_limit, groupStateHelper, well_flux, deferred_logger);
358 template<
typename TypeTag>
363 std::vector<Scalar>& well_flux,
366 const int np = this->number_of_phases_;
368 well_flux.resize(np, 0.0);
369 const bool allow_cf = this->getAllowCrossFlow();
370 const int nseg = this->numberOfSegments();
371 const WellStateType& well_state = simulator.problem().wellModel().wellState();
372 const auto& ws = well_state.
well(this->indexOfWell());
373 auto segments_copy = ws.segments;
374 segments_copy.scale_pressure(bhp);
375 const auto& segment_pressure = segments_copy.pressure;
376 for (
int seg = 0; seg < nseg; ++seg) {
377 for (
const int perf : this->segments_.perforations()[seg]) {
378 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
379 if (local_perf_index < 0)
381 const int cell_idx = this->well_cells_[local_perf_index];
382 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
384 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.);
385 getMobility(simulator, local_perf_index, mob, deferred_logger);
387 getTransMult(trans_mult, simulator, cell_idx);
388 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
389 std::vector<Scalar> Tw(this->num_conservation_quantities_,
390 this->well_index_[local_perf_index] * trans_mult);
391 this->getTw(Tw, local_perf_index, intQuants, trans_mult, wellstate_nupcol);
392 const Scalar seg_pressure = segment_pressure[seg];
393 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.);
396 computePerfRate(intQuants, mob, Tw, seg, perf, seg_pressure,
397 allow_cf, cq_s, perf_press, perf_rates, deferred_logger);
399 for(
int p = 0; p < np; ++p) {
400 well_flux[FluidSystem::activeCompToActivePhaseIdx(p)] += cq_s[p];
404 this->parallel_well_info_.communication().sum(well_flux.data(), well_flux.size());
408 template<
typename TypeTag>
414 std::vector<Scalar>& well_flux,
428 auto guard = groupStateHelper_copy.
pushWellState(well_state_copy);
429 auto& ws = well_state_copy.
well(this->index_of_well_);
432 const auto& summary_state = simulator.vanguard().summaryState();
433 auto inj_controls = well_copy.
well_ecl_.isInjector()
434 ? well_copy.
well_ecl_.injectionControls(summary_state)
435 : Well::InjectionControls(0);
436 auto prod_controls = well_copy.
well_ecl_.isProducer()
437 ? well_copy.
well_ecl_.productionControls(summary_state) :
438 Well::ProductionControls(0);
442 inj_controls.bhp_limit = bhp;
443 ws.injection_cmode = Well::InjectorCMode::BHP;
445 prod_controls.bhp_limit = bhp;
446 ws.production_cmode = Well::ProducerCMode::BHP;
452 const int np = this->number_of_phases_;
454 for (
int phase = 0; phase < np; ++phase){
455 trivial = trivial && (ws.well_potentials[phase] == 0.0) ;
459 for (
int phase = 0; phase < np; ++phase) {
460 ws.surface_rates[phase] = sign * ws.well_potentials[phase];
464 this->segments_.perforations(),
468 const double dt = simulator.timeStepSize();
471 well_state_copy, deferred_logger);
475 well_flux.resize(np, 0.0);
476 for (
int compIdx = 0; compIdx < this->num_conservation_quantities_; ++compIdx) {
478 well_flux[FluidSystem::activeCompToActivePhaseIdx(compIdx)] = rate.value();
485 template<
typename TypeTag>
486 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
493 std::vector<Scalar> potentials(this->number_of_phases_, 0.0);
494 const auto& summary_state = simulator.vanguard().summaryState();
496 const auto& well = this->well_ecl_;
497 if (well.isInjector()) {
498 auto bhp_at_thp_limit = computeBhpAtThpLimitInj(simulator, groupStateHelper, summary_state, deferred_logger);
499 if (bhp_at_thp_limit) {
500 const auto& controls = well.injectionControls(summary_state);
501 const Scalar bhp = std::min(*bhp_at_thp_limit,
502 static_cast<Scalar>(controls.bhp_limit));
503 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, potentials, deferred_logger);
504 deferred_logger.
debug(
"Converged thp based potential calculation for well "
507 deferred_logger.
warning(
"FAILURE_GETTING_CONVERGED_POTENTIAL",
508 "Failed in getting converged thp based potential calculation for well "
509 + this->name() +
". Instead the bhp based value is used");
510 const auto& controls = well.injectionControls(summary_state);
511 const Scalar bhp = controls.bhp_limit;
512 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, potentials, deferred_logger);
515 auto bhp_at_thp_limit = computeBhpAtThpLimitProd(
516 well_state, simulator, groupStateHelper, summary_state, deferred_logger);
517 if (bhp_at_thp_limit) {
518 const auto& controls = well.productionControls(summary_state);
519 const Scalar bhp = std::max(*bhp_at_thp_limit,
520 static_cast<Scalar>(controls.bhp_limit));
521 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, potentials, deferred_logger);
522 deferred_logger.
debug(
"Converged thp based potential calculation for well "
525 deferred_logger.
warning(
"FAILURE_GETTING_CONVERGED_POTENTIAL",
526 "Failed in getting converged thp based potential calculation for well "
527 + this->name() +
". Instead the bhp based value is used");
528 const auto& controls = well.productionControls(summary_state);
529 const Scalar bhp = controls.bhp_limit;
530 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, potentials, deferred_logger);
537 template<
typename TypeTag>
542 std::vector<Scalar>& well_potentials,
554 auto guard = groupStateHelper_copy.
pushWellState(well_state_copy);
555 auto& ws = well_state_copy.
well(this->index_of_well_);
558 const auto& summary_state = simulator.vanguard().summaryState();
559 auto inj_controls = well_copy.
well_ecl_.isInjector()
560 ? well_copy.
well_ecl_.injectionControls(summary_state)
561 : Well::InjectionControls(0);
562 auto prod_controls = well_copy.
well_ecl_.isProducer()
563 ? well_copy.
well_ecl_.productionControls(summary_state)
564 : Well::ProductionControls(0);
572 const int np = this->number_of_phases_;
574 for (
int phase = 0; phase < np; ++phase){
575 trivial = trivial && (ws.well_potentials[phase] == 0.0) ;
579 for (
int phase = 0; phase < np; ++phase) {
580 ws.surface_rates[phase] = sign * ws.well_potentials[phase];
584 this->segments_.perforations(),
588 const double dt = simulator.timeStepSize();
590 bool converged =
false;
591 if (this->well_ecl_.isProducer()) {
593 simulator, dt, inj_controls, prod_controls, groupStateHelper_copy, well_state_copy, deferred_logger
597 simulator, dt, inj_controls, prod_controls, groupStateHelper_copy, well_state_copy, deferred_logger,
605 well_potentials.clear();
606 well_potentials.resize(np, 0.0);
607 for (
int compIdx = 0; compIdx < this->num_conservation_quantities_; ++compIdx) {
609 well_potentials[FluidSystem::activeCompToActivePhaseIdx(compIdx)] = rate.value();
615 template <
typename TypeTag>
623 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
628 const BVectorWell dx_well = this->linSys_.solve();
629 updateWellState(simulator, dx_well, groupStateHelper, well_state, deferred_logger);
631 catch(
const NumericalProblem& exp) {
635 deferred_logger.
problem(
"In MultisegmentWell::solveEqAndUpdateWellState for well "
636 + this->name() +
": "+exp.what());
645 template <
typename TypeTag>
652 for (
int local_perf_index = 0; local_perf_index < this->number_of_local_perforations_; ++local_perf_index) {
655 std::vector<Scalar> kr(this->number_of_phases_, 0.0);
656 std::vector<Scalar> density(this->number_of_phases_, 0.0);
658 const int cell_idx = this->well_cells_[local_perf_index];
659 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
660 const auto& fs = intQuants.fluidState();
664 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
665 const int water_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::waterPhaseIdx);
666 kr[water_pos] = intQuants.relativePermeability(FluidSystem::waterPhaseIdx).value();
667 sum_kr += kr[water_pos];
668 density[water_pos] = fs.density(FluidSystem::waterPhaseIdx).value();
671 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
672 const int oil_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::oilPhaseIdx);
673 kr[oil_pos] = intQuants.relativePermeability(FluidSystem::oilPhaseIdx).value();
674 sum_kr += kr[oil_pos];
675 density[oil_pos] = fs.density(FluidSystem::oilPhaseIdx).value();
678 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
679 const int gas_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
680 kr[gas_pos] = intQuants.relativePermeability(FluidSystem::gasPhaseIdx).value();
681 sum_kr += kr[gas_pos];
682 density[gas_pos] = fs.density(FluidSystem::gasPhaseIdx).value();
685 assert(sum_kr != 0.);
688 Scalar average_density = 0.;
689 for (
int p = 0; p < this->number_of_phases_; ++p) {
690 average_density += kr[p] * density[p];
692 average_density /= sum_kr;
694 this->cell_perforation_pressure_diffs_[local_perf_index] = this->gravity_ * average_density * this->cell_perforation_depth_diffs_[local_perf_index];
702 template <
typename TypeTag>
708 for (
int seg = 0; seg < this->numberOfSegments(); ++seg) {
710 const Scalar surface_volume = getSegmentSurfaceVolume(simulator, seg, deferred_logger).value();
711 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
712 segment_fluid_initial_[seg][comp_idx] = surface_volume * this->primary_variables_.surfaceVolumeFraction(seg, comp_idx).value();
721 template <
typename TypeTag>
725 const BVectorWell& dwells,
729 const Scalar relaxation_factor)
731 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
733 const Scalar dFLimit = this->param_.dwell_fraction_max_;
734 const Scalar max_pressure_change = this->param_.max_pressure_change_ms_wells_;
735 const bool stop_or_zero_rate_target =
736 this->stoppedOrZeroRateTarget(groupStateHelper, deferred_logger);
737 this->primary_variables_.updateNewton(dwells,
740 stop_or_zero_rate_target,
741 max_pressure_change);
743 const auto& summary_state = simulator.vanguard().summaryState();
744 this->primary_variables_.copyToWellState(*
this, getRefDensity(),
750 auto& ws = well_state.
well(this->index_of_well_);
751 this->segments_.copyPhaseDensities(ws.segments);
754 Base::calculateReservoirRates(simulator.vanguard().eclState().runspec().co2Storage(), well_state.
well(this->index_of_well_));
761 template <
typename TypeTag>
768 updatePrimaryVariables(groupStateHelper, deferred_logger);
769 computePerfCellPressDiffs(simulator);
770 computeInitialSegmentFluids(simulator, deferred_logger);
777 template<
typename TypeTag>
785 auto fluidState = [&simulator,
this](
const int local_perf_index)
787 const auto cell_idx = this->well_cells_[local_perf_index];
788 return simulator.model()
789 .intensiveQuantities(cell_idx, 0).fluidState();
792 const int np = this->number_of_phases_;
793 auto setToZero = [np](
Scalar* x) ->
void
795 std::fill_n(x, np, 0.0);
798 auto addVector = [np](
const Scalar* src,
Scalar* dest) ->
void
800 std::transform(src, src + np, dest, dest, std::plus<>{});
803 auto& ws = well_state.
well(this->index_of_well_);
804 auto& perf_data = ws.perf_data;
805 auto* connPI = perf_data.prod_index.data();
806 auto* wellPI = ws.productivity_index.data();
810 const auto preferred_phase = this->well_ecl_.getPreferredPhase();
811 auto subsetPerfID = 0;
813 for (
const auto& perf : *this->perf_data_){
814 auto allPerfID = perf.ecl_index;
816 auto connPICalc = [&wellPICalc, allPerfID](
const Scalar mobility) ->
Scalar
821 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
826 getMobility(simulator,
static_cast<int>(subsetPerfID), mob, deferred_logger);
828 const auto& fs = fluidState(subsetPerfID);
831 if (this->isInjector()) {
832 this->computeConnLevelInjInd(fs, preferred_phase, connPICalc,
833 mob, connPI, deferred_logger);
836 this->computeConnLevelProdInd(fs, connPICalc, mob, connPI);
839 addVector(connPI, wellPI);
846 const auto& comm = this->parallel_well_info_.communication();
847 if (comm.size() > 1) {
848 comm.sum(wellPI, np);
851 assert (
static_cast<int>(subsetPerfID) == this->number_of_local_perforations_ &&
852 "Internal logic error in processing connections for PI/II");
859 template<
typename TypeTag>
863 [[maybe_unused]]
const int openConnIdx)
const
868 const auto segNum = this->wellEcl()
869 .getConnections()[globalConnIdx].segment();
871 const auto segIdx = this->wellEcl()
872 .getSegments().segmentNumberToIndex(segNum);
874 return this->segments_.density(segIdx).value();
881 template<
typename TypeTag>
886 if (this->number_of_local_perforations_ == 0) {
890 this->linSys_.extract(jacobian);
894 template<
typename TypeTag>
899 const int pressureVarIndex,
900 const bool use_well_weights,
903 if (this->number_of_local_perforations_ == 0) {
908 this->linSys_.extractCPRPressureMatrix(jacobian,
918 template<
typename TypeTag>
919 template<
class Value>
925 const std::vector<Value>& b_perfcells,
926 const std::vector<Value>& mob_perfcells,
927 const std::vector<Value>& Tw,
929 const Value& segment_pressure,
930 const Value& segment_density,
931 const bool& allow_cf,
932 const std::vector<Value>& cmix_s,
933 std::vector<Value>& cq_s,
938 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
939 if (local_perf_index < 0)
943 const Value perf_seg_press_diff = this->gravity() * segment_density *
944 this->segments_.local_perforation_depth_diff(local_perf_index);
946 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[local_perf_index];
950 perf_press = segment_pressure + perf_seg_press_diff;
953 const Value cell_press_at_perf = pressure_cell - cell_perf_press_diff;
956 const Value drawdown = cell_press_at_perf - perf_press;
959 if (drawdown > 0.0) {
961 if (!allow_cf && this->isInjector()) {
966 for (
int comp_idx = 0; comp_idx < this->numConservationQuantities(); ++comp_idx) {
967 const Value cq_p = - Tw[comp_idx] * (mob_perfcells[comp_idx] * drawdown);
968 cq_s[comp_idx] = b_perfcells[comp_idx] * cq_p;
971 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
972 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
973 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
974 const Value cq_s_oil = cq_s[oilCompIdx];
975 const Value cq_s_gas = cq_s[gasCompIdx];
976 cq_s[gasCompIdx] += rs * cq_s_oil;
977 cq_s[oilCompIdx] += rv * cq_s_gas;
981 if (!allow_cf && this->isProducer()) {
986 Value total_mob = mob_perfcells[0];
987 for (
int comp_idx = 1; comp_idx < this->numConservationQuantities(); ++comp_idx) {
988 total_mob += mob_perfcells[comp_idx];
992 Value volume_ratio = 0.0;
993 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
994 const unsigned waterCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
995 volume_ratio += cmix_s[waterCompIdx] / b_perfcells[waterCompIdx];
998 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
999 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
1000 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
1005 const Value d = 1.0 - rv * rs;
1007 if (getValue(d) == 0.0) {
1009 fmt::format(
"Zero d value obtained for well {} "
1010 "during flux calculation with rs {} and rv {}",
1011 this->name(), rs, rv),
1015 const Value tmp_oil = (cmix_s[oilCompIdx] - rv * cmix_s[gasCompIdx]) / d;
1016 volume_ratio += tmp_oil / b_perfcells[oilCompIdx];
1018 const Value tmp_gas = (cmix_s[gasCompIdx] - rs * cmix_s[oilCompIdx]) / d;
1019 volume_ratio += tmp_gas / b_perfcells[gasCompIdx];
1021 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
1022 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
1023 volume_ratio += cmix_s[oilCompIdx] / b_perfcells[oilCompIdx];
1025 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1026 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
1027 volume_ratio += cmix_s[gasCompIdx] / b_perfcells[gasCompIdx];
1031 for (
int componentIdx = 0; componentIdx < this->numConservationQuantities(); ++componentIdx) {
1032 const Value cqt_i = - Tw[componentIdx] * (total_mob * drawdown);
1033 Value cqt_is = cqt_i / volume_ratio;
1034 cq_s[componentIdx] = cmix_s[componentIdx] * cqt_is;
1039 if (this->isProducer()) {
1040 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1041 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
1042 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
1051 const Scalar d = 1.0 - getValue(rv) * getValue(rs);
1054 perf_rates.
vap_oil = getValue(rv) * (getValue(cq_s[gasCompIdx]) - getValue(rs) * getValue(cq_s[oilCompIdx])) / d;
1057 perf_rates.
dis_gas = getValue(rs) * (getValue(cq_s[oilCompIdx]) - getValue(rv) * getValue(cq_s[gasCompIdx])) / d;
1062 template <
typename TypeTag>
1063 template<
class Value>
1067 const std::vector<Value>& mob_perfcells,
1068 const std::vector<Value>& Tw,
1071 const Value& segment_pressure,
1072 const bool& allow_cf,
1073 std::vector<Value>& cq_s,
1079 auto obtain = [
this](
const Eval& value)
1081 if constexpr (std::is_same_v<Value, Scalar>) {
1082 static_cast<void>(
this);
1083 return getValue(value);
1085 return this->extendEval(value);
1088 auto obtainN = [](
const auto& value)
1090 if constexpr (std::is_same_v<Value, Scalar>) {
1091 return getValue(value);
1096 const auto& fs = int_quants.fluidState();
1098 const Value pressure_cell = obtain(this->getPerfCellPressure(fs));
1099 const Value rs = obtain(fs.Rs());
1100 const Value rv = obtain(fs.Rv());
1103 std::vector<Value> b_perfcells(this->num_conservation_quantities_, 0.0);
1105 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
1106 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1110 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
1111 b_perfcells[compIdx] = obtain(fs.invB(phaseIdx));
1114 std::vector<Value> cmix_s(this->numConservationQuantities(), 0.0);
1115 for (
int comp_idx = 0; comp_idx < this->numConservationQuantities(); ++comp_idx) {
1116 cmix_s[comp_idx] = obtainN(this->primary_variables_.surfaceVolumeFraction(seg, comp_idx));
1119 this->computePerfRate(pressure_cell,
1127 obtainN(this->segments_.density(seg)),
1136 template <
typename TypeTag>
1156 auto info = this->getFirstPerforationFluidStateInfo(simulator);
1157 temperature.setValue(std::get<0>(info));
1158 saltConcentration = this->extendEval(std::get<1>(info));
1160 this->segments_.computeFluidProperties(temperature,
1162 this->primary_variables_,
1166 template<
typename TypeTag>
1167 template<
class Value>
1172 const int cell_idx)
const
1174 auto obtain = [
this](
const Eval& value)
1176 if constexpr (std::is_same_v<Value, Scalar>) {
1177 static_cast<void>(
this);
1178 return getValue(value);
1180 return this->extendEval(value);
1186 template <
typename TypeTag>
1187 template<
class Value>
1191 const int local_perf_index,
1192 std::vector<Value>& mob,
1195 auto obtain = [
this](
const Eval& value)
1197 if constexpr (std::is_same_v<Value, Scalar>) {
1198 static_cast<void>(
this);
1199 return getValue(value);
1201 return this->extendEval(value);
1208 const auto perf_ecl_index = this->perforationData()[local_perf_index].ecl_index;
1209 const Connection& con = this->well_ecl_.getConnections()[perf_ecl_index];
1210 const int seg = this->segmentNumberToIndex(con.segment());
1214 const Scalar segment_pres = this->primary_variables_.getSegmentPressure(seg).value();
1215 const Scalar perf_seg_press_diff = this->gravity() * this->segments_.density(seg).value()
1216 * this->segments_.local_perforation_depth_diff(local_perf_index);
1217 const Scalar perf_press = segment_pres + perf_seg_press_diff;
1218 const Scalar multiplier = this->getInjMult(local_perf_index, segment_pres, perf_press, deferred_logger);
1219 for (std::size_t i = 0; i < mob.size(); ++i) {
1220 mob[i] *= multiplier;
1227 template<
typename TypeTag>
1232 return this->segments_.getRefDensity();
1235 template<
typename TypeTag>
1242 const auto& summaryState = simulator.vanguard().summaryState();
1246 const bool bhp_limit_not_defaulted = bhp_limit > 1.5 * unit::barsa;
1247 if ( bhp_limit_not_defaulted || !this->wellHasTHPConstraints(summaryState) ) {
1250 Scalar total_ipr_mass_rate = 0.0;
1251 for (
unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
1253 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1257 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
1258 const Scalar ipr_rate = this->ipr_a_[compIdx] - this->ipr_b_[compIdx] * bhp_limit;
1260 const Scalar rho = FluidSystem::referenceDensity( phaseIdx, Base::pvtRegionIdx() );
1261 total_ipr_mass_rate += ipr_rate * rho;
1263 if ( (this->isProducer() && total_ipr_mass_rate < 0.) || (this->isInjector() && total_ipr_mass_rate > 0.) ) {
1264 this->operability_status_.operable_under_only_bhp_limit =
false;
1268 if (this->operability_status_.operable_under_only_bhp_limit && this->wellHasTHPConstraints(summaryState)) {
1272 std::vector<Scalar> well_rates_bhp_limit;
1273 computeWellRatesWithBhp(simulator, bhp_limit, well_rates_bhp_limit, deferred_logger);
1275 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1278 this->getRefDensity(),
1279 this->wellEcl().alq_value(summaryState),
1282 if ( (this->isProducer() && thp < thp_limit) || (this->isInjector() && thp > thp_limit) ) {
1283 this->operability_status_.obey_thp_limit_under_bhp_limit =
false;
1294 this->operability_status_.operable_under_only_bhp_limit =
true;
1295 this->operability_status_.obey_thp_limit_under_bhp_limit =
false;
1301 template<
typename TypeTag>
1309 std::fill(this->ipr_a_.begin(), this->ipr_a_.end(), 0.);
1310 std::fill(this->ipr_b_.begin(), this->ipr_b_.end(), 0.);
1312 const int nseg = this->numberOfSegments();
1313 std::vector<Scalar> seg_dp(nseg, 0.0);
1314 for (
int seg = 0; seg < nseg; ++seg) {
1316 const Scalar dp = this->getSegmentDp(seg,
1317 this->segments_.density(seg).value(),
1320 for (
const int perf : this->segments_.perforations()[seg]) {
1321 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
1322 if (local_perf_index < 0)
1324 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
1327 getMobility(simulator, local_perf_index, mob, deferred_logger);
1329 const int cell_idx = this->well_cells_[local_perf_index];
1330 const auto& int_quantities = simulator.model().intensiveQuantities(cell_idx, 0);
1331 const auto& fs = int_quantities.fluidState();
1333 const Scalar perf_seg_press_diff = this->segments_.getPressureDiffSegLocalPerf(seg, local_perf_index);
1335 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[local_perf_index];
1336 const Scalar pressure_cell = this->getPerfCellPressure(fs).value();
1339 std::vector<Scalar> b_perf(this->num_conservation_quantities_);
1340 for (std::size_t phase = 0; phase < FluidSystem::numPhases; ++phase) {
1341 if (!FluidSystem::phaseIsActive(phase)) {
1344 const unsigned comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phase));
1345 b_perf[comp_idx] = fs.invB(phase).value();
1349 const Scalar h_perf = cell_perf_press_diff + perf_seg_press_diff + dp;
1350 const Scalar pressure_diff = pressure_cell - h_perf;
1353 if ( (this->isProducer() && pressure_diff < 0.) || (this->isInjector() && pressure_diff > 0.) ) {
1354 deferred_logger.
debug(
"CROSSFLOW_IPR",
1355 "cross flow found when updateIPR for well " + this->name());
1360 getTransMult(trans_mult, simulator, cell_idx);
1361 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1362 std::vector<Scalar> tw_perf(this->num_conservation_quantities_, this->well_index_[perf] * trans_mult);
1363 this->getTw(tw_perf, local_perf_index, int_quantities, trans_mult, wellstate_nupcol);
1364 std::vector<Scalar> ipr_a_perf(this->ipr_a_.size());
1365 std::vector<Scalar> ipr_b_perf(this->ipr_b_.size());
1366 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1367 const Scalar tw_mob = tw_perf[comp_idx] * mob[comp_idx] * b_perf[comp_idx];
1368 ipr_a_perf[comp_idx] += tw_mob * pressure_diff;
1369 ipr_b_perf[comp_idx] += tw_mob;
1373 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1374 const unsigned oil_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
1375 const unsigned gas_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
1376 const Scalar rs = (fs.Rs()).value();
1377 const Scalar rv = (fs.Rv()).value();
1379 const Scalar dis_gas_a = rs * ipr_a_perf[oil_comp_idx];
1380 const Scalar vap_oil_a = rv * ipr_a_perf[gas_comp_idx];
1382 ipr_a_perf[gas_comp_idx] += dis_gas_a;
1383 ipr_a_perf[oil_comp_idx] += vap_oil_a;
1385 const Scalar dis_gas_b = rs * ipr_b_perf[oil_comp_idx];
1386 const Scalar vap_oil_b = rv * ipr_b_perf[gas_comp_idx];
1388 ipr_b_perf[gas_comp_idx] += dis_gas_b;
1389 ipr_b_perf[oil_comp_idx] += vap_oil_b;
1392 for (std::size_t comp_idx = 0; comp_idx < ipr_a_perf.size(); ++comp_idx) {
1393 this->ipr_a_[comp_idx] += ipr_a_perf[comp_idx];
1394 this->ipr_b_[comp_idx] += ipr_b_perf[comp_idx];
1398 this->parallel_well_info_.communication().sum(this->ipr_a_.data(), this->ipr_a_.size());
1399 this->parallel_well_info_.communication().sum(this->ipr_b_.data(), this->ipr_b_.size());
1402 template<
typename TypeTag>
1417 auto rates = well_state.
well(this->index_of_well_).surface_rates;
1419 for (std::size_t p = 0; p < rates.size(); ++p) {
1420 zero_rates &= rates[p] == 0.0;
1422 auto& ws = well_state.
well(this->index_of_well_);
1424 const auto msg = fmt::format(
"updateIPRImplicit: Well {} has zero rate, IPRs might be problematic", this->name());
1425 deferred_logger.
debug(msg);
1438 std::fill(ws.implicit_ipr_a.begin(), ws.implicit_ipr_a.end(), 0.);
1439 std::fill(ws.implicit_ipr_b.begin(), ws.implicit_ipr_b.end(), 0.);
1441 auto inj_controls = Well::InjectionControls(0);
1442 auto prod_controls = Well::ProductionControls(0);
1443 prod_controls.addControl(Well::ProducerCMode::BHP);
1444 prod_controls.bhp_limit = well_state.
well(this->index_of_well_).bhp;
1447 const auto cmode = ws.production_cmode;
1448 ws.production_cmode = Well::ProducerCMode::BHP;
1449 const double dt = simulator.timeStepSize();
1450 assembleWellEqWithoutIteration(simulator, groupStateHelper, dt, inj_controls, prod_controls, well_state, deferred_logger,
1453 BVectorWell rhs(this->numberOfSegments());
1455 rhs[0][SPres] = -1.0;
1457 const BVectorWell x_well = this->linSys_.solve(rhs);
1458 constexpr int num_eq = MSWEval::numWellEq;
1459 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx){
1460 const EvalWell comp_rate = this->primary_variables_.getQs(comp_idx);
1461 const int idx = FluidSystem::activeCompToActivePhaseIdx(comp_idx);
1462 for (
size_t pvIdx = 0; pvIdx < num_eq; ++pvIdx) {
1464 ws.implicit_ipr_b[idx] -= x_well[0][pvIdx]*comp_rate.derivative(pvIdx+Indices::numEq);
1466 ws.implicit_ipr_a[idx] = ws.implicit_ipr_b[idx]*ws.bhp - comp_rate.value();
1469 ws.production_cmode = cmode;
1472 template<
typename TypeTag>
1480 const auto& summaryState = simulator.vanguard().summaryState();
1481 const auto obtain_bhp = this->isProducer()
1482 ? computeBhpAtThpLimitProd(
1483 well_state, simulator, groupStateHelper, summaryState, deferred_logger)
1484 : computeBhpAtThpLimitInj(simulator, groupStateHelper, summaryState, deferred_logger);
1487 this->operability_status_.can_obtain_bhp_with_thp_limit =
true;
1490 this->operability_status_.obey_bhp_limit_with_thp_limit = (*obtain_bhp >= bhp_limit);
1492 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1493 if (this->isProducer() && *obtain_bhp < thp_limit) {
1494 const std::string msg =
" obtained bhp " +
std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1495 +
" bars is SMALLER than thp limit "
1497 +
" bars as a producer for well " + this->name();
1498 deferred_logger.
debug(msg);
1500 else if (this->isInjector() && *obtain_bhp > thp_limit) {
1501 const std::string msg =
" obtained bhp " +
std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1502 +
" bars is LARGER than thp limit "
1504 +
" bars as a injector for well " + this->name();
1505 deferred_logger.
debug(msg);
1510 this->operability_status_.can_obtain_bhp_with_thp_limit =
false;
1511 this->operability_status_.obey_bhp_limit_with_thp_limit =
false;
1512 if (!this->wellIsStopped()) {
1513 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1514 deferred_logger.
debug(
" could not find bhp value at thp limit "
1516 +
" bar for well " + this->name() +
", the well might need to be closed ");
1525 template<
typename TypeTag>
1530 const Well::InjectionControls& inj_controls,
1531 const Well::ProductionControls& prod_controls,
1536 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return true;
1538 const int max_iter_number = this->param_.max_inner_iter_ms_wells_;
1542 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1547 updatePrimaryVariables(groupStateHelper, deferred_logger);
1549 std::vector<std::vector<Scalar> > residual_history;
1550 std::vector<Scalar> measure_history;
1553 Scalar relaxation_factor = 1.;
1554 bool converged =
false;
1555 bool relax_convergence =
false;
1556 this->regularize_ =
false;
1557 for (; it < max_iter_number; ++it, ++debug_cost_counter_) {
1559 if (it > this->param_.strict_inner_iter_wells_) {
1560 relax_convergence =
true;
1561 this->regularize_ =
true;
1564 assembleWellEqWithoutIteration(simulator, groupStateHelper, dt, inj_controls, prod_controls,
1565 well_state, deferred_logger,
1568 const auto report = getWellConvergence(groupStateHelper, Base::B_avg_, deferred_logger, relax_convergence);
1569 if (report.converged()) {
1576 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1580 residual_history.push_back(residuals);
1581 measure_history.push_back(this->getResidualMeasureValue(well_state,
1582 residual_history[it],
1583 this->param_.tolerance_wells_,
1584 this->param_.tolerance_pressure_ms_wells_,
1587 bool min_relaxation_reached = this->update_relaxation_factor(measure_history, relaxation_factor, this->regularize_, deferred_logger);
1588 if (min_relaxation_reached || this->repeatedStagnation(measure_history, this->regularize_, deferred_logger)) {
1590 const auto reportStag = getWellConvergence(groupStateHelper, Base::B_avg_, deferred_logger,
true);
1591 if (reportStag.converged()) {
1593 std::string message = fmt::format(
"Well stagnates/oscillates but {} manages to get converged with relaxed tolerances in {} inner iterations."
1595 deferred_logger.
debug(message);
1602 BVectorWell dx_well;
1604 dx_well = this->linSys_.solve();
1605 updateWellState(simulator, dx_well, groupStateHelper, well_state, deferred_logger, relaxation_factor);
1607 catch(
const NumericalProblem& exp) {
1611 deferred_logger.
problem(
"In MultisegmentWell::iterateWellEqWithControl for well "
1612 + this->name() +
": "+exp.what());
1619 std::ostringstream sstr;
1620 sstr <<
" Well " << this->name() <<
" converged in " << it <<
" inner iterations.";
1621 if (relax_convergence)
1622 sstr <<
" (A relaxed tolerance was used after "<< this->param_.strict_inner_iter_wells_ <<
" iterations)";
1626 deferred_logger.
debug(sstr.str(), OpmLog::defaultDebugVerbosityLevel + (it == 0));
1628 std::ostringstream sstr;
1629 sstr <<
" Well " << this->name() <<
" did not converge in " << it <<
" inner iterations.";
1630#define EXTRA_DEBUG_MSW 0
1632 sstr <<
"***** Outputting the residual history for well " << this->name() <<
" during inner iterations:";
1633 for (
int i = 0; i < it; ++i) {
1634 const auto& residual = residual_history[i];
1635 sstr <<
" residual at " << i <<
"th iteration ";
1636 for (
const auto& res : residual) {
1639 sstr <<
" " << measure_history[i] <<
" \n";
1642#undef EXTRA_DEBUG_MSW
1643 deferred_logger.
debug(sstr.str());
1650 template<
typename TypeTag>
1655 const Well::InjectionControls& inj_controls,
1656 const Well::ProductionControls& prod_controls,
1660 const bool fixed_control ,
1661 const bool fixed_status ,
1662 const bool solving_with_zero_rate )
1664 const int max_iter_number = this->param_.max_inner_iter_ms_wells_;
1668 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1673 updatePrimaryVariables(groupStateHelper, deferred_logger);
1675 std::vector<std::vector<Scalar> > residual_history;
1676 std::vector<Scalar> measure_history;
1679 Scalar relaxation_factor = 1.;
1680 bool converged =
false;
1681 bool relax_convergence =
false;
1682 this->regularize_ =
false;
1683 const auto& summary_state = groupStateHelper.
summaryState();
1688 const int min_its_after_switch = 3;
1690 const int max_status_switch = this->param_.max_well_status_switch_inner_iter_;
1691 int its_since_last_switch = min_its_after_switch;
1692 int switch_count= 0;
1693 int status_switch_count = 0;
1695 const auto well_status_orig = this->wellStatus_;
1696 const auto operability_orig = this->operability_status_;
1697 auto well_status_cur = well_status_orig;
1699 const bool allow_open = well_state.
well(this->index_of_well_).status == WellStatus::OPEN;
1701 const bool allow_switching = !this->wellUnderZeroRateTarget(groupStateHelper, deferred_logger) &&
1702 (!fixed_control || !fixed_status) && allow_open;
1703 bool final_check =
false;
1705 this->operability_status_.resetOperability();
1706 this->operability_status_.solvable =
true;
1708 for (; it < max_iter_number; ++it, ++debug_cost_counter_) {
1709 ++its_since_last_switch;
1710 if (allow_switching && its_since_last_switch >= min_its_after_switch && status_switch_count < max_status_switch){
1711 const Scalar wqTotal = this->primary_variables_.getWQTotal().value();
1712 bool changed = this->updateWellControlAndStatusLocalIteration(
1713 simulator, groupStateHelper, inj_controls, prod_controls, wqTotal,
1714 well_state, deferred_logger, fixed_control, fixed_status,
1715 solving_with_zero_rate
1718 its_since_last_switch = 0;
1720 if (well_status_cur != this->wellStatus_) {
1721 well_status_cur = this->wellStatus_;
1722 status_switch_count++;
1725 if (!changed && final_check) {
1728 final_check =
false;
1730 if (status_switch_count == max_status_switch) {
1731 this->wellStatus_ = well_status_orig;
1735 if (it > this->param_.strict_inner_iter_wells_) {
1736 relax_convergence =
true;
1737 this->regularize_ =
true;
1740 assembleWellEqWithoutIteration(simulator, groupStateHelper, dt, inj_controls, prod_controls,
1741 well_state, deferred_logger, solving_with_zero_rate);
1744 const auto report = getWellConvergence(groupStateHelper, Base::B_avg_, deferred_logger, relax_convergence);
1745 converged = report.converged();
1746 if (this->parallel_well_info_.communication().size() > 1 &&
1747 this->parallel_well_info_.communication().max(converged) != this->parallel_well_info_.communication().min(converged)) {
1748 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()));
1753 if (switch_count > 0 && its_since_last_switch < min_its_after_switch) {
1755 its_since_last_switch = min_its_after_switch;
1763 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1769 residual_history.push_back(residuals);
1773 measure_history.push_back(this->getResidualMeasureValue(well_state,
1774 residual_history[it],
1775 this->param_.tolerance_wells_,
1776 this->param_.tolerance_pressure_ms_wells_,
1778 bool min_relaxation_reached = this->update_relaxation_factor(measure_history, relaxation_factor, this->regularize_, deferred_logger);
1779 if (min_relaxation_reached || this->repeatedStagnation(measure_history, this->regularize_, deferred_logger)) {
1785 const BVectorWell dx_well = this->linSys_.solve();
1786 updateWellState(simulator, dx_well, groupStateHelper, well_state, deferred_logger, relaxation_factor);
1788 catch(
const NumericalProblem& exp) {
1792 deferred_logger.
problem(
"In MultisegmentWell::iterateWellEqWithSwitching for well "
1793 + this->name() +
": "+exp.what());
1799 if (allow_switching){
1801 const bool is_stopped = this->wellIsStopped();
1802 if (this->wellHasTHPConstraints(summary_state)){
1803 this->operability_status_.can_obtain_bhp_with_thp_limit = !is_stopped;
1804 this->operability_status_.obey_thp_limit_under_bhp_limit = !is_stopped;
1806 this->operability_status_.operable_under_only_bhp_limit = !is_stopped;
1809 std::string message = fmt::format(
" Well {} converged in {} inner iterations ("
1810 "{} control/status switches).", this->name(), it, switch_count);
1811 if (relax_convergence) {
1812 message.append(fmt::format(
" (A relaxed tolerance was used after {} iterations)",
1813 this->param_.strict_inner_iter_wells_));
1815 deferred_logger.
debug(message, OpmLog::defaultDebugVerbosityLevel + ((it == 0) && (switch_count == 0)));
1817 this->wellStatus_ = well_status_orig;
1818 this->operability_status_ = operability_orig;
1819 const std::string message = fmt::format(
" Well {} did not converge in {} inner iterations ("
1820 "{} switches, {} status changes).", this->name(), it, switch_count, status_switch_count);
1821 deferred_logger.
debug(message);
1822 this->primary_variables_.outputLowLimitPressureSegments(deferred_logger);
1829 template<
typename TypeTag>
1835 const Well::InjectionControls& inj_controls,
1836 const Well::ProductionControls& prod_controls,
1839 const bool solving_with_zero_rate)
1841 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
return;
1845 this->segments_.updateUpwindingSegments(this->primary_variables_);
1848 computeSegmentFluidProperties(simulator, deferred_logger);
1851 this->linSys_.clear();
1853 auto& ws = well_state.
well(this->index_of_well_);
1854 ws.phase_mixing_rates.fill(0.0);
1861 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
1863 const int nseg = this->numberOfSegments();
1865 const Scalar rhow = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) ?
1866 FluidSystem::referenceDensity( FluidSystem::waterPhaseIdx, Base::pvtRegionIdx() ) : 0.0;
1867 const unsigned watCompIdx = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) ?
1868 FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx) : 0;
1870 for (
int seg = 0; seg < nseg; ++seg) {
1872 const EvalWell seg_pressure = this->primary_variables_.getSegmentPressure(seg);
1873 auto& perf_data = ws.perf_data;
1874 auto& perf_rates = perf_data.phase_rates;
1875 auto& perf_press_state = perf_data.pressure;
1876 for (
const int perf : this->segments_.perforations()[seg]) {
1877 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
1878 if (local_perf_index < 0)
1880 const int cell_idx = this->well_cells_[local_perf_index];
1881 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
1882 std::vector<EvalWell> mob(this->num_conservation_quantities_, 0.0);
1883 getMobility(simulator, local_perf_index, mob, deferred_logger);
1885 getTransMult(trans_mult, simulator, cell_idx);
1886 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1887 std::vector<EvalWell> Tw(this->num_conservation_quantities_, this->well_index_[local_perf_index] * trans_mult);
1888 this->getTw(Tw, local_perf_index, int_quants, trans_mult, wellstate_nupcol);
1889 std::vector<EvalWell> cq_s(this->num_conservation_quantities_, 0.0);
1892 computePerfRate(int_quants, mob, Tw, seg, perf, seg_pressure,
1893 allow_cf, cq_s, perf_press, perfRates, deferred_logger);
1896 if (this->isProducer()) {
1897 ws.phase_mixing_rates[ws.dissolved_gas] += perfRates.
dis_gas;
1898 ws.phase_mixing_rates[ws.vaporized_oil] += perfRates.
vap_oil;
1899 perf_data.phase_mixing_rates[local_perf_index][ws.dissolved_gas] = perfRates.
dis_gas;
1900 perf_data.phase_mixing_rates[local_perf_index][ws.vaporized_oil] = perfRates.
vap_oil;
1904 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1905 perf_rates[local_perf_index*this->number_of_phases_ + FluidSystem::activeCompToActivePhaseIdx(comp_idx)] = cq_s[comp_idx].value();
1907 perf_press_state[local_perf_index] = perf_press.value();
1910 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
1911 perf_data.wat_mass_rates[local_perf_index] = cq_s[watCompIdx].value() * rhow;
1914 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1916 const EvalWell cq_s_effective = cq_s[comp_idx] * this->well_efficiency_factor_;
1918 this->connectionRates_[local_perf_index][comp_idx] = Base::restrictEval(cq_s_effective);
1921 assemblePerforationEq(seg, local_perf_index, comp_idx, cq_s_effective, this->linSys_);
1927 const auto& comm = this->parallel_well_info_.communication();
1928 comm.sum(ws.phase_mixing_rates.data(), ws.phase_mixing_rates.size());
1931 if (this->parallel_well_info_.communication().size() > 1) {
1933 this->linSys_.sumDistributed(this->parallel_well_info_.communication());
1935 for (
int seg = 0; seg < nseg; ++seg) {
1939 const EvalWell segment_surface_volume = getSegmentSurfaceVolume(simulator, seg, deferred_logger);
1944 const Scalar regularization_factor = this->regularize_? this->param_.regularization_factor_wells_ : 1.0;
1946 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1947 const EvalWell accumulation_term = regularization_factor * (segment_surface_volume * this->primary_variables_.surfaceVolumeFraction(seg, comp_idx)
1948 - segment_fluid_initial_[seg][comp_idx]) / dt;
1950 assembleAccumulationTerm(seg, comp_idx, accumulation_term, this->linSys_);
1955 const int seg_upwind = this->segments_.upwinding_segment(seg);
1956 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1958 this->primary_variables_.getSegmentRateUpwinding(seg,
1961 this->well_efficiency_factor_;
1963 assembleOutflowTerm(seg, seg_upwind, comp_idx, segment_rate, this->linSys_);
1969 for (
const int inlet : this->segments_.inlets()[seg]) {
1970 const int inlet_upwind = this->segments_.upwinding_segment(inlet);
1971 for (
int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1973 this->primary_variables_.getSegmentRateUpwinding(inlet,
1976 this->well_efficiency_factor_;
1978 assembleInflowTerm(seg, inlet, inlet_upwind, comp_idx, inlet_rate, this->linSys_);
1985 const bool stopped_or_zero_target = this->stoppedOrZeroRateTarget(groupStateHelper, deferred_logger);
1993 auto& group_state = solving_with_zero_rate
1997 auto group_guard = groupStateHelper_copy.
pushGroupState(group_state);
1999 assembleControlEq(groupStateHelper_copy,
2002 this->getRefDensity(),
2003 this->primary_variables_,
2005 stopped_or_zero_target,
2008 const UnitSystem& unit_system = simulator.vanguard().eclState().getDeckUnitSystem();
2009 const auto& summary_state = simulator.vanguard().summaryState();
2010 this->assemblePressureEq(seg, unit_system, well_state, summary_state, this->param_.use_average_density_ms_wells_, deferred_logger);
2014 this->parallel_well_info_.communication().sum(this->ipr_a_.data(), this->ipr_a_.size());
2015 this->linSys_.createSolver();
2021 template<
typename TypeTag>
2026 return !this->getAllowCrossFlow() && allDrawDownWrongDirection(simulator);
2030 template<
typename TypeTag>
2035 bool all_drawdown_wrong_direction =
true;
2036 const int nseg = this->numberOfSegments();
2038 for (
int seg = 0; seg < nseg; ++seg) {
2039 const EvalWell segment_pressure = this->primary_variables_.getSegmentPressure(seg);
2040 for (
const int perf : this->segments_.perforations()[seg]) {
2041 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
2042 if (local_perf_index < 0)
2045 const int cell_idx = this->well_cells_[local_perf_index];
2046 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
2047 const auto& fs = intQuants.fluidState();
2050 const EvalWell perf_seg_press_diff = this->segments_.getPressureDiffSegLocalPerf(seg, local_perf_index);
2052 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[local_perf_index];
2054 const Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2055 const Scalar perf_press = pressure_cell - cell_perf_press_diff;
2058 const EvalWell drawdown = perf_press - (segment_pressure + perf_seg_press_diff);
2063 if ( (drawdown < 0. && this->isInjector()) ||
2064 (drawdown > 0. && this->isProducer()) ) {
2065 all_drawdown_wrong_direction =
false;
2070 const auto& comm = this->parallel_well_info_.communication();
2071 if (comm.size() > 1)
2073 all_drawdown_wrong_direction =
2074 (comm.min(all_drawdown_wrong_direction ? 1 : 0) == 1);
2077 return all_drawdown_wrong_direction;
2083 template<
typename TypeTag>
2094 template<
typename TypeTag>
2104 auto info = this->getFirstPerforationFluidStateInfo(simulator);
2105 temperature.setValue(std::get<0>(info));
2106 saltConcentration = this->extendEval(std::get<1>(info));
2108 return this->segments_.getSurfaceVolume(temperature,
2110 this->primary_variables_,
2116 template<
typename TypeTag>
2117 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2122 const SummaryState& summary_state,
2129 this->getALQ(well_state),
2136 template<
typename TypeTag>
2137 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2141 const SummaryState& summary_state,
2144 bool iterate_if_no_solution)
const
2148 auto frates = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2154 std::vector<Scalar> rates(3);
2155 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2160 computeBhpAtThpLimitProd(frates,
2162 maxPerfPress(simulator),
2163 this->getRefDensity(),
2165 this->getTHPConstraint(summary_state),
2171 if (!iterate_if_no_solution)
2172 return std::nullopt;
2174 auto fratesIter = [
this, &simulator, &groupStateHelper, &deferred_logger](
const Scalar bhp) {
2178 std::vector<Scalar> rates(3);
2179 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, rates, deferred_logger);
2184 computeBhpAtThpLimitProd(fratesIter,
2186 maxPerfPress(simulator),
2187 this->getRefDensity(),
2189 this->getTHPConstraint(summary_state),
2193 template<
typename TypeTag>
2194 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2198 const SummaryState& summary_state,
2202 auto frates = [
this, &simulator, &deferred_logger](
const Scalar bhp) {
2208 std::vector<Scalar> rates(3);
2209 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2214 computeBhpAtThpLimitInj(frates,
2216 this->getRefDensity(),
2225 auto fratesIter = [
this, &simulator, &groupStateHelper, &deferred_logger](
const Scalar bhp) {
2229 std::vector<Scalar> rates(3);
2230 computeWellRatesWithBhpIterations(simulator, bhp, groupStateHelper, rates, deferred_logger);
2235 computeBhpAtThpLimitInj(fratesIter,
2237 this->getRefDensity(),
2248 template<
typename TypeTag>
2253 Scalar max_pressure = 0.0;
2254 const int nseg = this->numberOfSegments();
2255 for (
int seg = 0; seg < nseg; ++seg) {
2256 for (
const int perf : this->segments_.perforations()[seg]) {
2257 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
2258 if (local_perf_index < 0)
2261 const int cell_idx = this->well_cells_[local_perf_index];
2262 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2263 const auto& fs = int_quants.fluidState();
2264 Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2265 max_pressure = std::max(max_pressure, pressure_cell);
2268 max_pressure = this->parallel_well_info_.communication().max(max_pressure);
2269 return max_pressure;
2276 template<
typename TypeTag>
2277 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
2283 std::vector<Scalar> well_q_s(this->num_conservation_quantities_, 0.0);
2284 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
2285 const int nseg = this->numberOfSegments();
2286 for (
int seg = 0; seg < nseg; ++seg) {
2288 const Scalar seg_pressure = getValue(this->primary_variables_.getSegmentPressure(seg));
2289 for (
const int perf : this->segments_.perforations()[seg]) {
2290 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
2291 if (local_perf_index < 0)
2294 const int cell_idx = this->well_cells_[local_perf_index];
2295 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, 0);
2296 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
2297 getMobility(simulator, local_perf_index, mob, deferred_logger);
2299 getTransMult(trans_mult, simulator, cell_idx);
2300 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
2301 std::vector<Scalar> Tw(this->num_conservation_quantities_, this->well_index_[local_perf_index] * trans_mult);
2302 this->getTw(Tw, local_perf_index, int_quants, trans_mult, wellstate_nupcol);
2303 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.0);
2306 computePerfRate(int_quants, mob, Tw, seg, perf, seg_pressure,
2307 allow_cf, cq_s, perf_press, perf_rates, deferred_logger);
2308 for (
int comp = 0; comp < this->num_conservation_quantities_; ++comp) {
2309 well_q_s[comp] += cq_s[comp];
2313 const auto& comm = this->parallel_well_info_.communication();
2314 if (comm.size() > 1)
2316 comm.sum(well_q_s.data(), well_q_s.size());
2322 template <
typename TypeTag>
2323 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
2327 const int num_seg = this->numberOfSegments();
2328 constexpr int num_eq = MSWEval::numWellEq;
2329 std::vector<Scalar> retval(num_seg * num_eq);
2330 for (
int ii = 0; ii < num_seg; ++ii) {
2331 const auto& pv = this->primary_variables_.value(ii);
2332 std::copy(pv.begin(), pv.end(), retval.begin() + ii * num_eq);
2340 template <
typename TypeTag>
2345 const int num_seg = this->numberOfSegments();
2346 constexpr int num_eq = MSWEval::numWellEq;
2347 std::array<Scalar, num_eq> tmp;
2348 for (
int ii = 0; ii < num_seg; ++ii) {
2349 const auto start = it + ii * num_eq;
2350 std::copy(start, start + num_eq, tmp.begin());
2351 this->primary_variables_.setValue(ii, tmp);
2353 return num_seg * num_eq;
2356 template <
typename TypeTag>
2361 Scalar fsTemperature = 0.0;
2362 using SaltConcType =
typename std::decay<
decltype(std::declval<
decltype(simulator.model().intensiveQuantities(0, 0).fluidState())>().saltConcentration())>::type;
2363 SaltConcType fsSaltConcentration{};
2366 if (this->well_cells_.size() > 0) {
2369 const int cell_idx = this->well_cells_[0];
2370 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, 0);
2371 const auto& fs = intQuants.fluidState();
2373 fsTemperature = fs.temperature(FluidSystem::oilPhaseIdx).value();
2374 fsSaltConcentration = fs.saltConcentration();
2377 auto info = std::make_tuple(fsTemperature, fsSaltConcentration);
2381 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: GroupStateHelper.hpp:53
GroupState< Scalar > & groupState() const
Definition: GroupStateHelper.hpp:180
const SummaryState & summaryState() const
Definition: GroupStateHelper.hpp:257
const WellState< Scalar, IndexTraits > & wellState() const
Definition: GroupStateHelper.hpp:306
WellStateGuard pushWellState(WellState< Scalar, IndexTraits > &well_state)
Definition: GroupStateHelper.hpp:195
GroupStateGuard pushGroupState(GroupState< Scalar > &group_state)
Definition: GroupStateHelper.hpp:190
Definition: GroupState.hpp:41
Class 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
void updateWaterThroughput(const double dt, WellStateType &well_state) const override
Definition: MultisegmentWell_impl.hpp:2086
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:897
void computeWellRatesWithBhpIterations(const Simulator &simulator, const Scalar &bhp, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:411
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: MultisegmentWell_impl.hpp:862
void addWellContributions(SparseMatrixAdapter &jacobian) const override
Definition: MultisegmentWell_impl.hpp:884
void getTransMult(Value &trans_mult, const Simulator &simulator, const int cell_indx) const
Definition: MultisegmentWell_impl.hpp:1170
void updateWellStateWithTarget(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger) const override
updating the well state based the current control mode
Definition: MultisegmentWell_impl.hpp:184
void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) override
computing the well potentials for group control
Definition: MultisegmentWell_impl.hpp:298
void assembleWellEqWithoutIteration(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellStateType &well_state, DeferredLogger &deferred_logger, const bool solving_with_zero_rate) override
Definition: MultisegmentWell_impl.hpp:1832
Scalar getRefDensity() const override
Definition: MultisegmentWell_impl.hpp:1230
std::vector< Scalar > getPrimaryVars() const override
Definition: MultisegmentWell_impl.hpp:2325
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2196
std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger, bool iterate_if_no_solution) const override
Definition: MultisegmentWell_impl.hpp:2139
void solveEqAndUpdateWellState(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:618
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:265
std::vector< Scalar > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:2279
void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: MultisegmentWell_impl.hpp:230
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger, const Scalar relaxation_factor=1.0)
Definition: MultisegmentWell_impl.hpp:724
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:1238
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:1190
bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger, const bool fixed_control, const bool fixed_status, const bool solving_with_zero_rate) override
Definition: MultisegmentWell_impl.hpp:1653
ConvergenceReport getWellConvergence(const GroupStateHelperType &groupStateHelper, const std::vector< Scalar > &B_avg, DeferredLogger &deferred_logger, const bool relax_tolerance) const override
check whether the well equations get converged for this well
Definition: MultisegmentWell_impl.hpp:204
std::vector< Scalar > computeWellPotentialWithTHP(const WellStateType &well_state, const Simulator &simulator, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:488
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2024
void computeSegmentFluidProperties(const Simulator &simulator, DeferredLogger &deferred_logger)
Definition: MultisegmentWell_impl.hpp:1139
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: MultisegmentWell_impl.hpp:2343
void computeWellRatesWithBhp(const Simulator &simulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:361
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:2033
int debug_cost_counter_
Definition: MultisegmentWell.hpp:183
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:780
static constexpr bool has_polymer
Definition: WellInterface.hpp:115
bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1528
Scalar maxPerfPress(const Simulator &simulator) const override
Definition: MultisegmentWell_impl.hpp:2251
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:1066
void calculateExplicitQuantities(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:764
static constexpr bool has_solvent
Definition: WellInterface.hpp:113
void computeWellRatesAtBhpLimit(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:344
bool computeWellPotentialsImplicit(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:540
void computePerfCellPressDiffs(const Simulator &simulator)
Definition: MultisegmentWell_impl.hpp:648
void updateIPRImplicit(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1405
void updateIPR(const Simulator &ebos_simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:1304
void checkOperabilityUnderTHPLimit(const Simulator &ebos_simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1475
std::optional< Scalar > computeBhpAtThpLimitProd(const WellStateType &well_state, const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2119
FSInfo getFirstPerforationFluidStateInfo(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2359
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:705
EvalWell getSegmentSurfaceVolume(const Simulator &simulator, const int seg_idx, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2097
void updatePrimaryVariables(const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:157
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, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:589
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:2060
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:98
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:2073
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:87
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
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:391
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