flash/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 EWOMS_FLASH_PRIMARY_VARIABLES_HH
29#define EWOMS_FLASH_PRIMARY_VARIABLES_HH
30
31#include <opm/material/common/MathToolbox.hpp>
32#include <opm/material/common/Valgrind.hpp>
33
36
37#include <iostream>
38
39namespace Opm {
40
50template <class TypeTag>
52{
53 using ParentType = FvBasePrimaryVariables<TypeTag>;
54
60
61 // primary variable indices
62 enum { cTot0Idx = Indices::cTot0Idx };
63
64 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
65 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
66
67 using Toolbox = MathToolbox<Evaluation>;
69
70public:
72 : ParentType()
73 { Valgrind::SetDefined(*this); }
74
85
86 using ParentType::operator=;
87
91 template <class FluidState>
92 void assignMassConservative(const FluidState& fluidState,
93 const MaterialLawParams&,
94 bool = false)
95 {
96 // there is no difference between naive and mass conservative
97 // assignment in the flash model. (we only need the total
98 // concentrations.)
99 assignNaive(fluidState);
100 }
101
105 template <class FluidState>
106 void assignNaive(const FluidState& fluidState)
107 {
108 // reset everything
109 (*this) = 0.0;
110
111 // assign the phase temperatures. this is out-sourced to
112 // the energy module
113 EnergyModule::setPriVarTemperatures(*this, fluidState);
114
115 // determine the phase presence.
116 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
117 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
118 (*this)[cTot0Idx + compIdx] +=
119 fluidState.molarity(phaseIdx, compIdx) * fluidState.saturation(phaseIdx);
120 }
121 }
122 }
123
129 void print(std::ostream& os = std::cout) const
130 {
131 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
132 os << "(c_tot," << FluidSystem::componentName(compIdx) << " = "
133 << (*this)[cTot0Idx + compIdx];
134 }
135 os << ")" << std::flush;
136 }
137};
138
139} // namespace Opm
140
141#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 print(std::ostream &os=std::cout) const
Prints the names of the primary variables and their values.
Definition: flash/flashprimaryvariables.hh:129
void assignNaive(const FluidState &fluidState)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: flash/flashprimaryvariables.hh:106
FlashPrimaryVariables()
Definition: flash/flashprimaryvariables.hh:71
FlashPrimaryVariables(const FlashPrimaryVariables &value)=default
void assignMassConservative(const FluidState &fluidState, const MaterialLawParams &, bool=false)
< Import base class assignment operators.
Definition: flash/flashprimaryvariables.hh:92
FlashPrimaryVariables & operator=(const FlashPrimaryVariables &value)=default
Represents the primary variables used by the a model.
Definition: fvbaseprimaryvariables.hh:53
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