Buildings.Airflow.Multizone

Package with models for multizone airflow and contaminant transport

Information

This package provides models to compute the airflow and contaminant transport between different rooms and between a room and the exterior environment. See the User's Guide for more information.

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

Package Content

Name Description
Buildings.Airflow.Multizone.UsersGuide UsersGuide User's Guide
Buildings.Airflow.Multizone.DoorDiscretizedOpen DoorDiscretizedOpen Door model using discretization along height coordinate
Buildings.Airflow.Multizone.DoorDiscretizedOperable DoorDiscretizedOperable Door model using discretization along height coordinate
Buildings.Airflow.Multizone.DoorOpen DoorOpen Door model for bi-directional air flow between rooms
Buildings.Airflow.Multizone.DoorOperable DoorOperable Door model for bi-directional air flow between rooms that can be open or closed
Buildings.Airflow.Multizone.EffectiveAirLeakageArea EffectiveAirLeakageArea Effective air leakage area
Buildings.Airflow.Multizone.MediumColumn MediumColumn Vertical shaft with no friction and no storage of heat and mass
Buildings.Airflow.Multizone.MediumColumnDynamic MediumColumnDynamic Vertical shaft with no friction and storage of heat and mass
Buildings.Airflow.Multizone.Orifice Orifice Orifice
Buildings.Airflow.Multizone.ZonalFlow_ACS ZonalFlow_ACS Zonal flow with input air change per second
Buildings.Airflow.Multizone.ZonalFlow_m_flow ZonalFlow_m_flow Zonal flow with input air change per second
Buildings.Airflow.Multizone.Types Types Package with type definitions
Buildings.Airflow.Multizone.Examples Examples Collection of models that illustrate model use and test models
Buildings.Airflow.Multizone.Validation Validation Collection of validation models
Buildings.Airflow.Multizone.BaseClasses BaseClasses Package with base classes for Buildings.Airflow.Multizone

Buildings.Airflow.Multizone.DoorDiscretizedOpen Buildings.Airflow.Multizone.DoorDiscretizedOpen

Door model using discretization along height coordinate

Buildings.Airflow.Multizone.DoorDiscretizedOpen

Information

This model describes the bi-directional air flow through an open door.

To compute the bi-directional flow, the door is discretize along the height coordinate. An orifice equation is used to compute the flow for each compartment.

In this model, the door is always open. Use the model Buildings.Airflow.Multizone.DoorDiscretizedOperable for a door that can either be open or closed.

Extends from Buildings.Airflow.Multizone.BaseClasses.DoorDiscretized (Door model using discretization along height coordinate).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
VelocityvZer0.001Minimum velocity to prevent zero flow. Recommended: 0.001 [m/s]
IntegernCom10Number of compartments for the discretization
PressureDifferencedp_turbulent0.01Pressure difference where laminar and turbulent flow relation coincide. Recommended: 0.01 [Pa]
Geometry
LengthwOpe0.9Width of opening [m]
LengthhOpe2.1Height of opening [m]
LengthhA2.7/2Height of reference pressure zone A [m]
LengthhB2.7/2Height of reference pressure zone B [m]
Orifice characteristics
RealCD0.65Discharge coefficient
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
BooleanforceErrorControlOnFlowtrueFlag to force error control on m_flow. Set to true if interested in flow rate

Connectors

TypeNameDescription
FluidPort_aport_a1Fluid connector a1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_bport_b1Fluid connector b1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_aport_a2Fluid connector a2 (positive design flow direction is from port_a2 to port_b2)
FluidPort_bport_b2Fluid connector b2 (positive design flow direction is from port_a2 to port_b2)

Modelica definition

model DoorDiscretizedOpen "Door model using discretization along height coordinate" extends Buildings.Airflow.Multizone.BaseClasses.DoorDiscretized; parameter Real CD=0.65 "Discharge coefficient"; protected constant Real mFixed = 0.5 "Fixed value for flow coefficient"; constant Real gamma(min=1) = 1.5 "Normalized flow rate where dphi(0)/dpi intersects phi(1)"; constant Real a = gamma "Polynomial coefficient for regularized implementation of flow resistance"; constant Real b = 1/8*mFixed^2 - 3*gamma - 3/2*mFixed + 35.0/8 "Polynomial coefficient for regularized implementation of flow resistance"; constant Real c = -1/4*mFixed^2 + 3*gamma + 5/2*mFixed - 21.0/4 "Polynomial coefficient for regularized implementation of flow resistance"; constant Real d = 1/8*mFixed^2 - gamma - mFixed + 15.0/8 "Polynomial coefficient for regularized implementation of flow resistance"; equation m=mFixed; A = wOpe*hOpe; kVal = CD*dA*sqrt(2/rho_default); // orifice equation for i in 1:nCom loop dV_flow[i] = Buildings.Airflow.Multizone.BaseClasses.powerLawFixedM( k=kVal, dp=dpAB[i], m=mFixed, a=a, b=b, c=c, d=d, dp_turbulent=dp_turbulent); end for; end DoorDiscretizedOpen;

Buildings.Airflow.Multizone.DoorDiscretizedOperable Buildings.Airflow.Multizone.DoorDiscretizedOperable

