3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
compositionfromfugacities.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 *****************************************************************************/
25#ifndef DUMUX_COMPOSITION_FROM_FUGACITIES_HH
26#define DUMUX_COMPOSITION_FROM_FUGACITIES_HH
27
28#include <dune/common/fvector.hh>
29#include <dune/common/fmatrix.hh>
30
32
33namespace Dumux {
34
51template <class Scalar, class FluidSystem>
53{
54 enum { numComponents = FluidSystem::numComponents };
55
56 using ParameterCache = typename FluidSystem::ParameterCache;
57
58public:
59 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
60
68 template <class FluidState>
69 static void guessInitial(FluidState &fluidState,
70 ParameterCache &paramCache,
71 int phaseIdx,
72 const ComponentVector &fugVec)
73 {
74 if (!FluidSystem::isIdealMixture(phaseIdx))
75 {
76 // Pure component fugacities
77 for (unsigned int i = 0; i < numComponents; ++ i)
78 {
79 fluidState.setMoleFraction(phaseIdx,
80 i,
81 1.0/numComponents);
82 }
83 }
84 }
85
96 template <class FluidState>
97 static void solve(FluidState &fluidState,
98 ParameterCache &paramCache,
99 int phaseIdx,
100 const ComponentVector &targetFug)
101 {
102 // use a much more efficient method in case the phase is an
103 // ideal mixture
104 if (FluidSystem::isIdealMixture(phaseIdx)) {
105 solveIdealMix_(fluidState, paramCache, phaseIdx, targetFug);
106 return;
107 }
108
109 // save initial composition in case something goes wrong
110 Dune::FieldVector<Scalar, numComponents> xInit;
111 for (int i = 0; i < numComponents; ++i) {
112 xInit[i] = fluidState.moleFraction(phaseIdx, i);
113 }
114
116 // Newton method
118
119 // Jacobian matrix
120 Dune::FieldMatrix<Scalar, numComponents, numComponents> J;
121 // solution, i.e. phase composition
122 Dune::FieldVector<Scalar, numComponents> x;
123 // right hand side
124 Dune::FieldVector<Scalar, numComponents> b;
125
126 paramCache.updatePhase(fluidState, phaseIdx);
127
128 // maximum number of iterations
129 const int nMax = 25;
130 for (int nIdx = 0; nIdx < nMax; ++nIdx) {
131 // calculate Jacobian matrix and right hand side
132 linearize_(J, b, fluidState, paramCache, phaseIdx, targetFug);
133
134 // Solve J*x = b
135 x = 0;
136 try { J.solve(x, b); }
137 catch (Dune::FMatrixError &e)
138 { throw NumericalProblem(e.what()); }
139
140 // update the fluid composition. b is also used to store
141 // the defect for the next iteration.
142 Scalar relError = update_(fluidState, paramCache, x, b, phaseIdx, targetFug);
143
144 if (relError < 1e-9) {
145 Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
146 Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
147 fluidState.setDensity(phaseIdx, rho);
148 fluidState.setMolarDensity(phaseIdx, rhoMolar);
149
150 //std::cout << "num iterations: " << nIdx << "\n";
151 return;
152 }
153 }
154
155 DUNE_THROW(NumericalProblem,
156 "Calculating the " << FluidSystem::phaseName(phaseIdx)
157 << "Phase composition failed. Initial {x} = {"
158 << xInit
159 << "}, {fug_t} = {" << targetFug << "}, p = " << fluidState.pressure(phaseIdx)
160 << ", T = " << fluidState.temperature(phaseIdx));
161 }
162
163
164protected:
165 // update the phase composition in case the phase is an ideal
166 // mixture, i.e. the component's fugacity coefficients are
167 // independent of the phase's composition.
168 template <class FluidState>
169 static void solveIdealMix_(FluidState &fluidState,
170 ParameterCache &paramCache,
171 int phaseIdx,
172 const ComponentVector &fugacities)
173 {
174 for (int i = 0; i < numComponents; ++ i) {
175 Scalar phi = FluidSystem::fugacityCoefficient(fluidState,
176 paramCache,
177 phaseIdx,
178 i);
179 Scalar gamma = phi * fluidState.pressure(phaseIdx);
180 fluidState.setFugacityCoefficient(phaseIdx, i, phi);
181 fluidState.setMoleFraction(phaseIdx, i, fugacities[i]/gamma);
182 }
183
184 paramCache.updatePhase(fluidState, phaseIdx);
185
186 Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
187 fluidState.setDensity(phaseIdx, rho);
188 Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
189 fluidState.setMolarDensity(phaseIdx, rhoMolar);
190 }
191
192 template <class FluidState>
193 static void linearize_(Dune::FieldMatrix<Scalar, numComponents, numComponents> &J,
194 Dune::FieldVector<Scalar, numComponents> &defect,
195 FluidState &fluidState,
196 ParameterCache &paramCache,
197 int phaseIdx,
198 const ComponentVector &targetFug)
199 {
200 // reset jacobian
201 J = 0;
202
203 // calculate the defect (deviation of the current fugacities
204 // from the target fugacities)
205 for (int i = 0; i < numComponents; ++ i) {
206 Scalar phi = FluidSystem::fugacityCoefficient(fluidState,
207 paramCache,
208 phaseIdx,
209 i);
210 Scalar f = phi*fluidState.pressure(phaseIdx)*fluidState.moleFraction(phaseIdx, i);
211 fluidState.setFugacityCoefficient(phaseIdx, i, phi);
212
213 defect[i] = targetFug[i] - f;
214 }
215
216 // assemble jacobian matrix of the constraints for the composition
217 for (int i = 0; i < numComponents; ++ i) {
218 const Scalar eps = 1e-11; //std::max(1e-16, std::abs(x_i)*1e-9);
219
221 // approximately calculate partial derivatives of the
222 // fugacity defect of all components in regard to the mole
223 // fraction of the i-th component. This is done via
224 // forward differences
225
226 // deviate the mole fraction of the i-th component
227 Scalar x_i = fluidState.moleFraction(phaseIdx, i);
228 fluidState.setMoleFraction(phaseIdx, i, x_i + eps);
229 paramCache.updateSingleMoleFraction(fluidState, phaseIdx, i);
230
231 // compute new defect and derivative for all component
232 // fugacities
233 for (int j = 0; j < numComponents; ++j) {
234 // compute the j-th component's fugacity coefficient ...
235 Scalar phi = FluidSystem::fugacityCoefficient(fluidState,
236 paramCache,
237 phaseIdx,
238 j);
239 // ... and its fugacity ...
240 Scalar f =
241 phi *
242 fluidState.pressure(phaseIdx) *
243 fluidState.moleFraction(phaseIdx, j);
244 // as well as the defect for this fugacity
245 Scalar defJPlusEps = targetFug[j] - f;
246
247 // use forward differences to calculate the defect's
248 // derivative
249 J[j][i] = (defJPlusEps - defect[j])/eps;
250 }
251
252 // reset composition to original value
253 fluidState.setMoleFraction(phaseIdx, i, x_i);
254 paramCache.updateSingleMoleFraction(fluidState, phaseIdx, i);
255
256 // end forward differences
258 }
259 }
260
261 template <class FluidState>
262 static Scalar update_(FluidState &fluidState,
263 ParameterCache &paramCache,
264 Dune::FieldVector<Scalar, numComponents> &x,
265 Dune::FieldVector<Scalar, numComponents> &b,
266 int phaseIdx,
267 const Dune::FieldVector<Scalar, numComponents> &targetFug)
268 {
269 // store original composition and calculate relative error
270 Dune::FieldVector<Scalar, numComponents> origComp;
271 Scalar relError = 0;
272 Scalar sumDelta = 0;
273 using std::max;
274 using std::abs;
275 using std::min;
276 for (unsigned int i = 0; i < numComponents; ++i)
277 {
278 origComp[i] = fluidState.moleFraction(phaseIdx, i);
279 relError = max(relError, abs(x[i]));
280 sumDelta += abs(x[i]);
281 }
282
283 // chop update to at most 20% change in composition
284 const Scalar maxDelta = 0.2;
285 if (sumDelta > maxDelta)
286 x /= (sumDelta/maxDelta);
287
288 // change composition
289 for (unsigned int i = 0; i < numComponents; ++i)
290 {
291 Scalar newx = origComp[i] - x[i];
292
293 // only allow negative mole fractions if the target fugacity is negative
294 if (targetFug[i] > 0)
295 newx = max(0.0, newx);
296 // only allow positive mole fractions if the target fugacity is positive
297 else if (targetFug[i] < 0)
298 newx = min(0.0, newx);
299 // if the target fugacity is zero, the mole fraction must also be zero
300 else
301 newx = 0;
302
303 fluidState.setMoleFraction(phaseIdx, i, newx);
304 //sumMoleFrac += abs(newx);
305 }
306
307 paramCache.updateComposition(fluidState, phaseIdx);
308
309 return relError;
310 }
311
312 template <class FluidState>
313 static Scalar calculateDefect_(const FluidState &params,
314 int phaseIdx,
315 const ComponentVector &targetFug)
316 {
317 Scalar result = 0.0;
318 using std::abs;
319 for (int i = 0; i < numComponents; ++i) {
320 // sum of the fugacity defect weighted by the inverse
321 // fugacity coefficient
322 result += abs(
323 (targetFug[i] - params.fugacity(phaseIdx, i))
324 /
325 params.fugacityCoeff(phaseIdx, i) );
326 }
327 return result;
328 }
329};
330
331} // end namespace Dumux
332
333#endif
Some exceptions thrown in DuMux
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
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
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 Scalar calculateDefect_(const FluidState &params, int phaseIdx, const ComponentVector &targetFug)
Definition: compositionfromfugacities.hh:313
static Scalar update_(FluidState &fluidState, ParameterCache &paramCache, Dune::FieldVector< Scalar, numComponents > &x, Dune::FieldVector< Scalar, numComponents > &b, int phaseIdx, const Dune::FieldVector< Scalar, numComponents > &targetFug)
Definition: compositionfromfugacities.hh:262
Dune::FieldVector< Scalar, numComponents > ComponentVector
Definition: compositionfromfugacities.hh:59
static void linearize_(Dune::FieldMatrix< Scalar, numComponents, numComponents > &J, Dune::FieldVector< Scalar, numComponents > &defect, FluidState &fluidState, ParameterCache &paramCache, int phaseIdx, const ComponentVector &targetFug)
Definition: compositionfromfugacities.hh:193
static void solveIdealMix_(FluidState &fluidState, ParameterCache &paramCache, int phaseIdx, const ComponentVector &fugacities)
Definition: compositionfromfugacities.hh:169
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