| Name | Description | 
|---|---|
|  der_equalPercentage | Derivative of valve opening characteristics for equal percentage valve | 
|  der_exponentialDamper | Derivative of damper opening characteristics for an exponential damper | 
|  equalPercentage | Valve opening characteristics for equal percentage valve | 
|  Examples | Collection of models that illustrate model use and test models | 
|  exponentialDamper | Damper opening characteristics for an exponential damper | 
|  PartialActuator | Partial model of an actuator | 
|  PartialDamperExponential | Partial model for air dampers with exponential opening characteristics | 
|  PartialThreeWayValve | Partial three way valve | 
|  PartialTwoWayValve | Partial model for a two way valve | 
| ValveParameters | Model with parameters for valves | 
This function computes the derivative of the opening characteristics of an equal percentage valve.
The function is the derivative of TwoWayValveEqualPercentage.
| Type | Name | Default | Description | 
|---|---|---|---|
| Real | y | Valve opening signal, y=1 is fully open | |
| Real | R | Rangeability, R=50...100 typically | |
| Real | l | Valve leakage, l=Cv(y=0)/Cvs | |
| Real | delta | Range of significant deviation from equal percentage law | |
| Real | der_y | Derivative of valve opening signal | 
| Type | Name | Description | 
|---|---|---|
| Real | der_phi | Derivative of ratio actual to nominal mass flow rate, dphi/dy | 
function der_equalPercentage 
  "Derivative of valve opening characteristics for equal percentage valve"
  input Real y "Valve opening signal, y=1 is fully open";
  input Real R "Rangeability, R=50...100 typically";
  input Real l(min=0, max=1) "Valve leakage, l=Cv(y=0)/Cvs";
  input Real delta "Range of significant deviation from equal percentage law";
  input Real der_y "Derivative of valve opening signal";
  output Real der_phi 
    "Derivative of ratio actual to nominal mass flow rate, dphi/dy";
protected 
   Real a "Polynomial coefficient";
   Real b "Polynomial coefficient";
   Real c "Polynomial coefficient";
   Real logR "=log(R)";
   Real z "Auxiliary variable";
   Real q "Auxiliary variable";
   Real p "Auxiliary variable";
algorithm 
  if y < delta/2 then
    der_phi := (R^(delta-1) - l) / delta;
  else
    if (y > (3/2 * delta)) then
      der_phi := R^(y-1)*ln(R);
    else
      logR := Modelica.Math.log(R);
      z := (3*delta/2);
      q := delta*R^z*logR;
      p := R^z;
      a := (q - 2*p + 2*R^delta)/(delta^3*R);
      b := (-5*q + 12*p - 13*R^delta + l*R)/(2*delta^2*R);
      c := (7*q - 18*p + 24*R^delta - 6*l*R)/(4*delta*R);
      der_phi  := c + y * ( 2*b + 3*a*y);
    end if;
  end if;
end der_equalPercentage;
This function computes the derivative of the opening characteristics of an exponential damper.
The function is used by the model Dampers.Exponential.
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.
| Type | Name | Default | Description | 
|---|---|---|---|
| Real | y | Control signal, y=0 is closed, y=1 is open | |
| Real | a | Coefficient a for damper characteristics | |
| Real | b | Coefficient b for damper characteristics | |
| Real | cL[3] | Polynomial coefficients for curve fit for y < yl | |
| Real | cU[3] | Polynomial coefficients for curve fit for y > yu | |
| Real | yL | Lower value for damper curve | |
| Real | yU | Upper value for damper curve | |
| Real | der_y | Derivative of control signal | 
| Type | Name | Description | 
|---|---|---|
| Real | der_kTheta | Derivative of flow coefficient, der_kTheta=dkTheta/dy | 
function der_exponentialDamper 
  "Derivative of damper opening characteristics for an exponential damper"
  input Real y(unit="") "Control signal, y=0 is closed, y=1 is open";
  input Real a(unit="") "Coefficient a for damper characteristics";
  input Real b(unit="") "Coefficient b for damper characteristics";
  input Real[3] cL "Polynomial coefficients for curve fit for y < yl";
  input Real[3] cU "Polynomial coefficients for curve fit for y > yu";
  input Real yL "Lower value for damper curve";
  input Real yU "Upper value for damper curve";
  input Real der_y(unit="") "Derivative of control signal";
  output Real der_kTheta(min=0) 
    "Derivative of flow coefficient, der_kTheta=dkTheta/dy";