Door model using discretization along height coordinate

Buildings.Airflow.Multizone.DoorDiscretizedOperable

Information

This model describes the bi-directional air flow through an open door.

To compute the bi-directional flow, the door is discretize along the height coordinate, and uses an orifice equation to compute the flow for each compartment.

The door can be either open or closed, depending on the input signal y. Set y=0 if the door is closed, and y=1 if the door is open. Use the model Buildings.Airflow.Multizone.DoorDiscretizedOpen for a door that is always closed.

Extends from Buildings.Airflow.Multizone.BaseClasses.DoorDiscretized (Door model using discretization along height coordinate).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
VelocityvZer0.001Minimum velocity to prevent zero flow. Recommended: 0.001 [m/s]
IntegernCom10Number of compartments for the discretization
PressureDifferencedp_turbulent0.01Pressure difference where laminar and turbulent flow relation coincide. Recommended: 0.01 [Pa]
Geometry
LengthwOpe0.9Width of opening [m]
LengthhOpe2.1Height of opening [m]
LengthhA2.7/2Height of reference pressure zone A [m]
LengthhB2.7/2Height of reference pressure zone B [m]
Rating conditions
PressureDifferencedpCloRat4Pressure drop at rating condition of closed door [Pa]
RealCDCloRat1Discharge coefficient at rating conditions of closed door
Closed door
AreaLClo Effective leakage area of closed door [m2]
RealCDClo0.65Discharge coefficient of closed door
RealmClo0.65Flow exponent for crack of closed door
Open door
RealCDOpe0.65Discharge coefficient of open door
RealmOpe0.5Flow exponent for door of open door
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
BooleanforceErrorControlOnFlowtrueFlag to force error control on m_flow. Set to true if interested in flow rate

Connectors

TypeNameDescription
FluidPort_aport_a1Fluid connector a1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_bport_b1Fluid connector b1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_aport_a2Fluid connector a2 (positive design flow direction is from port_a2 to port_b2)
FluidPort_bport_b2Fluid connector b2 (positive design flow direction is from port_a2 to port_b2)
input RealInputyOpening signal, 0=closed, 1=open [1]

Modelica definition

model DoorDiscretizedOperable "Door model using discretization along height coordinate" extends Buildings.Airflow.Multizone.BaseClasses.DoorDiscretized; parameter Modelica.SIunits.PressureDifference dpCloRat(min=0, displayUnit="Pa") = 4 "Pressure drop at rating condition of closed door"; parameter Real CDCloRat(min=0, max=1)=1 "Discharge coefficient at rating conditions of closed door"; parameter Modelica.SIunits.Area LClo(min=0) "Effective leakage area of closed door"; parameter Real CDOpe=0.65 "Discharge coefficient of open door"; parameter Real CDClo=0.65 "Discharge coefficient of closed door"; parameter Real mOpe = 0.5 "Flow exponent for door of open door"; parameter Real mClo= 0.65 "Flow exponent for crack of closed door"; Modelica.Blocks.Interfaces.RealInput y(min=0, max=1, unit="1") "Opening signal, 0=closed, 1=open"; protected parameter Modelica.SIunits.Area AOpe=wOpe*hOpe "Open aperture area"; parameter Modelica.SIunits.Area AClo(fixed=false) "Closed aperture area"; Real kOpe "Open aperture flow coefficient, k = V_flow/ dp^m"; Real kClo "Closed aperture flow coefficient, k = V_flow/ dp^m"; Real fraOpe "Fraction of aperture that is open"; initial equation AClo=CDClo/CDCloRat * LClo * dpCloRat^(0.5-mClo); equation fraOpe =y; kClo = CDClo * AClo/nCom * sqrt(2/rho_default); kOpe = CDOpe * AOpe/nCom * sqrt(2/rho_default); // flow exponent m = fraOpe*mOpe + (1-fraOpe)*mClo; // opening area A = fraOpe*AOpe + (1-fraOpe)*AClo; // friction coefficient for power law kVal = fraOpe*kOpe + (1-fraOpe)*kClo; // orifice equation for i in 1:nCom loop dV_flow[i] = Buildings.Airflow.Multizone.BaseClasses.powerLaw( k=kVal, dp=dpAB[i], m=m, dp_turbulent=dp_turbulent); end for; end DoorDiscretizedOperable;

Buildings.Airflow.Multizone.DoorOpen Buildings.Airflow.Multizone.DoorOpen

Door model for bi-directional air flow between rooms

Buildings.Airflow.Multizone.DoorOpen

Information

Model for bi-directional air flow through a large opening such as a door.

In this model, the air flow is composed of two components, a one-directional bulk air flow due to static pressure difference in the adjoining two thermal zones, and a two-directional airflow due to temperature-induced differences in density of the air in the two thermal zones. Although turbulent air flow is a nonlinear phenomenon, the model is based on the simplifying assumption that these two air flow rates can be superposed. (Superposition is only exact for laminar flow.) This assumption is made because it leads to a simple model and because there is significant uncertainty and assumptions anyway in such simplified a model for bidirectional flow through a door.

Main equations

The air flow rate due to static pressure difference is

ab,p = CD w h (2/ρ0)0.5 Δpm,

