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 "flashindices.hh"
32
35
36#include <opm/material/constraintsolvers/NcpFlash.hpp>
37#include <opm/material/fluidstates/CompositionalFluidState.hpp>
38#include <opm/material/common/Valgrind.hpp>
39
40#include <dune/common/fvector.hh>
41
42#include <iostream>
43
44namespace Opm {
45
55template <class TypeTag>
57{
58 using ParentType = FvBasePrimaryVariables<TypeTag>;
59
65
66 // primary variable indices
67 enum { cTot0Idx = Indices::cTot0Idx };
68
69 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
70 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
71
72 using Toolbox = typename Opm::MathToolbox<Evaluation>;
74
75public:
76 FlashPrimaryVariables() : ParentType()
77 { Opm::Valgrind::SetDefined(*this); }
78
82 FlashPrimaryVariables(Scalar value) : ParentType(value)
83 {
84 Opm::Valgrind::CheckDefined(value);
85 Opm::Valgrind::SetDefined(*this);
86 }
87
94
98 template <class FluidState>
99 void assignMassConservative(const FluidState& fluidState,
100 const MaterialLawParams&,
101 bool = false)
102 {
103 // there is no difference between naive and mass conservative
104 // assignment in the flash model. (we only need the total
105 // concentrations.)
106 assignNaive(fluidState);
107 }
108
112 template <class FluidState>
113 void assignNaive(const FluidState& fluidState)
114 {
115 // reset everything
116 (*this) = 0.0;
117
118 // assign the phase temperatures. this is out-sourced to
119 // the energy module
120 EnergyModule::setPriVarTemperatures(*this, fluidState);
121
122 // determine the phase presence.
123 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
124 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
125 this->operator[](cTot0Idx + compIdx) +=
126 fluidState.molarity(phaseIdx, compIdx) * fluidState.saturation(phaseIdx);
127 }
128 }
129 }
130
136 void print(std::ostream& os = std::cout) const
137 {
138 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
139 os << "(c_tot," << FluidSystem::componentName(compIdx) << " = "
140 << this->operator[](cTot0Idx + compIdx);
141 }
142 os << ")" << std::flush;
143 }
144};
145
146} // namespace Opm
147
148#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:50
Represents the primary variables used by the compositional flow model based on flash calculations.
Definition: flash/flashprimaryvariables.hh:57
void print(std::ostream &os=std::cout) const
Prints the names of the primary variables and their values.
Definition: flash/flashprimaryvariables.hh:136
void assignNaive(const FluidState &fluidState)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: flash/flashprimaryvariables.hh:113
FlashPrimaryVariables(Scalar value)
Constructor with assignment from scalar.
Definition: flash/flashprimaryvariables.hh:82
FlashPrimaryVariables()
Definition: flash/flashprimaryvariables.hh:76
FlashPrimaryVariables(const FlashPrimaryVariables &value)=default
void assignMassConservative(const FluidState &fluidState, const MaterialLawParams &, bool=false)
Set the primary variables from an arbitrary fluid state in a mass conservative way.
Definition: flash/flashprimaryvariables.hh:99
FlashPrimaryVariables & operator=(const FlashPrimaryVariables &value)=default
Represents the primary variables used by the a model.
Definition: fvbaseprimaryvariables.hh:52
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