version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
17#ifndef DUMUX_COMPUTE_FROM_REFERENCE_PHASE_HH
18#define DUMUX_COMPUTE_FROM_REFERENCE_PHASE_HH
19
21
22#include <dune/common/fvector.hh>
23#include <dune/common/fmatrix.hh>
24
26
27namespace Dumux {
28
62template <class Scalar, class FluidSystem>
64{
65 enum { numPhases = FluidSystem::numPhases };
66 enum { numComponents = FluidSystem::numComponents };
68 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
69
70public:
98 template <class FluidState, class ParameterCache>
99 static void solve(FluidState &fluidState,
100 ParameterCache &paramCache,
101 int refPhaseIdx)
102 {
103 ComponentVector fugVec;
104
105 // compute the density and enthalpy of the
106 // reference phase
107 paramCache.updatePhase(fluidState, refPhaseIdx);
108 fluidState.setDensity(refPhaseIdx,
109 FluidSystem::density(fluidState,
110 paramCache,
111 refPhaseIdx));
112 fluidState.setMolarDensity(refPhaseIdx,
113 FluidSystem::molarDensity(fluidState,
114 paramCache,
115 refPhaseIdx));
116
117 // compute the fugacities of all components in the reference phase
118 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
119 {
120 fluidState.setFugacityCoefficient(refPhaseIdx,
121 compIdx,
122 FluidSystem::fugacityCoefficient(fluidState,
123 paramCache,
124 refPhaseIdx,
125 compIdx));
126 fugVec[compIdx] = fluidState.fugacity(refPhaseIdx, compIdx);
127 }
128
129 // compute all quantities for the non-reference phases
130 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
131 {
132 if (phaseIdx == refPhaseIdx)
133 continue; // reference phase is already calculated
134
135 CompositionFromFugacities::guessInitial(fluidState, paramCache, phaseIdx, fugVec);
136 CompositionFromFugacities::solve(fluidState, paramCache, phaseIdx, fugVec);
137 }
138 }
139};
140
141} // end namespace Dumux
142
143#endif
Calculates the chemical equilibrium from the component fugacities in the phase .
Definition: compositionfromfugacities.hh:41
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:56
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:84
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: computefromreferencephase.hh:64
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:99
Determines the fluid composition given the component fugacities and an arbitrary equation of state.
Some exceptions thrown in DuMux
std::string molarDensity(int phaseIdx) noexcept
I/O name of molar density for multiphase systems.
Definition: name.hh:71
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:53
Definition: adapt.hh:17