mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Report Solution GOR and OGR at Segment Level
This enables calculating the SOFRF and SGFRF summary vectors.
This commit is contained in:
@@ -1419,6 +1419,8 @@ updateWellStateFromPrimaryVariables(WellState& well_state,
|
|||||||
auto& ws = well_state.well(baseif_.indexOfWell());
|
auto& ws = well_state.well(baseif_.indexOfWell());
|
||||||
auto& segments = ws.segments;
|
auto& segments = ws.segments;
|
||||||
auto& segment_rates = segments.rates;
|
auto& segment_rates = segments.rates;
|
||||||
|
auto& disgas = segments.dissolved_gas_rate;
|
||||||
|
auto& vapoil = segments.vaporized_oil_rate;
|
||||||
auto& segment_pressure = segments.pressure;
|
auto& segment_pressure = segments.pressure;
|
||||||
for (int seg = 0; seg < this->numberOfSegments(); ++seg) {
|
for (int seg = 0; seg < this->numberOfSegments(); ++seg) {
|
||||||
std::vector<double> fractions(baseif_.numPhases(), 0.0);
|
std::vector<double> fractions(baseif_.numPhases(), 0.0);
|
||||||
@@ -1482,10 +1484,23 @@ updateWellStateFromPrimaryVariables(WellState& well_state,
|
|||||||
.saturatedOilVaporizationFactor(pvtReg, temperature, segment_pressure[seg]);
|
.saturatedOilVaporizationFactor(pvtReg, temperature, segment_pressure[seg]);
|
||||||
}
|
}
|
||||||
|
|
||||||
// 1) Local condition volume flow rates
|
// 1) Infer phase splitting for oil/gas.
|
||||||
const auto& [Rs, Rv] = this->baseif_.rateConverter().inferDissolvedVaporisedRatio
|
const auto& [Rs, Rv] = this->baseif_.rateConverter().inferDissolvedVaporisedRatio
|
||||||
(rsMax, rvMax, segment_rates.begin() + (seg + 0)*this->baseif_.numPhases());
|
(rsMax, rvMax, segment_rates.begin() + (seg + 0)*this->baseif_.numPhases());
|
||||||
|
|
||||||
|
if (! FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
|
||||||
|
vapoil[seg] = disgas[seg] = 0.0;
|
||||||
|
}
|
||||||
|
else {
|
||||||
|
const auto* qs = &segment_rates[seg*this->baseif_.numPhases() + 0];
|
||||||
|
const auto denom = 1.0 - (Rs * Rv);
|
||||||
|
const auto io = pu.phase_pos[Oil];
|
||||||
|
const auto ig = pu.phase_pos[Gas];
|
||||||
|
disgas[seg] = Rs * (qs[io] - Rv*qs[ig]) / denom;
|
||||||
|
vapoil[seg] = Rv * (qs[ig] - Rs*qs[io]) / denom;
|
||||||
|
}
|
||||||
|
|
||||||
|
// 2) Local condition volume flow rates
|
||||||
{
|
{
|
||||||
// Use std::span<> in C++20 and beyond.
|
// Use std::span<> in C++20 and beyond.
|
||||||
const auto rate_start = (seg + 0) * this->baseif_.numPhases();
|
const auto rate_start = (seg + 0) * this->baseif_.numPhases();
|
||||||
@@ -1499,7 +1514,7 @@ updateWellStateFromPrimaryVariables(WellState& well_state,
|
|||||||
temperature, saltConc, surf_rates, resv_rates);
|
temperature, saltConc, surf_rates, resv_rates);
|
||||||
}
|
}
|
||||||
|
|
||||||
// 2) Local condition holdup fractions.
|
// 3) Local condition holdup fractions.
|
||||||
const auto tot_resv =
|
const auto tot_resv =
|
||||||
std::accumulate(segments.phase_resv_rates.begin() + (seg + 0)*this->baseif_.numPhases(),
|
std::accumulate(segments.phase_resv_rates.begin() + (seg + 0)*this->baseif_.numPhases(),
|
||||||
segments.phase_resv_rates.begin() + (seg + 1)*this->baseif_.numPhases(),
|
segments.phase_resv_rates.begin() + (seg + 1)*this->baseif_.numPhases(),
|
||||||
@@ -1510,7 +1525,7 @@ updateWellStateFromPrimaryVariables(WellState& well_state,
|
|||||||
segments.phase_holdup.begin() + (seg + 0)*this->baseif_.numPhases(),
|
segments.phase_holdup.begin() + (seg + 0)*this->baseif_.numPhases(),
|
||||||
[tot_resv](const auto qr) { return std::clamp(qr / tot_resv, 0.0, 1.0); });
|
[tot_resv](const auto qr) { return std::clamp(qr / tot_resv, 0.0, 1.0); });
|
||||||
|
|
||||||
// 3) Local condition flow velocities for segments other than top segment.
|
// 4) Local condition flow velocities for segments other than top segment.
|
||||||
if (seg > 0) {
|
if (seg > 0) {
|
||||||
// Possibly poor approximation
|
// Possibly poor approximation
|
||||||
// Velocity = Flow rate / cross-sectional area.
|
// Velocity = Flow rate / cross-sectional area.
|
||||||
@@ -1525,7 +1540,7 @@ updateWellStateFromPrimaryVariables(WellState& well_state,
|
|||||||
[velocity](const auto hf) { return (hf > 0.0) ? velocity : 0.0; });
|
[velocity](const auto hf) { return (hf > 0.0) ? velocity : 0.0; });
|
||||||
}
|
}
|
||||||
|
|
||||||
// 4) Local condition phase viscosities.
|
// 5) Local condition phase viscosities.
|
||||||
segments.phase_viscosity[seg*this->baseif_.numPhases() + pu.phase_pos[Oil]] =
|
segments.phase_viscosity[seg*this->baseif_.numPhases() + pu.phase_pos[Oil]] =
|
||||||
FluidSystem::oilPvt().viscosity(pvtReg, temperature, segment_pressure[seg], Rs);
|
FluidSystem::oilPvt().viscosity(pvtReg, temperature, segment_pressure[seg], Rs);
|
||||||
|
|
||||||
|
|||||||
@@ -1,5 +1,5 @@
|
|||||||
/*
|
/*
|
||||||
Copyright Equinor ASA 2021
|
Copyright Equinor ASA 2021, 2022.
|
||||||
|
|
||||||
This file is part of the Open Porous Media project (OPM).
|
This file is part of the Open Porous Media project (OPM).
|
||||||
|
|
||||||
@@ -55,6 +55,8 @@ namespace Opm
|
|||||||
|
|
||||||
SegmentState::SegmentState(int num_phases, const WellSegments& segments)
|
SegmentState::SegmentState(int num_phases, const WellSegments& segments)
|
||||||
: rates (segments.size() * num_phases)
|
: rates (segments.size() * num_phases)
|
||||||
|
, dissolved_gas_rate (segments.size())
|
||||||
|
, vaporized_oil_rate (segments.size())
|
||||||
, phase_resv_rates (segments.size() * num_phases)
|
, phase_resv_rates (segments.size() * num_phases)
|
||||||
, phase_velocity (segments.size() * num_phases)
|
, phase_velocity (segments.size() * num_phases)
|
||||||
, phase_holdup (segments.size() * num_phases)
|
, phase_holdup (segments.size() * num_phases)
|
||||||
|
|||||||
@@ -1,5 +1,5 @@
|
|||||||
/*
|
/*
|
||||||
Copyright Equinor ASA 2021
|
Copyright Equinor ASA 2021, 2022.
|
||||||
|
|
||||||
This file is part of the Open Porous Media project (OPM).
|
This file is part of the Open Porous Media project (OPM).
|
||||||
|
|
||||||
@@ -45,6 +45,8 @@ public:
|
|||||||
std::size_t size() const;
|
std::size_t size() const;
|
||||||
|
|
||||||
std::vector<double> rates;
|
std::vector<double> rates;
|
||||||
|
std::vector<double> dissolved_gas_rate;
|
||||||
|
std::vector<double> vaporized_oil_rate;
|
||||||
|
|
||||||
/// Segment condition volume flow rates through segment (per phase)
|
/// Segment condition volume flow rates through segment (per phase)
|
||||||
std::vector<double> phase_resv_rates;
|
std::vector<double> phase_resv_rates;
|
||||||
|
|||||||
@@ -873,6 +873,7 @@ WellState::reportSegmentResults(const int well_id,
|
|||||||
const auto io = pu.phase_pos[Oil];
|
const auto io = pu.phase_pos[Oil];
|
||||||
|
|
||||||
seg_res.rates.set(data::Rates::opt::oil, rate[io]);
|
seg_res.rates.set(data::Rates::opt::oil, rate[io]);
|
||||||
|
seg_res.rates.set(data::Rates::opt::vaporized_oil, segments.vaporized_oil_rate[seg_ix]);
|
||||||
seg_res.rates.set(data::Rates::opt::reservoir_oil, resv[io]);
|
seg_res.rates.set(data::Rates::opt::reservoir_oil, resv[io]);
|
||||||
seg_res.velocity.set(PhaseQuant::Oil, velocity[io]);
|
seg_res.velocity.set(PhaseQuant::Oil, velocity[io]);
|
||||||
seg_res.holdup.set(PhaseQuant::Oil, holdup[io]);
|
seg_res.holdup.set(PhaseQuant::Oil, holdup[io]);
|
||||||
@@ -883,6 +884,7 @@ WellState::reportSegmentResults(const int well_id,
|
|||||||
const auto ig = pu.phase_pos[Gas];
|
const auto ig = pu.phase_pos[Gas];
|
||||||
|
|
||||||
seg_res.rates.set(data::Rates::opt::gas, rate[ig]);
|
seg_res.rates.set(data::Rates::opt::gas, rate[ig]);
|
||||||
|
seg_res.rates.set(data::Rates::opt::dissolved_gas, segments.dissolved_gas_rate[seg_ix]);
|
||||||
seg_res.rates.set(data::Rates::opt::reservoir_gas, resv[ig]);
|
seg_res.rates.set(data::Rates::opt::reservoir_gas, resv[ig]);
|
||||||
seg_res.velocity.set(PhaseQuant::Gas, velocity[ig]);
|
seg_res.velocity.set(PhaseQuant::Gas, velocity[ig]);
|
||||||
seg_res.holdup.set(PhaseQuant::Gas, holdup[ig]);
|
seg_res.holdup.set(PhaseQuant::Gas, holdup[ig]);
|
||||||
|
|||||||
Reference in New Issue
Block a user