 
Extends from Modelica.Icons.VariantsPackage (Icon for package containing variants).
| Name | Description | 
|---|---|
|  Exponential | Air damper with exponential opening characteristics | 
|  VAVBoxExponential | VAV box with a fixed resistance plus a damper model withe exponential characteristics | 
|  MixingBox | Outside air mixing box with interlocked air dampers | 
|  MixingBoxMinimumFlow | Outside air mixing box with parallel damper for minimum outside air flow rate | 
 Buildings.Fluid.Actuators.Dampers.Exponential
Buildings.Fluid.Actuators.Dampers.Exponential
 
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 aty=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 blades | single blades | |
| yL | 15/90 | 15/90 | 
| yU | 55/90 | 65/90 | 
| k0 | 1E6 | 1E6 | 
| k1 | 0.2 to 0.5 | 0.2 to 0.5 | 
| a | -1.51 | -1.51 | 
| b | 0.105*90 | 0.0842*90 | 
Extends from Buildings.Fluid.Actuators.BaseClasses.PartialDamperExponential (Partial model for air dampers with exponential opening characteristics).
| Type | Name | Default | Description | 
|---|---|---|---|
| replaceable package Medium | PartialMedium | Medium in the component | |
| MassFlowRate | m_flow_turbulent | if use_deltaM then deltaM*m_... | Turbulent flow if |m_flow| >= m_flow_turbulent [kg/s] | 
| Boolean | use_deltaM | true | Set to true to use deltaM for turbulent transition, else ReC is used | 
| Real | deltaM | 0.3 | Fraction of nominal mass flow rate where transition to turbulent occurs | 
| Boolean | use_v_nominal | true | Set to true to use face velocity to compute area | 
| Velocity | v_nominal | 1 | Nominal face velocity [m/s] | 
| Area | A | m_flow_nominal/rho_nominal/v... | Face area [m2] | 
| Boolean | roundDuct | false | Set to true for round duct, false for square cross section | 
| Real | ReC | 4000 | Reynolds number where transition to turbulent starts | 
| Real | kFixed | 0 | 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 | |||
| MassFlowRate | m_flow_nominal | Nominal mass flow rate [kg/s] | |
| Pressure | dp_nominal | (m_flow_nominal/kDam_nominal... | Pressure drop at nominal mass flow rate [Pa] | 
| Initialization | |||
| MassFlowRate | m_flow.start | 0 | Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction) [kg/s] | 
| Assumptions | |||
| Boolean | allowFlowReversal | system.allowFlowReversal | = true to allow flow reversal, false restricts to design direction (port_a -> port_b) | 
| Advanced | |||
| Boolean | homotopyInitialization | true | = true, use homotopy method | 
| Boolean | from_dp | false | = true, use m_flow = f(dp) else dp = f(m_flow) | 
| Boolean | linearized | false | = true, use linear relation between m_flow and dp for any flow rate | 
| Boolean | use_constant_density | true | Set to true to use constant density for flow friction | 
| Diagnostics | |||
| Boolean | show_V_flow | false | = true, if volume flow rate at inflowing port is computed | 
| Boolean | show_T | false | = true, if actual temperature at port is computed (may lead to events) | 
| Dynamics | |||
| Filtered opening | |||
| Boolean | filteredOpening | true | = true, if opening is filtered with a 2nd order CriticalDamping filter | 
| Time | riseTime | 120 | Rise time of the filter (time to reach 99.6 % of an opening step) [s] | 
| Init | init | Modelica.Blocks.Types.Init.I... | Type of initialization (no init/steady state/initial state/initial output) | 
| Real | y_start | 1 | Initial value of output | 
| Damper coefficients | |||
| Real | a | -1.51 | Coefficient a for damper characteristics | 
| Real | b | 0.105*90 | Coefficient b for damper characteristics | 
| Real | yL | 15/90 | Lower value for damper curve | 
| Real | yU | 55/90 | Upper value for damper curve | 
| Real | k0 | 1E6 | Flow coefficient for y=0, k0 = pressure drop divided by dynamic pressure | 
| Real | k1 | 0.45 | Flow coefficient for y=1, k1 = pressure drop divided by dynamic pressure | 
| 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) | 
| input RealInput | y | Actuator position (0: closed, 1: open) | 
| output RealOutput | y_actual | Actual valve position | 
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),
  final kFixed=0);
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)";
end Exponential;
 
 Buildings.Fluid.Actuators.Dampers.VAVBoxExponential
