StandardWellPrimaryVariables: remove unnecessary Scalar template parameter

use the Scalar type from the FluidSystem
This commit is contained in:
Arne Morten Kvarving 2024-02-22 15:17:09 +01:00
parent 9997cde07a
commit 5affbf4bd5
5 changed files with 41 additions and 40 deletions

View File

@ -32,7 +32,7 @@ class DeferredLogger;
class GroupState; class GroupState;
class Schedule; class Schedule;
template<class Scalar, int numEq> class StandardWellEquations; template<class Scalar, int numEq> class StandardWellEquations;
template<class FluidSystem, class Indices, class Scalar> class StandardWellPrimaryVariables; template<class FluidSystem, class Indices> class StandardWellPrimaryVariables;
class SummaryState; class SummaryState;
template<class FluidSystem> class WellInterfaceFluidSystem; template<class FluidSystem> class WellInterfaceFluidSystem;
class WellState; class WellState;
@ -43,7 +43,7 @@ class StandardWellAssemble
{ {
public: public:
using Scalar = typename FluidSystem::Scalar; using Scalar = typename FluidSystem::Scalar;
using PrimaryVariables = StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>; using PrimaryVariables = StandardWellPrimaryVariables<FluidSystem,Indices>;
using EvalWell = typename PrimaryVariables::EvalWell; using EvalWell = typename PrimaryVariables::EvalWell;
//! \brief Constructor initializes reference to well. //! \brief Constructor initializes reference to well.

View File

@ -94,7 +94,7 @@ public:
{ return perf_pressure_diffs_[perf]; } { return perf_pressure_diffs_[perf]; }
using Eval = typename WellInterfaceIndices<FluidSystem,Indices,Scalar>::Eval; using Eval = typename WellInterfaceIndices<FluidSystem,Indices,Scalar>::Eval;
using EvalWell = typename StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>::EvalWell; using EvalWell = typename StandardWellPrimaryVariables<FluidSystem,Indices>::EvalWell;
Eval connectionRateBrine(double& rate, Eval connectionRateBrine(double& rate,
const double vap_wat_rate, const double vap_wat_rate,

View File

@ -48,7 +48,7 @@ class StandardWellEval
{ {
protected: protected:
using Scalar = typename FluidSystem::Scalar; using Scalar = typename FluidSystem::Scalar;
using PrimaryVariables = StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>; using PrimaryVariables = StandardWellPrimaryVariables<FluidSystem,Indices>;
using StdWellConnections = StandardWellConnections<FluidSystem,Indices>; using StdWellConnections = StandardWellConnections<FluidSystem,Indices>;
static constexpr int Bhp = PrimaryVariables::Bhp; static constexpr int Bhp = PrimaryVariables::Bhp;
static constexpr int WQTotal= PrimaryVariables::WQTotal; static constexpr int WQTotal= PrimaryVariables::WQTotal;

View File

@ -95,8 +95,8 @@ Scalar relaxationFactorFraction(const Scalar old_value,
namespace Opm { namespace Opm {
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
void StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>:: void StandardWellPrimaryVariables<FluidSystem,Indices>::
init() init()
{ {
for (int eqIdx = 0; eqIdx < numWellEq_; ++eqIdx) { for (int eqIdx = 0; eqIdx < numWellEq_; ++eqIdx) {
@ -108,8 +108,8 @@ init()
} }
} }
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
void StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>:: void StandardWellPrimaryVariables<FluidSystem,Indices>::
resize(const int numWellEq) resize(const int numWellEq)
{ {
value_.resize(numWellEq, 0.0); value_.resize(numWellEq, 0.0);
@ -117,8 +117,8 @@ resize(const int numWellEq)
numWellEq_ = numWellEq; numWellEq_ = numWellEq;
} }
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
void StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>:: void StandardWellPrimaryVariables<FluidSystem,Indices>::
update(const WellState& well_state, update(const WellState& well_state,
const bool stop_or_zero_rate_target, const bool stop_or_zero_rate_target,
DeferredLogger& deferred_logger) DeferredLogger& deferred_logger)
@ -225,8 +225,8 @@ update(const WellState& well_state,
value_[Bhp] = ws.bhp; value_[Bhp] = ws.bhp;
} }
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
void StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>:: void StandardWellPrimaryVariables<FluidSystem,Indices>::
updatePolyMW(const WellState& well_state) updatePolyMW(const WellState& well_state)
{ {
if (well_.isInjector()) { if (well_.isInjector()) {
@ -241,8 +241,8 @@ updatePolyMW(const WellState& well_state)
} }
} }
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
void StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>:: void StandardWellPrimaryVariables<FluidSystem,Indices>::
updateNewton(const BVectorWell& dwells, updateNewton(const BVectorWell& dwells,
const bool stop_or_zero_rate_target, const bool stop_or_zero_rate_target,
const double dFLimit, const double dFLimit,
@ -303,8 +303,8 @@ updateNewton(const BVectorWell& dwells,
value_[Bhp] = std::max(value_[Bhp] - dx1_limited, bhp_lower_limit); value_[Bhp] = std::max(value_[Bhp] - dx1_limited, bhp_lower_limit);
} }
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
void StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>:: void StandardWellPrimaryVariables<FluidSystem,Indices>::
updateNewtonPolyMW(const BVectorWell& dwells) updateNewtonPolyMW(const BVectorWell& dwells)
{ {
if (well_.isInjector()) { if (well_.isInjector()) {
@ -322,8 +322,8 @@ updateNewtonPolyMW(const BVectorWell& dwells)
} }
} }
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
void StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>:: void StandardWellPrimaryVariables<FluidSystem,Indices>::
copyToWellState(WellState& well_state, copyToWellState(WellState& well_state,
DeferredLogger& deferred_logger) const DeferredLogger& deferred_logger) const
{ {
@ -427,8 +427,8 @@ copyToWellState(WellState& well_state,
} }
} }
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
void StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>:: void StandardWellPrimaryVariables<FluidSystem,Indices>::
copyToWellStatePolyMW(WellState& well_state) const copyToWellStatePolyMW(WellState& well_state) const
{ {
if (well_.isInjector()) { if (well_.isInjector()) {
@ -443,9 +443,9 @@ copyToWellStatePolyMW(WellState& well_state) const
} }
} }
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
typename StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>::EvalWell typename StandardWellPrimaryVariables<FluidSystem,Indices>::EvalWell
StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>:: StandardWellPrimaryVariables<FluidSystem,Indices>::
volumeFraction(const unsigned compIdx) const volumeFraction(const unsigned compIdx) const
{ {
if (FluidSystem::numActivePhases() == 1) { if (FluidSystem::numActivePhases() == 1) {
@ -480,9 +480,9 @@ volumeFraction(const unsigned compIdx) const
return well_fraction; return well_fraction;
} }
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
typename StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>::EvalWell typename StandardWellPrimaryVariables<FluidSystem,Indices>::EvalWell
StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>:: StandardWellPrimaryVariables<FluidSystem,Indices>::
volumeFractionScaled(const int compIdx) const volumeFractionScaled(const int compIdx) const
{ {
const int legacyCompIdx = well_.modelCompIdxToFlowCompIdx(compIdx); const int legacyCompIdx = well_.modelCompIdxToFlowCompIdx(compIdx);
@ -494,9 +494,9 @@ volumeFractionScaled(const int compIdx) const
return this->volumeFraction(compIdx); return this->volumeFraction(compIdx);
} }
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
typename StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>::EvalWell typename StandardWellPrimaryVariables<FluidSystem,Indices>::EvalWell
StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>:: StandardWellPrimaryVariables<FluidSystem,Indices>::
surfaceVolumeFraction(const int compIdx) const surfaceVolumeFraction(const int compIdx) const
{ {
EvalWell sum_volume_fraction_scaled(numWellEq_ + Indices::numEq, 0.); EvalWell sum_volume_fraction_scaled(numWellEq_ + Indices::numEq, 0.);
@ -509,9 +509,9 @@ surfaceVolumeFraction(const int compIdx) const
return this->volumeFractionScaled(compIdx) / sum_volume_fraction_scaled; return this->volumeFractionScaled(compIdx) / sum_volume_fraction_scaled;
} }
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
typename StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>::EvalWell typename StandardWellPrimaryVariables<FluidSystem,Indices>::EvalWell
StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>:: StandardWellPrimaryVariables<FluidSystem,Indices>::
getQs(const int comp_idx) const getQs(const int comp_idx) const
{ {
// Note: currently, the WQTotal definition is still depends on Injector/Producer. // Note: currently, the WQTotal definition is still depends on Injector/Producer.
@ -556,8 +556,8 @@ getQs(const int comp_idx) const
} }
} }
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
void StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>:: void StandardWellPrimaryVariables<FluidSystem,Indices>::
processFractions() processFractions()
{ {
static constexpr int Water = BlackoilPhases::Aqua; static constexpr int Water = BlackoilPhases::Aqua;
@ -659,8 +659,8 @@ processFractions()
} }
} }
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
double StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>:: double StandardWellPrimaryVariables<FluidSystem,Indices>::
relaxationFactorFractionsProducer(const BVectorWell& dwells, DeferredLogger& deferred_logger) const relaxationFactorFractionsProducer(const BVectorWell& dwells, DeferredLogger& deferred_logger) const
{ {
// TODO: not considering solvent yet // TODO: not considering solvent yet
@ -713,8 +713,8 @@ relaxationFactorFractionsProducer(const BVectorWell& dwells, DeferredLogger& def
return relaxation_factor; return relaxation_factor;
} }
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
void StandardWellPrimaryVariables<FluidSystem,Indices,Scalar>:: void StandardWellPrimaryVariables<FluidSystem,Indices>::
checkFinite(DeferredLogger& deferred_logger) const checkFinite(DeferredLogger& deferred_logger) const
{ {
for (const Scalar v : value_) { for (const Scalar v : value_) {
@ -726,7 +726,7 @@ checkFinite(DeferredLogger& deferred_logger) const
} }
#define INSTANCE(...) \ #define INSTANCE(...) \
template class StandardWellPrimaryVariables<BlackOilFluidSystem<double,BlackOilDefaultIndexTraits>,__VA_ARGS__,double>; template class StandardWellPrimaryVariables<BlackOilFluidSystem<double,BlackOilDefaultIndexTraits>,__VA_ARGS__>;
// 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>)

View File

@ -37,7 +37,7 @@ template<class FluidSystem, class Indices, class Scalar> class WellInterfaceIndi
class WellState; class WellState;
//! \brief Class holding primary variables for StandardWell. //! \brief Class holding primary variables for StandardWell.
template<class FluidSystem, class Indices, class Scalar> template<class FluidSystem, class Indices>
class StandardWellPrimaryVariables { class StandardWellPrimaryVariables {
protected: protected:
// the positions of the primary variables for StandardWell // the positions of the primary variables for StandardWell
@ -82,6 +82,7 @@ public:
static constexpr int GFrac = has_gfrac_variable ? has_wfrac_variable + 1 : -1000; static constexpr int GFrac = has_gfrac_variable ? has_wfrac_variable + 1 : -1000;
static constexpr int SFrac = !Indices::enableSolvent ? -1000 : has_wfrac_variable+has_gfrac_variable+1; static constexpr int SFrac = !Indices::enableSolvent ? -1000 : has_wfrac_variable+has_gfrac_variable+1;
using Scalar = typename FluidSystem::Scalar;
//! \brief Evaluation for the well equations. //! \brief Evaluation for the well equations.
using EvalWell = DenseAd::DynamicEvaluation<Scalar, numStaticWellEq + Indices::numEq + 1>; using EvalWell = DenseAd::DynamicEvaluation<Scalar, numStaticWellEq + Indices::numEq + 1>;
using BVectorWell = typename StandardWellEquations<Scalar,Indices::numEq>::BVectorWell; using BVectorWell = typename StandardWellEquations<Scalar,Indices::numEq>::BVectorWell;