28#ifndef EWOMS_NCP_PRIMARY_VARIABLES_HH 
   29#define EWOMS_NCP_PRIMARY_VARIABLES_HH 
   36#include <opm/material/constraintsolvers/NcpFlash.hpp> 
   37#include <opm/material/fluidstates/CompositionalFluidState.hpp> 
   38#include <opm/material/densead/Math.hpp> 
   40#include <dune/common/fvector.hh> 
   53template <
class TypeTag>
 
   65    enum { pressure0Idx = Indices::pressure0Idx };
 
   66    enum { saturation0Idx = Indices::saturation0Idx };
 
   67    enum { fugacity0Idx = Indices::fugacity0Idx };
 
   69    enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
 
   70    enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
 
   71    using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
 
   73    enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
 
   76    using NcpFlash = Opm::NcpFlash<Scalar, FluidSystem>;
 
   77    using Toolbox = Opm::MathToolbox<Evaluation>;
 
   99    template <
class Flu
idState>
 
  101                                const MaterialLawParams& matParams,
 
  102                                bool isInEquilibrium = 
false)
 
  104        using FsToolbox = Opm::MathToolbox<typename FluidState::Scalar>;
 
  108        for (
unsigned phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
 
  109            assert(fluidState.temperature(0) == fluidState.temperature(phaseIdx));
 
  115        if (isInEquilibrium) {
 
  122        typename FluidSystem::template ParameterCache<Scalar> paramCache;
 
  123        Opm::CompositionalFluidState<Scalar, FluidSystem> fsFlash;
 
  127        fsFlash.assign(fluidState);
 
  130        paramCache.updateAll(fsFlash);
 
  131        for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
 
  132            Scalar rho = FluidSystem::density(fsFlash, paramCache, phaseIdx);
 
  133            fsFlash.setDensity(phaseIdx, rho);
 
  137        ComponentVector globalMolarities(0.0);
 
  138        for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
 
  139            for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
 
  140                globalMolarities[compIdx] +=
 
  141                    FsToolbox::value(fsFlash.saturation(phaseIdx))
 
  142                    * FsToolbox::value(fsFlash.molarity(phaseIdx, compIdx));
 
  147        NcpFlash::template solve<MaterialLaw>(fsFlash, matParams, paramCache, globalMolarities);
 
  156    template <
class Flu
idState>
 
  157    void assignNaive(
const FluidState& fluidState, 
unsigned refPhaseIdx = 0)
 
  159        using FsToolbox = Opm::MathToolbox<typename FluidState::Scalar>;
 
  163        EnergyModule::setPriVarTemperatures(*
this, fluidState);
 
  166        typename FluidSystem::template ParameterCache<Scalar> paramCache;
 
  167        paramCache.updatePhase(fluidState, refPhaseIdx);
 
  168        Scalar pRef = FsToolbox::value(fluidState.pressure(refPhaseIdx));
 
  169        for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
 
  173                FluidSystem::template fugacityCoefficient<FluidState, Scalar>(fluidState,
 
  177            (*this)[fugacity0Idx + compIdx] =
 
  178                fugCoeff*fluidState.moleFraction(refPhaseIdx, compIdx)*pRef;
 
  182        (*this)[pressure0Idx] = FsToolbox::value(fluidState.pressure(0));
 
  185        for (
unsigned phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
 
  186            (*
this)[saturation0Idx + phaseIdx] = FsToolbox::value(fluidState.saturation(phaseIdx));
 
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 compositional multi-phase NCP model.
Definition: ncpprimaryvariables.hh:55
NcpPrimaryVariables()
Definition: ncpprimaryvariables.hh:80
void assignNaive(const FluidState &fluidState, unsigned refPhaseIdx=0)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: ncpprimaryvariables.hh:157
NcpPrimaryVariables & operator=(const NcpPrimaryVariables &value)=default
NcpPrimaryVariables(Scalar value)
Constructor with assignment from scalar.
Definition: ncpprimaryvariables.hh:86
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: ncpprimaryvariables.hh:100
NcpPrimaryVariables(const NcpPrimaryVariables &value)=default
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
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
Declares the properties required for the NCP compositional multi-phase model.