StandardWell: move connectionRateBrine to Connections class

This commit is contained in:
Arne Morten Kvarving 2023-05-22 11:10:02 +02:00
parent 41a38cc9d6
commit 73ece47d97
4 changed files with 46 additions and 33 deletions

View File

@ -438,11 +438,6 @@ namespace Opm
DeferredLogger& deferred_logger) const; DeferredLogger& deferred_logger) const;
private: private:
Eval connectionRateBrine(double& rate,
const double vap_wat_rate,
const std::vector<EvalWell>& cq_s,
const IntensiveQuantities& intQuants) const;
Eval connectionRateEnergy(const double maxOilSaturation, Eval connectionRateEnergy(const double maxOilSaturation,
const std::vector<EvalWell>& cq_s, const std::vector<EvalWell>& cq_s,
const IntensiveQuantities& intQuants, const IntensiveQuantities& intQuants,

View File

@ -514,6 +514,31 @@ computeProperties(const WellState& well_state,
this->computePressureDelta(); this->computePressureDelta();
} }
template<class FluidSystem, class Indices, class Scalar>
typename StandardWellConnections<FluidSystem,Indices,Scalar>::Eval
StandardWellConnections<FluidSystem,Indices,Scalar>::
connectionRateBrine(double& rate,
const double vap_wat_rate,
const std::vector<EvalWell>& cq_s,
const std::variant<Scalar,EvalWell>& saltConcentration) const
{
// TODO: the application of well efficiency factor has not been tested with an example yet
const unsigned waterCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
// Correction salt rate; evaporated water does not contain salt
EvalWell cq_s_sm = cq_s[waterCompIdx] - vap_wat_rate;
if (well_.isInjector()) {
cq_s_sm *= std::get<Scalar>(saltConcentration);
} else {
cq_s_sm *= std::get<EvalWell>(saltConcentration);
}
// Note. Efficiency factor is handled in the output layer
rate = cq_s_sm.value();
cq_s_sm *= well_.wellEfficiencyFactor();
return well_.restrictEval(cq_s_sm);
}
#define INSTANCE(...) \ #define INSTANCE(...) \
template class StandardWellConnections<BlackOilFluidSystem<double,BlackOilDefaultIndexTraits>,__VA_ARGS__,double>; template class StandardWellConnections<BlackOilFluidSystem<double,BlackOilDefaultIndexTraits>,__VA_ARGS__,double>;

View File

@ -23,7 +23,10 @@
#ifndef OPM_STANDARDWELL_CONNECTIONS_HEADER_INCLUDED #ifndef OPM_STANDARDWELL_CONNECTIONS_HEADER_INCLUDED
#define OPM_STANDARDWELL_CONNECTIONS_HEADER_INCLUDED #define OPM_STANDARDWELL_CONNECTIONS_HEADER_INCLUDED
#include <opm/simulators/wells/StandardWellPrimaryVariables.hpp>
#include <functional> #include <functional>
#include <variant>
#include <vector> #include <vector>
namespace Opm namespace Opm
@ -76,6 +79,14 @@ public:
Scalar pressure_diff(const unsigned perf) const Scalar pressure_diff(const unsigned perf) const
{ return perf_pressure_diffs_[perf]; } { return perf_pressure_diffs_[perf]; }
using Eval = typename WellInterfaceIndices<FluidSystem,Indices,Scalar>::Eval;
using EvalWell = typename StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>::EvalWell;
Eval connectionRateBrine(double& rate,
const double vap_wat_rate,
const std::vector<EvalWell>& cq_s,
const std::variant<Scalar,EvalWell>& saltConcentration) const;
private: private:
void computePressureDelta(); void computePressureDelta();

View File

@ -549,9 +549,17 @@ namespace Opm
} }
if constexpr (has_brine) { if constexpr (has_brine) {
std::variant<Scalar,EvalWell> saltConcentration;
if (this->isInjector()) {
saltConcentration = this->wsalt();
} else {
saltConcentration = this->extendEval(intQuants.fluidState().saltConcentration());
}
connectionRates[perf][Indices::contiBrineEqIdx] = connectionRates[perf][Indices::contiBrineEqIdx] =
connectionRateBrine(perf_data.brine_rates[perf], this->connections_.connectionRateBrine(perf_data.brine_rates[perf],
perf_rates.vap_wat, cq_s, intQuants); perf_rates.vap_wat, cq_s,
saltConcentration);
} }
if constexpr (has_micp) { if constexpr (has_micp) {
@ -2136,32 +2144,6 @@ namespace Opm
} }
template <typename TypeTag>
typename StandardWell<TypeTag>::Eval
StandardWell<TypeTag>::
connectionRateBrine(double& rate,
const double vap_wat_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 waterCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
// Correction salt rate; evaporated water does not contain salt
EvalWell cq_s_sm = cq_s[waterCompIdx] - vap_wat_rate;
if (this->isInjector()) {
cq_s_sm *= this->wsalt();
} else {
cq_s_sm *= this->extendEval(intQuants.fluidState().saltConcentration());
}
// Note. Efficiency factor is handled in the output layer
rate = cq_s_sm.value();
cq_s_sm *= this->well_efficiency_factor_;
return Base::restrictEval(cq_s_sm);
}
template <typename TypeTag> template <typename TypeTag>
typename StandardWell<TypeTag>::Eval typename StandardWell<TypeTag>::Eval
StandardWell<TypeTag>:: StandardWell<TypeTag>::