version 3.10-dev
ncpflash.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_NCP_FLASH_HH
14#define DUMUX_NCP_FLASH_HH
15
16#include <dune/common/fvector.hh>
17#include <dune/common/fmatrix.hh>
18
21
22namespace Dumux {
23
65template <class Scalar, class FluidSystem>
67{
68 enum { numPhases = FluidSystem::numPhases };
69 enum { numComponents = FluidSystem::numComponents };
70
71 using ParameterCache = typename FluidSystem::ParameterCache;
72
73 static constexpr int numEq = numPhases*(numComponents + 1);
74
75 using Matrix = Dune::FieldMatrix<Scalar, numEq, numEq>;
76 using Vector = Dune::FieldVector<Scalar, numEq>;
77
78public:
79 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
80
85 explicit NcpFlash(int wettingPhaseIdx = 0)
86 : wettingPhaseIdx_(wettingPhaseIdx)
87 {}
88
95 template <class FluidState>
96 void guessInitial(FluidState &fluidState,
97 ParameterCache &paramCache,
98 const ComponentVector &globalMolarities)
99 {
100 // the sum of all molarities
101 Scalar sumMoles = 0;
102 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
103 sumMoles += globalMolarities[compIdx];
104
105 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
106 // composition
107 for (int compIdx = 0; compIdx < numComponents; ++ compIdx)
108 fluidState.setMoleFraction(phaseIdx,
109 compIdx,
110 globalMolarities[compIdx]/sumMoles);
111
112 // pressure. use atmospheric pressure as initial guess
113 fluidState.setPressure(phaseIdx, 1.0135e5);
114
115 // saturation. assume all fluids to be equally distributed
116 fluidState.setSaturation(phaseIdx, 1.0/numPhases);
117 }
118
119 // set the fugacity coefficients of all components in all phases
120 paramCache.updateAll(fluidState);
121 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
122 for (int compIdx = 0; compIdx < numComponents; ++ compIdx) {
123 Scalar phi = FluidSystem::fugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx);
124 fluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
125 }
126 }
127 }
128
139 template <class MaterialLaw, class FluidState>
140 void solve(FluidState &fluidState,
141 ParameterCache &paramCache,
142 const MaterialLaw& material,
143 const ComponentVector &globalMolarities)
144 {
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 // make the fluid state consistent with the fluid system.
157 completeFluidState_(fluidState, paramCache, material);
158
159 /*
160 std::cout << "--------------------\n";
161 std::cout << "globalMolarities: ";
162 for (int compIdx = 0; compIdx < numComponents; ++ compIdx)
163 std::cout << globalMolarities[compIdx] << " ";
164 std::cout << "\n";
165 */
166
167 const int nMax = 50; // <- maximum number of newton iterations
168 for (int nIdx = 0; nIdx < nMax; ++nIdx) {
169 // calculate Jacobian matrix and right hand side
170 linearize_(J, b, fluidState, paramCache, material, globalMolarities);
171
172 // Solve J*x = b
173 deltaX = 0;
174
175 try {
176
177 // simple preconditioning of the linear system to reduce condition number
178 int eqIdx = 0;
179 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
180 {
181 for (int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
182 {
183 J[eqIdx] *= 10e-7; // roughly the magnitude of the fugacities
184 b[eqIdx] *= 10e-7;
185 ++eqIdx;
186 }
187 }
188
189 J.solve(deltaX, b);
190
191 }
192 catch (Dune::FMatrixError &e)
193 {
194 /*
195 printFluidState_(fluidState);
196 std::cout << "error: " << e << "\n";
197 std::cout << "b: " << b << "\n";
198 std::cout << "J: " << J << "\n";
199 */
200
201 throw NumericalProblem(e.what());
202 }
203
204 /*
205 printFluidState_(fluidState);
206 std::cout << "J:\n";
207 for (int i = 0; i < numEq; ++i) {
208 for (int j = 0; j < numEq; ++j) {
209 std::ostringstream os;
210 os << J[i][j];
211
212 std::string s(os.str());
213 do {
214 s += " ";
215 } while (s.size() < 20);
216 std::cout << s;
217 }
218 std::cout << "\n";
219 }
220
221 std::cout << "deltaX: " << deltaX << "\n";
222 std::cout << "---------------\n";
223 */
224
225 // update the fluid quantities.
226 Scalar relError = update_(fluidState, paramCache, material, deltaX);
227
228 if (relError < 1e-9)
229 return;
230 }
231
232 /*
233 printFluidState_(fluidState);
234 std::cout << "globalMolarities: ";
235 for (int compIdx = 0; compIdx < numComponents; ++ compIdx)
236 std::cout << globalMolarities[compIdx] << " ";
237 std::cout << "\n";
238 */
239
240 DUNE_THROW(NumericalProblem,
241 "Flash calculation failed."
242 " {c_alpha^kappa} = {" << globalMolarities << "}, T = "
243 << fluidState.temperature(/*phaseIdx=*/0));
244 }
245
246
247protected:
248 template <class FluidState>
249 void printFluidState_(const FluidState &fs)
250 {
251 std::cout << "saturations: ";
252 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
253 std::cout << fs.saturation(phaseIdx) << " ";
254 std::cout << "\n";
255
256 std::cout << "pressures: ";
257 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
258 std::cout << fs.pressure(phaseIdx) << " ";
259 std::cout << "\n";
260
261 std::cout << "densities: ";
262 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
263 std::cout << fs.density(phaseIdx) << " ";
264 std::cout << "\n";
265
266 std::cout << "molar densities: ";
267 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
268 std::cout << fs.molarDensity(phaseIdx) << " ";
269 std::cout << "\n";
270
271 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
272 std::cout << "composition " << FluidSystem::phaseName(phaseIdx) << "Phase: ";
273 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
274 std::cout << fs.moleFraction(phaseIdx, compIdx) << " ";
275 }
276 std::cout << "\n";
277 }
278
279 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
280 std::cout << "fugacities " << FluidSystem::phaseName(phaseIdx) << "Phase: ";
281 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
282 std::cout << fs.fugacity(phaseIdx, compIdx) << " ";
283 }
284 std::cout << "\n";
285 }
286
287 std::cout << "global component molarities: ";
288 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
289 Scalar sum = 0;
290 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
291 sum += fs.saturation(phaseIdx)*fs.molarity(phaseIdx, compIdx);
292 }
293 std::cout << sum << " ";
294 }
295 std::cout << "\n";
296 }
297
298 template <class MaterialLaw, class FluidState>
299 void linearize_(Matrix &J,
300 Vector &b,
301 FluidState &fluidState,
302 ParameterCache &paramCache,
303 const MaterialLaw& material,
304 const ComponentVector &globalMolarities)
305 {
306 FluidState origFluidState(fluidState);
307 ParameterCache origParamCache(paramCache);
308
309 Vector tmp;
310
311 // reset jacobian
312 J = 0;
313
314 calculateDefect_(b, fluidState, fluidState, globalMolarities);
315
316 // assemble jacobian matrix
317 for (int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
319 // approximately calculate partial derivatives of the
320 // fugacity defect of all components in regard to the mole
321 // fraction of the i-th component. This is done via
322 // forward differences
323
324 // deviate the mole fraction of the i-th component
325 Scalar x_i = getQuantity_(fluidState, pvIdx);
326 const Scalar eps = 1e-8/quantityWeight_(fluidState, pvIdx);
327 setQuantity_(fluidState, paramCache, material, pvIdx, x_i + eps);
328
329 // compute derivative of the defect
330 calculateDefect_(tmp, origFluidState, fluidState, globalMolarities);
331 tmp -= b;
332 tmp /= eps;
333
334 // store derivative in jacobian matrix
335 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
336 J[eqIdx][pvIdx] = tmp[eqIdx];
337
338 // fluid state and parameter cache to their original values
339 fluidState = origFluidState;
340 paramCache = origParamCache;
341
342 // end forward differences
344 }
345 }
346
347 template <class FluidState>
348 void calculateDefect_(Vector &b,
349 const FluidState &fluidStateEval,
350 const FluidState &fluidState,
351 const ComponentVector &globalMolarities)
352 {
353 int eqIdx = 0;
354
355 // fugacity of any component must be equal in all phases
356 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
357 for (int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
358 b[eqIdx] =
359 fluidState.fugacity(/*phaseIdx=*/0, compIdx) -
360 fluidState.fugacity(phaseIdx, compIdx);
361 ++eqIdx;
362 }
363 }
364
365 assert(eqIdx == numComponents*(numPhases - 1));
366
367 // the fact saturations must sum up to 1 is included explicitly!
368
369 // capillary pressures are explicitly included!
370
371 // global molarities are given
372 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
373 b[eqIdx] = 0.0;
374 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
375 b[eqIdx] +=
376 fluidState.saturation(phaseIdx)
377 * fluidState.molarity(phaseIdx, compIdx);
378 }
379
380 b[eqIdx] -= globalMolarities[compIdx];
381 ++eqIdx;
382 }
383
384 // model assumptions (-> non-linear complementarity functions)
385 // must be adhered
386 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
387 Scalar sumMoleFracEval = 0.0;
388 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
389 sumMoleFracEval += fluidStateEval.moleFraction(phaseIdx, compIdx);
390
391 if (1.0 - sumMoleFracEval > fluidStateEval.saturation(phaseIdx)) {
392 b[eqIdx] = fluidState.saturation(phaseIdx);
393 }
394 else {
395 Scalar sumMoleFrac = 0.0;
396 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
397 sumMoleFrac += fluidState.moleFraction(phaseIdx, compIdx);
398 b[eqIdx] = 1.0 - sumMoleFrac;
399 }
400
401 ++eqIdx;
402 }
403 }
404
405 template <class MaterialLaw, class FluidState>
406 Scalar update_(FluidState &fluidState,
407 ParameterCache &paramCache,
408 const MaterialLaw& material,
409 const Vector &deltaX)
410 {
411 Scalar relError = 0;
412 for (int pvIdx = 0; pvIdx < numEq; ++ pvIdx) {
413 Scalar tmp = getQuantity_(fluidState, pvIdx);
414 Scalar delta = deltaX[pvIdx];
415 using std::max;
416 using std::abs;
417 using std::min;
418 relError = max(relError, abs(delta)*quantityWeight_(fluidState, pvIdx));
419
420 if (isSaturationIdx_(pvIdx)) {
421 // dampen to at most 20% change in saturation per
422 // iteration
423 delta = min(0.2, max(-0.2, delta));
424 }
425 else if (isMoleFracIdx_(pvIdx)) {
426 // dampen to at most 15% change in mole fraction per
427 // iteration
428 delta = min(0.15, max(-0.15, delta));
429 }
430 else if (isPressureIdx_(pvIdx)) {
431 // dampen to at most 15% change in pressure per
432 // iteration
433 delta = min(0.15*fluidState.pressure(0), max(-0.15*fluidState.pressure(0), delta));
434 }
435
436 setQuantityRaw_(fluidState, pvIdx, tmp - delta);
437 }
438
439 // make sure all saturations, pressures and mole fractions are non-negative
440 Scalar sumSat = 0;
441 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
442 Scalar value = fluidState.saturation(phaseIdx);
443 if (value < -0.05) {
444 value = -0.05;
445 fluidState.setSaturation(phaseIdx, value);
446 }
447 sumSat += value;
448
449 value = fluidState.pressure(phaseIdx);
450 if (value < 0)
451 fluidState.setPressure(phaseIdx, 0.0);
452
453 bool allMoleFractionsAreZero = true;
454 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
455 value = fluidState.moleFraction(phaseIdx, compIdx);
456 if (value > 0)
457 allMoleFractionsAreZero = false;
458 else if (value < 0)
459 fluidState.setMoleFraction(phaseIdx, compIdx, 0.0);
460 }
461
462 // set at least one mole fraction to a positive value
463 // to prevent breakdown of the linear solver in completeFluidState_
464 if (allMoleFractionsAreZero)
465 fluidState.setMoleFraction(phaseIdx, /*compIdx=*/0, 1e-10);
466 }
467
468 // last saturation
469 if (sumSat > 1.05) {
470 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
471 Scalar value = fluidState.saturation(phaseIdx)/(0.95*sumSat);
472 fluidState.setSaturation(phaseIdx, value);
473 }
474 }
475
476 completeFluidState_(fluidState, paramCache, material);
477
478 return relError;
479 }
480
481 template <class MaterialLaw, class FluidState>
482 void completeFluidState_(FluidState &fluidState,
483 ParameterCache &paramCache,
484 const MaterialLaw& material)
485 {
486 // calculate the saturation of the last phase as a function of
487 // the other saturations
488 Scalar sumSat = 0.0;
489 for (int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
490 sumSat += fluidState.saturation(phaseIdx);
491 fluidState.setSaturation(/*phaseIdx=*/numPhases - 1, 1.0 - sumSat);
492
493 // update the pressures using the material law (saturations
494 // and first pressure are already set because it is implicitly
495 // solved for.)
496 const auto pc = material.capillaryPressures(fluidState, wettingPhaseIdx_);
497 for (int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
498 fluidState.setPressure(phaseIdx,
499 fluidState.pressure(0)
500 + (pc[phaseIdx] - pc[0]));
501
502 // update the parameter cache
503 paramCache.updateAll(fluidState, /*except=*/ParameterCache::Temperature);
504
505 // update all densities and fugacity coefficients
506 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
507
508 Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
509 Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, phaseIdx);
510 fluidState.setDensity(phaseIdx, rho);
511 fluidState.setMolarDensity(phaseIdx, rhoMolar);
512
513 for (int compIdx = 0; compIdx < numComponents; ++ compIdx) {
514 Scalar phi = FluidSystem::fugacityCoefficient( fluidState, paramCache, phaseIdx, compIdx);
515 fluidState.setFugacityCoefficient(phaseIdx, compIdx, phi);
516 }
517 }
518 }
519
520 bool isPressureIdx_(int pvIdx)
521 { return pvIdx == 0; }
522
523 bool isSaturationIdx_(int pvIdx)
524 { return 1 <= pvIdx && pvIdx < numPhases; }
525
526 bool isMoleFracIdx_(int pvIdx)
527 { return numPhases <= pvIdx && pvIdx < numPhases + numPhases*numComponents; }
528
529 // retrieves a quantity from the fluid state
530 template <class FluidState>
531 Scalar getQuantity_(const FluidState &fs, int pvIdx)
532 {
533 assert(pvIdx < numEq);
534
535 // first pressure
536 if (pvIdx < 1) {
537 int phaseIdx = 0;
538 return fs.pressure(phaseIdx);
539 }
540 // first M - 1 saturations
541 else if (pvIdx < numPhases) {
542 int phaseIdx = pvIdx - 1;
543 return fs.saturation(phaseIdx);
544 }
545 // mole fractions
546 else // if (pvIdx < numPhases + numPhases*numComponents)
547 {
548 int phaseIdx = (pvIdx - numPhases)/numComponents;
549 int compIdx = (pvIdx - numPhases)%numComponents;
550 return fs.moleFraction(phaseIdx, compIdx);
551 }
552 }
553
554 // set a quantity in the fluid state
555 template <class MaterialLaw, class FluidState>
556 void setQuantity_(FluidState &fs,
557 ParameterCache &paramCache,
558 const MaterialLaw& material,
559 int pvIdx,
560 Scalar value)
561 {
562 assert(0 <= pvIdx && pvIdx < numEq);
563
564 if (pvIdx < 1) { // <- first pressure
565 Scalar delta = value - fs.pressure(0);
566
567 // set all pressures. here we assume that the capillary
568 // pressure does not depend on absolute pressure.
569 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
570 fs.setPressure(phaseIdx, fs.pressure(phaseIdx) + delta);
571 paramCache.updateAllPressures(fs);
572
573 // update all densities and fugacity coefficients
574 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
575 Scalar rho = FluidSystem::density(fs, paramCache, phaseIdx);
576 Scalar rhoMolar = FluidSystem::molarDensity(fs, paramCache, phaseIdx);
577 fs.setDensity(phaseIdx, rho);
578 fs.setMolarDensity(phaseIdx, rhoMolar);
579
580 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
581 Scalar phi = FluidSystem::fugacityCoefficient(fs, paramCache, phaseIdx, compIdx);
582 fs.setFugacityCoefficient(phaseIdx, compIdx, phi);
583 }
584 }
585 }
586 else if (pvIdx < numPhases) { // <- first M - 1 saturations
587 Scalar delta = value - fs.saturation(/*phaseIdx=*/pvIdx - 1);
588 fs.setSaturation(/*phaseIdx=*/pvIdx - 1, value);
589
590 // set last saturation (-> minus the change of the
591 // satuation of the other phase)
592 fs.setSaturation(/*phaseIdx=*/numPhases - 1,
593 fs.saturation(numPhases - 1) - delta);
594
595 // update all fluid pressures using the capillary pressure
596 // law
597 const auto pc = material.capillaryPressures(fs, wettingPhaseIdx_);
598 for (int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx)
599 fs.setPressure(phaseIdx,
600 fs.pressure(0)
601 + (pc[phaseIdx] - pc[0]));
602 paramCache.updateAllPressures(fs);
603
604 // update all densities and fugacity coefficients
605 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
606 Scalar rho = FluidSystem::density(fs, paramCache, phaseIdx);
607 Scalar rhoMolar = FluidSystem::molarDensity(fs, paramCache, phaseIdx);
608 fs.setDensity(phaseIdx, rho);
609 fs.setMolarDensity(phaseIdx, rhoMolar);
610
611 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
612 Scalar phi = FluidSystem::fugacityCoefficient(fs, paramCache, phaseIdx, compIdx);
613 fs.setFugacityCoefficient(phaseIdx, compIdx, phi);
614 }
615 }
616 }
617 else if (pvIdx < numPhases + numPhases*numComponents) // <- mole fractions
618 {
619 int phaseIdx = (pvIdx - numPhases)/numComponents;
620 int compIdx = (pvIdx - numPhases)%numComponents;
621
622 fs.setMoleFraction(phaseIdx, compIdx, value);
623 paramCache.updateSingleMoleFraction(fs, phaseIdx, compIdx);
624
625 // update the density of the phase
626 Scalar rho = FluidSystem::density(fs, paramCache, phaseIdx);
627 Scalar rhoMolar = FluidSystem::molarDensity(fs, paramCache, phaseIdx);
628 fs.setDensity(phaseIdx, rho);
629 fs.setMolarDensity(phaseIdx, rhoMolar);
630
631 // if the phase's fugacity coefficients are composition
632 // dependent, update them as well.
633 if (!FluidSystem::isIdealMixture(phaseIdx)) {
634 for (int compIdx2 = 0; compIdx2 < numComponents; ++compIdx2) {
635 Scalar phi = FluidSystem::fugacityCoefficient(fs, paramCache, phaseIdx, compIdx2);
636 fs.setFugacityCoefficient(phaseIdx, compIdx2, phi);
637 }
638 }
639 }
640 else {
641 assert(false);
642 }
643 }
644
645 // set a quantity in the fluid state
646 template <class FluidState>
647 void setQuantityRaw_(FluidState &fs, int pvIdx, Scalar value)
648 {
649 assert(pvIdx < numEq);
650
651 // first pressure
652 if (pvIdx < 1) {
653 int phaseIdx = 0;
654 fs.setPressure(phaseIdx, value);
655 }
656 // first M - 1 saturations
657 else if (pvIdx < numPhases) {
658 int phaseIdx = pvIdx - 1;
659 fs.setSaturation(phaseIdx, value);
660 }
661 // mole fractions
662 else // if (pvIdx < numPhases + numPhases*numComponents)
663 {
664 int phaseIdx = (pvIdx - numPhases)/numComponents;
665 int compIdx = (pvIdx - numPhases)%numComponents;
666 fs.setMoleFraction(phaseIdx, compIdx, value);
667 }
668 }
669
670 template <class FluidState>
671 Scalar quantityWeight_(const FluidState &fs, int pvIdx)
672 {
673 // first pressure
674 if (pvIdx < 1)
675 return 1.0/1e5;
676 // first M - 1 saturations
677 else if (pvIdx < numPhases)
678 return 1.0;
679 // mole fractions
680 else // if (pvIdx < numPhases + numPhases*numComponents)
681 return 1.0;
682 }
683
684private:
685 int wettingPhaseIdx_ = 0;
686};
687
688} // end namespace Dumux
689
690#endif
Determines the phase compositions, pressures and saturations given the total mass of all components.
Definition: ncpflash.hh:67
Dune::FieldVector< Scalar, numComponents > ComponentVector
Definition: ncpflash.hh:79
NcpFlash(int wettingPhaseIdx=0)
Construct a new flash.
Definition: ncpflash.hh:85
bool isSaturationIdx_(int pvIdx)
Definition: ncpflash.hh:523
Scalar update_(FluidState &fluidState, ParameterCache &paramCache, const MaterialLaw &material, const Vector &deltaX)
Definition: ncpflash.hh:406
void completeFluidState_(FluidState &fluidState, ParameterCache &paramCache, const MaterialLaw &material)
Definition: ncpflash.hh:482
bool isPressureIdx_(int pvIdx)
Definition: ncpflash.hh:520
Scalar getQuantity_(const FluidState &fs, int pvIdx)
Definition: ncpflash.hh:531
void linearize_(Matrix &J, Vector &b, FluidState &fluidState, ParameterCache &paramCache, const MaterialLaw &material, const ComponentVector &globalMolarities)
Definition: ncpflash.hh:299
void guessInitial(FluidState &fluidState, ParameterCache &paramCache, const ComponentVector &globalMolarities)
Guess initial values for all quantities.
Definition: ncpflash.hh:96
void solve(FluidState &fluidState, ParameterCache &paramCache, const MaterialLaw &material, const ComponentVector &globalMolarities)
Calculates the chemical equilibrium from the component fugacities in a phase.
Definition: ncpflash.hh:140
void setQuantityRaw_(FluidState &fs, int pvIdx, Scalar value)
Definition: ncpflash.hh:647
Scalar quantityWeight_(const FluidState &fs, int pvIdx)
Definition: ncpflash.hh:671
bool isMoleFracIdx_(int pvIdx)
Definition: ncpflash.hh:526
void setQuantity_(FluidState &fs, ParameterCache &paramCache, const MaterialLaw &material, int pvIdx, Scalar value)
Definition: ncpflash.hh:556
void printFluidState_(const FluidState &fs)
Definition: ncpflash.hh:249
void calculateDefect_(Vector &b, const FluidState &fluidStateEval, const FluidState &fluidState, const ComponentVector &globalMolarities)
Definition: ncpflash.hh:348
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 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