/* Copyright 2017 SINTEF Digital, Mathematics and Cybernetics. Copyright 2017 Statoil ASA. This file is part of the Open Porous Media project (OPM). OPM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OPM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OPM. If not, see . */ #include #include #include #include #include #include #if HAVE_CUDA || HAVE_OPENCL #include #endif namespace Opm { template MultisegmentWell:: MultisegmentWell(const Well& well, const ParallelWellInfo& pw_info, const int time_step, const ModelParameters& param, const RateConverterType& rate_converter, const int pvtRegionIdx, const int num_components, const int num_phases, const int index_of_well, const std::vector& perf_data) : Base(well, pw_info, time_step, param, rate_converter, pvtRegionIdx, num_components, num_phases, index_of_well, perf_data) , MSWEval(static_cast&>(*this)) , segment_fluid_initial_(this->numberOfSegments(), std::vector(num_components_, 0.0)) { // not handling solvent or polymer for now with multisegment well if constexpr (has_solvent) { OPM_THROW(std::runtime_error, "solvent is not supported by multisegment well yet"); } if constexpr (has_polymer) { OPM_THROW(std::runtime_error, "polymer is not supported by multisegment well yet"); } if constexpr (Base::has_energy) { OPM_THROW(std::runtime_error, "energy is not supported by multisegment well yet"); } if constexpr (Base::has_foam) { OPM_THROW(std::runtime_error, "foam is not supported by multisegment well yet"); } if constexpr (Base::has_brine) { OPM_THROW(std::runtime_error, "brine is not supported by multisegment well yet"); } } template void MultisegmentWell:: init(const PhaseUsage* phase_usage_arg, const std::vector& depth_arg, const double gravity_arg, const int num_cells, const std::vector< Scalar >& B_avg) { Base::init(phase_usage_arg, depth_arg, gravity_arg, num_cells, B_avg); // TODO: for StandardWell, we need to update the perf depth here using depth_arg. // for MultisegmentWell, it is much more complicated. // It can be specified directly, it can be calculated from the segment depth, // it can also use the cell center, which is the same for StandardWell. // For the last case, should we update the depth with the depth_arg? For the // future, it can be a source of wrong result with Multisegment well. // An indicator from the opm-parser should indicate what kind of depth we should use here. // \Note: we do not update the depth here. And it looks like for now, we only have the option to use // specified perforation depth this->initMatrixAndVectors(num_cells); // calcuate the depth difference between the perforations and the perforated grid block for (int perf = 0; perf < number_of_perforations_; ++perf) { const int cell_idx = well_cells_[perf]; this->cell_perforation_depth_diffs_[perf] = depth_arg[cell_idx] - perf_depth_[perf]; } } template void MultisegmentWell:: initPrimaryVariablesEvaluation() const { this->MSWEval::initPrimaryVariablesEvaluation(); } template void MultisegmentWell:: updatePrimaryVariables(const WellState& well_state, DeferredLogger& /* deferred_logger */) const { this->MSWEval::updatePrimaryVariables(well_state); } template void MultisegmentWell:: updateWellStateWithTarget(const Simulator& ebos_simulator, WellState& well_state, DeferredLogger& deferred_logger) const { Base::updateWellStateWithTarget(ebos_simulator, well_state, deferred_logger); // scale segment rates based on the wellRates // and segment pressure based on bhp this->scaleSegmentRatesWithWellRates(well_state); this->scaleSegmentPressuresWithBhp(well_state); } template ConvergenceReport MultisegmentWell:: getWellConvergence(const WellState& well_state, const std::vector& B_avg, DeferredLogger& deferred_logger, const bool relax_tolerance) const { return this->MSWEval::getWellConvergence(well_state, B_avg, deferred_logger, param_.max_residual_allowed_, param_.tolerance_wells_, param_.relaxed_inner_tolerance_flow_ms_well_, param_.tolerance_pressure_ms_wells_, param_.relaxed_inner_tolerance_pressure_ms_well_, relax_tolerance); } template void MultisegmentWell:: apply(const BVector& x, BVector& Ax) const { if (!this->isOperable() && !this->wellIsStopped()) return; if ( param_.matrix_add_well_contributions_ ) { // Contributions are already in the matrix itself return; } BVectorWell Bx(this->duneB_.N()); this->duneB_.mv(x, Bx); // invDBx = duneD^-1 * Bx_ const BVectorWell invDBx = mswellhelpers::applyUMFPack(this->duneD_, this->duneDSolver_, Bx); // Ax = Ax - duneC_^T * invDBx this->duneC_.mmtv(invDBx,Ax); } template void MultisegmentWell:: apply(BVector& r) const { if (!this->isOperable() && !this->wellIsStopped()) return; // invDrw_ = duneD^-1 * resWell_ const BVectorWell invDrw = mswellhelpers::applyUMFPack(this->duneD_, this->duneDSolver_, this->resWell_); // r = r - duneC_^T * invDrw this->duneC_.mmtv(invDrw, r); } template void MultisegmentWell:: recoverWellSolutionAndUpdateWellState(const BVector& x, WellState& well_state, DeferredLogger& deferred_logger) const { if (!this->isOperable() && !this->wellIsStopped()) return; BVectorWell xw(1); this->recoverSolutionWell(x, xw); updateWellState(xw, well_state, deferred_logger); } template void MultisegmentWell:: computeWellPotentials(const Simulator& ebosSimulator, const WellState& well_state, std::vector& well_potentials, DeferredLogger& deferred_logger) { const int np = number_of_phases_; well_potentials.resize(np, 0.0); // Stopped wells have zero potential. if (this->wellIsStopped()) { return; } // If the well is pressure controlled the potential equals the rate. bool pressure_controlled_well = false; if (this->isInjector()) { const Well::InjectorCMode& current = well_state.currentInjectionControl(index_of_well_); if (current == Well::InjectorCMode::BHP || current == Well::InjectorCMode::THP) { pressure_controlled_well = true; } } else { const Well::ProducerCMode& current = well_state.currentProductionControl(index_of_well_); if (current == Well::ProducerCMode::BHP || current == Well::ProducerCMode::THP) { pressure_controlled_well = true; } } if (pressure_controlled_well) { // initialized the well rates with the potentials i.e. the well rates based on bhp const double sign = this->well_ecl_.isInjector() ? 1.0 : -1.0; for (int phase = 0; phase < np; ++phase){ well_potentials[phase] = sign * well_state.wellRates(index_of_well_)[phase]; } return; } debug_cost_counter_ = 0; // does the well have a THP related constraint? const auto& summaryState = ebosSimulator.vanguard().summaryState(); if (!Base::wellHasTHPConstraints(summaryState)) { computeWellRatesAtBhpLimit(ebosSimulator, well_potentials, deferred_logger); } else { well_potentials = computeWellPotentialWithTHP(ebosSimulator, deferred_logger); } deferred_logger.debug("Cost in iterations of finding well potential for well " + name() + ": " + std::to_string(debug_cost_counter_)); } template void MultisegmentWell:: computeWellRatesAtBhpLimit(const Simulator& ebosSimulator, std::vector& well_flux, DeferredLogger& deferred_logger) const { if (well_ecl_.isInjector()) { const auto controls = well_ecl_.injectionControls(ebosSimulator.vanguard().summaryState()); computeWellRatesWithBhp(ebosSimulator, controls.bhp_limit, well_flux, deferred_logger); } else { const auto controls = well_ecl_.productionControls(ebosSimulator.vanguard().summaryState()); computeWellRatesWithBhp(ebosSimulator, controls.bhp_limit, well_flux, deferred_logger); } } template void MultisegmentWell:: computeWellRatesWithBhp(const Simulator& ebosSimulator, const Scalar bhp, std::vector& well_flux, DeferredLogger& deferred_logger) const { // creating a copy of the well itself, to avoid messing up the explicit informations // during this copy, the only information not copied properly is the well controls MultisegmentWell well_copy(*this); well_copy.debug_cost_counter_ = 0; // store a copy of the well state, we don't want to update the real well state WellState well_state_copy = ebosSimulator.problem().wellModel().wellState(); const auto& group_state = ebosSimulator.problem().wellModel().groupState(); // Get the current controls. const auto& summary_state = ebosSimulator.vanguard().summaryState(); auto inj_controls = well_copy.well_ecl_.isInjector() ? well_copy.well_ecl_.injectionControls(summary_state) : Well::InjectionControls(0); auto prod_controls = well_copy.well_ecl_.isProducer() ? well_copy.well_ecl_.productionControls(summary_state) : Well::ProductionControls(0); // Set current control to bhp, and bhp value in state, modify bhp limit in control object. if (well_copy.well_ecl_.isInjector()) { inj_controls.bhp_limit = bhp; well_state_copy.currentInjectionControl(index_of_well_, Well::InjectorCMode::BHP); } else { prod_controls.bhp_limit = bhp; well_state_copy.currentProductionControl(index_of_well_, Well::ProducerCMode::BHP); } well_state_copy.update_bhp(well_copy.index_of_well_, bhp); well_copy.scaleSegmentPressuresWithBhp(well_state_copy); // initialized the well rates with the potentials i.e. the well rates based on bhp const int np = number_of_phases_; const double sign = well_copy.well_ecl_.isInjector() ? 1.0 : -1.0; for (int phase = 0; phase < np; ++phase){ well_state_copy.wellRates(well_copy.index_of_well_)[phase] = sign * well_state_copy.wellPotentials(well_copy.index_of_well_)[phase]; } well_copy.scaleSegmentRatesWithWellRates(well_state_copy); well_copy.calculateExplicitQuantities(ebosSimulator, well_state_copy, deferred_logger); const double dt = ebosSimulator.timeStepSize(); // iterate to get a solution at the given bhp. well_copy.iterateWellEqWithControl(ebosSimulator, dt, inj_controls, prod_controls, well_state_copy, group_state, deferred_logger); // compute the potential and store in the flux vector. well_flux.clear(); well_flux.resize(np, 0.0); for (int compIdx = 0; compIdx < num_components_; ++compIdx) { const EvalWell rate = well_copy.getQs(compIdx); well_flux[ebosCompIdxToFlowCompIdx(compIdx)] = rate.value(); } debug_cost_counter_ += well_copy.debug_cost_counter_; } template std::vector MultisegmentWell:: computeWellPotentialWithTHP(const Simulator& ebos_simulator, DeferredLogger& deferred_logger) const { std::vector potentials(number_of_phases_, 0.0); const auto& summary_state = ebos_simulator.vanguard().summaryState(); const auto& well = well_ecl_; if (well.isInjector()){ auto bhp_at_thp_limit = computeBhpAtThpLimitInj(ebos_simulator, summary_state, deferred_logger); if (bhp_at_thp_limit) { const auto& controls = well_ecl_.injectionControls(summary_state); const double bhp = std::min(*bhp_at_thp_limit, controls.bhp_limit); computeWellRatesWithBhp(ebos_simulator, bhp, potentials, deferred_logger); deferred_logger.debug("Converged thp based potential calculation for well " + name() + ", at bhp = " + std::to_string(bhp)); } else { deferred_logger.warning("FAILURE_GETTING_CONVERGED_POTENTIAL", "Failed in getting converged thp based potential calculation for well " + name() + ". Instead the bhp based value is used"); const auto& controls = well_ecl_.injectionControls(summary_state); const double bhp = controls.bhp_limit; computeWellRatesWithBhp(ebos_simulator, bhp, potentials, deferred_logger); } } else { auto bhp_at_thp_limit = computeBhpAtThpLimitProd(ebos_simulator, summary_state, deferred_logger); if (bhp_at_thp_limit) { const auto& controls = well_ecl_.productionControls(summary_state); const double bhp = std::max(*bhp_at_thp_limit, controls.bhp_limit); computeWellRatesWithBhp(ebos_simulator, bhp, potentials, deferred_logger); deferred_logger.debug("Converged thp based potential calculation for well " + name() + ", at bhp = " + std::to_string(bhp)); } else { deferred_logger.warning("FAILURE_GETTING_CONVERGED_POTENTIAL", "Failed in getting converged thp based potential calculation for well " + name() + ". Instead the bhp based value is used"); const auto& controls = well_ecl_.productionControls(summary_state); const double bhp = controls.bhp_limit; computeWellRatesWithBhp(ebos_simulator, bhp, potentials, deferred_logger); } } return potentials; } template void MultisegmentWell:: solveEqAndUpdateWellState(WellState& well_state, DeferredLogger& deferred_logger) { if (!this->isOperable() && !this->wellIsStopped()) return; // We assemble the well equations, then we check the convergence, // which is why we do not put the assembleWellEq here. const BVectorWell dx_well = mswellhelpers::applyUMFPack(this->duneD_, this->duneDSolver_, this->resWell_); updateWellState(dx_well, well_state, deferred_logger); } template void MultisegmentWell:: computePerfCellPressDiffs(const Simulator& ebosSimulator) { for (int perf = 0; perf < number_of_perforations_; ++perf) { std::vector kr(number_of_phases_, 0.0); std::vector density(number_of_phases_, 0.0); const int cell_idx = well_cells_[perf]; const auto& intQuants = *(ebosSimulator.model().cachedIntensiveQuantities(cell_idx, /*timeIdx=*/ 0)); const auto& fs = intQuants.fluidState(); double sum_kr = 0.; const PhaseUsage& pu = phaseUsage(); if (pu.phase_used[Water]) { const int water_pos = pu.phase_pos[Water]; kr[water_pos] = intQuants.relativePermeability(FluidSystem::waterPhaseIdx).value(); sum_kr += kr[water_pos]; density[water_pos] = fs.density(FluidSystem::waterPhaseIdx).value(); } if (pu.phase_used[Oil]) { const int oil_pos = pu.phase_pos[Oil]; kr[oil_pos] = intQuants.relativePermeability(FluidSystem::oilPhaseIdx).value(); sum_kr += kr[oil_pos]; density[oil_pos] = fs.density(FluidSystem::oilPhaseIdx).value(); } if (pu.phase_used[Gas]) { const int gas_pos = pu.phase_pos[Gas]; kr[gas_pos] = intQuants.relativePermeability(FluidSystem::gasPhaseIdx).value(); sum_kr += kr[gas_pos]; density[gas_pos] = fs.density(FluidSystem::gasPhaseIdx).value(); } assert(sum_kr != 0.); // calculate the average density double average_density = 0.; for (int p = 0; p < number_of_phases_; ++p) { average_density += kr[p] * density[p]; } average_density /= sum_kr; this->cell_perforation_pressure_diffs_[perf] = gravity_ * average_density * this->cell_perforation_depth_diffs_[perf]; } } template void MultisegmentWell:: computeInitialSegmentFluids(const Simulator& ebos_simulator) { for (int seg = 0; seg < this->numberOfSegments(); ++seg) { // TODO: trying to reduce the times for the surfaceVolumeFraction calculation const double surface_volume = getSegmentSurfaceVolume(ebos_simulator, seg).value(); for (int comp_idx = 0; comp_idx < num_components_; ++comp_idx) { segment_fluid_initial_[seg][comp_idx] = surface_volume * this->surfaceVolumeFraction(seg, comp_idx).value(); } } } template void MultisegmentWell:: updateWellState(const BVectorWell& dwells, WellState& well_state, DeferredLogger& deferred_logger, const double relaxation_factor) const { if (!this->isOperable() && !this->wellIsStopped()) return; const double dFLimit = param_.dwell_fraction_max_; const double max_pressure_change = param_.max_pressure_change_ms_wells_; this->MSWEval::updateWellState(dwells, relaxation_factor, dFLimit, max_pressure_change); this->updateWellStateFromPrimaryVariables(well_state, getRefDensity(), deferred_logger); Base::calculateReservoirRates(well_state); } template void MultisegmentWell:: calculateExplicitQuantities(const Simulator& ebosSimulator, const WellState& well_state, DeferredLogger& deferred_logger) { updatePrimaryVariables(well_state, deferred_logger); initPrimaryVariablesEvaluation(); computePerfCellPressDiffs(ebosSimulator); computeInitialSegmentFluids(ebosSimulator); } template void MultisegmentWell:: updateProductivityIndex(const Simulator& ebosSimulator, const WellProdIndexCalculator& wellPICalc, WellState& well_state, DeferredLogger& deferred_logger) const { auto fluidState = [&ebosSimulator, this](const int perf) { const auto cell_idx = this->well_cells_[perf]; return ebosSimulator.model() .cachedIntensiveQuantities(cell_idx, /*timeIdx=*/ 0)->fluidState(); }; const int np = this->number_of_phases_; auto setToZero = [np](double* x) -> void { std::fill_n(x, np, 0.0); }; auto addVector = [np](const double* src, double* dest) -> void { std::transform(src, src + np, dest, dest, std::plus<>{}); }; auto* wellPI = well_state.productivityIndex(this->index_of_well_).data(); auto* connPI = well_state.connectionProductivityIndex(this->index_of_well_).data(); setToZero(wellPI); const auto preferred_phase = this->well_ecl_.getPreferredPhase(); auto subsetPerfID = 0; for ( const auto& perf : *this->perf_data_){ auto allPerfID = perf.ecl_index; auto connPICalc = [&wellPICalc, allPerfID](const double mobility) -> double { return wellPICalc.connectionProdIndStandard(allPerfID, mobility); }; std::vector mob(this->num_components_, 0.0); getMobility(ebosSimulator, static_cast(subsetPerfID), mob); const auto& fs = fluidState(subsetPerfID); setToZero(connPI); if (this->isInjector()) { this->computeConnLevelInjInd(fs, preferred_phase, connPICalc, mob, connPI, deferred_logger); } else { // Production or zero flow rate this->computeConnLevelProdInd(fs, connPICalc, mob, connPI); } addVector(connPI, wellPI); ++subsetPerfID; connPI += np; } assert (static_cast(subsetPerfID) == this->number_of_perforations_ && "Internal logic error in processing connections for PI/II"); } template void MultisegmentWell:: addWellContributions(SparseMatrixAdapter& jacobian) const { const auto invDuneD = mswellhelpers::invertWithUMFPack(this->duneD_, this->duneDSolver_); // We need to change matrix A as follows // A -= C^T D^-1 B // D is a (nseg x nseg) block matrix with (4 x 4) blocks. // B and C are (nseg x ncells) block matrices with (4 x 4 blocks). // They have nonzeros at (i, j) only if this well has a // perforation at cell j connected to segment i. The code // assumes that no cell is connected to more than one segment, // i.e. the columns of B/C have no more than one nonzero. for (size_t rowC = 0; rowC < this->duneC_.N(); ++rowC) { for (auto colC = this->duneC_[rowC].begin(), endC = this->duneC_[rowC].end(); colC != endC; ++colC) { const auto row_index = colC.index(); for (size_t rowB = 0; rowB < this->duneB_.N(); ++rowB) { for (auto colB = this->duneB_[rowB].begin(), endB = this->duneB_[rowB].end(); colB != endB; ++colB) { const auto col_index = colB.index(); OffDiagMatrixBlockWellType tmp1; Detail::multMatrixImpl(invDuneD[rowC][rowB], (*colB), tmp1, std::true_type()); typename SparseMatrixAdapter::MatrixBlock tmp2; Detail::multMatrixTransposedImpl((*colC), tmp1, tmp2, std::false_type()); jacobian.addToBlock(row_index, col_index, tmp2); } } } } } template void MultisegmentWell:: computePerfRatePressure(const IntensiveQuantities& int_quants, const std::vector& mob_perfcells, const double Tw, const int seg, const int perf, const EvalWell& segment_pressure, const bool& allow_cf, std::vector& cq_s, EvalWell& perf_press, double& perf_dis_gas_rate, double& perf_vap_oil_rate, DeferredLogger& deferred_logger) const { const auto& fs = int_quants.fluidState(); const EvalWell pressure_cell = this->extendEval(fs.pressure(FluidSystem::oilPhaseIdx)); const EvalWell rs = this->extendEval(fs.Rs()); const EvalWell rv = this->extendEval(fs.Rv()); // not using number_of_phases_ because of solvent std::vector b_perfcells(num_components_, 0.0); for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) { if (!FluidSystem::phaseIsActive(phaseIdx)) { continue; } const unsigned compIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx)); b_perfcells[compIdx] = this->extendEval(fs.invB(phaseIdx)); } this->MSWEval::computePerfRatePressure(pressure_cell, rs, rv, b_perfcells, mob_perfcells, Tw, seg, perf, segment_pressure, allow_cf, cq_s, perf_press, perf_dis_gas_rate, perf_vap_oil_rate, deferred_logger); } template void MultisegmentWell:: computeSegmentFluidProperties(const Simulator& ebosSimulator) { // TODO: the concept of phases and components are rather confusing in this function. // needs to be addressed sooner or later. // get the temperature for later use. It is only useful when we are not handling // thermal related simulation // basically, it is a single value for all the segments EvalWell temperature; EvalWell saltConcentration; // not sure how to handle the pvt region related to segment // for the current approach, we use the pvt region of the first perforated cell // although there are some text indicating using the pvt region of the lowest // perforated cell // TODO: later to investigate how to handle the pvt region int pvt_region_index; { // using the first perforated cell const int cell_idx = well_cells_[0]; const auto& intQuants = *(ebosSimulator.model().cachedIntensiveQuantities(cell_idx, /*timeIdx=*/0)); const auto& fs = intQuants.fluidState(); temperature.setValue(fs.temperature(FluidSystem::oilPhaseIdx).value()); saltConcentration = this->extendEval(fs.saltConcentration()); pvt_region_index = fs.pvtRegionIndex(); } this->MSWEval::computeSegmentFluidProperties(temperature, saltConcentration, pvt_region_index); } template void MultisegmentWell:: getMobility(const Simulator& ebosSimulator, const int perf, std::vector& mob) const { // TODO: most of this function, if not the whole function, can be moved to the base class const int cell_idx = well_cells_[perf]; assert (int(mob.size()) == num_components_); const auto& intQuants = *(ebosSimulator.model().cachedIntensiveQuantities(cell_idx, /*timeIdx=*/0)); const auto& materialLawManager = ebosSimulator.problem().materialLawManager(); // either use mobility of the perforation cell or calcualte its own // based on passing the saturation table index const int satid = saturation_table_number_[perf] - 1; const int satid_elem = materialLawManager->satnumRegionIdx(cell_idx); if( satid == satid_elem ) { // the same saturation number is used. i.e. just use the mobilty from the cell for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) { if (!FluidSystem::phaseIsActive(phaseIdx)) { continue; } const unsigned activeCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx)); mob[activeCompIdx] = this->extendEval(intQuants.mobility(phaseIdx)); } // if (has_solvent) { // mob[contiSolventEqIdx] = extendEval(intQuants.solventMobility()); // } } else { const auto& paramsCell = materialLawManager->connectionMaterialLawParams(satid, cell_idx); Eval relativePerms[3] = { 0.0, 0.0, 0.0 }; MaterialLaw::relativePermeabilities(relativePerms, paramsCell, intQuants.fluidState()); // reset the satnumvalue back to original materialLawManager->connectionMaterialLawParams(satid_elem, cell_idx); // compute the mobility for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) { if (!FluidSystem::phaseIsActive(phaseIdx)) { continue; } const unsigned activeCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx)); mob[activeCompIdx] = this->extendEval(relativePerms[phaseIdx] / intQuants.fluidState().viscosity(phaseIdx)); } } } template double MultisegmentWell:: getRefDensity() const { return this->segment_densities_[0].value(); } template void MultisegmentWell:: checkOperabilityUnderBHPLimitProducer(const WellState& /*well_state*/, const Simulator& ebos_simulator, DeferredLogger& deferred_logger) { const auto& summaryState = ebos_simulator.vanguard().summaryState(); const double bhp_limit = Base::mostStrictBhpFromBhpLimits(summaryState); // Crude but works: default is one atmosphere. // TODO: a better way to detect whether the BHP is defaulted or not const bool bhp_limit_not_defaulted = bhp_limit > 1.5 * unit::barsa; if ( bhp_limit_not_defaulted || !this->wellHasTHPConstraints(summaryState) ) { // if the BHP limit is not defaulted or the well does not have a THP limit // we need to check the BHP limit double temp = 0; for (int p = 0; p < number_of_phases_; ++p) { temp += ipr_a_[p] - ipr_b_[p] * bhp_limit; } if (temp < 0.) { this->operability_status_.operable_under_only_bhp_limit = false; } // checking whether running under BHP limit will violate THP limit if (this->operability_status_.operable_under_only_bhp_limit && this->wellHasTHPConstraints(summaryState)) { // option 1: calculate well rates based on the BHP limit. // option 2: stick with the above IPR curve // we use IPR here std::vector well_rates_bhp_limit; computeWellRatesWithBhp(ebos_simulator, bhp_limit, well_rates_bhp_limit, deferred_logger); const double thp = this->calculateThpFromBhp(well_rates_bhp_limit, bhp_limit, getRefDensity(), deferred_logger); const double thp_limit = this->getTHPConstraint(summaryState); if (thp < thp_limit) { this->operability_status_.obey_thp_limit_under_bhp_limit = false; } } } else { // defaulted BHP and there is a THP constraint // default BHP limit is about 1 atm. // when applied the hydrostatic pressure correction dp, // most likely we get a negative value (bhp + dp)to search in the VFP table, // which is not desirable. // we assume we can operate under defaulted BHP limit and will violate the THP limit // when operating under defaulted BHP limit. this->operability_status_.operable_under_only_bhp_limit = true; this->operability_status_.obey_thp_limit_under_bhp_limit = false; } } template void MultisegmentWell:: updateIPR(const Simulator& ebos_simulator, DeferredLogger& deferred_logger) const { // TODO: not handling solvent related here for now // TODO: it only handles the producers for now // the formular for the injectors are not formulated yet if (this->isInjector()) { return; } // initialize all the values to be zero to begin with std::fill(ipr_a_.begin(), ipr_a_.end(), 0.); std::fill(ipr_b_.begin(), ipr_b_.end(), 0.); const int nseg = this->numberOfSegments(); double seg_bhp_press_diff = 0; double ref_depth = ref_depth_; for (int seg = 0; seg < nseg; ++seg) { // calculating the perforation rate for each perforation that belongs to this segment const double segment_depth = this->segmentSet()[seg].depth(); const double dp = wellhelpers::computeHydrostaticCorrection(ref_depth, segment_depth, this->segment_densities_[seg].value(), gravity_); ref_depth = segment_depth; seg_bhp_press_diff += dp; for (const int perf : this->segment_perforations_[seg]) { //std::vector mob(num_components_, {numWellEq_ + numEq, 0.0}); std::vector mob(num_components_, 0.0); // TODO: mabye we should store the mobility somewhere, so that we only need to calculate it one per iteration getMobility(ebos_simulator, perf, mob); const int cell_idx = well_cells_[perf]; const auto& int_quantities = *(ebos_simulator.model().cachedIntensiveQuantities(cell_idx, /*timeIdx=*/ 0)); const auto& fs = int_quantities.fluidState(); // the pressure of the reservoir grid block the well connection is in // pressure difference between the segment and the perforation const double perf_seg_press_diff = gravity_ * this->segment_densities_[seg].value() * this->perforation_segment_depth_diffs_[perf]; // pressure difference between the perforation and the grid cell const double cell_perf_press_diff = this->cell_perforation_pressure_diffs_[perf]; const double pressure_cell = fs.pressure(FluidSystem::oilPhaseIdx).value(); // calculating the b for the connection std::vector b_perf(num_components_); for (size_t phase = 0; phase < FluidSystem::numPhases; ++phase) { if (!FluidSystem::phaseIsActive(phase)) { continue; } const unsigned comp_idx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phase)); b_perf[comp_idx] = fs.invB(phase).value(); } // the pressure difference between the connection and BHP const double h_perf = cell_perf_press_diff + perf_seg_press_diff + seg_bhp_press_diff; const double pressure_diff = pressure_cell - h_perf; // Let us add a check, since the pressure is calculated based on zero value BHP // it should not be negative anyway. If it is negative, we might need to re-formulate // to taking into consideration the crossflow here. if (pressure_diff <= 0.) { deferred_logger.warning("NON_POSITIVE_DRAWDOWN_IPR", "non-positive drawdown found when updateIPR for well " + name()); } // the well index associated with the connection const double tw_perf = well_index_[perf]*ebos_simulator.problem().template rockCompTransMultiplier(int_quantities, cell_idx); // TODO: there might be some indices related problems here // phases vs components // ipr values for the perforation std::vector ipr_a_perf(ipr_a_.size()); std::vector ipr_b_perf(ipr_b_.size()); for (int p = 0; p < number_of_phases_; ++p) { const double tw_mob = tw_perf * mob[p].value() * b_perf[p]; ipr_a_perf[p] += tw_mob * pressure_diff; ipr_b_perf[p] += tw_mob; } // we need to handle the rs and rv when both oil and gas are present if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) { const unsigned oil_comp_idx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx); const unsigned gas_comp_idx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx); const double rs = (fs.Rs()).value(); const double rv = (fs.Rv()).value(); const double dis_gas_a = rs * ipr_a_perf[oil_comp_idx]; const double vap_oil_a = rv * ipr_a_perf[gas_comp_idx]; ipr_a_perf[gas_comp_idx] += dis_gas_a; ipr_a_perf[oil_comp_idx] += vap_oil_a; const double dis_gas_b = rs * ipr_b_perf[oil_comp_idx]; const double vap_oil_b = rv * ipr_b_perf[gas_comp_idx]; ipr_b_perf[gas_comp_idx] += dis_gas_b; ipr_b_perf[oil_comp_idx] += vap_oil_b; } for (int p = 0; p < number_of_phases_; ++p) { // TODO: double check the indices here ipr_a_[ebosCompIdxToFlowCompIdx(p)] += ipr_a_perf[p]; ipr_b_[ebosCompIdxToFlowCompIdx(p)] += ipr_b_perf[p]; } } } } template void MultisegmentWell:: checkOperabilityUnderTHPLimitProducer(const Simulator& ebos_simulator, const WellState& /*well_state*/, DeferredLogger& deferred_logger) { const auto& summaryState = ebos_simulator.vanguard().summaryState(); const auto obtain_bhp = computeBhpAtThpLimitProd(ebos_simulator, summaryState, deferred_logger); if (obtain_bhp) { this->operability_status_.can_obtain_bhp_with_thp_limit = true; const double bhp_limit = Base::mostStrictBhpFromBhpLimits(summaryState); this->operability_status_.obey_bhp_limit_with_thp_limit = (*obtain_bhp >= bhp_limit); const double thp_limit = this->getTHPConstraint(summaryState); if (*obtain_bhp < thp_limit) { const std::string msg = " obtained bhp " + std::to_string(unit::convert::to(*obtain_bhp, unit::barsa)) + " bars is SMALLER than thp limit " + std::to_string(unit::convert::to(thp_limit, unit::barsa)) + " bars as a producer for well " + name(); deferred_logger.debug(msg); } } else { // Shutting wells that can not find bhp value from thp // when under THP control this->operability_status_.can_obtain_bhp_with_thp_limit = false; this->operability_status_.obey_bhp_limit_with_thp_limit = false; if (!this->wellIsStopped()) { const double thp_limit = this->getTHPConstraint(summaryState); deferred_logger.debug(" could not find bhp value at thp limit " + std::to_string(unit::convert::to(thp_limit, unit::barsa)) + " bar for well " + name() + ", the well might need to be closed "); } } } template bool MultisegmentWell:: iterateWellEqWithControl(const Simulator& ebosSimulator, const double dt, const Well::InjectionControls& inj_controls, const Well::ProductionControls& prod_controls, WellState& well_state, const GroupState& group_state, DeferredLogger& deferred_logger) { if (!this->isOperable() && !this->wellIsStopped()) return true; const int max_iter_number = param_.max_inner_iter_ms_wells_; const WellState well_state0 = well_state; const std::vector residuals0 = this->getWellResiduals(Base::B_avg_, deferred_logger); std::vector > residual_history; std::vector measure_history; int it = 0; // relaxation factor double relaxation_factor = 1.; const double min_relaxation_factor = 0.6; bool converged = false; int stagnate_count = 0; bool relax_convergence = false; for (; it < max_iter_number; ++it, ++debug_cost_counter_) { assembleWellEqWithoutIteration(ebosSimulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger); const BVectorWell dx_well = mswellhelpers::applyUMFPack(this->duneD_, this->duneDSolver_, this->resWell_); if (it > param_.strict_inner_iter_ms_wells_) relax_convergence = true; const auto report = getWellConvergence(well_state, Base::B_avg_, deferred_logger, relax_convergence); if (report.converged()) { converged = true; break; } residual_history.push_back(this->getWellResiduals(Base::B_avg_, deferred_logger)); measure_history.push_back(this->getResidualMeasureValue(well_state, residual_history[it], param_.tolerance_wells_, param_.tolerance_pressure_ms_wells_, deferred_logger) ); bool is_oscillate = false; bool is_stagnate = false; this->detectOscillations(measure_history, it, is_oscillate, is_stagnate); // TODO: maybe we should have more sophiscated strategy to recover the relaxation factor, // for example, to recover it to be bigger if (is_oscillate || is_stagnate) { // HACK! std::ostringstream sstr; if (relaxation_factor == min_relaxation_factor) { // Still stagnating, terminate iterations if 5 iterations pass. ++stagnate_count; if (stagnate_count == 6) { sstr << " well " << name() << " observes severe stagnation and/or oscillation. We relax the tolerance and check for convergence. \n"; const auto reportStag = getWellConvergence(well_state, Base::B_avg_, deferred_logger, true); if (reportStag.converged()) { converged = true; sstr << " well " << name() << " manages to get converged with relaxed tolerances in " << it << " inner iterations"; deferred_logger.debug(sstr.str()); return converged; } } } // a factor value to reduce the relaxation_factor const double reduction_mutliplier = 0.9; relaxation_factor = std::max(relaxation_factor * reduction_mutliplier, min_relaxation_factor); // debug output if (is_stagnate) { sstr << " well " << name() << " observes stagnation in inner iteration " << it << "\n"; } if (is_oscillate) { sstr << " well " << name() << " observes oscillation in inner iteration " << it << "\n"; } sstr << " relaxation_factor is " << relaxation_factor << " now\n"; deferred_logger.debug(sstr.str()); } updateWellState(dx_well, well_state, deferred_logger, relaxation_factor); initPrimaryVariablesEvaluation(); } // TODO: we should decide whether to keep the updated well_state, or recover to use the old well_state if (converged) { std::ostringstream sstr; sstr << " Well " << name() << " converged in " << it << " inner iterations."; if (relax_convergence) sstr << " (A relaxed tolerance was used after "<< param_.strict_inner_iter_ms_wells_ << " iterations)"; deferred_logger.debug(sstr.str()); } else { std::ostringstream sstr; sstr << " Well " << name() << " did not converge in " << it << " inner iterations."; #define EXTRA_DEBUG_MSW 0 #if EXTRA_DEBUG_MSW sstr << "***** Outputting the residual history for well " << name() << " during inner iterations:"; for (int i = 0; i < it; ++i) { const auto& residual = residual_history[i]; sstr << " residual at " << i << "th iteration "; for (const auto& res : residual) { sstr << " " << res; } sstr << " " << measure_history[i] << " \n"; } #endif deferred_logger.debug(sstr.str()); } return converged; } template void MultisegmentWell:: assembleWellEqWithoutIteration(const Simulator& ebosSimulator, const double dt, const Well::InjectionControls& inj_controls, const Well::ProductionControls& prod_controls, WellState& well_state, const GroupState& group_state, DeferredLogger& deferred_logger) { if (!this->isOperable() && !this->wellIsStopped()) return; // update the upwinding segments this->updateUpwindingSegments(); // calculate the fluid properties needed. computeSegmentFluidProperties(ebosSimulator); // clear all entries this->duneB_ = 0.0; this->duneC_ = 0.0; this->duneD_ = 0.0; this->resWell_ = 0.0; this->duneDSolver_.reset(); well_state.wellVaporizedOilRates(index_of_well_) = 0.; well_state.wellDissolvedGasRates(index_of_well_) = 0.; // for the black oil cases, there will be four equations, // the first three of them are the mass balance equations, the last one is the pressure equations. // // but for the top segment, the pressure equation will be the well control equation, and the other three will be the same. const bool allow_cf = getAllowCrossFlow() || openCrossFlowAvoidSingularity(ebosSimulator); const int nseg = this->numberOfSegments(); for (int seg = 0; seg < nseg; ++seg) { // calculating the accumulation term // TODO: without considering the efficiencty factor for now { const EvalWell segment_surface_volume = getSegmentSurfaceVolume(ebosSimulator, seg); // Add a regularization_factor to increase the accumulation term // This will make the system less stiff and help convergence for // difficult cases const Scalar regularization_factor = param_.regularization_factor_ms_wells_; // for each component for (int comp_idx = 0; comp_idx < num_components_; ++comp_idx) { const EvalWell accumulation_term = regularization_factor * (segment_surface_volume * this->surfaceVolumeFraction(seg, comp_idx) - segment_fluid_initial_[seg][comp_idx]) / dt; this->resWell_[seg][comp_idx] += accumulation_term.value(); for (int pv_idx = 0; pv_idx < numWellEq; ++pv_idx) { this->duneD_[seg][seg][comp_idx][pv_idx] += accumulation_term.derivative(pv_idx + numEq); } } } // considering the contributions due to flowing out from the segment { for (int comp_idx = 0; comp_idx < num_components_; ++comp_idx) { const EvalWell segment_rate = this->getSegmentRateUpwinding(seg, comp_idx) * well_efficiency_factor_; const int seg_upwind = this->upwinding_segments_[seg]; // segment_rate contains the derivatives with respect to GTotal in seg, // and WFrac and GFrac in seg_upwind this->resWell_[seg][comp_idx] -= segment_rate.value(); this->duneD_[seg][seg][comp_idx][GTotal] -= segment_rate.derivative(GTotal + numEq); if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) { this->duneD_[seg][seg_upwind][comp_idx][WFrac] -= segment_rate.derivative(WFrac + numEq); } if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) { this->duneD_[seg][seg_upwind][comp_idx][GFrac] -= segment_rate.derivative(GFrac + numEq); } // pressure derivative should be zero } } // considering the contributions from the inlet segments { for (const int inlet : this->segment_inlets_[seg]) { for (int comp_idx = 0; comp_idx < num_components_; ++comp_idx) { const EvalWell inlet_rate = this->getSegmentRateUpwinding(inlet, comp_idx) * well_efficiency_factor_; const int inlet_upwind = this->upwinding_segments_[inlet]; // inlet_rate contains the derivatives with respect to GTotal in inlet, // and WFrac and GFrac in inlet_upwind this->resWell_[seg][comp_idx] += inlet_rate.value(); this->duneD_[seg][inlet][comp_idx][GTotal] += inlet_rate.derivative(GTotal + numEq); if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) { this->duneD_[seg][inlet_upwind][comp_idx][WFrac] += inlet_rate.derivative(WFrac + numEq); } if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) { this->duneD_[seg][inlet_upwind][comp_idx][GFrac] += inlet_rate.derivative(GFrac + numEq); } // pressure derivative should be zero } } } // calculating the perforation rate for each perforation that belongs to this segment const EvalWell seg_pressure = this->getSegmentPressure(seg); auto& perf_rates = well_state.perfPhaseRates(this->index_of_well_); auto& perf_press_state = well_state.perfPress(this->index_of_well_); for (const int perf : this->segment_perforations_[seg]) { const int cell_idx = well_cells_[perf]; const auto& int_quants = *(ebosSimulator.model().cachedIntensiveQuantities(cell_idx, /*timeIdx=*/ 0)); std::vector mob(num_components_, 0.0); getMobility(ebosSimulator, perf, mob); const double trans_mult = ebosSimulator.problem().template rockCompTransMultiplier(int_quants, cell_idx); const double Tw = well_index_[perf] * trans_mult; std::vector cq_s(num_components_, 0.0); EvalWell perf_press; double perf_dis_gas_rate = 0.; double perf_vap_oil_rate = 0.; computePerfRatePressure(int_quants, mob, Tw, seg, perf, seg_pressure, allow_cf, cq_s, perf_press, perf_dis_gas_rate, perf_vap_oil_rate, deferred_logger); // updating the solution gas rate and solution oil rate if (this->isProducer()) { well_state.wellDissolvedGasRates(index_of_well_) += perf_dis_gas_rate; well_state.wellVaporizedOilRates(index_of_well_) += perf_vap_oil_rate; } // store the perf pressure and rates for (int comp_idx = 0; comp_idx < num_components_; ++comp_idx) { perf_rates[perf*number_of_phases_ + ebosCompIdxToFlowCompIdx(comp_idx)] = cq_s[comp_idx].value(); } perf_press_state[perf] = perf_press.value(); for (int comp_idx = 0; comp_idx < num_components_; ++comp_idx) { // the cq_s entering mass balance equations need to consider the efficiency factors. const EvalWell cq_s_effective = cq_s[comp_idx] * well_efficiency_factor_; connectionRates_[perf][comp_idx] = Base::restrictEval(cq_s_effective); // subtract sum of phase fluxes in the well equations. this->resWell_[seg][comp_idx] += cq_s_effective.value(); // assemble the jacobians for (int pv_idx = 0; pv_idx < numWellEq; ++pv_idx) { // also need to consider the efficiency factor when manipulating the jacobians. this->duneC_[seg][cell_idx][pv_idx][comp_idx] -= cq_s_effective.derivative(pv_idx + numEq); // intput in transformed matrix // the index name for the D should be eq_idx / pv_idx this->duneD_[seg][seg][comp_idx][pv_idx] += cq_s_effective.derivative(pv_idx + numEq); } for (int pv_idx = 0; pv_idx < numEq; ++pv_idx) { // also need to consider the efficiency factor when manipulating the jacobians. this->duneB_[seg][cell_idx][comp_idx][pv_idx] += cq_s_effective.derivative(pv_idx); } } } // the fourth dequation, the pressure drop equation if (seg == 0) { // top segment, pressure equation is the control equation const auto& summaryState = ebosSimulator.vanguard().summaryState(); const Schedule& schedule = ebosSimulator.vanguard().schedule(); this->assembleControlEq(well_state, group_state, schedule, summaryState, inj_controls, prod_controls, getRefDensity(), deferred_logger); } else { const UnitSystem& unit_system = ebosSimulator.vanguard().eclState().getDeckUnitSystem(); this->assemblePressureEq(seg, unit_system, well_state, deferred_logger); } } } template bool MultisegmentWell:: openCrossFlowAvoidSingularity(const Simulator& ebos_simulator) const { return !getAllowCrossFlow() && allDrawDownWrongDirection(ebos_simulator); } template bool MultisegmentWell:: allDrawDownWrongDirection(const Simulator& ebos_simulator) const { bool all_drawdown_wrong_direction = true; const int nseg = this->numberOfSegments(); for (int seg = 0; seg < nseg; ++seg) { const EvalWell segment_pressure = this->getSegmentPressure(seg); for (const int perf : this->segment_perforations_[seg]) { const int cell_idx = well_cells_[perf]; const auto& intQuants = *(ebos_simulator.model().cachedIntensiveQuantities(cell_idx, /*timeIdx=*/ 0)); const auto& fs = intQuants.fluidState(); // pressure difference between the segment and the perforation const EvalWell perf_seg_press_diff = gravity_ * this->segment_densities_[seg] * this->perforation_segment_depth_diffs_[perf]; // pressure difference between the perforation and the grid cell const double cell_perf_press_diff = this->cell_perforation_pressure_diffs_[perf]; const double pressure_cell = (fs.pressure(FluidSystem::oilPhaseIdx)).value(); const double perf_press = pressure_cell - cell_perf_press_diff; // Pressure drawdown (also used to determine direction of flow) // TODO: not 100% sure about the sign of the seg_perf_press_diff const EvalWell drawdown = perf_press - (segment_pressure + perf_seg_press_diff); // for now, if there is one perforation can produce/inject in the correct // direction, we consider this well can still produce/inject. // TODO: it can be more complicated than this to cause wrong-signed rates if ( (drawdown < 0. && this->isInjector()) || (drawdown > 0. && this->isProducer()) ) { all_drawdown_wrong_direction = false; break; } } } return all_drawdown_wrong_direction; } template void MultisegmentWell:: updateWaterThroughput(const double /*dt*/, WellState& /*well_state*/) const { } template typename MultisegmentWell::EvalWell MultisegmentWell:: getSegmentSurfaceVolume(const Simulator& ebos_simulator, const int seg_idx) const { EvalWell temperature; EvalWell saltConcentration; int pvt_region_index; { // using the pvt region of first perforated cell // TODO: it should be a member of the WellInterface, initialized properly const int cell_idx = well_cells_[0]; const auto& intQuants = *(ebos_simulator.model().cachedIntensiveQuantities(cell_idx, /*timeIdx=*/0)); const auto& fs = intQuants.fluidState(); temperature.setValue(fs.temperature(FluidSystem::oilPhaseIdx).value()); saltConcentration = this->extendEval(fs.saltConcentration()); pvt_region_index = fs.pvtRegionIndex(); } return this->MSWEval::getSegmentSurfaceVolume(temperature, saltConcentration, pvt_region_index, seg_idx); } template std::optional MultisegmentWell:: computeBhpAtThpLimitProd(const Simulator& ebos_simulator, const SummaryState& summary_state, DeferredLogger& deferred_logger) const { // Make the frates() function. auto frates = [this, &ebos_simulator, &deferred_logger](const double bhp) { // Not solving the well equations here, which means we are // calculating at the current Fg/Fw values of the // well. This does not matter unless the well is // crossflowing, and then it is likely still a good // approximation. std::vector rates(3); computeWellRatesWithBhp(ebos_simulator, bhp, rates, deferred_logger); return rates; }; return this->MultisegmentWellGeneric:: computeBhpAtThpLimitProd(frates, summary_state, maxPerfPress(ebos_simulator), getRefDensity(), deferred_logger); } template std::optional MultisegmentWell:: computeBhpAtThpLimitInj(const Simulator& ebos_simulator, const SummaryState& summary_state, DeferredLogger& deferred_logger) const { // Make the frates() function. auto frates = [this, &ebos_simulator, &deferred_logger](const double bhp) { // Not solving the well equations here, which means we are // calculating at the current Fg/Fw values of the // well. This does not matter unless the well is // crossflowing, and then it is likely still a good // approximation. std::vector rates(3); computeWellRatesWithBhp(ebos_simulator, bhp, rates, deferred_logger); return rates; }; return this->MultisegmentWellGeneric:: computeBhpAtThpLimitInj(frates, summary_state, getRefDensity(), deferred_logger); } template double MultisegmentWell:: maxPerfPress(const Simulator& ebos_simulator) const { double max_pressure = 0.0; const int nseg = this->numberOfSegments(); for (int seg = 0; seg < nseg; ++seg) { for (const int perf : this->segment_perforations_[seg]) { const int cell_idx = well_cells_[perf]; const auto& int_quants = *(ebos_simulator.model().cachedIntensiveQuantities(cell_idx, /*timeIdx=*/ 0)); const auto& fs = int_quants.fluidState(); double pressure_cell = fs.pressure(FluidSystem::oilPhaseIdx).value(); max_pressure = std::max(max_pressure, pressure_cell); } } return max_pressure; } template std::vector MultisegmentWell:: computeCurrentWellRates(const Simulator& ebosSimulator, DeferredLogger& deferred_logger) const { // Calculate the rates that follow from the current primary variables. std::vector well_q_s(num_components_, 0.0); const bool allow_cf = getAllowCrossFlow() || openCrossFlowAvoidSingularity(ebosSimulator); const int nseg = this->numberOfSegments(); for (int seg = 0; seg < nseg; ++seg) { // calculating the perforation rate for each perforation that belongs to this segment const EvalWell seg_pressure = this->getSegmentPressure(seg); for (const int perf : this->segment_perforations_[seg]) { const int cell_idx = well_cells_[perf]; const auto& int_quants = *(ebosSimulator.model().cachedIntensiveQuantities(cell_idx, /*timeIdx=*/ 0)); std::vector mob(num_components_, 0.0); getMobility(ebosSimulator, perf, mob); const double trans_mult = ebosSimulator.problem().template rockCompTransMultiplier(int_quants, cell_idx); const double Tw = well_index_[perf] * trans_mult; std::vector cq_s(num_components_, 0.0); EvalWell perf_press; double perf_dis_gas_rate = 0.; double perf_vap_oil_rate = 0.; computePerfRatePressure(int_quants, mob, Tw, seg, perf, seg_pressure, allow_cf, cq_s, perf_press, perf_dis_gas_rate, perf_vap_oil_rate, deferred_logger); for (int comp = 0; comp < num_components_; ++comp) { well_q_s[comp] += cq_s[comp]; } } } std::vector well_q_s_noderiv(well_q_s.size()); for (int comp = 0; comp < num_components_; ++comp) { well_q_s_noderiv[comp] = well_q_s[comp].value(); } return well_q_s_noderiv; } template void MultisegmentWell:: computeConnLevelProdInd(const typename MultisegmentWell::FluidState& fs, const std::function& connPICalc, const std::vector& mobility, double* connPI) const { const auto& pu = this->phaseUsage(); const int np = this->number_of_phases_; for (int p = 0; p < np; ++p) { // Note: E100's notion of PI value phase mobility includes // the reciprocal FVF. const auto connMob = mobility[ flowPhaseToEbosCompIdx(p) ].value() * fs.invB(flowPhaseToEbosPhaseIdx(p)).value(); connPI[p] = connPICalc(connMob); } if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) { const auto io = pu.phase_pos[Oil]; const auto ig = pu.phase_pos[Gas]; const auto vapoil = connPI[ig] * fs.Rv().value(); const auto disgas = connPI[io] * fs.Rs().value(); connPI[io] += vapoil; connPI[ig] += disgas; } } template void MultisegmentWell:: computeConnLevelInjInd(const typename MultisegmentWell::FluidState& fs, const Phase preferred_phase, const std::function& connIICalc, const std::vector& mobility, double* connII, DeferredLogger& deferred_logger) const { // Assumes single phase injection const auto& pu = this->phaseUsage(); auto phase_pos = 0; if (preferred_phase == Phase::GAS) { phase_pos = pu.phase_pos[Gas]; } else if (preferred_phase == Phase::OIL) { phase_pos = pu.phase_pos[Oil]; } else if (preferred_phase == Phase::WATER) { phase_pos = pu.phase_pos[Water]; } else { OPM_DEFLOG_THROW(NotImplemented, "Unsupported Injector Type (" << static_cast(preferred_phase) << ") for well " << this->name() << " during connection I.I. calculation", deferred_logger); } const auto zero = EvalWell { 0.0 }; const auto mt = std::accumulate(mobility.begin(), mobility.end(), zero); connII[phase_pos] = connIICalc(mt.value() * fs.invB(flowPhaseToEbosPhaseIdx(phase_pos)).value()); } } // namespace Opm