28#ifndef EWOMS_IMMISCIBLE_PRIMARY_VARIABLES_HH
29#define EWOMS_IMMISCIBLE_PRIMARY_VARIABLES_HH
36#include <opm/material/constraintsolvers/ImmiscibleFlash.hpp>
37#include <opm/material/fluidstates/ImmiscibleFluidState.hpp>
38#include <opm/material/common/Valgrind.hpp>
40#include <dune/common/fvector.hh>
53template <
class TypeTag>
68 enum { pressure0Idx = Indices::pressure0Idx };
69 enum { saturation0Idx = Indices::saturation0Idx };
71 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
72 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
74 using Toolbox =
typename Opm::MathToolbox<Evaluation>;
75 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
76 using ImmiscibleFlash = Opm::ImmiscibleFlash<Scalar, FluidSystem>;
84 { Opm::Valgrind::SetUndefined(*
this); }
128 template <
class Flu
idState>
130 const MaterialLawParams& matParams,
131 bool isInEquilibrium =
false)
135 for (
unsigned phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
136 assert(std::abs(fluidState.temperature(0) - fluidState.temperature(phaseIdx)) < 1e-30);
142 if (isInEquilibrium) {
149 typename FluidSystem::template ParameterCache<Scalar> paramCache;
150 Opm::ImmiscibleFluidState<Scalar, FluidSystem> fsFlash;
154 fsFlash.assign(fluidState);
157 paramCache.updateAll(fsFlash);
158 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
159 Scalar rho = FluidSystem::density(fsFlash, paramCache, phaseIdx);
160 fsFlash.setDensity(phaseIdx, rho);
164 ComponentVector globalMolarities(0.0);
165 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
166 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
167 globalMolarities[compIdx] +=
168 fsFlash.saturation(phaseIdx) * fsFlash.molarity(phaseIdx, compIdx);
173 ImmiscibleFlash::template solve<MaterialLaw>(fsFlash, matParams, paramCache, globalMolarities);
195 template <
class Flu
idState>
200 EnergyModule::setPriVarTemperatures(asImp_(), fluidState);
202 (*this)[pressure0Idx] = fluidState.pressure(0);
203 for (
unsigned phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
204 (*
this)[saturation0Idx + phaseIdx] = fluidState.saturation(phaseIdx);
208 Implementation& asImp_()
209 {
return *
static_cast<Implementation *
>(
this); }
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:50
Represents the primary variables used by the a model.
Definition: fvbaseprimaryvariables.hh:52
Represents the primary variables used by the immiscible multi-phase, model.
Definition: immiscibleprimaryvariables.hh:55
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:196
ImmisciblePrimaryVariables(const ImmisciblePrimaryVariables &value)=default
Copy constructor.
ImmisciblePrimaryVariables()
Default constructor.
Definition: immiscibleprimaryvariables.hh:83
ImmisciblePrimaryVariables(Scalar value)
Constructor with assignment from scalar.
Definition: immiscibleprimaryvariables.hh:91
void assignMassConservative(const FluidState &fluidState, const MaterialLawParams &matParams, bool isInEquilibrium=false)
Set the primary variables from an arbitrary fluid state in a mass conservative way.
Definition: immiscibleprimaryvariables.hh:129
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:37
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:235