28#ifndef OPM_PTFLASH_PRIMARY_VARIABLES_HH
29#define OPM_PTFLASH_PRIMARY_VARIABLES_HH
36#include <opm/material/constraintsolvers/NcpFlash.hpp>
37#include <opm/material/fluidstates/CompositionalFluidState.hpp>
38#include <opm/material/common/Valgrind.hpp>
40#include <dune/common/fvector.hh>
55template <
class TypeTag>
56class FlashPrimaryVariables :
public FvBasePrimaryVariables<TypeTag>
58 using ParentType = FvBasePrimaryVariables<TypeTag>;
60 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
61 using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
62 using MaterialLawParams = GetPropType<TypeTag, Properties::MaterialLawParams>;
63 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
64 using Indices = GetPropType<TypeTag, Properties::Indices>;
67 enum { z0Idx = Indices::z0Idx };
68 enum { pressure0Idx = Indices::pressure0Idx };
70 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
71 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
73 using Toolbox =
typename Opm::MathToolbox<Evaluation>;
78 { Opm::Valgrind::SetDefined(*
this); }
85 Opm::Valgrind::CheckDefined(value);
86 Opm::Valgrind::SetDefined(*
this);
99 template <
class Flu
idState>
101 const MaterialLawParams&,
113 template <
class Flu
idState>
121 EnergyModule::setPriVarTemperatures(*
this, fluidState);
124 Dune::FieldVector<Scalar, numComponents> z(0.0);
125 Scalar sumMoles = 0.0;
126 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
127 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
128 Scalar tmp = Opm::getValue(fluidState.molarity(phaseIdx, compIdx) * fluidState.saturation(phaseIdx));
129 z[compIdx] += Opm::max(tmp, 1e-8);
135 for (
int i = 0; i < numComponents - 1; ++i)
136 (*
this)[z0Idx + i] = z[i];
138 (*this)[pressure0Idx] = Opm::getValue(fluidState.pressure(0));
146 void print(std::ostream& os = std::cout)
const
148 os <<
"(p_" << FluidSystem::phaseName(0) <<
" = "
149 << this->operator[](pressure0Idx);
150 for (
unsigned compIdx = 0; compIdx < numComponents - 2; ++compIdx) {
151 os <<
", z_" << FluidSystem::componentName(compIdx) <<
" = "
152 << this->operator[](z0Idx + compIdx);
154 os <<
")" << std::flush;
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:50
Represents the primary variables used by the compositional flow model based on flash calculations.
Definition: flash/flashprimaryvariables.hh:57
void print(std::ostream &os=std::cout) const
Prints the names of the primary variables and their values.
Definition: ptflash/flashprimaryvariables.hh:146
void assignNaive(const FluidState &fluidState)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: flash/flashprimaryvariables.hh:113
FlashPrimaryVariables(Scalar value)
Constructor with assignment from scalar.
Definition: ptflash/flashprimaryvariables.hh:83
FlashPrimaryVariables()
Definition: ptflash/flashprimaryvariables.hh:77
FlashPrimaryVariables(const FlashPrimaryVariables &value)=default
void assignMassConservative(const FluidState &fluidState, const MaterialLawParams &, bool=false)
Set the primary variables from an arbitrary fluid state in a mass conservative way.
Definition: ptflash/flashprimaryvariables.hh:100
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:37