Buildings.Fluid.Actuators.BaseClasses

Package with base classes for actuator models

Information

Extends from Modelica.Fluid.Icons.BaseClassLibrary (Icon for library).

Package Content

NameDescription
Buildings.Fluid.Actuators.BaseClasses.der_equalPercentage der_equalPercentage Derivative of valve opening characteristics for equal percentage valve
Buildings.Fluid.Actuators.BaseClasses.der_exponentialDamper der_exponentialDamper Derivative of damper opening characteristics for an exponential damper
Buildings.Fluid.Actuators.BaseClasses.equalPercentage equalPercentage Valve opening characteristics for equal percentage valve
Buildings.Fluid.Actuators.BaseClasses.exponentialDamper exponentialDamper Damper opening characteristics for an exponential damper
Buildings.Fluid.Actuators.BaseClasses.PartialActuator PartialActuator Partial model of an actuator
Buildings.Fluid.Actuators.BaseClasses.PartialDamperExponential PartialDamperExponential Partial model for air dampers with exponential opening characteristics
Buildings.Fluid.Actuators.BaseClasses.PartialThreeWayValve PartialThreeWayValve Partial three way valve
Buildings.Fluid.Actuators.BaseClasses.PartialTwoWayValve PartialTwoWayValve Partial model for a two way valve
ValveParameters Model with parameters for valves
Buildings.Fluid.Actuators.BaseClasses.Examples Examples Collection of models that illustrate model use and test models


Buildings.Fluid.Actuators.BaseClasses.der_equalPercentage

Derivative of valve opening characteristics for equal percentage valve

Information


This function computes the derivative of the opening characteristics of an equal percentage valve.

The function is the derivative of TwoWayValveEqualPercentage.

Inputs

TypeNameDefaultDescription
Realy Valve opening signal, y=1 is fully open
RealR Rangeability, R=50...100 typically
Reall Valve leakage, l=Cv(y=0)/Cvs
Realdelta Range of significant deviation from equal percentage law
Realder_y Derivative of valve opening signal
Realder_R  
Realder_l  
Realder_delta  

Outputs

TypeNameDescription
Realder_phiDerivative of ratio actual to nominal mass flow rate, dphi/dy

Modelica definition

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";
  input Real der_R;
  input Real der_l;
  input Real der_delta;
  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 * der_y;
  else
    if (y > (3/2 * delta)) then
      der_phi := R^(y-1)*ln(R) * der_y;
    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)) * der_y;
    end if;
  end if;
end der_equalPercentage;

Buildings.Fluid.Actuators.BaseClasses.der_exponentialDamper

Derivative of damper opening characteristics for an exponential damper

Information


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.

Inputs

TypeNameDefaultDescription
Realy Control signal, y=0 is closed, y=1 is open
Reala Coefficient a for damper characteristics
Realb Coefficient b for damper characteristics
RealcL[3] Polynomial coefficients for curve fit for y < yl
RealcU[3] Polynomial coefficients for curve fit for y > yu
RealyL Lower value for damper curve
RealyU Upper value for damper curve
Realder_y Derivative of control signal
Realder_a  
Realder_b  
Realder_cL[3]  
Realder_cU[3]  
Realder_yL  
Realder_yU  

Outputs

TypeNameDescription
Realder_kThetaDerivative of flow coefficient, der_kTheta=dkTheta/dy

Modelica definition

function der_exponentialDamper 
  "Derivative of damper opening characteristics for an exponential damper"

  input Real y "Control signal, y=0 is closed, y=1 is open";
  input Real a "Coefficient a for damper characteristics";
  input Real b "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 "Derivative of control signal";
  input Real der_a;
  input Real der_b;
  input Real[3] der_cL;
  input Real[3] der_cU;
  input Real der_yL;
  input Real der_yU;

  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]) * der_y;
  else
    if (y > yU) then
      der_kTheta := exp(cU[3] + y * (cU[2] + y * cU[1]))*(2 * cU[1] * y + cU[2]) * der_y;
    else
      der_kTheta := -b*exp(a+b*(1-y)) * der_y 
        "y=0 is closed, but theta=1 is closed in ASHRAE-825";
    end if;
  end if;
end der_exponentialDamper;

Buildings.Fluid.Actuators.BaseClasses.equalPercentage

Valve opening characteristics for equal percentage valve

