Buildings.Fluid.Actuators.Dampers

Package with air damper models

Information


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

Package Content

NameDescription
Buildings.Fluid.Actuators.Dampers.Exponential Exponential Air damper with exponential opening characteristics
Buildings.Fluid.Actuators.Dampers.VAVBoxExponential VAVBoxExponential VAV box with a fixed resistance plus a damper model withe exponential 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.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

  k = 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 k(y) is differentiable in y and the derivative is continuous.

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
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
Booleanuse_v_nominaltrueSet to true to use face velocity to compute area
Velocityv_nominal1Nominal face velocity [m/s]
AreaAm_flow_nominal/rho_nominal/v...Face area [m2]
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]
Pressuredp_nominal(m_flow_nominal/kDam_nominal...Pressure drop at nominal mass flow rate [Pa]
Initialization
MassFlowRatem_flow.start0Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction) [kg/s]
RealkDam.start1Flow coefficient for damper, kDam=m_flow/sqrt(dp)
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_small1E-4*m_flow_nominalSmall mass flow rate for regularization of zero flow [kg/s]
Booleanfrom_dpfalse= 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
Diagnostics
Booleanshow_V_flowfalse= true, if volume flow rate at inflowing port is computed
Booleanshow_Tfalse= true, if actual temperature at port is computed (may lead to events)
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 RealInputyDamper position (0: closed, 1: open)

Modelica definition

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

equation 
  k = kDam "flow coefficient for resistance base model";
end Exponential;

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 resistance 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
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
Booleanuse_v_nominaltrueSet to true to use face velocity to compute area
Velocityv_nominal1Nominal face velocity [m/s]
AreaAm_flow_nominal/rho_nominal/v...Face area [m2]
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]
Pressuredp_nominal Pressure drop at nominal mass flow rate [Pa]
Booleandp_nominalIncludesDampertrueset to true if dp_nominal includes the pressure loss of the open damper
Initialization
MassFlowRatem_flow.start0Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction) [kg/s]
Pressuredp.start0Pressure difference between port_a and port_b [Pa]
RealkDam.start1Flow coefficient for damper, kDam=m_flow/sqrt(dp)
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_small1E-4*m_flow_nominalSmall mass flow rate for regularization of zero flow [kg/s]
Booleanfrom_dpfalse= 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
Diagnostics
Booleanshow_V_flowfalse= true, if volume flow rate at inflowing port is computed
Booleanshow_Tfalse= true, if actual temperature at port is computed (may lead to events)
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 RealInputyDamper position (0: closed, 1: open)

Modelica definition

model VAVBoxExponential 
  "VAV box with a fixed resistance plus a damper model withe exponential characteristics"
  extends Buildings.Fluid.Actuators.BaseClasses.PartialDamperExponential;
  import SI = Modelica.SIunits;

  parameter SI.MassFlowRate m_flow_nominal "Mass flow rate";
  parameter Boolean dp_nominalIncludesDamper = true 
    "set to true if dp_nominal includes the pressure loss of the open damper";

protected 
  parameter SI.Pressure dpDamOpe_nominal = k1*m_flow_nominal^2/2/Medium.density(sta0)/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 = "       + realString(dp_nominal)
          + "\n  dpDamOpe_nominal = " + realString(dpDamOpe_nominal));
equation 
   k = if noEvent(kDam>Modelica.Constants.eps) then sqrt(1/(1/kResSqu + 1/kDam^2)) else 0 
    "flow coefficient for resistance base model";

end VAVBoxExponential;

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.BaseClasses.BaseIconLow (Base icon with model name below the icon).

Parameters

TypeNameDefaultDescription
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
Booleanuse_v_nominaltrueSet to true to use face velocity to compute area
Velocityv_nominal1Nominal face velocity [m/s]
BooleanroundDuctfalseSet to true for round duct, false for square cross section
RealReC4000Reynolds number where transition to turbulent starts
AreaAOutmOut_flow_nominal/rho_nomina...Face area outside air damper [m2]
AreaAExhmExh_flow_nominal/rho_nomina...Face area exhaust air damper [m2]
AreaARecmRec_flow_nominal/rho_nomina...Face area recirculation air damper [m2]
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]
PressuredpOut_nominal Pressure drop outside air leg [Pa]
MassFlowRatemRec_flow_nominal Mass flow rate recirculation air damper [kg/s]
PressuredpRec_nominal Pressure drop recirculation air leg [Pa]
MassFlowRatemExh_flow_nominal Mass flow rate exhaust air damper [kg/s]
PressuredpExh_nominal Pressure drop exhaust air leg [Pa]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_small1E-4*mOut_flow_nominalSmall mass flow rate for regularization of zero flow [kg/s]
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
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)
input RealInputyDamper position (0: closed, 1: open)

Modelica definition

