Buildings.Fluid.Storage

Package with thermal energy storage models

Information

This package contains thermal energy storage models.

Extends from Modelica.Icons.VariantsPackage (Icon for package containing variants).

Package Content

Name Description
Buildings.Fluid.Storage.UsersGuide UsersGuide User's Guide
Buildings.Fluid.Storage.ExpansionVessel ExpansionVessel Expansion vessel with fixed pressure
Buildings.Fluid.Storage.Stratified Stratified Model of a stratified tank for thermal energy storage
Buildings.Fluid.Storage.StratifiedEnhanced StratifiedEnhanced Stratified tank model with enhanced discretization
Buildings.Fluid.Storage.StratifiedEnhancedInternalHex StratifiedEnhancedInternalHex A model of a water storage tank with a secondary loop and intenral heat exchanger
Buildings.Fluid.Storage.Examples Examples Collection of models that illustrate model use and test models
Buildings.Fluid.Storage.Validation Validation Collection of models that validate the storage models
Buildings.Fluid.Storage.BaseClasses BaseClasses Package with base classes for Buildings.Fluid.Storage

Buildings.Fluid.Storage.ExpansionVessel Buildings.Fluid.Storage.ExpansionVessel

Expansion vessel with fixed pressure

Buildings.Fluid.Storage.ExpansionVessel

Information

This is a model of a pressure expansion vessel. The vessel has a constant pressure that is equal to the value of the parameter p_start. The model takes into account the energy and mass balance of the medium. It has no heat exchange with the ambient.

The expansion vessel needs to be used in closed loops that contain water to set a reference pressure and, for liquids where the density is modeled as a function of temperature, to allow for the thermal expansion of the liquid.

Note that alternatively, the model Buildings.Fluid.Sources.Boundary_pT may be used to set a reference pressure. The main difference between these two models is that in this model, there is an energy and mass balance for the volume. In contrast, for Buildings.Fluid.Sources.Boundary_pT, any mass flow rate that flows out of the model will be at a user-specified temperature. Therefore, Buildings.Fluid.Sources.Boundary_pT leads to smaller systems of equations, which may result in faster simulation.

Extends from Buildings.Fluid.Interfaces.LumpedVolumeDeclarations (Declarations for lumped volumes).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
VolumeV_start Volume of liquid stored in the vessel at the start of the simulation [m3]
PressurepMedium.p_defaultConstant pressure of the expansion vessel [Pa]
Dynamics
Equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Type of energy balance: dynamic (3 initialization options) or steady state
DynamicsmassDynamicsModelica.Fluid.Types.Dynamic...Type of mass balance: dynamic (3 initialization options) or steady state
RealmSenFac1Factor for scaling the sensible thermal mass of the volume
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.)

Connectors

TypeNameDescription
FluidPort_aport_aFluid port

Modelica definition

