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 "ThermalFunctionSolidProperties.h" 11 : #include "Function.h" 12 : 13 : registerMooseObject("SolidPropertiesApp", ThermalFunctionSolidProperties); 14 : 15 : InputParameters 16 31 : ThermalFunctionSolidProperties::validParams() 17 : { 18 31 : InputParameters params = ThermalSolidProperties::validParams(); 19 62 : params.addRequiredParam<FunctionName>("k", "Thermal conductivity"); 20 62 : params.addRequiredParam<FunctionName>("cp", "Isobaric specific heat"); 21 62 : params.addRequiredParam<FunctionName>("rho", "Density"); 22 31 : params.addClassDescription("Function-based thermal properties."); 23 31 : return params; 24 0 : } 25 : 26 16 : ThermalFunctionSolidProperties::ThermalFunctionSolidProperties(const InputParameters & parameters) 27 : : ThermalSolidProperties(parameters), 28 16 : _k_function(getFunction("k")), 29 16 : _cp_function(getFunction("cp")), 30 32 : _rho_function(getFunction("rho")) 31 : { 32 16 : } 33 : 34 : Real 35 607 : ThermalFunctionSolidProperties::cp_from_T(const Real & T) const 36 : { 37 607 : return _cp_function.value(T); 38 : } 39 : 40 : void 41 1 : ThermalFunctionSolidProperties::cp_from_T(const Real & T, Real & cp, Real & dcp_dT) const 42 : { 43 1 : cp = cp_from_T(T); 44 1 : dcp_dT = _cp_function.timeDerivative(T); 45 1 : } 46 : 47 : Real 48 14 : ThermalFunctionSolidProperties::cp_integral(const Real & T) const 49 : { 50 14 : return _cp_function.timeIntegral(_T_zero_e, T, Point(0, 0, 0)); 51 : } 52 : 53 : Real 54 605 : ThermalFunctionSolidProperties::k_from_T(const Real & T) const 55 : { 56 605 : return _k_function.value(T); 57 : } 58 : 59 : void 60 1 : ThermalFunctionSolidProperties::k_from_T(const Real & T, Real & k, Real & dk_dT) const 61 : { 62 1 : k = k_from_T(T); 63 1 : dk_dT = _k_function.timeDerivative(T); 64 1 : } 65 : 66 : Real 67 605 : ThermalFunctionSolidProperties::rho_from_T(const Real & T) const 68 : { 69 605 : return _rho_function.value(T); 70 : } 71 : 72 : void 73 1 : ThermalFunctionSolidProperties::rho_from_T(const Real & T, Real & rho, Real & drho_dT) const 74 : { 75 1 : rho = rho_from_T(T); 76 1 : drho_dT = _rho_function.timeDerivative(T); 77 1 : }