Line data Source code
1 : //* This file is part of the MOOSE framework
2 : //* https://mooseframework.inl.gov
3 : //*
4 : //* All rights reserved, see COPYRIGHT for full restrictions
5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 : //*
7 : //* Licensed under LGPL 2.1, please see LICENSE for details
8 : //* https://www.gnu.org/licenses/lgpl-2.1.html
9 :
10 : #include "IdealGasFluidProperties.h"
11 : #include "Conversion.h"
12 :
13 : #include "metaphysicl/raw_type.h"
14 :
15 : registerMooseObject("FluidPropertiesApp", IdealGasFluidProperties);
16 :
17 : InputParameters
18 1706 : IdealGasFluidProperties::validParams()
19 : {
20 1706 : InputParameters params = SinglePhaseFluidProperties::validParams();
21 1706 : params += NaNInterface::validParams();
22 :
23 3412 : params.addRangeCheckedParam<Real>("gamma", 1.4, "gamma > 1", "gamma value (cp/cv)");
24 3412 : params.addParam<Real>("molar_mass", 29.0e-3, "Constant molar mass of the fluid (kg/mol)");
25 3412 : params.addParam<Real>("e_ref", 0, "Reference specific internal energy [J/kg]");
26 3412 : params.addParam<Real>("mu", 18.23e-6, "Dynamic viscosity, Pa.s");
27 3412 : params.addParam<Real>("k", 25.68e-3, "Thermal conductivity, W/(m-K)");
28 3412 : params.addParam<Real>("T_c", 0, "Critical temperature, K");
29 3412 : params.addParam<Real>("rho_c", 0, "Critical density, kg/m3");
30 3412 : params.addParam<Real>("e_c", 0, "Internal energy at the critical point, J/kg");
31 :
32 1706 : params.addClassDescription("Fluid properties for an ideal gas");
33 :
34 1706 : return params;
35 0 : }
36 :
37 944 : IdealGasFluidProperties::IdealGasFluidProperties(const InputParameters & parameters)
38 : : SinglePhaseFluidProperties(parameters),
39 : NaNInterface(this),
40 :
41 944 : _gamma(getParam<Real>("gamma")),
42 1888 : _molar_mass(getParam<Real>("molar_mass")),
43 1888 : _e_ref(getParam<Real>("e_ref")),
44 :
45 944 : _R_specific(_R / _molar_mass),
46 944 : _cp(_gamma * _R_specific / (_gamma - 1.0)),
47 944 : _cv(_cp / _gamma),
48 :
49 1888 : _mu(getParam<Real>("mu")),
50 1888 : _k(getParam<Real>("k")),
51 :
52 1888 : _T_c(getParam<Real>("T_c")),
53 1888 : _rho_c(getParam<Real>("rho_c")),
54 2832 : _e_c(getParam<Real>("e_c"))
55 : {
56 944 : }
57 :
58 1816 : IdealGasFluidProperties::~IdealGasFluidProperties() {}
59 :
60 : std::string
61 1 : IdealGasFluidProperties::fluidName() const
62 : {
63 1 : return "ideal_gas";
64 : }
65 :
66 : Real
67 61581 : IdealGasFluidProperties::p_from_v_e(Real v, Real e) const
68 : {
69 61581 : if (v == 0.0)
70 0 : return getNaN("Invalid value of specific volume detected (v = " + Moose::stringify(v) + ").");
71 :
72 61581 : return (_gamma - 1.0) * (e - _e_ref) / v;
73 : }
74 :
75 : ADReal
76 1 : IdealGasFluidProperties::p_from_v_e(const ADReal & v, const ADReal & e) const
77 : {
78 1 : if (v.value() == 0.0)
79 0 : return getNaN("Invalid value of specific volume detected (v = " + Moose::stringify(v.value()) +
80 0 : ").");
81 :
82 2 : return (_gamma - 1.0) * (e - _e_ref) / v;
83 : }
84 :
85 : void
86 4 : IdealGasFluidProperties::p_from_v_e(Real v, Real e, Real & p, Real & dp_dv, Real & dp_de) const
87 : {
88 4 : p = p_from_v_e(v, e);
89 4 : dp_dv = -(_gamma - 1.0) * (e - _e_ref) / v / v;
90 4 : dp_de = (_gamma - 1.0) / v;
91 4 : }
92 :
93 : void
94 0 : IdealGasFluidProperties::p_from_v_e(
95 : const ADReal & v, const ADReal & e, ADReal & p, ADReal & dp_dv, ADReal & dp_de) const
96 : {
97 0 : p = p_from_v_e(v, e);
98 0 : dp_dv = -(_gamma - 1.0) * (e - _e_ref) / v / v;
99 0 : dp_de = (_gamma - 1.0) / v;
100 0 : }
101 :
102 : Real
103 62251 : IdealGasFluidProperties::T_from_v_e(Real /*v*/, Real e) const
104 : {
105 62251 : return (e - _e_ref) / _cv;
106 : }
107 :
108 : ADReal
109 0 : IdealGasFluidProperties::T_from_v_e(const ADReal & /*v*/, const ADReal & e) const
110 : {
111 0 : return (e - _e_ref) / _cv;
112 : }
113 :
114 : void
115 5 : IdealGasFluidProperties::T_from_v_e(Real v, Real e, Real & T, Real & dT_dv, Real & dT_de) const
116 : {
117 5 : T = T_from_v_e(v, e);
118 5 : dT_dv = 0.0;
119 5 : dT_de = 1.0 / _cv;
120 5 : }
121 :
122 : void
123 0 : IdealGasFluidProperties::T_from_v_e(
124 : const ADReal & v, const ADReal & e, ADReal & T, ADReal & dT_dv, ADReal & dT_de) const
125 : {
126 0 : T = T_from_v_e(v, e);
127 0 : dT_dv = 0.0;
128 0 : dT_de = 1.0 / _cv;
129 0 : }
130 :
131 : Real
132 20604 : IdealGasFluidProperties::c_from_v_e(Real v, Real e) const
133 : {
134 20604 : Real T = T_from_v_e(v, e);
135 :
136 20604 : Real c2 = _gamma * _R_specific * T;
137 20604 : if (c2 < 0)
138 : {
139 0 : c2 = 0;
140 0 : flagInvalidSolution(
141 : "Sound speed squared (gamma * R * T) is negative: c2 = " + Moose::stringify(c2) + ".");
142 : }
143 :
144 20604 : return std::sqrt(c2);
145 : }
146 :
147 : ADReal
148 0 : IdealGasFluidProperties::c_from_v_e(const ADReal & v, const ADReal & e) const
149 : {
150 0 : const auto T = T_from_v_e(v, e);
151 :
152 0 : auto c2 = _gamma * _R_specific * T;
153 0 : if (MetaPhysicL::raw_value(c2) < 0)
154 : {
155 0 : c2 = 0;
156 0 : flagInvalidSolution(
157 : "Sound speed squared (gamma * R * T) is negative: c2 = " + Moose::stringify(c2) + ".");
158 : }
159 :
160 0 : return std::sqrt(c2);
161 : }
162 :
163 : void
164 2 : IdealGasFluidProperties::c_from_v_e(Real v, Real e, Real & c, Real & dc_dv, Real & dc_de) const
165 : {
166 : Real T, dT_dv, dT_de;
167 2 : T_from_v_e(v, e, T, dT_dv, dT_de);
168 :
169 2 : c = std::sqrt(_gamma * _R_specific * T);
170 :
171 2 : const Real dc_dT = 0.5 / c * _gamma * _R_specific;
172 2 : dc_dv = dc_dT * dT_dv;
173 2 : dc_de = dc_dT * dT_de;
174 2 : }
175 :
176 : Real
177 120115 : IdealGasFluidProperties::c_from_p_T(Real /*p*/, Real T) const
178 : {
179 120115 : return std::sqrt(_cp * _R * T / (_cv * _molar_mass));
180 : }
181 :
182 : ADReal
183 0 : IdealGasFluidProperties::c_from_p_T(const ADReal & /*p*/, const ADReal & T) const
184 : {
185 0 : return std::sqrt(_cp * _R * T / (_cv * _molar_mass));
186 : }
187 :
188 : void
189 1 : IdealGasFluidProperties::c_from_p_T(
190 : const Real /*p*/, const Real T, Real & c, Real & dc_dp, Real & dc_dT) const
191 : {
192 1 : c = std::sqrt(_cp * _R * T / (_cv * _molar_mass));
193 1 : dc_dp = 0;
194 1 : dc_dT = 0.5 / c * _cp * _R / (_cv * _molar_mass);
195 1 : }
196 :
197 : Real
198 134 : IdealGasFluidProperties::beta_from_p_T(Real /*p*/, Real T) const
199 : {
200 134 : return 1.0 / T;
201 : }
202 :
203 : ADReal
204 0 : IdealGasFluidProperties::beta_from_p_T(const ADReal & /*p*/, const ADReal & T) const
205 : {
206 0 : return 1.0 / T;
207 : }
208 :
209 : void
210 1 : IdealGasFluidProperties::beta_from_p_T(
211 : const Real /*p*/, const Real T, Real & beta, Real & dbeta_dp, Real & dbeta_dT) const
212 : {
213 1 : beta = 1.0 / T;
214 1 : dbeta_dp = 0;
215 1 : dbeta_dT = -1.0 / Utility::pow<2>(T);
216 1 : }
217 :
218 20605 : Real IdealGasFluidProperties::cp_from_v_e(Real, Real) const { return _cp; }
219 :
220 : void
221 1 : IdealGasFluidProperties::cp_from_v_e(Real v, Real e, Real & cp, Real & dcp_dv, Real & dcp_de) const
222 : {
223 1 : cp = cp_from_v_e(v, e);
224 1 : dcp_dv = 0.0;
225 1 : dcp_de = 0.0;
226 1 : }
227 :
228 20605 : Real IdealGasFluidProperties::cv_from_v_e(Real, Real) const { return _cv; }
229 :
230 : void
231 1 : IdealGasFluidProperties::cv_from_v_e(Real v, Real e, Real & cv, Real & dcv_dv, Real & dcv_de) const
232 : {
233 1 : cv = cv_from_v_e(v, e);
234 1 : dcv_dv = 0.0;
235 1 : dcv_de = 0.0;
236 1 : }
237 :
238 1 : Real IdealGasFluidProperties::gamma_from_v_e(Real, Real) const { return _gamma; }
239 :
240 1 : Real IdealGasFluidProperties::gamma_from_p_T(Real, Real) const { return _gamma; }
241 :
242 20686 : Real IdealGasFluidProperties::mu_from_v_e(Real, Real) const { return _mu; }
243 :
244 : void
245 1 : IdealGasFluidProperties::mu_from_v_e(Real v, Real e, Real & mu, Real & dmu_dv, Real & dmu_de) const
246 : {
247 1 : mu = this->mu_from_v_e(v, e);
248 1 : dmu_dv = 0.0;
249 1 : dmu_de = 0.0;
250 1 : }
251 :
252 20685 : Real IdealGasFluidProperties::k_from_v_e(Real, Real) const { return _k; }
253 :
254 : void
255 1 : IdealGasFluidProperties::k_from_v_e(
256 : Real /*v*/, Real /*e*/, Real & k, Real & dk_dv, Real & dk_de) const
257 : {
258 1 : k = _k;
259 1 : dk_dv = 0;
260 1 : dk_de = 0;
261 1 : }
262 :
263 : Real
264 20602 : IdealGasFluidProperties::s_from_v_e(Real v, Real e) const
265 : {
266 20602 : const Real T = T_from_v_e(v, e);
267 20602 : const Real p = p_from_v_e(v, e);
268 20602 : const Real n = std::pow(T, _gamma) / std::pow(p, _gamma - 1.0);
269 20602 : if (n <= 0.0)
270 0 : return getNaN("Negative argument in the ln() function.");
271 20602 : return _cv * std::log(n);
272 : }
273 :
274 : void
275 1 : IdealGasFluidProperties::s_from_v_e(Real v, Real e, Real & s, Real & ds_dv, Real & ds_de) const
276 : {
277 : Real T, dT_dv, dT_de;
278 1 : T_from_v_e(v, e, T, dT_dv, dT_de);
279 :
280 : Real p, dp_dv, dp_de;
281 1 : p_from_v_e(v, e, p, dp_dv, dp_de);
282 :
283 1 : const Real n = std::pow(T, _gamma) / std::pow(p, _gamma - 1.0);
284 1 : if (n <= 0.0)
285 : {
286 0 : s = getNaN("Negative argument in the ln() function.");
287 0 : ds_dv = getNaN();
288 0 : ds_de = getNaN();
289 : }
290 : else
291 : {
292 1 : s = _cv * std::log(n);
293 :
294 1 : const Real dn_dT = _gamma * std::pow(T, _gamma - 1.0) / std::pow(p, _gamma - 1.0);
295 1 : const Real dn_dp = std::pow(T, _gamma) * (1.0 - _gamma) * std::pow(p, -_gamma);
296 :
297 1 : const Real dn_dv = dn_dT * dT_dv + dn_dp * dp_dv;
298 1 : const Real dn_de = dn_dT * dT_de + dn_dp * dp_de;
299 :
300 1 : ds_dv = _cv / n * dn_dv;
301 1 : ds_de = _cv / n * dn_de;
302 : }
303 1 : }
304 :
305 : Real
306 120226 : IdealGasFluidProperties::s_from_p_T(Real p, Real T) const
307 : {
308 120226 : const Real n = std::pow(T, _gamma) / std::pow(p, _gamma - 1.0);
309 120226 : if (n <= 0.0)
310 0 : return getNaN("Negative argument in the ln() function.");
311 120226 : return _cv * std::log(n);
312 : }
313 :
314 : void
315 93 : IdealGasFluidProperties::s_from_p_T(Real p, Real T, Real & s, Real & ds_dp, Real & ds_dT) const
316 : {
317 93 : const Real n = std::pow(T, _gamma) / std::pow(p, _gamma - 1.0);
318 93 : if (n <= 0.0)
319 : {
320 0 : s = getNaN("Negative argument in the ln() function.");
321 0 : ds_dp = getNaN();
322 0 : ds_dT = getNaN();
323 : }
324 : else
325 : {
326 93 : s = _cv * std::log(n);
327 :
328 93 : const Real dn_dT = _gamma * std::pow(T, _gamma - 1.0) / std::pow(p, _gamma - 1.0);
329 93 : const Real dn_dp = std::pow(T, _gamma) * (1.0 - _gamma) * std::pow(p, -_gamma);
330 :
331 93 : ds_dp = _cv / n * dn_dp;
332 93 : ds_dT = _cv / n * dn_dT;
333 : }
334 93 : }
335 :
336 : Real
337 43 : IdealGasFluidProperties::s_from_h_p(Real h, Real p) const
338 : {
339 43 : const Real aux = p * std::pow((h - _e_ref) / (_gamma * _cv), -_gamma / (_gamma - 1));
340 43 : if (aux <= 0.0)
341 0 : return getNaN("Non-positive argument in the ln() function.");
342 43 : return -(_gamma - 1) * _cv * std::log(aux);
343 : }
344 :
345 : void
346 1 : IdealGasFluidProperties::s_from_h_p(Real h, Real p, Real & s, Real & ds_dh, Real & ds_dp) const
347 : {
348 1 : s = s_from_h_p(h, p);
349 :
350 1 : const Real aux = p * std::pow((h - _e_ref) / (_gamma * _cv), -_gamma / (_gamma - 1));
351 1 : const Real daux_dh = p * std::pow((h - _e_ref) / (_gamma * _cv), -_gamma / (_gamma - 1) - 1) *
352 1 : (-_gamma / (_gamma - 1)) / (_gamma * _cv);
353 1 : const Real daux_dp = std::pow((h - _e_ref) / (_gamma * _cv), -_gamma / (_gamma - 1));
354 1 : ds_dh = -(_gamma - 1) * _cv / aux * daux_dh;
355 1 : ds_dp = -(_gamma - 1) * _cv / aux * daux_dp;
356 1 : }
357 :
358 : Real
359 22 : IdealGasFluidProperties::rho_from_p_s(Real p, Real s) const
360 : {
361 22 : const Real aux = (s + _cv * std::log(std::pow(p, _gamma - 1.0))) / _cv;
362 22 : const Real T = std::pow(std::exp(aux), 1.0 / _gamma);
363 22 : return rho_from_p_T(p, T);
364 : }
365 :
366 : void
367 1 : IdealGasFluidProperties::rho_from_p_s(
368 : Real p, Real s, Real & rho, Real & drho_dp, Real & drho_ds) const
369 : {
370 : // T(p,s)
371 1 : const Real aux = (s + _cv * std::log(std::pow(p, _gamma - 1.0))) / _cv;
372 1 : const Real T = std::pow(std::exp(aux), 1 / _gamma);
373 :
374 : // dT/dp
375 1 : const Real dT_dp = 1.0 / _gamma * std::pow(std::exp(aux), 1.0 / _gamma - 1.0) * std::exp(aux) /
376 1 : std::pow(p, _gamma - 1.0) * (_gamma - 1.0) * std::pow(p, _gamma - 2.0);
377 :
378 : // dT/ds
379 : const Real dT_ds =
380 1 : 1.0 / _gamma * std::pow(std::exp(aux), 1.0 / _gamma - 1.0) * std::exp(aux) / _cv;
381 :
382 : // Drho/Dp = d/dp[rho(p, T(p,s))] = drho/dp + drho/dT * dT/dp
383 : Real drho_dp_partial, drho_dT;
384 1 : rho_from_p_T(p, T, rho, drho_dp_partial, drho_dT);
385 1 : drho_dp = drho_dp_partial + drho_dT * dT_dp;
386 :
387 : // Drho/Ds = d/ds[rho(p, T(p,s))] = drho/dT * dT/ds
388 1 : drho_ds = drho_dT * dT_ds;
389 1 : }
390 :
391 : Real
392 9 : IdealGasFluidProperties::e_from_v_h(Real /*v*/, Real h) const
393 : {
394 9 : return (h + (_gamma - 1.0) * _e_ref) / _gamma;
395 : }
396 :
397 : void
398 2 : IdealGasFluidProperties::e_from_v_h(Real v, Real h, Real & e, Real & de_dv, Real & de_dh) const
399 : {
400 2 : e = e_from_v_h(v, h);
401 2 : de_dv = 0.0;
402 2 : de_dh = 1.0 / _gamma;
403 2 : }
404 :
405 : Real
406 1505606 : IdealGasFluidProperties::rho_from_p_T(Real p, Real T) const
407 : {
408 1505606 : return p * _molar_mass / (_R * T);
409 : }
410 :
411 : ADReal
412 0 : IdealGasFluidProperties::rho_from_p_T(const ADReal & p, const ADReal & T) const
413 : {
414 0 : return p * _molar_mass / (_R * T);
415 : }
416 :
417 : void
418 0 : IdealGasFluidProperties::rho_from_p_T(
419 : const ADReal & p, const ADReal & T, ADReal & rho, ADReal & drho_dp, ADReal & drho_dT) const
420 : {
421 0 : rho = rho_from_p_T(p, T);
422 0 : drho_dp = _molar_mass / (_R * T);
423 0 : drho_dT = -p * _molar_mass / (_R * T * T);
424 0 : }
425 :
426 : void
427 1385247 : IdealGasFluidProperties::rho_from_p_T(
428 : Real p, Real T, Real & rho, Real & drho_dp, Real & drho_dT) const
429 : {
430 1385247 : rho = rho_from_p_T(p, T);
431 1385247 : drho_dp = _molar_mass / (_R * T);
432 1385247 : drho_dT = -p * _molar_mass / (_R * T * T);
433 1385247 : }
434 :
435 : Real
436 104 : IdealGasFluidProperties::e_from_p_rho(Real p, Real rho) const
437 : {
438 104 : return p / (_gamma - 1.0) / rho + _e_ref;
439 : }
440 :
441 : ADReal
442 0 : IdealGasFluidProperties::e_from_p_rho(const ADReal & p, const ADReal & rho) const
443 : {
444 0 : return p / (_gamma - 1.0) / rho + _e_ref;
445 : }
446 :
447 : void
448 1 : IdealGasFluidProperties::e_from_p_rho(
449 : Real p, Real rho, Real & e, Real & de_dp, Real & de_drho) const
450 : {
451 1 : e = e_from_p_rho(p, rho);
452 1 : de_dp = 1.0 / (_gamma - 1.0) / rho;
453 1 : de_drho = -p / (_gamma - 1.0) / rho / rho;
454 1 : }
455 :
456 : void
457 0 : IdealGasFluidProperties::e_from_p_rho(
458 : const ADReal & p, const ADReal & rho, ADReal & e, ADReal & de_dp, ADReal & de_drho) const
459 : {
460 0 : e = e_from_p_rho(p, rho);
461 0 : de_dp = 1.0 / (_gamma - 1.0) / rho;
462 0 : de_drho = -p / (_gamma - 1.0) / rho / rho;
463 0 : }
464 :
465 : Real
466 756 : IdealGasFluidProperties::e_from_T_v(Real T, Real /*v*/) const
467 : {
468 756 : return _cv * T + _e_ref;
469 : }
470 :
471 : void
472 7 : IdealGasFluidProperties::e_from_T_v(Real T, Real /*v*/, Real & e, Real & de_dT, Real & de_dv) const
473 : {
474 7 : e = _cv * T + _e_ref;
475 7 : de_dT = _cv;
476 7 : de_dv = 0.0;
477 7 : }
478 :
479 : ADReal
480 0 : IdealGasFluidProperties::e_from_T_v(const ADReal & T, const ADReal & /*v*/) const
481 : {
482 0 : return _cv * T + _e_ref;
483 : }
484 :
485 : void
486 0 : IdealGasFluidProperties::e_from_T_v(
487 : const ADReal & T, const ADReal & /*v*/, ADReal & e, ADReal & de_dT, ADReal & de_dv) const
488 : {
489 0 : e = _cv * T + _e_ref;
490 0 : de_dT = _cv;
491 0 : de_dv = 0.0;
492 0 : }
493 :
494 : Real
495 8812 : IdealGasFluidProperties::p_from_T_v(Real T, Real v) const
496 : {
497 8812 : return (_gamma - 1.0) * _cv * T / v;
498 : }
499 :
500 : void
501 278 : IdealGasFluidProperties::p_from_T_v(Real T, Real v, Real & p, Real & dp_dT, Real & dp_dv) const
502 : {
503 278 : p = (_gamma - 1.0) * _cv * T / v;
504 278 : dp_dT = (_gamma - 1.0) * _cv / v;
505 278 : dp_dv = -(_gamma - 1.0) * _cv * T / (v * v);
506 278 : }
507 :
508 : Real
509 6 : IdealGasFluidProperties::h_from_T_v(Real T, Real /*v*/) const
510 : {
511 6 : return _gamma * _cv * T + _e_ref;
512 : }
513 :
514 : void
515 84 : IdealGasFluidProperties::h_from_T_v(Real T, Real /*v*/, Real & h, Real & dh_dT, Real & dh_dv) const
516 : {
517 84 : h = _gamma * _cv * T + _e_ref;
518 84 : dh_dT = _gamma * _cv;
519 84 : dh_dv = 0.0;
520 84 : }
521 :
522 : Real
523 86 : IdealGasFluidProperties::s_from_T_v(Real T, Real v) const
524 : {
525 86 : Real p = p_from_T_v(T, v);
526 86 : return s_from_p_T(p, T);
527 : }
528 :
529 : void
530 92 : IdealGasFluidProperties::s_from_T_v(Real T, Real v, Real & s, Real & ds_dT, Real & ds_dv) const
531 : {
532 : Real p, dp_dT_v, dp_dv_T;
533 : Real ds_dp_T, ds_dT_p;
534 92 : p_from_T_v(T, v, p, dp_dT_v, dp_dv_T);
535 92 : s_from_p_T(p, T, s, ds_dp_T, ds_dT_p);
536 92 : ds_dT = ds_dT_p + ds_dp_T * dp_dT_v;
537 92 : ds_dv = ds_dp_T * dp_dv_T;
538 92 : }
539 :
540 84 : Real IdealGasFluidProperties::cv_from_T_v(Real /*T*/, Real /*v*/) const { return _cv; }
541 :
542 32 : Real IdealGasFluidProperties::e_spndl_from_v(Real /*v*/) const { return _e_c; }
543 :
544 : void
545 0 : IdealGasFluidProperties::v_e_spndl_from_T(Real /*T*/, Real & v, Real & e) const
546 : {
547 0 : v = 1. / _rho_c;
548 0 : e = _e_c;
549 0 : }
550 :
551 : Real
552 812507 : IdealGasFluidProperties::h_from_p_T(Real /*p*/, Real T) const
553 : {
554 812507 : return _cp * T + _e_ref;
555 : }
556 :
557 : void
558 692125 : IdealGasFluidProperties::h_from_p_T(Real p, Real T, Real & h, Real & dh_dp, Real & dh_dT) const
559 : {
560 692125 : h = h_from_p_T(p, T);
561 692125 : dh_dp = 0.0;
562 692125 : dh_dT = _cp;
563 692125 : }
564 :
565 : Real
566 813269 : IdealGasFluidProperties::e_from_p_T(Real /*p*/, Real T) const
567 : {
568 813269 : return _cv * T + _e_ref;
569 : }
570 :
571 : void
572 693121 : IdealGasFluidProperties::e_from_p_T(Real p, Real T, Real & e, Real & de_dp, Real & de_dT) const
573 : {
574 693121 : e = e_from_p_T(p, T);
575 693121 : de_dp = 0.0;
576 693121 : de_dT = _cv;
577 693121 : }
578 :
579 : Real
580 23 : IdealGasFluidProperties::p_from_h_s(Real h, Real s) const
581 : {
582 23 : return std::pow((h - _e_ref) / (_gamma * _cv), _gamma / (_gamma - 1.0)) *
583 23 : std::exp(-s / ((_gamma - 1.0) * _cv));
584 : }
585 :
586 : void
587 1 : IdealGasFluidProperties::p_from_h_s(Real h, Real s, Real & p, Real & dp_dh, Real & dp_ds) const
588 : {
589 1 : p = p_from_h_s(h, s);
590 1 : dp_dh = _gamma / (_gamma - 1.0) / (_gamma * _cv) *
591 1 : std::pow((h - _e_ref) / (_gamma * _cv), 1.0 / (_gamma - 1.0)) *
592 1 : std::exp(-s / ((_gamma - 1.0) * _cv));
593 1 : dp_ds = std::pow((h - _e_ref) / (_gamma * _cv), _gamma / (_gamma - 1)) *
594 1 : std::exp(-s / ((_gamma - 1) * _cv)) / ((1 - _gamma) * _cv);
595 1 : }
596 :
597 : Real
598 469 : IdealGasFluidProperties::g_from_v_e(Real v, Real e) const
599 : {
600 : // g(p,T) for SGEOS is given by Equation (37) in the following reference:
601 : //
602 : // Ray A. Berry, Richard Saurel, Olivier LeMetayer
603 : // The discrete equation method (DEM) for fully compressible, two-phase flows in
604 : // ducts of spatially varying cross-section
605 : // Nuclear Engineering and Design 240 (2010) p. 3797-3818
606 : //
607 469 : const Real p = p_from_v_e(v, e);
608 469 : const Real T = T_from_v_e(v, e);
609 :
610 469 : return _gamma * _cv * T - _cv * T * std::log(std::pow(T, _gamma) / std::pow(p, _gamma - 1.0));
611 : }
612 :
613 : Real
614 1002 : IdealGasFluidProperties::molarMass() const
615 : {
616 1002 : return _molar_mass;
617 : }
618 :
619 : Real
620 32 : IdealGasFluidProperties::criticalTemperature() const
621 : {
622 32 : return _T_c;
623 : }
624 :
625 : Real
626 48 : IdealGasFluidProperties::criticalDensity() const
627 : {
628 48 : return _rho_c;
629 : }
630 :
631 : Real
632 16 : IdealGasFluidProperties::criticalInternalEnergy() const
633 : {
634 16 : return _e_c;
635 : }
636 :
637 : Real
638 42 : IdealGasFluidProperties::T_from_p_h(Real, Real h) const
639 : {
640 42 : return (h - _e_ref) / _gamma / _cv;
641 : }
642 :
643 : void
644 1 : IdealGasFluidProperties::T_from_p_h(Real /*p*/, Real h, Real & T, Real & dT_dp, Real & dT_dh) const
645 : {
646 1 : T = (h - _e_ref) / (_gamma * _cv);
647 1 : dT_dp = 0;
648 1 : dT_dh = 1.0 / (_gamma * _cv);
649 1 : }
650 :
651 120160 : Real IdealGasFluidProperties::cv_from_p_T(Real /* pressure */, Real /* temperature */) const
652 : {
653 120160 : return _cv;
654 : }
655 :
656 : void
657 1 : IdealGasFluidProperties::cv_from_p_T(Real p, Real T, Real & cv, Real & dcv_dp, Real & dcv_dT) const
658 : {
659 1 : cv = cv_from_p_T(p, T);
660 1 : dcv_dp = 0.0;
661 1 : dcv_dT = 0.0;
662 1 : }
663 :
664 120154 : Real IdealGasFluidProperties::cp_from_p_T(Real /* pressure */, Real /* temperature */) const
665 : {
666 120154 : return _cp;
667 : }
668 :
669 : void
670 1 : IdealGasFluidProperties::cp_from_p_T(Real p, Real T, Real & cp, Real & dcp_dp, Real & dcp_dT) const
671 : {
672 1 : cp = cp_from_p_T(p, T);
673 1 : dcp_dp = 0.0;
674 1 : dcp_dT = 0.0;
675 1 : }
676 :
677 120133 : Real IdealGasFluidProperties::mu_from_p_T(Real /* pressure */, Real /* temperature */) const
678 : {
679 120133 : return _mu;
680 : }
681 :
682 : void
683 1 : IdealGasFluidProperties::mu_from_p_T(Real p, Real T, Real & mu, Real & dmu_dp, Real & dmu_dT) const
684 : {
685 1 : mu = this->mu_from_p_T(p, T);
686 1 : dmu_dp = 0.0;
687 1 : dmu_dT = 0.0;
688 1 : }
689 :
690 120133 : Real IdealGasFluidProperties::k_from_p_T(Real /* pressure */, Real /* temperature */) const
691 : {
692 120133 : return _k;
693 : }
694 :
695 : void
696 1 : IdealGasFluidProperties::k_from_p_T(Real p, Real T, Real & k, Real & dk_dp, Real & dk_dT) const
697 : {
698 1 : k = k_from_p_T(p, T);
699 1 : dk_dp = 0.0;
700 1 : dk_dT = 0.0;
701 1 : }
702 :
703 0 : Real IdealGasFluidProperties::pp_sat_from_p_T(Real /*p*/, Real /*T*/) const
704 : {
705 0 : mooseError(__PRETTY_FUNCTION__, " not implemented. Use a real fluid property class!");
706 : }
|