model ExpansionVessel "Expansion vessel with fixed pressure" extends Buildings.Fluid.Interfaces.LumpedVolumeDeclarations( final energyDynamics=Modelica.Fluid.Types.Dynamics.FixedInitial, final massDynamics=Modelica.Fluid.Types.Dynamics.FixedInitial, final mSenFac=1); parameter Modelica.SIunits.Volume V_start(start=1) "Volume of liquid stored in the vessel at the start of the simulation"; parameter Modelica.SIunits.Pressure p = Medium.p_default "Constant pressure of the expansion vessel"; Modelica.Fluid.Interfaces.FluidPort_a port_a( redeclare package Medium = Medium) "Fluid port"; Modelica.SIunits.Mass m "Mass of liquid in the vessel"; protected final parameter Medium.ThermodynamicState state_start = Medium.setState_pTX( T=T_start, p=p_start, X=X_start[1:Medium.nXi]) "Medium state at start values"; final parameter Modelica.SIunits.Density rho_start=Medium.density( state=state_start) "Density, used to compute start and guess values"; Modelica.SIunits.Energy H "Internal energy of fluid"; Modelica.SIunits.Mass[Medium.nXi] mXi "Masses of independent components in the fluid"; Modelica.SIunits.Mass[Medium.nC] mC "Masses of trace substances in the fluid"; Medium.ExtraProperty C[Medium.nC](nominal=C_nominal) "Trace substance mixture content"; initial equation m = V_start * rho_start; H = m*Medium.specificInternalEnergy( Medium.setState_pTX(p=p_start, T=T_start, X= X_start[1:Medium.nXi])); mXi = m*X_start[1:Medium.nXi]; mC = m*C_start[1:Medium.nC]; equation assert(m > 1.0E-8, "Expansion vessel is undersized. You need to increase the value of the parameter V_start."); // Conservation equations der(m) = port_a.m_flow; der(H) = port_a.m_flow * actualStream(port_a.h_outflow); der(mXi) = port_a.m_flow * actualStream(port_a.Xi_outflow); der(mC) = port_a.m_flow * actualStream(port_a.C_outflow); // Properties of outgoing flow. // The port pressure is set to a constant value. port_a.p = p_start; m*port_a.h_outflow = H; m*port_a.Xi_outflow = mXi; m*port_a.C_outflow = mC; end ExpansionVessel;

Buildings.Fluid.Storage.Stratified Buildings.Fluid.Storage.Stratified

Model of a stratified tank for thermal energy storage

Buildings.Fluid.Storage.Stratified

Information

This is a model of a stratified storage tank.

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

For a model with enhanced stratification, use Buildings.Fluid.Storage.StratifiedEnhanced.

Extends from Buildings.Fluid.Storage.BaseClasses.PartialStratified (Partial model of a stratified tank for thermal energy storage).

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
Equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Formulation of energy balance
DynamicsmassDynamicsenergyDynamicsFormulation of mass 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
FluidPort_afluPorVol[nSeg]Fluid port that connects to the control volumes of the tank

Modelica definition

model Stratified "Model of a stratified tank for thermal energy storage" extends Buildings.Fluid.Storage.BaseClasses.PartialStratified(vol(each nPorts=3)); Modelica.Fluid.Interfaces.FluidPort_a fluPorVol[nSeg]( redeclare each final package Medium = Medium) "Fluid port that connects to the control volumes of the tank"; equation connect(port_a, vol[1].ports[1]); connect(vol[nSeg].ports[2], port_b); for i in 1:(nSeg-1) loop connect(vol[i].ports[2], vol[i + 1].ports[1]); end for; for i in 1:nSeg loop connect(fluPorVol[i], vol[i].ports[3]); end for; end Stratified;

Buildings.Fluid.Storage.StratifiedEnhanced Buildings.Fluid.Storage.StratifiedEnhanced

Stratified tank model with enhanced discretization

Buildings.Fluid.Storage.StratifiedEnhanced

Information

This is a model of a stratified storage tank for thermal energy storage.

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

Limitations

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