algorithm 
  if y < yL then
    der_kTheta := exp(cL[3] + y * (cL[2] + y * cL[1]))*(2 * cL[1] * y + cL[2]);
  else
    if (y > yU) then
      der_kTheta := exp(cU[3] + y * (cU[2] + y * cU[1]))*(2 * cU[1] * y + cU[2]);
    else
      der_kTheta := -b*exp(a+b*(1-y)) 
        "y=0 is closed, but theta=1 is closed in ASHRAE-825";
    end if;
  end if;
end der_exponentialDamper;
This function computes the opening characteristics of an equal percentage valve.
The function is used by the model TwoWayValveEqualPercentage.
For y < delta/2, the valve characteristics is linear. For y > 3*delta/2 the valve characteristics is equal percentage. In between, a cubic spline is used to ensure that the valve characteristics is once continuously differentiable with respect to y.
| Type | Name | Default | Description | 
|---|---|---|---|
| Real | y | Valve opening signal, y=1 is fully open | |
| Real | R | Rangeability, R=50...100 typically | |
| Real | l | Valve leakage, l=Cv(y=0)/Cvs | |
| Real | delta | Range of significant deviation from equal percentage law | 
| Type | Name | Description | 
|---|---|---|
| Real | phi | Ratio actual to nominal mass flow rate, phi=Cv(y)/Cv(y=1) | 
function equalPercentage 
  "Valve opening characteristics for equal percentage valve"
  annotation(derivative=der_equalPercentage);
  input Real y "Valve opening signal, y=1 is fully open";
  input Real R "Rangeability, R=50...100 typically";
  input Real l(min=0, max=1) "Valve leakage, l=Cv(y=0)/Cvs";
  input Real delta "Range of significant deviation from equal percentage law";
  output Real phi "Ratio actual to nominal mass flow rate, phi=Cv(y)/Cv(y=1)";
protected 
   Real a "Polynomial coefficient";
   Real b "Polynomial coefficient";
   Real c "Polynomial coefficient";
   Real d "Polynomial coefficient";
   Real logR "=log(R)";
   Real z "Auxiliary variable";
   Real q "Auxiliary variable";
   Real p "Auxiliary variable";
algorithm 
  if y < delta/2 then
    phi := l + y * (R^(delta-1) - l) / delta;
  else
    if (y > (3/2 * delta)) then
      phi := R^(y-1);
    else
      logR := Modelica.Math.log(R);
      z := (3*delta/2);
      q := delta*R^z*logR;
      p := R^z;
      a := (q - 2*p + 2*R^delta)/(delta^3*R);
      b := (-5*q + 12*p - 13*R^delta + l*R)/(2*delta^2*R);
      c := (7*q - 18*p + 24*R^delta - 6*l*R)/(4*delta*R);
      d := (-3*q + 8*p - 9*R^delta + 9*l*R)/(8*R);
      phi  := d + y * ( c + y * ( b + y * a));
    end if;
  end if;
end equalPercentage;
This function computes the opening characteristics of an exponential damper.
The function is used by the model Dampers.Exponential.
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.
| Type | Name | Default | Description | 
|---|---|---|---|
| Real | y | Control signal, y=0 is closed, y=1 is open | |
| Real | a | Coefficient a for damper characteristics | |
| Real | b | Coefficient b for damper characteristics | |
| Real | cL[3] | Polynomial coefficients for curve fit for y < yl | |
| Real | cU[3] | Polynomial coefficients for curve fit for y > yu | |
| Real | yL | Lower value for damper curve | |
| Real | yU | Upper value for damper curve | 
| Type | Name | Description | 
|---|---|---|
| Real | kTheta | Flow coefficient, kTheta = pressure drop divided by dynamic pressure | 
function exponentialDamper 
  "Damper opening characteristics for an exponential damper"
  annotation(derivative=der_exponentialDamper);
  input Real y(unit="") "Control signal, y=0 is closed, y=1 is open";
  input Real a(unit="") "Coefficient a for damper characteristics";
  input Real b(unit="") "Coefficient b for damper characteristics";
  input Real[3] cL "Polynomial coefficients for curve fit for y < yl";
  input Real[3] cU "Polynomial coefficients for curve fit for y > yu";
  input Real yL "Lower value for damper curve";
  input Real yU "Upper value for damper curve";
  output Real kTheta(min=0) 
    "Flow coefficient, kTheta = pressure drop divided by dynamic pressure";
