ComputeFromReferencePhase.hpp
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*/
27#ifndef OPM_COMPUTE_FROM_REFERENCE_PHASE_HPP
28#define OPM_COMPUTE_FROM_REFERENCE_PHASE_HPP
29
31
33
34#include <dune/common/fvector.hh>
35
36namespace Opm {
37
63template <class Scalar, class FluidSystem, class Evaluation = Scalar>
65{
66 enum { numPhases = FluidSystem::numPhases };
67 enum { numComponents = FluidSystem::numComponents };
68 typedef ::Opm::CompositionFromFugacities<Scalar, FluidSystem, Evaluation> CompositionFromFugacities;
69 typedef Dune::FieldVector<Evaluation, numComponents> ComponentVector;
70
71public:
106 template <class FluidState>
107 static void solve(FluidState& fluidState,
108 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>& paramCache,
109 unsigned refPhaseIdx,
110 bool setViscosity,
111 bool setEnthalpy)
112 {
113 // compute the density and enthalpy of the
114 // reference phase
115 paramCache.updatePhase(fluidState, refPhaseIdx);
116 fluidState.setDensity(refPhaseIdx,
117 FluidSystem::density(fluidState,
118 paramCache,
119 refPhaseIdx));
120
121 if (setEnthalpy)
122 fluidState.setEnthalpy(refPhaseIdx,
123 FluidSystem::enthalpy(fluidState,
124 paramCache,
125 refPhaseIdx));
126
127 if (setViscosity)
128 fluidState.setViscosity(refPhaseIdx,
129 FluidSystem::viscosity(fluidState,
130 paramCache,
131 refPhaseIdx));
132
133 // compute the fugacities of all components in the reference phase
134 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
135 fluidState.setFugacityCoefficient(refPhaseIdx,
136 compIdx,
137 FluidSystem::fugacityCoefficient(fluidState,
138 paramCache,
139 refPhaseIdx,
140 compIdx));
141 }
142
143 // compute all quantities for the non-reference phases
144 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
145 if (phaseIdx == refPhaseIdx)
146 continue; // reference phase is already calculated
147
148 ComponentVector fugVec;
149 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
150 const auto& fug = fluidState.fugacity(refPhaseIdx, compIdx);
151 fugVec[compIdx] = decay<Evaluation>(fug);
152 }
153
154 CompositionFromFugacities::solve(fluidState, paramCache, phaseIdx, fugVec);
155
156 if (setViscosity)
157 fluidState.setViscosity(phaseIdx,
158 FluidSystem::viscosity(fluidState,
159 paramCache,
160 phaseIdx));
161
162 if (setEnthalpy)
163 fluidState.setEnthalpy(phaseIdx,
164 FluidSystem::enthalpy(fluidState,
165 paramCache,
166 phaseIdx));
167 }
168 }
169};
170
171} // namespace Opm
172
173#endif
Some templates to wrap the valgrind client request macros.
Calculates the chemical equilibrium from the component fugacities in a phase.
Definition: CompositionFromFugacities.hpp:47
static void solve(FluidState &fluidState, typename FluidSystem::template ParameterCache< typename FluidState::Scalar > &paramCache, unsigned phaseIdx, const ComponentVector &targetFug)
Calculates the chemical equilibrium from the component fugacities in a phase.
Definition: CompositionFromFugacities.hpp:80
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: ComputeFromReferencePhase.hpp:65
static void solve(FluidState &fluidState, typename FluidSystem::template ParameterCache< typename FluidState::Scalar > &paramCache, unsigned refPhaseIdx, bool setViscosity, bool setEnthalpy)
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: ComputeFromReferencePhase.hpp:107
Definition: Air_Mesitylene.hpp:34