version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_IMMISCIBLE_FLASH_HH
14#define DUMUX_IMMISCIBLE_FLASH_HH
15
16#include <dune/common/fvector.hh>
17#include <dune/common/fmatrix.hh>
18
21
22namespace Dumux {
23
51template <class Scalar, class FluidSystem>
53{
54 static constexpr int numPhases = FluidSystem::numPhases;
55 static constexpr int numComponents = FluidSystem::numComponents;
56 static_assert(numPhases == numComponents,
57 "Immiscibility assumes that the number of phases"
58 " is equal to the number of components");
59
60 using ParameterCache = typename FluidSystem::ParameterCache;
61
62 static constexpr int numEq = numPhases;
63
64 using Matrix = Dune::FieldMatrix<Scalar, numEq, numEq>;
65 using Vector = Dune::FieldVector<Scalar, numEq>;
66
67public:
68 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
69
74 explicit ImmiscibleFlash(int wettingPhaseIdx = 0)
75 : wettingPhaseIdx_(wettingPhaseIdx)
76 {}
77
84 template <class FluidState>
85 void guessInitial(FluidState &fluidState,
86 ParameterCache &paramCache,
87 const ComponentVector &globalMolarities)
88 {
89 // the sum of all molarities
90 Scalar sumMoles = 0;
91 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
92 sumMoles += globalMolarities[compIdx];
93
94 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
95 // pressure. use atmospheric pressure as initial guess
96 fluidState.setPressure(phaseIdx, 2e5);
97
98 // saturation. assume all fluids to be equally distributed
99 fluidState.setSaturation(phaseIdx, 1.0/numPhases);
100 }
101 }
102
113 template <class MaterialLaw, class FluidState>
114 void solve(FluidState& fluidState,
115 ParameterCache& paramCache,
116 const MaterialLaw& material,
117 const ComponentVector& globalMolarities)
118 {
120 // Check if all fluid phases are incompressible
122 bool allIncompressible = true;
123 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
124 if (FluidSystem::isCompressible(phaseIdx)) {
125 allIncompressible = false;
126 break;
127 }
128 }
129
130 if (allIncompressible) {
131 paramCache.updateAll(fluidState);
132 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
133 Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
134 Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
135 fluidState.setDensity(phaseIdx, rho);
136 fluidState.setMolarDensity(phaseIdx, rhoMolar);
137
138 Scalar saturation =
139 globalMolarities[/*compIdx=*/phaseIdx]
140 / fluidState.molarDensity(phaseIdx);
141 fluidState.setSaturation(phaseIdx, saturation);
142 }
143 }
144
145// /////////////////////////
146 // Newton method
148
149 // Jacobian matrix
150 Matrix J;
151 // solution, i.e. phase composition
152 Vector deltaX;
153 // right hand side
154 Vector b;
155
156 completeFluidState_(fluidState, paramCache, material);
157
158 const int nMax = 50; // <- maximum number of newton iterations
159 for (int nIdx = 0; nIdx < nMax; ++nIdx) {
160 // calculate Jacobian matrix and right hand side
161 linearize_(J, b, fluidState, paramCache, material, globalMolarities);
162
163 // Solve J*x = b
164 deltaX = 0;
165
166 try { J.solve(deltaX, b); }
167 catch (Dune::FMatrixError& e)
168 {
169 /*
170 std::cout << "error: " << e << "\n";
171 std::cout << "b: " << b << "\n";
172 */
173
174 throw NumericalProblem(e.what());
175 }
176
177 /*
178 printFluidState_(fluidState);
179 std::cout << "J:\n";
180 for (int i = 0; i < numEq; ++i) {
181 for (int j = 0; j < numEq; ++j) {
182 std::ostringstream os;
183 os << J[i][j];
184
185 std::string s(os.str());
186 do {
187 s += " ";
188 } while (s.size() < 20);
189 std::cout << s;
190 }
191 std::cout << "\n";
192 }
193 std::cout << "residual: " << b << "\n";
194 std::cout << "deltaX: " << deltaX << "\n";
195 std::cout << "---------------\n";
196 */
197
198 // update the fluid quantities.
199 Scalar relError = update_(fluidState, paramCache, material, deltaX);
200
201 if (relError < 1e-9)
202 return;
203 }
204
205 /*
206 printFluidState_(fluidState);
207 std::cout << "globalMolarities: ";
208 for (int compIdx = 0; compIdx < numComponents; ++ compIdx)
209 std::cout << globalMolarities[compIdx] << " ";
210 std::cout << "\n";
211 */
212
213 DUNE_THROW(NumericalProblem,
214 "Flash calculation failed."
215 " {c_alpha^kappa} = {" << globalMolarities << "}, T = "
216 << fluidState.temperature(/*phaseIdx=*/0));
217 }
218
219
220protected:
221 template <class FluidState>
222 void printFluidState_(const FluidState &fs)
223 {
224 std::cout << "saturations: ";
225 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
226 std::cout << fs.saturation(phaseIdx) << " ";
227 std::cout << "\n";
228
229 std::cout << "pressures: ";
230 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
231 std::cout << fs.pressure(phaseIdx) << " ";
232 std::cout << "\n";
233
234 std::cout << "densities: ";
235 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
236 std::cout << fs.density(phaseIdx) << " ";
237 std::cout << "\n";
238
239 std::cout << "global component molarities: ";
240 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
241 Scalar sum = 0;
242 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
243 sum += fs.saturation(phaseIdx)*fs.molarity(phaseIdx, compIdx);
244 }
245 std::cout << sum << " ";
246 }
247 std::cout << "\n";
248 }
249
250 template <class MaterialLaw, class FluidState>
251 void linearize_(Matrix &J,
252 Vector &b,
253 FluidState &fluidState,
254 ParameterCache &paramCache,
255 const MaterialLaw& material,
256 const ComponentVector &globalMolarities)
257 {
258 FluidState origFluidState(fluidState);
259 ParameterCache origParamCache(paramCache);
260
261 Vector tmp;
262
263 // reset jacobian
264 J = 0;
265
266 calculateDefect_(b, fluidState, fluidState, globalMolarities);
267
268 // assemble jacobian matrix
269 for (int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
271 // approximately calculate partial derivatives of the
272 // fugacity defect of all components in regard to the mole
273 // fraction of the i-th component. This is done via
274 // forward differences
275
276 // deviate the mole fraction of the i-th component
277 Scalar x_i = getQuantity_(fluidState, pvIdx);
278 const Scalar eps = 1e-10/quantityWeight_(fluidState, pvIdx);
279 setQuantity_<MaterialLaw>(fluidState, paramCache, material, pvIdx, x_i + eps);
280
281 // compute derivative of the defect
282 calculateDefect_(tmp, origFluidState, fluidState, globalMolarities);
283 tmp -= b;
284 tmp /= eps;
285 // store derivative in jacobian matrix
286 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
287 J[eqIdx][pvIdx] = tmp[eqIdx];
288
289 // fluid state and parameter cache to their original values
290 fluidState = origFluidState;
291 paramCache = origParamCache;
292
293 // end forward differences
295 }
296 }
297
298 template <class FluidState>
299 void calculateDefect_(Vector &b,
300 const FluidState &fluidStateEval,
301 const FluidState &fluidState,
302 const ComponentVector &globalMolarities)
303 {
304 // global molarities are given
305 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
306 b[phaseIdx] =
307 fluidState.saturation(phaseIdx)
308 * fluidState.molarity(phaseIdx, /*compIdx=*/phaseIdx);
309
310 b[phaseIdx] -= globalMolarities[/*compIdx=*/phaseIdx];
311 }
312 }
313
314 template <class MaterialLaw, class FluidState>
315 Scalar update_(FluidState &fluidState,
316 ParameterCache &paramCache,
317 const MaterialLaw& material,
318 const Vector &deltaX)
319 {
320 Scalar relError = 0;
321 for (int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
322 Scalar tmp = getQuantity_(fluidState, pvIdx);
323 Scalar delta = deltaX[pvIdx];
324 using std::max;
325 using std::abs;
326 using std::min;
327 relError = max(relError, abs(delta)*quantityWeight_(fluidState, pvIdx));
328
329 if (isSaturationIdx_(pvIdx)) {
330 // dampen to at most 20% change in saturation per
331 // iteration
332 delta = min(0.2, max(-0.2, delta));
333 }
334 else if (isPressureIdx_(pvIdx)) {
335 // dampen to at most 30% change in pressure per
336 // iteration
337 delta = min(0.30*fluidState.pressure(0), max(-0.30*fluidState.pressure(0), delta));
338 }
339
340 setQuantityRaw_(fluidState, pvIdx, tmp - delta);
341 }
342
343 completeFluidState_<MaterialLaw>(fluidState, paramCache, material);
344
345 return relError;
346 }
347
348 template <class MaterialLaw, class FluidState>
349 void completeFluidState_(FluidState &fluidState,
350 ParameterCache &paramCache,
351 const MaterialLaw& material)
352 {
353 // calculate the saturation of the last phase as a function of
354 // the other saturations
355 Scalar sumSat = 0.0;
356 for (int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
357 sumSat += fluidState.saturation(phaseIdx);
358
359 if (sumSat > 1.0) {
360 // make sure that the last saturation does not become
361 // negative
362 for (int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
363 {
364 Scalar S = fluidState.saturation(phaseIdx);
365 fluidState.setSaturation(phaseIdx, S/sumSat);
366 }
367 sumSat = 1;
368 }
369
370 // set the last saturation
371 fluidState.setSaturation(/*phaseIdx=*/numPhases - 1, 1.0 - sumSat);
372
373 // update the pressures using the material law (saturations
374 // and first pressure are already set because it is implicitly
375 // solved for.)
376 const auto pc = material.capillaryPressures(fluidState, wettingPhaseIdx_);
377 for (int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
378 fluidState.setPressure(phaseIdx,
379 fluidState.pressure(0)
380 + (pc[phaseIdx] - pc[0]));
381
382 // update the parameter cache
383 paramCache.updateAll(fluidState, /*except=*/ParameterCache::Temperature|ParameterCache::Composition);
384
385 // update all densities
386 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
387 Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
388 Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
389 fluidState.setDensity(phaseIdx, rho);
390 fluidState.setMolarDensity(phaseIdx, rhoMolar);
391 }
392 }
393
394 bool isPressureIdx_(int pvIdx)
395 { return pvIdx == 0; }
396
397 bool isSaturationIdx_(int pvIdx)
398 { return 1 <= pvIdx && pvIdx < numPhases; }
399
400 // retrieves a quantity from the fluid state
401 template <class FluidState>
402 Scalar getQuantity_(const FluidState &fs, int pvIdx)
403 {
404 assert(pvIdx < numEq);
405
406 // first pressure
407 if (pvIdx < 1) {
408 int phaseIdx = 0;
409 return fs.pressure(phaseIdx);
410 }
411 // saturations
412 else {
413 assert(pvIdx < numPhases);
414 int phaseIdx = pvIdx - 1;
415 return fs.saturation(phaseIdx);
416 }
417 }
418
419 // set a quantity in the fluid state
420 template <class MaterialLaw, class FluidState>
421 void setQuantity_(FluidState &fs,
422 ParameterCache &paramCache,
423 const MaterialLaw& material,
424 int pvIdx,
425 Scalar value)
426 {
427 assert(0 <= pvIdx && pvIdx < numEq);
428
429 if (pvIdx < 1) {
430 // -> first pressure
431 Scalar delta = value - fs.pressure(0);
432
433 // set all pressures. here we assume that the capillary
434 // pressure does not depend on absolute pressure.
435 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
436 fs.setPressure(phaseIdx, fs.pressure(phaseIdx) + delta);
437 paramCache.updateAllPressures(fs);
438
439 // update all densities
440 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
441 Scalar rho = FluidSystem::density(fs, paramCache, phaseIdx);
442 Scalar rhoMolar = FluidSystem::molarDensity(fs, paramCache, phaseIdx);
443 fs.setDensity(phaseIdx, rho);
444 fs.setMolarDensity(phaseIdx, rhoMolar);
445 }
446 }
447 else {
448 assert(pvIdx < numPhases);
449
450 // -> first M-1 saturations
451 Scalar delta = value - fs.saturation(/*phaseIdx=*/pvIdx - 1);
452 fs.setSaturation(/*phaseIdx=*/pvIdx - 1, value);
453
454 // set last saturation (-> minus the change of the
455 // saturation of the other phases)
456 fs.setSaturation(/*phaseIdx=*/numPhases - 1,
457 fs.saturation(numPhases - 1) - delta);
458
459 // update all fluid pressures using the capillary pressure
460 // law
461 const auto pc = material.capillaryPressures(fs, wettingPhaseIdx_);
462 for (int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
463 fs.setPressure(phaseIdx,
464 fs.pressure(0)
465 + (pc[phaseIdx] - pc[0]));
466 paramCache.updateAllPressures(fs);
467
468 // update all densities
469 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
470 Scalar rho = FluidSystem::density(fs, paramCache, phaseIdx);
471 Scalar rhoMolar = FluidSystem::molarDensity(fs, paramCache, phaseIdx);
472 fs.setDensity(phaseIdx, rho);
473 fs.setMolarDensity(phaseIdx, rhoMolar);
474
475 }
476 }
477 }
478
479 // set a quantity in the fluid state
480 template <class FluidState>
481 void setQuantityRaw_(FluidState &fs, int pvIdx, Scalar value)
482 {
483 assert(pvIdx < numEq);
484
485 // first pressure
486 if (pvIdx < 1) {
487 int phaseIdx = 0;
488 fs.setPressure(phaseIdx, value);
489 }
490 // saturations
491 else {
492 assert(pvIdx < numPhases);
493 int phaseIdx = pvIdx - 1;
494
495 // make sure that the first M-1 saturations does not get
496 // negative
497 using std::max;
498 value = max(0.0, value);
499 fs.setSaturation(phaseIdx, value);
500 }
501 }
502
503 template <class FluidState>
504 Scalar quantityWeight_(const FluidState &fs, int pvIdx)
505 {
506 // first pressure
507 if (pvIdx < 1)
508 return 1.0/1e5;
509 // first M - 1 saturations
510 else {
511 assert(pvIdx < numPhases);
512 return 1.0;
513 }
514 }
515
516private:
517 int wettingPhaseIdx_ = 0;
518
519};
520
521} // end namespace Dumux
522
523#endif
Determines the pressures and saturations of all fluid phases given the total mass of all components.
Definition: immiscibleflash.hh:53
void guessInitial(FluidState &fluidState, ParameterCache &paramCache, const ComponentVector &globalMolarities)
Guess initial values for all quantities.
Definition: immiscibleflash.hh:85
void printFluidState_(const FluidState &fs)
Definition: immiscibleflash.hh:222
void solve(FluidState &fluidState, ParameterCache &paramCache, const MaterialLaw &material, const ComponentVector &globalMolarities)
Calculates the chemical equilibrium from the component fugacities in a phase.
Definition: immiscibleflash.hh:114
bool isPressureIdx_(int pvIdx)
Definition: immiscibleflash.hh:394
void setQuantityRaw_(FluidState &fs, int pvIdx, Scalar value)
Definition: immiscibleflash.hh:481
void calculateDefect_(Vector &b, const FluidState &fluidStateEval, const FluidState &fluidState, const ComponentVector &globalMolarities)
Definition: immiscibleflash.hh:299
Scalar getQuantity_(const FluidState &fs, int pvIdx)
Definition: immiscibleflash.hh:402
void completeFluidState_(FluidState &fluidState, ParameterCache &paramCache, const MaterialLaw &material)
Definition: immiscibleflash.hh:349
void linearize_(Matrix &J, Vector &b, FluidState &fluidState, ParameterCache &paramCache, const MaterialLaw &material, const ComponentVector &globalMolarities)
Definition: immiscibleflash.hh:251
Scalar update_(FluidState &fluidState, ParameterCache &paramCache, const MaterialLaw &material, const Vector &deltaX)
Definition: immiscibleflash.hh:315
void setQuantity_(FluidState &fs, ParameterCache &paramCache, const MaterialLaw &material, int pvIdx, Scalar value)
Definition: immiscibleflash.hh:421
bool isSaturationIdx_(int pvIdx)
Definition: immiscibleflash.hh:397
Scalar quantityWeight_(const FluidState &fs, int pvIdx)
Definition: immiscibleflash.hh:504
ImmiscibleFlash(int wettingPhaseIdx=0)
Construct a new flash.
Definition: immiscibleflash.hh:74
Dune::FieldVector< Scalar, numComponents > ComponentVector
Definition: immiscibleflash.hh:68
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:27
Some exceptions thrown in DuMux
Makes the capillary pressure-saturation relations available under the M-phase API for material laws.
std::string saturation(int phaseIdx) noexcept
I/O name of saturation for multiphase systems.
Definition: name.hh:31
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