algorithm 
  if y < yL then
    kTheta := exp(cL[3] + y * (cL[2] + y * cL[1]));
  else
    if (y > yU) then
      kTheta := exp(cU[3] + y * (cU[2] + y * cU[1]));
    else
      kTheta := exp(a+b*(1-y)) "y=0 is closed";
    end if;
  end if;
end exponentialDamper;
 Buildings.Fluids.Actuators.BaseClasses.PartialActuator
Buildings.Fluids.Actuators.BaseClasses.PartialActuator
 
| Type | Name | Default | Description | 
|---|---|---|---|
| replaceable package Medium | Modelica.Media.Interfaces.Pa... | Medium in the component | |
| Nominal condition | |||
| MassFlowRate | m0_flow | Nominal mass flow rate [kg/s] | |
| Assumptions | |||
| Boolean | allowFlowReversal | system.allowFlowReversal | = true to allow flow reversal, false restricts to design direction (port_a -> port_b) | 
| Advanced | |||
| AbsolutePressure | dp_start | 0.01*system.p_start | Guess value of dp = port_a.p - port_b.p [Pa] | 
| MassFlowRate | m_flow_start | system.m_flow_start | Guess value of m_flow = port_a.m_flow [kg/s] | 
| MassFlowRate | m_flow_small | 1E-4*m0_flow | Small mass flow rate for regularization of zero flow [kg/s] | 
| 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 | 
| Diagnostics | |||
| Boolean | show_T | true | = true, if temperatures at port_a and port_b are computed | 
| Boolean | show_V_flow | true | = true, if volume flow rate at inflowing 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) | 
| input RealInput | y | Damper position (0: closed, 1: open) | 
partial model PartialActuator "Partial model of an actuator"
  extends Buildings.Fluids.BaseClasses.PartialResistance;
  import SI = Modelica.SIunits;
public 
  Modelica.Blocks.Interfaces.RealInput y(min=0, max=1) 
    "Damper position (0: closed, 1: open)";
end PartialActuator;
 Buildings.Fluids.Actuators.BaseClasses.PartialDamperExponential
Buildings.Fluids.Actuators.BaseClasses.PartialDamperExponential
 
Partial model for air dampers with exponential opening characteristics. This is the base model for air dampers. The model defines the flow rate where the linearization near the origin occurs. The model also defines parameters that are used by different air damper models.
This model does not assign k=kDam because the model VAVBoxExponential consists of a fixed resistance and a resistance due to the air damper. If k would be assigned here, then this partial model could not be used as a base class for VAVBoxExponential.
For a description of the opening characteristics and typical parameter values, see the damper model Exponential.
Extends from PartialActuator (Partial model of an actuator).
| Type | Name | Default | Description | 
|---|---|---|---|
| replaceable package Medium | Modelica.Media.Interfaces.Pa... | Medium in the component | |
| Area | A | 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 laminar starts | 
| 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 | 
| Nominal condition | |||
| MassFlowRate | m0_flow | Nominal mass flow rate [kg/s] | |
| Assumptions | |||
| Boolean | allowFlowReversal | system.allowFlowReversal | = true to allow flow reversal, false restricts to design direction (port_a -> port_b) | 
| Advanced | |||
| AbsolutePressure | dp_start | 0.01*system.p_start | Guess value of dp = port_a.p - port_b.p [Pa] | 
| MassFlowRate | m_flow_start | system.m_flow_start | Guess value of m_flow = port_a.m_flow [kg/s] | 
| MassFlowRate | m_flow_small | 1E-4*m0_flow | Small mass flow rate for regularization of zero flow [kg/s] | 
| 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 | 
| Diagnostics | |||
| Boolean | show_T | true | = true, if temperatures at port_a and port_b are computed | 
| Boolean | show_V_flow | true | = true, if volume flow rate at inflowing 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) | 
| input RealInput | y | Damper position (0: closed, 1: open) | 
partial model PartialDamperExponential 
  "Partial model for air dampers with exponential opening characteristics"
  extends PartialActuator;
  parameter Modelica.SIunits.Area A "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 laminar starts";
  parameter Real a(unit="")=-1.51 "Coefficient a for damper characteristics";
  parameter Real b(unit="")=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(min=0) = 1E6 
    "Flow coefficient for y=0, k0 = pressure drop divided by dynamic pressure";
  parameter Real k1(min=0) = 0.45 
    "Flow coefficient for y=1, k1 = pressure drop divided by dynamic pressure";
  Real kDam(start=1) 
    "Flow coefficient for damper, k=m_flow/sqrt(dp), with unit=(kg*m)^(1/2)";
