49 typedef ::Opm::IdealGas<Scalar> IdealGas;
50 static const int liquidPhaseIdx = 0;
51 static const int gasPhaseIdx = 1;
64 template <
class Evaluation,
class CO2Params>
65 OPM_HOST_DEVICE
static Evaluation
gasDiffCoeff(
const CO2Params& params,
66 const Evaluation& temperature,
67 const Evaluation& pressure,
68 bool extrapolate =
false)
71 Scalar k = 1.3806504e-23;
73 Scalar R_h = 1.72e-10;
74 const Evaluation& mu =
CO2::gasViscosity(params, temperature, pressure, extrapolate);
75 return k / (c * std::numbers::pi * R_h) * (temperature / mu);
84 template <
class Evaluation>
85 OPM_HOST_DEVICE
static Evaluation
liquidDiffCoeff(
const Evaluation& ,
const Evaluation& )
111 template <
class Evaluation,
class CO2Params>
112 OPM_HOST_DEVICE
static void
114 const Evaluation& temperature,
115 const Evaluation& pg,
117 const int knownPhaseIdx,
120 const int& activityModel,
121 bool extrapolate =
false)
123 OPM_TIMEFUNCTION_LOCAL(Subsystem::PvtProps);
126 bool iterate =
false;
127 if ((activityModel == 1 &&
salinity > 0.0) || (activityModel == 2 && temperature > 372.15)) {
133 if (knownPhaseIdx < 0) {
134 Evaluation molalityNaCl = massFracToMolality_(
salinity);
140 if (activityModel == 3) {
141 auto [xCO2, yH2O] = mutualSolubilitySpycherPruess2005_(params, temperature, pg, molalityNaCl, extrapolate);
149 auto [xCO2, yH2O] = fixPointIterSolubility_(params, temperature, pg, molalityNaCl, activityModel, extrapolate);
156 auto [xCO2, yH2O] = nonIterSolubility_(params, temperature, pg, molalityNaCl, activityModel, extrapolate);
166 else if (knownPhaseIdx == liquidPhaseIdx && activityModel == 3) {
167 Evaluation x_NaCl = salinityToMolFrac_(
salinity);
168 const Evaluation& A = computeA_(params, temperature, pg, Evaluation(0.0), Evaluation(0.0),
false, extrapolate,
true);
169 ygH2O = A * (1 - xlCO2 - x_NaCl);
175 else if (knownPhaseIdx == gasPhaseIdx && activityModel == 3) {
177 Evaluation x_NaCl = salinityToMolFrac_(
salinity);
178 const Evaluation& A = computeA_(params, temperature, pg, Evaluation(0.0), Evaluation(0.0),
false, extrapolate,
true);
179 xlCO2 = 1 - x_NaCl - ygH2O / A;
186 template <
class Evaluation>
187 OPM_HOST_DEVICE
static Evaluation
henry(
const Evaluation& temperature,
bool extrapolate =
false)
203 template <
class Evaluation,
class CO2Params>
205 const Evaluation& temperature,
206 const Evaluation& pg,
207 const Evaluation& yH2O,
209 bool extrapolate =
false,
210 bool spycherPruess2005 =
false)
212 OPM_TIMEFUNCTION_LOCAL(Subsystem::PvtProps);
217 const Evaluation pg_bar = pg / 1.e5;
221 const Evaluation a_CO2 = aCO2_(temperature, highTemp);
222 const Evaluation a_CO2_H2O = aCO2_H2O_(temperature, yH2O, highTemp);
223 const Evaluation a_mix = aMix_(temperature, yH2O, highTemp);
224 const Scalar b_CO2 = bCO2_(highTemp);
225 const Evaluation b_mix = bMix_(yH2O, highTemp);
226 const Evaluation Rt15 = R * pow(temperature, 1.5);
229 if (spycherPruess2005) {
230 const Evaluation logVpb_V = log((V + b_CO2) / V);
231 lnPhiCO2 = log(V / (V - b_CO2));
232 lnPhiCO2 += b_CO2 / (V - b_CO2);
233 lnPhiCO2 -= 2 * a_CO2 / (Rt15 * b_CO2) * logVpb_V;
240 - b_CO2 / (V + b_CO2));
241 lnPhiCO2 -= log(pg_bar * V / (R * temperature));
244 lnPhiCO2 = (b_CO2 / b_mix) * (pg_bar * V / (R * temperature) - 1);
245 lnPhiCO2 -= log(pg_bar * (V - b_mix) / (R * temperature));
246 lnPhiCO2 += (2 * (yH2O * a_CO2_H2O + (1 - yH2O) * a_CO2) / a_mix - (b_CO2 / b_mix)) *
247 a_mix / (b_mix * Rt15) * log(V / (V + b_mix));
249 return exp(lnPhiCO2);
265 template <
class Evaluation,
class CO2Params>
267 const Evaluation& temperature,
268 const Evaluation& pg,
269 const Evaluation& yH2O,
271 bool extrapolate =
false,
272 bool spycherPruess2005 =
false)
274 OPM_TIMEFUNCTION_LOCAL(Subsystem::PvtProps);
279 const Evaluation& pg_bar = pg / 1.e5;
283 const Evaluation a_H2O = aH2O_(temperature, highTemp);
284 const Evaluation a_CO2_H2O = aCO2_H2O_(temperature, yH2O, highTemp);
285 const Evaluation a_mix = aMix_(temperature, yH2O, highTemp);
286 const Scalar b_H2O = bH2O_(highTemp);
287 const Evaluation b_mix = bMix_(yH2O, highTemp);
288 const Evaluation Rt15 = R * pow(temperature, 1.5);
291 if (spycherPruess2005) {
292 const Evaluation logVpb_V = log((V + b_mix) / V);
295 + b_H2O/(V - b_mix) - 2*a_CO2_H2O
296 / (Rt15*b_mix)*logVpb_V
297 + a_mix*b_H2O/(Rt15*b_mix*b_mix)
298 *(logVpb_V - b_mix/(V + b_mix))
299 - log(pg_bar*V/(R*temperature));
302 lnPhiH2O = (b_H2O / b_mix) * (pg_bar * V / (R * temperature) - 1);
303 lnPhiH2O -= log(pg_bar * (V - b_mix) / (R * temperature));
304 lnPhiH2O += (2 * (yH2O * a_H2O + (1 - yH2O) * a_CO2_H2O) / a_mix - (b_H2O / b_mix)) *
305 a_mix / (b_mix * Rt15) * log(V / (V + b_mix));
307 return exp(lnPhiH2O);
314 template <
class Evaluation>
315 OPM_HOST_DEVICE
static Evaluation aCO2_(
const Evaluation& temperature,
const bool& highTemp)
318 return 8.008e7 - 4.984e4 * temperature;
321 return 7.54e7 - 4.13e4 * temperature;
328 template <
class Evaluation>
329 OPM_HOST_DEVICE
static Evaluation aH2O_(
const Evaluation& temperature,
const bool& highTemp)
332 return 1.337e8 - 1.4e4 * temperature;
342 template <
class Evaluation>
343 OPM_HOST_DEVICE
static Evaluation aCO2_H2O_(
const Evaluation& temperature,
const Evaluation& yH2O,
const bool& highTemp)
347 Evaluation aCO2 = aCO2_(temperature, highTemp);
348 Evaluation aH2O = aH2O_(temperature, highTemp);
351 Evaluation K_CO2_H2O = 0.4228 - 7.422e-4 * temperature;
352 Evaluation K_H2O_CO2 = 1.427e-2 - 4.037e-4 * temperature;
353 Evaluation k_CO2_H2O = yH2O * K_H2O_CO2 + (1 - yH2O) * K_CO2_H2O;
356 return sqrt(aCO2 * aH2O) * (1 - k_CO2_H2O);
366 template <
class Evaluation>
367 OPM_HOST_DEVICE
static Evaluation aMix_(
const Evaluation& temperature,
const Evaluation& yH2O,
const bool& highTemp)
371 Evaluation aCO2 = aCO2_(temperature, highTemp);
372 Evaluation aH2O = aH2O_(temperature, highTemp);
373 Evaluation a_CO2_H2O = aCO2_H2O_(temperature, yH2O, highTemp);
375 return yH2O * yH2O * aH2O + 2 * yH2O * (1 - yH2O) * a_CO2_H2O + (1 - yH2O) * (1 - yH2O) * aCO2;
378 return aCO2_(temperature, highTemp);
385 OPM_HOST_DEVICE
static Scalar bCO2_(
const bool& highTemp)
398 OPM_HOST_DEVICE
static Scalar bH2O_(
const bool& highTemp)
411 template <
class Evaluation>
412 OPM_HOST_DEVICE
static Evaluation bMix_(
const Evaluation& yH2O,
const bool& highTemp)
416 Scalar bCO2 = bCO2_(highTemp);
417 Scalar bH2O = bH2O_(highTemp);
419 return yH2O * bH2O + (1 - yH2O) * bCO2;
422 return bCO2_(highTemp);
429 template <
class Evaluation>
430 OPM_HOST_DEVICE
static Evaluation V_avg_CO2_(
const Evaluation& temperature,
const bool& highTemp)
432 if (highTemp && (temperature > 373.15)) {
433 return 32.6 + 3.413e-2 * (temperature - 373.15);
443 template <
class Evaluation>
444 OPM_HOST_DEVICE
static Evaluation V_avg_H2O_(
const Evaluation& temperature,
const bool& highTemp)
446 if (highTemp && (temperature > 373.15)) {
447 return 18.1 + 3.137e-2 * (temperature - 373.15);
457 template <
class Evaluation>
458 OPM_HOST_DEVICE
static Evaluation AM_(
const Evaluation& temperature,
const bool& highTemp)
460 if (highTemp && temperature > 373.15) {
461 Evaluation deltaTk = temperature - 373.15;
462 return deltaTk * (-3.084e-2 + 1.927e-5 * deltaTk);
472 template <
class Evaluation>
473 OPM_HOST_DEVICE
static Evaluation Pref_(
const Evaluation& temperature,
const bool& highTemp)
475 if (highTemp && temperature > 373.15) {
476 const Evaluation& temperatureCelcius = temperature - 273.15;
477 static const Scalar c[5] = { -1.9906e-1, 2.0471e-3, 1.0152e-4, -1.4234e-6, 1.4168e-8 };
478 return c[0] + temperatureCelcius * (c[1] + temperatureCelcius * (c[2] +
479 temperatureCelcius * (c[3] + temperatureCelcius * c[4])));
489 template <
class Evaluation>
490 OPM_HOST_DEVICE
static Evaluation activityCoefficientCO2_(
const Evaluation& temperature,
491 const Evaluation& xCO2,
492 const bool& highTemp)
496 Evaluation AM = AM_(temperature, highTemp);
497 Evaluation lnGammaCO2 = 2 * AM * xCO2 * (1 - xCO2) * (1 - xCO2);
498 return exp(lnGammaCO2);
508 template <
class Evaluation>
509 OPM_HOST_DEVICE
static Evaluation activityCoefficientH2O_(
const Evaluation& temperature,
510 const Evaluation& xCO2,
511 const bool& highTemp)
515 Evaluation AM = AM_(temperature, highTemp);
516 Evaluation lnGammaH2O = (1 - 2 * (1 - xCO2)) * AM * xCO2 * xCO2;
517 return exp(lnGammaH2O);
529 template <
class Evaluation>
530 OPM_HOST_DEVICE
static Evaluation salinityToMolFrac_(
const Evaluation&
salinity) {
531 OPM_TIMEFUNCTION_LOCAL(Subsystem::PvtProps);
533 const Scalar Ms = 58.44e-3;
537 const Evaluation x_NaCl = -Mw * X_NaCl / ((Ms - Mw) * X_NaCl - Ms);
547 template <
class Evaluation>
548 OPM_HOST_DEVICE
static Evaluation moleFracToMolality_(
const Evaluation& x_NaCl)
551 return 55.508 * x_NaCl / (1 - x_NaCl);
555 template <
class Evaluation>
556 OPM_HOST_DEVICE
static Evaluation massFracToMolality_(
const Evaluation& X_NaCl)
558 const Scalar MmNaCl = 58.44e-3;
559 return X_NaCl / (MmNaCl * (1 - X_NaCl));
567 template <
class Evaluation>
568 OPM_HOST_DEVICE
static Evaluation molalityToMoleFrac_(
const Evaluation& m_NaCl)
571 return m_NaCl / (55.508 + m_NaCl);
577 template <
class Evaluation,
class CO2Parameters>
578 OPM_HOST_DEVICE
static std::pair<Evaluation, Evaluation> fixPointIterSolubility_(
const CO2Parameters& params,
579 const Evaluation& temperature,
580 const Evaluation& pg,
581 const Evaluation& m_NaCl,
582 const int& activityModel,
583 bool extrapolate =
false)
585 OPM_TIMEFUNCTION_LOCAL(Subsystem::PvtProps);
588 Evaluation xCO2 = 0.009;
589 Evaluation gammaNaCl = 1.0;
592 if (m_NaCl > 0.0 && activityModel == 2) {
593 gammaNaCl = activityCoefficientSalt_(temperature, pg, m_NaCl, Evaluation(0.0), activityModel);
599 bool highTemp =
true;
600 if (activityModel == 1) {
603 const bool iterate =
true;
606 for (
int i = 0; i < max_iter; ++i) {
608 if (m_NaCl > 0.0 && activityModel == 1) {
609 gammaNaCl = activityCoefficientSalt_(temperature, pg, m_NaCl, xCO2, activityModel);
613 auto [xCO2_new, yH2O_new] = mutualSolubility_(params, temperature, pg, xCO2, yH2O, m_NaCl, gammaNaCl, highTemp,
614 iterate, extrapolate);
617 if (abs(xCO2_new - xCO2) < tol && abs(yH2O_new - yH2O) < tol) {
636 template <
class Evaluation,
class CO2Parameters>
637 OPM_HOST_DEVICE
static std::pair<Evaluation, Evaluation> nonIterSolubility_(
const CO2Parameters& params,
638 const Evaluation& temperature,
639 const Evaluation& pg,
640 const Evaluation& m_NaCl,
641 const int& activityModel,
642 bool extrapolate =
false)
645 Evaluation gammaNaCl = 1.0;
646 if (m_NaCl > 0.0 && activityModel > 0 && activityModel < 3) {
647 gammaNaCl = activityCoefficientSalt_(temperature, pg, m_NaCl, Evaluation(0.0), activityModel);
652 const bool highTemp =
false;
653 const bool iterate =
false;
654 auto [xCO2, yH2O] = mutualSolubility_(params, temperature, pg, Evaluation(0.0), Evaluation(0.0), m_NaCl, gammaNaCl,
655 highTemp, iterate, extrapolate);
663 template <
class Evaluation,
class CO2Parameters>
664 OPM_HOST_DEVICE
static std::pair<Evaluation, Evaluation> mutualSolubility_(
const CO2Parameters& params,
665 const Evaluation& temperature,
666 const Evaluation& pg,
667 const Evaluation& xCO2,
668 const Evaluation& yH2O,
669 const Evaluation& m_NaCl,
670 const Evaluation& gammaNaCl,
671 const bool& highTemp,
673 bool extrapolate =
false)
676 const Evaluation& A = computeA_(params, temperature, pg, yH2O, xCO2, highTemp, extrapolate);
677 Evaluation B = computeB_(params, temperature, pg, yH2O, xCO2, highTemp, extrapolate);
683 Evaluation yH2O_new = (1. - B) * 55.508 / ((1. / A - B) * (2 * m_NaCl + 55.508) + 2 * m_NaCl * B);
686 xCO2_new = B * (1 - yH2O);
689 xCO2_new = B * (1 - yH2O_new);
692 return {xCO2_new, yH2O_new};
698 template <
class Evaluation,
class CO2Parameters>
699 OPM_HOST_DEVICE
static std::pair<Evaluation, Evaluation> mutualSolubilitySpycherPruess2005_(
const CO2Parameters& params,
700 const Evaluation& temperature,
701 const Evaluation& pg,
702 const Evaluation& m_NaCl,
703 bool extrapolate =
false)
706 const Evaluation& A = computeA_(params, temperature, pg, Evaluation(0.0), Evaluation(0.0),
false, extrapolate,
true);
707 const Evaluation& B = computeB_(params, temperature, pg, Evaluation(0.0), Evaluation(0.0),
false, extrapolate,
true);
710 Evaluation yH2O = (1 - B) / (1. / A - B);
711 Evaluation xCO2 = B * (1 - yH2O);
715 const Evaluation& gammaNaCl = activityCoefficientSalt_(temperature, pg, m_NaCl, Evaluation(0.0), 3);
718 Evaluation mCO2 = (xCO2 * 55.508) / (1 - xCO2);
720 xCO2 = mCO2 / (m_NaCl + 55.508 + mCO2);
723 const Evaluation& xNaCl = molalityToMoleFrac_(m_NaCl);
724 yH2O = A * (1 - xCO2 - xNaCl);
738 template <
class Evaluation,
class CO2Params>
739 OPM_HOST_DEVICE
static Evaluation computeA_(
const CO2Params& params,
740 const Evaluation& temperature,
741 const Evaluation& pg,
742 const Evaluation& yH2O,
743 const Evaluation& xCO2,
744 const bool& highTemp,
745 bool extrapolate =
false,
746 bool spycherPruess2005 =
false)
748 OPM_TIMEFUNCTION_LOCAL(Subsystem::PvtProps);
750 const Evaluation& deltaP = pg / 1e5 - Pref_(temperature, highTemp);
751 Evaluation v_av_H2O = V_avg_H2O_(temperature, highTemp);
752 Evaluation k0_H2O = equilibriumConstantH2O_(temperature, highTemp);
753 Evaluation phi_H2O =
fugacityCoefficientH2O(params, temperature, pg, yH2O, highTemp, extrapolate, spycherPruess2005);
754 Evaluation gammaH2O = activityCoefficientH2O_(temperature, xCO2, highTemp);
758 if ( temperature > 372.15 && temperature < 382.15 && !spycherPruess2005) {
759 const Evaluation weight = (382.15 - temperature) / 10.;
760 const Evaluation& k0_H2O_low = equilibriumConstantH2O_(temperature,
false);
761 const Evaluation& phi_H2O_low =
fugacityCoefficientH2O(params, temperature, pg, Evaluation(0.0),
false, extrapolate);
762 k0_H2O = k0_H2O * (1 - weight) + k0_H2O_low * weight;
763 phi_H2O = phi_H2O * (1 - weight) + phi_H2O_low * weight;
767 const Evaluation& pg_bar = pg / 1.e5;
769 return k0_H2O * gammaH2O / (phi_H2O * pg_bar) * exp(deltaP * v_av_H2O / (R * temperature));
780 template <
class Evaluation,
class CO2Parameters>
781 OPM_HOST_DEVICE
static Evaluation computeB_(
const CO2Parameters& params,
782 const Evaluation& temperature,
783 const Evaluation& pg,
784 const Evaluation& yH2O,
785 const Evaluation& xCO2,
786 const bool& highTemp,
787 bool extrapolate =
false,
788 bool spycherPruess2005 =
false)
790 OPM_TIMEFUNCTION_LOCAL(Subsystem::PvtProps);
792 const Evaluation& deltaP = pg / 1e5 - Pref_(temperature, highTemp);
793 Evaluation v_av_CO2 = V_avg_CO2_(temperature, highTemp);
794 Evaluation k0_CO2 = equilibriumConstantCO2_(temperature, pg, highTemp, spycherPruess2005);
795 Evaluation phi_CO2 =
fugacityCoefficientCO2(params, temperature, pg, yH2O, highTemp, extrapolate, spycherPruess2005);
796 Evaluation gammaCO2 = activityCoefficientCO2_(temperature, xCO2, highTemp);
800 if ( temperature > 372.15 && temperature < 382.15 && !spycherPruess2005) {
801 const Evaluation weight = (382.15 - temperature) / 10.;
802 const Evaluation& k0_CO2_low = equilibriumConstantCO2_(temperature, pg,
false, spycherPruess2005);
803 const Evaluation& phi_CO2_low =
fugacityCoefficientCO2(params, temperature, pg, Evaluation(0.0),
false, extrapolate);
804 k0_CO2 = k0_CO2 * (1 - weight) + k0_CO2_low * weight;
805 phi_CO2 = phi_CO2 * (1 - weight) + phi_CO2_low * weight;
809 const Evaluation& pg_bar = pg / 1.e5;
811 return phi_CO2 * pg_bar / (55.508 * k0_CO2 * gammaCO2) * exp(-deltaP * v_av_CO2 / (R * temperature));
817 template <
class Evaluation>
818 OPM_HOST_DEVICE
static Evaluation activityCoefficientSalt_(
const Evaluation& temperature,
819 const Evaluation& pg,
820 const Evaluation& m_NaCl,
821 const Evaluation& xCO2,
822 const int& activityModel)
824 OPM_TIMEFUNCTION_LOCAL(Subsystem::PvtProps);
831 if (activityModel == 1) {
832 lambda = computeLambdaRumpfetal_(temperature);
834 Evaluation m_CO2 = xCO2 * (2 * m_NaCl + 55.508) / (1 - xCO2);
835 convTerm = (1 + (m_CO2 + 2 * m_NaCl) / 55.508) / (1 + m_CO2 / 55.508);
837 else if (activityModel == 2) {
838 lambda = computeLambdaSpycherPruess2009_(temperature);
839 xi = computeXiSpycherPruess2009_(temperature);
840 convTerm = 1 + 2 * m_NaCl / 55.508;
842 else if (activityModel == 3) {
843 lambda = computeLambdaDuanSun_(temperature, pg);
844 xi = computeXiDuanSun_(temperature, pg);
848 OPM_THROW(std::invalid_argument,
"Activity model for salt-out effect has not been implemented!.");
852 const Evaluation& lnGamma = 2 * lambda * m_NaCl + xi * m_NaCl * m_NaCl;
855 return convTerm * exp(lnGamma);
861 template <
class Evaluation>
862 OPM_HOST_DEVICE
static Evaluation computeLambdaSpycherPruess2009_(
const Evaluation& temperature)
865 static const Scalar c[3] = { 2.217e-4, 1.074, 2648. };
868 return c[0] * temperature + c[1] / temperature + c[2] / (temperature * temperature);
874 template <
class Evaluation>
875 OPM_HOST_DEVICE
static Evaluation computeXiSpycherPruess2009_(
const Evaluation& temperature)
878 static const Scalar c[3] = { 1.3e-5, -20.12, 5259. };
881 return c[0] * temperature + c[1] / temperature + c[2] / (temperature * temperature);
887 template <
class Evaluation>
888 OPM_HOST_DEVICE
static Evaluation computeLambdaRumpfetal_(
const Evaluation& temperature)
891 static const Scalar c[4] = { 0.254, -76.82, -10656, 6312e3 };
893 return c[0] + c[1] / temperature + c[2] / (temperature * temperature) +
894 c[3] / (temperature * temperature * temperature);
904 template <
class Evaluation>
905 OPM_HOST_DEVICE
static Evaluation computeLambdaDuanSun_(
const Evaluation& temperature,
const Evaluation& pg)
907 static const Scalar c[6] =
908 { -0.411370585, 6.07632013E-4, 97.5347708, -0.0237622469, 0.0170656236, 1.41335834E-5 };
910 Evaluation pg_bar = pg / 1.0E5;
911 return c[0] + c[1]*temperature + c[2]/temperature + c[3]*pg_bar/temperature + c[4]*pg_bar/(630.0 - temperature)
912 + c[5]*temperature*log(pg_bar);
922 template <
class Evaluation>
923 OPM_HOST_DEVICE
static Evaluation computeXiDuanSun_(
const Evaluation& temperature,
const Evaluation& pg)
925 static const Scalar c[4] =
926 { 3.36389723E-4, -1.98298980E-5, 2.12220830E-3, -5.24873303E-3 };
928 Evaluation pg_bar = pg / 1.0E5;
929 return c[0] + c[1]*temperature + c[2]*pg_bar/temperature + c[3]*pg_bar/(630.0 - temperature);
938 template <
class Evaluation>
939 OPM_HOST_DEVICE
static Evaluation equilibriumConstantCO2_(
const Evaluation& temperature,
940 const Evaluation& pg,
941 const bool& highTemp,
942 bool spycherPruess2005 =
false)
944 OPM_TIMEFUNCTION_LOCAL(Subsystem::PvtProps);
945 Evaluation temperatureCelcius = temperature - 273.15;
946 std::array<Scalar, 4> c;
948 c = { 1.668, 3.992e-3, -1.156e-5, 1.593e-9 };
959 c = { 1.169, 1.368e-2, -5.38e-5, 0.0 };
962 c = { 1.189, 1.304e-2, -5.446e-5, 0.0 };
965 Evaluation logk0_CO2 = c[0] + temperatureCelcius * (c[1] + temperatureCelcius *
966 (c[2] + temperatureCelcius * c[3]));
967 Evaluation k0_CO2 = pow(10.0, logk0_CO2);
977 template <
class Evaluation>
978 OPM_HOST_DEVICE
static Evaluation equilibriumConstantH2O_(
const Evaluation& temperature,
const bool& highTemp)
980 Evaluation temperatureCelcius = temperature - 273.15;
981 std::array<Scalar, 5> c;
983 c = { -2.1077, 2.8127e-2, -8.4298e-5, 1.4969e-7, -1.1812e-10 };
986 c = { -2.209, 3.097e-2, -1.098e-4, 2.048e-7, 0.0 };
988 Evaluation logk0_H2O = c[0] + temperatureCelcius * (c[1] + temperatureCelcius * (c[2] +
989 temperatureCelcius * (c[3] + temperatureCelcius * c[4])));
990 return pow(10.0, logk0_H2O);