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