protected 
  Real kTheta(min=0) 
    "Flow coefficient, kTheta = pressure drop divided by dynamic pressure";
  parameter Real[3] cL(fixed=false) 
    "Polynomial coefficients for curve fit for y < yl";
  parameter Real[3] cU(fixed=false) 
    "Polynomial coefficients for curve fit for y > yu";
protected 
  parameter Real facRouDuc= if roundDuct then sqrt(Modelica.Constants.pi)/2 else 1;
initial equation 
 cL[1] = (ln(k0) - b - a)/yL^2;
 cL[2] = (-b*yL - 2*ln(k0) + 2*b + 2*a)/yL;
 cL[3] = ln(k0);
 cU[1] = (ln(k1) - a)/(yU^2 - 2*yU + 1);
 cU[2] = (-b*yU^2 - 2*ln(k1)*yU - (-2*b - 2*a)*yU - b)/(yU^2 - 2*yU + 1);
 cU[3] = (ln(k1)*yU^2 + b*yU^2 + (-2*b - 2*a)*yU + b + a)/(yU^2 - 2*yU + 1);
 assert(k0 > k1, "k0 must be bigger than k1.");
equation 
   m_flow_laminar=eta0*ReC*sqrt(A)*facRouDuc;
   kTheta = exponentialDamper(y=y, a=a, b=b, cL=cL, cU=cU, yL=yL, yU=yU) 
    "y=0 is closed";
   assert(kTheta>=0, "Flow coefficient must not be negative");
   kDam = sqrt(2*Medium.density(state_a)/kTheta) * A 
    "flow coefficient for resistance base model, kDam=k=m_flow/sqrt(dp)";
end PartialDamperExponential;
 Buildings.Fluids.Actuators.BaseClasses.PartialThreeWayValve
Buildings.Fluids.Actuators.BaseClasses.PartialThreeWayValve
 
Partial model of a three way valve. This is the base model for valves with different opening characteristics, such as linear, equal percentage or quick opening. The three way valve model consists of a mixer where valves are placed in two of the flow legs. The third flow leg has no friction. The flow coefficient Kv_SI for flow from port_1 -> port_2 is a parameter and the flow coefficient for flow from port_3 -> port_2 is computed as
         Kv_SI(port_1 -> port_2)
  fraK = ----------------------
         Kv_SI(port_3 -> port_2)
 
