LBL logo

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

NameDescription
Buildings.Fluid.Storage.ExpansionVessel ExpansionVessel Pressure expansion vessel with fixed gas cushion
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.Examples Examples Collection of models that illustrate model use and test models
Buildings.Fluid.Storage.BaseClasses BaseClasses Package with base classes for Buildings.Fluid.Storage


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

Pressure expansion vessel with fixed gas cushion

Buildings.Fluid.Storage.ExpansionVessel

Information

This is a model of a pressure expansion vessel. The vessel has a fixed total volume. A fraction of the volume is occupied by a fixed mass of gas, and the other fraction is occupied by the liquid that flows through the port. The pressure p in the vessel is

 VGas0 * p_start = (VTot-VLiquid) * p
where VGas0 is the initial volume occupied by the gas, p_start is the initial pressure, VTot is the total volume of the vessel and VLiquid is the amount of liquid in the vessel.

Optionally, a heat port can be activated by setting use_HeatTransfer=true. This heat port connects directly to the liquid. The gas does not participate in the energy balance.

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.FixedBoundary 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.FixedBoundary, any mass flow rate that flows out of the model will be at a user-specified temperature. Therefore, Buildings.Fluid.Sources.FixedBoundary 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
VolumeVTot Total volume of vessel (water and gas side) [m3]
VolumeVGas0VTot/2Initial volume of gas [m3]
PressurepMax5E5Maximum pressure before simulation stops with an error [Pa]
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]
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 
  "Pressure expansion vessel with fixed gas cushion"
 extends Buildings.Fluid.Interfaces.LumpedVolumeDeclarations;
 parameter Modelica.SIunits.Volume VTot 
    "Total volume of vessel (water and gas side)";
 parameter Modelica.SIunits.Volume VGas0 = VTot/2 "Initial volume of gas";

//////////////////////////////////////////////////////////////
  Modelica.Fluid.Interfaces.FluidPort_a port_a(redeclare package Medium =
        Medium) "Fluid port";

  parameter Modelica.SIunits.Pressure pMax = 5E5 
    "Maximum pressure before simulation stops with an error";

 Modelica.SIunits.Volume VLiq "Volume of liquid in the vessel";
  // We set m(start=(VTot-VGas0)*1000, stateSelect=StateSelect.always)
  // since the mass accumulated in the volume should be a state.
  // This often leads to smaller systems of equations.
protected 
  Buildings.Fluid.Interfaces.ConservationEquation vol(
    redeclare final package Medium = Medium,
    m(start=(VTot-VGas0)*1000, stateSelect=StateSelect.always),
    final nPorts=1,
    final energyDynamics=energyDynamics,
    final massDynamics=massDynamics,
    final p_start=p_start,
    final T_start=T_start,
    final X_start=X_start,
    final C_start=C_start,
    final C_nominal=C_nominal,
    final fluidVolume = VLiq) 
    "Model for mass and energy balance of water in expansion vessel";
  Modelica.Blocks.Sources.Constant heaInp(final k=0) 
    "Block to set heat input into volume";

  Modelica.Blocks.Sources.Constant       masExc[Medium.nXi](final k=zeros(Medium.nXi)) if 
        Medium.nXi > 0 "Block to set mass exchange in volume";
equation 
  assert(port_a.p < pMax, "Pressure exceeds maximum pressure.\n" +
       "You may need to increase VTot of the ExpansionVessel.");

  // Water content and pressure
  port_a.p * (VTot-VLiq) = p_start * VGas0;

  connect(heaInp.y, vol.Q_flow);
  connect(masExc.y, vol.mXi_flow);
  connect(port_a, vol.ports[1]);
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. The tank uses several volumes to model the stratification. Heat conduction is modeled between the volumes through the fluid, and between the volumes and the ambient. The port heaPorVol may be used to connect a temperature sensor that measures the fluid temperature of an individual volume. It may also be used to add heat to individual volumes.

The tank has nSeg fluid volumes. The top volume has the index 1. Thus, to add a heating element to the bottom element, connect a heat input to heaPorVol[nSeg].

