Buildings.Fluid.Actuators.Dampers

Package with air damper models

Information

This package contains component models for air dampers. For motor models, see Buildings.Fluid.Actuators.Motors.

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

Package Content

Name Description
Buildings.Fluid.Actuators.Dampers.Exponential Exponential Air damper with exponential opening characteristics
Buildings.Fluid.Actuators.Dampers.MixingBox MixingBox Outside air mixing box with interlocked air dampers
Buildings.Fluid.Actuators.Dampers.MixingBoxMinimumFlow MixingBoxMinimumFlow Outside air mixing box with parallel damper for minimum outside air flow rate
Buildings.Fluid.Actuators.Dampers.PressureIndependent PressureIndependent Model for an air damper whose mass flow is proportional to the input signal
Buildings.Fluid.Actuators.Dampers.VAVBoxExponential VAVBoxExponential VAV box with a fixed resistance plus a damper model withe exponential characteristics
Buildings.Fluid.Actuators.Dampers.Examples Examples Collection of models that illustrate model use and test models

Buildings.Fluid.Actuators.Dampers.Exponential Buildings.Fluid.Actuators.Dampers.Exponential

Air damper with exponential opening characteristics

Buildings.Fluid.Actuators.Dampers.Exponential

Information

This model is an air damper with flow coefficient that is an exponential function of the opening angle. The model is as in ASHRAE 825-RP. A control signal of y=0 means the damper is closed, and y=1 means the damper is open. This is opposite of the implementation of ASHRAE 825-RP, but used here for consistency within this library.

For yL < y < yU, the damper characteristics is

kd(y) = exp(a+b (1-y)).

Outside this range, the damper characteristic is defined by a quadratic polynomial that matches the damper resistance at y=0 and y=yL or y=yU and y=1, respectively. In addition, the polynomials are such that kd(y) is differentiable in y and the derivative is continuous.

The damper characteristics kd(y) is then used to compute the flow coefficient k(y) as

k(y) = (2 ρ ⁄ kd(y))1/2 A,

where A is the face area, which is computed using the nominal mass flow rate m_flow_nominal, the nominal velocity v_nominal and the density of the medium. The flow coefficient k(y) is used to compute the mass flow rate versus pressure drop relation as

m = sign(Δp) k(y) √ Δp  

with regularization near the origin.

ASHRAE 825-RP lists the following parameter values as typical:

opposed bladessingle blades
yL15/9015/90
yU55/9065/90
k01E61E6
k10.2 to 0.50.2 to 0.5
a-1.51-1.51
b0.105*900.0842*90

References

P. Haves, L. K. Norford, M. DeSimone and L. Mei, A Standard Simulation Testbed for the Evaluation of Control Algorithms & Strategies, ASHRAE Final Report 825-RP, Atlanta, GA.

Extends from Buildings.Fluid.Actuators.BaseClasses.PartialDamperExponential (Partial model for air dampers with exponential opening characteristics).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
MassFlowRatem_flow_turbulentif use_deltaM then deltaM*m_...Turbulent flow if |m_flow| >= m_flow_turbulent [kg/s]
Booleanuse_deltaMtrueSet to true to use deltaM for turbulent transition, else ReC is used
RealdeltaM0.3Fraction of nominal mass flow rate where transition to turbulent occurs
Velocityv_nominal1Nominal face velocity [m/s]
BooleanroundDuctfalseSet to true for round duct, false for square cross section
RealReC4000Reynolds number where transition to turbulent starts
RealkFixed0Flow coefficient of fixed resistance that may be in series with damper, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2).
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
PressureDifferencedp_nominal(m_flow_nominal/kDam_default...Pressure drop at nominal mass flow rate [Pa]
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
Booleanfrom_dpfalse= true, use m_flow = f(dp) else dp = f(m_flow)
BooleanhomotopyInitializationtrue= true, use homotopy method
Booleanlinearizedfalse= true, use linear relation between m_flow and dp for any flow rate
Booleanuse_constant_densitytrueSet to true to use constant density for flow friction
Dynamics
Filtered opening
Booleanuse_inputFiltertrue= true, if opening is filtered with a 2nd order CriticalDamping filter
TimeriseTime120Rise time of the filter (time to reach 99.6 % of an opening step) [s]
Integerorder2Order of filter
InitinitModelica.Blocks.Types.Init.I...Type of initialization (no init/steady state/initial state/initial output)
Realy_start1Initial value of output
Damper coefficients
Reala-1.51Coefficient a for damper characteristics
Realb0.105*90Coefficient b for damper characteristics
RealyL15/90Lower value for damper curve
RealyU55/90Upper value for damper curve
Realk01E6Flow coefficient for y=0, k0 = pressure drop divided by dynamic pressure
Realk10.45Flow coefficient for y=1, k1 = pressure drop divided by dynamic pressure

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)
input RealInputyActuator position (0: closed, 1: open)
output RealOutputy_actualActual valve position