Buildings.Fluid.Actuators.Dampers.VAVBoxExponential
 
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).
| Type | Name | Default | Description | 
|---|---|---|---|
| replaceable package Medium | PartialMedium | Medium in the component | |
| MassFlowRate | m_flow_turbulent | if use_deltaM then deltaM*m_... | Turbulent flow if |m_flow| >= m_flow_turbulent [kg/s] | 
| Boolean | use_deltaM | true | Set to true to use deltaM for turbulent transition, else ReC is used | 
| Real | deltaM | 0.3 | Fraction of nominal mass flow rate where transition to turbulent occurs | 
| Boolean | use_v_nominal | true | Set to true to use face velocity to compute area | 
| Velocity | v_nominal | 1 | Nominal face velocity [m/s] | 
| Area | A | m_flow_nominal/rho_nominal/v... | Face area [m2] | 
| Boolean | roundDuct | false | Set to true for round duct, false for square cross section | 
| Real | ReC | 4000 | Reynolds number where transition to turbulent starts | 
| Real | kFixed | sqrt(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 | |||
| MassFlowRate | m_flow_nominal | Nominal mass flow rate [kg/s] | |
| Pressure | dp_nominal | Pressure drop at nominal mass flow rate [Pa] | |
| Boolean | dp_nominalIncludesDamper | true | set to true if dp_nominal includes the pressure loss of the open damper | 
| Initialization | |||
| MassFlowRate | m_flow.start | 0 | Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction) [kg/s] | 
| Assumptions | |||
| Boolean | allowFlowReversal | system.allowFlowReversal | = true to allow flow reversal, false restricts to design direction (port_a -> port_b) | 
| Advanced | |||
| Boolean | homotopyInitialization | true | = true, use homotopy method | 
| Boolean | from_dp | false | = true, use m_flow = f(dp) else dp = f(m_flow) | 
| Boolean | linearized | false | = true, use linear relation between m_flow and dp for any flow rate | 
| Boolean | use_constant_density | true | Set to true to use constant density for flow friction | 
| Diagnostics | |||
| Boolean | show_V_flow | false | = true, if volume flow rate at inflowing port is computed | 
| Boolean | show_T | false | = true, if actual temperature at port is computed (may lead to events) | 
| Dynamics | |||
| Filtered opening | |||
| Boolean | filteredOpening | true | = true, if opening is filtered with a 2nd order CriticalDamping filter | 
| Time | riseTime | 120 | Rise time of the filter (time to reach 99.6 % of an opening step) [s] | 
| Init | init | Modelica.Blocks.Types.Init.I... | Type of initialization (no init/steady state/initial state/initial output) | 
| Real | y_start | 1 | Initial value of output | 
| Damper coefficients | |||
| Real | a | -1.51 | Coefficient a for damper characteristics | 
| Real | b | 0.105*90 | Coefficient b for damper characteristics | 
| Real | yL | 15/90 | Lower value for damper curve | 
| Real | yU | 55/90 | Upper value for damper curve | 
| Real | k0 | 1E6 | Flow coefficient for y=0, k0 = pressure drop divided by dynamic pressure | 
| Real | k1 | 0.45 | Flow coefficient for y=1, k1 = pressure drop divided by dynamic pressure | 
| 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) | 
| input RealInput | y | Actuator position (0: closed, 1: open) | 
| output RealOutput | y_actual | Actual valve position | 
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.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 = "       + String(dp_nominal)
          + "\n  dpDamOpe_nominal = " + String(dpDamOpe_nominal));
end VAVBoxExponential;
 
 Buildings.Fluid.Actuators.Dampers.MixingBox
