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