Extends from BaseClasses.PartialStratified (Partial model of a stratified tank for thermal energy storage).

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)]
IntegernSeg4Number 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
Equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Formulation of energy balance
DynamicsmassDynamicsenergyDynamicsFormulation of mass 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 StratifiedEnhanced "Stratified tank model with enhanced discretization" extends BaseClasses.PartialStratified( nSeg=4, vol(each prescribedHeatFlowRate=true, each nPorts=3)); protected Buildings.Fluid.Sensors.EnthalpyFlowRate H_a_flow( redeclare package Medium = Medium, final m_flow_nominal=m_flow_nominal, final tau=0, final allowFlowReversal=allowFlowReversal, final m_flow_small=m_flow_small) "Enthalpy flow rate at port a"; Buildings.Fluid.Sensors.EnthalpyFlowRate[nSeg - 1] H_vol_flow( redeclare package Medium = Medium, each final m_flow_nominal=m_flow_nominal, each final tau=0, each final allowFlowReversal=allowFlowReversal, each final m_flow_small=m_flow_small) "Enthalpy flow rate between the volumes"; Buildings.Fluid.Sensors.EnthalpyFlowRate H_b_flow( redeclare package Medium = Medium, final m_flow_nominal=m_flow_nominal, final tau=0, final allowFlowReversal=allowFlowReversal, final m_flow_small=m_flow_small) "Enthalpy flow rate at port b"; BaseClasses.ThirdOrderStratifier str( redeclare package Medium = Medium, nSeg=nSeg, m_flow_small=m_flow_small) "Model to reduce numerical dissipation"; Modelica.Blocks.Sources.RealExpression mTan_flow(y=port_a.m_flow) "Mass flow rate at port a"; equation connect(H_a_flow.port_b, vol[1].ports[1]); connect(vol[nSeg].ports[2], H_b_flow.port_a); connect(H_b_flow.port_b, port_b); for i in 1:(nSeg-1) loop connect(vol[i].ports[2], H_vol_flow[i].port_a); connect(H_vol_flow[i].port_b, vol[i + 1].ports[1]); end for; connect(port_a, H_a_flow.port_a); connect(vol[1:nSeg].ports[3], str.fluidPort[2:nSeg+1]); connect(H_a_flow.H_flow, str.H_flow[1]); connect(H_vol_flow[1:nSeg-1].H_flow, str.H_flow[2:nSeg]); connect(H_b_flow.H_flow, str.H_flow[nSeg + 1]); connect(str.heatPort, vol.heatPort); connect(port_a, str.fluidPort[1]); connect(port_b, str.fluidPort[nSeg + 2]); connect(mTan_flow.y, str.m_flow); end StratifiedEnhanced;

Buildings.Fluid.Storage.StratifiedEnhancedInternalHex Buildings.Fluid.Storage.StratifiedEnhancedInternalHex

A model of a water storage tank with a secondary loop and intenral heat exchanger

Buildings.Fluid.Storage.StratifiedEnhancedInternalHex

Information

This is a model of a stratified storage tank for thermal energy storage with built-in heat exchanger.

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

Limitations

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

