3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
immiscibleflash.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_IMMISCIBLE_FLASH_HH
26#define DUMUX_IMMISCIBLE_FLASH_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
64template <class Scalar, class FluidSystem>
66{
67 static constexpr int numPhases = FluidSystem::numPhases;
68 static constexpr int numComponents = FluidSystem::numComponents;
69 static_assert(numPhases == numComponents,
70 "Immiscibility assumes that the number of phases"
71 " is equal to the number of components");
72
73 using ParameterCache = typename FluidSystem::ParameterCache;
74
75 static constexpr int numEq = numPhases;
76
77 using Matrix = Dune::FieldMatrix<Scalar, numEq, numEq>;
78 using Vector = Dune::FieldVector<Scalar, numEq>;
79
80public:
81 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
82
87 explicit ImmiscibleFlash(int wettingPhaseIdx = 0)
88 : wettingPhaseIdx_(wettingPhaseIdx)
89 {}
90
97 template <class FluidState>
98 void guessInitial(FluidState &fluidState,
99 ParameterCache &paramCache,
100 const ComponentVector &globalMolarities)
101 {
102 // the sum of all molarities
103 Scalar sumMoles = 0;
104 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
105 sumMoles += globalMolarities[compIdx];
106
107 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
108 // pressure. use atmospheric pressure as initial guess
109 fluidState.setPressure(phaseIdx, 2e5);
110
111 // saturation. assume all fluids to be equally distributed
112 fluidState.setSaturation(phaseIdx, 1.0/numPhases);
113 }
114 }
115
126 template <class MaterialLaw, class FluidState>
127 void solve(FluidState &fluidState,
128 ParameterCache &paramCache,
129 const typename MaterialLaw::Params &matParams,
130 const ComponentVector &globalMolarities)
131 {
132#if DUNE_VERSION_LT(DUNE_COMMON, 2, 7)
133 Dune::FMatrixPrecision<Scalar>::set_singular_limit(1e-25);
134#endif
135
137 // Check if all fluid phases are incompressible
139 bool allIncompressible = true;
140 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
141 if (FluidSystem::isCompressible(phaseIdx)) {
142 allIncompressible = false;
143 break;
144 }
145 }
146
147 if (allIncompressible) {
148 paramCache.updateAll(fluidState);
149 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
150 Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
151 Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
152 fluidState.setDensity(phaseIdx, rho);
153 fluidState.setMolarDensity(phaseIdx, rhoMolar);
154
155 Scalar saturation =
156 globalMolarities[/*compIdx=*/phaseIdx]
157 / fluidState.molarDensity(phaseIdx);
158 fluidState.setSaturation(phaseIdx, saturation);
159 }
160 }
161
162// /////////////////////////
163 // Newton method
165
166 // Jacobian matrix
167 Matrix J;
168 // solution, i.e. phase composition
169 Vector deltaX;
170 // right hand side
171 Vector b;
172
176
177 completeFluidState_<MaterialLaw>(fluidState, paramCache, matParams);
178
179 const int nMax = 50; // <- maximum number of newton iterations
180 for (int nIdx = 0; nIdx < nMax; ++nIdx) {
181 // calculate Jacobian matrix and right hand side
182 linearize_<MaterialLaw>(J, b, fluidState, paramCache, matParams, globalMolarities);
185
186 // Solve J*x = b
187 deltaX = 0;
188
189 try { J.solve(deltaX, b); }
190 catch (Dune::FMatrixError& e)
191 {
192 /*
193 std::cout << "error: " << e << "\n";
194 std::cout << "b: " << b << "\n";
195 */
196
197 throw NumericalProblem(e.what());
198 }
200
201 /*
202 printFluidState_(fluidState);
203 std::cout << "J:\n";
204 for (int i = 0; i < numEq; ++i) {
205 for (int j = 0; j < numEq; ++j) {
206 std::ostringstream os;
207 os << J[i][j];
208
209 std::string s(os.str());
210 do {
211 s += " ";
212 } while (s.size() < 20);
213 std::cout << s;
214 }
215 std::cout << "\n";
216 }
217 std::cout << "residual: " << b << "\n";
218 std::cout << "deltaX: " << deltaX << "\n";
219 std::cout << "---------------\n";
220 */
221
222 // update the fluid quantities.
223 Scalar relError = update_<MaterialLaw>(fluidState, paramCache, matParams, deltaX);
224
225 if (relError < 1e-9)
226 return;
227 }
228
229 /*
230 printFluidState_(fluidState);
231 std::cout << "globalMolarities: ";
232 for (int compIdx = 0; compIdx < numComponents; ++ compIdx)
233 std::cout << globalMolarities[compIdx] << " ";
234 std::cout << "\n";
235 */
236
237 DUNE_THROW(NumericalProblem,
238 "Flash calculation failed."
239 " {c_alpha^kappa} = {" << globalMolarities << "}, T = "
240 << fluidState.temperature(/*phaseIdx=*/0));
241 }
242
243
244protected:
245 template <class FluidState>
246 void printFluidState_(const FluidState &fs)
247 {
248 std::cout << "saturations: ";
249 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
250 std::cout << fs.saturation(phaseIdx) << " ";
251 std::cout << "\n";
252
253 std::cout << "pressures: ";
254 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
255 std::cout << fs.pressure(phaseIdx) << " ";
256 std::cout << "\n";
257
258 std::cout << "densities: ";
259 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
260 std::cout << fs.density(phaseIdx) << " ";
261 std::cout << "\n";
262
263 std::cout << "global component molarities: ";
264 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
265 Scalar sum = 0;
266 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
267 sum += fs.saturation(phaseIdx)*fs.molarity(phaseIdx, compIdx);
268 }
269 std::cout << sum << " ";
270 }
271 std::cout << "\n";
272 }
273
274 template <class MaterialLaw, class FluidState>
275 void linearize_(Matrix &J,
276 Vector &b,
277 FluidState &fluidState,
278 ParameterCache &paramCache,
279 const typename MaterialLaw::Params &matParams,
280 const ComponentVector &globalMolarities)
281 {
282 FluidState origFluidState(fluidState);
283 ParameterCache origParamCache(paramCache);
284
285 Vector tmp;
286
287 // reset jacobian
288 J = 0;
289
291 calculateDefect_(b, fluidState, fluidState, globalMolarities);
293
294 // assemble jacobian matrix
295 for (int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
297 // approximately calculate partial derivatives of the
298 // fugacity defect of all components in regard to the mole
299 // fraction of the i-th component. This is done via
300 // forward differences
301
302 // deviate the mole fraction of the i-th component
303 Scalar x_i = getQuantity_(fluidState, pvIdx);
304 const Scalar eps = 1e-10/quantityWeight_(fluidState, pvIdx);
305 setQuantity_<MaterialLaw>(fluidState, paramCache, matParams, pvIdx, x_i + eps);
306
307 // compute derivative of the defect
308 calculateDefect_(tmp, origFluidState, fluidState, globalMolarities);
309 tmp -= b;
310 tmp /= eps;
311 // store derivative in jacobian matrix
312 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
313 J[eqIdx][pvIdx] = tmp[eqIdx];
314
315 // fluid state and parameter cache to their original values
316 fluidState = origFluidState;
317 paramCache = origParamCache;
318
319 // end forward differences
321 }
322 }
323
324 template <class FluidState>
325 void calculateDefect_(Vector &b,
326 const FluidState &fluidStateEval,
327 const FluidState &fluidState,
328 const ComponentVector &globalMolarities)
329 {
330 // global molarities are given
331 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
332 b[phaseIdx] =
333 fluidState.saturation(phaseIdx)
334 * fluidState.molarity(phaseIdx, /*compIdx=*/phaseIdx);
335
336 b[phaseIdx] -= globalMolarities[/*compIdx=*/phaseIdx];
337 }
338 }
339
340 template <class MaterialLaw, class FluidState>
341 Scalar update_(FluidState &fluidState,
342 ParameterCache &paramCache,
343 const typename MaterialLaw::Params &matParams,
344 const Vector &deltaX)
345 {
346 Scalar relError = 0;
347 for (int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
348 Scalar tmp = getQuantity_(fluidState, pvIdx);
349 Scalar delta = deltaX[pvIdx];
350 using std::max;
351 using std::abs;
352 using std::min;
353 relError = max(relError, abs(delta)*quantityWeight_(fluidState, pvIdx));
354
355 if (isSaturationIdx_(pvIdx)) {
356 // dampen to at most 20% change in saturation per
357 // iteration
358 delta = min(0.2, max(-0.2, delta));
359 }
360 else if (isPressureIdx_(pvIdx)) {
361 // dampen to at most 30% change in pressure per
362 // iteration
363 delta = min(0.30*fluidState.pressure(0), max(-0.30*fluidState.pressure(0), delta));
364 }
365
366 setQuantityRaw_(fluidState, pvIdx, tmp - delta);
367 }
368
369 completeFluidState_<MaterialLaw>(fluidState, paramCache, matParams);
370
371 return relError;
372 }
373
374 template <class MaterialLaw, class FluidState>
375 void completeFluidState_(FluidState &fluidState,
376 ParameterCache &paramCache,
377 const typename MaterialLaw::Params &matParams)
378 {
379 // calculate the saturation of the last phase as a function of
380 // the other saturations
381 Scalar sumSat = 0.0;
382 for (int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
383 sumSat += fluidState.saturation(phaseIdx);
384
385 if (sumSat > 1.0) {
386 // make sure that the last saturation does not become
387 // negative
388 for (int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
389 {
390 Scalar S = fluidState.saturation(phaseIdx);
391 fluidState.setSaturation(phaseIdx, S/sumSat);
392 }
393 sumSat = 1;
394 }
395
396 // set the last saturation
397 fluidState.setSaturation(/*phaseIdx=*/numPhases - 1, 1.0 - sumSat);
398
399 // update the pressures using the material law (saturations
400 // and first pressure are already set because it is implicitly
401 // solved for.)
403 MaterialLaw::capillaryPressures(pc, matParams, fluidState, wettingPhaseIdx_);
404 for (int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
405 fluidState.setPressure(phaseIdx,
406 fluidState.pressure(0)
407 + (pc[phaseIdx] - pc[0]));
408
409 // update the parameter cache
410 paramCache.updateAll(fluidState, /*except=*/ParameterCache::Temperature|ParameterCache::Composition);
411
412 // update all densities
413 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
414 Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
415 Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
416 fluidState.setDensity(phaseIdx, rho);
417 fluidState.setMolarDensity(phaseIdx, rhoMolar);
418 }
419 }
420
421 bool isPressureIdx_(int pvIdx)
422 { return pvIdx == 0; }
423
424 bool isSaturationIdx_(int pvIdx)
425 { return 1 <= pvIdx && pvIdx < numPhases; }
426
427 // retrieves a quantity from the fluid state
428 template <class FluidState>
429 Scalar getQuantity_(const FluidState &fs, int pvIdx)
430 {
431 assert(pvIdx < numEq);
432
433 // first pressure
434 if (pvIdx < 1) {
435 int phaseIdx = 0;
436 return fs.pressure(phaseIdx);
437 }
438 // saturations
439 else {
440 assert(pvIdx < numPhases);
441 int phaseIdx = pvIdx - 1;
442 return fs.saturation(phaseIdx);
443 }
444 }
445
446 // set a quantity in the fluid state
447 template <class MaterialLaw, class FluidState>
448 void setQuantity_(FluidState &fs,
449 ParameterCache &paramCache,
450 const typename MaterialLaw::Params &matParams,
451 int pvIdx,
452 Scalar value)
453 {
454 assert(0 <= pvIdx && pvIdx < numEq);
455
456 if (pvIdx < 1) {
457 // -> first pressure
458 Scalar delta = value - fs.pressure(0);
459
460 // set all pressures. here we assume that the capillary
461 // pressure does not depend on absolute pressure.
462 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
463 fs.setPressure(phaseIdx, fs.pressure(phaseIdx) + delta);
464 paramCache.updateAllPressures(fs);
465
466 // update all densities
467 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
468 Scalar rho = FluidSystem::density(fs, paramCache, phaseIdx);
469 Scalar rhoMolar = FluidSystem::molarDensity(fs, paramCache, phaseIdx);
470 fs.setDensity(phaseIdx, rho);
471 fs.setMolarDensity(phaseIdx, rhoMolar);
472 }
473 }
474 else {
475 assert(pvIdx < numPhases);
476
477 // -> first M-1 saturations
478 Scalar delta = value - fs.saturation(/*phaseIdx=*/pvIdx - 1);
479 fs.setSaturation(/*phaseIdx=*/pvIdx - 1, value);
480
481 // set last saturation (-> minus the change of the
482 // saturation of the other phases)
483 fs.setSaturation(/*phaseIdx=*/numPhases - 1,
484 fs.saturation(numPhases - 1) - delta);
485
486 // update all fluid pressures using the capillary pressure
487 // law
489 MaterialLaw::capillaryPressures(pc, matParams, fs, wettingPhaseIdx_);
490 for (int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
491 fs.setPressure(phaseIdx,
492 fs.pressure(0)
493 + (pc[phaseIdx] - pc[0]));
494 paramCache.updateAllPressures(fs);
495
496 // update all densities
497 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
498 Scalar rho = FluidSystem::density(fs, paramCache, phaseIdx);
499 Scalar rhoMolar = FluidSystem::molarDensity(fs, paramCache, phaseIdx);
500 fs.setDensity(phaseIdx, rho);
501 fs.setMolarDensity(phaseIdx, rhoMolar);
502
503 }
504 }
505 }
506
507 // set a quantity in the fluid state
508 template <class FluidState>
509 void setQuantityRaw_(FluidState &fs, int pvIdx, Scalar value)
510 {
511 assert(pvIdx < numEq);
512
513 // first pressure
514 if (pvIdx < 1) {
515 int phaseIdx = 0;
516 fs.setPressure(phaseIdx, value);
517 }
518 // saturations
519 else {
520 assert(pvIdx < numPhases);
521 int phaseIdx = pvIdx - 1;
522
523 // make sure that the first M-1 saturations does not get
524 // negative
525 using std::max;
526 value = max(0.0, value);
527 fs.setSaturation(phaseIdx, value);
528 }
529 }
530
531 template <class FluidState>
532 Scalar quantityWeight_(const FluidState &fs, int pvIdx)
533 {
534 // first pressure
535 if (pvIdx < 1)
536 return 1.0/1e5;
537 // first M - 1 saturations
538 else {
539 assert(pvIdx < numPhases);
540 return 1.0;
541 }
542 }
543
544private:
545 int wettingPhaseIdx_ = 0;
546
547};
548
549} // end namespace Dumux
550
551#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
void SetUndefined(const T &value)
Make the memory on which an object resides undefined.
Definition: valgrind.hh:102
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
std::string saturation(int phaseIdx) noexcept
I/O name of saturation for multiphase systems.
Definition: name.hh:43
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
Determines the pressures and saturations of all fluid phases given the total mass of all components.
Definition: immiscibleflash.hh:66
void guessInitial(FluidState &fluidState, ParameterCache &paramCache, const ComponentVector &globalMolarities)
Guess initial values for all quantities.
Definition: immiscibleflash.hh:98
void printFluidState_(const FluidState &fs)
Definition: immiscibleflash.hh:246
bool isPressureIdx_(int pvIdx)
Definition: immiscibleflash.hh:421
void setQuantityRaw_(FluidState &fs, int pvIdx, Scalar value)
Definition: immiscibleflash.hh:509
void calculateDefect_(Vector &b, const FluidState &fluidStateEval, const FluidState &fluidState, const ComponentVector &globalMolarities)
Definition: immiscibleflash.hh:325
Scalar getQuantity_(const FluidState &fs, int pvIdx)
Definition: immiscibleflash.hh:429
Scalar update_(FluidState &fluidState, ParameterCache &paramCache, const typename MaterialLaw::Params &matParams, const Vector &deltaX)
Definition: immiscibleflash.hh:341
void solve(FluidState &fluidState, ParameterCache &paramCache, const typename MaterialLaw::Params &matParams, const ComponentVector &globalMolarities)
Calculates the chemical equilibrium from the component fugacities in a phase.
Definition: immiscibleflash.hh:127
void linearize_(Matrix &J, Vector &b, FluidState &fluidState, ParameterCache &paramCache, const typename MaterialLaw::Params &matParams, const ComponentVector &globalMolarities)
Definition: immiscibleflash.hh:275
void setQuantity_(FluidState &fs, ParameterCache &paramCache, const typename MaterialLaw::Params &matParams, int pvIdx, Scalar value)
Definition: immiscibleflash.hh:448
bool isSaturationIdx_(int pvIdx)
Definition: immiscibleflash.hh:424
void completeFluidState_(FluidState &fluidState, ParameterCache &paramCache, const typename MaterialLaw::Params &matParams)
Definition: immiscibleflash.hh:375
Scalar quantityWeight_(const FluidState &fs, int pvIdx)
Definition: immiscibleflash.hh:532
ImmiscibleFlash(int wettingPhaseIdx=0)
Contruct a new flash.
Definition: immiscibleflash.hh:87
Dune::FieldVector< Scalar, numComponents > ComponentVector
Definition: immiscibleflash.hh:81