Buildings.Fluid.Storage.BaseClasses

Package with base classes for Buildings.Fluid.Storage

Information

This package contains base classes that are used to construct the models in Buildings.Fluid.Storage.

Extends from Modelica.Icons.BasesPackage (Icon for packages containing base classes).

Package Content

Name Description
Buildings.Fluid.Storage.BaseClasses.Buoyancy Buoyancy Model to add buoyancy if there is a temperature inversion in the tank
Buildings.Fluid.Storage.BaseClasses.IndirectTankHeatExchanger IndirectTankHeatExchanger Heat exchanger typically submerged in a fluid with a second fluid circulating through it
Buildings.Fluid.Storage.BaseClasses.PartialStratified PartialStratified Partial model of a stratified tank for thermal energy storage
Buildings.Fluid.Storage.BaseClasses.PartialTwoPortInterface PartialTwoPortInterface Partial model transporting fluid between two ports without storing mass or energy
Buildings.Fluid.Storage.BaseClasses.ThirdOrderStratifier ThirdOrderStratifier Model to reduce the numerical dissipation in a tank
Buildings.Fluid.Storage.BaseClasses.Examples Examples Examples for BaseClasses models

Buildings.Fluid.Storage.BaseClasses.Buoyancy Buildings.Fluid.Storage.BaseClasses.Buoyancy

Model to add buoyancy if there is a temperature inversion in the tank

Buildings.Fluid.Storage.BaseClasses.Buoyancy

Information

This model outputs a heat flow rate that can be added to fluid volumes in order to emulate buoyancy during a temperature inversion. For simplicity, this model does not compute a buoyancy induced mass flow rate, but rather a heat flow that has the same magnitude as the enthalpy flow associated with the buoyancy induced mass flow rate.

Extends from Modelica.Blocks.Icons.Block (Basic graphical layout of input/output block).

Parameters

TypeNameDefaultDescription
replaceable package MediumModelica.Media.Interfaces.Pa...Medium model
VolumeV Volume [m3]
IntegernSeg2Number of volume segments
Timetau Time constant for mixing [s]

Connectors

TypeNameDescription
replaceable package MediumMedium model
HeatPort_aheatPort[nSeg]Heat input into the volumes

Modelica definition

model Buoyancy "Model to add buoyancy if there is a temperature inversion in the tank" extends Modelica.Blocks.Icons.Block; replaceable package Medium = Modelica.Media.Interfaces.PartialMedium "Medium model"; parameter Modelica.Units.SI.Volume V "Volume"; parameter Integer nSeg(min=2) = 2 "Number of volume segments"; parameter Modelica.Units.SI.Time tau(min=0) "Time constant for mixing"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a[nSeg] heatPort "Heat input into the volumes"; Modelica.Units.SI.HeatFlowRate[nSeg - 1] Q_flow "Heat flow rate from segment i+1 to i"; protected parameter Medium.ThermodynamicState sta_default = Medium.setState_pTX(T=Medium.T_default, p=Medium.p_default, X=Medium.X_default[1:Medium.nXi]) "Medium state at default properties"; parameter Modelica.Units.SI.Density rho_default=Medium.density(sta_default) "Density, used to compute fluid mass"; parameter Modelica.Units.SI.SpecificHeatCapacity cp_default= Medium.specificHeatCapacityCp(sta_default) "Specific heat capacity"; parameter Real k(unit="W/K") = V*rho_default*cp_default/tau/nSeg "Proportionality constant, since we use dT instead of dH"; Modelica.Units.SI.TemperatureDifference dT[nSeg - 1] "Temperature difference between adjoining volumes"; equation for i in 1:nSeg-1 loop dT[i] = heatPort[i+1].T-heatPort[i].T; Q_flow[i] = k*noEvent(smooth(1, if dT[i]>0 then dT[i]^2 else 0)); end for; heatPort[1].Q_flow = -Q_flow[1]; for i in 2:nSeg-1 loop heatPort[i].Q_flow = -Q_flow[i]+Q_flow[i-1]; end for; heatPort[nSeg].Q_flow = Q_flow[nSeg-1]; end Buoyancy;

Buildings.Fluid.Storage.BaseClasses.IndirectTankHeatExchanger Buildings.Fluid.Storage.BaseClasses.IndirectTankHeatExchanger

Heat exchanger typically submerged in a fluid with a second fluid circulating through it