where is the volumetric air flow rate, CD is the discharge coefficient, w and h are the width and height of the opening, ρ0 is the mass density at the medium default pressure, temperature and humidity, m is the flow exponent and Δp = pa - pb is the static pressure difference between the thermal zones. For this model explanation, we will assume pa > pb. For turbulent flow, m=1/2 and for laminar flow m=1.

The air flow rate due to temperature difference in the thermal zones is ab,t for flow from thermal zone a to b, and ba,t for air flow rate from thermal zone b to a. The model has two air flow paths to allow bi-directional air flow. The mass flow rates at these two air flow paths are

a1 = ρ0   (+V̇ab,p/2 +   V̇ab,t),

and, similarly,

ba = ρ0   (-V̇ab,p/2 +   V̇ba,t),

where we simplified the calculation by using the density ρ0. To calculate ba,t, we again use the density ρ0 and because of this simplification, we can write

ab,t = -ṁba,t = ρ0   V̇ab,t = -ρ0   V̇ba,t,

from which follows that the neutral height, e.g., the height where the air flow rate due to flow induced by temperature difference is zero, is at h/2. Hence,

ab,t = CD0h/2 w v(z) dz,

where v(z) is the velocity at height z. From the Bernoulli equation, we obtain

v(z) = (2 g z Δρ ⁄ ρ0)1/2.

The density difference can be written as

Δρ = ρab ≈ ρ0 (Tb - Ta) ⁄ T0,

where we used ρa = p0 /(R Ta) and Ta Tb ≈ T02. Substituting this expression into the integral and integrating from 0 to z yields

ab,t = 1⁄3 CD w h (g h ⁄ (R T0 ρ0))1/2 Δp1/2.

The above equation is equivalent to (6) in Brown and Solvason (1962).

Main assumptions

The main assumptions are as follows:

From these assumptions follows that the neutral height for buoyancy-driven air flow is at half of the height of the opening.

Notes

For a more detailed model, use Buildings.Airflow.Multizone.DoorDiscretizedOpen.

References

Extends from Buildings.Airflow.Multizone.BaseClasses.Door (Partial door model for bi-directional flow).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Geometry
LengthwOpe0.9Width of opening [m]
LengthhOpe2.1Height of opening [m]
Custom Parameters
VelocityvABVAB_flow/AOpeAverage velocity from A to B [m/s]
VelocityvBAVBA_flow/AOpeAverage velocity from B to A [m/s]
Orifice characteristics
RealCD0.65Discharge coefficient
Realm0.5Flow coefficient
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
PressureDifferencedp_turbulent0.01Pressure difference where laminar and turbulent flow relation coincide [Pa]

Connectors

TypeNameDescription
FluidPort_aport_a1Fluid connector a1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_bport_b1Fluid connector b1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_aport_a2Fluid connector a2 (positive design flow direction is from port_a2 to port_b2)
FluidPort_bport_b2Fluid connector b2 (positive design flow direction is from port_a2 to port_b2)

Modelica definition

model DoorOpen "Door model for bi-directional air flow between rooms" extends Buildings.Airflow.Multizone.BaseClasses.Door( final vAB = VAB_flow/AOpe, final vBA = VBA_flow/AOpe); parameter Real CD=0.65 "Discharge coefficient"; parameter Real m = 0.5 "Flow coefficient"; protected constant Real gamma(min=1) = 1.5 "Normalized flow rate where dphi(0)/dpi intersects phi(1)"; constant Real a = gamma "Polynomial coefficient for regularized implementation of flow resistance"; parameter Real b = 1/8*m^2 - 3*gamma - 3/2*m + 35.0/8 "Polynomial coefficient for regularized implementation of flow resistance"; parameter Real c = -1/4*m^2 + 3*gamma + 5/2*m - 21.0/4 "Polynomial coefficient for regularized implementation of flow resistance"; parameter Real d = 1/8*m^2 - gamma - m + 15.0/8 "Polynomial coefficient for regularized implementation of flow resistance"; parameter Real kVal=CD*AOpe*sqrt(2/rho_default) "Flow coefficient, k = V_flow/ dp^m"; parameter Real kT = rho_default * CD * AOpe/3 * sqrt(Modelica.Constants.g_n /(Medium.T_default*conTP) * hOpe) "Constant coefficient for buoyancy driven air flow rate"; parameter Modelica.SIunits.MassFlowRate m_flow_turbulent= kVal * rho_default * sqrt(dp_turbulent) "Mass flow rate where regularization to laminar flow occurs for temperature-driven flow"; equation // Air flow rate due to static pressure difference VABp_flow = Buildings.Airflow.Multizone.BaseClasses.powerLawFixedM( k=kVal, dp=port_a1.p-port_a2.p, m=m, a=a, b=b, c=c, d=d, dp_turbulent=dp_turbulent); // Air flow rate due to buoyancy // Because powerLawFixedM requires as an input a pressure difference pa-pb, // we convert Ta-Tb by multiplying it with rho*R, and we divide // above the constant expression by (rho*R)^m on the right hand-side of kT. // Note that here, k is for mass flow rate, whereas in powerLawFixedM, it is for volume flow rate. // We can use m=0.5 as this is from Bernoulli. mABt_flow = Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_dp( k=kT, dp=conTP*(Medium.temperature(state_a1_inflow)-Medium.temperature(state_a2_inflow)), m_flow_turbulent=m_flow_turbulent); end DoorOpen;

