CoolProp 6.8.1dev
An open-source fluid property and humid air property database
IncompressibleFluid.cpp
Go to the documentation of this file.
1
2#include "DataStructures.h"
5#include "math.h"
6#include "MatrixMath.h"
7#include "PolyMath.h"
8#include <Eigen/Core>
9
10constexpr double INCOMP_EPSILON = DBL_EPSILON * 100.0;
11constexpr double INCOMP_DELTA = INCOMP_EPSILON * 10.0;
12
13namespace CoolProp {
14
16
21//IncompressibleFluid::IncompressibleFluid();
22
24 return;
25 // TODO: Implement validation function
26
27 // u and s have to be of the polynomial type!
28 //throw NotImplementedError("TODO");
29}
30
32 if (density.coeffs.cols() == 1) return true;
33 return false;
34}
35
37double IncompressibleFluid::baseExponential(IncompressibleData data, double y, double ybase) {
38 Eigen::VectorXd coeffs = makeVector(data.coeffs);
39 size_t r = coeffs.rows(), c = coeffs.cols();
40 if (strict && (r != 3 || c != 1)) {
41 throw ValueError(format("%s (%d): You have to provide a 3,1 matrix of coefficients, not (%d,%d).", __FILE__, __LINE__, r, c));
42 }
43
44 // Guard the function against zero denominators in exp((double)(coeffs[0] / ((y - ybase) + coeffs[1]) - coeffs[2]))
45 auto fnc = [&](double x) { return exp((double)(coeffs[0] / (x) - coeffs[2])); } ;
46 double x_den = (y - ybase) + coeffs[1];
47 double x_lo = -INCOMP_EPSILON;
48 double x_hi = +INCOMP_EPSILON;
49 if (x_den < x_lo || x_den > x_hi) {
50 return fnc(x_den);
51 }
52 // ... now we know that we are in the danger zone
53 // step away from the zero-crossing
54 x_lo -= INCOMP_DELTA;
55 x_hi += INCOMP_DELTA;
56 const double f_lo = fnc(x_lo);
57 const double f_hi = fnc(x_hi);
58 // Linearize around the zero-crossing
59 return (f_hi - f_lo) / (x_hi - x_lo) * (x_den - x_lo) + f_lo;
60}
61
63double IncompressibleFluid::baseLogexponential(IncompressibleData data, double y, double ybase) {
64 Eigen::VectorXd coeffs = makeVector(data.coeffs);
65 size_t r = coeffs.rows(), c = coeffs.cols();
66 if (strict && (r != 3 || c != 1)) {
67 throw ValueError(format("%s (%d): You have to provide a 3,1 matrix of coefficients, not (%d,%d).", __FILE__, __LINE__, r, c));
68 }
69
70 // Guard the function against zero denominators in exp((double)(log((double)(1.0 / ((y - ybase) + coeffs[0]) + 1.0 / ((y - ybase) + coeffs[0]) / ((y - ybase) + coeffs[0]))) * coeffs[1] + coeffs[2]))
71 auto fnc = [&](double x) { return exp((double)(log((double)(1.0 / (x) + 1.0 / (x) / (x))) * coeffs[1] + coeffs[2])); } ;
72 double x_den = (y - ybase) + coeffs[0];
73 double x_lo = -INCOMP_EPSILON;
74 double x_hi = +INCOMP_EPSILON;
75 if (x_den < x_lo || x_den > x_hi) {
76 return fnc(x_den);
77 }
78 // ... now we know that we are in the danger zone
79 // step away from the zero-crossing
80 x_lo -= INCOMP_DELTA;
81 x_hi += INCOMP_DELTA;
82 const double f_lo = fnc(x_lo);
83 const double f_hi = fnc(x_hi);
84 // Linearize around the zero-crossing
85 return (f_hi - f_lo) / (x_hi - x_lo) * (x_den - x_lo) + f_lo;
86}
87
89 size_t r = data.coeffs.rows(), c = data.coeffs.cols();
90 double offset = 0.0;
91 double in = 0.0;
92 Eigen::MatrixXd coeffs;
93 if (r > 0 && c > 0) {
94 offset = data.coeffs(0, 0);
95 if (r == 1 && c > 1) { // row vector -> function of z
96 coeffs = Eigen::MatrixXd(data.coeffs.block(0, 1, r, c - 1));
97 in = z;
98 } else if (r > 1 && c == 1) { // column vector -> function of y
99 coeffs = Eigen::MatrixXd(data.coeffs.block(1, 0, r - 1, c));
100 in = y;
101 } else {
102 throw ValueError(format("%s (%d): You have to provide a vector (1D matrix) of coefficients, not (%d,%d).", __FILE__, __LINE__, r, c));
103 }
104 return poly.evaluate(coeffs, in, 0, offset);
105 }
106 throw ValueError(format("%s (%d): You have to provide a vector (1D matrix) of coefficients, not (%d,%d).", __FILE__, __LINE__, r, c));
107}
108
110double IncompressibleFluid::rho(double T, double p, double x) {
111 switch (density.type) {
113 return poly.evaluate(density.coeffs, T, x, 0, 0, Tbase, xbase);
115 return baseExponential(density, T, 0.0);
117 return baseLogexponential(density, T, 0.0);
119 return exp(poly.evaluate(density.coeffs, T, x, 0, 0, Tbase, xbase));
121 return basePolyOffset(density, T, x);
123 throw ValueError(format("%s (%d): The function type is not specified (\"[%d]\"), are you sure the coefficients have been set?", __FILE__,
124 __LINE__, density.type));
125 default:
126 throw ValueError(format("%s (%d): Your function type \"[%d]\" is unknown.", __FILE__, __LINE__, density.type));
127 }
128}
129
131double IncompressibleFluid::c(double T, double p, double x) {
132 switch (specific_heat.type) {
134 //throw NotImplementedError("Here you should implement the polynomial.");
135 return poly.evaluate(specific_heat.coeffs, T, x, 0, 0, Tbase, xbase);
137 throw ValueError(format("%s (%d): The function type is not specified (\"[%d]\"), are you sure the coefficients have been set?", __FILE__,
138 __LINE__, specific_heat.type));
139 default:
140 throw ValueError(format("%s (%d): There is no predefined way to use this function type \"[%d]\" for specific heat.", __FILE__, __LINE__,
142 }
143}
144
146double IncompressibleFluid::visc(double T, double p, double x) {
147 switch (viscosity.type) {
149 return poly.evaluate(viscosity.coeffs, T, x, 0, 0, Tbase, xbase);
151 return baseExponential(viscosity, T, 0.0);
153 return baseLogexponential(viscosity, T, 0.0);
155 return exp(poly.evaluate(viscosity.coeffs, T, x, 0, 0, Tbase, xbase));
157 return basePolyOffset(viscosity, T, x);
159 throw ValueError(format("%s (%d): The function type is not specified (\"[%d]\"), are you sure the coefficients have been set?", __FILE__,
160 __LINE__, viscosity.type));
161 default:
162 throw ValueError(format("%s (%d): Your function type \"[%d]\" is unknown.", __FILE__, __LINE__, viscosity.type));
163 }
164}
166double IncompressibleFluid::cond(double T, double p, double x) {
167 switch (conductivity.type) {
169 return poly.evaluate(conductivity.coeffs, T, x, 0, 0, Tbase, xbase);
171 return baseExponential(conductivity, T, 0.0);
173 return baseLogexponential(conductivity, T, 0.0);
175 return exp(poly.evaluate(conductivity.coeffs, T, x, 0, 0, Tbase, xbase));
177 return basePolyOffset(conductivity, T, x);
179 throw ValueError(format("%s (%d): The function type is not specified (\"[%d]\"), are you sure the coefficients have been set?", __FILE__,
180 __LINE__, conductivity.type));
181 default:
182 throw ValueError(format("%s (%d): Your function type \"[%d]\" is unknown.", __FILE__, __LINE__, conductivity.type));
183 }
184}
186double IncompressibleFluid::psat(double T, double x) {
187 if (T <= this->TminPsat) return 0.0;
188 switch (p_sat.type) {
190 return poly.evaluate(p_sat.coeffs, T, x, 0, 0, Tbase, xbase);
192 return baseExponential(p_sat, T, 0.0);
194 return baseLogexponential(p_sat, T, 0.0);
196 return exp(poly.evaluate(p_sat.coeffs, T, x, 0, 0, Tbase, xbase));
198 return basePolyOffset(p_sat, T, x);
200 throw ValueError(format("%s (%d): The function type is not specified (\"[%d]\"), are you sure the coefficients have been set?", __FILE__,
201 __LINE__, p_sat.type));
202 default:
203 throw ValueError(format("%s (%d): Your function type \"[%d]\" is unknown.", __FILE__, __LINE__, p_sat.type));
204 }
205}
207double IncompressibleFluid::Tfreeze(double p, double x) {
208 switch (T_freeze.type) {
210 return poly.evaluate(T_freeze.coeffs, p, x, 0, 0, 0.0, xbase);
212 return baseExponential(T_freeze, x, 0.0);
214 return baseLogexponential(T_freeze, x, 0.0);
216 return exp(poly.evaluate(T_freeze.coeffs, p, x, 0, 0, 0.0, xbase));
218 return basePolyOffset(T_freeze, p, x);
220 throw ValueError(format("%s (%d): The function type is not specified (\"[%d]\"), are you sure the coefficients have been set?", __FILE__,
221 __LINE__, T_freeze.type));
222 default:
223 throw ValueError(format("%s (%d): Your function type \"[%d]\" is unknown.", __FILE__, __LINE__, T_freeze.type));
224 }
225}
226
227/* Below are direct calculations of the derivatives. Nothing
228 * special is going on, we simply use the polynomial class to
229 * derive the different functions with respect to temperature.
230 */
232double IncompressibleFluid::drhodTatPx(double T, double p, double x) {
233 switch (density.type) {
235 return poly.derivative(density.coeffs, T, x, 0, 0, 0, Tbase, xbase);
237 throw ValueError(format("%s (%d): The function type is not specified (\"[%d]\"), are you sure the coefficients have been set?", __FILE__,
238 __LINE__, density.type));
239 default:
240 throw ValueError(
241 format("%s (%d): There is no predefined way to use this function type \"[%d]\" for density.", __FILE__, __LINE__, density.type));
242 }
243}
245// with respect to temperature at constant pressure and composition
246// integrated in temperature
247double IncompressibleFluid::dsdTatPxdT(double T, double p, double x) {
248 switch (specific_heat.type) {
250 return poly.integral(specific_heat.coeffs, T, x, 0, -1, 0, Tbase, xbase);
252 throw ValueError(format("%s (%d): The function type is not specified (\"[%d]\"), are you sure the coefficients have been set?", __FILE__,
253 __LINE__, specific_heat.type));
254 default:
255 throw ValueError(
256 format("%s (%d): There is no predefined way to use this function type \"[%d]\" for entropy.", __FILE__, __LINE__, specific_heat.type));
257 }
258}
260// with respect to temperature at constant pressure and composition
261// integrated in temperature
262double IncompressibleFluid::dhdTatPxdT(double T, double p, double x) {
263 switch (specific_heat.type) {
265 return poly.integral(specific_heat.coeffs, T, x, 0, 0, 0, Tbase, xbase);
267 throw ValueError(format("%s (%d): The function type is not specified (\"[%d]\"), are you sure the coefficients have been set?", __FILE__,
268 __LINE__, specific_heat.type));
269 default:
270 throw ValueError(
271 format("%s (%d): There is no predefined way to use this function type \"[%d]\" for entropy.", __FILE__, __LINE__, specific_heat.type));
272 }
273}
274
276
278double IncompressibleFluid::inputFromMass(double T, double x) {
279 if (this->xid == IFRAC_PURE) {
280 return _HUGE;
281 } else if (this->xid == IFRAC_MASS) {
282 return x;
283 } else {
284 throw NotImplementedError("Mass composition conversion has not been implemented.");
285 //switch (mass2input.type) {
286 // case IncompressibleData::INCOMPRESSIBLE_POLYNOMIAL:
287 // return poly.evaluate(mass2input.coeffs, T, x, 0, 0, 0.0, 0.0); // TODO: make sure Tbase and xbase are defined in the correct way
288 // break;
289 // case IncompressibleData::INCOMPRESSIBLE_EXPONENTIAL:
290 // return baseExponential(mass2input, x, 0.0);
291 // break;
292 // case IncompressibleData::INCOMPRESSIBLE_LOGEXPONENTIAL:
293 // return baseLogexponential(mass2input, x, 0.0);
294 // break;
295 // case IncompressibleData::INCOMPRESSIBLE_EXPPOLYNOMIAL:
296 // return exp(poly.evaluate(mass2input.coeffs, T, x, 0, 0, 0.0, 0.0)); // TODO: make sure Tbase and xbase are defined in the correct way
297 // break;
298 // case IncompressibleData::INCOMPRESSIBLE_POLYOFFSET:
299 // return basePolyOffset(mass2input, T, x);
300 // break;
301 // case IncompressibleData::INCOMPRESSIBLE_NOT_SET:
302 // throw ValueError(format("%s (%d): The function type is not specified (\"[%d]\"), are you sure the coefficients have been set?",__FILE__,__LINE__,mass2input.type));
303 // break;
304 // default:
305 // throw ValueError(format("%s (%d): Your function type \"[%d]\" is unknown.",__FILE__,__LINE__,mass2input.type));
306 // break;
307 //}
308 //return _HUGE;
309 }
310}
311
313
315double IncompressibleFluid::inputFromVolume(double T, double x) {
316 if (this->xid == IFRAC_PURE) {
317 return _HUGE;
318 } else if (this->xid == IFRAC_VOLUME) {
319 return x;
320 } else {
321 throw NotImplementedError("Volume composition conversion has not been implemented.");
322 //switch (volume2input.type) {
323 // case IncompressibleData::INCOMPRESSIBLE_POLYNOMIAL:
324 // return poly.evaluate(volume2input.coeffs, T, x, 0, 0, 0.0, 0.0); // TODO: make sure Tbase and xbase are defined in the correct way
325 // break;
326 // case IncompressibleData::INCOMPRESSIBLE_EXPONENTIAL:
327 // return baseExponential(volume2input, x, 0.0);
328 // break;
329 // case IncompressibleData::INCOMPRESSIBLE_LOGEXPONENTIAL:
330 // return baseLogexponential(volume2input, x, 0.0);
331 // break;
332 // case IncompressibleData::INCOMPRESSIBLE_EXPPOLYNOMIAL:
333 // return exp(poly.evaluate(volume2input.coeffs, T, x, 0, 0, 0.0, 0.0)); // TODO: make sure Tbase and xbase are defined in the correct way
334 // break;
335 // case IncompressibleData::INCOMPRESSIBLE_POLYOFFSET:
336 // return basePolyOffset(volume2input, T, x);
337 // break;
338 // case IncompressibleData::INCOMPRESSIBLE_NOT_SET:
339 // throw ValueError(format("%s (%d): The function type is not specified (\"[%d]\"), are you sure the coefficients have been set?",__FILE__,__LINE__,volume2input.type));
340 // break;
341 // default:
342 // throw ValueError(format("%s (%d): Your function type \"[%d]\" is unknown.",__FILE__,__LINE__,volume2input.type));
343 // break;
344 //}
345 //return _HUGE;
346 }
347}
348
350
352double IncompressibleFluid::inputFromMole(double T, double x) {
353 if (this->xid == IFRAC_PURE) {
354 return _HUGE;
355 } else if (this->xid == IFRAC_MOLE) {
356 return x;
357 } else {
358 throw NotImplementedError("Mole composition conversion has not been implemented.");
359 /*
360 switch (mole2input.type) {
361 case IncompressibleData::INCOMPRESSIBLE_POLYNOMIAL:
362 return poly.evaluate(mole2input.coeffs, T, x, 0, 0, 0.0, 0.0); // TODO: make sure Tbase and xbase are defined in the correct way
363 break;
364 case IncompressibleData::INCOMPRESSIBLE_EXPONENTIAL:
365 return baseExponential(mole2input, x, 0.0);
366 break;
367 case IncompressibleData::INCOMPRESSIBLE_LOGEXPONENTIAL:
368 return baseLogexponential(mole2input, x, 0.0);
369 break;
370 case IncompressibleData::INCOMPRESSIBLE_EXPPOLYNOMIAL:
371 return exp(poly.evaluate(mole2input.coeffs, T, x, 0, 0, 0.0, 0.0)); // TODO: make sure Tbase and xbase are defined in the correct way
372 break;
373 case IncompressibleData::INCOMPRESSIBLE_POLYOFFSET:
374 return basePolyOffset(mole2input, T, x);
375 break;
376 case IncompressibleData::INCOMPRESSIBLE_NOT_SET:
377 throw ValueError(format("%s (%d): The function type is not specified (\"[%d]\"), are you sure the coefficients have been set?",__FILE__,__LINE__,mole2input.type));
378 break;
379 default:
380 throw ValueError(format("%s (%d): Your function type \"[%d]\" is unknown.",__FILE__,__LINE__,mole2input.type));
381 break;
382 }
383 return _HUGE;
384 */
385 }
386}
387
388/* Some functions can be inverted directly, those are listed
389 * here. It is also possible to solve for other quantities, but
390 * that involves some more sophisticated processing and is not
391 * done here, but in the backend, T(h,p) for example.
392 */
394double IncompressibleFluid::T_rho(double Dmass, double p, double x) {
395 double d_raw = Dmass; // No changes needed, no reference values...
396 switch (density.type) {
398 return poly.solve_limits(density.coeffs, x, d_raw, Tmin, Tmax, 0, 0, 0, Tbase, xbase);
400 throw ValueError(format("%s (%d): The function type is not specified (\"[%d]\"), are you sure the coefficients have been set?", __FILE__,
401 __LINE__, specific_heat.type));
402 default:
403 throw ValueError(format("%s (%d): There is no predefined way to use this function type \"[%d]\" for inverse density.", __FILE__, __LINE__,
405 }
406}
408double IncompressibleFluid::T_c(double Cmass, double p, double x) {
409 double c_raw = Cmass; // No changes needed, no reference values...
410 switch (specific_heat.type) {
412 return poly.solve_limits(specific_heat.coeffs, x, c_raw, Tmin, Tmax, 0, 0, 0, Tbase, xbase);
414 throw ValueError(format("%s (%d): The function type is not specified (\"[%d]\"), are you sure the coefficients have been set?", __FILE__,
415 __LINE__, specific_heat.type));
416 default:
417 throw ValueError(format("%s (%d): There is no predefined way to use this function type \"[%d]\" for inverse specific heat.", __FILE__,
418 __LINE__, specific_heat.type));
419 }
420}
421
422/*
423 * Some more functions to provide a single implementation
424 * of important routines.
425 * We start with the check functions that can validate input
426 * in terms of pressure p, temperature T and composition x.
427 */
429
432bool IncompressibleFluid::checkT(double T, double p, double x) {
433 if (Tmin <= 0.) throw ValueError("Please specify the minimum temperature.");
434 if (Tmax <= 0.) throw ValueError("Please specify the maximum temperature.");
435 if ((Tmin > T) || (T > Tmax)) throw ValueError(format("Your temperature %f is not between %f and %f.", T, Tmin, Tmax));
436 double TF = 0.0;
438 if (T < TF) throw ValueError(format("Your temperature %f is below the freezing point of %f.", T, TF));
439 return true;
440}
441
443
449bool IncompressibleFluid::checkP(double T, double p, double x) {
450 double ps = 0.0;
452 if (p < 0.0) throw ValueError(format("You cannot use negative pressures: %f < %f. ", p, 0.0));
453 if (ps > 0.0 && p < ps) throw ValueError(format("Equations are valid for liquid phase only: %f < %f (psat). ", p, ps));
454 return true;
455}
456
458
462 if (xmin < 0.0 || xmin > 1.0) throw ValueError("Please specify the minimum concentration between 0 and 1.");
463 if (xmax < 0.0 || xmax > 1.0) throw ValueError("Please specify the maximum concentration between 0 and 1.");
464 if ((x < xmin * (1 - INCOMP_EPSILON)) || (x > xmax * (1 + INCOMP_EPSILON))) {
465 throw ValueError(format("Your composition %g is not between %g and %g.", x, xmin, xmax));
466 }
467 return true;
468}
469
470} /* namespace CoolProp */
471
472// Testing still needs to be enhanced.
473/* Below, I try to carry out some basic tests for both 2D and 1D
474 * polynomials as well as the exponential functions for vapour
475 * pressure etc.
476 */
477#ifdef ENABLE_CATCH
478# include <math.h>
479# include <iostream>
480# include <catch2/catch_all.hpp>
481# include "TestObjects.h"
482
483Eigen::MatrixXd makeMatrix(const std::vector<double>& coefficients) {
484 //IncompressibleClass::checkCoefficients(coefficients,18);
485 std::vector<std::vector<double>> matrix;
486 std::vector<double> tmpVector;
487
488 tmpVector.clear();
489 tmpVector.push_back(coefficients[0]);
490 tmpVector.push_back(coefficients[6]);
491 tmpVector.push_back(coefficients[11]);
492 tmpVector.push_back(coefficients[15]);
493 matrix.push_back(tmpVector);
494
495 tmpVector.clear();
496 tmpVector.push_back(coefficients[1] * 100.0);
497 tmpVector.push_back(coefficients[7] * 100.0);
498 tmpVector.push_back(coefficients[12] * 100.0);
499 tmpVector.push_back(coefficients[16] * 100.0);
500 matrix.push_back(tmpVector);
501
502 tmpVector.clear();
503 tmpVector.push_back(coefficients[2] * 100.0 * 100.0);
504 tmpVector.push_back(coefficients[8] * 100.0 * 100.0);
505 tmpVector.push_back(coefficients[13] * 100.0 * 100.0);
506 tmpVector.push_back(coefficients[17] * 100.0 * 100.0);
507 matrix.push_back(tmpVector);
508
509 tmpVector.clear();
510 tmpVector.push_back(coefficients[3] * 100.0 * 100.0 * 100.0);
511 tmpVector.push_back(coefficients[9] * 100.0 * 100.0 * 100.0);
512 tmpVector.push_back(coefficients[14] * 100.0 * 100.0 * 100.0);
513 tmpVector.push_back(0.0);
514 matrix.push_back(tmpVector);
515
516 tmpVector.clear();
517 tmpVector.push_back(coefficients[4] * 100.0 * 100.0 * 100.0 * 100.0);
518 tmpVector.push_back(coefficients[10] * 100.0 * 100.0 * 100.0 * 100.0);
519 tmpVector.push_back(0.0);
520 tmpVector.push_back(0.0);
521 matrix.push_back(tmpVector);
522
523 tmpVector.clear();
524 tmpVector.push_back(coefficients[5] * 100.0 * 100.0 * 100.0 * 100.0 * 100.0);
525 tmpVector.push_back(0.0);
526 tmpVector.push_back(0.0);
527 tmpVector.push_back(0.0);
528 matrix.push_back(tmpVector);
529
530 tmpVector.clear();
531 return CoolProp::vec_to_eigen(matrix).transpose();
532}
533
534TEST_CASE("Internal consistency checks and example use cases for the incompressible fluids", "[IncompressibleFluids]") {
535 bool PRINT = false;
536 std::string tmpStr;
537 std::vector<double> tmpVector;
538 std::vector<std::vector<double>> tmpMatrix;
539
540 SECTION("Test case for \"SylthermXLT\" by Dow Chemicals") {
541
542 std::vector<double> cRho;
543 cRho.push_back(+1.1563685145E+03);
544 cRho.push_back(-1.0269048032E+00);
545 cRho.push_back(-9.3506079577E-07);
546 cRho.push_back(+1.0368116627E-09);
549 density.coeffs = CoolProp::vec_to_eigen(cRho);
550
551 std::vector<double> cHeat;
552 cHeat.push_back(+1.1562261074E+03);
553 cHeat.push_back(+2.0994549103E+00);
554 cHeat.push_back(+7.7175381057E-07);
555 cHeat.push_back(-3.7008444051E-20);
556 CoolProp::IncompressibleData specific_heat;
558 specific_heat.coeffs = CoolProp::vec_to_eigen(cHeat);
559
560 std::vector<double> cCond;
561 cCond.push_back(+1.6121957379E-01);
562 cCond.push_back(-1.3023781944E-04);
563 cCond.push_back(-1.4395238766E-07);
564 CoolProp::IncompressibleData conductivity;
566 conductivity.coeffs = CoolProp::vec_to_eigen(cCond);
567
568 std::vector<double> cVisc;
569 cVisc.push_back(+1.0337654989E+03);
570 cVisc.push_back(-4.3322764383E+01);
571 cVisc.push_back(+1.0715062356E+01);
574 viscosity.coeffs = CoolProp::vec_to_eigen(cVisc);
575
577 XLT.setName("XLT");
578 XLT.setDescription("SylthermXLT");
579 XLT.setReference("Dow Chemicals data sheet");
580 XLT.setTmax(533.15);
581 XLT.setTmin(173.15);
582 XLT.setxmax(0.0);
583 XLT.setxmin(0.0);
584 XLT.setTminPsat(533.15);
585
586 XLT.setTbase(0.0);
587 XLT.setxbase(0.0);
588
590 XLT.setDensity(density);
591 XLT.setSpecificHeat(specific_heat);
592 XLT.setViscosity(viscosity);
593 XLT.setConductivity(conductivity);
594 //XLT.setPsat(parse_coefficients(fluid_json, "saturation_pressure", false));
595 //XLT.setTfreeze(parse_coefficients(fluid_json, "T_freeze", false));
596 //XLT.setVolToMass(parse_coefficients(fluid_json, "volume2mass", false));
597 //XLT.setMassToMole(parse_coefficients(fluid_json, "mass2mole", false));
598
600 //XLT.validate();
601 double acc = 0.0001;
602 double val = 0;
603 double res = 0;
604
605 // Prepare the results and compare them to the calculated values
606 double T = 273.15 + 50;
607 double p = 10e5;
608 double x = 0.0;
609
610 // Compare density
611 val = 824.4615702148608;
612 res = XLT.rho(T, p, x);
613 {
614 CAPTURE(T);
615 CAPTURE(val);
616 CAPTURE(res);
617 CHECK(check_abs(val, res, acc));
618 }
619
620 // Compare cp
621 val = 1834.7455527670554;
622 res = XLT.c(T, p, x);
623 {
624 CAPTURE(T);
625 CAPTURE(val);
626 CAPTURE(res);
627 CHECK(check_abs(val, res, acc));
628 }
629
630 // Check property functions
631 CHECK_THROWS(XLT.s(T, p, x));
632 CHECK_THROWS(XLT.h(T, p, x));
633 CHECK_THROWS(XLT.u(T, p, x));
634
635 // Compare v
636 val = 0.0008931435169681835;
637 res = XLT.visc(T, p, x);
638 {
639 CAPTURE(T);
640 CAPTURE(val);
641 CAPTURE(res);
642 CHECK(check_abs(val, res, acc));
643 }
644
645 // Compare l
646 val = 0.10410086156049088;
647 res = XLT.cond(T, p, x);
648 {
649 CAPTURE(T);
650 CAPTURE(val);
651 CAPTURE(res);
652 CHECK(check_abs(val, res, acc));
653 }
654 }
655
656 SECTION("Test case for Methanol from SecCool") {
657
659
660 // Prepare the results and compare them to the calculated values
661 double acc = 0.0001;
662 double T = 273.15 + 10;
663 double p = 10e5;
664 double x = 0.25;
665 double expected = 0;
666 double actual = 0;
667
668 // Compare density
669 expected = 963.2886528091547;
670 actual = CH3OH.rho(T, p, x);
671 {
672 CAPTURE(T);
673 CAPTURE(p);
674 CAPTURE(x);
675 CAPTURE(expected);
676 CAPTURE(actual);
677 CHECK(check_abs(expected, actual, acc));
678 }
679
680 // Compare cp
681 expected = 3993.9748117022423;
682 actual = CH3OH.c(T, p, x);
683 {
684 CAPTURE(T);
685 CAPTURE(p);
686 CAPTURE(x);
687 CAPTURE(expected);
688 CAPTURE(actual);
689 CHECK(check_abs(expected, actual, acc));
690 }
691
692 // Check property functions
693 CHECK_THROWS(CH3OH.s(T, p, x));
694 CHECK_THROWS(CH3OH.h(T, p, x));
695 CHECK_THROWS(CH3OH.u(T, p, x));
696
697 // Compare viscosity
698 expected = 0.0023970245009602097; // Pa-s // Melinder 3.48 mPas @ 24.87 wt%, 0C
699 actual = CH3OH.visc(T, p, x);
700 {
701 CAPTURE(T);
702 CAPTURE(p);
703 CAPTURE(x);
704 CAPTURE(expected);
705 CAPTURE(actual);
706 std::string errmsg = CoolProp::get_global_param_string("errstring");
707 CAPTURE(errmsg);
708 CHECK(check_abs(expected, actual, acc));
709 }
710
711 // Compare conductivity
712 expected = 0.44791148414693727;
713 actual = CH3OH.cond(T, p, x);
714 {
715 CAPTURE(T);
716 CAPTURE(p);
717 CAPTURE(x);
718 CAPTURE(expected);
719 CAPTURE(actual);
720 std::string errmsg = CoolProp::get_global_param_string("errstring");
721 CAPTURE(errmsg);
722 CHECK(check_abs(expected, actual, acc));
723 }
724
725 // Compare Tfreeze
726 expected = -20.02 + 273.15; // 253.1293105454671;
727 actual = CH3OH.Tfreeze(p, x);
728 {
729 CAPTURE(T);
730 CAPTURE(p);
731 CAPTURE(x);
732 CAPTURE(expected);
733 CAPTURE(actual);
734 std::string errmsg = CoolProp::get_global_param_string("errstring");
735 CAPTURE(errmsg);
736 CHECK(check_abs(expected, actual, acc));
737 }
738 }
739}
740
741#endif /* ENABLE_CATCH */