46 std::vector<Scalar> vaporPressure_;
48 std::vector<Scalar> minLiquidDensity__;
49 std::vector<Scalar> maxLiquidDensity__;
51 std::vector<Scalar> minGasDensity__;
52 std::vector<Scalar> maxGasDensity__;
56 std::vector<Scalar> gasEnthalpy_;
57 std::vector<Scalar> liquidEnthalpy_;
59 std::vector<Scalar> gasHeatCapacity_;
60 std::vector<Scalar> liquidHeatCapacity_;
62 std::vector<Scalar> gasDensity_;
63 std::vector<Scalar> liquidDensity_;
65 std::vector<Scalar> gasViscosity_;
66 std::vector<Scalar> liquidViscosity_;
68 std::vector<Scalar> gasThermalConductivity_;
69 std::vector<Scalar> liquidThermalConductivity_;
73 std::vector<Scalar> gasPressure_;
74 std::vector<Scalar> liquidPressure_;
89 void init(Scalar tempMin, Scalar tempMax,
unsigned nTemp,
90 Scalar pressMin, Scalar pressMax,
unsigned nPress)
101 vaporPressure_.resize(nTemp_);
102 minGasDensity__.resize(nTemp_);
103 maxGasDensity__.resize(nTemp_);
104 minLiquidDensity__.resize(nTemp_);
105 maxLiquidDensity__.resize(nTemp_);
107 gasEnthalpy_.resize(nTemp_*nPress_);
108 liquidEnthalpy_.resize(nTemp_*nPress_);
109 gasHeatCapacity_.resize(nTemp_*nPress_);
110 liquidHeatCapacity_.resize(nTemp_*nPress_);
111 gasDensity_.resize(nTemp_*nPress_);
112 liquidDensity_.resize(nTemp_*nPress_);
113 gasViscosity_.resize(nTemp_*nPress_);
114 liquidViscosity_.resize(nTemp_*nPress_);
115 gasThermalConductivity_.resize(nTemp_*nPress_);
116 liquidThermalConductivity_.resize(nTemp_*nPress_);
117 gasPressure_.resize(nTemp_*nDensity_);
118 liquidPressure_.resize(nTemp_*nDensity_);
141 using Scalar = ScalarT;
143 static constexpr bool isTabulated =
true;
155 static void init(Scalar tempMin, Scalar tempMax,
unsigned nTemp,
156 Scalar pressMin, Scalar pressMax,
unsigned nPress)
158 data_.init(tempMin, tempMax, nTemp, pressMin, pressMax, nPress);
160 assert(std::numeric_limits<Scalar>::has_quiet_NaN);
161 constexpr Scalar NaN = std::numeric_limits<Scalar>::quiet_NaN();
164 for (
unsigned iT = 0; iT < data_.nTemp_; ++ iT) {
165 const Scalar temperature = iT * (data_.tempMax_ - data_.tempMin_) / (data_.nTemp_ - 1) + data_.tempMin_;
168 data_.vaporPressure_[iT] = RawComponent::vaporPressure(temperature);
170 catch (
const std::exception&) {
171 data_.vaporPressure_[iT] = NaN;
174 const Scalar pgMax = maxGasPressure_(iT);
175 const Scalar pgMin = minGasPressure_(iT);
178 for (
unsigned iP = 0; iP < data_.nPress_; ++ iP) {
179 const Scalar pressure = iP * (pgMax - pgMin) / (data_.nPress_ - 1) + pgMin;
181 const unsigned i = iT + iP * data_.nTemp_;
184 data_.gasEnthalpy_[i] = RawComponent::gasEnthalpy(temperature, pressure);
186 catch (
const std::exception&) {
187 data_.gasEnthalpy_[i] = NaN;
191 data_.gasHeatCapacity_[i] = RawComponent::gasHeatCapacity(temperature, pressure);
193 catch (
const std::exception&) {
194 data_.gasHeatCapacity_[i] = NaN;
198 data_.gasDensity_[i] = RawComponent::gasDensity(temperature, pressure);
200 catch (
const std::exception&) {
201 data_.gasDensity_[i] = NaN;
205 data_.gasViscosity_[i] = RawComponent::gasViscosity(temperature, pressure);
207 catch (
const std::exception&) {
208 data_.gasViscosity_[i] = NaN;
212 data_.gasThermalConductivity_[i] = RawComponent::gasThermalConductivity(temperature, pressure);
214 catch (
const std::exception&) {
215 data_.gasThermalConductivity_[i] = NaN;
219 const Scalar plMin = minLiquidPressure_(iT);
220 const Scalar plMax = maxLiquidPressure_(iT);
221 for (
unsigned iP = 0; iP < data_.nPress_; ++ iP) {
222 Scalar pressure = iP * (plMax - plMin) / (data_.nPress_ - 1) + plMin;
224 const unsigned i = iT + iP*data_.nTemp_;
227 data_.liquidEnthalpy_[i] = RawComponent::liquidEnthalpy(temperature, pressure);
229 catch (
const std::exception&) {
230 data_.liquidEnthalpy_[i] = NaN;
234 data_.liquidHeatCapacity_[i] = RawComponent::liquidHeatCapacity(temperature, pressure);
236 catch (
const std::exception&) {
237 data_.liquidHeatCapacity_[i] = NaN;
241 data_.liquidDensity_[i] = RawComponent::liquidDensity(temperature, pressure);
243 catch (
const std::exception&) {
244 data_.liquidDensity_[i] = NaN;
248 data_.liquidViscosity_[i] = RawComponent::liquidViscosity(temperature, pressure);
250 catch (
const std::exception&) {
251 data_.liquidViscosity_[i] = NaN;
255 data_.liquidThermalConductivity_[i] = RawComponent::liquidThermalConductivity(temperature, pressure);
257 catch (
const std::exception&) {
258 data_.liquidThermalConductivity_[i] = NaN;
264 for (
unsigned iT = 0; iT < data_.nTemp_; ++ iT) {
265 const Scalar temperature = iT * (data_.tempMax_ - data_.tempMin_) / (data_.nTemp_ - 1) + data_.tempMin_;
269 data_.minGasDensity__[iT] = RawComponent::gasDensity(temperature, minGasPressure_(iT));
270 if (iT < data_.nTemp_ - 1) {
271 data_.maxGasDensity__[iT] = RawComponent::gasDensity(temperature, maxGasPressure_(iT + 1));
274 data_.maxGasDensity__[iT] = RawComponent::gasDensity(temperature, maxGasPressure_(iT));
278 for (
unsigned iRho = 0; iRho < data_.nDensity_; ++ iRho) {
279 const Scalar density =
280 Scalar(iRho) / (data_.nDensity_ - 1) *
281 (data_.maxGasDensity__[iT] - data_.minGasDensity__[iT])
283 data_.minGasDensity__[iT];
285 const unsigned i = iT + iRho * data_.nTemp_;
288 data_.gasPressure_[i] = RawComponent::gasPressure(temperature, density);
290 catch (
const std::exception&) {
291 data_.gasPressure_[i] = NaN;
297 data_.minLiquidDensity__[iT] = RawComponent::liquidDensity(temperature, minLiquidPressure_(iT));
298 if (iT < data_.nTemp_ - 1) {
299 data_.maxLiquidDensity__[iT] = RawComponent::liquidDensity(temperature, maxLiquidPressure_(iT + 1));
302 data_.maxLiquidDensity__[iT] = RawComponent::liquidDensity(temperature, maxLiquidPressure_(iT));
306 for (
unsigned iRho = 0; iRho < data_.nDensity_; ++ iRho) {
307 const Scalar density =
308 Scalar(iRho) / (data_.nDensity_ - 1) *
309 (data_.maxLiquidDensity__[iT] - data_.minLiquidDensity__[iT])
311 data_.minLiquidDensity__[iT];
313 const unsigned i = iT + iRho * data_.nTemp_;
316 data_.liquidPressure_[i] = RawComponent::liquidPressure(temperature, density);
318 catch (
const std::exception&) {
319 data_.liquidPressure_[i] = NaN;
329 {
return RawComponent::name(); }
335 {
return RawComponent::molarMass(); }
341 {
return RawComponent::criticalTemperature(); }
347 {
return RawComponent::criticalPressure(); }
353 {
throw std::runtime_error(
"Not implemented: acentricFactor of the component"); }
359 {
throw std::runtime_error(
"Not implemented: criticalVolume of the compoenent"); }
365 {
return RawComponent::tripleTemperature(); }
371 {
return RawComponent::triplePressure(); }
379 template <
class Evaluation>
382 const Evaluation& result = interpolateT_(data_.vaporPressure_, temperature);
383 if (std::isnan(scalarValue(result))) {
384 return RawComponent::vaporPressure(temperature);
395 template <
class Evaluation>
396 static Evaluation
gasEnthalpy(
const Evaluation& temperature,
const Evaluation& pressure)
398 const Evaluation& result = interpolateGasTP_(data_.gasEnthalpy_,
401 if (std::isnan(scalarValue(result))) {
402 return RawComponent::gasEnthalpy(temperature, pressure);
413 template <
class Evaluation>
414 static Evaluation
liquidEnthalpy(
const Evaluation& temperature,
const Evaluation& pressure)
416 const Evaluation& result = interpolateLiquidTP_(data_.liquidEnthalpy_,
419 if (std::isnan(scalarValue(result))) {
420 return RawComponent::liquidEnthalpy(temperature, pressure);
431 template <
class Evaluation>
432 static Evaluation
gasHeatCapacity(
const Evaluation& temperature,
const Evaluation& pressure)
434 const Evaluation& result = interpolateGasTP_(data_.gasHeatCapacity_,
437 if (std::isnan(scalarValue(result))) {
438 return RawComponent::gasHeatCapacity(temperature, pressure);
449 template <
class Evaluation>
452 const Evaluation& result = interpolateLiquidTP_(data_.liquidHeatCapacity_,
455 if (std::isnan(scalarValue(result))) {
456 return RawComponent::liquidHeatCapacity(temperature, pressure);
467 template <
class Evaluation>
477 template <
class Evaluation>
487 template <
class Evaluation>
488 static Evaluation
gasPressure(
const Evaluation& temperature, Scalar density)
490 const Evaluation& result = interpolateGasTRho_(data_.gasPressure_,
493 if (std::isnan(scalarValue(result))) {
494 return RawComponent::gasPressure(temperature,
506 template <
class Evaluation>
509 const Evaluation& result = interpolateLiquidTRho_(data_.liquidPressure_,
512 if (std::isnan(scalarValue(result))) {
513 return RawComponent::liquidPressure(temperature,
523 {
return RawComponent::gasIsCompressible(); }
529 {
return RawComponent::liquidIsCompressible(); }
535 {
return RawComponent::gasIsIdeal(); }
544 template <
class Evaluation>
545 static Evaluation
gasDensity(
const Evaluation& temperature,
const Evaluation& pressure)
547 const Evaluation& result = interpolateGasTP_(data_.gasDensity_,
550 if (std::isnan(scalarValue(result))) {
551 return RawComponent::gasDensity(temperature, pressure);
563 template <
class Evaluation>
564 static Evaluation
liquidDensity(
const Evaluation& temperature,
const Evaluation& pressure)
566 const Evaluation& result = interpolateLiquidTP_(data_.liquidDensity_,
569 if (std::isnan(scalarValue(result))) {
570 return RawComponent::liquidDensity(temperature, pressure);
581 template <
class Evaluation>
582 static Evaluation
gasViscosity(
const Evaluation& temperature,
const Evaluation& pressure)
584 const Evaluation& result = interpolateGasTP_(data_.gasViscosity_,
587 if (std::isnan(scalarValue(result))) {
588 return RawComponent::gasViscosity(temperature, pressure);
599 template <
class Evaluation>
600 static Evaluation
liquidViscosity(
const Evaluation& temperature,
const Evaluation& pressure)
602 const Evaluation& result = interpolateLiquidTP_(data_.liquidViscosity_,
605 if (std::isnan(scalarValue(result))) {
606 return RawComponent::liquidViscosity(temperature, pressure);
617 template <
class Evaluation>
620 const Evaluation& result = interpolateGasTP_(data_.gasThermalConductivity_,
623 if (std::isnan(scalarValue(result))) {
624 return RawComponent::gasThermalConductivity(temperature, pressure);
635 template <
class Evaluation>
638 const Evaluation& result = interpolateLiquidTP_(data_.liquidThermalConductivity_,
641 if (std::isnan(scalarValue(result))) {
642 return RawComponent::liquidThermalConductivity(temperature, pressure);
649 template <
class Evaluation>
650 static Evaluation interpolateT_(
const std::vector<Scalar>& values,
const Evaluation& T)
652 Evaluation alphaT = tempIdx_(T);
653 if (alphaT < 0 || alphaT >= data_.nTemp_ - 1) {
654 return std::numeric_limits<Scalar>::quiet_NaN();
657 const std::size_t iT =
static_cast<std::size_t
>(scalarValue(alphaT));
661 values[iT ]*(1 - alphaT) +
662 values[iT + 1]*( alphaT);
667 template <
class Evaluation>
668 static Evaluation interpolateLiquidTP_(
const std::vector<Scalar>& values,
672 Evaluation alphaT = tempIdx_(T);
673 if (alphaT < 0 || alphaT >= data_.nTemp_ - 1) {
674 return std::numeric_limits<Scalar>::quiet_NaN();
677 const std::size_t iT =
static_cast<std::size_t
>(scalarValue(alphaT));
680 Evaluation alphaP1 = pressLiquidIdx_(p, iT);
681 Evaluation alphaP2 = pressLiquidIdx_(p, iT + 1);
683 const std::size_t iP1 =
684 static_cast<std::size_t
>(
685 std::max<int>(0, std::min(
static_cast<int>(data_.nPress_) - 2,
686 static_cast<int>(scalarValue(alphaP1)))));
687 const std::size_t iP2 =
688 static_cast<std::size_t
>(
689 std::max(0, std::min(
static_cast<int>(data_.nPress_) - 2,
690 static_cast<int>(scalarValue(alphaP2)))));
695 values[(iT ) + (iP1 ) * data_.nTemp_] * (1 - alphaT) * (1 - alphaP1) +
696 values[(iT ) + (iP1 + 1) * data_.nTemp_] * (1 - alphaT) * ( alphaP1) +
697 values[(iT + 1) + (iP2 ) * data_.nTemp_] * ( alphaT) * (1 - alphaP2) +
698 values[(iT + 1) + (iP2 + 1) * data_.nTemp_] * ( alphaT) * ( alphaP2);
703 template <
class Evaluation>
704 static Evaluation interpolateGasTP_(
const std::vector<Scalar>& values,
708 Evaluation alphaT = tempIdx_(T);
709 if (alphaT < 0 || alphaT >= data_.nTemp_ - 1) {
710 return std::numeric_limits<Scalar>::quiet_NaN();
713 const std::size_t iT =
714 static_cast<std::size_t
>(
715 std::max(0, std::min(
static_cast<int>(data_.nTemp_) - 2,
716 static_cast<int>(scalarValue(alphaT)))));
719 Evaluation alphaP1 = pressGasIdx_(p, iT);
720 Evaluation alphaP2 = pressGasIdx_(p, iT + 1);
721 const std::size_t iP1 =
722 static_cast<std::size_t
>(
723 std::max(0, std::min(
static_cast<int>(data_.nPress_) - 2,
724 static_cast<int>(scalarValue(alphaP1)))));
725 const std::size_t iP2 =
726 static_cast<std::size_t
>(
727 std::max(0, std::min(
static_cast<int>(data_.nPress_) - 2,
728 static_cast<int>(scalarValue(alphaP2)))));
733 values[(iT ) + (iP1 ) * data_.nTemp_] * (1 - alphaT) * (1 - alphaP1) +
734 values[(iT ) + (iP1 + 1) * data_.nTemp_] * (1 - alphaT) * ( alphaP1) +
735 values[(iT + 1) + (iP2 ) * data_.nTemp_] * ( alphaT) * (1 - alphaP2) +
736 values[(iT + 1) + (iP2 + 1) * data_.nTemp_] * ( alphaT) * ( alphaP2);
741 template <
class Evaluation>
742 static Evaluation interpolateGasTRho_(
const std::vector<Scalar>& values,
744 const Evaluation& rho)
746 Evaluation alphaT = tempIdx_(T);
749 std::min(
static_cast<int>(data_.nTemp_ - 2),
750 static_cast<int>(alphaT)));
753 Evaluation alphaP1 = densityGasIdx_(rho, iT);
754 Evaluation alphaP2 = densityGasIdx_(rho, iT + 1);
757 std::min(
static_cast<int>(data_.nDensity_ - 2),
758 static_cast<int>(alphaP1)));
761 std::min(
static_cast<int>(data_.nDensity_ - 2),
762 static_cast<int>(alphaP2)));
767 values[(iT ) + (iP1 ) * data_.nTemp_] * (1 - alphaT) * (1 - alphaP1) +
768 values[(iT ) + (iP1 + 1) * data_.nTemp_] * (1 - alphaT) * ( alphaP1) +
769 values[(iT + 1) + (iP2 ) * data_.nTemp_] * ( alphaT) * (1 - alphaP2) +
770 values[(iT + 1) + (iP2 + 1) * data_.nTemp_] * ( alphaT) * ( alphaP2);
775 template <
class Evaluation>
776 static Evaluation interpolateLiquidTRho_(
const std::vector<Scalar>& values,
778 const Evaluation& rho)
780 Evaluation alphaT = tempIdx_(T);
781 const unsigned iT = std::max<int>(0, std::min<int>(data_.nTemp_ - 2,
static_cast<int>(alphaT)));
784 Evaluation alphaP1 = densityLiquidIdx_(rho, iT);
785 Evaluation alphaP2 = densityLiquidIdx_(rho, iT + 1);
786 const unsigned iP1 = std::max<int>(0, std::min<int>(data_.nDensity_ - 2,
static_cast<int>(alphaP1)));
787 const unsigned iP2 = std::max<int>(0, std::min<int>(data_.nDensity_ - 2,
static_cast<int>(alphaP2)));
792 values[(iT ) + (iP1 ) * data_.nTemp_] * (1 - alphaT) * (1 - alphaP1) +
793 values[(iT ) + (iP1 + 1) * data_.nTemp_] * (1 - alphaT) * ( alphaP1) +
794 values[(iT + 1) + (iP2 ) * data_.nTemp_] * ( alphaT) * (1 - alphaP2) +
795 values[(iT + 1) + (iP2 + 1) * data_.nTemp_] * ( alphaT) * ( alphaP2);
799 template <
class Evaluation>
800 static Evaluation tempIdx_(
const Evaluation& temperature)
802 return (data_.nTemp_ - 1) * (temperature - data_.tempMin_) / (data_.tempMax_ - data_.tempMin_);
806 template <
class Evaluation>
807 static Evaluation pressLiquidIdx_(
const Evaluation& pressure, std::size_t tempIdx)
809 const Scalar plMin = minLiquidPressure_(tempIdx);
810 const Scalar plMax = maxLiquidPressure_(tempIdx);
812 return (data_.nPress_ - 1) * (pressure - plMin) / (plMax - plMin);
816 template <
class Evaluation>
817 static Evaluation pressGasIdx_(
const Evaluation& pressure, std::size_t tempIdx)
819 const Scalar pgMin = minGasPressure_(tempIdx);
820 const Scalar pgMax = maxGasPressure_(tempIdx);
822 return (data_.nPress_ - 1) * (pressure - pgMin) / (pgMax - pgMin);
826 template <
class Evaluation>
827 static Evaluation densityLiquidIdx_(
const Evaluation& density, std::size_t tempIdx)
829 const Scalar densityMin = minLiquidDensity_(tempIdx);
830 const Scalar densityMax = maxLiquidDensity_(tempIdx);
831 return (data_.nDensity_ - 1) * (density - densityMin) / (densityMax - densityMin);
835 template <
class Evaluation>
836 static Evaluation densityGasIdx_(
const Evaluation& density, std::size_t tempIdx)
838 const Scalar densityMin = minGasDensity_(tempIdx);
839 const Scalar densityMax = maxGasDensity_(tempIdx);
840 return (data_.nDensity_ - 1) * (density - densityMin) / (densityMax - densityMin);
845 static Scalar minLiquidPressure_(std::size_t tempIdx)
847 if (!useVaporPressure) {
848 return data_.pressMin_;
851 return std::max<Scalar>(data_.pressMin_, data_.vaporPressure_[tempIdx] / 1.1);
857 static Scalar maxLiquidPressure_(std::size_t tempIdx)
859 if (!useVaporPressure) {
860 return data_.pressMax_;
863 return std::max<Scalar>(data_.pressMax_, data_.vaporPressure_[tempIdx] * 1.1);
869 static Scalar minGasPressure_(std::size_t tempIdx)
871 if (!useVaporPressure) {
872 return data_.pressMin_;
875 return std::min<Scalar>(data_.pressMin_, data_.vaporPressure_[tempIdx] / 1.1);
881 static Scalar maxGasPressure_(std::size_t tempIdx)
883 if (!useVaporPressure) {
884 return data_.pressMax_;
887 return std::min<Scalar>(data_.pressMax_, data_.vaporPressure_[tempIdx] * 1.1);
893 static Scalar minLiquidDensity_(std::size_t tempIdx)
894 {
return data_.minLiquidDensity__[tempIdx]; }
898 static Scalar maxLiquidDensity_(std::size_t tempIdx)
899 {
return data_.maxLiquidDensity__[tempIdx]; }
903 static Scalar minGasDensity_(std::size_t tempIdx)
904 {
return data_.minGasDensity__[tempIdx]; }
908 static Scalar maxGasDensity_(std::size_t tempIdx)
909 {
return data_.maxGasDensity__[tempIdx]; }
911 static TabulatedComponentData<Scalar> data_;