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 "SinglePhaseFluidProperties.h"
11 :
12 : InputParameters
13 2787 : SinglePhaseFluidProperties::validParams()
14 : {
15 2787 : InputParameters params = FluidProperties::validParams();
16 2787 : params.set<std::string>("fp_type") = "single-phase-fp";
17 :
18 : // Variable set conversion parameters
19 8361 : params.addRangeCheckedParam<Real>(
20 5574 : "tolerance", 1e-8, "tolerance > 0", "Tolerance for 2D Newton variable set conversion");
21 5574 : params.addRangeCheckedParam<Real>(
22 : "T_initial_guess",
23 : 400,
24 : "T_initial_guess > 0",
25 : "Temperature initial guess for Newton Method variable set conversion");
26 8361 : params.addRangeCheckedParam<Real>(
27 : "p_initial_guess",
28 5574 : 2e5,
29 : "p_initial_guess > 0",
30 : "Pressure initial guess for Newton Method variable set conversion");
31 5574 : params.addParam<unsigned int>(
32 5574 : "max_newton_its", 100, "Maximum number of Newton iterations for variable set conversions");
33 5574 : params.addParam<bool>(
34 5574 : "verbose_newton", false, "Whether to output Newton inversion iterations to console");
35 5574 : params.addParamNamesToGroup(
36 : "tolerance T_initial_guess p_initial_guess max_newton_its verbose_newton",
37 : "Variable set conversions Newton solve");
38 :
39 2787 : return params;
40 0 : }
41 :
42 1467 : SinglePhaseFluidProperties::SinglePhaseFluidProperties(const InputParameters & parameters)
43 : : FluidProperties(parameters),
44 : // downstream apps are creating fluid properties without their parameters, hence the workaround
45 1467 : _tolerance(getParam<Real>("tolerance")),
46 2934 : _T_initial_guess(getParam<Real>("T_initial_guess")),
47 2934 : _p_initial_guess(getParam<Real>("p_initial_guess")),
48 2934 : _max_newton_its(getParam<unsigned int>("max_newton_its")),
49 4401 : _verbose_newton(getParam<bool>("verbose_newton"))
50 : {
51 1467 : }
52 :
53 1377 : SinglePhaseFluidProperties::~SinglePhaseFluidProperties() {}
54 :
55 : #pragma GCC diagnostic push
56 : #pragma GCC diagnostic ignored "-Woverloaded-virtual"
57 :
58 : Real
59 126 : SinglePhaseFluidProperties::s_from_p_T(const Real pressure, const Real temperature) const
60 : {
61 : Real v, e;
62 126 : v_e_from_p_T(pressure, temperature, v, e);
63 126 : return s_from_v_e(v, e);
64 : }
65 :
66 : void
67 1 : SinglePhaseFluidProperties::s_from_p_T(
68 : const Real pressure, const Real temperature, Real & s, Real & ds_dp, Real & ds_dT) const
69 : {
70 : Real v, e, dv_dp, dv_dT, de_dp, de_dT;
71 1 : v_e_from_p_T(pressure, temperature, v, dv_dp, dv_dT, e, de_dp, de_dT);
72 :
73 : Real ds_dv, ds_de;
74 1 : s_from_v_e(v, e, s, ds_dv, ds_de);
75 1 : ds_dp = ds_dv * dv_dp + ds_de * de_dp;
76 1 : ds_dT = ds_dv * dv_dT + ds_de * de_dT;
77 1 : }
78 :
79 : Real
80 0 : SinglePhaseFluidProperties::s_from_v_e(const Real v, const Real e) const
81 : {
82 0 : const Real p0 = _p_initial_guess;
83 0 : const Real T0 = _T_initial_guess;
84 : Real p, T;
85 0 : bool conversion_succeeded = true;
86 0 : p_T_from_v_e(v, e, p0, T0, p, T, conversion_succeeded);
87 0 : const Real s = s_from_p_T(p, T);
88 0 : return s;
89 : }
90 :
91 : void
92 0 : SinglePhaseFluidProperties::s_from_v_e(
93 : const Real v, const Real e, Real & s, Real & ds_dv, Real & ds_de) const
94 : {
95 0 : const Real p0 = _p_initial_guess;
96 0 : const Real T0 = _T_initial_guess;
97 : Real p, T;
98 0 : bool conversion_succeeded = true;
99 0 : p_T_from_v_e(v, e, p0, T0, p, T, conversion_succeeded);
100 0 : s = s_from_p_T(p, T);
101 0 : ds_dv = p / T;
102 0 : ds_de = 1 / T;
103 0 : }
104 :
105 : Real
106 150 : SinglePhaseFluidProperties::c_from_p_T(Real p, Real T) const
107 : {
108 : Real v, e;
109 150 : v_e_from_p_T(p, T, v, e);
110 150 : return c_from_v_e(v, e);
111 : }
112 :
113 : void
114 1 : SinglePhaseFluidProperties::c_from_p_T(Real p, Real T, Real & c, Real & dc_dp, Real & dc_dT) const
115 : {
116 : Real v, e, dv_dp, dv_dT, de_dp, de_dT;
117 1 : v_e_from_p_T(p, T, v, dv_dp, dv_dT, e, de_dp, de_dT);
118 :
119 : Real dc_dv, dc_de;
120 1 : c_from_v_e(v, e, c, dc_dv, dc_de);
121 1 : dc_dp = dc_dv * dv_dp + dc_de * de_dp;
122 1 : dc_dT = dc_dv * dv_dT + dc_de * de_dT;
123 1 : }
124 :
125 : Real
126 30 : SinglePhaseFluidProperties::mu_from_p_T(Real p, Real T) const
127 : {
128 : Real v, e;
129 30 : v_e_from_p_T(p, T, v, e);
130 30 : return mu_from_v_e(v, e);
131 : }
132 :
133 : void
134 1 : SinglePhaseFluidProperties::mu_from_p_T(
135 : Real p, Real T, Real & mu, Real & dmu_dp, Real & dmu_dT) const
136 : {
137 : Real v, e, dv_dp, dv_dT, de_dp, de_dT;
138 1 : v_e_from_p_T(p, T, v, dv_dp, dv_dT, e, de_dp, de_dT);
139 :
140 : Real dmu_dv, dmu_de;
141 1 : mu_from_v_e(v, e, mu, dmu_dv, dmu_de);
142 1 : dmu_dp = dmu_dv * dv_dp + dmu_de * de_dp;
143 1 : dmu_dT = dmu_dv * dv_dT + dmu_de * de_dT;
144 1 : }
145 :
146 : Real
147 78 : SinglePhaseFluidProperties::cv_from_p_T(Real p, Real T) const
148 : {
149 : Real v, e;
150 78 : v_e_from_p_T(p, T, v, e);
151 78 : return cv_from_v_e(v, e);
152 : }
153 :
154 : void
155 1 : SinglePhaseFluidProperties::cv_from_p_T(
156 : Real p, Real T, Real & cv, Real & dcv_dp, Real & dcv_dT) const
157 : {
158 : Real v, e, dv_dp, dv_dT, de_dp, de_dT;
159 1 : v_e_from_p_T(p, T, v, dv_dp, dv_dT, e, de_dp, de_dT);
160 :
161 : Real dcv_dv, dcv_de;
162 1 : cv_from_v_e(v, e, cv, dcv_dv, dcv_de);
163 1 : dcv_dp = dcv_dv * dv_dp + dcv_de * de_dp;
164 1 : dcv_dT = dcv_dv * dv_dT + dcv_de * de_dT;
165 1 : }
166 :
167 : Real
168 78 : SinglePhaseFluidProperties::cp_from_p_T(Real p, Real T) const
169 : {
170 : Real v, e;
171 78 : v_e_from_p_T(p, T, v, e);
172 78 : return cp_from_v_e(v, e);
173 : }
174 :
175 : void
176 1 : SinglePhaseFluidProperties::cp_from_p_T(
177 : Real p, Real T, Real & cp, Real & dcp_dp, Real & dcp_dT) const
178 : {
179 : Real v, e, dv_dp, dv_dT, de_dp, de_dT;
180 1 : v_e_from_p_T(p, T, v, dv_dp, dv_dT, e, de_dp, de_dT);
181 :
182 : Real dcp_dv, dcp_de;
183 1 : cp_from_v_e(v, e, cp, dcp_dv, dcp_de);
184 1 : dcp_dp = dcp_dv * dv_dp + dcp_de * de_dp;
185 1 : dcp_dT = dcp_dv * dv_dT + dcp_de * de_dT;
186 1 : }
187 :
188 : Real
189 78 : SinglePhaseFluidProperties::k_from_p_T(Real p, Real T) const
190 : {
191 : Real v, e;
192 78 : v_e_from_p_T(p, T, v, e);
193 78 : return k_from_v_e(v, e);
194 : }
195 :
196 : void
197 1 : SinglePhaseFluidProperties::k_from_p_T(Real p, Real T, Real & k, Real & dk_dp, Real & dk_dT) const
198 : {
199 : Real v, e, dv_dp, dv_dT, de_dp, de_dT;
200 1 : v_e_from_p_T(p, T, v, dv_dp, dv_dT, e, de_dp, de_dT);
201 :
202 : Real dk_dv, dk_de;
203 1 : k_from_v_e(v, e, k, dk_dv, dk_de);
204 1 : dk_dp = dk_dv * dv_dp + dk_de * de_dp;
205 1 : dk_dT = dk_dv * dv_dT + dk_de * de_dT;
206 1 : }
207 :
208 : Real
209 30605 : SinglePhaseFluidProperties::h_from_v_e(Real v, Real e) const
210 : {
211 30605 : return e + v * p_from_v_e(v, e);
212 : }
213 :
214 : void
215 30 : SinglePhaseFluidProperties::h_from_v_e(Real v, Real e, Real & h, Real & dh_dv, Real & dh_de) const
216 : {
217 : Real p, dp_dv, dp_de;
218 30 : p_from_v_e(v, e, p, dp_dv, dp_de);
219 30 : h = e + v * p;
220 30 : dh_dv = p + v * dp_dv;
221 30 : dh_de = 1 + v * dp_de;
222 30 : }
223 :
224 : Real
225 685 : SinglePhaseFluidProperties::e_from_p_T(Real p, Real T) const
226 : {
227 685 : const Real rho = rho_from_p_T(p, T);
228 685 : return e_from_p_rho(p, rho);
229 : }
230 :
231 : void
232 1 : SinglePhaseFluidProperties::e_from_p_T(Real p, Real T, Real & e, Real & de_dp, Real & de_dT) const
233 : {
234 : // From rho(p,T), compute: drho(p,T)/dp, drho(p,T)/dT
235 1 : Real rho = 0., drho_dp = 0., drho_dT = 0.;
236 1 : rho_from_p_T(p, T, rho, drho_dp, drho_dT);
237 :
238 : // From e(p, rho), compute: de(p,rho)/dp, de(p,rho)/drho
239 1 : Real depr_dp = 0., depr_drho = 0.;
240 1 : e_from_p_rho(p, rho, e, depr_dp, depr_drho);
241 : // Using partial derivative rules, we have:
242 : // de(p,T)/dp = de(p,rho)/dp * dp/dp + de(p,rho)/drho * drho(p,T)/dp, (dp/dp == 1)
243 : // de(p,T)/dT = de(p,rho)/dp * dp/dT + de(p,rho)/drho * drho(p,T)/dT, (dp/dT == 0)
244 1 : de_dp = depr_dp + depr_drho * drho_dp;
245 1 : de_dT = depr_drho * drho_dT;
246 1 : }
247 :
248 : Real
249 72 : SinglePhaseFluidProperties::v_from_p_T(Real p, Real T) const
250 : {
251 72 : const Real rho = rho_from_p_T(p, T);
252 72 : return 1.0 / rho;
253 : }
254 :
255 : void
256 928839 : SinglePhaseFluidProperties::v_from_p_T(Real p, Real T, Real & v, Real & dv_dp, Real & dv_dT) const
257 : {
258 : Real rho, drho_dp, drho_dT;
259 928839 : rho_from_p_T(p, T, rho, drho_dp, drho_dT);
260 :
261 928839 : v = 1.0 / rho;
262 928839 : const Real dv_drho = -1.0 / (rho * rho);
263 :
264 928839 : dv_dp = dv_drho * drho_dp;
265 928839 : dv_dT = dv_drho * drho_dT;
266 928839 : }
267 :
268 : void
269 0 : SinglePhaseFluidProperties::beta_from_p_T(Real, Real, Real &, Real &, Real &) const
270 : {
271 0 : mooseError(__PRETTY_FUNCTION__, " is not implemented.");
272 : }
273 :
274 : Real
275 0 : SinglePhaseFluidProperties::beta_from_p_T(Real p, Real T) const
276 : {
277 : // The volumetric thermal expansion coefficient is defined as
278 : // 1/v dv/dT)_p
279 : // It is the fractional change rate of volume with respect to temperature change
280 : // at constant pressure. Here it is coded as
281 : // - 1/rho drho/dT)_p
282 : // using chain rule with v = v(rho)
283 :
284 : Real rho, drho_dp, drho_dT;
285 0 : rho_from_p_T(p, T, rho, drho_dp, drho_dT);
286 0 : return -drho_dT / rho;
287 : }
288 :
289 : Real
290 0 : SinglePhaseFluidProperties::molarMass() const
291 : {
292 0 : mooseError(__PRETTY_FUNCTION__, " not implemented.");
293 : }
294 :
295 : std::string
296 0 : SinglePhaseFluidProperties::fluidName() const
297 : {
298 0 : return std::string("");
299 : }
300 :
301 : Real
302 0 : SinglePhaseFluidProperties::criticalPressure() const
303 : {
304 0 : mooseError(__PRETTY_FUNCTION__, " not implemented.");
305 : }
306 :
307 : Real
308 0 : SinglePhaseFluidProperties::criticalTemperature() const
309 : {
310 0 : mooseError(__PRETTY_FUNCTION__, " not implemented.");
311 : }
312 :
313 : Real
314 2 : SinglePhaseFluidProperties::criticalDensity() const
315 : {
316 2 : return rho_from_p_T(criticalPressure(), criticalTemperature());
317 : }
318 :
319 : Real
320 1 : SinglePhaseFluidProperties::criticalInternalEnergy() const
321 : {
322 1 : return e_from_p_rho(criticalPressure(), criticalDensity());
323 : }
324 :
325 : Real
326 0 : SinglePhaseFluidProperties::triplePointPressure() const
327 : {
328 0 : mooseError(__PRETTY_FUNCTION__, " not implemented.");
329 : }
330 :
331 : Real
332 0 : SinglePhaseFluidProperties::triplePointTemperature() const
333 : {
334 0 : mooseError(__PRETTY_FUNCTION__, " not implemented.");
335 : }
336 :
337 : Real
338 0 : SinglePhaseFluidProperties::gamma_from_v_e(Real v, Real e) const
339 : {
340 0 : return cp_from_v_e(v, e) / cv_from_v_e(v, e);
341 : }
342 :
343 : void
344 0 : SinglePhaseFluidProperties::gamma_from_v_e(
345 : Real v, Real e, Real & gamma, Real & dgamma_dv, Real & dgamma_de) const
346 : {
347 0 : unimplementedDerivativeMethod(__PRETTY_FUNCTION__);
348 :
349 0 : dgamma_dv = 0.0;
350 0 : dgamma_de = 0.0;
351 0 : gamma = gamma_from_v_e(v, e);
352 0 : }
353 :
354 : Real
355 0 : SinglePhaseFluidProperties::gamma_from_p_T(Real p, Real T) const
356 : {
357 0 : return cp_from_p_T(p, T) / cv_from_p_T(p, T);
358 : }
359 :
360 : void
361 0 : SinglePhaseFluidProperties::gamma_from_p_T(
362 : Real p, Real T, Real & gamma, Real & dgamma_dp, Real & dgamma_dT) const
363 : {
364 0 : unimplementedDerivativeMethod(__PRETTY_FUNCTION__);
365 :
366 0 : dgamma_dp = 0.0;
367 0 : dgamma_dT = 0.0;
368 0 : gamma = gamma_from_p_T(p, T);
369 0 : }
370 :
371 : Real
372 0 : SinglePhaseFluidProperties::vaporPressure(Real) const
373 : {
374 0 : mooseError(__PRETTY_FUNCTION__, " not implemented.");
375 : }
376 :
377 : std::vector<Real>
378 0 : SinglePhaseFluidProperties::henryCoefficients() const
379 : {
380 0 : mooseError(__PRETTY_FUNCTION__, " not implemented.");
381 : }
382 :
383 : void
384 0 : SinglePhaseFluidProperties::vaporPressure(Real T, Real & p, Real & dp_dT) const
385 : {
386 0 : unimplementedDerivativeMethod(__PRETTY_FUNCTION__);
387 :
388 0 : dp_dT = 0.0;
389 0 : p = vaporPressure(T);
390 0 : }
391 :
392 : ADReal
393 0 : SinglePhaseFluidProperties::vaporPressure(const ADReal & T) const
394 : {
395 0 : Real p = 0.0;
396 0 : Real temperature = T.value();
397 0 : Real dpdT = 0.0;
398 :
399 0 : vaporPressure(temperature, p, dpdT);
400 :
401 0 : ADReal result = p;
402 0 : result.derivatives() = T.derivatives() * dpdT;
403 :
404 0 : return result;
405 : }
406 :
407 : Real
408 0 : SinglePhaseFluidProperties::vaporTemperature(Real) const
409 : {
410 0 : mooseError(__PRETTY_FUNCTION__, " not implemented.");
411 : }
412 :
413 : void
414 0 : SinglePhaseFluidProperties::vaporTemperature(Real p, Real & T, Real & dT_dp) const
415 : {
416 0 : unimplementedDerivativeMethod(__PRETTY_FUNCTION__);
417 :
418 0 : dT_dp = 0.0;
419 0 : T = vaporTemperature(p);
420 0 : }
421 :
422 : ADReal
423 0 : SinglePhaseFluidProperties::vaporTemperature(const ADReal & p) const
424 : {
425 0 : Real T = 0.0;
426 0 : Real pressure = p.value();
427 0 : Real dTdp = 0.0;
428 :
429 0 : vaporTemperature(pressure, T, dTdp);
430 :
431 0 : ADReal result = T;
432 0 : result.derivatives() = p.derivatives() * dTdp;
433 :
434 0 : return result;
435 : }
436 :
437 : void
438 9 : SinglePhaseFluidProperties::rho_e_from_p_T(Real p,
439 : Real T,
440 : Real & rho,
441 : Real & drho_dp,
442 : Real & drho_dT,
443 : Real & e,
444 : Real & de_dp,
445 : Real & de_dT) const
446 : {
447 9 : rho_from_p_T(p, T, rho, drho_dp, drho_dT);
448 9 : e_from_p_T(p, T, e, de_dp, de_dT);
449 9 : }
450 :
451 : void
452 4 : SinglePhaseFluidProperties::rho_mu_from_p_T(Real p, Real T, Real & rho, Real & mu) const
453 : {
454 4 : rho = rho_from_p_T(p, T);
455 4 : mu = mu_from_p_T(p, T);
456 4 : }
457 :
458 : void
459 4 : SinglePhaseFluidProperties::rho_mu_from_p_T(Real p,
460 : Real T,
461 : Real & rho,
462 : Real & drho_dp,
463 : Real & drho_dT,
464 : Real & mu,
465 : Real & dmu_dp,
466 : Real & dmu_dT) const
467 : {
468 4 : rho_from_p_T(p, T, rho, drho_dp, drho_dT);
469 4 : mu_from_p_T(p, T, mu, dmu_dp, dmu_dT);
470 4 : }
471 :
472 : void
473 0 : SinglePhaseFluidProperties::rho_mu_from_p_T(const ADReal & p,
474 : const ADReal & T,
475 : ADReal & rho,
476 : ADReal & mu) const
477 : {
478 0 : rho = rho_from_p_T(p, T);
479 0 : mu = mu_from_p_T(p, T);
480 0 : }
481 :
482 : Real
483 0 : SinglePhaseFluidProperties::e_spndl_from_v(Real) const
484 : {
485 0 : mooseError(__PRETTY_FUNCTION__, " not implemented.");
486 : }
487 :
488 : void
489 0 : SinglePhaseFluidProperties::v_e_spndl_from_T(Real, Real &, Real &) const
490 : {
491 0 : mooseError(__PRETTY_FUNCTION__, " not implemented.");
492 : }
493 :
494 : Real
495 0 : SinglePhaseFluidProperties::T_from_p_h(Real p, Real h) const
496 : {
497 0 : const Real s = s_from_h_p(h, p);
498 0 : const Real rho = rho_from_p_s(p, s);
499 0 : const Real v = 1. / rho;
500 0 : const Real e = e_from_v_h(v, h);
501 0 : return T_from_v_e(v, e);
502 : }
503 :
504 : Real
505 0 : SinglePhaseFluidProperties::p_from_h_s(Real h, Real s) const
506 : {
507 0 : Real p0 = _p_initial_guess;
508 0 : Real T0 = _T_initial_guess;
509 : Real p, T;
510 0 : bool conversion_succeeded = true;
511 0 : p_T_from_h_s(h, s, p0, T0, p, T, conversion_succeeded);
512 0 : return p;
513 : }
514 :
515 : void
516 0 : SinglePhaseFluidProperties::p_from_h_s(Real h, Real s, Real & p, Real & dp_dh, Real & dp_ds) const
517 : {
518 0 : Real p0 = _p_initial_guess;
519 0 : Real T0 = _T_initial_guess;
520 : Real T;
521 0 : bool conversion_succeeded = true;
522 0 : p_T_from_h_s(h, s, p0, T0, p, T, conversion_succeeded);
523 0 : dp_dh = rho_from_p_T(p, T);
524 0 : dp_ds = -T * rho_from_p_T(p, T);
525 0 : }
526 :
527 : void
528 0 : SinglePhaseFluidProperties::T_from_p_h(Real p, Real h, Real & T, Real & dT_dp, Real & dT_dh) const
529 : {
530 : Real s, ds_dh, ds_dp;
531 0 : s_from_h_p(h, p, s, ds_dh, ds_dp);
532 :
533 : Real rho, drho_dp_partial, drho_ds;
534 0 : rho_from_p_s(p, s, rho, drho_dp_partial, drho_ds);
535 0 : const Real drho_dp = drho_dp_partial + drho_ds * ds_dp;
536 0 : const Real drho_dh = drho_ds * ds_dh;
537 :
538 0 : const Real v = 1.0 / rho;
539 0 : const Real dv_drho = -1.0 / (rho * rho);
540 0 : const Real dv_dp = dv_drho * drho_dp;
541 0 : const Real dv_dh = dv_drho * drho_dh;
542 :
543 : Real e, de_dv, de_dh_partial;
544 0 : e_from_v_h(v, h, e, de_dv, de_dh_partial);
545 0 : const Real de_dp = de_dv * dv_dp;
546 0 : const Real de_dh = de_dh_partial + de_dv * dv_dh;
547 :
548 : Real dT_dv, dT_de;
549 0 : T_from_v_e(v, e, T, dT_dv, dT_de);
550 0 : dT_dp = dT_dv * dv_dp + dT_de * de_dp;
551 0 : dT_dh = dT_dv * dv_dh + dT_de * de_dh;
552 0 : }
553 :
554 : #pragma GCC diagnostic pop
|