28#ifndef EWOMS_PVS_PRIMARY_VARIABLES_HH
29#define EWOMS_PVS_PRIMARY_VARIABLES_HH
34#include <opm/common/Exceptions.hpp>
39#include <opm/material/constraintsolvers/NcpFlash.hpp>
40#include <opm/material/fluidstates/CompositionalFluidState.hpp>
41#include <opm/material/common/Valgrind.hpp>
43#include <dune/common/fvector.hh>
58template <
class TypeTag>
73 enum { pressure0Idx = Indices::pressure0Idx };
74 enum { switch0Idx = Indices::switch0Idx };
76 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
77 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
78 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
80 using Toolbox = MathToolbox<Evaluation>;
81 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
83 using NcpFlash = Opm::NcpFlash<Scalar, FluidSystem>;
87 { Valgrind::SetDefined(*
this); }
94 Valgrind::CheckDefined(value);
95 Valgrind::SetDefined(*
this);
106 Valgrind::SetDefined(*
this);
108 phasePresence_ = value.phasePresence_;
114 template <
class Flu
idState>
116 const MaterialLawParams& matParams,
117 bool isInEquilibrium =
false)
121 for (
unsigned phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
122 assert(std::abs(fluidState.temperature(0) - fluidState.temperature(phaseIdx)) < 1e-30);
128 if (isInEquilibrium) {
135 typename FluidSystem::template ParameterCache<Scalar> paramCache;
136 CompositionalFluidState<Scalar, FluidSystem> fsFlash;
140 fsFlash.assign(fluidState);
143 paramCache.updateAll(fsFlash);
144 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
145 Scalar rho = FluidSystem::density(fsFlash, paramCache, phaseIdx);
146 fsFlash.setDensity(phaseIdx, rho);
149 ComponentVector globalMolarities(0.0);
150 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
151 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
152 globalMolarities[compIdx] +=
153 fsFlash.saturation(phaseIdx) * fsFlash.molarity(phaseIdx, compIdx);
158 NcpFlash::template solve<MaterialLaw>(fsFlash, matParams, paramCache, globalMolarities);
169 {
return phasePresence_; }
178 { phasePresence_ = value; }
220 {
return phasePresence_& (1 << phaseIdx); }
226 {
return static_cast<unsigned>(ffs(phasePresence_) - 1); }
234 phasePresence_ = value.phasePresence_;
263 unsigned varIdx = switch0Idx + phaseIdx - 1;
264 if (std::is_same<Evaluation, Scalar>::value)
265 return (*
this)[varIdx];
269 Toolbox::createConstant((*
this)[varIdx]);
270 return Toolbox::createVariable((*
this)[varIdx], varIdx);
277 template <
class Flu
idState>
280 using FsToolbox = MathToolbox<typename FluidState::Scalar>;
284 EnergyModule::setPriVarTemperatures(*
this, fluidState);
287 (*this)[pressure0Idx] = FsToolbox::value(fluidState.pressure(0));
288 Valgrind::CheckDefined((*
this)[pressure0Idx]);
292 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
296 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
297 a -= FsToolbox::value(fluidState.moleFraction(phaseIdx, compIdx));
299 Scalar b = FsToolbox::value(fluidState.saturation(phaseIdx));
302 phasePresence_ |= (1 << phaseIdx);
306 if (phasePresence_ == 0)
307 throw NumericalProblem(
"Phase state was 0, i.e., no fluid is present");
312 for (
unsigned switchIdx = 0; switchIdx < numPhases - 1; ++switchIdx) {
313 unsigned phaseIdx = switchIdx;
314 unsigned compIdx = switchIdx + 1;
315 if (switchIdx >= lowestPhaseIdx)
319 (*this)[switch0Idx + switchIdx] = FsToolbox::value(fluidState.saturation(phaseIdx));
320 Valgrind::CheckDefined((*
this)[switch0Idx + switchIdx]);
323 (*this)[switch0Idx + switchIdx] =
324 FsToolbox::value(fluidState.moleFraction(lowestPhaseIdx, compIdx));
325 Valgrind::CheckDefined((*
this)[switch0Idx + switchIdx]);
331 for (
unsigned compIdx = numPhases - 1; compIdx < numComponents - 1; ++compIdx) {
332 (*this)[switch0Idx + compIdx] =
333 FsToolbox::value(fluidState.moleFraction(lowestPhaseIdx, compIdx + 1));
334 Valgrind::CheckDefined((*
this)[switch0Idx + compIdx]);
341 void print(std::ostream& os = std::cout)
const
343 os <<
"(p_" << FluidSystem::phaseName(0) <<
" = "
344 << this->operator[](pressure0Idx);
346 for (
unsigned switchIdx = 0; switchIdx < numPhases - 1; ++switchIdx) {
347 unsigned phaseIdx = switchIdx;
348 unsigned compIdx = switchIdx + 1;
349 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;
365 os <<
", x_" << FluidSystem::phaseName(lowestPhaseIdx) <<
"^"
366 << FluidSystem::componentName(compIdx + 1) <<
" = "
367 << (*this)[switch0Idx + compIdx];
370 os <<
", phase presence: " <<
static_cast<int>(phasePresence_) << std::flush;
374 short phasePresence_;
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
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:60
void print(std::ostream &os=std::cout) const
Prints the names of the primary variables and their values.
Definition: pvsprimaryvariables.hh:341
unsigned implicitSaturationIdx() const
Returns the index of the phase with's its saturation is determined by the closure condition of satura...
Definition: pvsprimaryvariables.hh:199
void setPhasePresence(short value)
Set which fluid phases are present in a given control volume.
Definition: pvsprimaryvariables.hh:177
static bool phaseIsPresent(unsigned phaseIdx, short phasePresence)
Returns true iff a phase is present for a given phase presence.
Definition: pvsprimaryvariables.hh:210
void setPhasePresent(unsigned phaseIdx, bool yesno=true)
Set whether a given indivividual phase should be present or not.
Definition: pvsprimaryvariables.hh:187
ThisType & operator=(const Implementation &value)
Assignment operator from an other primary variables object.
Definition: pvsprimaryvariables.hh:231
void assignNaive(const FluidState &fluidState)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: pvsprimaryvariables.hh:278
PvsPrimaryVariables(const PvsPrimaryVariables &value)
Definition: pvsprimaryvariables.hh:104
PvsPrimaryVariables(Scalar value)
Constructor with assignment from scalar.
Definition: pvsprimaryvariables.hh:92
unsigned lowestPresentPhaseIdx() const
Returns the phase with the lowest index that is present.
Definition: pvsprimaryvariables.hh:225
ThisType & operator=(Scalar value)
Assignment operator from a scalar value.
Definition: pvsprimaryvariables.hh:242
Evaluation explicitSaturationValue(unsigned phaseIdx, unsigned timeIdx) const
Returns an explcitly stored saturation for a given phase.
Definition: pvsprimaryvariables.hh:257
bool phaseIsPresent(unsigned phaseIdx) const
Returns true iff a phase is present for the current phase presence.
Definition: pvsprimaryvariables.hh:219
short phasePresence() const
Return the fluid phases which are present in a given control volume.
Definition: pvsprimaryvariables.hh:168
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: pvsprimaryvariables.hh:115
PvsPrimaryVariables()
Definition: pvsprimaryvariables.hh:86
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 compositional multi-phase primary variable switching model.