The heat ports outside the tank insulation can be used to specify an ambient temperature. Leave these ports unconnected to force adiabatic boundary conditions. Note, however, that all heat conduction elements through the tank wall (but not the top and bottom) are connected to the heat port heaPorSid. Thus, not connecting heaPorSid means an adiabatic boundary condition in the sense that heaPorSid.Q_flow = 0. This, however, still allows heat to flow through the tank walls, modelled by conWal, from one fluid volume to another one.

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

Extends from Buildings.Fluid.Interfaces.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 due to temperature inversion [s]
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
Initialization
MassFlowRatem_flow.start0Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction) [kg/s]
Pressuredp.start0Pressure difference between port_a and port_b [Pa]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_small1E-4*abs(m_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
BooleanhomotopyInitializationtrue= true, use homotopy method
Diagnostics
Booleanshow_V_flowfalse= true, if volume flow rate at inflowing port is computed
Booleanshow_Ttrue= true, if actual temperature at port is computed (may lead to events)
Dynamics
Equations
DynamicsenergyDynamicssystem.energyDynamicsFormulation of energy balance
DynamicsmassDynamicsenergyDynamicsFormulation of mass balance
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

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)
replaceable package Medium 
HeatPort_aheaPorVol[nSeg]Heat port of fluid volumes
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
output RealOutputQl_flowHeat loss of tank (positive if heat flows from tank to ambient)

Modelica definition

model Stratified 
  "Model of a stratified tank for thermal energy storage"
  extends Buildings.Fluid.Interfaces.PartialTwoPortInterface(show_T=true);

  replaceable package Medium =
      Modelica.Media.Interfaces.PartialSimpleMedium;
  import Modelica.Fluid.Types;
  import Modelica.Fluid.Types.Dynamics;
  parameter Modelica.SIunits.Volume VTan "Tank volume";
  parameter Modelica.SIunits.Length hTan "Height of tank (without insulation)";
  parameter Modelica.SIunits.Length dIns "Thickness of insulation";
  parameter Modelica.SIunits.ThermalConductivity kIns = 0.04 
    "Specific heat conductivity of insulation";
  parameter Integer nSeg(min=2) = 2 "Number of volume segments";

  ////////////////////////////////////////////////////////////////////
  // Assumptions
  parameter Types.Dynamics energyDynamics=system.energyDynamics 
    "Formulation of energy balance";
  parameter Types.Dynamics massDynamics=energyDynamics 
    "Formulation of mass 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 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";

  ////////////////////////////////////////////////////////////////////

  MixingVolumes.MixingVolume[nSeg] vol(
    redeclare each package Medium = Medium,
    each energyDynamics=energyDynamics,
    each massDynamics=massDynamics,
    each p_start=p_start,
    each T_start=T_start,
    each X_start=X_start,
    each C_start=C_start,
    each V=VTan/nSeg,
    each nPorts=nPorts,
    each m_flow_nominal = m_flow_nominal) "Tank segment";
  Sensors.EnthalpyFlowRate hA_flow(redeclare package Medium = Medium,
      m_flow_nominal=m_flow_nominal) "Enthalpy flow rate at port a";
  Sensors.EnthalpyFlowRate[nSeg-1] hVol_flow(redeclare package Medium = Medium,
      each m_flow_nominal=m_flow_nominal);
  Sensors.EnthalpyFlowRate hB_flow(redeclare package Medium = Medium,
      m_flow_nominal=m_flow_nominal) "Enthalpy flow rate at port b";
  BaseClasses.Buoyancy buo(
    redeclare package Medium = Medium,
    V=VTan,
    nSeg=nSeg,
    tau=tau) "Model to prevent unstable tank stratification";
  parameter Modelica.SIunits.Time tau=1 
    "Time constant for mixing due to temperature inversion";
  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.Interfaces.HeatPort_a[nSeg] heaPorVol 
    "Heat port of fluid volumes";
  Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a heaPorSid 
    "Heat port tank side (outside insulation)";
  Modelica.Thermal.HeatTransfer.Sensors.HeatFlowSensor heaFloSid[
                                                         nSeg] 
    "Heat flow at wall of tank (outside insulation)";
  Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a heaPorTop 
    "Heat port tank top (outside insulation)";
  Modelica.Thermal.HeatTransfer.Sensors.HeatFlowSensor heaFloTop 
    "Heat flow at top of tank (outside insulation)";
  Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a heaPorBot 
    "Heat port tank bottom (outside insulation). Leave unconnected for adiabatic condition";
  Modelica.Thermal.HeatTransfer.Sensors.HeatFlowSensor heaFloBot 
    "Heat flow at bottom of tank (outside insulation)";
  Modelica.Blocks.Interfaces.RealOutput Ql_flow 
    "Heat loss of tank (positive if heat flows from tank to ambient)";

