3.5-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
38
39namespace Dumux {
40
74template <class Scalar, class FluidSystem>
76{
77 enum { numPhases = FluidSystem::numPhases };
78 enum { numComponents = FluidSystem::numComponents };
80 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
81
82public:
110 template <class FluidState, class ParameterCache>
111 static void solve(FluidState &fluidState,
112 ParameterCache &paramCache,
113 int refPhaseIdx)
114 {
115 ComponentVector fugVec;
116
117 // compute the density and enthalpy of the
118 // reference phase
119 paramCache.updatePhase(fluidState, refPhaseIdx);
120 fluidState.setDensity(refPhaseIdx,
121 FluidSystem::density(fluidState,
122 paramCache,
123 refPhaseIdx));
124 fluidState.setMolarDensity(refPhaseIdx,
125 FluidSystem::molarDensity(fluidState,
126 paramCache,
127 refPhaseIdx));
128
129 // compute the fugacities of all components in the reference phase
130 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
131 {
132 fluidState.setFugacityCoefficient(refPhaseIdx,
133 compIdx,
134 FluidSystem::fugacityCoefficient(fluidState,
135 paramCache,
136 refPhaseIdx,
137 compIdx));
138 fugVec[compIdx] = fluidState.fugacity(refPhaseIdx, compIdx);
139 }
140
141 // compute all quantities for the non-reference phases
142 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
143 {
144 if (phaseIdx == refPhaseIdx)
145 continue; // reference phase is already calculated
146
147 CompositionFromFugacities::guessInitial(fluidState, paramCache, phaseIdx, fugVec);
148 CompositionFromFugacities::solve(fluidState, paramCache, phaseIdx, fugVec);
149 }
150 }
151};
152
153} // end namespace Dumux
154
155#endif
Some exceptions thrown in DuMux
Determines the fluid composition given the component fugacities and an arbitary equation of state.
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:53
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:69
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:97
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: computefromreferencephase.hh:76
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:111