mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
StandardWell: move connectionRatezFraction to Connections class
This commit is contained in:
parent
56014ccff9
commit
62a2ee1713
@ -443,12 +443,6 @@ namespace Opm
|
||||
const IntensiveQuantities& intQuants,
|
||||
DeferredLogger& deferred_logger) const;
|
||||
|
||||
std::tuple<Eval,EvalWell>
|
||||
connectionRatezFraction(double& rate,
|
||||
const double dis_gas_rate,
|
||||
const std::vector<EvalWell>& cq_s,
|
||||
const IntensiveQuantities& intQuants) const;
|
||||
|
||||
template<class Value>
|
||||
void gasOilPerfRateInj(const std::vector<Value>& cq_s,
|
||||
PerforationRates& perf_rates,
|
||||
|
@ -639,6 +639,32 @@ connectionRatePolymer(double& rate,
|
||||
return {well_.restrictEval(cq_s_poly), cq_s_poly};
|
||||
}
|
||||
|
||||
template<class FluidSystem, class Indices, class Scalar>
|
||||
std::tuple<typename StandardWellConnections<FluidSystem,Indices,Scalar>::Eval,
|
||||
typename StandardWellConnections<FluidSystem,Indices,Scalar>::EvalWell>
|
||||
StandardWellConnections<FluidSystem,Indices,Scalar>::
|
||||
connectionRatezFraction(double& rate,
|
||||
const double dis_gas_rate,
|
||||
const std::vector<EvalWell>& cq_s,
|
||||
const std::variant<Scalar, std::array<EvalWell,2>>& solventConcentration) const
|
||||
{
|
||||
// TODO: the application of well efficiency factor has not been tested with an example yet
|
||||
const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
|
||||
EvalWell cq_s_zfrac_effective = cq_s[gasCompIdx];
|
||||
if (well_.isInjector()) {
|
||||
cq_s_zfrac_effective *= std::get<Scalar>(solventConcentration);
|
||||
} else if (cq_s_zfrac_effective.value() != 0.0) {
|
||||
const double dis_gas_frac = dis_gas_rate / cq_s_zfrac_effective.value();
|
||||
const auto& vol = std::get<std::array<EvalWell,2>>(solventConcentration);
|
||||
cq_s_zfrac_effective *= dis_gas_frac * vol[0] + (1.0 - dis_gas_frac) * vol[1];
|
||||
}
|
||||
|
||||
rate = cq_s_zfrac_effective.value();
|
||||
|
||||
cq_s_zfrac_effective *= well_.wellEfficiencyFactor();
|
||||
return {well_.restrictEval(cq_s_zfrac_effective), cq_s_zfrac_effective};
|
||||
}
|
||||
|
||||
#define INSTANCE(...) \
|
||||
template class StandardWellConnections<BlackOilFluidSystem<double,BlackOilDefaultIndexTraits>, \
|
||||
__VA_ARGS__,double>;
|
||||
|
@ -104,6 +104,12 @@ public:
|
||||
const std::variant<Scalar,EvalWell>& oxygenConcentration,
|
||||
const std::variant<Scalar,EvalWell>& ureaConcentration) const;
|
||||
|
||||
std::tuple<Eval,EvalWell>
|
||||
connectionRatezFraction(double& rate,
|
||||
const double dis_gas_rate,
|
||||
const std::vector<EvalWell>& cq_s,
|
||||
const std::variant<Scalar, std::array<EvalWell,2>>& solventConcentration) const;
|
||||
|
||||
private:
|
||||
void computePressureDelta();
|
||||
|
||||
|
@ -558,10 +558,18 @@ namespace Opm
|
||||
}
|
||||
|
||||
if constexpr (has_zFraction) {
|
||||
std::variant<Scalar,std::array<EvalWell,2>> solventConcentration;
|
||||
if (this->isInjector()) {
|
||||
solventConcentration = this->wsolvent();
|
||||
} else {
|
||||
solventConcentration = std::array{this->extendEval(intQuants.xVolume()),
|
||||
this->extendEval(intQuants.yVolume())};
|
||||
}
|
||||
std::tie(connectionRates[perf][Indices::contiZfracEqIdx],
|
||||
cq_s_zfrac_effective) =
|
||||
connectionRatezFraction(perf_data.solvent_rates[perf],
|
||||
perf_rates.dis_gas, cq_s, intQuants);
|
||||
this->connections_.connectionRatezFraction(perf_data.solvent_rates[perf],
|
||||
perf_rates.dis_gas, cq_s,
|
||||
solventConcentration);
|
||||
}
|
||||
|
||||
if constexpr (has_brine) {
|
||||
@ -2256,32 +2264,6 @@ namespace Opm
|
||||
}
|
||||
|
||||
|
||||
template <typename TypeTag>
|
||||
std::tuple<typename StandardWell<TypeTag>::Eval,
|
||||
typename StandardWell<TypeTag>::EvalWell>
|
||||
StandardWell<TypeTag>::
|
||||
connectionRatezFraction(double& rate,
|
||||
const double dis_gas_rate,
|
||||
const std::vector<EvalWell>& cq_s,
|
||||
const IntensiveQuantities& intQuants) const
|
||||
{
|
||||
// TODO: the application of well efficiency factor has not been tested with an example yet
|
||||
const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
|
||||
EvalWell cq_s_zfrac_effective = cq_s[gasCompIdx];
|
||||
if (this->isInjector()) {
|
||||
cq_s_zfrac_effective *= this->wsolvent();
|
||||
} else if (cq_s_zfrac_effective.value() != 0.0) {
|
||||
const double dis_gas_frac = dis_gas_rate / cq_s_zfrac_effective.value();
|
||||
cq_s_zfrac_effective *= this->extendEval(dis_gas_frac*intQuants.xVolume() + (1.0-dis_gas_frac)*intQuants.yVolume());
|
||||
}
|
||||
|
||||
rate = cq_s_zfrac_effective.value();
|
||||
|
||||
cq_s_zfrac_effective *= this->well_efficiency_factor_;
|
||||
return {Base::restrictEval(cq_s_zfrac_effective), cq_s_zfrac_effective};
|
||||
}
|
||||
|
||||
|
||||
template <typename TypeTag>
|
||||
template<class Value>
|
||||
void
|
||||
|
Loading…
Reference in New Issue
Block a user