28#ifndef EWOMS_PVS_PRIMARY_VARIABLES_HH
29#define EWOMS_PVS_PRIMARY_VARIABLES_HH
31#include <dune/common/fvector.hh>
33#include <opm/common/Exceptions.hpp>
35#include <opm/material/common/MathToolbox.hpp>
36#include <opm/material/common/Valgrind.hpp>
37#include <opm/material/constraintsolvers/NcpFlash.hpp>
38#include <opm/material/fluidstates/CompositionalFluidState.hpp>
60template <
class TypeTag>
75 enum { pressure0Idx = Indices::pressure0Idx };
76 enum { switch0Idx = Indices::switch0Idx };
78 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
79 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
80 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
82 using Toolbox = MathToolbox<Evaluation>;
83 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
85 using NcpFlash = ::Opm::NcpFlash<Scalar, FluidSystem>;
89 { Valgrind::SetDefined(*
this); }
97 Valgrind::SetDefined(*
this);
99 phasePresence_ = value.phasePresence_;
105 template <
class Flu
idState>
107 const MaterialLawParams& matParams,
108 bool isInEquilibrium =
false)
112 for (
unsigned phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
113 assert(std::abs(fluidState.temperature(0) - fluidState.temperature(phaseIdx)) < 1e-30);
119 if (isInEquilibrium) {
126 typename FluidSystem::template ParameterCache<Scalar> paramCache;
127 CompositionalFluidState<Scalar, FluidSystem> fsFlash;
131 fsFlash.assign(fluidState);
134 paramCache.updateAll(fsFlash);
135 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
136 const Scalar rho = FluidSystem::density(fsFlash, paramCache, phaseIdx);
137 fsFlash.setDensity(phaseIdx, rho);
140 ComponentVector globalMolarities(0.0);
141 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
142 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
143 globalMolarities[compIdx] +=
144 fsFlash.saturation(phaseIdx) * fsFlash.molarity(phaseIdx, compIdx);
149 NcpFlash::template solve<MaterialLaw>(fsFlash, matParams, paramCache, globalMolarities);
160 {
return phasePresence_; }
169 { phasePresence_ = value; }
213 {
return phasePresence_ & (1 << phaseIdx); }
219 {
return static_cast<unsigned>(ffs(phasePresence_) - 1); }
227 phasePresence_ = value.phasePresence_;
257 const unsigned varIdx = switch0Idx + phaseIdx - 1;
258 if constexpr (std::is_same_v<Evaluation, Scalar>) {
259 return (*
this)[varIdx];
264 Toolbox::createConstant((*
this)[varIdx]);
266 return Toolbox::createVariable((*
this)[varIdx], varIdx);
273 template <
class Flu
idState>
276 using FsToolbox = MathToolbox<typename FluidState::Scalar>;
280 EnergyModule::setPriVarTemperatures(*
this, fluidState);
283 (*this)[pressure0Idx] = FsToolbox::value(fluidState.pressure(0));
284 Valgrind::CheckDefined((*
this)[pressure0Idx]);
288 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
292 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
293 a -= FsToolbox::value(fluidState.moleFraction(phaseIdx, compIdx));
295 const Scalar b = FsToolbox::value(fluidState.saturation(phaseIdx));
298 phasePresence_ |= (1 << phaseIdx);
303 if (phasePresence_ == 0) {
304 throw NumericalProblem(
"Phase state was 0, i.e., no fluid is present");
310 for (
unsigned switchIdx = 0; switchIdx < numPhases - 1; ++switchIdx) {
311 unsigned phaseIdx = switchIdx;
312 const unsigned compIdx = switchIdx + 1;
313 if (switchIdx >= lowestPhaseIdx) {
318 (*this)[switch0Idx + switchIdx] = FsToolbox::value(fluidState.saturation(phaseIdx));
319 Valgrind::CheckDefined((*
this)[switch0Idx + switchIdx]);
322 (*this)[switch0Idx + switchIdx] =
323 FsToolbox::value(fluidState.moleFraction(lowestPhaseIdx, compIdx));
324 Valgrind::CheckDefined((*
this)[switch0Idx + switchIdx]);
330 for (
unsigned compIdx = numPhases - 1; compIdx < numComponents - 1; ++compIdx) {
331 (*this)[switch0Idx + compIdx] =
332 FsToolbox::value(fluidState.moleFraction(lowestPhaseIdx, compIdx + 1));
333 Valgrind::CheckDefined((*
this)[switch0Idx + compIdx]);
342 os <<
"(p_" << FluidSystem::phaseName(0) <<
" = "
343 << (*this)[pressure0Idx];
345 for (
unsigned switchIdx = 0; switchIdx < numPhases - 1; ++switchIdx) {
346 unsigned phaseIdx = switchIdx;
347 const unsigned compIdx = switchIdx + 1;
348 if (phaseIdx >= lowestPhaseIdx) {
354 os <<
", S_" << FluidSystem::phaseName(phaseIdx) <<
" = "
355 << (*this)[switch0Idx + switchIdx];
358 os <<
", x_" << FluidSystem::phaseName(lowestPhaseIdx) <<
"^"
359 << FluidSystem::componentName(compIdx) <<
" = "
360 << (*this)[switch0Idx + switchIdx];
363 for (
unsigned compIdx = numPhases - 1; compIdx < numComponents - 1; ++compIdx) {
364 os <<
", x_" << FluidSystem::phaseName(lowestPhaseIdx) <<
"^"
365 << FluidSystem::componentName(compIdx + 1) <<
" = "
366 << (*this)[switch0Idx + compIdx];
369 os <<
", phase presence: " <<
static_cast<int>(phasePresence_) << std::flush;
373 short phasePresence_{};
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
FvBasePrimaryVariables & operator=(const FvBasePrimaryVariables &value)=default
Assignment from another primary variables object.
Represents the primary variables used in the primary variable switching compositional model.
Definition: pvsprimaryvariables.hh:62
unsigned implicitSaturationIdx() const
Returns the index of the phase with's its saturation is determined by the closure condition of satura...
Definition: pvsprimaryvariables.hh:192
void setPhasePresence(short value)
Set which fluid phases are present in a given control volume.
Definition: pvsprimaryvariables.hh:168
static bool phaseIsPresent(unsigned phaseIdx, short phasePresence)
Returns true iff a phase is present for a given phase presence.
Definition: pvsprimaryvariables.hh:203
ThisType & operator=(const Implementation &value)
Assignment operator from an other primary variables object.
Definition: pvsprimaryvariables.hh:224
void assignNaive(const FluidState &fluidState)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: pvsprimaryvariables.hh:274
PvsPrimaryVariables(const PvsPrimaryVariables &value)
Definition: pvsprimaryvariables.hh:95
unsigned lowestPresentPhaseIdx() const
Returns the phase with the lowest index that is present.
Definition: pvsprimaryvariables.hh:218
ThisType & operator=(Scalar value)
Assignment operator from a scalar value.
Definition: pvsprimaryvariables.hh:235
void print(std::ostream &os) const
Prints the names of the primary variables and their values.
Definition: pvsprimaryvariables.hh:340
Evaluation explicitSaturationValue(unsigned phaseIdx, unsigned timeIdx) const
Returns an explcitly stored saturation for a given phase.
Definition: pvsprimaryvariables.hh:250
bool phaseIsPresent(unsigned phaseIdx) const
Returns true iff a phase is present for the current phase presence.
Definition: pvsprimaryvariables.hh:212
short phasePresence() const
Return the fluid phases which are present in a given control volume.
Definition: pvsprimaryvariables.hh:159
void assignMassConservative(const FluidState &fluidState, const MaterialLawParams &matParams, bool isInEquilibrium=false)
< Import base class assignment operators.
Definition: pvsprimaryvariables.hh:106
void setPhasePresent(unsigned phaseIdx, bool yesno)
Set whether a given indivividual phase should be present or not.
Definition: pvsprimaryvariables.hh:178
PvsPrimaryVariables()
Definition: pvsprimaryvariables.hh:88
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 compositional multi-phase primary variable switching model.