Extends from StratifiedEnhanced (Stratified tank model with enhanced discretization).

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)]
IntegernSeg4Number of volume segments
Timetau1Time constant for mixing [s]
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
Heat exchanger
replaceable package MediumHexModelica.Media.Interfaces.Pa...Medium in the heat exchanger
HeighthHex_a Height of portHex_a of the heat exchanger, measured from tank bottom [m]
HeighthHex_b Height of portHex_b of the heat exchanger, measured from tank bottom [m]
IntegerhexSegMult2Number of heat exchanger segments in each tank segment
DiameterdExtHex0.025Exterior diameter of the heat exchanger pipe [m]
HeatFlowRateQ_flow_nominal Heat transfer at nominal conditions [W]
TemperatureTTan_nominal Temperature of fluid inside the tank at nominal heat transfer conditions [K]
TemperatureTHex_nominal Temperature of fluid inside the heat exchanger at nominal heat transfer conditions [K]
Realr_nominal0.5Ratio between coil inside and outside convective heat transfer at nominal heat transfer conditions
MassFlowRatemHex_flow_nominal Nominal mass flow rate through the heat exchanger [kg/s]
PressureDifferencedpHex_nominal2500Pressure drop across the heat exchanger at nominal conditions [Pa]
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Heat exchanger
BooleanallowFlowReversalHextrue= true to allow flow reversal in heat exchanger, false restricts to design direction (portHex_a -> portHex_b)
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
Equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Formulation of energy balance
DynamicsmassDynamicsenergyDynamicsFormulation of mass 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
Flow resistance heat exchanger
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 heat exchanger
Equations
DynamicsenergyDynamicsHexModelica.Fluid.Types.Dynamic...Formulation of energy balance for heat exchanger internal fluid mass
DynamicsmassDynamicsHexenergyDynamicsHexFormulation of mass balance for heat exchanger
DynamicsenergyDynamicsHexSolidenergyDynamicsHexFormulation of energy balance for heat exchanger solid mass
LengthlHexrTan*abs(segHex_a - segHex_b...Approximate length of the heat exchanger [m]
AreaACroHex(dExtHex^2 - (0.8*dExtHex)^2...Cross sectional area of the heat exchanger [m2]
SpecificHeatCapacitycHex490Specific heat capacity of the heat exchanger material [J/(kg.K)]
DensitydHex8000Density of the heat exchanger material [kg/m3]
HeatCapacityCHexACroHex*lHex*dHex*cHexCapacitance of the heat exchanger without the fluid [J/K]

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
FluidPort_aportHex_aHeat exchanger inlet
FluidPort_bportHex_bHeat exchanger outlet
Heat exchanger
replaceable package MediumHexMedium in the heat exchanger

Modelica definition

model StratifiedEnhancedInternalHex "A model of a water storage tank with a secondary loop and intenral heat exchanger" extends StratifiedEnhanced; replaceable package MediumHex = Modelica.Media.Interfaces.PartialMedium "Medium in the heat exchanger"; parameter Modelica.SIunits.Height hHex_a "Height of portHex_a of the heat exchanger, measured from tank bottom"; parameter Modelica.SIunits.Height hHex_b "Height of portHex_b of the heat exchanger, measured from tank bottom"; parameter Integer hexSegMult(min=1) = 2 "Number of heat exchanger segments in each tank segment"; parameter Modelica.SIunits.Diameter dExtHex = 0.025 "Exterior diameter of the heat exchanger pipe"; parameter Modelica.SIunits.HeatFlowRate Q_flow_nominal "Heat transfer at nominal conditions"; parameter Modelica.SIunits.Temperature TTan_nominal "Temperature of fluid inside the tank at nominal heat transfer conditions"; parameter Modelica.SIunits.Temperature THex_nominal "Temperature of fluid inside the heat exchanger at nominal heat transfer conditions"; parameter Real r_nominal(min=0, max=1)=0.5 "Ratio between coil inside and outside convective heat transfer at nominal heat transfer conditions"; parameter Modelica.SIunits.MassFlowRate mHex_flow_nominal "Nominal mass flow rate through the heat exchanger"; parameter Modelica.SIunits.PressureDifference dpHex_nominal(displayUnit="Pa") = 2500 "Pressure drop across the heat exchanger at nominal conditions"; parameter Boolean computeFlowResistance=true "=true, compute flow resistance. Set to false to assume no friction"; parameter Boolean from_dp=false "= true, use m_flow = f(dp) else dp = f(m_flow)"; parameter Boolean linearizeFlowResistance=false "= true, use linear relation between m_flow and dp for any flow rate"; parameter Real deltaM=0.1 "Fraction of nominal flow rate where flow transitions to laminar"; parameter Modelica.Fluid.Types.Dynamics energyDynamicsHex= Modelica.Fluid.Types.Dynamics.DynamicFreeInitial "Formulation of energy balance for heat exchanger internal fluid mass"; parameter Modelica.Fluid.Types.Dynamics massDynamicsHex= energyDynamicsHex "Formulation of mass balance for heat exchanger"; parameter Modelica.Fluid.Types.Dynamics energyDynamicsHexSolid=energyDynamicsHex "Formulation of energy balance for heat exchanger solid mass"; parameter Modelica.SIunits.Length lHex= rTan*abs(segHex_a-segHex_b)*Modelica.Constants.pi "Approximate length of the heat exchanger"; parameter Modelica.SIunits.Area ACroHex= (dExtHex^2-(0.8*dExtHex)^2)*Modelica.Constants.pi/4 "Cross sectional area of the heat exchanger"; parameter Modelica.SIunits.SpecificHeatCapacity cHex=490 "Specific heat capacity of the heat exchanger material"; parameter Modelica.SIunits.Density dHex=8000 "Density of the heat exchanger material"; parameter Modelica.SIunits.HeatCapacity CHex= ACroHex*lHex*dHex*cHex "Capacitance of the heat exchanger without the fluid"; parameter Boolean allowFlowReversalHex = true "= true to allow flow reversal in heat exchanger, false restricts to design direction (portHex_a -> portHex_b)"; Modelica.Fluid.Interfaces.FluidPort_a portHex_a( redeclare final package Medium =MediumHex, m_flow(min=if allowFlowReversalHex then -Modelica.Constants.inf else 0)) "Heat exchanger inlet"; Modelica.Fluid.Interfaces.FluidPort_b portHex_b( redeclare final package Medium = MediumHex, m_flow(max=if allowFlowReversalHex then Modelica.Constants.inf else 0)) "Heat exchanger outlet"; BaseClasses.IndirectTankHeatExchanger indTanHex( redeclare final package MediumTan = Medium, redeclare final package MediumHex = MediumHex, final nSeg=nSegHex, final CHex=CHex, final volHexFlu=volHexFlu, final Q_flow_nominal=Q_flow_nominal, final TTan_nominal=TTan_nominal, final THex_nominal=THex_nominal, final r_nominal=r_nominal, final dExtHex=dExtHex, final dp_nominal=dpHex_nominal, final m_flow_nominal=mHex_flow_nominal, final energyDynamics=energyDynamicsHex, final energyDynamicsSolid=energyDynamicsHexSolid, final massDynamics=massDynamicsHex, final computeFlowResistance=computeFlowResistance, from_dp=from_dp, final linearizeFlowResistance=linearizeFlowResistance, final deltaM=deltaM, final allowFlowReversal=allowFlowReversalHex, final m_flow_small=1e-4*abs(mHex_flow_nominal)) "Heat exchanger inside the tank"; Modelica.SIunits.HeatFlowRate QHex_flow = -sum(indTanHex.port.Q_flow) "Heat transferred from the heat exchanger to the tank"; protected final parameter Integer segHex_a = nSeg-integer(hHex_a/segHeight) "Tank segment in which port a1 of the heat exchanger is located in"; final parameter Integer segHex_b = nSeg-integer(hHex_b/segHeight) "Tank segment in which port b1 of the heat exchanger is located in"; final parameter Modelica.SIunits.Height segHeight = hTan/nSeg "Height of each tank segment (relative to bottom of same segment)"; final parameter Modelica.SIunits.Length dHHex = abs(hHex_a-hHex_b) "Vertical distance between the heat exchanger inlet and outlet"; final parameter Modelica.SIunits.Volume volHexFlu= Modelica.Constants.pi * (0.8*dExtHex)^2/4 *lHex "Volume of the heat exchanger"; final parameter Integer nSegHexTan= if segHex_a > segHex_b then segHex_a-segHex_b + 1 else segHex_b-segHex_a + 1 "Number of tank segments the heat exchanger resides in"; final parameter Integer nSegHex = nSegHexTan*hexSegMult "Number of heat exchanger segments"; initial equation assert(hHex_a >= 0 and hHex_a <= hTan, "The parameter hHex_a is outside its valid range."); assert(hHex_b >= 0 and hHex_b <= hTan, "The parameter hHex_b is outside its valid range."); assert(dHHex > 0, "The parameters hHex_a and hHex_b must not be equal."); equation for j in 0:nSegHexTan-1 loop for i in 1:hexSegMult loop connect(indTanHex.port[j*hexSegMult+i], heaPorVol[segHex_a + (if hHex_a > hHex_b then j else -j)]); end for; end for; connect(portHex_a, indTanHex.port_a); connect(indTanHex.port_b, portHex_b); end StratifiedEnhancedInternalHex;