13 #include "MathUtils.h"
14 #include "Conversion.h"
23 params.addRequiredParam<UserObjectName>(
"brine_fp",
"The name of the user object for brine");
24 params.addRequiredParam<UserObjectName>(
"co2_fp",
"The name of the user object for CO2");
25 params.addParam<
unsigned int>(
"salt_component", 2,
"The component number of salt");
26 params.addClassDescription(
"Fluid state class for brine and CO2");
32 _salt_component(getParam<unsigned int>(
"salt_component")),
35 _Mh2o(_brine_fp.molarMassH2O()),
36 _invMh2o(1.0 / _Mh2o),
37 _Mco2(_co2_fp.molarMass()),
38 _Mnacl(_brine_fp.molarMassNaCl()),
43 _co2_henry(_co2_fp.henryCoefficients())
46 if (
_co2_fp.fluidName() !=
"co2")
47 paramError(
"co2_fp",
"A valid CO2 FluidProperties UserObject must be provided");
50 paramError(
"brine_fp",
"A valid Brine FluidProperties UserObject must be provided");
60 paramError(
"liquid_phase_number",
61 "This value is larger than the possible number of phases ",
66 paramError(
"liquid_fluid_component",
67 "This value is larger than the possible number of fluid components",
74 "The value provided must be different from the value entered in liquid_fluid_component");
78 paramError(
"salt_component",
79 "The value provided is larger than the possible number of fluid components",
97 std::vector<FluidStateProperties> & fsp)
const
107 Moose::derivInsert(p.derivatives(),
_pidx, 1.0);
109 Moose::derivInsert(T.derivatives(),
_Tidx, 1.0);
111 Moose::derivInsert(Zco2.derivatives(),
_Zidx, 1.0);
113 Moose::derivInsert(X.derivatives(),
_Xidx, 1.0);
163 const DualReal & Xnacl,
166 std::vector<FluidStateProperties> & fsp)
const
188 phaseState(Z.value(), Xco2.value(), Yco2.value(), phase_state);
204 Moose::derivInsert(Xco2.derivatives(),
_pidx, 0.0);
205 Moose::derivInsert(Xco2.derivatives(),
_Tidx, 0.0);
206 Moose::derivInsert(Xco2.derivatives(),
_Xidx, 0.0);
207 Moose::derivInsert(Xco2.derivatives(),
_Zidx, 1.0);
208 Moose::derivInsert(Yco2.derivatives(),
_pidx, 0.0);
209 Moose::derivInsert(Yco2.derivatives(),
_Tidx, 0.0);
210 Moose::derivInsert(Yco2.derivatives(),
_Xidx, 0.0);
219 Moose::derivInsert(Xco2.derivatives(),
_pidx, 0.0);
220 Moose::derivInsert(Xco2.derivatives(),
_Tidx, 0.0);
221 Moose::derivInsert(Xco2.derivatives(),
_Xidx, 0.0);
222 Moose::derivInsert(Yco2.derivatives(),
_pidx, 0.0);
223 Moose::derivInsert(Yco2.derivatives(),
_Tidx, 0.0);
224 Moose::derivInsert(Yco2.derivatives(),
_Xidx, 0.0);
225 Moose::derivInsert(Yco2.derivatives(),
_Zidx, 1.0);
248 std::vector<FluidStateProperties> & fsp)
const
254 DualReal co2_density, co2_viscosity;
264 mooseAssert(gas.
density.value() > 0.0,
"Gas density must be greater than zero");
271 const DualReal & Xnacl,
272 std::vector<FluidStateProperties> & fsp)
const
285 const DualReal liquid_density = 1.0 / (Xco2 / co2_partial_density + (1.0 - Xco2) / brine_density);
299 const DualReal liquid_enthalpy = (1.0 - Xco2) * brine_enthalpy + Xco2 * (co2_enthalpy + hdis);
302 liquid.
density = liquid_density;
306 mooseAssert(liquid.
density.value() > 0.0,
"Liquid density must be greater than zero");
313 const DualReal & Xnacl,
315 std::vector<FluidStateProperties> & fsp)
const
338 const DualReal liquid_density = 1.0 / (Xco2 / co2_partial_density + (1.0 - Xco2) / brine_density);
343 const DualReal K0 = Yco2 / Xco2;
344 const DualReal K1 = (1.0 - Yco2) / (1.0 - Xco2);
348 const DualReal
saturation = vapor_mass_fraction * liquid_density /
349 (gas_density + vapor_mass_fraction * (liquid_density - gas_density));
357 const DualReal & Xnacl,
360 std::vector<FluidStateProperties> & fsp)
const
378 const DualReal & Xnacl,
380 DualReal & Yh2o)
const
391 const DualReal mnacl = Xnacl / (1.0 - Xnacl) /
_Mnacl;
394 const DualReal mco2 = xco2 * (2.0 * mnacl +
_invMh2o) / (1.0 - xco2);
396 const DualReal denominator = (1.0 + mnacl *
_Mnacl + mco2 *
_Mco2);
397 Xco2 = mco2 *
_Mco2 / denominator;
403 const DualReal & co2_density,
405 DualReal & fh2o)
const
409 ": fugacityCoefficientsLowTemp() is not valid for T > 373.15K. Use "
410 "fugacityCoefficientsHighTemp() instead");
413 const DualReal pbar =
pressure * 1.0e-5;
416 const DualReal V =
_Mco2 / co2_density * 1.0e6;
419 const DualReal aCO2 = 7.54e7 - 4.13e4 *
temperature;
420 const Real bCO2 = 27.8;
421 const Real aCO2H2O = 7.89e7;
422 const Real bH2O = 18.18;
427 auto lnPhi = [V, aCO2, bCO2, t15,
this](DualReal a, DualReal b) {
428 return std::log(V / (V - bCO2)) + b / (V - bCO2) -
429 2.0 * a / (
_Rbar * t15 * bCO2) * std::log((V + bCO2) / V) +
430 aCO2 * b / (
_Rbar * t15 * bCO2 * bCO2) * (std::log((V + bCO2) / V) - bCO2 / (V + bCO2));
433 const DualReal lnPhiH2O = lnPhi(aCO2H2O, bH2O) - std::log(pbar * V / (
_Rbar *
temperature));
434 const DualReal lnPhiCO2 = lnPhi(aCO2, bCO2) - std::log(pbar * V / (
_Rbar *
temperature));
436 fh2o = std::exp(lnPhiH2O);
437 fco2 = std::exp(lnPhiCO2);
443 const DualReal & co2_density,
444 const DualReal & xco2,
445 const DualReal & yh2o,
447 DualReal & fh2o)
const
451 ": fugacityCoefficientsHighTemp() is not valid for T <= 373.15K. Use "
452 "fugacityCoefficientsLowTemp() instead");
461 const DualReal & co2_density,
462 const DualReal & xco2,
463 const DualReal & yh2o)
const
466 const DualReal pbar =
pressure * 1.0e-5;
468 const DualReal V =
_Mco2 / co2_density * 1.0e6;
471 const DualReal yco2 = 1.0 - yh2o;
472 const DualReal xh2o = 1.0 - xco2;
474 const DualReal aCO2 = 8.008e7 - 4.984e4 *
temperature;
475 const DualReal aH2O = 1.337e8 - 1.4e4 *
temperature;
476 const Real bCO2 = 28.25;
477 const Real bH2O = 15.7;
478 const DualReal KH2OCO2 = 1.427e-2 - 4.037e-4 *
temperature;
479 const DualReal KCO2H2O = 0.4228 - 7.422e-4 *
temperature;
480 const DualReal kH2OCO2 = KH2OCO2 * yh2o + KCO2H2O * yco2;
481 const DualReal kCO2H2O = KCO2H2O * yh2o + KH2OCO2 * yco2;
483 const DualReal aH2OCO2 = std::sqrt(aCO2 * aH2O) * (1.0 - kH2OCO2);
484 const DualReal aCO2H2O = std::sqrt(aCO2 * aH2O) * (1.0 - kCO2H2O);
486 const DualReal amix = yh2o * yh2o * aH2O + yh2o * yco2 * (aH2OCO2 + aCO2H2O) + yco2 * yco2 * aCO2;
487 const DualReal bmix = yh2o * bH2O + yco2 * bCO2;
491 DualReal lnPhiH2O = bH2O / bmix * (pbar * V / (
_Rbar *
temperature) - 1.0) -
493 DualReal term3 = (2.0 * yh2o * aH2O + yco2 * (aH2OCO2 + aCO2H2O) -
494 yh2o * yco2 * std::sqrt(aH2O * aCO2) * (kH2OCO2 - kCO2H2O) * (yh2o - yco2) +
495 xh2o * xco2 * std::sqrt(aH2O * aCO2) * (kH2OCO2 - kCO2H2O)) /
497 term3 -= bH2O / bmix;
498 term3 *= amix / (bmix *
_Rbar * t15) * std::log(V / (V + bmix));
501 return std::exp(lnPhiH2O);
507 const DualReal & co2_density,
508 const DualReal & xco2,
509 const DualReal & yh2o)
const
512 const DualReal pbar =
pressure * 1.0e-5;
514 const DualReal V =
_Mco2 / co2_density * 1.0e6;
517 const DualReal yco2 = 1.0 - yh2o;
518 const DualReal xh2o = 1.0 - xco2;
520 const DualReal aCO2 = 8.008e7 - 4.984e4 *
temperature;
521 const DualReal aH2O = 1.337e8 - 1.4e4 *
temperature;
522 const Real bCO2 = 28.25;
523 const Real bH2O = 15.7;
524 const DualReal KH2OCO2 = 1.427e-2 - 4.037e-4 *
temperature;
525 const DualReal KCO2H2O = 0.4228 - 7.422e-4 *
temperature;
526 const DualReal kH2OCO2 = KH2OCO2 * yh2o + KCO2H2O * yco2;
527 const DualReal kCO2H2O = KCO2H2O * yh2o + KH2OCO2 * yco2;
529 const DualReal aH2OCO2 = std::sqrt(aCO2 * aH2O) * (1.0 - kH2OCO2);
530 const DualReal aCO2H2O = std::sqrt(aCO2 * aH2O) * (1.0 - kCO2H2O);
532 const DualReal amix = yh2o * yh2o * aH2O + yh2o * yco2 * (aH2OCO2 + aCO2H2O) + yco2 * yco2 * aCO2;
533 const DualReal bmix = yh2o * bH2O + yco2 * bCO2;
537 DualReal lnPhiCO2 = bCO2 / bmix * (pbar * V / (
_Rbar *
temperature) - 1.0) -
540 DualReal term3 = (2.0 * yco2 * aCO2 + yh2o * (aH2OCO2 + aCO2H2O) -
541 yh2o * yco2 * std::sqrt(aH2O * aCO2) * (kH2OCO2 - kCO2H2O) * (yh2o - yco2) +
542 xh2o * xco2 * std::sqrt(aH2O * aCO2) * (kCO2H2O - kH2OCO2)) /
545 lnPhiCO2 += (term3 - bCO2 / bmix) * amix / (bmix *
_Rbar * t15) * std::log(V / (V + bmix));
547 return std::exp(lnPhiCO2);
552 const DualReal & xco2)
const
559 const DualReal xh2o = 1.0 - xco2;
560 const DualReal Am = -3.084e-2 * Tref + 1.927e-5 * Tref * Tref;
562 return std::exp((Am - 2.0 * Am * xh2o) * xco2 * xco2);
568 const DualReal & xco2)
const
575 const DualReal xh2o = 1.0 - xco2;
576 const DualReal Am = -3.084e-2 * Tref + 1.927e-5 * Tref * Tref;
578 return std::exp(2.0 * Am * xco2 * xh2o * xh2o);
585 const DualReal & Xnacl)
const
588 const DualReal pbar =
pressure * 1.0e-5;
590 const DualReal mnacl = Xnacl / (1.0 - Xnacl) /
_Mnacl;
597 const DualReal xi = 3.36389723e-4 - 1.9829898e-5 *
temperature +
601 return std::exp(2.0 * lambda * mnacl + xi * mnacl * mnacl);
606 const DualReal & Xnacl)
const
609 const DualReal mnacl = Xnacl / (1.0 - Xnacl) /
_Mnacl;
617 return (1.0 - mnacl /
_invMh2o) * std::exp(2.0 * lambda * mnacl + xi * mnacl * mnacl);
625 const DualReal Tc2 = Tc * Tc;
626 const DualReal Tc3 = Tc2 * Tc;
627 const DualReal Tc4 = Tc3 * Tc;
632 logK0H2O = -2.209 + 3.097e-2 * Tc - 1.098e-4 * Tc2 + 2.048e-7 * Tc3;
634 else if (Tc > 99.0 && Tc < 109.0)
636 const DualReal Tint = (Tc - 99.0) / 10.0;
637 const DualReal Tint2 = Tint * Tint;
638 logK0H2O = -0.0204026 + 0.0152513 * Tint + 0.417565 * Tint2 - 0.278636 * Tint * Tint2;
642 logK0H2O = -2.1077 + 2.8127e-2 * Tc - 8.4298e-5 * Tc2 + 1.4969e-7 * Tc3 - 1.1812e-10 * Tc4;
652 const DualReal Tc2 = Tc * Tc;
653 const DualReal Tc3 = Tc2 * Tc;
658 logK0CO2 = 1.189 + 1.304e-2 * Tc - 5.446e-5 * Tc2;
660 else if (Tc > 99.0 && Tc < 109.0)
662 const DualReal Tint = (Tc - 99.0) / 10.0;
663 const DualReal Tint2 = Tint * Tint;
664 logK0CO2 = 1.9462 + 2.25692e-2 * Tint - 9.49577e-3 * Tint2 - 6.77721e-3 * Tint * Tint2;
668 logK0CO2 = 1.668 + 3.992e-3 * Tc - 1.156e-5 * Tc2 + 1.593e-9 * Tc3;
676 const DualReal & Xnacl,
678 DualReal & yh2o)
const
691 Moose::derivInsert(Tlower.derivatives(),
_Tidx, 1.0);
693 DualReal xco2_lower, yh2o_lower;
696 const Real dxco2_dT_lower = xco2_lower.derivatives()[
_Tidx];
697 const Real dyh2o_dT_lower = yh2o_lower.derivatives()[
_Tidx];
700 Real xco2_upper, yh2o_upper;
704 pressure.value(),
_Tupper, Xnacl.value(), co2_density_upper, xco2_upper, yh2o_upper);
706 Real A, dA_dp, dA_dT, B, dB_dp, dB_dT, dB_dX;
721 const Real dyh2o_dT_upper =
722 ((1.0 - B) * dA_dT + (A - 1.0) * A * dB_dT) / (1.0 - A * B) / (1.0 - A * B);
723 const Real dxco2_dT_upper = dB_dT * (1.0 - yh2o_upper) - B * dyh2o_dT_upper;
726 Real xco2r, yh2or, dxco2_dT, dyh2o_dT;
728 Tint, xco2_lower.value(), dxco2_dT_lower, xco2_upper, dxco2_dT_upper, xco2r, dxco2_dT);
730 Tint, yh2o_lower.value(), dyh2o_dT_lower, yh2o_upper, dyh2o_dT_upper, yh2or, dyh2o_dT);
733 Moose::derivInsert(xco2.derivatives(),
_pidx, xco2_lower.derivatives()[
_pidx]);
734 Moose::derivInsert(xco2.derivatives(),
_Tidx, dxco2_dT);
735 Moose::derivInsert(xco2.derivatives(),
_Xidx, xco2_lower.derivatives()[
_Xidx]);
738 Moose::derivInsert(yh2o.derivatives(),
_pidx, yh2o_lower.derivatives()[
_pidx]);
739 Moose::derivInsert(yh2o.derivatives(),
_Tidx, dyh2o_dT);
740 Moose::derivInsert(yh2o.derivatives(),
_Xidx, yh2o_lower.derivatives()[
_Xidx]);
753 Real A, dA_dp, dA_dT, B, dB_dp, dB_dT, dB_dX;
768 const Real dyh2o_dp =
769 ((1.0 - B) * dA_dp + (A - 1.0) * A * dB_dp) / (1.0 - A * B) / (1.0 - A * B);
770 const Real dxco2_dp = dB_dp * (1.0 - yh2or) - B * dyh2o_dp;
772 const Real dyh2o_dT =
773 ((1.0 - B) * dA_dT + (A - 1.0) * A * dB_dT) / (1.0 - A * B) / (1.0 - A * B);
774 const Real dxco2_dT = dB_dT * (1.0 - yh2or) - B * dyh2o_dT;
776 const Real dyh2o_dX = ((A - 1.0) * A * dB_dX) / (1.0 - A * B) / (1.0 - A * B);
777 const Real dxco2_dX = dB_dX * (1.0 - yh2or) - B * dyh2o_dX;
780 Moose::derivInsert(xco2.derivatives(),
_pidx, dxco2_dp);
781 Moose::derivInsert(xco2.derivatives(),
_Tidx, dxco2_dT);
782 Moose::derivInsert(xco2.derivatives(),
_Xidx, dxco2_dX);
785 Moose::derivInsert(yh2o.derivatives(),
_pidx, dyh2o_dp);
786 Moose::derivInsert(yh2o.derivatives(),
_Tidx, dyh2o_dT);
787 Moose::derivInsert(yh2o.derivatives(),
_Xidx, dyh2o_dX);
794 const DualReal & Xnacl,
796 DualReal & yh2o)
const
800 ": equilibriumMoleFractionsLowTemp() is not valid for T > 373.15K. Use "
801 "equilibriumMoleFractions() instead");
814 const DualReal yh2ow = (1.0 - B) / (1.0 / A - B);
815 const DualReal xco2w = B * (1.0 - yh2ow);
818 const DualReal mco2w = xco2w *
_invMh2o / (1.0 - xco2w);
821 const DualReal mco2 = mco2w / gamma;
824 const DualReal mnacl = Xnacl / (1.0 - Xnacl) /
_Mnacl;
827 const DualReal total_moles = 2.0 * mnacl +
_invMh2o + mco2;
828 xco2 = mco2 / total_moles;
829 yh2o = A * (1.0 - xco2 - 2.0 * mnacl / total_moles);
835 const DualReal & co2_density,
841 ": funcABLowTemp() is not valid for T > 373.15K. Use funcABHighTemp() instead");
844 const DualReal pbar =
pressure * 1.0e-5;
847 const Real pref = 1.0;
848 const Real vCO2 = 32.6;
849 const Real vH2O = 18.1;
851 const DualReal delta_pbar = pbar - pref;
859 DualReal phiH2O, phiCO2;
862 A = K0H2O / (phiH2O * pbar) * std::exp(delta_pbar * vH2O / Rt);
863 B = phiCO2 * pbar / (
_invMh2o * K0CO2) * std::exp(-delta_pbar * vCO2 / Rt);
878 ": funcABHighTemp() is not valid for T <= 373.15K. Use funcABLowTemp() instead");
881 const Real pbar =
pressure * 1.0e-5;
886 const Real pref = -1.9906e-1 + 2.0471e-3 * Tc + 1.0152e-4 * Tc * Tc - 1.4234e-6 * Tc * Tc * Tc +
887 1.4168e-8 * Tc * Tc * Tc * Tc;
888 const Real vCO2 = 32.6 + 3.413e-2 * (Tc - 100.0);
889 const Real vH2O = 18.1 + 3.137e-2 * (Tc - 100.0);
891 const Real delta_pbar = pbar - pref;
903 const DualReal rhoco2 = co2_density;
904 const DualReal x = xco2;
905 const DualReal y = yh2o;
907 DualReal phiH2O, phiCO2;
916 const DualReal X = Xnacl;
919 A = K0H2O * gammaH2O / (phiH2O.value() * pbar) * std::exp(delta_pbar * vH2O / Rt);
920 B = phiCO2.value() * pbar / (
_invMh2o * K0CO2 * gamma * gammaCO2) *
921 std::exp(-delta_pbar * vCO2 / Rt);
942 const Real dp = 1.0e-2;
943 const Real dT = 1.0e-6;
944 const Real dX = 1.0e-8;
948 dA_dp = (A2 - A) / dp;
949 dB_dp = (B2 - B) / dp;
952 dA_dT = (A2 - A) / dT;
953 dB_dT = (B2 - B) / dT;
956 dB_dX = (B2 - B) / dX;
961 Real
pressure, Real
temperature, Real Xnacl, Real co2_density, Real & xco2, Real & yh2o)
const
968 const Real mnacl = Xnacl / (1.0 - Xnacl) /
_Mnacl;
979 auto fy = [mnacl,
this](Real y, Real A, Real B) {
981 (1.0 - B) *
_invMh2o / ((1.0 / A - B) * (2.0 * mnacl +
_invMh2o) + 2.0 * mnacl * B);
985 auto dfy = [mnacl,
this](Real A, Real B, Real dA, Real dB) {
986 const Real denominator = (1.0 / A - B) * (2.0 * mnacl +
_invMh2o) + 2.0 * mnacl * B;
987 return 1.0 +
_invMh2o * dB / denominator +
989 (2.0 * mnacl * dB - (2.0 * mnacl +
_invMh2o) * (dB + dA / A / A)) / denominator /
995 const Real dy = 1.0e-8;
998 unsigned int iter = 0;
999 const Real
tol = 1.0e-12;
1000 const unsigned int max_its = 10;
1003 while (std::abs(fy(y, A, B)) >
tol)
1011 y = y - fy(y, A, B) / dfy(A, B, dA, dB);
1031 const DualReal V = 37.51 - 9.585e-2 * Tc + 8.74e-4 * Tc * Tc - 5.044e-7 * Tc * Tc * Tc;
1033 return 1.0e6 *
_Mco2 / V;
1046 const DualReal X = Xnacl;
1054 DualReal Xco2, Yh2o;
1058 const DualReal Yco2 = 1.0 - Yh2o;
1068 const DualReal liquid_saturation = 1.0 -
saturation;
1086 const std::vector<Real> b{1.19784e-1, -7.17823e-4, 4.93854e-6, -1.03826e-8, 1.08233e-11};
1092 for (
unsigned int i = 0; i < b.size(); ++i)
1096 const DualReal xmol = Xnacl / (1.0 - Xnacl) /
_Mnacl;
1099 return Kh_h2o *
std::pow(10.0, xmol * kb);
1104 const DualReal & Xnacl)
const
1117 const Real dhdis_dT =
1118 (-
_R * T2 * T2 * Kh2.derivatives()[
_Tidx] / Kh2 /
_Mco2 - hdis).value() / dT;
1120 const Real dX = Xnacl.value() * 1.0e-8;
1121 const DualReal X2 = Xnacl + dX;
1124 const Real dhdis_dX =
1128 hdis.derivatives() =
temperature.derivatives() * dhdis_dT + Xnacl.derivatives() * dhdis_dX;
1137 const DualReal delta_h = -58.3533 + 0.134519 *
temperature;
1140 return delta_h * 1000.0 /
_Mco2;
1145 Real
temperature, Real f0, Real df0, Real f1, Real df1, Real & value, Real & deriv)
const
1151 const Real b = df0 * dT;
1152 const Real c = 3.0 * (f1 - f0) - (2.0 * df0 + df1) * dT;
1153 const Real d = 2.0 * (f0 - f1) + (df0 + df1) * dT;
1167 if (temperature < 285.15 || temperature > 573.15)
1168 mooseException(
name() +
": temperature " + Moose::stringify(
temperature) +
1169 " is outside range 285.15 K <= T <= 573.15 K");
1172 mooseException(
name() +
": pressure " + Moose::stringify(
pressure) +
1173 " must be less than 60 MPa");