where fraK is a parameter.
Since this model uses two way valves to construct a three way valve, see PartialTwoWayValve for details regarding the valve implementation.
Extends from Buildings.Fluids.BaseClasses.PartialThreeWayResistance (Flow splitter with partial resistance model at each port), Buildings.Fluids.Actuators.BaseClasses.ValveParameters (Model with parameters for valves).
| Type | Name | Default | Description | 
|---|---|---|---|
| replaceable package Medium | Modelica.Media.Interfaces.Pa... | Fluid medium model | |
| PartialTwoPortTransport | res1 | redeclare Modelica_Fluid.Int... | Partial model, to be replaced with a fluid component | 
| PartialTwoPortTransport | res3 | redeclare Modelica_Fluid.Int... | Partial model, to be replaced with a fluid component | 
| Real | fraK | 0.7 | Fraction Kv_SI(port_1->port_2)/Kv_SI(port_3->port_2) | 
| Real | l[2] | {0,0} | Valve leakage, l=Cv(y=0)/Cvs | 
| Flow Coefficient | |||
| CvTypes | CvData | Buildings.Fluids.Types.CvTyp... | Selection of flow coefficient | 
| Real | Kv | 0 | Kv (metric) flow coefficient [m3/h/(bar)^(1/2)] | 
| Real | Cv | 0 | Cv (US) flow coefficient [USG/min/(psi)^(1/2)] | 
| Area | Av | 0 | Av (metric) flow coefficient [m2] | 
| Real | Kv_SI | m0_flow/sqrt(dp0) | Flow coefficient for fully open valve in SI units, Kv=m_flow/sqrt(dp) [kg/s/(Pa)^(1/2)] | 
| Pressure-flow linearization | |||
| Real | deltaM | 0.02 | Fraction of nominal flow rate where linearization starts, if y=1 | 
| Nominal condition | |||
| MassFlowRate | m0_flow | Nominal mass flow rate [kg/s] | |
| Pressure | dp0 | 6000 | Nominal pressure drop [Pa] | 
| Advanced | |||
| Boolean | from_dp | true | = true, use m_flow = f(dp) else dp = f(m_flow) | 
| Boolean | linearized[2] | {false,false} | = true, use linear relation between m_flow and dp for any flow rate | 
| Nominal condition | |||
| Density | rhoStd | Medium.density_pTX(101325, 2... | Inlet density for which valve coefficients are defined [kg/m3] | 
| Type | Name | Description | 
|---|---|---|
| FluidPort_a | port_1 | |
| FluidPort_b | port_2 | |
| FluidPort_a | port_3 | |
| input RealInput | y | Damper position (0: closed, 1: open) | 
partial model PartialThreeWayValve "Partial three way valve"
    extends Buildings.Fluids.BaseClasses.PartialThreeWayResistance(
        redeclare FixedResistances.LosslessPipe res2(
            redeclare package Medium = Medium));
    extends Buildings.Fluids.Actuators.BaseClasses.ValveParameters(
      rhoStd=Medium.density_pTX(101325, 273.15+4, Medium.X_default));
  parameter Real fraK(min=0, max=1) = 0.7 
    "Fraction Kv_SI(port_1->port_2)/Kv_SI(port_3->port_2)";
  parameter Real[2] l(min=0, max=1) = {0, 0} "Valve leakage, l=Cv(y=0)/Cvs";
  parameter Real deltaM = 0.02 
    "Fraction of nominal flow rate where linearization starts, if y=1";
  parameter Medium.MassFlowRate m0_flow(min=0) "Nominal mass flow rate";
  parameter Modelica.SIunits.Pressure dp0 = 6000 "Nominal pressure drop";
  parameter Boolean[2] linearized = {false, false} 
    "= true, use linear relation between m_flow and dp for any flow rate";
  Modelica.Blocks.Interfaces.RealInput y "Damper position (0: closed, 1: open)";
protected 
  Modelica.Blocks.Math.Feedback inv "Inversion of control signal";
  Modelica.Blocks.Sources.Constant uni(final k=1) 
    "Outputs one for bypass valve";
equation 
  connect(uni.y, inv.u1);
end PartialThreeWayValve;
 Buildings.Fluids.Actuators.BaseClasses.PartialTwoWayValve
Buildings.Fluids.Actuators.BaseClasses.PartialTwoWayValve
 
Partial model for a two way valve. This is the base model for valves with different opening characteristics, such as linear, equal percentage or quick opening.
Modelling options
The following options have been adapted from the valve implementation in Modelica_Fluid and are described in Buildings.Fluids.Actuators.BaseClasses.ValveParameters.
In contrast to the model in Modelica_Fluid, this model uses the parameter Kv_SI, which is the flow coefficient in SI units, i.e., it is the ratio between mass flow rate in kg/s and square root of pressure drop in Pa.
To prevent the derivative d/dP (m_flow) to be infinite near the origin, this model linearizes the pressure drop vs. flow relation ship. The region in which it is linearized is parameterized by
m_small_flow = deltaM * Kv_SI * sqrt(dp0)Because the parameterization contains Kv_SI, the values for deltaM and dp0 need not be changed if the valve size changes.
The two way valve models are implemented using this partial model, as opposed to using different functions for the valve opening characteristics, because each valve opening characteristics has different parameters.
Extends from Buildings.Fluids.Actuators.BaseClasses.PartialActuator (Partial model of an actuator), Buildings.Fluids.Actuators.BaseClasses.ValveParameters (Model with parameters for valves).
| Type | Name | Default | Description | 
|---|---|---|---|
| replaceable package Medium | Modelica.Media.Interfaces.Pa... | Medium in the component | |
| Real | l | 0.005 | Valve leakage, l=Cv(y=0)/Cvs | 
| Nominal condition | |||
| MassFlowRate | m0_flow | Nominal mass flow rate [kg/s] | |
| Pressure | dp0 | 6000 | Nominal pressure drop [Pa] | 
| Flow Coefficient | |||
| CvTypes | CvData | Buildings.Fluids.Types.CvTyp... | Selection of flow coefficient | 
| Real | Kv | 0 | Kv (metric) flow coefficient [m3/h/(bar)^(1/2)] | 
| Real | Cv | 0 | Cv (US) flow coefficient [USG/min/(psi)^(1/2)] | 
| Area | Av | 0 | Av (metric) flow coefficient [m2] | 
| Real | Kv_SI | m0_flow/sqrt(dp0) | Flow coefficient for fully open valve in SI units, Kv=m_flow/sqrt(dp) [kg/s/(Pa)^(1/2)] | 
| Pressure-flow linearization | |||
| Real | deltaM | 0.02 | Fraction of nominal flow rate where linearization starts, if y=1 | 
| Assumptions | |||
| Boolean | allowFlowReversal | system.allowFlowReversal | = true to allow flow reversal, false restricts to design direction (port_a -> port_b) | 
| Advanced | |||
| AbsolutePressure | dp_start | 0.01*system.p_start | Guess value of dp = port_a.p - port_b.p [Pa] | 
| MassFlowRate | m_flow_start | system.m_flow_start | Guess value of m_flow = port_a.m_flow [kg/s] | 
| MassFlowRate | m_flow_small | 1E-4*m0_flow | Small mass flow rate for regularization of zero flow [kg/s] | 
| 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 | 
| Diagnostics | |||
| Boolean | show_T | true | = true, if temperatures at port_a and port_b are computed | 
| Boolean | show_V_flow | true | = true, if volume flow rate at inflowing port is computed | 
| Nominal condition | |||
| Density | rhoStd | Medium.density_pTX(101325, 2... | Inlet density for which valve coefficients are defined [kg/m3] | 
| 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 | Damper position (0: closed, 1: open) | 
partial model PartialTwoWayValve "Partial model for a two way valve"
  extends Buildings.Fluids.Actuators.BaseClasses.PartialActuator;
  extends Buildings.Fluids.Actuators.BaseClasses.ValveParameters(
      rhoStd=Medium.density_pTX(101325, 273.15+4, Medium.X_default));
  parameter Real l(min=0, max=1) = 0.005 "Valve leakage, l=Cv(y=0)/Cvs";
  Real phi "Ratio actual to nominal mass flow rate, phi=Cv(y)/Cv(y=1)";
equation 
 m_flow_laminar = deltaM * m0_flow;
 k = phi * Kv_SI;
end PartialTwoWayValve;
Model that computes the flow coefficients of valves. This base class allows the following modeling options, which have been adapted from the valve implementation in Modelica_Fluid to specify the valve flow coefficient in fully open conditions:
The treatment of parameters Kv and Cv is explained in detail in the Users Guide.
In contrast to the model in Modelica_Fluid, this model uses the parameter Kv_SI, which is the flow coefficient in SI units, i.e., it is the ratio between mass flow rate in kg/s and square root of pressure drop in Pa.
| Type | Name | Default | Description | 
|---|---|---|---|
| Flow Coefficient | |||
| CvTypes | CvData | Buildings.Fluids.Types.CvTyp... | Selection of flow coefficient | 
| Real | Kv | 0 | Kv (metric) flow coefficient [m3/h/(bar)^(1/2)] | 
| Real | Cv | 0 | Cv (US) flow coefficient [USG/min/(psi)^(1/2)] | 
| Area | Av | 0 | Av (metric) flow coefficient [m2] | 
| Real | Kv_SI | m0_flow/sqrt(dp0) | Flow coefficient for fully open valve in SI units, Kv=m_flow/sqrt(dp) [kg/s/(Pa)^(1/2)] | 
| Pressure-flow linearization | |||
| Real | deltaM | 0.02 | Fraction of nominal flow rate where linearization starts, if y=1 | 
| Nominal condition | |||
| MassFlowRate | m0_flow | Nominal mass flow rate [kg/s] | |
| Pressure | dp0 | 6000 | Nominal pressure drop [Pa] | 
| Advanced | |||
| Nominal condition | |||
| Density | rhoStd | Inlet density for which valve coefficients are defined [kg/m3] | |
partial model ValveParameters "Model with parameters for valves"
  parameter Buildings.Fluids.Types.CvTypes CvData=Buildings.Fluids.Types.CvTypes.OpPoint 
    "Selection of flow coefficient";
  parameter Real Kv(
    fixed= if CvData==Buildings.Fluids.Types.CvTypes.Kv then true else false) = 0 
    "Kv (metric) flow coefficient [m3/h/(bar)^(1/2)]";
  parameter Real Cv(
    fixed= if CvData==Buildings.Fluids.Types.CvTypes.Cv then true else false) = 0 
    "Cv (US) flow coefficient [USG/min/(psi)^(1/2)]";
  parameter Modelica.SIunits.Area Av(
    fixed= if CvData==Buildings.Fluids.Types.CvTypes.Av then true else false) = 0 
    "Av (metric) flow coefficient";
  parameter Real deltaM = 0.02 
    "Fraction of nominal flow rate where linearization starts, if y=1";
  parameter Modelica.SIunits.MassFlowRate m0_flow(min=0) 
    "Nominal mass flow rate";
  parameter Modelica.SIunits.Pressure dp0 = 6000 "Nominal pressure drop";
  parameter Real Kv_SI(
    min=0,
    fixed= if CvData==Buildings.Fluids.Types.CvTypes.OpPoint then true else false,
    start=m0_flow/sqrt(dp0)) = m0_flow/sqrt(dp0) 
    "Flow coefficient for fully open valve in SI units, Kv=m_flow/sqrt(dp) [kg/s/(Pa)^(1/2)]";
  parameter Modelica.SIunits.Density rhoStd 
    "Inlet density for which valve coefficients are defined";
initial equation 
/*
  if CvData == Buildings.Fluids.Types.CvTypes.Av then
    Kv_SI = Av * sqrt(rhoStd);
  elseif CvData == Buildings.Fluids.Types.CvTypes.Kv then
    Kv_SI = Kv*rhoStd/3600/sqrt(1E5) 
      "Unit conversion m3/(h*sqrt(bar)) to kg/(s*sqrt(Pa))";
  elseif CvData == Buildings.Fluids.Types.CvTypes.Cv then
    Kv_SI = Cv*rhoStd*0.0631/1000/sqrt(6895) 
      "Unit conversion USG/(min*sqrt(psi)) to kg/(s*sqrt(Pa))";
 end if;
*/
    Kv_SI = Av * sqrt(rhoStd);
    Kv_SI = Kv*rhoStd/3600/sqrt(1E5) 
    "Unit conversion m3/(h*sqrt(bar)) to kg/(s*sqrt(Pa))";
    Kv_SI = Cv*rhoStd*0.0631/1000/sqrt(6895) 
    "Unit conversion USG/(min*sqrt(psi)) to kg/(s*sqrt(Pa))";
end ValveParameters;