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