mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
StandardWell: move connectionRateFoam to Connections class
This commit is contained in:
parent
73ece47d97
commit
2f6623993a
@ -443,10 +443,6 @@ namespace Opm
|
|||||||
const IntensiveQuantities& intQuants,
|
const IntensiveQuantities& intQuants,
|
||||||
DeferredLogger& deferred_logger) const;
|
DeferredLogger& deferred_logger) const;
|
||||||
|
|
||||||
Eval connectionRateFoam(const std::vector<EvalWell>& cq_s,
|
|
||||||
const IntensiveQuantities& intQuants,
|
|
||||||
DeferredLogger& deferred_logger) const;
|
|
||||||
|
|
||||||
std::tuple<Eval,Eval,Eval>
|
std::tuple<Eval,Eval,Eval>
|
||||||
connectionRatesMICP(const std::vector<EvalWell>& cq_s,
|
connectionRatesMICP(const std::vector<EvalWell>& cq_s,
|
||||||
const IntensiveQuantities& intQuants) const;
|
const IntensiveQuantities& intQuants) const;
|
||||||
|
@ -30,7 +30,7 @@
|
|||||||
#include <opm/models/blackoil/blackoilonephaseindices.hh>
|
#include <opm/models/blackoil/blackoilonephaseindices.hh>
|
||||||
#include <opm/models/blackoil/blackoiltwophaseindices.hh>
|
#include <opm/models/blackoil/blackoiltwophaseindices.hh>
|
||||||
|
|
||||||
#include <opm/simulators/utils/DeferredLogger.hpp>
|
#include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>
|
||||||
#include <opm/simulators/wells/ParallelWellInfo.hpp>
|
#include <opm/simulators/wells/ParallelWellInfo.hpp>
|
||||||
#include <opm/simulators/wells/WellInterfaceIndices.hpp>
|
#include <opm/simulators/wells/WellInterfaceIndices.hpp>
|
||||||
#include <opm/simulators/wells/WellState.hpp>
|
#include <opm/simulators/wells/WellState.hpp>
|
||||||
@ -539,8 +539,48 @@ connectionRateBrine(double& rate,
|
|||||||
return well_.restrictEval(cq_s_sm);
|
return well_.restrictEval(cq_s_sm);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
template<class FluidSystem, class Indices, class Scalar>
|
||||||
|
typename StandardWellConnections<FluidSystem,Indices,Scalar>::Eval
|
||||||
|
StandardWellConnections<FluidSystem,Indices,Scalar>::
|
||||||
|
connectionRateFoam(const std::vector<EvalWell>& cq_s,
|
||||||
|
const std::variant<Scalar,EvalWell>& foamConcentration,
|
||||||
|
const Phase transportPhase,
|
||||||
|
DeferredLogger& deferred_logger) const
|
||||||
|
{
|
||||||
|
// TODO: the application of well efficiency factor has not been tested with an example yet
|
||||||
|
auto getFoamTransportIdx = [&deferred_logger,transportPhase] {
|
||||||
|
switch (transportPhase) {
|
||||||
|
case Phase::WATER: {
|
||||||
|
return Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
|
||||||
|
}
|
||||||
|
case Phase::GAS: {
|
||||||
|
return Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
|
||||||
|
}
|
||||||
|
case Phase::SOLVENT: {
|
||||||
|
if constexpr (Indices::enableSolvent)
|
||||||
|
return static_cast<unsigned>(Indices::contiSolventEqIdx);
|
||||||
|
else
|
||||||
|
OPM_DEFLOG_THROW(std::runtime_error, "Foam transport phase is SOLVENT but SOLVENT is not activated.", deferred_logger);
|
||||||
|
}
|
||||||
|
default: {
|
||||||
|
OPM_DEFLOG_THROW(std::runtime_error, "Foam transport phase must be GAS/WATER/SOLVENT.", deferred_logger);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
};
|
||||||
|
|
||||||
|
EvalWell cq_s_foam = cq_s[getFoamTransportIdx()] * well_.wellEfficiencyFactor();;
|
||||||
|
if (well_.isInjector()) {
|
||||||
|
cq_s_foam *= std::get<Scalar>(foamConcentration);
|
||||||
|
} else {
|
||||||
|
cq_s_foam *= std::get<EvalWell>(foamConcentration);
|
||||||
|
}
|
||||||
|
|
||||||
|
return well_.restrictEval(cq_s_foam);
|
||||||
|
}
|
||||||
|
|
||||||
#define INSTANCE(...) \
|
#define INSTANCE(...) \
|
||||||
template class StandardWellConnections<BlackOilFluidSystem<double,BlackOilDefaultIndexTraits>,__VA_ARGS__,double>;
|
template class StandardWellConnections<BlackOilFluidSystem<double,BlackOilDefaultIndexTraits>, \
|
||||||
|
__VA_ARGS__,double>;
|
||||||
|
|
||||||
// One phase
|
// One phase
|
||||||
INSTANCE(BlackOilOnePhaseIndices<0u,0u,0u,0u,false,false,0u,1u,0u>)
|
INSTANCE(BlackOilOnePhaseIndices<0u,0u,0u,0u,false,false,0u,1u,0u>)
|
||||||
|
@ -33,6 +33,7 @@ namespace Opm
|
|||||||
{
|
{
|
||||||
|
|
||||||
class DeferredLogger;
|
class DeferredLogger;
|
||||||
|
enum class Phase;
|
||||||
template<class FluidSystem, class Indices, class Scalar> class WellInterfaceIndices;
|
template<class FluidSystem, class Indices, class Scalar> class WellInterfaceIndices;
|
||||||
class WellState;
|
class WellState;
|
||||||
|
|
||||||
@ -87,6 +88,11 @@ public:
|
|||||||
const std::vector<EvalWell>& cq_s,
|
const std::vector<EvalWell>& cq_s,
|
||||||
const std::variant<Scalar,EvalWell>& saltConcentration) const;
|
const std::variant<Scalar,EvalWell>& saltConcentration) const;
|
||||||
|
|
||||||
|
Eval connectionRateFoam(const std::vector<EvalWell>& cq_s,
|
||||||
|
const std::variant<Scalar,EvalWell>& foamConcentration,
|
||||||
|
const Phase transportPhase,
|
||||||
|
DeferredLogger& deferred_logger) const;
|
||||||
|
|
||||||
private:
|
private:
|
||||||
void computePressureDelta();
|
void computePressureDelta();
|
||||||
|
|
||||||
|
@ -537,8 +537,16 @@ namespace Opm
|
|||||||
}
|
}
|
||||||
|
|
||||||
if constexpr (has_foam) {
|
if constexpr (has_foam) {
|
||||||
|
std::variant<Scalar,EvalWell> foamConcentration;
|
||||||
|
if (this->isInjector()) {
|
||||||
|
foamConcentration = this->wfoam();
|
||||||
|
} else {
|
||||||
|
foamConcentration = this->extendEval(intQuants.foamConcentration());
|
||||||
|
}
|
||||||
connectionRates[perf][Indices::contiFoamEqIdx] =
|
connectionRates[perf][Indices::contiFoamEqIdx] =
|
||||||
connectionRateFoam(cq_s, intQuants, deferred_logger);
|
this->connections_.connectionRateFoam(cq_s, foamConcentration,
|
||||||
|
FoamModule::transportPhase(),
|
||||||
|
deferred_logger);
|
||||||
}
|
}
|
||||||
|
|
||||||
if constexpr (has_zFraction) {
|
if constexpr (has_zFraction) {
|
||||||
@ -2225,43 +2233,6 @@ namespace Opm
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
template <typename TypeTag>
|
|
||||||
typename StandardWell<TypeTag>::Eval
|
|
||||||
StandardWell<TypeTag>::
|
|
||||||
connectionRateFoam(const std::vector<EvalWell>& cq_s,
|
|
||||||
const IntensiveQuantities& intQuants,
|
|
||||||
DeferredLogger& deferred_logger) const
|
|
||||||
{
|
|
||||||
// TODO: the application of well efficiency factor has not been tested with an example yet
|
|
||||||
auto getFoamTransportIdx = [&deferred_logger] {
|
|
||||||
switch (FoamModule::transportPhase()) {
|
|
||||||
case Phase::WATER: {
|
|
||||||
return Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
|
|
||||||
}
|
|
||||||
case Phase::GAS: {
|
|
||||||
return Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
|
|
||||||
}
|
|
||||||
case Phase::SOLVENT: {
|
|
||||||
if constexpr (has_solvent)
|
|
||||||
return static_cast<unsigned>(Indices::contiSolventEqIdx);
|
|
||||||
else
|
|
||||||
OPM_DEFLOG_THROW(std::runtime_error, "Foam transport phase is SOLVENT but SOLVENT is not activated.", deferred_logger);
|
|
||||||
}
|
|
||||||
default: {
|
|
||||||
OPM_DEFLOG_THROW(std::runtime_error, "Foam transport phase must be GAS/WATER/SOLVENT.", deferred_logger);
|
|
||||||
}
|
|
||||||
}
|
|
||||||
};
|
|
||||||
EvalWell cq_s_foam = cq_s[getFoamTransportIdx()] * this->well_efficiency_factor_;
|
|
||||||
if (this->isInjector()) {
|
|
||||||
cq_s_foam *= this->wfoam();
|
|
||||||
} else {
|
|
||||||
cq_s_foam *= this->extendEval(intQuants.foamConcentration());
|
|
||||||
}
|
|
||||||
return Base::restrictEval(cq_s_foam);
|
|
||||||
}
|
|
||||||
|
|
||||||
|
|
||||||
template <typename TypeTag>
|
template <typename TypeTag>
|
||||||
std::tuple<typename StandardWell<TypeTag>::Eval,
|
std::tuple<typename StandardWell<TypeTag>::Eval,
|
||||||
typename StandardWell<TypeTag>::Eval,
|
typename StandardWell<TypeTag>::Eval,
|
||||||
|
Loading…
Reference in New Issue
Block a user