Buildings.Fluid.HeatExchangers.RadiantSlabs.BaseClasses.Functions

Functions for convective heat transfer

Information

This package contains functions that are used to construct the models in Buildings.Fluid.HeatExchangers.RadiantSlabs.

Package Content

Name Description
Buildings.Fluid.HeatExchangers.RadiantSlabs.BaseClasses.Functions.AverageResistance AverageResistance Average fictitious resistance for plane that contains the pipes
Buildings.Fluid.HeatExchangers.RadiantSlabs.BaseClasses.Functions.heatFlowRate heatFlowRate Heat flow rate for epsilon-NTU model

Buildings.Fluid.HeatExchangers.RadiantSlabs.BaseClasses.Functions.AverageResistance Buildings.Fluid.HeatExchangers.RadiantSlabs.BaseClasses.Functions.AverageResistance

Average fictitious resistance for plane that contains the pipes

Information

This function computes a fictitious thermal resistance between the pipe outer wall and a fictitious, average temperature of the plane that contains the pipes. The equation is the same as is implemented in TRNSYS 17 Type 56 active layer component, manual page 197-201. Different equations are used for

Limitations

The resistance Rx is based on a steady-state heat transfer analysis. Therefore, it is only valid during steady-state. For a fully dynamic model, a finite element method for the radiant slab would need to be implemented.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
DistancedisPip pipe distance [m]
DiameterdPipOut pipe outside diameter [m]
ThermalConductivityk pipe level construction element thermal conductivity [W/(m.K)]
SystemTypesysTyp Type of radiant system
ThermalConductivitykIns floor slab insulation thermal conductivity [W/(m.K)]
ThicknessdIns floor slab insulation thickness [m]

Outputs

TypeNameDescription
ThermalInsulanceRxThermal insulance [m2.K/W]

Modelica definition

function AverageResistance "Average fictitious resistance for plane that contains the pipes" extends Modelica.Icons.Function; input Modelica.SIunits.Distance disPip "pipe distance"; input Modelica.SIunits.Diameter dPipOut "pipe outside diameter"; input Modelica.SIunits.ThermalConductivity k "pipe level construction element thermal conductivity"; input Buildings.Fluid.HeatExchangers.RadiantSlabs.Types.SystemType sysTyp "Type of radiant system"; input Modelica.SIunits.ThermalConductivity kIns "floor slab insulation thermal conductivity"; input Modelica.SIunits.Thickness dIns "floor slab insulation thickness"; output Modelica.SIunits.ThermalInsulance Rx "Thermal insulance"; protected Real cri(unit="1") "Criteria used to select formula for computation of resistance"; Real infSum "Approximation to infinite sum used to compute the thermal resistance"; Real alpha(unit="W/(m2.K)", min=0) "Criteria used to select formula for computation of resistance"; Real fac "Factor used for systems in wall or ceiling, or for capillary tubes"; algorithm if sysTyp == Buildings.Fluid.HeatExchangers.RadiantSlabs.Types.SystemType.Floor then alpha := kIns/dIns; assert(alpha < 1.212, "Warning: In RadiantAverageResistance, require alpha = kIns/dIns <= 1.212 W/(m2.K).\n" + " Obtained alpha = " + String(alpha) + " W/(m2.K)\n" + " kIns = " + String(kIns) + " W/(m.K)\n" + " dIns = " + String(dIns) + " m\n" + " For these values, the radiant slab model is outside its valid range.\n", level=AssertionLevel.warning); infSum := - sum(((alpha/k*disPip - 2*Modelica.Constants.pi*s)/ (alpha/k*disPip + 2*Modelica.Constants.pi*s)) *Modelica.Math.exp(-4*Modelica.Constants.pi*s*dIns/disPip)/s for s in 1:100); Rx := disPip*(Modelica.Math.log(disPip/Modelica.Constants.pi/dPipOut) + infSum) /(2*Modelica.Constants.pi*k); fac := 0; // not needed. elseif sysTyp == Buildings.Fluid.HeatExchangers.RadiantSlabs.Types.SystemType.Ceiling_Wall_or_Capillary then // Branch for radiant ceilings, radiant walls, and systems with capillary heat exchangers cri := disPip/dPipOut; fac := if (cri >= 5.8) then Modelica.Math.log(cri/Modelica.Constants.pi) else (cri/Modelica.Constants.pi/3); Rx := disPip/2/Modelica.Constants.pi/k * fac; else assert(sysTyp == Buildings.Fluid.HeatExchangers.RadiantSlabs.Types.SystemType.Floor or sysTyp == Buildings.Fluid.HeatExchangers.RadiantSlabs.Types.SystemType.Ceiling_Wall_or_Capillary, "Invalid value for sysTyp in \"Buildings.Fluid.HeatExchangers.RadiantSlabs.BaseClasses.Functions.AverageResistance\" Check parameters of the radiant slab model."); cri := 0; fac := 0; Rx := 1; end if; end AverageResistance;

