28#ifndef OPM_PTFLASH_PRIMARY_VARIABLES_HH
29#define OPM_PTFLASH_PRIMARY_VARIABLES_HH
31#include <opm/material/common/MathToolbox.hpp>
32#include <opm/material/common/Valgrind.hpp>
51template <
class TypeTag>
52class FlashPrimaryVariables :
public FvBasePrimaryVariables<TypeTag>
54 using ParentType = FvBasePrimaryVariables<TypeTag>;
56 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
57 using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
58 using MaterialLawParams = GetPropType<TypeTag, Properties::MaterialLawParams>;
59 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
60 using Indices = GetPropType<TypeTag, Properties::Indices>;
63 enum { z0Idx = Indices::z0Idx };
64 enum { pressure0Idx = Indices::pressure0Idx };
65 enum { water0Idx = Indices::water0Idx };
67 static constexpr bool waterEnabled = Indices::waterEnabled;
70 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
71 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
72 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
74 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
75 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
77 using Toolbox = MathToolbox<Evaluation>;
82 { Opm::Valgrind::SetDefined(*
this); }
91 using ParentType::operator=;
96 template <
class Flu
idState>
98 const MaterialLawParams&,
110 template <
class Flu
idState>
118 EnergyModule::setPriVarTemperatures(*
this, fluidState);
121 for (
int i = 0; i < numComponents - 1; ++i)
122 (*
this)[z0Idx + i] = getValue(fluidState.moleFraction(i));
125 (*this)[pressure0Idx] = getValue(fluidState.pressure(oilPhaseIdx));
128 if constexpr (waterEnabled) {
129 (*this)[water0Idx] = getValue(fluidState.saturation(waterPhaseIdx));
140 os <<
"(p_" << FluidSystem::phaseName(FluidSystem::oilPhaseIdx) <<
" = "
141 << (*this)[pressure0Idx];
142 for (
unsigned compIdx = 0; compIdx < numComponents - 2; ++compIdx) {
143 os <<
", z_" << FluidSystem::componentName(compIdx) <<
" = "
144 << (*this)[z0Idx + compIdx];
146 if constexpr (waterEnabled) {
147 os <<
", S_w = " << (*this)[water0Idx];
149 os <<
")" << std::flush;
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:54
Represents the primary variables used by the compositional flow model based on flash calculations.
Definition: flash/flashprimaryvariables.hh:52
void assignNaive(const FluidState &fluidState)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: flash/flashprimaryvariables.hh:106
FlashPrimaryVariables()
Definition: ptflash/flashprimaryvariables.hh:81
FlashPrimaryVariables(const FlashPrimaryVariables &value)=default
void print(std::ostream &os) const
Prints the names of the primary variables and their values.
Definition: ptflash/flashprimaryvariables.hh:138
void assignMassConservative(const FluidState &fluidState, const MaterialLawParams &, bool=false)
< Import base class assignment operators.
Definition: ptflash/flashprimaryvariables.hh:97
FlashPrimaryVariables & operator=(const FlashPrimaryVariables &value)=default
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Definition: blackoilboundaryratevector.hh:39