Buildings.Airflow.Multizone.DoorOperable Buildings.Airflow.Multizone.DoorOperable

Door model for bi-directional air flow between rooms that can be open or closed

Buildings.Airflow.Multizone.DoorOperable

Information

Model for bi-directional air flow through a large opening such as a door which can be opened or closed based on the control input signal y.

For the control input signal y=1, this model is identical to Buildings.Airflow.Multizone.DoorOpen, and for y=0, the door is assumed to be closed and the air flow rate is set to the air flow rate through the crack posed by the open door, clo.

The air flow rate for the closed door is computed as

clo = kclo ΔpmClo,

where clo is the volume flow rate, kclo is a flow coefficient and mClo is the flow exponent. The flow coefficient is

kclo = Lclo CDCloRat ΔpRat(0.5-mClo) (2/ρ0)0.5,

where Lclo is the effective air leakage area, CDCloRat is the discharge coefficient at the reference condition, ΔpRat is the pressure drop at the rating condition, and ρ0 is the mass density at the medium default pressure, temperature and humidity.

The effective air leakage area Lclo can be obtained, for example, from the ASHRAE fundamentals (ASHRAE, 1997, p. 25.18). In the ASHRAE fundamentals, the effective air leakage area is based on a reference pressure difference of ΔpRat = 4 Pa and a discharge coefficient of CDCloRat = 1. A similar model is also used in the CONTAM software (Dols and Walton, 2002). Dols and Walton (2002) recommend to use for the flow exponent mClo=0.6 to mClo=0.7 if the flow exponent is not reported with the test results.

For the open door, the air flow rate ope is computed as described in Buildings.Airflow.Multizone.DoorOpen with the parameters CDOpe and mOpe.

The actual air flow rate is computed as

clo = (y-1) V̇clo + y V̇ope,

where y ∈ [0, 1] is the control signal. Note that for values of y that are different from 0 and 1, the model simply interpolates the air flow rate between a fully open and a fully closed door. In practice, the air flow rate would likely increase quickly if the door is slightly opened, and hence we do not claim that the model is accurate for values other than y = 0 and y = 1.

References

Extends from Buildings.Airflow.Multizone.BaseClasses.Door (Partial door model for bi-directional flow).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Geometry
LengthwOpe0.9Width of opening [m]
LengthhOpe2.1Height of opening [m]
Custom Parameters
VelocityvABVAB_flow/AAverage velocity from A to B [m/s]
VelocityvBAVBA_flow/AAverage velocity from B to A [m/s]
Open door
RealCDOpe0.65Discharge coefficient of open door
RealmOpe0.5Flow exponent for door of open door
Closed door
AreaLClo Effective leakage area of closed door [m2]
RealmClo0.65Flow exponent for crack of closed door
Closed door rating conditions
PressureDifferencedpCloRat4Pressure drop at rating condition of closed door [Pa]
RealCDCloRat1Discharge coefficient at rating conditions of closed door
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
PressureDifferencedp_turbulent0.01Pressure difference where laminar and turbulent flow relation coincide [Pa]

Connectors

TypeNameDescription
FluidPort_aport_a1Fluid connector a1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_bport_b1Fluid connector b1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_aport_a2Fluid connector a2 (positive design flow direction is from port_a2 to port_b2)
FluidPort_bport_b2Fluid connector b2 (positive design flow direction is from port_a2 to port_b2)
input RealInputyOpening signal, 0=closed, 1=open [1]

Modelica definition

