28#ifndef EWOMS_IMMISCIBLE_PRIMARY_VARIABLES_HH
29#define EWOMS_IMMISCIBLE_PRIMARY_VARIABLES_HH
31#include <dune/common/fvector.hh>
37#include <opm/material/common/Valgrind.hpp>
38#include <opm/material/constraintsolvers/ImmiscibleFlash.hpp>
39#include <opm/material/fluidstates/ImmiscibleFluidState.hpp>
52template <
class TypeTag>
67 enum { pressure0Idx = Indices::pressure0Idx };
68 enum { saturation0Idx = Indices::saturation0Idx };
70 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
71 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
73 using Toolbox =
typename Opm::MathToolbox<Evaluation>;
74 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
75 using ImmiscibleFlash = Opm::ImmiscibleFlash<Scalar, FluidSystem>;
83 { Opm::Valgrind::SetUndefined(*
this); }
99 using ParentType::operator=;
121 template <
class Flu
idState>
123 const MaterialLawParams& matParams,
124 bool isInEquilibrium =
false)
128 for (
unsigned phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
129 assert(std::abs(fluidState.temperature(0) - fluidState.temperature(phaseIdx)) < 1e-30);
135 if (isInEquilibrium) {
142 typename FluidSystem::template ParameterCache<Scalar> paramCache;
143 Opm::ImmiscibleFluidState<Scalar, FluidSystem> fsFlash;
147 fsFlash.assign(fluidState);
150 paramCache.updateAll(fsFlash);
151 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
152 Scalar rho = FluidSystem::density(fsFlash, paramCache, phaseIdx);
153 fsFlash.setDensity(phaseIdx, rho);
157 ComponentVector globalMolarities(0.0);
158 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
159 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
160 globalMolarities[compIdx] +=
161 fsFlash.saturation(phaseIdx) * fsFlash.molarity(phaseIdx, compIdx);
166 ImmiscibleFlash::template solve<MaterialLaw>(fsFlash, matParams, paramCache, globalMolarities);
188 template <
class Flu
idState>
193 EnergyModule::setPriVarTemperatures(asImp_(), fluidState);
195 (*this)[pressure0Idx] = fluidState.pressure(0);
196 for (
unsigned phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
197 (*
this)[saturation0Idx + phaseIdx] = fluidState.saturation(phaseIdx);
201 Implementation& asImp_()
202 {
return *
static_cast<Implementation *
>(
this); }
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 immiscible multi-phase, model.
Definition: immiscibleprimaryvariables.hh:54
ImmisciblePrimaryVariables & operator=(const ImmisciblePrimaryVariables &value)=default
Assignment operator.
void assignNaive(const FluidState &fluidState)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: immiscibleprimaryvariables.hh:189
ImmisciblePrimaryVariables(const ImmisciblePrimaryVariables &value)=default
Copy constructor.
ImmisciblePrimaryVariables()
Default constructor.
Definition: immiscibleprimaryvariables.hh:82
void assignMassConservative(const FluidState &fluidState, const MaterialLawParams &matParams, bool isInEquilibrium=false)
< Import base class assignment operators.
Definition: immiscibleprimaryvariables.hh:122
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Defines the properties required for the immiscible multi-phase model.
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