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 : #pragma once
11 :
12 : #include "FluidProperties.h"
13 : #include "NewtonInversion.h"
14 : #include "metaphysicl/dualnumberarray.h"
15 :
16 : /**
17 : * Adds AD versions of each fluid property. These functions use the Real versions of these methods
18 : * to compute the AD variables complete with derivatives. Typically, these do not need to be
19 : * overriden in derived classes.
20 : */
21 : #define propfuncAD(want, prop1, prop2) \
22 : virtual ADReal want##_from_##prop1##_##prop2(const ADReal & p1, const ADReal & p2) const \
23 : { \
24 : Real x = 0; \
25 : Real raw1 = p1.value(); \
26 : Real raw2 = p2.value(); \
27 : Real dxd1 = 0; \
28 : Real dxd2 = 0; \
29 : want##_from_##prop1##_##prop2(raw1, raw2, x, dxd1, dxd2); \
30 : \
31 : ADReal result = x; \
32 : result.derivatives() = p1.derivatives() * dxd1 + p2.derivatives() * dxd2; \
33 : return result; \
34 : } \
35 : \
36 : virtual void want##_from_##prop1##_##prop2(const ADReal & prop1, \
37 : const ADReal & prop2, \
38 : ADReal & val, \
39 : ADReal & d##want##d1, \
40 : ADReal & d##want##d2) const \
41 : { \
42 : unimplementedDerivativeMethod(__PRETTY_FUNCTION__); \
43 : Real dummy, tmp1, tmp2; \
44 : val = want##_from_##prop1##_##prop2(prop1, prop2); \
45 : want##_from_##prop1##_##prop2(prop1.value(), prop2.value(), dummy, tmp1, tmp2); \
46 : d##want##d1 = tmp1; \
47 : d##want##d2 = tmp2; \
48 : }
49 :
50 : /**
51 : * Adds function definitions with not implemented error. These functions should be overriden in
52 : * derived classes where required. AD versions are constructed automatically using propfuncAD.
53 : */
54 : #define propfunc(want, prop1, prop2) \
55 : virtual Real want##_from_##prop1##_##prop2(Real, Real) const \
56 : { \
57 : mooseError( \
58 : "The fluid properties class '", \
59 : type(), \
60 : "' has not implemented the method below. If your application requires this method, you " \
61 : "must either implement it or use a different fluid properties class.\n\n", \
62 : __PRETTY_FUNCTION__); \
63 : } \
64 : \
65 : virtual void want##_from_##prop1##_##prop2( \
66 : Real prop1, Real prop2, Real & val, Real & d##want##d1, Real & d##want##d2) const \
67 : { \
68 : unimplementedDerivativeMethod(__PRETTY_FUNCTION__); \
69 : d##want##d1 = 0; \
70 : d##want##d2 = 0; \
71 : val = want##_from_##prop1##_##prop2(prop1, prop2); \
72 : } \
73 : \
74 : propfuncAD(want, prop1, prop2)
75 :
76 : /**
77 : * Adds Real declarations of functions that have a default implementation.
78 : * Important: properties declared using this macro must be defined in SinglePhaseFluidProperties.C.
79 : * AD versions are constructed automatically using propfuncAD.
80 : */
81 : #define propfuncWithDefault(want, prop1, prop2) \
82 : virtual Real want##_from_##prop1##_##prop2(Real, Real) const; \
83 : virtual void want##_from_##prop1##_##prop2( \
84 : Real prop1, Real prop2, Real & val, Real & d##want##d1, Real & d##want##d2) const; \
85 : \
86 : propfuncAD(want, prop1, prop2)
87 :
88 : /**
89 : * Adds Real and ADReal declarations of functions that have an implementation.
90 : */
91 : #define propfuncWithDefinitionOverride(want, prop1, prop2) \
92 : Real want##_from_##prop1##_##prop2(Real, Real) const override; \
93 : void want##_from_##prop1##_##prop2( \
94 : Real prop1, Real prop2, Real & val, Real & d##want##d1, Real & d##want##d2) const override; \
95 : ADReal want##_from_##prop1##_##prop2(const ADReal &, const ADReal &) const override; \
96 : void want##_from_##prop1##_##prop2(const ADReal & prop1, \
97 : const ADReal & prop2, \
98 : ADReal & val, \
99 : ADReal & d##want##d1, \
100 : ADReal & d##want##d2) const override; \
101 : template <typename CppType> \
102 : CppType want##_from_##prop1##_##prop2##_template(const CppType & prop1, const CppType & prop2) \
103 : const; \
104 : template <typename CppType> \
105 : void want##_from_##prop1##_##prop2##_template(const CppType & prop1, \
106 : const CppType & prop2, \
107 : CppType & val, \
108 : CppType & d##want##d1, \
109 : CppType & d##want##d2) const
110 :
111 : /**
112 : * Common class for single phase fluid properties
113 : */
114 : class SinglePhaseFluidProperties : public FluidProperties
115 : {
116 : public:
117 : static InputParameters validParams();
118 :
119 : SinglePhaseFluidProperties(const InputParameters & parameters);
120 : virtual ~SinglePhaseFluidProperties();
121 :
122 : #pragma GCC diagnostic push
123 : #pragma GCC diagnostic ignored "-Woverloaded-virtual"
124 : // clang-format off
125 :
126 : /**
127 : * @brief Compute a fluid property given for the state defined by two given properties.
128 : *
129 : * For all functions, the first two arguments are the given properties that define the fluid
130 : * state. For the two-argument variants, the desired property is the return value.
131 : * The five-argument variants also provide partial derivatives dx/da and dx/db where x is the
132 : * desired property being computed, a is the first given property, and b is the second given
133 : * property. The desired property, dx/da, and dx/db are stored into the 3rd, 4th, and 5th
134 : * arguments respectively.
135 : *
136 : * Properties/parameters used in these function are listed below with their units:
137 : *
138 : * @begincode
139 : * p pressure [Pa]
140 : * T temperature [K]
141 : * e specific internal energy [J/kg]
142 : * v specific volume [m^3/kg]
143 : * rho density [kg/m^3]
144 : * h specific enthalpy [J/kg]
145 : * s specific entropy [J/(kg*K)]
146 : * mu viscosity [Pa*s]
147 : * k thermal conductivity [W/(m*K)]
148 : * c speed of sound [m/s]
149 : * cp constant-pressure specific heat [J/K]
150 : * cv constant-volume specific heat [J/K]
151 : * beta volumetric thermal expansion coefficient [1/K]
152 : * g Gibbs free energy [J]
153 : * pp_sat partial pressure at saturation [Pa]
154 : * gamma Adiabatic ratio (cp/cv) [-]
155 : * @endcode
156 : *
157 : * As an example:
158 : *
159 : * @begincode
160 : * // calculate pressure given specific vol and energy:
161 : * auto pressure = your_fluid_properties_object.p_from_v_e(specific_vol, specific_energy);
162 : *
163 : * // or use the derivative variant:
164 : * Real dp_dv = 0; // derivative will be stored into here
165 : * Real dp_de = 0; // derivative will be stored into here
166 : * your_fluid_properties_object.p_from_v_e(specific_vol, specific_energy, pressure, dp_dv, dp_de);
167 : * @endcode
168 : *
169 : * Automatic differentiation (AD) support is provided through x_from_a_b(ADReal a, ADReal b) and
170 : * x_from_a_b(ADReal a, ADReal b, ADReal x, ADReal dx_da, ADReal dx_db) versions of the
171 : * functions where a and b must be ADReal/DualNumber's calculated using all AD-supporting values:
172 : *
173 : * @begincode
174 : * auto v = 1/rho; // rho must be an AD non-linear variable.
175 : * auto e = rhoE/rho - vel_energy; // rhoE and vel_energy must be AD variables/numbers also.
176 : * auto pressure = your_fluid_properties_object.p_from_v_e(v, e);
177 : * // pressure now contains partial derivatives w.r.t. all degrees of freedom
178 : * @endcode
179 : */
180 : ///@{
181 2 : propfunc(p, v, e)
182 2 : propfunc(T, v, e)
183 0 : propfunc(c, v, e)
184 0 : propfunc(cp, v, e)
185 0 : propfunc(cv, v, e)
186 0 : propfunc(mu, v, e)
187 0 : propfunc(k, v, e)
188 0 : propfuncWithDefault(s, v, e)
189 0 : propfunc(s, h, p)
190 0 : propfunc(rho, p, s)
191 0 : propfunc(e, v, h)
192 0 : propfuncWithDefault(s, p, T)
193 0 : propfunc(pp_sat, p, T)
194 0 : propfunc(mu, rho, T)
195 0 : propfunc(k, rho, T)
196 0 : propfuncWithDefault(c, p, T)
197 6 : propfuncWithDefault(cp, p, T)
198 0 : propfuncWithDefault(cv, p, T)
199 0 : propfuncWithDefault(mu, p, T)
200 0 : propfuncWithDefault(k, p, T)
201 20 : propfunc(rho, p, T)
202 0 : propfunc(e, p, rho)
203 0 : propfunc(e, T, v)
204 0 : propfunc(p, T, v)
205 0 : propfunc(h, T, v)
206 0 : propfunc(s, T, v)
207 0 : propfunc(cv, T, v)
208 2 : propfunc(h, p, T)
209 0 : propfuncWithDefault(h, v, e)
210 7 : propfunc(g, v, e)
211 0 : propfuncWithDefault(p, h, s)
212 0 : propfunc(T, h, p) // temporary, until uniformization
213 0 : propfuncWithDefault(T, p, h)
214 0 : propfuncWithDefault(beta, p, T)
215 3 : propfuncWithDefault(v, p, T)
216 6 : propfuncWithDefault(e, p, T)
217 0 : propfuncWithDefault(gamma, v, e)
218 0 : propfuncWithDefault(gamma, p, T)
219 : ///@}
220 :
221 : // clang-format on
222 :
223 : #undef propfunc
224 : #undef propfuncWithDefault
225 : #undef propfuncAD
226 :
227 : /**
228 : * Fluid name
229 : * @return string representing fluid name
230 : */
231 : virtual std::string fluidName() const;
232 :
233 : /**
234 : * Molar mass [kg/mol]
235 : * @return molar mass
236 : */
237 : virtual Real molarMass() const;
238 :
239 : /**
240 : * Critical pressure
241 : * @return critical pressure (Pa)
242 : */
243 : virtual Real criticalPressure() const;
244 :
245 : /**
246 : * Critical temperature
247 : * @return critical temperature (K)
248 : */
249 : virtual Real criticalTemperature() const;
250 :
251 : /**
252 : * Critical density
253 : * @return critical density (kg/m^3)
254 : */
255 : virtual Real criticalDensity() const;
256 :
257 : /**
258 : * Critical specific internal energy
259 : * @return specific internal energy (J/kg)
260 : */
261 : virtual Real criticalInternalEnergy() const;
262 :
263 : /**
264 : * Triple point pressure
265 : * @return triple point pressure (Pa)
266 : */
267 : virtual Real triplePointPressure() const;
268 :
269 : /**
270 : * Triple point temperature
271 : * @return triple point temperature (K)
272 : */
273 : virtual Real triplePointTemperature() const;
274 :
275 : /**
276 : * Specific internal energy from temperature and specific volume
277 : *
278 : * @param[in] T temperature
279 : * @param[in] v specific volume
280 : */
281 : virtual Real e_spndl_from_v(Real v) const;
282 :
283 : /**
284 : * Specific internal energy from temperature and specific volume
285 : *
286 : * @param[in] T temperature
287 : * @param[in] v specific volume
288 : */
289 : virtual void v_e_spndl_from_T(Real T, Real & v, Real & e) const;
290 :
291 : /**
292 : * Vapor pressure. Used to delineate liquid and gas phases.
293 : * Valid for temperatures between the triple point temperature
294 : * and the critical temperature
295 : *
296 : * @param T fluid temperature (K)
297 : * @param[out] saturation pressure (Pa)
298 : * @param[out] derivative of saturation pressure wrt temperature (Pa/K)
299 : */
300 : virtual Real vaporPressure(Real T) const;
301 : virtual void vaporPressure(Real T, Real & psat, Real & dpsat_dT) const;
302 : virtual ADReal vaporPressure(const ADReal & T) const;
303 :
304 : /**
305 : * Vapor temperature. Used to delineate liquid and gas phases.
306 : * Valid for pressures between the triple point pressure
307 : * and the critical pressure
308 : *
309 : * @param p fluid pressure (Pa)
310 : * @param[out] saturation temperature (K)
311 : * @param[out] derivative of saturation temperature wrt pressure
312 : */
313 : virtual Real vaporTemperature(Real p) const;
314 : virtual void vaporTemperature(Real p, Real & Tsat, Real & dTsat_dp) const;
315 : virtual ADReal vaporTemperature(const ADReal & p) const;
316 :
317 : /**
318 : * Henry's law coefficients for dissolution in water
319 : * @return Henry's constant coefficients
320 : */
321 : virtual std::vector<Real> henryCoefficients() const;
322 :
323 : template <typename CppType>
324 : void v_e_from_p_T(const CppType & p, const CppType & T, CppType & v, CppType & e) const;
325 : template <typename CppType>
326 : void v_e_from_p_T(const CppType & p,
327 : const CppType & T,
328 : CppType & v,
329 : CppType & dv_dp,
330 : CppType & dv_dT,
331 : CppType & e,
332 : CppType & de_dp,
333 : CppType & de_dT) const;
334 :
335 : /**
336 : * Combined methods. These methods are particularly useful for the PorousFlow
337 : * module, where density and viscosity are typically both computed everywhere.
338 : * The combined methods allow the most efficient means of calculating both
339 : * properties, especially where rho(p, T) and mu(rho, T). In this case, an
340 : * extra density calculation would be required to calculate mu(p, T). All
341 : * property names are described above.
342 : */
343 : virtual void rho_mu_from_p_T(Real p, Real T, Real & rho, Real & mu) const;
344 : virtual void rho_mu_from_p_T(Real p,
345 : Real T,
346 : Real & rho,
347 : Real & drho_dp,
348 : Real & drho_dT,
349 : Real & mu,
350 : Real & dmu_dp,
351 : Real & dmu_dT) const;
352 : virtual void rho_mu_from_p_T(const ADReal & p, const ADReal & T, ADReal & rho, ADReal & mu) const;
353 :
354 : virtual void rho_e_from_p_T(Real p,
355 : Real T,
356 : Real & rho,
357 : Real & drho_dp,
358 : Real & drho_dT,
359 : Real & e,
360 : Real & de_dp,
361 : Real & de_dT) const;
362 :
363 : /**
364 : * Determines (p,T) from (v,e) using Newton Solve in 2D
365 : * Useful for conversion between different sets of state variables
366 : *
367 : * @param[in] v specific volume (m^3 / kg)
368 : * @param[in] e specific internal energy (J / kg)
369 : * @param[in] p0 initial guess for pressure (Pa / kg)
370 : * @param[in] T0 initial guess for temperature (K)
371 : * @param[out] fluid pressure (Pa / kg)
372 : * @param[out] Temperature (K)
373 : */
374 : template <typename CppType>
375 : void p_T_from_v_e(const CppType & v,
376 : const CppType & e,
377 : Real p0,
378 : Real T0,
379 : CppType & p,
380 : CppType & T,
381 : bool & conversion_succeeded) const;
382 :
383 : /**
384 : * Determines (p,T) from (v,h) using Newton Solve in 2D
385 : * Useful for conversion between different sets of state variables
386 : *
387 : * @param[in] v specific volume (m^3 / kg)
388 : * @param[in] h specific enthalpy (J / kg)
389 : * @param[in] p0 initial guess for pressure (Pa / kg)
390 : * @param[in] T0 initial guess for temperature (K)
391 : * @param[out] fluid pressure (Pa / kg)
392 : * @param[out] Temperature (K)
393 : */
394 : template <typename T>
395 : void p_T_from_v_h(const T & v,
396 : const T & h,
397 : Real p0,
398 : Real T0,
399 : T & pressure,
400 : T & temperature,
401 : bool & conversion_succeeded) const;
402 : /**
403 : * Determines (p,T) from (h,s) using Newton Solve in 2D
404 : * Useful for conversion between different sets of state variables
405 : *
406 : * @param[in] h specific enthalpy (J / kg)
407 : * @param[in] s specific entropy (J/K*kg)
408 : * @param[in] p0 initial guess for pressure (Pa / kg)
409 : * @param[in] T0 initial guess for temperature (K)
410 : * @param[out] fluid pressure (Pa / kg)
411 : * @param[out] Temperature (K)
412 : */
413 : template <typename T>
414 : void p_T_from_h_s(const T & h,
415 : const T & s,
416 : Real p0,
417 : Real T0,
418 : T & pressure,
419 : T & temperature,
420 : bool & conversion_succeeded) const;
421 :
422 : protected:
423 : /**
424 : * Computes the dependent variable z and its derivatives with respect to the independent
425 : * variables x and y using the simple two parameter \p z_from_x_y functor. The derivatives are
426 : * computed using a compound automatic differentiation type
427 : */
428 : template <typename T, typename Functor>
429 : static void
430 : xyDerivatives(const T x, const T & y, T & z, T & dz_dx, T & dz_dy, const Functor & z_from_x_y);
431 :
432 : /**
433 : * Given a type example, this method returns zero and unity representations of that type (first
434 : * and second members of returned pair respectively)
435 : */
436 : template <typename T>
437 : static std::pair<T, T> makeZeroAndOne(const T &);
438 :
439 : /**
440 : * Newton's method may be used to convert between variable sets
441 : */
442 : /// Relative tolerance of the solves
443 : const Real _tolerance;
444 : /// Initial guess for temperature (or temperature used to compute the initial guess)
445 : const Real _T_initial_guess;
446 : /// Initial guess for pressure (or pressure used to compute the initial guess)
447 : const Real _p_initial_guess;
448 : /// Maximum number of iterations for the variable conversion newton solves
449 : const unsigned int _max_newton_its;
450 : /// Whether to output information about newton solves to console
451 : const bool _verbose_newton;
452 :
453 : private:
454 3 : void unimplementedDerivativeMethod(const std::string & property_function_name) const
455 : {
456 : const std::string message =
457 3 : "The fluid properties class '" + type() +
458 : "' has not implemented the method below, which computes derivatives of fluid properties "
459 : "with regards to the flow variables. If your application requires this "
460 : "method, you must either implement it or use a different fluid properties "
461 : " class.\n\n" +
462 3 : property_function_name;
463 :
464 3 : if (_allow_imperfect_jacobians)
465 2 : mooseDoOnce(mooseWarning(message + "\nThe unimplemented derivatives for this fluid property "
466 : "are currently neglected, set to 0."));
467 : else
468 3 : mooseError(message + "\n\nYou can avoid this error by neglecting the "
469 : "unimplemented derivatives of fluid properties by setting the "
470 : "'allow_imperfect_jacobians' parameter");
471 1 : }
472 : };
473 :
474 : #pragma GCC diagnostic pop
475 :
476 : template <typename T>
477 : std::pair<T, T>
478 317 : SinglePhaseFluidProperties::makeZeroAndOne(const T & /*ex*/)
479 : {
480 317 : return {T{0, 0}, T{1, 0}};
481 : }
482 :
483 : template <>
484 : inline std::pair<Real, Real>
485 : SinglePhaseFluidProperties::makeZeroAndOne(const Real & /*ex*/)
486 : {
487 : return {Real{0}, Real{1}};
488 : }
489 :
490 : template <typename T, typename Functor>
491 : void
492 418 : SinglePhaseFluidProperties::xyDerivatives(
493 : const T x, const T & y, T & z, T & dz_dx, T & dz_dy, const Functor & z_from_x_y)
494 : {
495 : typedef MetaPhysicL::DualNumber<T, MetaPhysicL::NumberArray<2, T>> CompoundType;
496 317 : const auto [zero, one] = makeZeroAndOne(x);
497 :
498 317 : CompoundType x_c(x, zero);
499 : auto & x_cd = x_c.derivatives();
500 358 : x_cd[0] = one;
501 317 : CompoundType y_c(y, zero);
502 : auto & y_cd = y_c.derivatives();
503 358 : y_cd[1] = one;
504 :
505 : const auto z_c = z_from_x_y(x_c, y_c);
506 358 : z = z_c.value();
507 358 : dz_dx = z_c.derivatives()[0];
508 358 : dz_dy = z_c.derivatives()[1];
509 418 : }
510 :
511 : template <typename CppType>
512 : void
513 260059 : SinglePhaseFluidProperties::p_T_from_v_e(const CppType & v, // v value
514 : const CppType & e, // e value
515 : const Real p0, // initial guess
516 : const Real T0, // initial guess
517 : CppType & p, // returned pressure
518 : CppType & T, // returned temperature
519 : bool & conversion_succeeded) const
520 : {
521 260059 : auto v_lambda = [&](const CppType & pressure,
522 : const CppType & temperature,
523 : CppType & new_v,
524 : CppType & dv_dp,
525 6738665 : CppType & dv_dT) { v_from_p_T(pressure, temperature, new_v, dv_dp, dv_dT); };
526 260059 : auto e_lambda = [&](const CppType & pressure,
527 : const CppType & temperature,
528 : CppType & new_e,
529 : CppType & de_dp,
530 6738665 : CppType & de_dT) { e_from_p_T(pressure, temperature, new_e, de_dp, de_dT); };
531 : try
532 : {
533 780177 : FluidPropertiesUtils::NewtonSolve2D(v,
534 : e,
535 : p0,
536 : T0,
537 : p,
538 : T,
539 : _tolerance,
540 260059 : _tolerance,
541 : v_lambda,
542 : e_lambda,
543 : "p_T_from_v_e",
544 260059 : _max_newton_its,
545 260059 : _verbose_newton);
546 200945 : conversion_succeeded = true;
547 : }
548 59114 : catch (MooseException &)
549 : {
550 59114 : conversion_succeeded = false;
551 : }
552 :
553 260059 : if (!conversion_succeeded)
554 59114 : mooseDoOnce(mooseWarning("Conversion from (v, e)=(", v, ", ", e, ") to (p, T) failed"));
555 260059 : }
556 :
557 : template <typename T>
558 : void
559 260004 : SinglePhaseFluidProperties::p_T_from_v_h(const T & v, // v value
560 : const T & h, // e value
561 : const Real p0, // initial guess
562 : const Real T0, // initial guess
563 : T & pressure, // returned pressure
564 : T & temperature, // returned temperature
565 : bool & conversion_succeeded) const
566 : {
567 260004 : auto v_lambda = [&](const T & pressure, const T & temperature, T & new_v, T & dv_dp, T & dv_dT)
568 6877620 : { v_from_p_T(pressure, temperature, new_v, dv_dp, dv_dT); };
569 260004 : auto h_lambda = [&](const T & pressure, const T & temperature, T & new_h, T & dh_dp, T & dh_dT)
570 6877620 : { h_from_p_T(pressure, temperature, new_h, dh_dp, dh_dT); };
571 : try
572 : {
573 780012 : FluidPropertiesUtils::NewtonSolve2D(v,
574 : h,
575 : p0,
576 : T0,
577 : pressure,
578 : temperature,
579 : _tolerance,
580 260004 : _tolerance,
581 : v_lambda,
582 : h_lambda,
583 : "p_T_from_v_h",
584 260004 : _max_newton_its,
585 260004 : _verbose_newton);
586 199444 : conversion_succeeded = true;
587 : }
588 60560 : catch (MooseException &)
589 : {
590 60560 : conversion_succeeded = false;
591 : }
592 :
593 260004 : if (!conversion_succeeded)
594 60560 : mooseDoOnce(mooseWarning("Conversion from (v, h)=(", v, ", ", h, ") to (p, T) failed"));
595 260004 : }
596 :
597 : template <typename T>
598 : void
599 2 : SinglePhaseFluidProperties::p_T_from_h_s(const T & h, // h value
600 : const T & s, // s value
601 : const Real p0, // initial guess
602 : const Real T0, // initial guess
603 : T & pressure, // returned pressure
604 : T & temperature, // returned temperature
605 : bool & conversion_succeeded) const
606 : {
607 2 : auto h_lambda = [&](const T & pressure, const T & temperature, T & new_h, T & dh_dp, T & dh_dT)
608 6 : { h_from_p_T(pressure, temperature, new_h, dh_dp, dh_dT); };
609 2 : auto s_lambda = [&](const T & pressure, const T & temperature, T & new_s, T & ds_dp, T & ds_dT)
610 6 : { s_from_p_T(pressure, temperature, new_s, ds_dp, ds_dT); };
611 : try
612 : {
613 6 : FluidPropertiesUtils::NewtonSolve2D(h,
614 : s,
615 : p0,
616 : T0,
617 : pressure,
618 : temperature,
619 : _tolerance,
620 2 : _tolerance,
621 : h_lambda,
622 : s_lambda,
623 : "p_T_from_h_s",
624 2 : _max_newton_its,
625 2 : _verbose_newton);
626 2 : conversion_succeeded = true;
627 : }
628 0 : catch (MooseException &)
629 : {
630 0 : conversion_succeeded = false;
631 : }
632 :
633 2 : if (!conversion_succeeded)
634 0 : mooseDoOnce(mooseWarning("Conversion from (h, s)=(", h, ", ", s, ") to (p, T) failed"));
635 2 : }
636 :
637 : template <typename CppType>
638 : void
639 543 : SinglePhaseFluidProperties::v_e_from_p_T(const CppType & p,
640 : const CppType & T,
641 : CppType & v,
642 : CppType & e) const
643 : {
644 543 : const CppType rho = rho_from_p_T(p, T);
645 543 : v = 1.0 / rho;
646 : try
647 : {
648 : // more likely to not involve a Newton search
649 543 : e = e_from_p_T(p, T);
650 : }
651 0 : catch (...)
652 : {
653 0 : e = e_from_p_rho(p, rho);
654 : }
655 543 : }
656 :
657 : template <typename CppType>
658 : void
659 31 : SinglePhaseFluidProperties::v_e_from_p_T(const CppType & p,
660 : const CppType & T,
661 : CppType & v,
662 : CppType & dv_dp,
663 : CppType & dv_dT,
664 : CppType & e,
665 : CppType & de_dp,
666 : CppType & de_dT) const
667 : {
668 : CppType rho, drho_dp, drho_dT;
669 31 : rho_from_p_T(p, T, rho, drho_dp, drho_dT);
670 :
671 31 : v = 1.0 / rho;
672 31 : const CppType dv_drho = -1.0 / (rho * rho);
673 31 : dv_dp = dv_drho * drho_dp;
674 31 : dv_dT = dv_drho * drho_dT;
675 :
676 : CppType de_dp_partial, de_drho;
677 31 : e_from_p_rho(p, rho, e, de_dp_partial, de_drho);
678 31 : de_dp = de_dp_partial + de_drho * drho_dp;
679 31 : de_dT = de_drho * drho_dT;
680 31 : }
|