3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
computefromreferencephase.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
29#ifndef DUMUX_COMPUTE_FROM_REFERENCE_PHASE_HH
30#define DUMUX_COMPUTE_FROM_REFERENCE_PHASE_HH
31
33
34#include <dune/common/fvector.hh>
35#include <dune/common/fmatrix.hh>
36
39
40namespace Dumux {
41
75template <class Scalar, class FluidSystem>
77{
78 enum { numPhases = FluidSystem::numPhases };
79 enum { numComponents = FluidSystem::numComponents };
81 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
82
83public:
111 template <class FluidState, class ParameterCache>
112 static void solve(FluidState &fluidState,
113 ParameterCache &paramCache,
114 int refPhaseIdx)
115 {
116 ComponentVector fugVec;
117
118 // compute the density and enthalpy of the
119 // reference phase
120 paramCache.updatePhase(fluidState, refPhaseIdx);
121 fluidState.setDensity(refPhaseIdx,
122 FluidSystem::density(fluidState,
123 paramCache,
124 refPhaseIdx));
125 fluidState.setMolarDensity(refPhaseIdx,
126 FluidSystem::molarDensity(fluidState,
127 paramCache,
128 refPhaseIdx));
129
130 // compute the fugacities of all components in the reference phase
131 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
132 {
133 fluidState.setFugacityCoefficient(refPhaseIdx,
134 compIdx,
135 FluidSystem::fugacityCoefficient(fluidState,
136 paramCache,
137 refPhaseIdx,
138 compIdx));
139 fugVec[compIdx] = fluidState.fugacity(refPhaseIdx, compIdx);
140 }
141
142 // compute all quantities for the non-reference phases
143 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
144 {
145 if (phaseIdx == refPhaseIdx)
146 continue; // reference phase is already calculated
147
148 CompositionFromFugacities::guessInitial(fluidState, paramCache, phaseIdx, fugVec);
149 CompositionFromFugacities::solve(fluidState, paramCache, phaseIdx, fugVec);
150 }
151 }
152};
153
154} // end namespace Dumux
155
156#endif
Some exceptions thrown in DuMux
Some templates to wrap the valgrind macros.
Determines the fluid composition given the component fugacities and an arbitary equation of state.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
std::string molarDensity(int phaseIdx) noexcept
I/O name of molar density for multiphase systems.
Definition: name.hh:83
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Calculates the chemical equilibrium from the component fugacities in the phase .
Definition: compositionfromfugacities.hh:55
static void guessInitial(FluidState &fluidState, ParameterCache &paramCache, int phaseIdx, const ComponentVector &fugVec)
Guess an initial value for the composition of the phase.
Definition: compositionfromfugacities.hh:71
static void solve(FluidState &fluidState, ParameterCache &paramCache, int phaseIdx, const ComponentVector &targetFug)
Calculates the chemical equilibrium from the component fugacities in a phase.
Definition: compositionfromfugacities.hh:99
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: computefromreferencephase.hh:77
static void solve(FluidState &fluidState, ParameterCache &paramCache, int refPhaseIdx)
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: computefromreferencephase.hh:112