66 static constexpr int numPhases = FluidSystem::numPhases;
67 static constexpr int numComponents = FluidSystem::numComponents;
68 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx};
69 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx};
70 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx};
71 static constexpr int numMiscibleComponents = FluidSystem::numMiscibleComponents;
72 static constexpr int numMisciblePhases = FluidSystem::numMisciblePhases;
73 static constexpr int numEq = numMisciblePhases + numMisciblePhases * numMiscibleComponents;
75 using EOSType = CompositionalConfig::EOSType;
82 template <
class Flu
idState>
83 static bool solve(FluidState& fluid_state,
84 const std::string& twoPhaseMethod,
85 Scalar flash_tolerance,
86 const EOSType& eos_type,
90 ScalarFluidState fluid_state_scalar;
92 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
93 fluid_state_scalar.setKvalue(compIdx, Opm::getValue(fluid_state.K(compIdx) ) );
94 fluid_state_scalar.setMoleFraction(compIdx, Opm::getValue(fluid_state.moleFraction(compIdx) ) );
97 fluid_state_scalar.setLvalue(Opm::getValue(fluid_state.L()));
99 fluid_state_scalar.setPressure(FluidSystem::oilPhaseIdx,
100 Opm::getValue(fluid_state.pressure(FluidSystem::oilPhaseIdx)));
101 fluid_state_scalar.setPressure(FluidSystem::gasPhaseIdx,
102 Opm::getValue(fluid_state.pressure(FluidSystem::gasPhaseIdx)));
104 fluid_state_scalar.setTemperature(Opm::getValue(fluid_state.temperature(0)));
106 const auto is_single_phase = flash_solve_scalar_(fluid_state_scalar, twoPhaseMethod, flash_tolerance, eos_type, verbosity);
110 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
111 const auto x_i = fluid_state_scalar.moleFraction(oilPhaseIdx, compIdx);
112 fluid_state.setMoleFraction(oilPhaseIdx, compIdx, x_i);
113 const auto y_i = fluid_state_scalar.moleFraction(gasPhaseIdx, compIdx);
114 fluid_state.setMoleFraction(gasPhaseIdx, compIdx, y_i);
118 updateDerivatives_(fluid_state_scalar, fluid_state, eos_type, is_single_phase);
120 return is_single_phase;
130 template <
class Flu
idState,
class ComponentVector>
131 static void solve(FluidState& fluid_state,
132 const ComponentVector& globalMolarities,
133 Scalar tolerance = 0.0)
137 using MaterialLawParams =
typename MaterialLaw::Params;
139 MaterialLawParams matParams;
140 solve<MaterialLaw>(fluid_state, matParams, globalMolarities, tolerance);
143 template <
class Vector>
144 static typename Vector::field_type solveRachfordRice_g_(
const Vector& K,
const Vector& z,
int verbosity)
148 using field_type =
typename Vector::field_type;
149 constexpr field_type tol = 1e-12;
150 constexpr int itmax = 10000;
151 field_type Kmin = K[0];
152 field_type Kmax = K[0];
153 for (
int compIdx = 1; compIdx < numComponents; ++compIdx){
154 if (K[compIdx] < Kmin)
156 else if (K[compIdx] >= Kmax)
160 auto Vmin = 1 / (1 - Kmax);
161 auto Vmax = 1 / (1 - Kmin);
163 auto V = (Vmin + Vmax)/2;
165 if (verbosity == 3 || verbosity == 4) {
166 std::cout <<
"Initial guess " << numComponents <<
"c : V = " << V <<
" and [Vmin, Vmax] = [" << Vmin <<
", " << Vmax <<
"]" << std::endl;
167 std::cout << std::setw(10) <<
"Iteration" << std::setw(16) <<
"abs(step)" << std::setw(16) <<
"V" << std::endl;
170 for (
int iteration = 1; iteration < itmax; ++iteration) {
172 field_type denum = 0.0;
174 for (
int compIdx = 0; compIdx < numComponents; ++compIdx){
175 auto dK = K[compIdx] - 1.0;
176 auto a = z[compIdx] * dK;
177 auto b = (1 + V * dK);
179 denum += z[compIdx] * (dK*dK) / (b*b);
181 auto delta = r / denum;
185 if (V < Vmin || V > Vmax)
188 if (verbosity == 3 || verbosity == 4) {
189 std::cout <<
"V = " << V <<
" is not within the the range [Vmin, Vmax], solve using Bisection method!" << std::endl;
197 decltype(Vmax) Lmin = 1.0;
198 decltype(Vmin) Lmax = 0.0;
199 auto L = bisection_g_(K, Lmin, Lmax, z, verbosity);
202 if (verbosity >= 1) {
203 std::cout <<
"Rachford-Rice (Bisection) converged to final solution L = " << L << std::endl;
209 if (verbosity == 3 || verbosity == 4) {
210 std::cout << std::setw(10) << iteration << std::setw(16) << Opm::abs(delta) << std::setw(16) << V << std::endl;
213 if ( Opm::abs(r) < tol ) {
218 if (verbosity >= 1) {
219 std::cout <<
"Rachford-Rice converged to final solution L = " << L << std::endl;
226 throw std::runtime_error(
" Rachford-Rice did not converge within maximum number of iterations" );
230 template <
typename Flu
idState>
231 static bool flash_solve_scalar_(FluidState& fluid_state,
232 const std::string& twoPhaseMethod,
233 const Scalar flash_tolerance,
234 const EOSType& eos_type,
235 const int verbosity = 0)
238 bool is_stable =
false;
239 auto L_scalar = fluid_state.L();
240 using ScalarVector = Dune::FieldVector<Scalar, numComponents>;
241 ScalarVector K_scalar, z_scalar;
242 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
243 K_scalar[compIdx] = fluid_state.K(compIdx);
244 z_scalar[compIdx] = fluid_state.moleFraction(compIdx);
247 if ( L_scalar <= 0 || L_scalar == 1 ) {
248 if (verbosity >= 1) {
249 std::cout <<
"Perform stability test (L <= 0 or L == 1)!" << std::endl;
251 phaseStabilityTest_(is_stable, K_scalar, fluid_state, z_scalar, eos_type, verbosity);
253 if (verbosity >= 1) {
254 std::cout <<
"Inputs after stability test are K = [" << K_scalar <<
"], L = [" << L_scalar <<
"], z = [" << z_scalar <<
"], P = " << fluid_state.pressure(0) <<
", and T = " << fluid_state.temperature(0) << std::endl;
258 const bool is_single_phase = is_stable;
261 if ( !is_single_phase ) {
263 L_scalar = solveRachfordRice_g_(K_scalar, z_scalar, verbosity);
264 flash_2ph(z_scalar, twoPhaseMethod, K_scalar, L_scalar, fluid_state, flash_tolerance, eos_type, verbosity);
267 L_scalar = li_single_phase_label_(fluid_state, z_scalar, verbosity);
269 fluid_state.setLvalue(L_scalar);
270 return is_single_phase;
273 template <
class Vector>
274 static typename Vector::field_type bisection_g_(
const Vector& K,
typename Vector::field_type Lmin,
275 typename Vector::field_type Lmax,
const Vector& z,
int verbosity)
278 typename Vector::field_type gLmin = rachfordRice_g_(K, Lmin, z);
281 if (verbosity >= 3) {
282 std::cout << std::setw(10) <<
"Iteration" << std::setw(16) <<
"g(Lmid)" << std::setw(16) <<
"L" << std::endl;
285 constexpr int max_it = 10000;
287 auto closeLmaxLmin = [](
double max_v,
double min_v) {
288 return Opm::abs(max_v - min_v) / 2. < 1e-10;
293 if (closeLmaxLmin(Lmax, Lmin) ){
294 throw std::runtime_error(fmt::format(
"Strange bisection with Lmax {} and Lmin {}?", Lmax, Lmin));
296 for (
int iteration = 0; iteration < max_it; ++iteration){
298 auto L = (Lmin + Lmax) / 2;
299 auto gMid = rachfordRice_g_(K, L, z);
301 if (verbosity == 3 || verbosity == 4) {
302 std::cout << std::setw(10) << iteration << std::setw(16) << gMid << std::setw(16) << L << std::endl;
306 if (Opm::abs(gMid) < 1e-16 || closeLmaxLmin(Lmax, Lmin)){
310 else if (Dune::sign(gMid) != Dune::sign(gLmin)) {
320 throw std::runtime_error(
" Rachford-Rice bisection failed with " + std::to_string(max_it) +
" iterations!");
323 template <
class Vector,
class FlashFlu
idState>
324 static typename Vector::field_type li_single_phase_label_(
const FlashFluidState& fluid_state,
const Vector& z,
int verbosity)
327 typename Vector::field_type sumVz = 0.0;
328 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
330 const auto& V_crit = FluidSystem::criticalVolume(compIdx);
333 sumVz += (V_crit * z[compIdx]);
337 typename Vector::field_type Tc_est = 0.0;
338 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
340 const auto& V_crit = FluidSystem::criticalVolume(compIdx);
341 const auto& T_crit = FluidSystem::criticalTemperature(compIdx);
344 Tc_est += (V_crit * T_crit * z[compIdx] / sumVz);
348 const auto& T = fluid_state.temperature(0);
351 typename Vector::field_type L;
357 if (verbosity >= 1) {
358 std::cout <<
"Cell is single-phase, liquid (L = 1.0) due to Li's phase labeling method giving T < Tc_est (" << T <<
" < " << Tc_est <<
")!" << std::endl;
366 if (verbosity >= 1) {
367 std::cout <<
"Cell is single-phase, vapor (L = 0.0) due to Li's phase labeling method giving T >= Tc_est (" << T <<
" >= " << Tc_est <<
")!" << std::endl;
375 template <
class FlashFlu
idState,
class ComponentVector>
376 static void phaseStabilityTest_(
bool& isStable, ComponentVector& K, FlashFluidState& fluid_state,
const ComponentVector& z,
const EOSType& eos_type,
int verbosity)
379 bool isTrivialL, isTrivialV;
380 ComponentVector x, y;
381 typename FlashFluidState::Scalar S_l, S_v;
382 ComponentVector K0 = K;
383 ComponentVector K1 = K;
386 if (verbosity == 3 || verbosity == 4) {
387 std::cout <<
"Stability test for vapor phase:" << std::endl;
389 checkStability_(fluid_state, isTrivialV, K0, y, S_v, z,
true, eos_type, verbosity);
390 bool V_unstable = (S_v < (1.0 + 1e-5)) || isTrivialV;
393 if (verbosity == 3 || verbosity == 4) {
394 std::cout <<
"Stability test for liquid phase:" << std::endl;
396 checkStability_(fluid_state, isTrivialL, K1, x, S_l, z,
false, eos_type, verbosity);
397 bool L_stable = (S_l < (1.0 + 1e-5)) || isTrivialL;
400 isStable = L_stable && V_unstable;
404 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
405 fluid_state.setMoleFraction(gasPhaseIdx, compIdx, z[compIdx]);
406 fluid_state.setMoleFraction(oilPhaseIdx, compIdx, z[compIdx]);
411 for (
int compIdx = 0; compIdx<numComponents; ++compIdx) {
412 K[compIdx] = y[compIdx] / x[compIdx];
419 template <
class FlashFlu
idState>
420 static typename FlashFluidState::Scalar wilsonK_(
const FlashFluidState& fluid_state,
int compIdx)
422 const auto& acf = FluidSystem::acentricFactor(compIdx);
423 const auto& T_crit = FluidSystem::criticalTemperature(compIdx);
424 const auto& T = fluid_state.temperature(0);
425 const auto& p_crit = FluidSystem::criticalPressure(compIdx);
426 const auto& p = fluid_state.pressure(0);
428 const auto& tmp = Opm::exp(5.3727 * (1+acf) * (1-T_crit/T)) * (p_crit/p);
432 template <
class Vector>
433 static typename Vector::field_type rachfordRice_g_(
const Vector& K,
typename Vector::field_type L,
const Vector& z)
435 typename Vector::field_type g=0;
436 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
437 g += (z[compIdx]*(K[compIdx]-1))/(K[compIdx]-L*(K[compIdx]-1));
442 template <
class Vector>
443 static typename Vector::field_type rachfordRice_dg_dL_(
const Vector& K,
const typename Vector::field_type L,
const Vector& z)
445 typename Vector::field_type dg=0;
446 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
447 dg += (z[compIdx]*(K[compIdx]-1)*(K[compIdx]-1))/((K[compIdx]-L*(K[compIdx]-1))*(K[compIdx]-L*(K[compIdx]-1)));
452 template <
class FlashFlu
idState,
class ComponentVector>
453 static void checkStability_(
const FlashFluidState& fluid_state,
bool& isTrivial, ComponentVector& K, ComponentVector& xy_loc,
454 typename FlashFluidState::Scalar& S_loc,
const ComponentVector& z,
bool isGas,
const EOSType& eos_type,
457 using FlashEval =
typename FlashFluidState::Scalar;
461 FlashFluidState fluid_state_fake = fluid_state;
462 FlashFluidState fluid_state_global = fluid_state;
465 if (verbosity >= 3) {
466 std::cout << std::setw(10) <<
"Iteration" << std::setw(16) <<
"K-Norm" << std::setw(16) <<
"R-Norm" << std::endl;
471 for (
int i = 0; i < 20000; ++i) {
474 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
475 xy_loc[compIdx] = K[compIdx] * z[compIdx];
476 S_loc += xy_loc[compIdx];
478 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
479 xy_loc[compIdx] /= S_loc;
480 fluid_state_fake.setMoleFraction(gasPhaseIdx, compIdx, xy_loc[compIdx]);
484 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
485 xy_loc[compIdx] = z[compIdx]/K[compIdx];
486 S_loc += xy_loc[compIdx];
488 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
489 xy_loc[compIdx] /= S_loc;
490 fluid_state_fake.setMoleFraction(oilPhaseIdx, compIdx, xy_loc[compIdx]);
494 int phaseIdx = (isGas ?
static_cast<int>(gasPhaseIdx) : static_cast<int>(oilPhaseIdx));
495 int phaseIdx2 = (isGas ?
static_cast<int>(oilPhaseIdx) : static_cast<int>(gasPhaseIdx));
497 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
498 fluid_state_global.setMoleFraction(phaseIdx2, compIdx, z[compIdx]);
501 typename FluidSystem::template ParameterCache<FlashEval> paramCache_fake(eos_type);
502 paramCache_fake.updatePhase(fluid_state_fake, phaseIdx);
504 typename FluidSystem::template ParameterCache<FlashEval> paramCache_global(eos_type);
505 paramCache_global.updatePhase(fluid_state_global, phaseIdx2);
508 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
509 auto phiFake = CubicEOS::computeFugacityCoefficient(fluid_state_fake, paramCache_fake, phaseIdx, compIdx);
510 auto phiGlobal = CubicEOS::computeFugacityCoefficient(fluid_state_global, paramCache_global, phaseIdx2, compIdx);
512 fluid_state_fake.setFugacityCoefficient(phaseIdx, compIdx, phiFake);
513 fluid_state_global.setFugacityCoefficient(phaseIdx2, compIdx, phiGlobal);
518 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
520 auto fug_fake = fluid_state_fake.fugacity(phaseIdx, compIdx);
521 auto fug_global = fluid_state_global.fugacity(phaseIdx2, compIdx);
522 auto fug_ratio = fug_global / fug_fake;
523 R[compIdx] = fug_ratio / S_loc;
526 auto fug_fake = fluid_state_fake.fugacity(phaseIdx, compIdx);
527 auto fug_global = fluid_state_global.fugacity(phaseIdx2, compIdx);
528 auto fug_ratio = fug_fake / fug_global;
529 R[compIdx] = fug_ratio * S_loc;
533 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
534 K[compIdx] *= R[compIdx];
538 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
539 auto a = Opm::getValue(R[compIdx]) - 1.0;
540 auto b = Opm::log(Opm::getValue(K[compIdx]));
546 if (verbosity >= 3) {
547 std::cout << std::setw(10) << i << std::setw(16) << K_norm << std::setw(16) << R_norm << std::endl;
551 isTrivial = (K_norm < 1e-5);
552 if (isTrivial || R_norm < 1e-10)
557 throw std::runtime_error(
" Stability test did not converge");
561 template <
class FlashFlu
idState,
class ComponentVector>
562 static void computeLiquidVapor_(FlashFluidState& fluid_state,
typename FlashFluidState::Scalar& L, ComponentVector& K,
const ComponentVector& z)
567 typename FlashFluidState::Scalar sumx=0;
568 typename FlashFluidState::Scalar sumy=0;
569 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
570 x[compIdx] = z[compIdx]/(L + (1-L)*K[compIdx]);
572 y[compIdx] = (K[compIdx]*z[compIdx])/(L + (1-L)*K[compIdx]);
578 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
579 fluid_state.setMoleFraction(oilPhaseIdx, compIdx, x[compIdx]);
580 fluid_state.setMoleFraction(gasPhaseIdx, compIdx, y[compIdx]);
584 template <
class Flu
idState,
class ComponentVector>
585 static void flash_2ph(
const ComponentVector& z_scalar,
586 const std::string& flash_2p_method,
587 ComponentVector& K_scalar,
588 typename FluidState::Scalar& L_scalar,
589 FluidState& fluid_state_scalar,
590 const Scalar flash_tolerance,
591 const EOSType& eos_type,
593 if (verbosity >= 1) {
594 std::cout <<
"Cell is two-phase! Solve Rachford-Rice with initial K = [" << K_scalar <<
"]" << std::endl;
599 bool converged =
false;
600 if (flash_2p_method ==
"newton") {
601 if (verbosity >= 1) {
602 std::cout <<
"Calculate composition using Newton." << std::endl;
604 converged = newtonComposition_(K_scalar, L_scalar, fluid_state_scalar, z_scalar, flash_tolerance, eos_type, verbosity);
605 }
else if (flash_2p_method ==
"ssi") {
607 if (verbosity >= 1) {
608 std::cout <<
"Calculate composition using Succcessive Substitution." << std::endl;
610 converged = successiveSubstitutionComposition_(K_scalar, L_scalar, fluid_state_scalar, z_scalar,
false, flash_tolerance, eos_type, verbosity);
611 }
else if (flash_2p_method ==
"ssi+newton") {
612 converged = successiveSubstitutionComposition_(K_scalar, L_scalar, fluid_state_scalar, z_scalar,
true, flash_tolerance, eos_type, verbosity);
614 converged = newtonComposition_(K_scalar, L_scalar, fluid_state_scalar, z_scalar, flash_tolerance, eos_type, verbosity);
617 throw std::logic_error(
"unknown two phase flash method " + flash_2p_method +
" is specified");
621 throw std::runtime_error(
"flash calculation did not get converged with " + flash_2p_method);
625 template <
class FlashFlu
idState,
class ComponentVector>
626 static bool newtonComposition_(ComponentVector& K,
627 typename FlashFluidState::Scalar& L,
628 FlashFluidState& fluid_state,
629 const ComponentVector& z,
630 const Scalar tolerance,
631 const EOSType& eos_type,
636 constexpr size_t num_equations = numMisciblePhases * numMiscibleComponents + 1;
637 constexpr size_t num_primary_variables = numMisciblePhases * numMiscibleComponents + 1;
638 using NewtonVector = Dune::FieldVector<Scalar, num_equations>;
639 using NewtonMatrix = Dune::FieldMatrix<Scalar, num_equations, num_primary_variables>;
641 NewtonVector soln(0.);
642 NewtonVector res(0.);
643 NewtonMatrix jac (0.);
646 computeLiquidVapor_(fluid_state, L, K, z);
649 if (verbosity >= 1) {
650 std::cout <<
" the current L is " << Opm::getValue(L) << std::endl;
651 std::cout <<
"Initial guess: x = [";
652 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
653 if (compIdx < numComponents - 1)
654 std::cout << fluid_state.moleFraction(oilPhaseIdx, compIdx) <<
" ";
656 std::cout << fluid_state.moleFraction(oilPhaseIdx, compIdx);
658 std::cout <<
"], y = [";
659 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
660 if (compIdx < numComponents - 1)
661 std::cout << fluid_state.moleFraction(gasPhaseIdx, compIdx) <<
" ";
663 std::cout << fluid_state.moleFraction(gasPhaseIdx, compIdx);
665 std::cout <<
"], and " <<
"L = " << L << std::endl;
669 if (verbosity == 2 || verbosity == 4) {
670 std::cout << std::setw(10) <<
"Iteration" << std::setw(16) <<
"Norm2(step)" << std::setw(16) <<
"Norm2(Residual)" << std::endl;
674 using Eval = DenseAd::Evaluation<Scalar, num_primary_variables>;
676 std::vector<Eval> x(numComponents), y(numComponents);
680 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
681 x[compIdx] = Eval(fluid_state.moleFraction(oilPhaseIdx, compIdx), compIdx);
682 const unsigned idx = compIdx + numComponents;
683 y[compIdx] = Eval(fluid_state.moleFraction(gasPhaseIdx, compIdx), idx);
685 l = Eval(L, num_primary_variables - 1);
688 CompositionalFluidState<Eval, FluidSystem> flash_fluid_state;
689 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
690 flash_fluid_state.setMoleFraction(FluidSystem::oilPhaseIdx, compIdx, x[compIdx]);
691 flash_fluid_state.setMoleFraction(FluidSystem::gasPhaseIdx, compIdx, y[compIdx]);
693 flash_fluid_state.setKvalue(compIdx, y[compIdx] / x[compIdx]);
695 flash_fluid_state.setLvalue(l);
697 flash_fluid_state.setPressure(FluidSystem::oilPhaseIdx,
698 fluid_state.pressure(FluidSystem::oilPhaseIdx));
699 flash_fluid_state.setPressure(FluidSystem::gasPhaseIdx,
700 fluid_state.pressure(FluidSystem::gasPhaseIdx));
703 flash_fluid_state.setSaturation(FluidSystem::gasPhaseIdx,
704 fluid_state.saturation(FluidSystem::gasPhaseIdx));
705 flash_fluid_state.setSaturation(FluidSystem::oilPhaseIdx,
706 fluid_state.saturation(FluidSystem::oilPhaseIdx));
707 flash_fluid_state.setTemperature(fluid_state.temperature(0));
709 using ParamCache =
typename FluidSystem::template ParameterCache<typename CompositionalFluidState<Eval, FluidSystem>::Scalar>;
710 ParamCache paramCache(eos_type);
712 for (
unsigned phaseIdx = 0; phaseIdx < numMisciblePhases; ++phaseIdx) {
713 paramCache.updatePhase(flash_fluid_state, phaseIdx);
714 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
716 Eval phi = FluidSystem::fugacityCoefficient(flash_fluid_state, paramCache, phaseIdx, compIdx);
717 flash_fluid_state.setFugacityCoefficient(phaseIdx, compIdx, phi);
720 bool converged =
false;
722 constexpr unsigned max_iter = 1000;
723 while (iter < max_iter) {
724 assembleNewton_<CompositionalFluidState<Eval, FluidSystem>, ComponentVector, num_primary_variables, num_equations>
725 (flash_fluid_state, z, jac, res);
726 if (verbosity >= 1) {
727 std::cout <<
" newton residual is " << res.two_norm() << std::endl;
729 converged = res.two_norm() < tolerance;
734 jac.solve(soln, res);
735 constexpr Scalar damping_factor = 1.0;
737 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
738 x[compIdx] -= soln[compIdx] * damping_factor;
739 y[compIdx] -= soln[compIdx + numComponents] * damping_factor;
741 l -= soln[num_equations - 1] * damping_factor;
743 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
744 flash_fluid_state.setMoleFraction(FluidSystem::oilPhaseIdx, compIdx, x[compIdx]);
745 flash_fluid_state.setMoleFraction(FluidSystem::gasPhaseIdx, compIdx, y[compIdx]);
747 flash_fluid_state.setKvalue(compIdx, y[compIdx] / x[compIdx]);
749 flash_fluid_state.setLvalue(l);
751 for (
unsigned phaseIdx = 0; phaseIdx < numMisciblePhases; ++phaseIdx) {
752 paramCache.updatePhase(flash_fluid_state, phaseIdx);
753 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
754 Eval phi = FluidSystem::fugacityCoefficient(flash_fluid_state, paramCache, phaseIdx, compIdx);
755 flash_fluid_state.setFugacityCoefficient(phaseIdx, compIdx, phi);
760 if (verbosity >= 1) {
761 for (
unsigned i = 0; i < num_equations; ++i) {
762 for (
unsigned j = 0; j < num_primary_variables; ++j) {
763 std::cout <<
" " << jac[i][j];
765 std::cout << std::endl;
767 std::cout << std::endl;
770 throw std::runtime_error(
" Newton composition update did not converge within maxIterations " + std::to_string(max_iter));
774 for (
unsigned idx = 0; idx < numComponents; ++idx) {
775 const auto x_i = Opm::getValue(flash_fluid_state.moleFraction(oilPhaseIdx, idx));
776 fluid_state.setMoleFraction(FluidSystem::oilPhaseIdx, idx, x_i);
777 const auto y_i = Opm::getValue(flash_fluid_state.moleFraction(gasPhaseIdx, idx));
778 fluid_state.setMoleFraction(FluidSystem::gasPhaseIdx, idx, y_i);
779 const auto K_i = Opm::getValue(flash_fluid_state.K(idx));
780 fluid_state.setKvalue(idx, K_i);
784 L = Opm::getValue(l);
785 fluid_state.setLvalue(L);
790 template<
typename FlashFlu
idState,
typename ComponentVector,
size_t num_primary,
size_t num_equation >
791 static void assembleNewton_(
const FlashFluidState& fluid_state,
792 const ComponentVector& global_composition,
793 Dune::FieldMatrix<double, num_equation, num_primary>& jac,
794 Dune::FieldVector<double, num_equation>& res)
796 using Eval = DenseAd::Evaluation<double, num_primary>;
797 std::vector<Eval> x(numComponents), y(numComponents);
798 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
799 x[compIdx] = fluid_state.moleFraction(oilPhaseIdx, compIdx);
800 y[compIdx] = fluid_state.moleFraction(gasPhaseIdx, compIdx);
802 const Eval& l = fluid_state.L();
806 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
809 auto local_res = -global_composition[compIdx] + l * x[compIdx] + (1 - l) * y[compIdx];
810 res[compIdx] = Opm::getValue(local_res);
811 for (
unsigned i = 0; i < num_primary; ++i) {
812 jac[compIdx][i] = local_res.derivative(i);
818 auto local_res = (fluid_state.fugacity(oilPhaseIdx, compIdx) -
819 fluid_state.fugacity(gasPhaseIdx, compIdx));
820 res[compIdx + numComponents] = Opm::getValue(local_res);
821 for (
unsigned i = 0; i < num_primary; ++i) {
822 jac[compIdx + numComponents][i] = local_res.derivative(i);
829 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
833 auto local_res = sumx - sumy;
834 res[num_equation - 1] = Opm::getValue(local_res);
835 for (
unsigned i = 0; i < num_primary; ++i) {
836 jac[num_equation - 1][i] = local_res.derivative(i);
840 template <
typename FlashFlu
idStateScalar,
typename Flu
idState>
841 static void updateDerivatives_(
const FlashFluidStateScalar& fluid_state_scalar,
842 FluidState& fluid_state,
843 const EOSType& eos_type,
844 bool is_single_phase)
847 updateDerivativesTwoPhase_(fluid_state_scalar, fluid_state, eos_type);
849 updateDerivativesSinglePhase_(fluid_state_scalar, fluid_state);
853 template <
typename FlashFlu
idStateScalar,
typename Flu
idState>
854 static void updateDerivativesTwoPhase_(
const FlashFluidStateScalar& fluid_state_scalar,
855 FluidState& fluid_state,
856 const EOSType& eos_type)
858 using InputEval =
typename FluidState::Scalar;
859 using ComponentVector = Dune::FieldVector<InputEval, numComponents>;
861 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
862 z[compIdx] = fluid_state.moleFraction(compIdx);
866 constexpr size_t num_equations = numMisciblePhases * numMiscibleComponents + 1;
867 constexpr size_t secondary_num_pv = numComponents + 1;
869 using SecondaryComponentVector = Dune::FieldVector<SecondaryEval, numComponents>;
872 SecondaryFlashFluidState secondary_fluid_state;
873 SecondaryComponentVector secondary_z;
876 const SecondaryEval sec_p = SecondaryEval(fluid_state_scalar.pressure(FluidSystem::oilPhaseIdx), 0);
877 secondary_fluid_state.setPressure(FluidSystem::oilPhaseIdx, sec_p);
878 secondary_fluid_state.setPressure(FluidSystem::gasPhaseIdx, sec_p);
881 secondary_fluid_state.setTemperature(Opm::getValue(fluid_state_scalar.temperature(0)));
883 for (
unsigned idx = 0; idx < numComponents; ++idx) {
884 secondary_z[idx] = SecondaryEval(Opm::getValue(z[idx]), idx + 1);
887 for (
unsigned idx = 0; idx < numComponents; ++idx) {
889 const auto x_i = fluid_state_scalar.moleFraction(oilPhaseIdx, idx);
890 secondary_fluid_state.setMoleFraction(FluidSystem::oilPhaseIdx, idx, x_i);
891 const auto y_i = fluid_state_scalar.moleFraction(gasPhaseIdx, idx);
892 secondary_fluid_state.setMoleFraction(FluidSystem::gasPhaseIdx, idx, y_i);
894 const auto K_i = fluid_state_scalar.K(idx);
895 secondary_fluid_state.setKvalue(idx, K_i);
897 const auto L = fluid_state_scalar.L();
898 secondary_fluid_state.setLvalue(L);
902 using SecondaryParamCache =
typename FluidSystem::template ParameterCache<typename SecondaryFlashFluidState::Scalar>;
903 SecondaryParamCache secondary_param_cache(eos_type);
904 for (
unsigned phase_idx = 0; phase_idx < numMisciblePhases; ++phase_idx) {
905 secondary_param_cache.updatePhase(secondary_fluid_state, phase_idx);
906 for (
unsigned comp_idx = 0; comp_idx < numComponents; ++comp_idx) {
907 SecondaryEval phi = FluidSystem::fugacityCoefficient(secondary_fluid_state, secondary_param_cache, phase_idx, comp_idx);
908 secondary_fluid_state.setFugacityCoefficient(phase_idx, comp_idx, phi);
912 using SecondaryNewtonVector = Dune::FieldVector<Scalar, num_equations>;
913 using SecondaryNewtonMatrix = Dune::FieldMatrix<Scalar, num_equations, secondary_num_pv>;
914 SecondaryNewtonMatrix sec_jac;
915 SecondaryNewtonVector sec_res;
919 assembleNewton_<SecondaryFlashFluidState, SecondaryComponentVector, secondary_num_pv, num_equations>
920 (secondary_fluid_state, secondary_z, sec_jac, sec_res);
925 constexpr size_t primary_num_pv = numMisciblePhases * numMiscibleComponents + 1;
927 using PrimaryComponentVector = Dune::FieldVector<double, numComponents>;
930 PrimaryFlashFluidState primary_fluid_state;
932 PrimaryComponentVector primary_z;
933 for (
unsigned comp_idx = 0; comp_idx < numComponents; ++comp_idx) {
934 primary_z[comp_idx] = Opm::getValue(z[comp_idx]);
936 for (
unsigned comp_idx = 0; comp_idx < numComponents; ++comp_idx) {
937 const auto x_ii = PrimaryEval(fluid_state_scalar.moleFraction(oilPhaseIdx, comp_idx), comp_idx);
938 primary_fluid_state.setMoleFraction(oilPhaseIdx, comp_idx, x_ii);
939 const unsigned idx = comp_idx + numComponents;
940 const auto y_ii = PrimaryEval(fluid_state_scalar.moleFraction(gasPhaseIdx, comp_idx), idx);
941 primary_fluid_state.setMoleFraction(gasPhaseIdx, comp_idx, y_ii);
942 primary_fluid_state.setKvalue(comp_idx, y_ii / x_ii);
945 l = PrimaryEval(fluid_state_scalar.L(), primary_num_pv - 1);
946 primary_fluid_state.setLvalue(l);
947 primary_fluid_state.setPressure(oilPhaseIdx, fluid_state_scalar.pressure(oilPhaseIdx));
948 primary_fluid_state.setPressure(gasPhaseIdx, fluid_state_scalar.pressure(gasPhaseIdx));
949 primary_fluid_state.setTemperature(fluid_state_scalar.temperature(0));
952 using PrimaryParamCache =
typename FluidSystem::template ParameterCache<typename PrimaryFlashFluidState::Scalar>;
953 PrimaryParamCache primary_param_cache(eos_type);
954 for (
unsigned phase_idx = 0; phase_idx < numMisciblePhases; ++phase_idx) {
955 primary_param_cache.updatePhase(primary_fluid_state, phase_idx);
956 for (
unsigned comp_idx = 0; comp_idx < numComponents; ++comp_idx) {
957 PrimaryEval phi = FluidSystem::fugacityCoefficient(primary_fluid_state, primary_param_cache, phase_idx, comp_idx);
958 primary_fluid_state.setFugacityCoefficient(phase_idx, comp_idx, phi);
962 using PrimaryNewtonVector = Dune::FieldVector<Scalar, num_equations>;
963 using PrimaryNewtonMatrix = Dune::FieldMatrix<Scalar, num_equations, primary_num_pv>;
964 PrimaryNewtonVector pri_res;
965 PrimaryNewtonMatrix pri_jac;
968 assembleNewton_<PrimaryFlashFluidState, PrimaryComponentVector, primary_num_pv, num_equations>
969 (primary_fluid_state, primary_z, pri_jac, pri_res);
975 sec_jac.template leftmultiply<PrimaryNewtonMatrix>(pri_jac);
977 ComponentVector x(numComponents), y(numComponents);
978 InputEval L_eval = L;
985 const auto p_l = fluid_state.pressure(FluidSystem::oilPhaseIdx);
986 const auto p_v = fluid_state.pressure(FluidSystem::gasPhaseIdx);
988 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
989 x[compIdx] = fluid_state_scalar.moleFraction(FluidSystem::oilPhaseIdx,compIdx);
990 y[compIdx] = fluid_state_scalar.moleFraction(FluidSystem::gasPhaseIdx,compIdx);
997 constexpr size_t num_deri = InputEval::numVars;
998 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
999 std::vector<double> deri(num_deri, 0.);
1001 for (
unsigned idx = 0; idx < num_deri; ++idx) {
1002 deri[idx] = -sec_jac[compIdx][0] * p_l.derivative(idx);
1005 for (
unsigned cIdx = 0; cIdx < numComponents; ++cIdx) {
1006 const double pz = -sec_jac[compIdx][cIdx + 1];
1007 const auto& zi = z[cIdx];
1008 for (
unsigned idx = 0; idx < num_deri; ++idx) {
1009 deri[idx] += pz * zi.derivative(idx);
1012 for (
unsigned idx = 0; idx < num_deri; ++idx) {
1013 x[compIdx].setDerivative(idx, deri[idx]);
1016 for (
unsigned idx = 0; idx < num_deri; ++idx) {
1017 deri[idx] = -sec_jac[compIdx + numComponents][0] * p_v.derivative(idx);
1019 for (
unsigned cIdx = 0; cIdx < numComponents; ++cIdx) {
1020 const double pz = -sec_jac[compIdx + numComponents][cIdx + 1];
1021 const auto& zi = z[cIdx];
1022 for (
unsigned idx = 0; idx < num_deri; ++idx) {
1023 deri[idx] += pz * zi.derivative(idx);
1026 for (
unsigned idx = 0; idx < num_deri; ++idx) {
1027 y[compIdx].setDerivative(idx, deri[idx]);
1031 std::vector<double> deriL(num_deri, 0.);
1032 for (
unsigned idx = 0; idx < num_deri; ++idx) {
1033 deriL[idx] = -sec_jac[2 * numComponents][0] * p_v.derivative(idx);
1035 for (
unsigned cIdx = 0; cIdx < numComponents; ++cIdx) {
1036 const double pz = -sec_jac[2 * numComponents][cIdx + 1];
1037 const auto& zi = z[cIdx];
1038 for (
unsigned idx = 0; idx < num_deri; ++idx) {
1039 deriL[idx] += pz * zi.derivative(idx);
1043 for (
unsigned idx = 0; idx < num_deri; ++idx) {
1044 L_eval.setDerivative(idx, deriL[idx]);
1049 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
1050 fluid_state.setMoleFraction(FluidSystem::oilPhaseIdx, compIdx, x[compIdx]);
1051 fluid_state.setMoleFraction(FluidSystem::gasPhaseIdx, compIdx, y[compIdx]);
1053 fluid_state.setLvalue(L_eval);
1056 template <
typename FlashFlu
idStateScalar,
typename Flu
idState>
1057 static void updateDerivativesSinglePhase_(
const FlashFluidStateScalar& fluid_state_scalar,
1058 FluidState& fluid_state)
1060 using InputEval =
typename FluidState::Scalar;
1062 InputEval L_eval = fluid_state_scalar.L();;
1066 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
1067 fluid_state.setMoleFraction(FluidSystem::oilPhaseIdx, compIdx, fluid_state.moleFraction(compIdx) );
1068 fluid_state.setMoleFraction(FluidSystem::gasPhaseIdx, compIdx, fluid_state.moleFraction(compIdx) );
1070 fluid_state.setLvalue(L_eval);
1075 template <
class FlashFlu
idState,
class ComponentVector>
1076 static bool successiveSubstitutionComposition_(ComponentVector& K,
1077 typename ComponentVector::field_type& L,
1078 FlashFluidState& fluid_state,
1079 const ComponentVector& z,
1080 const bool newton_afterwards,
1081 const Scalar flash_tolerance,
1082 const EOSType& eos_type,
1083 const int verbosity)
1086 const int maxIterations = newton_afterwards ? 5 : 100;
1089 std::ios_base::fmtflags f(std::cout.flags());
1093 std::cout <<
"Initial guess: K = [" << K <<
"] and L = " << L << std::endl;
1095 if (verbosity == 2 || verbosity == 4) {
1097 int fugWidth = (numComponents * 12)/2;
1098 int convWidth = fugWidth + 7;
1099 std::cout << std::setw(10) <<
"Iteration" << std::setw(fugWidth) <<
"fL/fV" << std::setw(convWidth) <<
"norm2(fL/fv-1)" << std::endl;
1104 for (
int i=0; i < maxIterations; ++i){
1106 computeLiquidVapor_(fluid_state, L, K, z);
1109 using ParamCache =
typename FluidSystem::template ParameterCache<typename FlashFluidState::Scalar>;
1110 ParamCache paramCache(eos_type);
1111 for (
int phaseIdx=0; phaseIdx<numMisciblePhases; ++phaseIdx){
1112 paramCache.updatePhase(fluid_state, phaseIdx);
1113 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
1114 auto phi = FluidSystem::fugacityCoefficient(fluid_state, paramCache, phaseIdx, compIdx);
1115 fluid_state.setFugacityCoefficient(phaseIdx, compIdx, phi);
1120 ComponentVector newFugRatio;
1121 ComponentVector convFugRatio;
1122 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
1123 newFugRatio[compIdx] = fluid_state.fugacity(oilPhaseIdx, compIdx)/fluid_state.fugacity(gasPhaseIdx, compIdx);
1124 convFugRatio[compIdx] = newFugRatio[compIdx] - 1.0;
1128 if (verbosity >= 2) {
1130 int fugWidth = (prec + 3);
1131 int convWidth = prec + 9;
1132 std::cout << std::defaultfloat;
1133 std::cout << std::fixed;
1134 std::cout << std::setw(5) << i;
1135 std::cout << std::setw(fugWidth);
1136 std::cout << std::setprecision(prec);
1137 std::cout << newFugRatio;
1138 std::cout << std::scientific;
1139 std::cout << std::setw(convWidth) << convFugRatio.two_norm() << std::endl;
1143 if (convFugRatio.two_norm() < flash_tolerance) {
1148 if (verbosity >= 1) {
1149 std::cout <<
"Solution converged to the following result :" << std::endl;
1150 std::cout <<
"x = [";
1151 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
1152 if (compIdx < numComponents - 1)
1153 std::cout << fluid_state.moleFraction(oilPhaseIdx, compIdx) <<
" ";
1155 std::cout << fluid_state.moleFraction(oilPhaseIdx, compIdx);
1157 std::cout <<
"]" << std::endl;
1158 std::cout <<
"y = [";
1159 for (
int compIdx = 0; compIdx < numComponents; ++compIdx) {
1160 if (compIdx < numComponents - 1)
1161 std::cout << fluid_state.moleFraction(gasPhaseIdx, compIdx) <<
" ";
1163 std::cout << fluid_state.moleFraction(gasPhaseIdx, compIdx);
1165 std::cout <<
"]" << std::endl;
1166 std::cout <<
"K = [" << K <<
"]" << std::endl;
1167 std::cout <<
"L = " << L << std::endl;
1175 for (
int compIdx=0; compIdx<numComponents; ++compIdx){
1176 K[compIdx] *= newFugRatio[compIdx];
1180 L = solveRachfordRice_g_(K, z, 0);
1185 const std::string msg = fmt::format(
"Successive substitution composition update did not converge within maxIterations {}.", maxIterations);
1186 if (!newton_afterwards) {
1187 throw std::runtime_error(msg);
1188 }
else if (verbosity > 0) {
1189 std::cout << msg << std::endl;