Buildings.Fluid.Storage.BaseClasses.IndirectTankHeatExchanger

Information

This model is a heat exchanger with a moving fluid on one side and a stagnant fluid on the other. It is intended for use when a heat exchanger is submerged in a stagnant fluid. For example, the heat exchanger in a storage tank which is part of a solar thermal system.

This component models the fluid in the heat exchanger, convection between the fluid and the heat exchanger, and convection from the heat exchanger to the surrounding fluid.

The model is based on Buildings.Fluid.HeatExchangers.BaseClasses.HACoilInside and Buildings.Fluid.HeatExchangers.BaseClasses.HANaturalCylinder.

The fluid ports are intended to be connected to a circulated heat transfer fluid while the heat port is intended to be connected to a stagnant fluid.

Extends from Buildings.Fluid.Interfaces.TwoPortFlowResistanceParameters (Parameters for flow resistance for models with two ports), Buildings.Fluid.Interfaces.LumpedVolumeDeclarations (Declarations for lumped volumes), Buildings.Fluid.Interfaces.PartialTwoPortInterface (Partial model with two ports and declaration of quantities that are used by many models).

Parameters

TypeNameDefaultDescription
replaceable package MediumHexModelica.Media.Interfaces.Pa...Heat transfer fluid flowing through the heat exchanger
replaceable package MediumTanModelica.Media.Interfaces.Pa...Heat transfer fluid inside the tank
replaceable package MediumPartialMediumMedium in the component
IntegernSeg Number of segments in the heat exchanger
HeatCapacityCHex Capacitance of the heat exchanger [J/K]
VolumevolHexFlu Volume of heat transfer fluid in the heat exchanger [m3]
DiameterdExtHex Exterior diameter of the heat exchanger pipe [m]
Nominal condition
PressureDifferencedp_nominal Pressure difference [Pa]
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
HeatFlowRateQ_flow_nominal Heat transfer at nominal conditions [W]
TemperatureTTan_nominal Temperature of fluid inside the tank at UA_nominal [K]
TemperatureTHex_nominal Temperature of fluid inside the heat exchanger at UA_nominal [K]
Realr_nominal0.5Ratio between coil inside and outside convective heat transfer
Flow resistance
BooleancomputeFlowResistancetrue=true, compute flow resistance. Set to false to assume no friction
Booleanfrom_dpfalse= true, use m_flow = f(dp) else dp = f(m_flow)
BooleanlinearizeFlowResistancefalse= true, use linear relation between m_flow and dp for any flow rate
RealdeltaM0.1Fraction of nominal flow rate where flow transitions to laminar
Dynamics
Conservation equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Type of energy balance: dynamic (3 initialization options) or steady state
DynamicsenergyDynamicsSolidenergyDynamicsFormulation of energy balance for heat exchanger solid mass
RealmSenFac1Factor for scaling the sensible thermal mass of the volume
Advanced
Dynamics
DynamicsmassDynamicsenergyDynamicsType of mass balance: dynamic (3 initialization options) or steady state, must be steady state if energyDynamics is steady state
MassFlowRatem_flow_small1E-4*abs(m_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
Modeling detail
BooleanhA_flowDependenttrueSet to false to make the convective heat coefficient calculation of the fluid inside the coil independent of mass flow rate
BooleanhA_temperatureDependenttrueSet to false to make the convective heat coefficient calculation of the fluid inside the coil independent of temperature
Initialization
AbsolutePressurep_startMedium.p_defaultStart value of pressure [Pa]
TemperatureT_startMedium.T_defaultStart value of temperature [K]
MassFractionX_start[Medium.nX]Medium.X_defaultStart value of mass fractions m_i/m [kg/kg]
ExtraPropertyC_start[Medium.nC]fill(0, Medium.nC)Start value of trace substances
ExtraPropertyC_nominal[Medium.nC]fill(1E-2, Medium.nC)Nominal value of trace substances. (Set to typical order of magnitude.)
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal

Connectors

TypeNameDescription
replaceable package MediumHexHeat transfer fluid flowing through the heat exchanger
replaceable package MediumTanHeat transfer fluid inside the tank
replaceable package MediumMedium in the component
FluidPort_aport_aFluid connector a (positive design flow direction is from port_a to port_b)
FluidPort_bport_bFluid connector b (positive design flow direction is from port_a to port_b)
HeatPort_aport[nSeg]Heat port connected to water inside the tank

Modelica definition

model IndirectTankHeatExchanger "Heat exchanger typically submerged in a fluid with a second fluid circulating through it" replaceable package MediumHex = Modelica.Media.Interfaces.PartialMedium "Heat transfer fluid flowing through the heat exchanger"; replaceable package MediumTan = Modelica.Media.Interfaces.PartialMedium "Heat transfer fluid inside the tank"; extends Buildings.Fluid.Interfaces.TwoPortFlowResistanceParameters; extends Buildings.Fluid.Interfaces.LumpedVolumeDeclarations( final massDynamics=energyDynamics, redeclare final package Medium = MediumHex); extends Buildings.Fluid.Interfaces.PartialTwoPortInterface( redeclare final package Medium = MediumHex, final show_T=false); constant Boolean homotopyInitialization = true "= true, use homotopy method"; parameter Integer nSeg(min=2) "Number of segments in the heat exchanger"; parameter Modelica.Units.SI.HeatCapacity CHex "Capacitance of the heat exchanger"; parameter Modelica.Units.SI.Volume volHexFlu "Volume of heat transfer fluid in the heat exchanger"; parameter Modelica.Units.SI.HeatFlowRate Q_flow_nominal "Heat transfer at nominal conditions"; final parameter Modelica.Units.SI.ThermalConductance UA_nominal=abs( Q_flow_nominal/(THex_nominal - TTan_nominal)) "Nominal UA value for the heat exchanger"; parameter Modelica.Units.SI.Temperature TTan_nominal "Temperature of fluid inside the tank at UA_nominal"; parameter Modelica.Units.SI.Temperature THex_nominal "Temperature of fluid inside the heat exchanger at UA_nominal"; parameter Real r_nominal(min=0, max=1)=0.5 "Ratio between coil inside and outside convective heat transfer"; parameter Modelica.Units.SI.Diameter dExtHex "Exterior diameter of the heat exchanger pipe"; parameter Modelica.Fluid.Types.Dynamics energyDynamics=Modelica.Fluid.Types.Dynamics.DynamicFreeInitial "Formulation of energy balance for heat exchanger internal fluid mass"; parameter Modelica.Fluid.Types.Dynamics energyDynamicsSolid=energyDynamics "Formulation of energy balance for heat exchanger solid mass"; parameter Boolean hA_flowDependent = true "Set to false to make the convective heat coefficient calculation of the fluid inside the coil independent of mass flow rate"; parameter Boolean hA_temperatureDependent = true "Set to false to make the convective heat coefficient calculation of the fluid inside the coil independent of temperature"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a port[nSeg] "Heat port connected to water inside the tank"; Buildings.Fluid.FixedResistances.PressureDrop res( redeclare final package Medium = MediumHex, final dp_nominal=dp_nominal, final m_flow_nominal=m_flow_nominal, final allowFlowReversal=allowFlowReversal, final homotopyInitialization=homotopyInitialization, final show_T=show_T, final from_dp=from_dp, final linearized=linearizeFlowResistance) "Calculates the flow resistance and pressure drop through the heat exchanger"; Buildings.Fluid.MixingVolumes.MixingVolume vol[nSeg]( redeclare each package Medium = MediumHex, each nPorts=2, each m_flow_nominal=m_flow_nominal, each V=volHexFlu/nSeg, each energyDynamics=energyDynamics, each massDynamics=energyDynamics, each p_start=p_start, each T_start=T_start, each X_start=X_start, each C_start=C_start, each C_nominal=C_nominal, each final prescribedHeatFlowRate=false, each final allowFlowReversal=allowFlowReversal) "Heat exchanger fluid"; Modelica.Thermal.HeatTransfer.Components.HeatCapacitor cap[nSeg]( each C=CHex/nSeg, each T(start=T_start, fixed=(energyDynamicsSolid == Modelica.Fluid.Types.Dynamics.FixedInitial)), each der_T( fixed=(energyDynamicsSolid == Modelica.Fluid.Types.Dynamics.SteadyStateInitial))) if not energyDynamicsSolid == Modelica.Fluid.Types.Dynamics.SteadyState "Thermal mass of the heat exchanger"; protected Buildings.Fluid.Sensors.MassFlowRate senMasFlo( redeclare package Medium = MediumHex, allowFlowReversal=allowFlowReversal) "Mass flow rate of the heat transfer fluid"; Modelica.Thermal.HeatTransfer.Components.Convection htfToHex[nSeg] "Convection coefficient between the heat transfer fluid and heat exchanger"; Modelica.Thermal.HeatTransfer.Components.Convection HexToTan[nSeg] "Convection coefficient between the heat exchanger and the surrounding medium"; Modelica.Thermal.HeatTransfer.Sensors.TemperatureSensor temSenHex[nSeg] "Temperature of the heat transfer fluid"; Modelica.Thermal.HeatTransfer.Sensors.TemperatureSensor temSenWat[nSeg] "Temperature sensor of the fluid surrounding the heat exchanger"; Modelica.Blocks.Routing.Replicator rep(nout=nSeg) "Replicates senMasFlo signal from 1 seg to nSeg"; Buildings.Fluid.HeatExchangers.BaseClasses.HACoilInside hAPipIns[nSeg]( each m_flow_nominal=m_flow_nominal, each hA_nominal=UA_nominal/nSeg*(r_nominal + 1)/r_nominal, each T_nominal=THex_nominal, each final flowDependent=hA_flowDependent, each final temperatureDependent=hA_temperatureDependent) "Computation of convection coefficients inside the coil"; Buildings.Fluid.HeatExchangers.BaseClasses.HANaturalCylinder hANatCyl[nSeg]( redeclare each final package Medium = Medium, each final ChaLen=dExtHex, each final hA_nominal=UA_nominal/nSeg*(1 + r_nominal), each final TFlu_nominal=TTan_nominal, each final TSur_nominal=TTan_nominal - (r_nominal/(1 + r_nominal))*( TTan_nominal - THex_nominal)) "Calculates an hA value for each side of the heat exchanger"; Modelica.Thermal.HeatTransfer.Sensors.TemperatureSensor temSenSur[nSeg] "Temperature at the external surface of the heat exchanger"; initial equation assert(homotopyInitialization, "In " + getInstanceName() + ": The constant homotopyInitialization has been modified from its default value. This constant will be removed in future releases.", level = AssertionLevel.warning); equation for i in 1:(nSeg - 1) loop connect(vol[i].ports[2], vol[i + 1].ports[1]); end for; connect(rep.u,senMasFlo. m_flow); connect(port,HexToTan. fluid); connect(vol[1].ports[1],senMasFlo.port_b); connect(cap.port,HexToTan.solid); connect(vol.heatPort,htfToHex. fluid); connect(htfToHex.solid,HexToTan. solid); connect(temSenHex.T, hAPipIns.T); connect(hAPipIns.hA, htfToHex.Gc); connect(HexToTan.solid,temSenSur. port); connect(temSenWat.port, port); connect(temSenSur.T, hANatCyl.TSur); connect(hANatCyl.TFlu, temSenWat.T); connect(port_a, senMasFlo.port_a); connect(vol[nSeg].ports[2], res.port_a); connect(res.port_b, port_b); connect(temSenHex.port, vol.heatPort); connect(rep.y, hAPipIns.m_flow); connect(hANatCyl.hA,HexToTan. Gc); end IndirectTankHeatExchanger;

Buildings.Fluid.Storage.BaseClasses.PartialStratified Buildings.Fluid.Storage.BaseClasses.PartialStratified

Partial model of a stratified tank for thermal energy storage

Buildings.Fluid.Storage.BaseClasses.PartialStratified

Information

This is a partial model of a stratified storage tank.

See the Buildings.Fluid.Storage.UsersGuide for more information.

Implementation

This model does not include the ports that connect to the fluid from the outside, because these ports cannot be used for the models that contain the Buildings.Fluid.Storage.BaseClasses.ThirdOrderStratifier.

Extends from Buildings.Fluid.Storage.BaseClasses.PartialTwoPortInterface (Partial model transporting fluid between two ports without storing mass or energy).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
VolumeVTan Tank volume [m3]
LengthhTan Height of tank (without insulation) [m]
LengthdIns Thickness of insulation [m]
ThermalConductivitykIns0.04Specific heat conductivity of insulation [W/(m.K)]
IntegernSeg2Number of volume segments
Timetau1Time constant for mixing [s]
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
MassFlowRatem_flow_small1E-4*abs(m_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
Dynamics
Conservation equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Formulation of energy balance
Initialization
AbsolutePressurep_startMedium.p_defaultStart value of pressure [Pa]
TemperatureT_startMedium.T_defaultStart value of temperature [K]
TemperatureTFlu_start[nSeg]T_start*ones(nSeg)Initial temperature of the tank segments, with TFlu_start[1] being the top segment [K]
MassFractionX_start[Medium.nX]Medium.X_defaultStart value of mass fractions m_i/m [kg/kg]
ExtraPropertyC_start[Medium.nC]fill(0, Medium.nC)Start value of trace substances

Connectors

TypeNameDescription
FluidPort_aport_aFluid connector a (positive design flow direction is from port_a to port_b)
FluidPort_bport_bFluid connector b (positive design flow direction is from port_a to port_b)
output RealOutputQl_flowHeat loss of tank (positive if heat flows from tank to ambient)
HeatPort_aheaPorVol[nSeg]Heat port that connects to the control volumes of the tank
HeatPort_aheaPorSidHeat port tank side (outside insulation)
HeatPort_aheaPorTopHeat port tank top (outside insulation)
HeatPort_aheaPorBotHeat port tank bottom (outside insulation). Leave unconnected for adiabatic condition

Modelica definition

model PartialStratified "Partial model of a stratified tank for thermal energy storage" extends Buildings.Fluid.Storage.BaseClasses.PartialTwoPortInterface; import Modelica.Fluid.Types; import Modelica.Fluid.Types.Dynamics; parameter Modelica.Units.SI.Volume VTan "Tank volume"; parameter Modelica.Units.SI.Length hTan "Height of tank (without insulation)"; parameter Modelica.Units.SI.Length dIns "Thickness of insulation"; parameter Modelica.Units.SI.ThermalConductivity kIns=0.04 "Specific heat conductivity of insulation"; parameter Integer nSeg(min=2) = 2 "Number of volume segments"; //////////////////////////////////////////////////////////////////// // Assumptions parameter Types.Dynamics energyDynamics=Modelica.Fluid.Types.Dynamics.FixedInitial "Formulation of energy balance"; // Initialization parameter Medium.AbsolutePressure p_start = Medium.p_default "Start value of pressure"; parameter Medium.Temperature T_start=Medium.T_default "Start value of temperature"; parameter Modelica.Units.SI.Temperature TFlu_start[nSeg]=T_start*ones(nSeg) "Initial temperature of the tank segments, with TFlu_start[1] being the top segment"; parameter Medium.MassFraction X_start[Medium.nX] = Medium.X_default "Start value of mass fractions m_i/m"; parameter Medium.ExtraProperty C_start[Medium.nC]( quantity=Medium.extraPropertiesNames)=fill(0, Medium.nC) "Start value of trace substances"; // Dynamics parameter Modelica.Units.SI.Time tau=1 "Time constant for mixing"; //////////////////////////////////////////////////////////////////// // Connectors Modelica.Blocks.Interfaces.RealOutput Ql_flow "Heat loss of tank (positive if heat flows from tank to ambient)"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a[nSeg] heaPorVol "Heat port that connects to the control volumes of the tank"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a heaPorSid "Heat port tank side (outside insulation)"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a heaPorTop "Heat port tank top (outside insulation)"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a heaPorBot "Heat port tank bottom (outside insulation). Leave unconnected for adiabatic condition"; // Models Buildings.Fluid.MixingVolumes.MixingVolume[nSeg] vol( redeclare each package Medium = Medium, each energyDynamics=energyDynamics, each massDynamics=energyDynamics, each p_start=p_start, T_start=TFlu_start, each X_start=X_start, each C_start=C_start, each V=VTan/nSeg, each m_flow_nominal=m_flow_nominal, each final mSenFac=1, each final m_flow_small=m_flow_small, each final allowFlowReversal=allowFlowReversal) "Tank segment"; protected parameter Medium.ThermodynamicState sta_default = Medium.setState_pTX( T=Medium.T_default, p=Medium.p_default, X=Medium.X_default[1:Medium.nXi]) "Medium state at default properties"; parameter Modelica.Units.SI.Length hSeg=hTan/nSeg "Height of a tank segment"; parameter Modelica.Units.SI.Area ATan=VTan/hTan "Tank cross-sectional area (without insulation)"; parameter Modelica.Units.SI.Length rTan=sqrt(ATan/Modelica.Constants.pi) "Tank diameter (without insulation)"; parameter Modelica.Units.SI.ThermalConductance conFluSeg=ATan* Medium.thermalConductivity(sta_default)/hSeg "Thermal conductance between fluid volumes"; parameter Modelica.Units.SI.ThermalConductance conTopSeg=ATan*kIns/dIns "Thermal conductance from center of top (or bottom) volume through tank insulation at top (or bottom)"; BaseClasses.Buoyancy buo( redeclare final package Medium = Medium, final V=VTan, final nSeg=nSeg, final tau=tau) "Model to prevent unstable tank stratification"; Modelica.Thermal.HeatTransfer.Components.ThermalConductor[nSeg - 1] conFlu( each G=conFluSeg) "Thermal conductance in fluid between the segments"; Modelica.Thermal.HeatTransfer.Components.ThermalConductor[nSeg] conWal( each G=2*Modelica.Constants.pi*kIns*hSeg/Modelica.Math.log((rTan+dIns)/rTan)) "Thermal conductance through tank wall"; Modelica.Thermal.HeatTransfer.Components.ThermalConductor conTop( G=conTopSeg) "Thermal conductance through tank top"; Modelica.Thermal.HeatTransfer.Components.ThermalConductor conBot( G=conTopSeg) "Thermal conductance through tank bottom"; Modelica.Thermal.HeatTransfer.Sensors.HeatFlowSensor heaFloTop "Heat flow at top of tank (outside insulation)"; Modelica.Thermal.HeatTransfer.Sensors.HeatFlowSensor heaFloBot "Heat flow at bottom of tank (outside insulation)"; Modelica.Thermal.HeatTransfer.Sensors.HeatFlowSensor heaFloSid[nSeg] "Heat flow at wall of tank (outside insulation)"; Modelica.Blocks.Routing.Multiplex3 mul( n1=1, n2=nSeg, n3=1) "Multiplex to collect heat flow rates"; Modelica.Blocks.Math.Sum sum1(nin=nSeg + 2); Modelica.Thermal.HeatTransfer.Components.ThermalCollector theCol(m=nSeg) "Connector to assign multiple heat ports to one heat port"; equation connect(buo.heatPort, vol.heatPort); for i in 1:nSeg-1 loop // heat conduction between fluid nodes connect(vol[i].heatPort, conFlu[i].port_a); connect(vol[i+1].heatPort, conFlu[i].port_b); end for; connect(vol[1].heatPort, conTop.port_a); connect(vol.heatPort, conWal.port_a); connect(conBot.port_a, vol[nSeg].heatPort); connect(vol.heatPort, heaPorVol); connect(conWal.port_b, heaFloSid.port_a); connect(conTop.port_b, heaFloTop.port_a); connect(conBot.port_b, heaFloBot.port_a); connect(heaFloTop.port_b, heaPorTop); connect(heaFloBot.port_b, heaPorBot); connect(heaFloTop.Q_flow, mul.u1[1]); connect(heaFloSid.Q_flow, mul.u2); connect(heaFloBot.Q_flow, mul.u3[1]); connect(mul.y, sum1.u); connect(sum1.y, Ql_flow); connect(heaFloSid.port_b, theCol.port_a); connect(theCol.port_b, heaPorSid); end PartialStratified;

Buildings.Fluid.Storage.BaseClasses.PartialTwoPortInterface Buildings.Fluid.Storage.BaseClasses.PartialTwoPortInterface

Partial model transporting fluid between two ports without storing mass or energy

Buildings.Fluid.Storage.BaseClasses.PartialTwoPortInterface

Information

This partial class implements the same functionality as Buildings.Fluid.Interfaces.PartialTwoPortInterface, except that port_a and port_b are placed at the top and bottom of the component.

Implementation

The implementation is done in this package as opposed to Buildings.Fluid.Interfaces as it is only used by the storage model, and may be removed when the tool limitations that are discussed in IBPSA, #1794. are removed.

Parameters

TypeNameDefaultDescription
replaceable package MediumModelica.Media.Interfaces.Pa...Medium in the component
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
MassFlowRatem_flow_small1E-4*abs(m_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed

Connectors

TypeNameDescription
replaceable package MediumMedium in the component
FluidPort_aport_aFluid connector a (positive design flow direction is from port_a to port_b)
FluidPort_bport_bFluid connector b (positive design flow direction is from port_a to port_b)

Modelica definition

partial model PartialTwoPortInterface "Partial model transporting fluid between two ports without storing mass or energy" replaceable package Medium = Modelica.Media.Interfaces.PartialMedium "Medium in the component"; parameter Boolean allowFlowReversal = true "= false to simplify equations, assuming, but not enforcing, no flow reversal"; Modelica.Fluid.Interfaces.FluidPort_a port_a( redeclare final package Medium = Medium, m_flow(min=if allowFlowReversal then -Modelica.Constants.inf else 0), p(start=Medium.p_default), h_outflow(start = Medium.h_default, nominal = Medium.h_default)) "Fluid connector a (positive design flow direction is from port_a to port_b)"; Modelica.Fluid.Interfaces.FluidPort_b port_b( redeclare final package Medium = Medium, m_flow(max=if allowFlowReversal then +Modelica.Constants.inf else 0), p(start=Medium.p_default), h_outflow(start = Medium.h_default, nominal = Medium.h_default)) "Fluid connector b (positive design flow direction is from port_a to port_b)"; parameter Modelica.Units.SI.MassFlowRate m_flow_nominal "Nominal mass flow rate"; parameter Modelica.Units.SI.MassFlowRate m_flow_small(min=0) = 1E-4*abs( m_flow_nominal) "Small mass flow rate for regularization of zero flow"; // Diagnostics parameter Boolean show_T = false "= true, if actual temperature at port is computed"; Modelica.Units.SI.MassFlowRate m_flow(start=_m_flow_start) = port_a.m_flow "Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction)"; Modelica.Units.SI.PressureDifference dp( start=_dp_start, displayUnit="Pa") = port_a.p - port_b.p "Pressure difference between port_a and port_b"; Medium.ThermodynamicState sta_a= if allowFlowReversal then Medium.setState_phX(port_a.p, noEvent(actualStream(port_a.h_outflow)), noEvent(actualStream(port_a.Xi_outflow))) else Medium.setState_phX(port_a.p, noEvent(inStream(port_a.h_outflow)), noEvent(inStream(port_a.Xi_outflow))) if show_T "Medium properties in port_a"; Medium.ThermodynamicState sta_b= if allowFlowReversal then Medium.setState_phX(port_b.p, noEvent(actualStream(port_b.h_outflow)), noEvent(actualStream(port_b.Xi_outflow))) else Medium.setState_phX(port_b.p, noEvent(port_b.h_outflow), noEvent(port_b.Xi_outflow)) if show_T "Medium properties in port_b"; protected final parameter Modelica.Units.SI.MassFlowRate _m_flow_start=0 "Start value for m_flow, used to avoid a warning if not set in m_flow, and to avoid m_flow.start in parameter window"; final parameter Modelica.Units.SI.PressureDifference _dp_start(displayUnit= "Pa") = 0 "Start value for dp, used to avoid a warning if not set in dp, and to avoid dp.start in parameter window"; end PartialTwoPortInterface;

Buildings.Fluid.Storage.BaseClasses.ThirdOrderStratifier Buildings.Fluid.Storage.BaseClasses.ThirdOrderStratifier

Model to reduce the numerical dissipation in a tank

Buildings.Fluid.Storage.BaseClasses.ThirdOrderStratifier

Information

This model reduces the numerical dissipation that is introduced by the standard first-order upwind discretization scheme which is created when connecting fluid volumes in series.

The model is used in conjunction with Buildings.Fluid.Storage.Stratified. It computes a heat flux that needs to be added to each volume of Buildings.Fluid.Storage.Stratified in order to give the results that a third-order upwind discretization scheme (QUICK) would give.

The QUICK method can cause oscillations in the tank temperatures since the high order method introduces numerical dispersion. There are two ways to reduce the oscillations:

Both approaches are implemented in the model.

The model is used by Buildings.Fluid.Storage.StratifiedEnhanced.

Limitations

The model requires at least 4 fluid segments. Hence, set nSeg to 4 or higher.

Extends from Modelica.Blocks.Icons.Block (Basic graphical layout of input/output block).

Parameters

TypeNameDefaultDescription
replaceable package MediumModelica.Media.Interfaces.Pa...Medium model
MassFlowRatem_flow_small Small mass flow rate for regularization of zero flow [kg/s]
IntegernSeg Number of volume segments
Realalpha0.5Under-relaxation coefficient (1: QUICK; 0: 1st order upwind)

Connectors

TypeNameDescription
replaceable package MediumMedium model
HeatPort_aheatPort[nSeg]Heat input into the volumes
input RealInputm_flowMass flow rate from port a to port b
input RealInputH_flow[nSeg + 1]Enthalpy flow between the volumes
FluidPort_afluidPort[nSeg + 2]Fluid port, needed to get pressure, temperature and species concentration

Modelica definition

model ThirdOrderStratifier "Model to reduce the numerical dissipation in a tank" extends Modelica.Blocks.Icons.Block; replaceable package Medium = Modelica.Media.Interfaces.PartialMedium "Medium model"; parameter Modelica.Units.SI.MassFlowRate m_flow_small(min=0) "Small mass flow rate for regularization of zero flow"; parameter Integer nSeg(min=4) "Number of volume segments"; parameter Real alpha( min=0, max=1) = 0.5 "Under-relaxation coefficient (1: QUICK; 0: 1st order upwind)"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a[nSeg] heatPort "Heat input into the volumes"; Modelica.Blocks.Interfaces.RealInput m_flow "Mass flow rate from port a to port b"; Modelica.Blocks.Interfaces.RealInput[nSeg + 1] H_flow "Enthalpy flow between the volumes"; Modelica.Fluid.Interfaces.FluidPort_a[nSeg + 2] fluidPort(redeclare each package Medium = Medium) "Fluid port, needed to get pressure, temperature and species concentration"; protected Modelica.Units.SI.SpecificEnthalpy[nSeg + 1] hOut "Extended vector with new outlet enthalpies to reduce numerical dissipation (at the boundary between two volumes)"; Modelica.Units.SI.SpecificEnthalpy[nSeg + 2] h "Extended vector with port enthalpies, needed to simplify loop"; Modelica.Units.SI.HeatFlowRate Q_flow[nSeg] "Heat exchange computed using upwind third order discretization scheme"; // Modelica.Units.SI.HeatFlowRate Q_flow_upWind // "Heat exchange computed using upwind third order discretization scheme"; //Used to test the energy conservation Real sig "Sign used to implement the third order upwind scheme without triggering a state event"; Real comSig "Sign used to implement the third order upwind scheme without triggering a state event"; equation assert(nSeg >= 4, "Number of segments of the enhanced stratified tank should be no less than 4 (nSeg>=4)."); // assign zero flow conditions at port fluidPort[:].m_flow = zeros(nSeg + 2); fluidPort[:].h_outflow = zeros(nSeg + 2); fluidPort[:].Xi_outflow = zeros(nSeg + 2, Medium.nXi); fluidPort[:].C_outflow = zeros(nSeg + 2, Medium.nC); // assign extended enthalpy vectors for i in 1:nSeg + 2 loop h[i] = inStream(fluidPort[i].h_outflow); end for; // Value that transitions between 0 and 1 as the flow reverses. sig = Modelica.Fluid.Utilities.regStep( m_flow, 1, 0, m_flow_small); // at surface between port_a and vol1 comSig = 1 - sig; // at surface between port_a and vol1 hOut[1] = sig*h[1] + comSig*h[2]; // at surface between vol[nSeg] and port_b hOut[nSeg + 1] = sig*h[nSeg + 1] + comSig*h[nSeg + 2]; // Pros: These two equations can further reduce the temperature overshoot by using the upwind // Cons: The minimum of nSeg hase to be 4 instead of 2. hOut[2] = sig*h[2] + comSig*h[3]; // at surface between vol1 and vol2 hOut[nSeg] = sig*h[nSeg] + comSig*h[nSeg + 1]; // at surface between vol[nSeg-1] and vol[nSeg] for i in 3:nSeg - 1 loop // at surface between vol[i-1] and vol[i] // QUICK method hOut[i] = 0.5*(h[i] + h[i + 1]) - comSig*0.125*(h[i + 2] + h[i] - 2*h[i + 1]) - sig*0.125*(h[i - 1] + h[i + 1] - 2*h[i]); // hOut[i] = 0.5*(h[i]+h[i+1]); // Central difference method end for; for i in 1:nSeg loop // difference between QUICK and UPWIND; index of H_flow is same as hOut Q_flow[i] = m_flow*(hOut[i + 1] - hOut[i]) - (H_flow[i + 1] - H_flow[i]); end for; // Q_flow_upWind = sum(Q_flow[i] for i in 1:nSeg); //Used to test the energy conservation for i in 1:nSeg loop // Add the difference back to the volume as heat flow. An under-relaxation is needed to reduce // oscillations caused by high order method heatPort[i].Q_flow = Q_flow[i]*alpha; end for; end ThirdOrderStratifier;