Information


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.

Inputs

TypeNameDefaultDescription
Realy Valve opening signal, y=1 is fully open
RealR Rangeability, R=50...100 typically
Reall Valve leakage, l=Cv(y=0)/Cvs
Realdelta Range of significant deviation from equal percentage law

Outputs

TypeNameDescription
RealphiRatio actual to nominal mass flow rate, phi=Cv(y)/Cv(y=1)

Modelica definition

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;

Buildings.Fluid.Actuators.BaseClasses.exponentialDamper

Damper opening characteristics for an exponential damper

Information


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.

Inputs

TypeNameDefaultDescription
Realy Control signal, y=0 is closed, y=1 is open
Reala Coefficient a for damper characteristics
Realb Coefficient b for damper characteristics
RealcL[3] Polynomial coefficients for curve fit for y < yl
RealcU[3] Polynomial coefficients for curve fit for y > yu
RealyL Lower value for damper curve
RealyU Upper value for damper curve

Outputs

TypeNameDescription
RealkThetaFlow coefficient, kTheta = pressure drop divided by dynamic pressure

Modelica definition

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.Fluid.Actuators.BaseClasses.PartialActuator Buildings.Fluid.Actuators.BaseClasses.PartialActuator

Partial model of an actuator

Buildings.Fluid.Actuators.BaseClasses.PartialActuator

Information


Partial actuator that is the base class for dampers and two way valves.

Extends from Buildings.Fluid.BaseClasses.PartialResistance (Partial model for a hydraulic resistance).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
Pressuredp_nominal Pressure [Pa]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
AbsolutePressuredp_start0.01*system.p_startGuess value of dp = port_a.p - port_b.p [Pa]
MassFlowRatem_flow_startsystem.m_flow_startGuess value of m_flow = port_a.m_flow [kg/s]
MassFlowRatem_flow_small1E-4*m_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
Diagnostics
Booleanshow_Ttrue= true, if temperatures at port_a and port_b are computed
Booleanshow_V_flowtrue= true, if volume flow rate at inflowing port is computed

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

partial model PartialActuator "Partial model of an actuator"
  extends Buildings.Fluid.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.Fluid.Actuators.BaseClasses.PartialDamperExponential Buildings.Fluid.Actuators.BaseClasses.PartialDamperExponential

Partial model for air dampers with exponential opening characteristics

Buildings.Fluid.Actuators.BaseClasses.PartialDamperExponential

Information


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).

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]
AreaA 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 [Pa]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
AbsolutePressuredp_start0.01*system.p_startGuess value of dp = port_a.p - port_b.p [Pa]
MassFlowRatem_flow_startsystem.m_flow_startGuess value of m_flow = port_a.m_flow [kg/s]
MassFlowRatem_flow_small1E-4*m_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
Diagnostics
Booleanshow_Ttrue= true, if temperatures at port_a and port_b are computed
Booleanshow_V_flowtrue= true, if volume flow rate at inflowing port is computed
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

partial model PartialDamperExponential 
  "Partial model for air dampers with exponential opening characteristics"
  extends PartialActuator;

  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 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 turbulent 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 kDamSqu(start=1, unit="kg.m") 
    "Flow coefficient for damper, kDam=k^2=m_flow^2/|dp|";

  parameter Boolean use_constant_density=true 
    "Set to true to use constant density for flow friction";
  Medium.Density rho "Medium density";
protected 
  parameter Medium.Density rho_nominal=Medium.density(sta0) 
    "Density, used to compute fluid volume";

  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;
  parameter Modelica.SIunits.Area area=
     if use_v_nominal then m_flow_nominal/rho_nominal/v_nominal else A 
    "Face velocity used in the computation";
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 
   rho = if use_constant_density then rho_nominal else Medium.density(state_a);
   m_flow_turbulent=if use_deltaM then deltaM * m_flow_nominal else eta_nominal*ReC*sqrt(area)*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");
   kDamSqu = 2*rho/kTheta * area * area 
    "flow coefficient for resistance base model, kDamSqu=k*k=m_flow*m_flow/dp";

end PartialDamperExponential;

Buildings.Fluid.Actuators.BaseClasses.PartialThreeWayValve Buildings.Fluid.Actuators.BaseClasses.PartialThreeWayValve

Partial three way valve

Buildings.Fluid.Actuators.BaseClasses.PartialThreeWayValve

