ptflash/flashprimaryvariables.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef OPM_PTFLASH_PRIMARY_VARIABLES_HH
29#define OPM_PTFLASH_PRIMARY_VARIABLES_HH
30
31#include <opm/material/common/MathToolbox.hpp>
32#include <opm/material/common/Valgrind.hpp>
33
37
38#include <ostream>
39
40namespace Opm {
41
51template <class TypeTag>
52class FlashPrimaryVariables : public FvBasePrimaryVariables<TypeTag>
53{
54 using ParentType = FvBasePrimaryVariables<TypeTag>;
55
56 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
57 using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
58 using MaterialLawParams = GetPropType<TypeTag, Properties::MaterialLawParams>;
59 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
60 using Indices = GetPropType<TypeTag, Properties::Indices>;
61
62 // primary variable indices
63 enum { z0Idx = Indices::z0Idx };
64 enum { pressure0Idx = Indices::pressure0Idx };
65 enum { water0Idx = Indices::water0Idx };
66
67 static constexpr bool waterEnabled = Indices::waterEnabled;
68
69 // phase indices
70 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
71 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
72 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
73
74 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
75 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
76
77 using Toolbox = MathToolbox<Evaluation>;
79
80public:
81 FlashPrimaryVariables() : ParentType()
82 { Opm::Valgrind::SetDefined(*this); }
83
90
91 using ParentType::operator=;
92
96 template <class FluidState>
97 void assignMassConservative(const FluidState& fluidState,
98 const MaterialLawParams&,
99 bool = false)
100 {
101 // there is no difference between naive and mass conservative
102 // assignment in the flash model. (we only need the total
103 // concentrations.)
104 assignNaive(fluidState);
105 }
106
110 template <class FluidState>
111 void assignNaive(const FluidState& fluidState)
112 {
113 // reset everything
114 (*this) = 0.0;
115
116 // assign the phase temperatures. this is out-sourced to
117 // the energy module
118 EnergyModule::setPriVarTemperatures(*this, fluidState);
119
120 // assign components total fraction
121 for (int i = 0; i < numComponents - 1; ++i)
122 (*this)[z0Idx + i] = getValue(fluidState.moleFraction(i));
123
124 // assign pressure
125 (*this)[pressure0Idx] = getValue(fluidState.pressure(oilPhaseIdx));
126
127 // assign water saturation
128 if constexpr (waterEnabled) {
129 (*this)[water0Idx] = getValue(fluidState.saturation(waterPhaseIdx));
130 }
131 }
132
138 void print(std::ostream& os) const
139 {
140 os << "(p_" << FluidSystem::phaseName(FluidSystem::oilPhaseIdx) << " = "
141 << (*this)[pressure0Idx];
142 for (unsigned compIdx = 0; compIdx < numComponents - 2; ++compIdx) {
143 os << ", z_" << FluidSystem::componentName(compIdx) << " = "
144 << (*this)[z0Idx + compIdx];
145 }
146 if constexpr (waterEnabled) {
147 os << ", S_w = " << (*this)[water0Idx];
148 }
149 os << ")" << std::flush;
150 }
151};
152
153} // namespace Opm
154
155#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:54
Represents the primary variables used by the compositional flow model based on flash calculations.
Definition: flash/flashprimaryvariables.hh:52
void assignNaive(const FluidState &fluidState)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: flash/flashprimaryvariables.hh:106
FlashPrimaryVariables()
Definition: ptflash/flashprimaryvariables.hh:81
FlashPrimaryVariables(const FlashPrimaryVariables &value)=default
void print(std::ostream &os) const
Prints the names of the primary variables and their values.
Definition: ptflash/flashprimaryvariables.hh:138
void assignMassConservative(const FluidState &fluidState, const MaterialLawParams &, bool=false)
< Import base class assignment operators.
Definition: ptflash/flashprimaryvariables.hh:97
FlashPrimaryVariables & operator=(const FlashPrimaryVariables &value)=default
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Definition: blackoilboundaryratevector.hh:39