Modelica definition

model Exponential "Air damper with exponential opening characteristics" extends Buildings.Fluid.Actuators.BaseClasses.PartialDamperExponential( final dp_nominal=(m_flow_nominal/kDam_default)^2, dp(nominal=10), final kFixed=0); protected parameter Real kDam_default = sqrt(2*rho_default)*A/kThetaSqRt_default "Flow coefficient for damper, k=m_flow/sqrt(dp), with unit=(kg*m)^(1/2)"; parameter Real kThetaSqRt_default(min=0)= Buildings.Fluid.Actuators.BaseClasses.exponentialDamper( y=1, a=a, b=b, cL=cL, cU=cU, yL=yL, yU=yU) "Flow coefficient, kThetaSqRt = sqrt(pressure drop divided by dynamic pressure)"; initial algorithm assert(kThetaSqRt_default>0, "Flow coefficient must be strictly positive."); end Exponential;

Buildings.Fluid.Actuators.Dampers.MixingBox Buildings.Fluid.Actuators.Dampers.MixingBox

Outside air mixing box with interlocked air dampers

Buildings.Fluid.Actuators.Dampers.MixingBox

Information

Model of an outside air mixing box with air dampers. Set y=0 to close the outside air and exhast air dampers.

If dp_nominalIncludesDamper=true, then the parameter dp_nominal is equal to the pressure drop of the damper plus the fixed flow resistance at the nominal flow rate. If dp_nominalIncludesDamper=false, then dp_nominal does not include the flow resistance of the air damper.

Extends from Buildings.Fluid.Actuators.BaseClasses.ActuatorSignal (Partial model that implements the filtered opening for valves and dampers).

Parameters