model DoorOperable "Door model for bi-directional air flow between rooms that can be open or closed" extends Buildings.Airflow.Multizone.BaseClasses.Door( final vAB = VAB_flow/A, final vBA = VBA_flow/A); parameter Real CDOpe=0.65 "Discharge coefficient of open door"; parameter Real mOpe = 0.5 "Flow exponent for door of open door"; parameter Modelica.SIunits.Area LClo(min=0) "Effective leakage area of closed door"; parameter Real mClo= 0.65 "Flow exponent for crack of closed door"; parameter Modelica.SIunits.PressureDifference dpCloRat(min=0, displayUnit="Pa") = 4 "Pressure drop at rating condition of closed door"; parameter Real CDCloRat(min=0, max=1)=1 "Discharge coefficient at rating conditions of closed door"; Modelica.Blocks.Interfaces.RealInput y(min=0, max=1, unit="1") "Opening signal, 0=closed, 1=open"; protected parameter Real gamma(min=1) = 1.5 "Normalized flow rate where dphi(0)/dpi intersects phi(1)"; parameter Real[2] a = {gamma, gamma} "Polynomial coefficient for regularized implementation of flow resistance"; parameter Real b[2] = {1/8*m^2 - 3*gamma - 3/2*m + 35.0/8 for m in {mOpe, mClo}} "Polynomial coefficient for regularized implementation of flow resistance"; parameter Real c[2] = {-1/4*m^2 + 3*gamma + 5/2*m - 21.0/4 for m in {mOpe, mClo}} "Polynomial coefficient for regularized implementation of flow resistance"; parameter Real d[2] = {1/8*m^2 - gamma - m + 15.0/8 for m in {mOpe, mClo}} "Polynomial coefficient for regularized implementation of flow resistance"; parameter Modelica.SIunits.Area AClo = LClo * dpCloRat^(0.5-mClo) "Closed area"; parameter Real kVal[2]={ CDOpe *AOpe*sqrt(2/rho_default), CDCloRat*AClo*sqrt(2/rho_default)} "Flow coefficient, k = V_flow/ dp^m"; parameter Real kT = rho_default * CDOpe * AOpe/3 * sqrt(Modelica.Constants.g_n /(Medium.T_default*conTP) * hOpe) "Constant coefficient for buoyancy driven air flow rate"; parameter Modelica.SIunits.MassFlowRate m_flow_turbulent= kVal[1] * rho_default * sqrt(dp_turbulent) "Mass flow rate where regularization to laminar flow occurs for temperature-driven flow"; Modelica.SIunits.VolumeFlowRate VABpOpeClo_flow[2](each nominal=0.001) "Volume flow rate from A to B if positive due to static pressure difference"; Modelica.SIunits.VolumeFlowRate VABp_flow(nominal=0.001) "Volume flow rate from A to B if positive due to static pressure difference"; Modelica.SIunits.Area A "Current opening area"; equation // Air flow rate due to static pressure difference VABpOpeClo_flow = Buildings.Airflow.Multizone.BaseClasses.powerLawFixedM( k=kVal, dp=port_a1.p-port_a2.p, m={mOpe, mClo}, a=a, b=b, c=c, d=d, dp_turbulent=dp_turbulent); VABp_flow = y*VABpOpeClo_flow[1] + (1-y)*VABpOpeClo_flow[2]; A = y*AOpe + (1-y)*AClo; // Air flow rate due to buoyancy // Because powerLawFixedM requires as an input a pressure difference pa-pb, // we convert Ta-Tb by multiplying it with rho*R, and we divide // above the constant expression by (rho*R)^m on the right hand-side of kT. mABt_flow = y*Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_dp( k=kT, dp=conTP*(Medium.temperature(state_a1_inflow)-Medium.temperature(state_a2_inflow)), m_flow_turbulent=m_flow_turbulent); end DoorOperable;

Buildings.Airflow.Multizone.EffectiveAirLeakageArea Buildings.Airflow.Multizone.EffectiveAirLeakageArea

Effective air leakage area

Buildings.Airflow.Multizone.EffectiveAirLeakageArea

Information

This model describes the one-directional pressure driven air flow through a crack-like opening, using the equation

V̇ = k Δpm,

where is the volume flow rate, k is a flow coefficient and m is the flow exponent. The flow coefficient is

k = L CD,Rat ΔpRat(0.5-m) (2/ρ0)0.5,

where L is the effective air leakage area, CD,Rat is the discharge coefficient at the reference condition, ΔpRat is the pressure drop at the rating condition, and ρ0 is the mass density at the medium default pressure, temperature and humidity.

The effective air leakage area L can be obtained, for example, from the ASHRAE fundamentals (ASHRAE, 1997, p. 25.18). In the ASHRAE fundamentals, the effective air leakage area is based on a reference pressure difference of ΔpRat = 4 Pa and a discharge coefficient of CD,Rat = 1. A similar model is also used in the CONTAM software (Dols and Walton, 2002). Dols and Walton (2002) recommend to use for the flow exponent m=0.6 to m=0.7 if the flow exponent is not reported with the test results.

References

Extends from Buildings.Airflow.Multizone.BaseClasses.PowerLawResistance (Flow resistance that uses the power law).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Realm0.65Flow exponent, m=0.5 for turbulent, m=1 for laminar
AreaL Effective leakage area [m2]
Rating conditions
PressureDifferencedpRat4Pressure drop [Pa]
RealCDRat1Discharge coefficient
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
BooleanforceErrorControlOnFlowtrueFlag to force error control on m_flow. Set to true if interested in flow rate
BooleanuseDefaultPropertiestrueSet to false to use density and viscosity based on actual medium state, rather than using default values
PressureDifferencedp_turbulent0.1Pressure difference where laminar and turbulent flow relation coincide. Recommended = 0.1 [Pa]

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)

Modelica definition

model EffectiveAirLeakageArea "Effective air leakage area" extends Buildings.Airflow.Multizone.BaseClasses.PowerLawResistance( m=0.65, final k=L * CDRat * sqrt(2.0/rho_default) * dpRat^(0.5-m)); parameter Modelica.SIunits.PressureDifference dpRat( min=0, displayUnit="Pa") = 4 "Pressure drop"; parameter Real CDRat( min=0, max=1) = 1 "Discharge coefficient"; parameter Modelica.SIunits.Area L(min=0) "Effective leakage area"; equation v = V_flow/L; end EffectiveAirLeakageArea;

Buildings.Airflow.Multizone.MediumColumn Buildings.Airflow.Multizone.MediumColumn

Vertical shaft with no friction and no storage of heat and mass