Buildings.Fluid.Actuators.Dampers.MixingBox
 
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), Buildings.Fluid.Actuators.BaseClasses.ActuatorSignal (Partial model that implements the filtered opening for valves and dampers).
| Type | Name | Default | Description | 
|---|---|---|---|
| Boolean | use_deltaM | true | Set to true to use deltaM for turbulent transition, else ReC is used | 
| Real | deltaM | 0.3 | Fraction of nominal mass flow rate where transition to turbulent occurs | 
| Boolean | use_v_nominal | true | Set to true to use face velocity to compute area | 
| Velocity | v_nominal | 1 | Nominal face velocity [m/s] | 
| Boolean | roundDuct | false | Set to true for round duct, false for square cross section | 
| Real | ReC | 4000 | Reynolds number where transition to turbulent starts | 
| Area | AOut | mOut_flow_nominal/rho_nomina... | Face area outside air damper [m2] | 
| Area | AExh | mExh_flow_nominal/rho_nomina... | Face area exhaust air damper [m2] | 
| Area | ARec | mRec_flow_nominal/rho_nomina... | Face area recirculation air damper [m2] | 
| Nominal condition | |||
| Boolean | dp_nominalIncludesDamper | false | set to true if dp_nominal includes the pressure loss of the open damper | 
| MassFlowRate | mOut_flow_nominal | Mass flow rate outside air damper [kg/s] | |
| Pressure | dpOut_nominal | Pressure drop outside air leg [Pa] | |
| MassFlowRate | mRec_flow_nominal | Mass flow rate recirculation air damper [kg/s] | |
| Pressure | dpRec_nominal | Pressure drop recirculation air leg [Pa] | |
| MassFlowRate | mExh_flow_nominal | Mass flow rate exhaust air damper [kg/s] | |
| Pressure | dpExh_nominal | Pressure drop exhaust air leg [Pa] | |
| Dynamics | |||
| Filtered opening | |||
| Boolean | filteredOpening | true | = true, if opening is filtered with a 2nd order CriticalDamping filter | 
| Time | riseTime | 120 | Rise time of the filter (time to reach 99.6 % of an opening step) [s] | 
| Init | init | Modelica.Blocks.Types.Init.I... | Type of initialization (no init/steady state/initial state/initial output) | 
| Real | y_start | 1 | Initial value of output | 
| Assumptions | |||
| Boolean | allowFlowReversal | system.allowFlowReversal | = true to allow flow reversal, false restricts to design direction (port_a -> port_b) | 
| Advanced | |||
| Boolean | from_dp | true | = true, use m_flow = f(dp) else dp = f(m_flow) | 
| Boolean | linearized | false | = true, use linear relation between m_flow and dp for any flow rate | 
| Boolean | use_constant_density | true | Set to true to use constant density for flow friction | 
| Damper coefficients | |||
| Real | a | -1.51 | Coefficient a for damper characteristics | 
| Real | b | 0.105*90 | Coefficient b for damper characteristics | 
| Real | yL | 15/90 | Lower value for damper curve | 
| Real | yU | 55/90 | Upper value for damper curve | 
| Real | k0 | 1E6 | Flow coefficient for y=0, k0 = pressure drop divided by dynamic pressure | 
| Real | k1 | 0.45 | Flow coefficient for y=1, k1 = pressure drop divided by dynamic pressure | 
| Type | Name | Description | 
|---|---|---|
| input RealInput | y | Actuator position (0: closed, 1: open) | 
| output RealOutput | y_actual | Actual valve position | 
| FluidPort_a | port_Out | Fluid connector a (positive design flow direction is from port_a to port_b) | 
| FluidPort_b | port_Exh | Fluid connector b (positive design flow direction is from port_a to port_b) | 
| FluidPort_a | port_Ret | Fluid connector a (positive design flow direction is from port_a to port_b) | 
| FluidPort_b | port_Sup | Fluid connector b (positive design flow direction is from port_a to port_b) | 
model MixingBox "Outside air mixing box with interlocked air dampers"
  extends Buildings.BaseClasses.BaseIconLow;
  extends Buildings.Fluid.Actuators.BaseClasses.ActuatorSignal;
  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,
    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 filteredOpening=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 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,
    a=a,
    b=b,
    yL=yL,
    yU=yU,
    k0=k0,
    k1=k1,
    use_constant_density=use_constant_density,
    allowFlowReversal=allowFlowReversal,
    final filteredOpening=false) "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,
    a=a,
    b=b,
    yL=yL,
    yU=yU,
    k0=k0,
    k1=k1,
    use_constant_density=use_constant_density,
    allowFlowReversal=allowFlowReversal,
    final filteredOpening=false) "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 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); 
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(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
 
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).
| Type | Name | Default | Description | 
|---|---|---|---|
| replaceable package Medium | PartialMedium | Medium in the component | |
| Boolean | use_deltaM | true | Set to true to use deltaM for turbulent transition, else ReC is used | 
| Real | deltaM | 0.3 | Fraction of nominal mass flow rate where transition to turbulent occurs | 
| Boolean | use_v_nominal | true | Set to true to use face velocity to compute area | 
| Velocity | v_nominal | 1 | Nominal face velocity [m/s] | 
| Boolean | roundDuct | false | Set to true for round duct, false for square cross section | 
| Real | ReC | 4000 | Reynolds number where transition to turbulent starts | 
| Area | AOut | mOut_flow_nominal/rho_nomina... | Face area outside air damper [m2] | 
| Area | AExh | mExh_flow_nominal/rho_nomina... | Face area exhaust air damper [m2] | 
| Area | ARec | mRec_flow_nominal/rho_nomina... | Face area recirculation air damper [m2] | 
| Area | AOutMin | Face area minimum outside air damper [m2] | |
| Nominal condition | |||
| Boolean | dp_nominalIncludesDamper | false | set to true if dp_nominal includes the pressure loss of the open damper | 
| MassFlowRate | mOut_flow_nominal | Mass flow rate outside air damper [kg/s] | |
| Pressure | dpOut_nominal | Pressure drop outside air leg [Pa] | |
| MassFlowRate | mRec_flow_nominal | Mass flow rate recirculation air damper [kg/s] | |
| Pressure | dpRec_nominal | Pressure drop recirculation air leg [Pa] | |
| MassFlowRate | mExh_flow_nominal | Mass flow rate exhaust air damper [kg/s] | |
| Pressure | dpExh_nominal | Pressure drop exhaust air leg [Pa] | |
| MassFlowRate | mOutMin_flow_nominal | Mass flow rate minimum outside air damper [kg/s] | |
| Pressure | dpOutMin_nominal | Pressure drop minimum outside air leg [Pa] | |
| Dynamics | |||
| Filtered opening | |||
| Boolean | filteredOpening | true | = true, if opening is filtered with a 2nd order CriticalDamping filter | 
| Time | riseTime | 120 | Rise time of the filter (time to reach 99.6 % of an opening step) [s] | 
| Init | init | Modelica.Blocks.Types.Init.I... | Type of initialization (no init/steady state/initial state/initial output) | 
| Real | y_start | 1 | Initial value of output | 
| Real | yOutMin_start | y_start | Initial value of signal for minimum outside air damper | 
| Assumptions | |||
| Boolean | allowFlowReversal | system.allowFlowReversal | = true to allow flow reversal, false restricts to design direction (port_a -> port_b) | 
| Advanced | |||
| Boolean | from_dp | true | = true, use m_flow = f(dp) else dp = f(m_flow) | 
| Boolean | linearized | false | = true, use linear relation between m_flow and dp for any flow rate | 
| Boolean | use_constant_density | true | Set to true to use constant density for flow friction | 
| Damper coefficients | |||
| Real | a | -1.51 | Coefficient a for damper characteristics | 
| Real | b | 0.105*90 | Coefficient b for damper characteristics | 
| Real | yL | 15/90 | Lower value for damper curve | 
| Real | yU | 55/90 | Upper value for damper curve | 
| Real | k0 | 1E6 | Flow coefficient for y=0, k0 = pressure drop divided by dynamic pressure | 
| Real | k1 | 0.45 | Flow coefficient for y=1, k1 = pressure drop divided by dynamic pressure | 
| Type | Name | Description | 
|---|---|---|
| input RealInput | y | Actuator position (0: closed, 1: open) | 
| output RealOutput | y_actual | Actual valve position | 
| FluidPort_a | port_Out | Fluid connector a (positive design flow direction is from port_a to port_b) | 
| FluidPort_b | port_Exh | Fluid connector b (positive design flow direction is from port_a to port_b) | 
| FluidPort_a | port_Ret | Fluid connector a (positive design flow direction is from port_a to port_b) | 
| FluidPort_b | port_Sup | Fluid connector b (positive design flow direction is from port_a to port_b) | 
| FluidPort_a | port_OutMin | Fluid connector a (positive design flow direction is from port_a to port_b) | 
| input RealInput | yOutMin | Damper position minimum outside air (0: closed, 1: open) | 
| output RealOutput | yOutMin_actual | Actual valve position | 
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";
  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,
    use_v_nominal=use_v_nominal,
    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,
    A=AOutMin,
    final filteredOpening=false) "Damper for minimum outside air intake";
protected 
  Modelica.Blocks.Interfaces.RealOutput yOutMin_filtered if filteredOpening 
    "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 
        filteredOpening 
    "Second order filter to approximate valve opening time, and to improve numerics"; 
equation 
 connect(filterOutMin.y, yOutMin_filtered); 
  if filteredOpening 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);
  connect(yOutMin_actual, yOutMin_actual); 
end MixingBoxMinimumFlow;