28#ifndef EWOMS_NCP_PRIMARY_VARIABLES_HH
29#define EWOMS_NCP_PRIMARY_VARIABLES_HH
31#include <dune/common/fvector.hh>
33#include <opm/material/constraintsolvers/NcpFlash.hpp>
34#include <opm/material/densead/Math.hpp>
35#include <opm/material/fluidstates/CompositionalFluidState.hpp>
54template <
class TypeTag>
66 enum { pressure0Idx = Indices::pressure0Idx };
67 enum { saturation0Idx = Indices::saturation0Idx };
68 enum { fugacity0Idx = Indices::fugacity0Idx };
70 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
71 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
72 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
74 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
77 using NcpFlash = ::Opm::NcpFlash<Scalar, FluidSystem>;
78 using Toolbox = MathToolbox<Evaluation>;
90 using ParentType::operator=;
95 template <
class Flu
idState>
97 const MaterialLawParams& matParams,
98 bool isInEquilibrium =
false)
100 using FsToolbox = MathToolbox<typename FluidState::Scalar>;
104 for (
unsigned phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
105 assert(fluidState.temperature(0) == fluidState.temperature(phaseIdx));
111 if (isInEquilibrium) {
118 typename FluidSystem::template ParameterCache<Scalar> paramCache;
119 CompositionalFluidState<Scalar, FluidSystem> fsFlash;
123 fsFlash.assign(fluidState);
126 paramCache.updateAll(fsFlash);
127 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
128 const Scalar rho = FluidSystem::density(fsFlash, paramCache, phaseIdx);
129 fsFlash.setDensity(phaseIdx, rho);
133 ComponentVector globalMolarities(0.0);
134 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
135 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
136 globalMolarities[compIdx] +=
137 FsToolbox::value(fsFlash.saturation(phaseIdx)) *
138 FsToolbox::value(fsFlash.molarity(phaseIdx, compIdx));
143 NcpFlash::template solve<MaterialLaw>(fsFlash, matParams, paramCache, globalMolarities);
152 template <
class Flu
idState>
153 void assignNaive(
const FluidState& fluidState,
unsigned refPhaseIdx = 0)
155 using FsToolbox = MathToolbox<typename FluidState::Scalar>;
159 EnergyModule::setPriVarTemperatures(*
this, fluidState);
162 typename FluidSystem::template ParameterCache<Scalar> paramCache;
163 paramCache.updatePhase(fluidState, refPhaseIdx);
164 const Scalar pRef = FsToolbox::value(fluidState.pressure(refPhaseIdx));
165 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
168 const Scalar fugCoeff =
169 FluidSystem::template fugacityCoefficient<FluidState, Scalar>(fluidState,
173 (*this)[fugacity0Idx + compIdx] =
174 fugCoeff * fluidState.moleFraction(refPhaseIdx, compIdx) * pRef;
178 (*this)[pressure0Idx] = FsToolbox::value(fluidState.pressure(0));
181 for (
unsigned phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx) {
182 (*this)[saturation0Idx + phaseIdx] = FsToolbox::value(fluidState.saturation(phaseIdx));
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:54
Represents the primary variables used by the a model.
Definition: fvbaseprimaryvariables.hh:53
Represents the primary variables used by the compositional multi-phase NCP model.
Definition: ncpprimaryvariables.hh:56
void assignNaive(const FluidState &fluidState, unsigned refPhaseIdx=0)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: ncpprimaryvariables.hh:153
NcpPrimaryVariables & operator=(const NcpPrimaryVariables &value)=default
NcpPrimaryVariables()=default
void assignMassConservative(const FluidState &fluidState, const MaterialLawParams &matParams, bool isInEquilibrium=false)
< Import base class assignment operators.
Definition: ncpprimaryvariables.hh:96
NcpPrimaryVariables(const NcpPrimaryVariables &value)=default
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Definition: blackoilboundaryratevector.hh:39
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
Declares the properties required for the NCP compositional multi-phase model.