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