Buildings.Airflow.Multizone.MediumColumn

Information

This model describes the pressure difference of a vertical medium column. It can be used to model the pressure difference caused by stack effect.

Typical use and important parameters

The model can be used with the following three configurations, which are controlled by the setting of the parameter densitySelection:

The settings top and bottom should be used when rooms or different floors of a building are connected since multizone airflow models assume that each floor is completely mixed. For these two seetings, this model will compute the pressure between the center of the room and an opening that is at height h relative to the center of the room. The setting actual may be used to model a chimney in which a column of air will change its density based on the flow direction.

In this model, the parameter h must always be positive, and the port port_a must be at the top of the column.

Dynamics

For a dynamic model, use Buildings.Airflow.Multizone.MediumColumnDynamic instead of this model.

Parameters

TypeNameDefaultDescription
replaceable package MediumModelica.Media.Interfaces.Pa...Medium in the component
Lengthh3Height of shaft [m]
densitySelectiondensitySelection Select how to pick density

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

model MediumColumn "Vertical shaft with no friction and no storage of heat and mass" replaceable package Medium = Modelica.Media.Interfaces.PartialMedium "Medium in the component"; parameter Modelica.SIunits.Length h(min=0) = 3 "Height of shaft"; parameter Buildings.Airflow.Multizone.Types.densitySelection densitySelection "Select how to pick density"; Modelica.Fluid.Interfaces.FluidPort_a port_a( redeclare package Medium = Medium, p(start=Medium.p_default)) "Fluid connector a (positive design flow direction is from port_a to port_b)"; Modelica.Fluid.Interfaces.FluidPort_b port_b( redeclare package Medium = Medium, p(start=Medium.p_default)) "Fluid connector b (positive design flow direction is from port_a to port_b)"; Modelica.SIunits.VolumeFlowRate V_flow "Volume flow rate at inflowing port (positive when flow from port_a to port_b)"; Modelica.SIunits.MassFlowRate m_flow "Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction)"; Modelica.SIunits.PressureDifference dp(displayUnit="Pa") "Pressure difference between port_a and port_b"; Modelica.SIunits.Density rho "Density in medium column"; protected Medium.ThermodynamicState sta_a=Medium.setState_phX( port_a.p, actualStream(port_a.h_outflow), actualStream(port_a.Xi_outflow)) "Medium properties in port_a"; Medium.MassFraction Xi[Medium.nXi] "Mass fraction used to compute density"; initial equation // The next assert tests for all allowed values of the enumeration. // Testing against densitySelection > 0 gives an error in OpenModelica as enumerations start with 1. assert(densitySelection == Buildings.Airflow.Multizone.Types.densitySelection.fromTop or densitySelection == Buildings.Airflow.Multizone.Types.densitySelection.fromBottom or densitySelection == Buildings.Airflow.Multizone.Types.densitySelection.actual, "You need to set the parameter \"densitySelection\" for the model MediumColumn."); equation // Design direction of mass flow rate m_flow = port_a.m_flow; // Pressure difference between ports // Xi is computed first as it is used in two expression, and in one // of them only one component is used. // We test for Medium.nXi == 0 as Modelica.Media.Air.SimpleAir has no // moisture and hence Xi[1] is an illegal statement. // We first compute temperature and then invoke a density function that // takes temperature as an argument. Simply calling a density function // of a medium that takes enthalpy as an argument would be dangerous // as different media can have different datum for the enthalpy. if (densitySelection == Buildings.Airflow.Multizone.Types.densitySelection.fromTop) then Xi = inStream(port_a.Xi_outflow); rho = Buildings.Utilities.Psychrometrics.Functions.density_pTX( p=Medium.p_default, T=Medium.temperature(Medium.setState_phX(port_a.p, inStream(port_a.h_outflow), Xi)), X_w=if Medium.nXi == 0 then 0 else Xi[1]); elseif (densitySelection == Buildings.Airflow.Multizone.Types.densitySelection.fromBottom) then Xi = inStream(port_b.Xi_outflow); rho = Buildings.Utilities.Psychrometrics.Functions.density_pTX( p=Medium.p_default, T=Medium.temperature(Medium.setState_phX(port_b.p, inStream(port_b.h_outflow), Xi)), X_w=if Medium.nXi == 0 then 0 else Xi[1]); else Xi = actualStream(port_a.Xi_outflow); rho = Buildings.Utilities.Psychrometrics.Functions.density_pTX( p=Medium.p_default, T=Medium.temperature(Medium.setState_phX(port_a.p, actualStream(port_a.h_outflow), Xi)), X_w=if Medium.nXi == 0 then 0 else Xi[1]); end if; V_flow = m_flow/Medium.density(sta_a); dp = port_a.p - port_b.p; dp = -h*rho*Modelica.Constants.g_n; // Isenthalpic state transformation (no storage and no loss of energy) port_a.h_outflow = inStream(port_b.h_outflow); port_b.h_outflow = inStream(port_a.h_outflow); // Mass balance (no storage) port_a.m_flow + port_b.m_flow = 0; // Transport of substances port_a.Xi_outflow = inStream(port_b.Xi_outflow); port_b.Xi_outflow = inStream(port_a.Xi_outflow); port_a.C_outflow = inStream(port_b.C_outflow); port_b.C_outflow = inStream(port_a.C_outflow); end MediumColumn;