TypeNameDefaultDescription
replaceable package MediumModelica.Media.Interfaces.Pa...Medium in the component
Booleanuse_deltaMtrueSet to true to use deltaM for turbulent transition, else ReC is used
RealdeltaM0.3Fraction of nominal mass flow rate where transition to turbulent occurs
Velocityv_nominal1Nominal face velocity [m/s]
BooleanroundDuctfalseSet to true for round duct, false for square cross section
RealReC4000Reynolds number where transition to turbulent starts
Nominal condition
Booleandp_nominalIncludesDamperfalseset to true if dp_nominal includes the pressure loss of the open damper
MassFlowRatemOut_flow_nominal Mass flow rate outside air damper [kg/s]
PressureDifferencedpOut_nominaldpOut_nominal( ...Pressure drop outside air leg [Pa]
MassFlowRatemRec_flow_nominal Mass flow rate recirculation air damper [kg/s]
PressureDifferencedpRec_nominaldpRec_nominal( ...Pressure drop recirculation air leg [Pa]
MassFlowRatemExh_flow_nominal Mass flow rate exhaust air damper [kg/s]
PressureDifferencedpExh_nominaldpExh_nominal( ...Pressure drop exhaust air leg [Pa]
Dynamics
Filtered opening
Booleanuse_inputFiltertrue= true, if opening is filtered with a 2nd order CriticalDamping filter
TimeriseTime120Rise time of the filter (time to reach 99.6 % of an opening step) [s]
Integerorder2Order of filter
InitinitModelica.Blocks.Types.Init.I...Type of initialization (no init/steady state/initial state/initial output)
Realy_start1Initial value of output
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
Booleanfrom_dptrue= true, use m_flow = f(dp) else dp = f(m_flow)
Booleanlinearizedfalse= true, use linear relation between m_flow and dp for any flow rate
Booleanuse_constant_densitytrueSet to true to use constant density for flow friction
Damper coefficients
Reala-1.51Coefficient a for damper characteristics
Realb0.105*90Coefficient b for damper characteristics
RealyL15/90Lower value for damper curve
RealyU55/90Upper value for damper curve
Realk01E6Flow coefficient for y=0, k0 = pressure drop divided by dynamic pressure
Realk10.45Flow coefficient for y=1, k1 = pressure drop divided by dynamic pressure

Connectors

TypeNameDescription
input RealInputyActuator position (0: closed, 1: open)
output RealOutputy_actualActual valve position
replaceable package MediumMedium in the component
FluidPort_aport_OutFluid connector a (positive design flow direction is from port_a to port_b)
FluidPort_bport_ExhFluid connector b (positive design flow direction is from port_a to port_b)
FluidPort_aport_RetFluid connector a (positive design flow direction is from port_a to port_b)
FluidPort_bport_SupFluid connector b (positive design flow direction is from port_a to port_b)

Modelica definition

model MixingBox "Outside air mixing box with interlocked air dampers" extends Buildings.Fluid.Actuators.BaseClasses.ActuatorSignal; replaceable package Medium = Modelica.Media.Interfaces.PartialMedium "Medium in the component"; import Modelica.Constants; parameter Boolean allowFlowReversal = true "= false to simplify equations, assuming, but not enforcing, no flow reversal"; VAVBoxExponential damOA( redeclare package Medium = Medium, dp_nominal=dpOut_nominal, dp_nominalIncludesDamper=dp_nominalIncludesDamper, from_dp=from_dp, linearized=linearized, use_deltaM=use_deltaM, deltaM=deltaM, roundDuct=roundDuct, ReC=ReC, a=a, b=b, yL=yL, yU=yU, k0=k0, k1=k1, use_constant_density=use_constant_density, allowFlowReversal=allowFlowReversal, m_flow_nominal=mOut_flow_nominal, final use_inputFilter=false); parameter Boolean use_deltaM = true "Set to true to use deltaM for turbulent transition, else ReC is used"; parameter Real deltaM = 0.3 "Fraction of nominal mass flow rate where transition to turbulent occurs"; parameter Modelica.SIunits.Velocity v_nominal=1 "Nominal face velocity"; parameter Boolean roundDuct = false "Set to true for round duct, false for square cross section"; parameter Real ReC=4000 "Reynolds number where transition to turbulent starts"; VAVBoxExponential damExh( redeclare package Medium = Medium, m_flow_nominal=mExh_flow_nominal, dp_nominal=dpExh_nominal, dp_nominalIncludesDamper=dp_nominalIncludesDamper, from_dp=from_dp, linearized=linearized, use_deltaM=use_deltaM, deltaM=deltaM, roundDuct=roundDuct, ReC=ReC, a=a, b=b, yL=yL, yU=yU, k0=k0, k1=k1, use_constant_density=use_constant_density, allowFlowReversal=allowFlowReversal, final use_inputFilter=false) "Exhaust air damper"; VAVBoxExponential damRec( redeclare package Medium = Medium, m_flow_nominal=mRec_flow_nominal, dp_nominal=dpRec_nominal, dp_nominalIncludesDamper=dp_nominalIncludesDamper, from_dp=from_dp, linearized=linearized, use_deltaM=use_deltaM, deltaM=deltaM, roundDuct=roundDuct, ReC=ReC, a=a, b=b, yL=yL, yU=yU, k0=k0, k1=k1, use_constant_density=use_constant_density, allowFlowReversal=allowFlowReversal, final use_inputFilter=false) "Recirculation air damper"; parameter Boolean dp_nominalIncludesDamper=false "set to true if dp_nominal includes the pressure loss of the open damper"; parameter Modelica.SIunits.MassFlowRate mOut_flow_nominal "Mass flow rate outside air damper"; parameter Modelica.SIunits.PressureDifference dpOut_nominal(min=0, displayUnit="Pa") "Pressure drop outside air leg"; parameter Modelica.SIunits.MassFlowRate mRec_flow_nominal "Mass flow rate recirculation air damper"; parameter Modelica.SIunits.PressureDifference dpRec_nominal(min=0, displayUnit="Pa") "Pressure drop recirculation air leg"; parameter Modelica.SIunits.MassFlowRate mExh_flow_nominal "Mass flow rate exhaust air damper"; parameter Modelica.SIunits.PressureDifference dpExh_nominal(min=0, displayUnit="Pa") "Pressure drop exhaust air leg"; parameter Boolean from_dp=true "= true, use m_flow = f(dp) else dp = f(m_flow)"; parameter Boolean linearized=false "= true, use linear relation between m_flow and dp for any flow rate"; parameter Boolean use_constant_density=true "Set to true to use constant density for flow friction"; parameter Real a=-1.51 "Coefficient a for damper characteristics"; parameter Real b=0.105*90 "Coefficient b for damper characteristics"; parameter Real yL=15/90 "Lower value for damper curve"; parameter Real yU=55/90 "Upper value for damper curve"; parameter Real k0=1E6 "Flow coefficient for y=0, k0 = pressure drop divided by dynamic pressure"; parameter Real k1=0.45 "Flow coefficient for y=1, k1 = pressure drop divided by dynamic pressure"; Modelica.Fluid.Interfaces.FluidPort_a port_Out(redeclare package Medium = Medium, m_flow(start=0, min=if allowFlowReversal then -Constants.inf else 0)) "Fluid connector a (positive design flow direction is from port_a to port_b)"; Modelica.Fluid.Interfaces.FluidPort_b port_Exh(redeclare package Medium = Medium, m_flow(start=0, max=if allowFlowReversal then +Constants.inf else 0)) "Fluid connector b (positive design flow direction is from port_a to port_b)"; Modelica.Fluid.Interfaces.FluidPort_a port_Ret(redeclare package Medium = Medium, m_flow(start=0, min=if allowFlowReversal then -Constants.inf else 0)) "Fluid connector a (positive design flow direction is from port_a to port_b)"; Modelica.Fluid.Interfaces.FluidPort_b port_Sup(redeclare package Medium = Medium, m_flow(start=0, max=if allowFlowReversal then +Constants.inf else 0)) "Fluid connector b (positive design flow direction is from port_a to port_b)"; Modelica.Blocks.Sources.Constant uni(k=1) "Unity signal"; Modelica.Blocks.Math.Add add(k2=-1) "Adder"; protected parameter Medium.Density rho_default=Medium.density(sta_default) "Density, used to compute fluid volume"; parameter Medium.ThermodynamicState sta_default= Medium.setState_pTX(T=Medium.T_default, p=Medium.p_default, X=Medium.X_default) "Default medium state"; equation connect(uni.y, add.u1); connect(add.y, damRec.y); connect(damOA.port_a, port_Out); connect(damExh.port_b, port_Exh); connect(port_Sup, damOA.port_b); connect(damRec.port_b, port_Sup); connect(port_Ret, damExh.port_a); connect(port_Ret, damRec.port_a); connect(y_actual, add.u2); connect(y_actual, damOA.y); connect(y_actual, damExh.y); end MixingBox;

Buildings.Fluid.Actuators.Dampers.MixingBoxMinimumFlow Buildings.Fluid.Actuators.Dampers.MixingBoxMinimumFlow

Outside air mixing box with parallel damper for minimum outside air flow rate

Buildings.Fluid.Actuators.Dampers.MixingBoxMinimumFlow

Information

Model of an outside air mixing box with air dampers and a flow path for the minimum outside air flow rate.

If dp_nominalIncludesDamper=true, then the parameter dp_nominal is equal to the pressure drop of the damper plus the fixed flow resistance at the nominal flow rate. If dp_nominalIncludesDamper=false, then dp_nominal does not include the flow resistance of the air damper.

Extends from Buildings.Fluid.Actuators.Dampers.MixingBox (Outside air mixing box with interlocked air dampers).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Booleanuse_deltaMtrueSet to true to use deltaM for turbulent transition, else ReC is used
RealdeltaM0.3Fraction of nominal mass flow rate where transition to turbulent occurs
Velocityv_nominal1Nominal face velocity [m/s]
BooleanroundDuctfalseSet to true for round duct, false for square cross section
RealReC4000Reynolds number where transition to turbulent starts
Nominal condition
Booleandp_nominalIncludesDamperfalseset to true if dp_nominal includes the pressure loss of the open damper
MassFlowRatemOut_flow_nominal Mass flow rate outside air damper [kg/s]
PressureDifferencedpOut_nominal Pressure drop outside air leg [Pa]
MassFlowRatemRec_flow_nominal Mass flow rate recirculation air damper [kg/s]
PressureDifferencedpRec_nominal Pressure drop recirculation air leg [Pa]
MassFlowRatemExh_flow_nominal Mass flow rate exhaust air damper [kg/s]
PressureDifferencedpExh_nominal Pressure drop exhaust air leg [Pa]
MassFlowRatemOutMin_flow_nominal Mass flow rate minimum outside air damper [kg/s]
PressureDifferencedpOutMin_nominaldpOutMin_nominal( ...Pressure drop minimum outside air leg [Pa]
Dynamics
Filtered opening
Booleanuse_inputFiltertrue= true, if opening is filtered with a 2nd order CriticalDamping filter
TimeriseTime120Rise time of the filter (time to reach 99.6 % of an opening step) [s]
Integerorder2Order of filter
InitinitModelica.Blocks.Types.Init.I...Type of initialization (no init/steady state/initial state/initial output)
Realy_start1Initial value of output
RealyOutMin_starty_startInitial value of signal for minimum outside air damper
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
Booleanfrom_dptrue= true, use m_flow = f(dp) else dp = f(m_flow)
Booleanlinearizedfalse= true, use linear relation between m_flow and dp for any flow rate
Booleanuse_constant_densitytrueSet to true to use constant density for flow friction
Damper coefficients
Reala-1.51Coefficient a for damper characteristics
Realb0.105*90Coefficient b for damper characteristics
RealyL15/90Lower value for damper curve
RealyU55/90Upper value for damper curve
Realk01E6Flow coefficient for y=0, k0 = pressure drop divided by dynamic pressure
Realk10.45Flow coefficient for y=1, k1 = pressure drop divided by dynamic pressure

Connectors

TypeNameDescription
input RealInputyActuator position (0: closed, 1: open)
output RealOutputy_actualActual valve position
FluidPort_aport_OutFluid connector a (positive design flow direction is from port_a to port_b)
FluidPort_bport_ExhFluid connector b (positive design flow direction is from port_a to port_b)
FluidPort_aport_RetFluid connector a (positive design flow direction is from port_a to port_b)
FluidPort_bport_SupFluid connector b (positive design flow direction is from port_a to port_b)
FluidPort_aport_OutMinFluid connector a (positive design flow direction is from port_a to port_b)
input RealInputyOutMinDamper position minimum outside air (0: closed, 1: open)
output RealOutputyOutMin_actualActual valve position

Modelica definition

model MixingBoxMinimumFlow "Outside air mixing box with parallel damper for minimum outside air flow rate" extends Buildings.Fluid.Actuators.Dampers.MixingBox; import Modelica.Constants; final parameter Modelica.SIunits.Area AOutMin = mOutMin_flow_nominal/rho_default/v_nominal "Face area minimum outside air damper"; parameter Modelica.SIunits.MassFlowRate mOutMin_flow_nominal "Mass flow rate minimum outside air damper"; parameter Modelica.SIunits.PressureDifference dpOutMin_nominal(min=0, displayUnit="Pa") "Pressure drop minimum outside air leg"; parameter Real yOutMin_start=y_start "Initial value of signal for minimum outside air damper"; Modelica.Fluid.Interfaces.FluidPort_a port_OutMin(redeclare package Medium = Medium, m_flow(start=0, min=if allowFlowReversal then -Constants.inf else 0)) "Fluid connector a (positive design flow direction is from port_a to port_b)"; Modelica.Blocks.Interfaces.RealInput yOutMin "Damper position minimum outside air (0: closed, 1: open)"; Modelica.Blocks.Interfaces.RealOutput yOutMin_actual "Actual valve position"; Buildings.Fluid.Actuators.Dampers.VAVBoxExponential damOAMin( redeclare package Medium = Medium, dp_nominalIncludesDamper=dp_nominalIncludesDamper, from_dp=from_dp, linearized=linearized, use_deltaM=use_deltaM, deltaM=deltaM, v_nominal=v_nominal, roundDuct=roundDuct, ReC=ReC, a=a, b=b, yL=yL, yU=yU, k0=k0, k1=k1, use_constant_density=use_constant_density, allowFlowReversal=allowFlowReversal, m_flow_nominal=mOutMin_flow_nominal, dp_nominal=dpOutMin_nominal, final use_inputFilter=false) "Damper for minimum outside air intake"; protected Modelica.Blocks.Interfaces.RealOutput yOutMin_filtered if use_inputFilter "Filtered damper position in the range 0..1"; Modelica.Blocks.Continuous.Filter filterOutMin( order=2, f_cut=5/(2*Modelica.Constants.pi*riseTime), final init=init, final y_start=yOutMin_start, final analogFilter=Modelica.Blocks.Types.AnalogFilter.CriticalDamping, final filterType=Modelica.Blocks.Types.FilterType.LowPass, x(each stateSelect=StateSelect.always)) if use_inputFilter "Second order filter to approximate valve opening time, and to improve numerics"; equation connect(filterOutMin.y, yOutMin_filtered); if use_inputFilter then connect(yOutMin, filterOutMin.u); connect(filterOutMin.y, yOutMin_actual); else connect(yOutMin, yOutMin_actual); end if; ////// connect(port_OutMin, damOAMin.port_a); connect(damOAMin.port_b, port_Sup); connect(yOutMin_actual, damOAMin.y); end MixingBoxMinimumFlow;

Buildings.Fluid.Actuators.Dampers.PressureIndependent Buildings.Fluid.Actuators.Dampers.PressureIndependent

Model for an air damper whose mass flow is proportional to the input signal

Buildings.Fluid.Actuators.Dampers.PressureIndependent

Information

Model for an air damper whose airflow is proportional to the input signal, assuming that at y = 1, m_flow = m_flow_nominal. This is unless the pressure difference dp is too low, in which case a kDam = m_flow_nominal/sqrt(dp_nominal) characteristic is used.

The model is similar to Buildings.Fluid.Actuators.Valves.TwoWayPressureIndependent, except for adaptations for damper parameters. Please see that documentation for more information.

Extends from Buildings.Fluid.BaseClasses.PartialResistance (Partial model for a hydraulic resistance), Buildings.Fluid.Actuators.BaseClasses.ActuatorSignal (Partial model that implements the filtered opening for valves and dampers).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
MassFlowRatem_flow_turbulentif use_deltaM then deltaM*m_...Turbulent flow if |m_flow| >= m_flow_turbulent [kg/s]
Booleanuse_deltaMtrueSet to true to use deltaM for turbulent transition, else ReC is used
RealdeltaM0.3Fraction of nominal mass flow rate where transition to turbulent occurs
Velocityv_nominal1Nominal face velocity [m/s]
BooleanroundDuctfalseSet to true for round duct, false for square cross section
RealReC4000Reynolds number where transition to turbulent starts
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
PressureDifferencedp_nominal Pressure drop at nominal mass flow rate [Pa]
PressureDifferencedpFixed_nominal0Pressure drop of duct and other resistances that are in series [Pa]
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
Booleanfrom_dptrue= true, use m_flow = f(dp) else dp = f(m_flow)
BooleanhomotopyInitializationtrue= true, use homotopy method
Booleanlinearizedfalse= true, use linear relation between m_flow and dp for any flow rate
Booleanuse_constant_densitytrueSet to true to use constant density for flow friction
Reall0.0001Damper leakage, l=kDam(y=0)/kDam(y=1)
Reall20.01Gain for mass flow increase if pressure is above nominal pressure [1]
Realdeltax0.02Transition interval for flow rate [1]
Dynamics
Filtered opening
Booleanuse_inputFiltertrue= true, if opening is filtered with a 2nd order CriticalDamping filter
TimeriseTime120Rise time of the filter (time to reach 99.6 % of an opening step) [s]
Integerorder2Order of filter
InitinitModelica.Blocks.Types.Init.I...Type of initialization (no init/steady state/initial state/initial output)
Realy_start1Initial value of output

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)
input RealInputyActuator position (0: closed, 1: open)
output RealOutputy_actualActual valve position

Modelica definition

model PressureIndependent "Model for an air damper whose mass flow is proportional to the input signal" extends Buildings.Fluid.BaseClasses.PartialResistance( m_flow_turbulent=if use_deltaM then deltaM * m_flow_nominal else eta_default*ReC*sqrt(A)*facRouDuc, final linearized = false, from_dp=true); extends Buildings.Fluid.Actuators.BaseClasses.ActuatorSignal; parameter Boolean use_deltaM = true "Set to true to use deltaM for turbulent transition, else ReC is used"; parameter Real deltaM = 0.3 "Fraction of nominal mass flow rate where transition to turbulent occurs"; parameter Modelica.SIunits.Velocity v_nominal = 1 "Nominal face velocity"; final parameter Modelica.SIunits.Area A=m_flow_nominal/rho_default/v_nominal "Face area"; parameter Boolean roundDuct = false "Set to true for round duct, false for square cross section"; parameter Real ReC=4000 "Reynolds number where transition to turbulent starts"; parameter Boolean use_constant_density=true "Set to true to use constant density for flow friction"; parameter Modelica.SIunits.PressureDifference dpFixed_nominal(displayUnit="Pa", min=0) = 0 "Pressure drop of duct and other resistances that are in series"; parameter Real l(min=1e-10, max=1) = 0.0001 "Damper leakage, l=kDam(y=0)/kDam(y=1)"; input Real phi = l + y_actual*(1 - l) "Ratio actual to nominal mass flow rate of damper, phi=kDam(y)/kDam(y=1)"; parameter Real l2(unit="1", min=1e-10) = 0.01 "Gain for mass flow increase if pressure is above nominal pressure"; parameter Real deltax(unit="1", min=1E-5) = 0.02 "Transition interval for flow rate"; Medium.Density rho "Medium density"; protected parameter Medium.Density rho_default=Medium.density(sta_default) "Density, used to compute fluid volume"; parameter Real facRouDuc= if roundDuct then sqrt(Modelica.Constants.pi)/2 else 1; parameter Real kDam(unit="") = m_flow_nominal/sqrt(dp_nominal_pos) "Flow coefficient of damper, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2)"; parameter Real kFixed(unit="") = if dpFixed_nominal > Modelica.Constants.eps then m_flow_nominal / sqrt(dpFixed_nominal) else 0 "Flow coefficient of fixed resistance in series with damper, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2)"; parameter Real k = if (dpFixed_nominal > Modelica.Constants.eps) then sqrt(1/(1/kFixed^2 + 1/kDam^2)) else kDam "Flow coefficient of damper plus fixed resistance"; parameter Real coeff1 = l2/dp_nominal*m_flow_nominal "Parameter for avoiding unnecessary computations"; parameter Real coeff2 = 1/coeff1 "Parameter for avoiding unnecessary computations"; constant Real y2dd = 0 "Second derivative at second support point"; Modelica.SIunits.MassFlowRate m_flow_set "Requested mass flow rate"; Modelica.SIunits.PressureDifference dp_min(displayUnit="Pa") "Minimum pressure difference required for delivering requested mass flow rate"; Modelica.SIunits.PressureDifference dp_x, dp_x1, dp_x2, dp_y2, dp_y1 "Support points for interpolation flow functions"; Modelica.SIunits.MassFlowRate m_flow_x, m_flow_x1, m_flow_x2, m_flow_y2, m_flow_y1 "Support points for interpolation flow functions"; Modelica.SIunits.MassFlowRate m_flow_smooth "Smooth interpolation result between two flow regimes"; Modelica.SIunits.PressureDifference dp_smooth "Smooth interpolation result between two flow regimes"; initial equation assert(m_flow_turbulent > 0, "m_flow_turbulent must be bigger than zero."); equation rho = if use_constant_density then rho_default else Medium.density(Medium.setState_phX(port_a.p, inStream(port_a.h_outflow), inStream(port_a.Xi_outflow))); // From TwoWayPressureIndependent valve model m_flow_set = m_flow_nominal*phi; dp_min = Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_m_flow( m_flow=m_flow_set, k=k, m_flow_turbulent=m_flow_turbulent); if from_dp then m_flow_x=0; m_flow_x1=0; m_flow_x2=0; dp_y1=0; dp_y2=0; dp_smooth=0; dp_x = dp-dp_min; dp_x1 = -dp_x2; dp_x2 = deltax*dp_min; // min function ensures that m_flow_y1 does not increase further for dp_x > dp_x1 m_flow_y1 = Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_dp( dp=min(dp, dp_min+dp_x1), k=k, m_flow_turbulent=m_flow_turbulent); // max function ensures that m_flow_y2 does not decrease further for dp_x < dp_x2 m_flow_y2 = m_flow_set + coeff1*max(dp_x,dp_x2); m_flow_smooth = noEvent(smooth(2, if dp_x <= dp_x1 then m_flow_y1 elseif dp_x >=dp_x2 then m_flow_y2 else Buildings.Utilities.Math.Functions.quinticHermite( x=dp_x, x1=dp_x1, x2=dp_x2, y1=m_flow_y1, y2=m_flow_y2, y1d= Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_dp_der( dp=dp_min + dp_x1, k=k, m_flow_turbulent=m_flow_turbulent, dp_der=1), y2d=coeff1, y1dd=Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_dp_der2( dp=dp_min + dp_x1, k=k, m_flow_turbulent=m_flow_turbulent, dp_der=1, dp_der2=0), y2dd=y2dd))); else dp_x=0; dp_x1=0; dp_x2=0; m_flow_y1=0; m_flow_y2=0; m_flow_smooth=0; m_flow_x = m_flow-m_flow_set; m_flow_x1 = -m_flow_x2; m_flow_x2 = deltax*m_flow_set; // min function ensures that dp_y1 does not increase further for m_flow_x > m_flow_x1 dp_y1 = Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_m_flow( m_flow=min(m_flow, m_flow_set + m_flow_x1), k=k, m_flow_turbulent=m_flow_turbulent); // max function ensures that dp_y2 does not decrease further for m_flow_x < m_flow_x2 dp_y2 = dp_min + coeff2*max(m_flow_x, m_flow_x2); dp_smooth = noEvent(smooth(2, if m_flow_x <= m_flow_x1 then dp_y1 elseif m_flow_x >=m_flow_x2 then dp_y2 else Buildings.Utilities.Math.Functions.quinticHermite( x=m_flow_x, x1=m_flow_x1, x2=m_flow_x2, y1=dp_y1, y2=dp_y2, y1d=Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_m_flow_der( m_flow=m_flow_set + m_flow_x1, k=k, m_flow_turbulent=m_flow_turbulent, m_flow_der=1), y2d=coeff2, y1dd=Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_m_flow_der2( m_flow=m_flow_set + m_flow_x1, k=k, m_flow_turbulent=m_flow_turbulent, m_flow_der=1, m_flow_der2=0), y2dd=y2dd))); end if; if homotopyInitialization then if from_dp then m_flow=homotopy(actual=m_flow_smooth, simplified=m_flow_nominal_pos*dp/dp_nominal_pos); else dp=homotopy(actual=dp_smooth, simplified=dp_nominal_pos*m_flow/m_flow_nominal_pos); end if; else if from_dp then m_flow=m_flow_smooth; else dp=dp_smooth; end if; end if; end PressureIndependent;

Buildings.Fluid.Actuators.Dampers.VAVBoxExponential Buildings.Fluid.Actuators.Dampers.VAVBoxExponential

VAV box with a fixed resistance plus a damper model withe exponential characteristics

Buildings.Fluid.Actuators.Dampers.VAVBoxExponential

Information

Model of two resistances in series. One resistance has a fixed flow coefficient, the other resistance is an air damper whose flow coefficient is an exponential function of the opening angle.

If dp_nominalIncludesDamper=true, then the parameter dp_nominal is equal to the pressure drop of the damper plus the fixed flow resistance at the nominal flow rate. If dp_nominalIncludesDamper=false, then dp_nominal does not include the flow resistance of the air damper.

Extends from Buildings.Fluid.Actuators.BaseClasses.PartialDamperExponential (Partial model for air dampers with exponential opening characteristics).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
MassFlowRatem_flow_turbulentif use_deltaM then deltaM*m_...Turbulent flow if |m_flow| >= m_flow_turbulent [kg/s]
Booleanuse_deltaMtrueSet to true to use deltaM for turbulent transition, else ReC is used
RealdeltaM0.3Fraction of nominal mass flow rate where transition to turbulent occurs
Velocityv_nominal1Nominal face velocity [m/s]
BooleanroundDuctfalseSet to true for round duct, false for square cross section
RealReC4000Reynolds number where transition to turbulent starts
RealkFixedsqrt(kResSqu)Flow coefficient of fixed resistance that may be in series with damper, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2).
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
PressureDifferencedp_nominal Pressure drop at nominal mass flow rate [Pa]
Booleandp_nominalIncludesDampertrueset to true if dp_nominal includes the pressure loss of the open damper
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
Booleanfrom_dpfalse= true, use m_flow = f(dp) else dp = f(m_flow)
BooleanhomotopyInitializationtrue= true, use homotopy method
Booleanlinearizedfalse= true, use linear relation between m_flow and dp for any flow rate
Booleanuse_constant_densitytrueSet to true to use constant density for flow friction
Dynamics
Filtered opening
Booleanuse_inputFiltertrue= true, if opening is filtered with a 2nd order CriticalDamping filter
TimeriseTime120Rise time of the filter (time to reach 99.6 % of an opening step) [s]
Integerorder2Order of filter
InitinitModelica.Blocks.Types.Init.I...Type of initialization (no init/steady state/initial state/initial output)
Realy_start1Initial value of output
Damper coefficients
Reala-1.51Coefficient a for damper characteristics
Realb0.105*90Coefficient b for damper characteristics
RealyL15/90Lower value for damper curve
RealyU55/90Upper value for damper curve
Realk01E6Flow coefficient for y=0, k0 = pressure drop divided by dynamic pressure
Realk10.45Flow coefficient for y=1, k1 = pressure drop divided by dynamic pressure

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)
input RealInputyActuator position (0: closed, 1: open)
output RealOutputy_actualActual valve position

Modelica definition

model VAVBoxExponential "VAV box with a fixed resistance plus a damper model withe exponential characteristics" extends Buildings.Fluid.Actuators.BaseClasses.PartialDamperExponential( dp(nominal=dp_nominal), final kFixed=sqrt(kResSqu)); parameter Boolean dp_nominalIncludesDamper = true "set to true if dp_nominal includes the pressure loss of the open damper"; protected parameter Modelica.SIunits.PressureDifference dpDamOpe_nominal(displayUnit="Pa")= k1*m_flow_nominal^2/2/Medium.density(sta_default)/A^2 "Pressure drop of fully open damper at nominal flow rate"; parameter Real kResSqu(unit="kg.m", fixed=false) "Resistance coefficient for fixed resistance element"; initial equation kResSqu = if dp_nominalIncludesDamper then m_flow_nominal^2 / (dp_nominal-dpDamOpe_nominal) else m_flow_nominal^2 / dp_nominal; assert(kResSqu > 0, "Wrong parameters in damper model: dp_nominal < dpDamOpe_nominal" + "\n dp_nominal = " + String(dp_nominal) + "\n dpDamOpe_nominal = " + String(dpDamOpe_nominal)); end VAVBoxExponential;