model MixingBox "Outside air mixing box with interlocked air dampers"
  extends Buildings.BaseClasses.BaseIconLow;
  outer Modelica.Fluid.System system "System wide properties";
  replaceable package Medium =
      Modelica.Media.Interfaces.PartialMedium "Medium in the component";
  import Modelica.Constants;

  parameter Boolean allowFlowReversal = system.allowFlowReversal 
    "= true to allow flow reversal, false restricts to design direction (port_a -> port_b)";

  VAVBoxExponential damOA(A=AOut,
    redeclare package Medium = Medium,
    dp_nominal=dpOut_nominal,
    dp_nominalIncludesDamper=dp_nominalIncludesDamper,
    from_dp=from_dp,
    linearized=linearized,
    use_deltaM=use_deltaM,
    deltaM=deltaM,
    use_v_nominal=use_v_nominal,
    v_nominal=v_nominal,
    roundDuct=roundDuct,
    ReC=ReC,
    m_flow_small=m_flow_small,
    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);
  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 Boolean use_v_nominal = true 
    "Set to true to use face velocity to compute area";
  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";

  parameter Modelica.SIunits.Area AOut=mOut_flow_nominal/rho_nominal/v_nominal 
    "Face area outside air damper";
  VAVBoxExponential damExh(                            A=AExh,
    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,
    use_v_nominal=use_v_nominal,
    v_nominal=v_nominal,
    roundDuct=roundDuct,
    ReC=ReC,
    m_flow_small=m_flow_small,
    a=a,
    b=b,
    yL=yL,
    yU=yU,
    k0=k0,
    k1=k1,
    use_constant_density=use_constant_density,
    allowFlowReversal=allowFlowReversal) "Exhaust air damper";
  parameter Modelica.SIunits.Area AExh=mExh_flow_nominal/rho_nominal/v_nominal 
    "Face area exhaust air damper";
  VAVBoxExponential damRec(                            A=ARec,
    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,
    use_v_nominal=use_v_nominal,
    v_nominal=v_nominal,
    roundDuct=roundDuct,
    ReC=ReC,
    m_flow_small=m_flow_small,
    a=a,
    b=b,
    yL=yL,
    yU=yU,
    k0=k0,
    k1=k1,
    use_constant_density=use_constant_density,
    allowFlowReversal=allowFlowReversal) "Recirculation air damper";
  parameter Modelica.SIunits.Area ARec=mRec_flow_nominal/rho_nominal/v_nominal 
    "Face area 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.Pressure 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.Pressure 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.Pressure dpExh_nominal(min=0, displayUnit="Pa") 
    "Pressure drop exhaust air leg";

  parameter Modelica.Media.Interfaces.PartialMedium.MassFlowRate m_flow_small=1E-4
      *mOut_flow_nominal "Small mass flow rate for regularization of zero flow";
  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.Interfaces.RealInput y "Damper position (0: closed, 1: open)";
  Modelica.Blocks.Sources.Constant uni(k=1) "Unity signal";

  Modelica.Blocks.Math.Add add(k2=-1);

protected 
  parameter Medium.Density rho_nominal=Medium.density(sta_nominal) 
    "Density, used to compute fluid volume";
  parameter Medium.ThermodynamicState sta_nominal=
     Medium.setState_pTX(T=Medium.T_default, p=Medium.p_default, X=Medium.X_default);
equation 
  connect(y, damOA.y);
  connect(y, damExh.y);
  connect(uni.y, add.u1);
  connect(y, add.u2);
  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);
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
Booleanuse_v_nominaltrueSet to true to use face velocity to compute area
Velocityv_nominal1Nominal face velocity [m/s]
BooleanroundDuctfalseSet to true for round duct, false for square cross section
RealReC4000Reynolds number where transition to turbulent starts
AreaAOutmOut_flow_nominal/rho_nomina...Face area outside air damper [m2]
AreaAExhmExh_flow_nominal/rho_nomina...Face area exhaust air damper [m2]
AreaARecmRec_flow_nominal/rho_nomina...Face area recirculation air damper [m2]
AreaAOutMin Face area minimum outside air damper [m2]
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]
PressuredpOut_nominal Pressure drop outside air leg [Pa]
MassFlowRatemRec_flow_nominal Mass flow rate recirculation air damper [kg/s]
PressuredpRec_nominal Pressure drop recirculation air leg [Pa]
MassFlowRatemExh_flow_nominal Mass flow rate exhaust air damper [kg/s]
PressuredpExh_nominal Pressure drop exhaust air leg [Pa]
MassFlowRatemOutMin_flow_nominal Mass flow rate minimum outside air damper [kg/s]
PressuredpOutMin_nominal Pressure drop minimum outside air leg [Pa]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_small1E-4*mOut_flow_nominalSmall mass flow rate for regularization of zero flow [kg/s]
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
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)
input RealInputyDamper position (0: closed, 1: open)
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)

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;

  parameter Modelica.SIunits.Area AOutMin 
    "Face area minimum outside air damper";

  parameter Modelica.SIunits.MassFlowRate mOutMin_flow_nominal 
    "Mass flow rate minimum outside air damper";
  parameter Modelica.SIunits.Pressure dpOutMin_nominal(min=0, displayUnit="Pa") 
    "Pressure drop minimum outside air leg";

  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)";

  VAVBoxExponential damOAMin(
    redeclare package Medium = Medium,
    dp_nominalIncludesDamper=dp_nominalIncludesDamper,
    from_dp=from_dp,
    linearized=linearized,
    use_deltaM=use_deltaM,
    deltaM=deltaM,
    use_v_nominal=use_v_nominal,
    v_nominal=v_nominal,
    roundDuct=roundDuct,
    ReC=ReC,
    m_flow_small=m_flow_small,
    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,
    A=AOutMin) "Damper for minimum outside air intake";
equation 
  connect(port_OutMin, damOAMin.port_a);
  connect(damOAMin.port_b, port_Sup);
  connect(yOutMin, damOAMin.y);
end MixingBoxMinimumFlow;

HTML-documentation generated by Dymola Fri Jul 30 18:05:53 2010.