Buildings.Airflow.Multizone.MediumColumnDynamic Buildings.Airflow.Multizone.MediumColumnDynamic

Vertical shaft with no friction and storage of heat and mass

Buildings.Airflow.Multizone.MediumColumnDynamic

Information

This model contains a completely mixed fluid volume and models that take into account the pressure difference of a medium column that is at the same temperature as the fluid volume. It can be used to model the pressure difference caused by a stack effect.

Typical use and important parameters

Set the parameter use_HeatTransfer=true to expose a heatPort. This heatPort can be used to add or subtract heat from the volume. This allows, for example, to use this model in conjunction with a model for heat transfer through walls to model a solar chimney that stores heat.

Dynamics

For a steady-state model, use Buildings.Airflow.Multizone.MediumColumn instead of this model.

In this model, the parameter h must always be positive, and the port port_a must be at the top of the column.

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

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Lengthh3Height of shaft [m]
VolumeV Volume in medium shaft [m3]
Dynamics
Equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Type of energy balance: dynamic (3 initialization options) or steady state
DynamicsmassDynamicsenergyDynamicsType 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.)
Assumptions
Heat transfer
Booleanuse_HeatTransferfalse= true to use the HeatTransfer model
replaceable model HeatTransferModelica.Fluid.Vessels.BaseC...Wall heat transfer

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)
HeatPort_aheatPortHeat port to exchange energy with the fluid volume
Assumptions
Heat transfer
replaceable model HeatTransferWall heat transfer

Modelica definition

model MediumColumnDynamic "Vertical shaft with no friction and storage of heat and mass" extends Buildings.Fluid.Interfaces.LumpedVolumeDeclarations; replaceable package Medium = Modelica.Media.Interfaces.PartialMedium "Medium in the component"; parameter Modelica.SIunits.Length h(min=0) = 3 "Height of shaft"; Modelica.Fluid.Interfaces.FluidPort_a port_a( redeclare final package Medium = Medium, p(start=Medium.p_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, p(start=Medium.p_default)) "Fluid connector b (positive design flow direction is from port_a to port_b)"; parameter Modelica.SIunits.Volume V "Volume in medium shaft"; // Heat transfer through boundary parameter Boolean use_HeatTransfer = false "= true to use the HeatTransfer model"; // m_flow_nominal is not used by vol, since this component // can only be configured as a dynamic model. // Therefore we set it to about 10 air changes per hour Buildings.Fluid.MixingVolumes.MixingVolume vol( final nPorts=2, redeclare final package Medium = Medium, final m_flow_nominal=V/360, final V=V, 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) "Air volume in the shaft"; MediumColumn colTop( redeclare final package Medium = Medium, final densitySelection=Buildings.Airflow.Multizone.Types.densitySelection.fromBottom, h=h/2) "Medium column that connects to top port"; MediumColumn colBot( redeclare final package Medium = Medium, final densitySelection=Buildings.Airflow.Multizone.Types.densitySelection.fromTop, h=h/2) "Medium colum that connects to bottom port"; replaceable model HeatTransfer = Modelica.Fluid.Vessels.BaseClasses.HeatTransfer.IdealHeatTransfer constrainedby Modelica.Fluid.Vessels.BaseClasses.HeatTransfer.PartialVesselHeatTransfer "Wall heat transfer"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a heatPort if use_HeatTransfer "Heat port to exchange energy with the fluid volume"; equation connect(colBot.port_a, vol.ports[1]); connect(vol.ports[2], colTop.port_b); connect(colTop.port_a, port_a); connect(colBot.port_b, port_b); connect(heatPort, vol.heatPort); end MediumColumnDynamic;

Buildings.Airflow.Multizone.Orifice Buildings.Airflow.Multizone.Orifice

Orifice

Buildings.Airflow.Multizone.Orifice

Information

This model describes the mass flow rate and pressure difference relation of an orifice in the form

V̇ = k Δpm,

where is the volume flow rate, k is a flow coefficient and m is the flow exponent. The flow coefficient is

k = CD A (2/ρ0)0.5,

where CD is the discharge coefficient, A is the cross section area and ρ0 is the mass density at the medium default pressure, temperature and humidity.

For turbulent flow, set m=1/2 and for laminar flow, set m=1. Large openings are characterized by values close to 0.5, while values near 0.65 have been found for small crack-like openings (Dols and Walton, 2002).

References

Extends from Buildings.Airflow.Multizone.BaseClasses.PowerLawResistance (Flow resistance that uses the power law).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Realm0.5Flow exponent, m=0.5 for turbulent, m=1 for laminar
Orifice characteristics
AreaA Area of orifice [m2]
RealCD0.65Discharge coefficient
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
BooleanforceErrorControlOnFlowtrueFlag to force error control on m_flow. Set to true if interested in flow rate
BooleanuseDefaultPropertiestrueSet to false to use density and viscosity based on actual medium state, rather than using default values
PressureDifferencedp_turbulent0.1Pressure difference where laminar and turbulent flow relation coincide. Recommended = 0.1 [Pa]

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)

