This package contains base classes that are used to construct the models in Buildings.Airflow.Multizone.
Extends from Modelica.Icons.BasesPackage (Icon for packages containing base classes).Name | Description |
---|---|
DoorDiscretized | Door model using discretization along height coordinate |
PowerLawResistance | Flow resistance that uses the power law |
TwoWayFlowElement | Flow resistance that uses the power law |
TwoWayFlowElementBuoyancy | Flow resistance that uses the power law |
ZonalFlow | Flow across zonal boundaries of a room |
ErrorControl | Interface that defines parameters for error control |
powerLaw | Power law used in orifice equations |
powerLawFixedM | Power law used in orifice equations when m is constant |
windPressureLowRise | Wind pressure coefficient for low-rise buildings |
Examples | Collection of models that illustrate model use and test models |
This is a partial model for the bi-directional air flow through a 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 compartment area dA
is a variable, which allows
using the model for a door that can be open or closed.
Type | Name | Default | Description |
---|---|---|---|
Boolean | forceErrorControlOnFlow | true | Flag to force error control on m_flow. Set to true if interested in flow rate |
replaceable package Medium | PartialMedium | ||
Velocity | vZer | 0.001 | Minimum velocity to prevent zero flow. Recommended: 0.001 [m/s] |
Integer | nCom | 10 | Number of compartments for the discretization |
Pressure | dp_turbulent | 0.01 | Pressure difference where laminar and turbulent flow relation coincide. Recommended: 0.01 [Pa] |
Initialization | |||
MassFlowRate | m1_flow.start | 0 | Mass flow rate from port_a1 to port_b1 (m1_flow > 0 is design flow direction) [kg/s] |
Pressure | dp1.start | 0 | Pressure difference between port_a1 and port_b1 [Pa] |
MassFlowRate | m2_flow.start | 0 | Mass flow rate from port_a2 to port_b2 (m2_flow > 0 is design flow direction) [kg/s] |
Pressure | dp2.start | 0 | Pressure difference between port_a2 and port_b2 [Pa] |
Density | rho_a1_inflow.start | 1.2 | Density of air flowing in from port_a1 [kg/m3] |
Density | rho_a2_inflow.start | 1.2 | Density of air flowing in from port_a2 [kg/m3] |
Geometry | |||
Length | wOpe | 0.9 | Width of opening [m] |
Length | hOpe | 2.1 | Height of opening [m] |
Length | hA | 2.7/2 | Height of reference pressure zone A [m] |
Length | hB | 2.7/2 | Height of reference pressure zone B [m] |
Orifice characteristics | |||
Real | CD | 0.65 | Discharge coefficient |
Advanced | |||
Initialization | |||
SpecificEnthalpy | h_outflow_a1_start | Medium1.h_default | Start value for enthalpy flowing out of port a1 [J/kg] |
SpecificEnthalpy | h_outflow_b1_start | Medium1.h_default | Start value for enthalpy flowing out of port b1 [J/kg] |
SpecificEnthalpy | h_outflow_a2_start | Medium2.h_default | Start value for enthalpy flowing out of port a2 [J/kg] |
SpecificEnthalpy | h_outflow_b2_start | Medium2.h_default | Start value for enthalpy flowing out of port b2 [J/kg] |
MassFlowRate | m1_flow_small | 1E-4*abs(m1_flow_nominal) | Small mass flow rate for regularization of zero flow [kg/s] |
MassFlowRate | m2_flow_small | 1E-4*abs(m2_flow_nominal) | Small mass flow rate for regularization of zero flow [kg/s] |
Boolean | homotopyInitialization | true | = true, use homotopy method |
Diagnostics | |||
Boolean | show_T | false | = true, if actual temperature at port is computed |
Type | Name | Description |
---|---|---|
FluidPort_a | port_a1 | Fluid connector a1 (positive design flow direction is from port_a1 to port_b1) |
FluidPort_b | port_b1 | Fluid connector b1 (positive design flow direction is from port_a1 to port_b1) |
FluidPort_a | port_a2 | Fluid connector a2 (positive design flow direction is from port_a2 to port_b2) |
FluidPort_b | port_b2 | Fluid connector b2 (positive design flow direction is from port_a2 to port_b2) |
partial model DoorDiscretized "Door model using discretization along height coordinate" extends Buildings.Airflow.Multizone.BaseClasses.TwoWayFlowElementBuoyancy; parameter Integer nCom=10 "Number of compartments for the discretization"; parameter Modelica.SIunits.Pressure dp_turbulent(min=0) = 0.01 "Pressure difference where laminar and turbulent flow relation coincide. Recommended: 0.01"; parameter Real CD=0.65 "|Orifice characteristics|Discharge coefficient"; Modelica.SIunits.Pressure dpAB[nCom](each nominal=1) "Pressure difference between compartments"; Modelica.SIunits.Velocity v[nCom](each nominal=0.01) "Velocity in compartment from A to B"; Modelica.SIunits.Velocity vTop "Velocity at top of opening from A to B"; Modelica.SIunits.Velocity vBot "Velocity at bottom of opening from A to B"; protected parameter Modelica.SIunits.Length dh=hOpe/nCom "Height of each compartment"; Modelica.SIunits.AbsolutePressure pA[nCom](each nominal=101325) "Pressure in compartments of room A"; Modelica.SIunits.AbsolutePressure pB[nCom](each nominal=101325) "Pressure in compartments of room B"; Modelica.SIunits.VolumeFlowRate dV_flow[nCom] "Volume flow rate through compartment from A to B"; Modelica.SIunits.VolumeFlowRate dVAB_flow[nCom] "Volume flow rate through compartment from A to B if positive"; Modelica.SIunits.VolumeFlowRate dVBA_flow[nCom] "Volume flow rate through compartment from B to A if positive"; Real m(min=0.5, max=1) "Flow exponent, m=0.5 for turbulent, m=1 for laminar"; Real kVal "Flow coefficient for each compartment, k = V_flow/ dp^m"; Modelica.SIunits.Area dA "Compartment area"; 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"; equation dA = A/nCom; for i in 1:nCom loop // pressure drop in each compartment pA[i] = port_a1.p + rho_a1_inflow*Modelica.Constants.g_n*(hA - (i - 0.5)*dh); pB[i] = port_a2.p + rho_a2_inflow*Modelica.Constants.g_n*(hB - (i - 0.5)*dh); dpAB[i] = pA[i] - pB[i]; v[i] = dV_flow[i]/dA; // assignment of net volume flows dVAB_flow[i] = dV_flow[i]* Buildings.Utilities.Math.Functions.smoothHeaviside(x=dV_flow[i], delta= VZer_flow/nCom) + VZer_flow/nCom; dVBA_flow[i] = -dV_flow[i] + dVAB_flow[i] + 2*VZer_flow/nCom; end for; // add positive and negative flows VAB_flow = ones(nCom)*dVAB_flow; VBA_flow = ones(nCom)*dVBA_flow; vTop = v[nCom]; vBot = v[1];end DoorDiscretized;
This model describes the mass flow rate and pressure difference relation of an orifice in the form
V_flow = k * dp^m,where
k
is a variable and
m
a parameter.
For turbulent flow, set m=1/2
and
for laminar flow, set m=1
.
The model is used as a base for the interzonal air flow models.
Extends from Buildings.Fluid.Interfaces.PartialTwoPortInterface (Partial model transporting fluid between two ports without storing mass or energy), Buildings.Airflow.Multizone.BaseClasses.ErrorControl (Interface that defines parameters for error control).
Type | Name | Default | Description |
---|---|---|---|
replaceable package Medium | PartialMedium | Medium in the component | |
Boolean | forceErrorControlOnFlow | true | Flag to force error control on m_flow. Set to true if interested in flow rate |
Real | m | Flow exponent, m=0.5 for turbulent, m=1 for laminar | |
Boolean | useDefaultProperties | true | Set to false to use density and viscosity based on actual medium state, rather than using default values |
Pressure | dp_turbulent | 0.1 | Pressure difference where laminar and turbulent flow relation coincide. Recommended = 0.1 [Pa] |
Length | lWet | sqrt(A) | Wetted perimeter used for Reynolds number calculation [m] |
Nominal condition | |||
MassFlowRate | m_flow_nominal | rho_default*k*dp_turbulent | Nominal mass flow rate [kg/s] |
Initialization | |||
MassFlowRate | m_flow.start | 0 | Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction) [kg/s] |
Pressure | dp.start | 0 | Pressure difference between port_a and port_b [Pa] |
Orifice characteristics | |||
Area | A | Area of orifice [m2] | |
Assumptions | |||
Boolean | allowFlowReversal | system.allowFlowReversal | = true to allow flow reversal, false restricts to design direction (port_a -> port_b) |
Advanced | |||
MassFlowRate | m_flow_small | 1E-4*abs(m_flow_nominal) | Small mass flow rate for regularization of zero flow [kg/s] |
Boolean | homotopyInitialization | true | = true, use homotopy method |
Diagnostics | |||
Boolean | show_T | false | = true, if actual temperature at port is computed |
Type | Name | Description |
---|---|---|
FluidPort_a | port_a | Fluid connector a (positive design flow direction is from port_a to port_b) |
FluidPort_b | port_b | Fluid connector b (positive design flow direction is from port_a to port_b) |
partial model PowerLawResistance "Flow resistance that uses the power law" extends Buildings.Fluid.Interfaces.PartialTwoPortInterface( final m_flow_nominal=rho_default*k*dp_turbulent, final showDesignFlowDirection=false); extends Buildings.Airflow.Multizone.BaseClasses.ErrorControl; parameter Modelica.SIunits.Area A "|Orifice characteristics|Area of orifice"; parameter Real m(min=0.5, max=1) "Flow exponent, m=0.5 for turbulent, m=1 for laminar"; parameter Boolean useDefaultProperties=true "Set to false to use density and viscosity based on actual medium state, rather than using default values"; parameter Modelica.SIunits.Pressure dp_turbulent(min=0, displayUnit="Pa") = 0.1 "Pressure difference where laminar and turbulent flow relation coincide. Recommended = 0.1"; parameter Modelica.SIunits.Length lWet=sqrt(A) "Wetted perimeter used for Reynolds number calculation"; Modelica.SIunits.VolumeFlowRate V_flow "Volume flow rate through the component"; Modelica.SIunits.Velocity v(nominal=1) "Average velocity"; Modelica.SIunits.Density rho "Fluid density at port_a"; Real Re "Reynolds number"; protected constant Real gamma(min=1) = 1.5 "Normalized flow rate where dphi(0)/dpi intersects phi(1)"; parameter Real k "Flow coefficient, k = V_flow/ dp^m"; parameter Medium.ThermodynamicState sta_default=Medium.setState_pTX( T=Medium.T_default, p=Medium.p_default, X=Medium.X_default) "State of the medium at the medium default properties"; parameter Modelica.SIunits.Density rho_default=Medium.density(sta_default) "Density at the medium default properties"; parameter Modelica.SIunits.DynamicViscosity dynVis_default= Medium.dynamicViscosity(sta_default) "Dynamic viscosity at the medium default properties"; parameter 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"; Medium.ThermodynamicState sta "State of the medium in the component"; Modelica.SIunits.DynamicViscosity dynVis "Dynamic viscosity"; Modelica.SIunits.Mass mExc "Air mass exchanged (for purpose of error control only)"; initial equation mExc=0; equation if forceErrorControlOnFlow then der(mExc) = port_a.m_flow; else der(mExc) = 0; end if; if useDefaultProperties then sta = sta_default; rho = rho_default; dynVis = dynVis_default; else sta = if homotopyInitialization then Medium.setState_phX(port_a.p, homotopy(actual=actualStream(port_a.h_outflow), simplified=inStream(port_a.h_outflow)), homotopy(actual=actualStream(port_a.Xi_outflow), simplified=inStream(port_a.Xi_outflow))) else Medium.setState_phX(port_a.p, actualStream(port_a.h_outflow), actualStream(port_a.Xi_outflow)); rho = Medium.density(sta); dynVis = Medium.dynamicViscosity(sta); end if; V_flow = Buildings.Airflow.Multizone.BaseClasses.powerLawFixedM( k=k, dp=dp, m=m, a=a, b=b, c=c, d=d, dp_turbulent=dp_turbulent); port_a.m_flow = rho*V_flow; v = V_flow/A; Re = v*lWet*rho/dynVis; // 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 PowerLawResistance;
This is a partial model for models that describe the bi-directional air flow through large openings.
Models that extend this model need to compute
mAB_flow
and mBA_flow
,
or alternatively VAB_flow
and VBA_flow
,
and the face area area
.
The face area is a variable to allow this partial model to be used
for doors that can be open or closed as a function of an input signal.
Type | Name | Default | Description |
---|---|---|---|
replaceable package Medium1 | PartialMedium | Medium 1 in the component | |
replaceable package Medium2 | PartialMedium | Medium 2 in the component | |
Boolean | forceErrorControlOnFlow | true | Flag to force error control on m_flow. Set to true if interested in flow rate |
replaceable package Medium | Modelica.Media.Interfaces.Pa... | ||
Velocity | vZer | 0.001 | Minimum velocity to prevent zero flow. Recommended: 0.001 [m/s] |
Nominal condition | |||
MassFlowRate | m1_flow_nominal | 10/3600*1.2 | Nominal mass flow rate [kg/s] |
MassFlowRate | m2_flow_nominal | m1_flow_nominal | Nominal mass flow rate [kg/s] |
Initialization | |||
MassFlowRate | m1_flow.start | 0 | Mass flow rate from port_a1 to port_b1 (m1_flow > 0 is design flow direction) [kg/s] |
Pressure | dp1.start | 0 | Pressure difference between port_a1 and port_b1 [Pa] |
MassFlowRate | m2_flow.start | 0 | Mass flow rate from port_a2 to port_b2 (m2_flow > 0 is design flow direction) [kg/s] |
Pressure | dp2.start | 0 | Pressure difference between port_a2 and port_b2 [Pa] |
Assumptions | |||
Boolean | allowFlowReversal1 | false | = true to allow flow reversal in medium 1, false restricts to design direction (port_a -> port_b) |
Boolean | allowFlowReversal2 | false | = true to allow flow reversal in medium 2, false restricts to design direction (port_a -> port_b) |
Advanced | |||
Initialization | |||
SpecificEnthalpy | h_outflow_a1_start | Medium1.h_default | Start value for enthalpy flowing out of port a1 [J/kg] |
SpecificEnthalpy | h_outflow_b1_start | Medium1.h_default | Start value for enthalpy flowing out of port b1 [J/kg] |
SpecificEnthalpy | h_outflow_a2_start | Medium2.h_default | Start value for enthalpy flowing out of port a2 [J/kg] |
SpecificEnthalpy | h_outflow_b2_start | Medium2.h_default | Start value for enthalpy flowing out of port b2 [J/kg] |
MassFlowRate | m1_flow_small | 1E-4*abs(m1_flow_nominal) | Small mass flow rate for regularization of zero flow [kg/s] |
MassFlowRate | m2_flow_small | 1E-4*abs(m2_flow_nominal) | Small mass flow rate for regularization of zero flow [kg/s] |
Boolean | homotopyInitialization | true | = true, use homotopy method |
Diagnostics | |||
Boolean | show_T | false | = true, if actual temperature at port is computed |
Type | Name | Description |
---|---|---|
replaceable package Medium1 | Medium 1 in the component | |
replaceable package Medium2 | Medium 2 in the component | |
FluidPort_a | port_a1 | Fluid connector a1 (positive design flow direction is from port_a1 to port_b1) |
FluidPort_b | port_b1 | Fluid connector b1 (positive design flow direction is from port_a1 to port_b1) |
FluidPort_a | port_a2 | Fluid connector a2 (positive design flow direction is from port_a2 to port_b2) |
FluidPort_b | port_b2 | Fluid connector b2 (positive design flow direction is from port_a2 to port_b2) |
replaceable package Medium |
partial model TwoWayFlowElement "Flow resistance that uses the power law" extends Buildings.Fluid.Interfaces.PartialFourPortInterface( redeclare final package Medium1 = Medium, redeclare final package Medium2 = Medium, final allowFlowReversal1=false, final allowFlowReversal2=false, final m1_flow_nominal=10/3600*1.2, final m2_flow_nominal=m1_flow_nominal); extends Buildings.Airflow.Multizone.BaseClasses.ErrorControl; replaceable package Medium = Modelica.Media.Interfaces.PartialMedium; Modelica.SIunits.VolumeFlowRate VAB_flow(nominal=0.001) "Volume flow rate from A to B if positive"; Modelica.SIunits.VolumeFlowRate VBA_flow(nominal=0.001) "Volume flow rate from B to A if positive"; Modelica.SIunits.MassFlowRate mAB_flow(nominal=0.001) "Mass flow rate from A to B if positive"; Modelica.SIunits.MassFlowRate mBA_flow(nominal=0.001) "Mass flow rate from B to A if positive"; Modelica.SIunits.Velocity vAB(nominal=0.01) "Average velocity from A to B"; Modelica.SIunits.Velocity vBA(nominal=0.01) "Average velocity from B to A"; Modelica.SIunits.Density rho_a1_inflow(start=1.2) "Density of air flowing in from port_a1"; Modelica.SIunits.Density rho_a2_inflow(start=1.2) "Density of air flowing in from port_a2"; Modelica.SIunits.Area A "Face area"; parameter Modelica.SIunits.Velocity vZer=0.001 "Minimum velocity to prevent zero flow. Recommended: 0.001"; protected Modelica.SIunits.VolumeFlowRate VZer_flow(fixed=false) "Minimum net volume flow rate to prevent zero flow"; protected Modelica.SIunits.Mass mExcAB(start=0) "Air mass exchanged (for purpose of error control only)"; Modelica.SIunits.Mass mExcBA(start=0) "Air mass exchanged (for purpose of error control only)"; equation // enforcing error control on both direction rather than on the sum only // gives higher robustness. The reason may be that for bi-directional flow, // (VAB_flow - VBA_flow) may be close to zero. if forceErrorControlOnFlow then der(mExcAB) = mAB_flow; der(mExcBA) = mBA_flow; else der(mExcAB) = 0; der(mExcBA) = 0; end if; rho_a1_inflow = Medium.density(state_a1_inflow); rho_a2_inflow = Medium.density(state_a2_inflow); VZer_flow = vZer*A; mAB_flow = rho_a1_inflow*VAB_flow; mBA_flow = rho_a2_inflow*VBA_flow; // Average velocity (using the whole orifice area) vAB = VAB_flow/A; vBA = VBA_flow/A; port_a1.m_flow = mAB_flow; port_a2.m_flow = mBA_flow; // Energy balance (no storage, no heat loss/gain) port_a1.h_outflow = inStream(port_b1.h_outflow); port_b1.h_outflow = inStream(port_a1.h_outflow); port_a2.h_outflow = inStream(port_b2.h_outflow); port_b2.h_outflow = inStream(port_a2.h_outflow); // Mass balance (no storage) port_a1.m_flow = -port_b1.m_flow; port_a2.m_flow = -port_b2.m_flow; port_a1.Xi_outflow = inStream(port_b1.Xi_outflow); port_b1.Xi_outflow = inStream(port_a1.Xi_outflow); port_a2.Xi_outflow = inStream(port_b2.Xi_outflow); port_b2.Xi_outflow = inStream(port_a2.Xi_outflow); // Transport of trace substances port_a1.C_outflow = inStream(port_b1.C_outflow); port_b1.C_outflow = inStream(port_a1.C_outflow); port_a2.C_outflow = inStream(port_b2.C_outflow); port_b2.C_outflow = inStream(port_a2.C_outflow);end TwoWayFlowElement;
This is a partial model for models that describe the bi-directional air flow through large openings.
Models that extend this model need to compute
mAB_flow
and mBA_flow
,
or alternatively VAB_flow
and VBA_flow
,
and the face area area
.
The face area is a variable to allow this partial model to be used
for doors that can be open or closed as a function of an input signal.
Type | Name | Default | Description |
---|---|---|---|
Boolean | forceErrorControlOnFlow | true | Flag to force error control on m_flow. Set to true if interested in flow rate |
replaceable package Medium | PartialMedium | ||
Velocity | vZer | 0.001 | Minimum velocity to prevent zero flow. Recommended: 0.001 [m/s] |
Initialization | |||
MassFlowRate | m1_flow.start | 0 | Mass flow rate from port_a1 to port_b1 (m1_flow > 0 is design flow direction) [kg/s] |
Pressure | dp1.start | 0 | Pressure difference between port_a1 and port_b1 [Pa] |
MassFlowRate | m2_flow.start | 0 | Mass flow rate from port_a2 to port_b2 (m2_flow > 0 is design flow direction) [kg/s] |
Pressure | dp2.start | 0 | Pressure difference between port_a2 and port_b2 [Pa] |
Density | rho_a1_inflow.start | 1.2 | Density of air flowing in from port_a1 [kg/m3] |
Density | rho_a2_inflow.start | 1.2 | Density of air flowing in from port_a2 [kg/m3] |
Geometry | |||
Length | wOpe | 0.9 | Width of opening [m] |
Length | hOpe | 2.1 | Height of opening [m] |
Length | hA | 2.7/2 | Height of reference pressure zone A [m] |
Length | hB | 2.7/2 | Height of reference pressure zone B [m] |
Advanced | |||
Initialization | |||
SpecificEnthalpy | h_outflow_a1_start | Medium1.h_default | Start value for enthalpy flowing out of port a1 [J/kg] |
SpecificEnthalpy | h_outflow_b1_start | Medium1.h_default | Start value for enthalpy flowing out of port b1 [J/kg] |
SpecificEnthalpy | h_outflow_a2_start | Medium2.h_default | Start value for enthalpy flowing out of port a2 [J/kg] |
SpecificEnthalpy | h_outflow_b2_start | Medium2.h_default | Start value for enthalpy flowing out of port b2 [J/kg] |
MassFlowRate | m1_flow_small | 1E-4*abs(m1_flow_nominal) | Small mass flow rate for regularization of zero flow [kg/s] |
MassFlowRate | m2_flow_small | 1E-4*abs(m2_flow_nominal) | Small mass flow rate for regularization of zero flow [kg/s] |
Boolean | homotopyInitialization | true | = true, use homotopy method |
Diagnostics | |||
Boolean | show_T | false | = true, if actual temperature at port is computed |
Type | Name | Description |
---|---|---|
FluidPort_a | port_a1 | Fluid connector a1 (positive design flow direction is from port_a1 to port_b1) |
FluidPort_b | port_b1 | Fluid connector b1 (positive design flow direction is from port_a1 to port_b1) |
FluidPort_a | port_a2 | Fluid connector a2 (positive design flow direction is from port_a2 to port_b2) |
FluidPort_b | port_b2 | Fluid connector b2 (positive design flow direction is from port_a2 to port_b2) |
partial model TwoWayFlowElementBuoyancy "Flow resistance that uses the power law" extends Buildings.Airflow.Multizone.BaseClasses.TwoWayFlowElement; parameter Modelica.SIunits.Length wOpe=0.9 "|Geometry|Width of opening"; parameter Modelica.SIunits.Length hOpe=2.1 "|Geometry|Height of opening"; parameter Modelica.SIunits.Length hA=2.7/2 "|Geometry|Height of reference pressure zone A"; parameter Modelica.SIunits.Length hB=2.7/2 "|Geometry|Height of reference pressure zone B";end TwoWayFlowElementBuoyancy;
This is a partial model for computing the air exchange between volumes.
Models that extend this model need to provide an equation for
port_a1.m_flow
and port_a2.m_flow
.
Type | Name | Default | Description |
---|---|---|---|
replaceable package Medium1 | PartialMedium | Medium 1 in the component | |
replaceable package Medium2 | PartialMedium | Medium 2 in the component | |
Boolean | forceErrorControlOnFlow | true | Flag to force error control on m_flow. Set to true if interested in flow rate |
replaceable package Medium | Modelica.Media.Interfaces.Pa... | ||
Nominal condition | |||
MassFlowRate | m1_flow_nominal | 10/3600*1.2 | Nominal mass flow rate [kg/s] |
MassFlowRate | m2_flow_nominal | m1_flow_nominal | Nominal mass flow rate [kg/s] |
Initialization | |||
MassFlowRate | m1_flow.start | 0 | Mass flow rate from port_a1 to port_b1 (m1_flow > 0 is design flow direction) [kg/s] |
Pressure | dp1.start | 0 | Pressure difference between port_a1 and port_b1 [Pa] |
MassFlowRate | m2_flow.start | 0 | Mass flow rate from port_a2 to port_b2 (m2_flow > 0 is design flow direction) [kg/s] |
Pressure | dp2.start | 0 | Pressure difference between port_a2 and port_b2 [Pa] |
Assumptions | |||
Boolean | allowFlowReversal1 | false | = true to allow flow reversal in medium 1, false restricts to design direction (port_a -> port_b) |
Boolean | allowFlowReversal2 | false | = true to allow flow reversal in medium 2, false restricts to design direction (port_a -> port_b) |
Advanced | |||
Initialization | |||
SpecificEnthalpy | h_outflow_a1_start | Medium1.h_default | Start value for enthalpy flowing out of port a1 [J/kg] |
SpecificEnthalpy | h_outflow_b1_start | Medium1.h_default | Start value for enthalpy flowing out of port b1 [J/kg] |
SpecificEnthalpy | h_outflow_a2_start | Medium2.h_default | Start value for enthalpy flowing out of port a2 [J/kg] |
SpecificEnthalpy | h_outflow_b2_start | Medium2.h_default | Start value for enthalpy flowing out of port b2 [J/kg] |
MassFlowRate | m1_flow_small | 1E-4*abs(m1_flow_nominal) | Small mass flow rate for regularization of zero flow [kg/s] |
MassFlowRate | m2_flow_small | 1E-4*abs(m2_flow_nominal) | Small mass flow rate for regularization of zero flow [kg/s] |
Boolean | homotopyInitialization | true | = true, use homotopy method |
Diagnostics | |||
Boolean | show_T | false | = true, if actual temperature at port is computed |
Type | Name | Description |
---|---|---|
replaceable package Medium1 | Medium 1 in the component | |
replaceable package Medium2 | Medium 2 in the component | |
FluidPort_a | port_a1 | Fluid connector a1 (positive design flow direction is from port_a1 to port_b1) |
FluidPort_b | port_b1 | Fluid connector b1 (positive design flow direction is from port_a1 to port_b1) |
FluidPort_a | port_a2 | Fluid connector a2 (positive design flow direction is from port_a2 to port_b2) |
FluidPort_b | port_b2 | Fluid connector b2 (positive design flow direction is from port_a2 to port_b2) |
replaceable package Medium |
partial model ZonalFlow "Flow across zonal boundaries of a room" extends Buildings.Fluid.Interfaces.PartialFourPortInterface( redeclare final package Medium1 = Medium, redeclare final package Medium2 = Medium, final allowFlowReversal1 = false, final allowFlowReversal2 = false, final m1_flow_nominal = 10/3600*1.2, final m2_flow_nominal = m1_flow_nominal); extends Buildings.Airflow.Multizone.BaseClasses.ErrorControl; replaceable package Medium = Modelica.Media.Interfaces.PartialMedium; equation // Energy balance (no storage, no heat loss/gain) port_a1.h_outflow = inStream(port_b1.h_outflow); port_b1.h_outflow = inStream(port_a1.h_outflow); port_a2.h_outflow = inStream(port_b2.h_outflow); port_b2.h_outflow = inStream(port_a2.h_outflow); // Mass balance (no storage) port_a1.m_flow = -port_b1.m_flow; port_a2.m_flow = -port_b2.m_flow; port_a1.Xi_outflow = inStream(port_b1.Xi_outflow); port_b1.Xi_outflow = inStream(port_a1.Xi_outflow); port_a2.Xi_outflow = inStream(port_b2.Xi_outflow); port_b2.Xi_outflow = inStream(port_a2.Xi_outflow); // Transport of trace substances port_a1.C_outflow = inStream(port_b1.C_outflow); port_b1.C_outflow = inStream(port_a1.C_outflow); port_a2.C_outflow = inStream(port_b2.C_outflow); port_b2.C_outflow = inStream(port_a2.C_outflow);end ZonalFlow;
This is an interface that defines parameters used for error control.
Dymola does error control on state variables, such as temperature, pressure and
species concentration.
Flow variables such as m_flow
are typically not checked during the error control.
This can give large errors in flow variables, as long as the error on the volume's state variables
that are coupled to the flow variables is small.
Obtaining accurate flow variables can be achieved by imposing an error control
on the exchanged mass, which can be defined as
dm/dt = m_flow.By setting
enforceErrorControlOnFlow = true
, such an equation is imposed
by models that extend this class.
Type | Name | Default | Description |
---|---|---|---|
Boolean | forceErrorControlOnFlow | true | Flag to force error control on m_flow. Set to true if interested in flow rate |
model ErrorControl "Interface that defines parameters for error control" parameter Boolean forceErrorControlOnFlow = true "Flag to force error control on m_flow. Set to true if interested in flow rate";end ErrorControl;
This model describes the mass flow rate and pressure difference relation of an orifice in the form
V = k sign(Δp) |Δp|m
where V is the volume flow rate, k > 0 is a flow coefficient Δ p is the pressure drop and m ∈ [0.5, 1] is a flow coefficient. The equation is regularized for |Δp| < Δpt, where Δpt is a parameter. For turbulent flow, set m=1 ⁄ 2 and for laminar flow, set m=1.
The model is used for the interzonal air flow models.
For |Δp| < Δpt, the equation is regularized so that it is twice continuously differentiable in Δp, and that it has an infinite number of continuous derivatives in m and in k.
If m is not a function of time, then a, b, c and d can be pre-computed. In this situation, use Buildings.Airflow.Multizone.BaseClasses.powerLawFixedM, which allows to compute these values outside of this function, for example as parameters of a model.
Type | Name | Default | Description |
---|---|---|---|
Real | k | Flow coefficient, k = V_flow/ dp^m | |
Pressure | dp | Pressure difference [Pa] | |
Real | m | Flow exponent, m=0.5 for turbulent, m=1 for laminar | |
Pressure | dp_turbulent | 0.001 | Pressure difference where regularization starts [Pa] |
Type | Name | Description |
---|---|---|
VolumeFlowRate | V_flow | Volume flow rate [m3/s] |
function powerLaw "Power law used in orifice equations" input Real k "Flow coefficient, k = V_flow/ dp^m"; input Modelica.SIunits.Pressure dp "Pressure difference"; input Real m(min=0.5, max=1) "Flow exponent, m=0.5 for turbulent, m=1 for laminar"; input Modelica.SIunits.Pressure dp_turbulent(min=0)=0.001 "Pressure difference where regularization starts"; output Modelica.SIunits.VolumeFlowRate V_flow "Volume flow rate"; protected constant Real gamma(min=1) = 1.5 "Normalized flow rate where dphi(0)/dpi intersects phi(1)"; Real a "Polynomial coefficient for regularized implementation of flow resistance"; Real b "Polynomial coefficient for regularized implementation of flow resistance"; Real c "Polynomial coefficient for regularized implementation of flow resistance"; Real d "Polynomial coefficient for regularized implementation of flow resistance"; Real pi "Normalized pressure"; Real pi2 "Square of normalized pressure"; algorithm if (dp >= dp_turbulent) then V_flow := k*dp^m; elseif (dp <= -dp_turbulent) then V_flow :=-k*(-dp)^m; else a := gamma; b := 1/8*m^2 - 3*gamma - 3/2*m + 35.0/8; c := -1/4*m^2 + 3*gamma + 5/2*m - 21.0/4; d := 1/8*m^2 - gamma - m + 15.0/8; pi := dp/dp_turbulent; pi2 := pi*pi; V_flow := k*dp_turbulent^m * pi * ( a + pi2 * ( b + pi2 * ( c + pi2 * d))); end if;end powerLaw;
This model describes the mass flow rate and pressure difference relation of an orifice in the form
V = k sign(Δp) |Δp|m
where V is the volume flow rate, k > 0 is a flow coefficient Δ p is the pressure drop and m ∈ [0.5, 1] is a flow coefficient. The equation is regularized for |Δp| < Δpt, where Δpt is a parameter. For turbulent flow, set m=1 ⁄ 2 and for laminar flow, set m=1.
The model is used for the interzonal air flow models. It is identical to Buildings.Airflow.Multizone.BaseClasses.powerLaw but it requires the polynomial coefficients as an input. This allows a more efficient simulation if m and therefore also a, b, c and d are constant.
For |Δp| < Δpt, the equation is regularized so that it is twice continuously differentiable in Δp, and that it has an infinite number of continuous derivatives in m and in k.
If m, and therefore also a, b, c and d, change with time, then it is more convenient and efficient to use Buildings.Airflow.Multizone.BaseClasses.powerLaw.
Type | Name | Default | Description |
---|---|---|---|
Real | k | Flow coefficient, k = V_flow/ dp^m | |
Pressure | dp | Pressure difference [Pa] | |
Real | m | Flow exponent, m=0.5 for turbulent, m=1 for laminar | |
Real | a | Polynomial coefficient | |
Real | b | Polynomial coefficient | |
Real | c | Polynomial coefficient | |
Real | d | Polynomial coefficient | |
Pressure | dp_turbulent | 0.001 | Pressure difference where regularization starts [Pa] |
Type | Name | Description |
---|---|---|
VolumeFlowRate | V_flow | Volume flow rate [m3/s] |
function powerLawFixedM "Power law used in orifice equations when m is constant" input Real k "Flow coefficient, k = V_flow/ dp^m"; input Modelica.SIunits.Pressure dp "Pressure difference"; input Real m(min=0.5, max=1) "Flow exponent, m=0.5 for turbulent, m=1 for laminar"; input Real a "Polynomial coefficient"; input Real b "Polynomial coefficient"; input Real c "Polynomial coefficient"; input Real d "Polynomial coefficient"; input Modelica.SIunits.Pressure dp_turbulent(min=0)=0.001 "Pressure difference where regularization starts"; output Modelica.SIunits.VolumeFlowRate V_flow "Volume flow rate"; protected constant Real gamma(min=1) = 1.5 "Normalized flow rate where dphi(0)/dpi intersects phi(1)"; Real pi "Normalized pressure"; Real pi2 "Square of normalized pressure"; algorithm if (dp >= dp_turbulent) then V_flow := k*dp^m; elseif (dp <= -dp_turbulent) then V_flow :=-k*(-dp)^m; else pi := dp/dp_turbulent; pi2 := pi*pi; V_flow := k*dp_turbulent^m * pi * ( a + pi2 * ( b + pi2 * ( c + pi2 * d))); end if;end powerLawFixedM;
This function computes the wind pressure coefficient for low-rise buildings with rectangular shape. The correlation is the data fit from Swami and Chandra (1987), who fitted a function to various wind pressure coefficients from the literature. The same correlation is also implemented in CONTAM (Persily and Ivy, 2001).
The wind pressure coefficient is computed based on the natural logarithm of the side ratio of the walls, which is defined as
G = ln(x ⁄ y)
where x is the length of the wall that will be connected to this model, and y is the length of the adjacent wall as shown in the figure below.
Based on the wind incidence angle α and the side ratio
of the walls, the model computes how much the wind pressure
is attenuated compared to the reference wind pressure Cp0
.
The reference wind pressure Cp0
is a user-defined parameter,
and must be equal to the wind pressure at zero wind incidence angle, i.e.,
α = 0.
Swami and Chandra (1987) recommend Cp0 = 0.6 for
all low-rise buildings as this represents the average of
various values reported in the literature.
The attenuation factor is
Cp ⁄ Cp0 = ln(1.248 - 0.703 sin(α ⁄ 2) - 1.175 sin2(α) - 0.131 sin3(2 α G) + 0.769 cos(α ⁄ 2) +0.071 G2 * sin2(α ⁄ 2) + 0.717 cos2(α ⁄ 2)),
where Cp is the wind pressure coefficient for the current angle of incidence.
This function is used in Buildings.Fluid.Sources.Outside_CpLowRise which can be used directly with components of this package.
Symmetry requires that the first derivative of the wind pressure coefficient with respect to the incidence angle is zero for incidence angles of zero and π. However, the correlation of Swami and Chandra has non-zero derivatives at these values. In this implementation, the original function is therefore slightly modified for incidence angles between 0 and 5 degree, and between 175 and 180 degree. This leads to a model that is differentiable in the incidence angle, which generally leads to better numeric performance.
Type | Name | Default | Description |
---|---|---|---|
Real | Cp0 | Wind pressure coefficient for normal wind incidence angle | |
Angle | incAng | Wind incidence angle (0: normal to wall) [rad] | |
Real | G | Natural logarithm of side ratio |
Type | Name | Description |
---|---|---|
Real | Cp | Wind pressure coefficient |
function windPressureLowRise "Wind pressure coefficient for low-rise buildings" input Real Cp0(min=0) "Wind pressure coefficient for normal wind incidence angle"; input Modelica.SIunits.Angle incAng "Wind incidence angle (0: normal to wall)"; input Real G "Natural logarithm of side ratio"; output Real Cp "Wind pressure coefficient"; protected Modelica.SIunits.Angle aR "alpha, restricted to 0...pi"; Modelica.SIunits.Angle incAng2 "0.5*wind incidence angle"; Real sinA2 "=sin(alpha/2)"; Real cosA2 "=cos(alpha/2)"; Real a "Attenuation factor"; constant Modelica.SIunits.Angle pi2 = 2*Modelica.Constants.pi; constant Modelica.SIunits.Angle aRDel = 5*Modelica.Constants.pi/180 "Lower bound where transition occurs"; constant Modelica.SIunits.Angle aRDel2 = aRDel/2 "Half-width of transition interval"; constant Modelica.SIunits.Angle aRMax = 175*Modelica.Constants.pi/180 "Upper bound where transition occurs"; constant Real a180 = Modelica.Math.log(1.248 - 0.703 + 0.131*Modelica.Math.sin(2*Modelica.Constants.pi*G)^3 + 0.071*G^2) "Attenuation factor at 180 degree incidence angle"; algorithm // Restrict incAng to [0...pi] // Change sign to positive aR := if incAng < 0 then -incAng else incAng; // Constrain to [0...2*pi] if aR > pi2 then aR := aR - integer(aR/pi2)*pi2; end if; // Constrain to [0...pi] if aR > Modelica.Constants.pi then aR := pi2-aR; end if; // Evaluate eqn. 2-1 from FSEC-CR-163-86 incAng2 :=aR/2; sinA2 :=Modelica.Math.sin(incAng2); cosA2 :=Modelica.Math.cos(incAng2); // Implementation of the wind pressure coefficient that is once // continuously differentiable for all incidence angles if aR < aRDel then Cp :=Cp0*Buildings.Utilities.Math.Functions.spliceFunction( pos= Modelica.Math.log(1.248 - 0.703*sinA2 - 1.175*Modelica.Math.sin(aR)^2 + 0.131*Modelica.Math.sin(2*aR*G)^3 + 0.769*cosA2 + 0.071*G^2*sinA2^2 + 0.717*cosA2^2), neg=1, x=aR-aRDel2, deltax=aRDel2); elseif aR > aRMax then Cp := Cp0*Buildings.Utilities.Math.Functions.spliceFunction( pos= a180, neg=Modelica.Math.log(1.248 - 0.703*sinA2 - 1.175*Modelica.Math.sin(aR)^2 + 0.131*Modelica.Math.sin(2*aR*G)^3 + 0.769*cosA2 + 0.071*G^2*sinA2^2 + 0.717*cosA2^2), x=aR+aRDel2-Modelica.Constants.pi, deltax=aRDel2); else Cp :=Cp0*Modelica.Math.log(1.248 - 0.703*sinA2 - 1.175*Modelica.Math.sin(aR)^2 + 0.131*Modelica.Math.sin(2*aR*G)^3 + 0.769*cosA2 + 0.071*G^2*sinA2^2 + 0.717*cosA2^2); end if;end windPressureLowRise;