protected 
  constant Integer nPorts = 2 "Number of ports of volume";
  parameter Modelica.SIunits.Length hSeg = hTan / nSeg 
    "Height of a tank segment";
  parameter Modelica.SIunits.Area ATan = VTan/hTan 
    "Tank cross-sectional area (without insulation)";
  parameter Modelica.SIunits.Length rTan = sqrt(ATan/Modelica.Constants.pi) 
    "Tank diameter (without insulation)";
  parameter Modelica.SIunits.ThermalConductance conFluSeg = ATan*Medium.lambda_const/hSeg 
    "Thermal conductance between fluid volumes";
  parameter Modelica.SIunits.ThermalConductance conTopSeg = ATan*kIns/dIns 
    "Thermal conductance from center of top (or bottom) volume through tank insulation at top (or bottom)";

protected 
  Modelica.Blocks.Routing.Multiplex3 mul(
    n1=1,
    n2=nSeg,
    n3=1);
  Modelica.Blocks.Math.Sum sum1(nin=nSeg + 2);
public 
  Modelica.Thermal.HeatTransfer.Components.ThermalCollector theCol(m=nSeg) 
    "Connector to assign multiple heat ports to one heat port";
equation 
  connect(hA_flow.port_b, vol[1].ports[1]);
  connect(vol[nSeg].ports[2], hB_flow.port_a);
  connect(hB_flow.port_b, port_b);
  for i in 1:(nSeg-1) loop

  connect(vol[i].ports[2], hVol_flow[i].port_a);
  connect(hVol_flow[i].port_b, vol[i+1].ports[1]);
  end for;
  connect(port_a, hA_flow.port_a);
  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);
  for i in 1:nSeg loop

  end for;

  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 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. The model is identical to Buildings.Fluid.Storage.Stratified, except that it adds a correction that reduces the numerical dissipation. The correction uses a third order upwind scheme to compute the outlet temperatures of the segments in the tank. This model is implemented in Buildings.Fluid.Storage.BaseClasses.ThirdOrderStratifier.

Limitations

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

Extends from Stratified (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 due to temperature inversion [s]
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
Initialization
MassFlowRatem_flow.start0Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction) [kg/s]
Pressuredp.start0Pressure difference between port_a and port_b [Pa]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_small1E-4*abs(m_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
BooleanhomotopyInitializationtrue= true, use homotopy method
Diagnostics
Booleanshow_V_flowfalse= true, if volume flow rate at inflowing port is computed
Booleanshow_Ttrue= true, if actual temperature at port is computed (may lead to events)
Dynamics
Equations
DynamicsenergyDynamicssystem.energyDynamicsFormulation of energy balance
DynamicsmassDynamicsenergyDynamicsFormulation of mass balance
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

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)
HeatPort_aheaPorVol[nSeg]Heat port of fluid volumes
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
output RealOutputQl_flowHeat loss of tank (positive if heat flows from tank to ambient)

Modelica definition

model StratifiedEnhanced 
  "Stratified tank model with enhanced discretization"
  extends Stratified(nSeg=4, nPorts=3, vol(each prescribedHeatFlowRate=true));
  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(vol[1:nSeg].ports[3], str.fluidPort[2:nSeg+1]);
  connect(hA_flow.H_flow, str.H_flow[1]);
  connect(hVol_flow[1:nSeg-1].H_flow, str.H_flow[2:nSeg]);
  connect(hB_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;

Automatically generated Tue Jan 8 08:30:45 2013.