Information


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.Fluid.BaseClasses.PartialThreeWayResistance (Flow splitter with partial resistance model at each port), Buildings.Fluid.Actuators.BaseClasses.ValveParameters (Model with parameters for valves).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumFluid medium model
PartialTwoPortTransportres1redeclare Modelica.Fluid.Int...Partial model, to be replaced with a fluid component
PartialTwoPortTransportres3redeclare Modelica.Fluid.Int...Partial model, to be replaced with a fluid component
RealfraK0.7Fraction Kv_SI(port_1->port_2)/Kv_SI(port_3->port_2)
Reall[2]{0,0}Valve leakage, l=Cv(y=0)/Cvs
Flow Coefficient
CvTypesCvDataBuildings.Fluid.Types.CvType...Selection of flow coefficient
RealKv0Kv (metric) flow coefficient [m3/h/(bar)^(1/2)]
RealCv0Cv (US) flow coefficient [USG/min/(psi)^(1/2)]
AreaAv0Av (metric) flow coefficient [m2]
RealKv_SIm_flow_nominal/sqrt(dpVal_no...Flow coefficient for fully open valve in SI units, Kv=m_flow/sqrt(dp) [kg/s/(Pa)^(1/2)]
Pressure-flow linearization
RealdeltaM0.02Fraction of nominal flow rate where linearization starts, if y=1
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
PressuredpVal_nominaldp_nominalNominal pressure drop [Pa]
Pressuredp_nominal6000Nominal pressure drop [Pa]
Advanced
Booleanfrom_dptrue= true, use m_flow = f(dp) else dp = f(m_flow)
Booleanlinearized[2]{false,false}= true, use linear relation between m_flow and dp for any flow rate
Nominal condition
DensityrhoStdMedium.density_pTX(101325, 2...Inlet density for which valve coefficients are defined [kg/m3]
Assumptions
Dynamics
BooleandynamicBalancetrueSet to true to use a dynamic balance, which often leads to smaller systems of equations
Timetau10Time constant at nominal flow for dynamic energy and momentum balance [s]
MassFlowRatemDyn_flow_nominalm_flow_nominalNominal mass flow rate for dynamic momentum and energy balance [kg/s]
DynamicsenergyDynamicssystem.energyDynamicsFormulation of energy balance
DynamicsmassDynamicssystem.massDynamicsFormulation of mass balance
Initialization
AbsolutePressurep_startsystem.p_startStart value of pressure [Pa]
Booleanuse_T_starttrue= true, use T_start, otherwise h_start
TemperatureT_startif use_T_start then system.T...Start value of temperature [K]
SpecificEnthalpyh_startif use_T_start then Medium.s...Start value of specific enthalpy [J/kg]
MassFractionX_start[Medium.nX]Medium.X_defaultStart value of mass fractions m_i/m [kg/kg]
ExtraPropertyC_start[Medium.nC]fill(0, Medium.nC)Start value of trace substances

Connectors

TypeNameDescription
FluidPort_aport_1 
FluidPort_bport_2 
FluidPort_aport_3 
input RealInputyValve position (0: closed, 1: open)

Modelica definition

partial model PartialThreeWayValve "Partial three way valve"
    extends Buildings.Fluid.BaseClasses.PartialThreeWayResistance(
      final mDyn_flow_nominal = m_flow_nominal,
        redeclare FixedResistances.LosslessPipe res2(
            redeclare package Medium = Medium));
    extends Buildings.Fluid.Actuators.BaseClasses.ValveParameters(
      rhoStd=Medium.density_pTX(101325, 273.15+4, Medium.X_default),
      final dpVal_nominal=dp_nominal);

  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 m_flow_nominal(min=0) "Nominal mass flow rate";
  parameter Modelica.SIunits.Pressure dp_nominal(displayUnit="Pa") = 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 "Valve 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.Fluid.Actuators.BaseClasses.PartialTwoWayValve Buildings.Fluid.Actuators.BaseClasses.PartialTwoWayValve

Partial model for a two way valve

Buildings.Fluid.Actuators.BaseClasses.PartialTwoWayValve

Information


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.Fluid.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(dp_nominal)
Because the parameterization contains Kv_SI, the values for deltaM and dp_nominal 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.Fluid.Actuators.BaseClasses.PartialActuator (Partial model of an actuator), Buildings.Fluid.Actuators.BaseClasses.ValveParameters (Model with parameters for valves).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Reall0.0001Valve leakage, l=Cv(y=0)/Cvs
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
Pressuredp_nominal6000Pressure [Pa]
PressuredpVal_nominaldp_nominalNominal pressure drop [Pa]
Flow Coefficient
CvTypesCvDataBuildings.Fluid.Types.CvType...Selection of flow coefficient
RealKv0Kv (metric) flow coefficient [m3/h/(bar)^(1/2)]
RealCv0Cv (US) flow coefficient [USG/min/(psi)^(1/2)]
AreaAv0Av (metric) flow coefficient [m2]
RealKv_SIm_flow_nominal/sqrt(dpVal_no...Flow coefficient for fully open valve in SI units, Kv=m_flow/sqrt(dp) [kg/s/(Pa)^(1/2)]
Pressure-flow linearization
RealdeltaM0.02Fraction of nominal flow rate where linearization starts, if y=1
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
AbsolutePressuredp_start0.01*system.p_startGuess value of dp = port_a.p - port_b.p [Pa]
MassFlowRatem_flow_startsystem.m_flow_startGuess value of m_flow = port_a.m_flow [kg/s]
MassFlowRatem_flow_small1E-4*m_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
Diagnostics
Booleanshow_Ttrue= true, if temperatures at port_a and port_b are computed
Booleanshow_V_flowtrue= true, if volume flow rate at inflowing port is computed
Nominal condition
DensityrhoStdMedium.density_pTX(101325, 2...Inlet density for which valve coefficients are defined [kg/m3]

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

partial model PartialTwoWayValve "Partial model for a two way valve"
  extends Buildings.Fluid.Actuators.BaseClasses.PartialActuator(
       dp_nominal=6000);
  extends Buildings.Fluid.Actuators.BaseClasses.ValveParameters(
      rhoStd=Medium.density_pTX(101325, 273.15+4, Medium.X_default),
      final dpVal_nominal=dp_nominal);


  parameter Real l(min=0, max=1) = 0.0001 "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_turbulent = deltaM * m_flow_nominal;
 k = phi * Kv_SI;
end PartialTwoWayValve;

Buildings.Fluid.Actuators.BaseClasses.ValveParameters

Model with parameters for valves

Information


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.

Parameters

TypeNameDefaultDescription
Flow Coefficient
CvTypesCvDataBuildings.Fluid.Types.CvType...Selection of flow coefficient
RealKv0Kv (metric) flow coefficient [m3/h/(bar)^(1/2)]
RealCv0Cv (US) flow coefficient [USG/min/(psi)^(1/2)]
AreaAv0Av (metric) flow coefficient [m2]
RealKv_SIm_flow_nominal/sqrt(dpVal_no...Flow coefficient for fully open valve in SI units, Kv=m_flow/sqrt(dp) [kg/s/(Pa)^(1/2)]
Pressure-flow linearization
RealdeltaM0.02Fraction of nominal flow rate where linearization starts, if y=1
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
PressuredpVal_nominal Nominal pressure drop [Pa]
Advanced
Nominal condition
DensityrhoStd Inlet density for which valve coefficients are defined [kg/m3]

Modelica definition

partial model ValveParameters "Model with parameters for valves"


  parameter Buildings.Fluid.Types.CvTypes CvData=Buildings.Fluid.Types.CvTypes.OpPoint 
    "Selection of flow coefficient";
  parameter Real Kv(
    fixed= if CvData==Buildings.Fluid.Types.CvTypes.Kv then true else false) = 0 
    "Kv (metric) flow coefficient [m3/h/(bar)^(1/2)]";
  parameter Real Cv(
    fixed= if CvData==Buildings.Fluid.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.Fluid.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 m_flow_nominal(min=0) 
    "Nominal mass flow rate";
  parameter Modelica.SIunits.Pressure dpVal_nominal "Nominal pressure drop";

  parameter Real Kv_SI(
    min=0,
    fixed= if CvData==Buildings.Fluid.Types.CvTypes.OpPoint then true else false,
    start=m_flow_nominal/sqrt(dpVal_nominal)) = m_flow_nominal/sqrt(dpVal_nominal) 
    "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 
    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;

HTML-documentation generated by Dymola Sat Feb 6 17:25:35 2010.