Buildings.Fluid.HeatExchangers.RadiantSlabs.BaseClasses.Functions.heatFlowRate Buildings.Fluid.HeatExchangers.RadiantSlabs.BaseClasses.Functions.heatFlowRate

Heat flow rate for epsilon-NTU model

Information

This function computes the heat flow rate for the radiant slab.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
TemperatureT_a Temperature at port_a [K]
TemperatureT_b Temperature at port_b [K]
TemperatureT_s Temperature of solid [K]
TemperatureT_f Temperature of fluid control volume [K]
SpecificHeatCapacityc_p Specific heat capacity [J/(kg.K)]
ThermalConductanceUA UA value [W/K]
MassFlowRatem_flow Mass flow rate from port_a to port_b [kg/s]
MassFlowRatem_flow_nominal Nominal mass flow rate from port_a to port_b [kg/s]

Outputs

TypeNameDescription
HeatFlowRateQ_flowHeat flow rate [W]

Modelica definition

function heatFlowRate "Heat flow rate for epsilon-NTU model" extends Modelica.Icons.Function; input Modelica.SIunits.Temperature T_a "Temperature at port_a"; input Modelica.SIunits.Temperature T_b "Temperature at port_b"; input Modelica.SIunits.Temperature T_s "Temperature of solid"; input Modelica.SIunits.Temperature T_f "Temperature of fluid control volume"; input Modelica.SIunits.SpecificHeatCapacity c_p "Specific heat capacity"; input Modelica.SIunits.ThermalConductance UA "UA value"; input Modelica.SIunits.MassFlowRate m_flow "Mass flow rate from port_a to port_b"; input Modelica.SIunits.MassFlowRate m_flow_nominal "Nominal mass flow rate from port_a to port_b"; output Modelica.SIunits.HeatFlowRate Q_flow "Heat flow rate"; protected Modelica.SIunits.MassFlowRate m_abs_flow "Absolute value of mass flow rate"; Modelica.SIunits.Temperature T_in "Inlet fluid temperature"; Real eps "Heat transfer effectiveness"; algorithm m_abs_flow :=noEvent(abs(m_flow)); T_in :=smooth(1, noEvent(if m_flow >= 0 then T_a else T_b)); if m_abs_flow > 0.15*m_flow_nominal then // High flow rate. Use epsilon-NTU formula. eps := 1-Modelica.Math.exp(-UA/m_abs_flow/c_p); Q_flow :=eps*(T_s-T_in)*m_abs_flow*c_p; elseif (m_abs_flow < 0.05*m_flow_nominal) then // Low flow rate. Use heat conduction to temperature of the control volume. Q_flow :=UA*(T_s-T_f); else // Transition. Use a once continuously differentiable transition between the // above regimes. eps := 1-Modelica.Math.exp(-UA/m_abs_flow/c_p); Q_flow := Buildings.Utilities.Math.Functions.spliceFunction( pos=eps*(T_s-T_in)*m_abs_flow*c_p, neg=UA*(T_s-T_f), x=m_abs_flow/m_flow_nominal-0.1, deltax=0.05); end if; end heatFlowRate;