Modelica definition

model Orifice "Orifice" extends Buildings.Airflow.Multizone.BaseClasses.PowerLawResistance( m=0.5, k=CD*A*sqrt(2.0/rho_default)); parameter Modelica.SIunits.Area A "Area of orifice"; parameter Real CD=0.65 "Discharge coefficient"; equation v = V_flow/A; end Orifice;

Buildings.Airflow.Multizone.ZonalFlow_ACS Buildings.Airflow.Multizone.ZonalFlow_ACS

Zonal flow with input air change per second

Buildings.Airflow.Multizone.ZonalFlow_ACS

Information

This model computes the air exchange between volumes.

Input is the air change per seconds. The volume flow rate is computed as

  V_flow = ACS * V

where ACS is an input and the volume V is a parameter.

Extends from Buildings.Airflow.Multizone.BaseClasses.ZonalFlow (Flow across zonal boundaries of a room).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMedium 
BooleanuseDefaultPropertiesfalseSet to true to use constant density
VolumeV Volume of room [m3]
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed

Connectors

TypeNameDescription
FluidPort_aport_a1Fluid connector a1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_bport_b1Fluid connector b1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_aport_a2Fluid connector a2 (positive design flow direction is from port_a2 to port_b2)
FluidPort_bport_b2Fluid connector b2 (positive design flow direction is from port_a2 to port_b2)
input RealInputACSAir change per seconds, relative to the smaller of the two volumes

Modelica definition

model ZonalFlow_ACS "Zonal flow with input air change per second" extends Buildings.Airflow.Multizone.BaseClasses.ZonalFlow; parameter Boolean useDefaultProperties = false "Set to true to use constant density"; parameter Modelica.SIunits.Volume V "Volume of room"; Modelica.Blocks.Interfaces.RealInput ACS "Air change per seconds, relative to the smaller of the two volumes"; protected Modelica.SIunits.VolumeFlowRate V_flow "Volume flow rate at standard pressure"; Modelica.SIunits.MassFlowRate m_flow "Mass flow rate"; parameter Medium.ThermodynamicState sta_default = Medium.setState_pTX(T=Medium.T_default, p=Medium.p_default, X=Medium.X_default); parameter Modelica.SIunits.Density rho_default=Medium.density(sta_default) "Density, used to compute fluid volume"; Medium.ThermodynamicState sta_a1_inflow= Medium1.setState_phX(port_a1.p, port_b1.h_outflow, port_b1.Xi_outflow) "Medium properties in port_a1"; Medium.ThermodynamicState sta_a2_inflow= Medium1.setState_phX(port_a2.p, port_b2.h_outflow, port_b2.Xi_outflow) "Medium properties in port_a2"; equation when useDefaultProperties and initial() then assert( abs(1-rho_default/((Medium.density(sta_a1_inflow) + Medium.density(sta_a2_inflow))/2)) < 0.2, "Wrong density. Densities need to match."); // The next three lines have been removed as they cause // a compilation error in OpenModelica. // + "\n Medium.density(sta_a1) = " + String(Medium.density(sta_a1_inflow)) // + "\n Medium.density(sta_a2) = " + String(Medium.density(sta_a2_inflow)) // + "\n rho_nominal = " + String(rho_default)); end when; V_flow = V * ACS; m_flow / V_flow = if useDefaultProperties then rho_default else (Medium.density(sta_a1_inflow) + Medium.density(sta_a2_inflow))/2; // assign variable in base class port_a1.m_flow = m_flow; port_a2.m_flow = m_flow; end ZonalFlow_ACS;

Buildings.Airflow.Multizone.ZonalFlow_m_flow Buildings.Airflow.Multizone.ZonalFlow_m_flow

Zonal flow with input air change per second

Buildings.Airflow.Multizone.ZonalFlow_m_flow

Information

This model computes the air exchange between volumes.

Input is the mass flow rate from port_a1 to port_b1 and from port_a2 to port_b2.

Extends from Buildings.Airflow.Multizone.BaseClasses.ZonalFlow (Flow across zonal boundaries of a room).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMedium 
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed

Connectors

TypeNameDescription
FluidPort_aport_a1Fluid connector a1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_bport_b1Fluid connector b1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_aport_a2Fluid connector a2 (positive design flow direction is from port_a2 to port_b2)
FluidPort_bport_b2Fluid connector b2 (positive design flow direction is from port_a2 to port_b2)
input RealInputmAB_flowMass flow rate from A to B, positive if flow from port_a1 to port_b1
input RealInputmBA_flowMass flow rate from B to A, positive if flow from port_a2 to port_b2

Modelica definition

model ZonalFlow_m_flow "Zonal flow with input air change per second" extends Buildings.Airflow.Multizone.BaseClasses.ZonalFlow; Modelica.Blocks.Interfaces.RealInput mAB_flow "Mass flow rate from A to B, positive if flow from port_a1 to port_b1"; Modelica.Blocks.Interfaces.RealInput mBA_flow "Mass flow rate from B to A, positive if flow from port_a2 to port_b2"; equation port_a1.m_flow = mAB_flow; port_a2.m_flow = mBA_flow; end ZonalFlow_m_flow;