Buildings.Fluid.Actuators.BaseClasses

Package with base classes for Buildings.Fluid.Actuators

Information

This package contains base classes that are used to construct the models in Buildings.Fluid.Actuators.

Extends from Modelica.Icons.BasesPackage (Icon for packages containing base classes).

Package Content

Name Description
Buildings.Fluid.Actuators.BaseClasses.ActuatorSignal ActuatorSignal Partial model that implements the filtered opening for valves and dampers
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
Buildings.Fluid.Actuators.BaseClasses.PartialTwoWayValveKv PartialTwoWayValveKv Partial model for a two way valve using a Kv characteristic
Buildings.Fluid.Actuators.BaseClasses.ValveParameters ValveParameters Model with parameters for valves
Buildings.Fluid.Actuators.BaseClasses.der_equalPercentage der_equalPercentage Derivative of valve opening characteristics for equal percentage valve
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.Examples Examples Collection of models that illustrate model use and test models

Buildings.Fluid.Actuators.BaseClasses.ActuatorSignal Buildings.Fluid.Actuators.BaseClasses.ActuatorSignal

Partial model that implements the filtered opening for valves and dampers

Buildings.Fluid.Actuators.BaseClasses.ActuatorSignal

Information

This model implements the filter that is used to approximate the travel time of the actuator. Models that extend this model use the signal y_actual to obtain the current position of the actuator.

See Buildings.Fluid.Actuators.UsersGuide for a description of the filter.

Parameters

TypeNameDefaultDescription
Dynamics
Filtered opening
Booleanuse_inputFiltertrue= true, if opening is filtered with a 2nd order CriticalDamping filter
TimeriseTime120Rise time of the filter (time to reach 99.6 % of an opening step) [s]
InitinitModelica.Blocks.Types.Init.I...Type of initialization (no init/steady state/initial state/initial output)
Realy_start1Initial position of actuator

Connectors

TypeNameDescription
input RealInputyActuator position (0: closed, 1: open)
output RealOutputy_actualActual actuator position

Modelica definition

model ActuatorSignal "Partial model that implements the filtered opening for valves and dampers" constant Integer order(min=1) = 2 "Order of filter"; parameter Boolean use_inputFilter=true "= true, if opening is filtered with a 2nd order CriticalDamping filter"; parameter Modelica.Units.SI.Time riseTime=120 "Rise time of the filter (time to reach 99.6 % of an opening step)"; parameter Modelica.Blocks.Types.Init init=Modelica.Blocks.Types.Init.InitialOutput "Type of initialization (no init/steady state/initial state/initial output)"; parameter Real y_start=1 "Initial position of actuator"; Modelica.Blocks.Interfaces.RealInput y(min=0, max=1) "Actuator position (0: closed, 1: open)"; Modelica.Blocks.Interfaces.RealOutput y_actual "Actual actuator position"; // Classes used to implement the filtered opening protected final parameter Modelica.Units.SI.Frequency fCut=5/(2*Modelica.Constants.pi* riseTime) "Cut-off frequency of filter"; parameter Boolean casePreInd = false "In case of PressureIndependent the model I/O is modified"; Modelica.Blocks.Interfaces.RealOutput y_internal(unit="1") "Output connector for internal use (= y_actual if not casePreInd)"; Modelica.Blocks.Interfaces.RealOutput y_filtered if use_inputFilter "Filtered valve position in the range 0..1"; Buildings.Fluid.BaseClasses.ActuatorFilter filter( final n=order, final f=fCut, final normalized=true, final initType=init, final y_start=y_start) if use_inputFilter "Second order filter to approximate actuator opening time, and to improve numerics"; equation connect(filter.y, y_filtered); if use_inputFilter then connect(y, filter.u); connect(filter.y, y_internal); else connect(y, y_internal); end if; if not casePreInd then connect(y_internal, y_actual); end if; end ActuatorSignal;

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 implements the functions that relate the opening signal and the flow coefficient. The model also defines parameters that are used by different air damper models.

The model is as in ASHRAE 825-RP except that a control signal of y=0 means the damper is closed, and y=1 means the damper is open. This is opposite of the implementation of ASHRAE 825-RP, but used here for consistency within this library.

For yL < y < yU, the damper characteristics is:

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

where kd is the loss coefficient (total pressure drop divided by dynamic pressure) and y is the fractional opening.

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

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

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

where A is the face area, which is computed using the nominal mass flow rate m_flow_nominal, the nominal velocity v_nominal and the density of the medium.

ASHRAE 825-RP lists the following parameter values as typical (note that the default values in the model correspond to opposed blades).

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

(The loss coefficient in fully closed position k0 is computed based on the leakage coefficient and the coefficient in fully open position.)

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.BaseClasses.PartialResistance (Partial model for a hydraulic resistance), Buildings.Fluid.Actuators.BaseClasses.ActuatorSignal (Partial model that implements the filtered opening for valves and dampers).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
MassFlowRatem_flow_turbulentif use_deltaM then deltaM*m_...Turbulent flow if |m_flow| >= m_flow_turbulent [kg/s]
Booleanuse_deltaMtrueSet to true to use deltaM for turbulent transition, else ReC is used
RealdeltaM0.3Fraction of nominal mass flow rate where transition to turbulent occurs
BooleanroundDuctfalseSet to true for round duct, false for square cross section
RealReC4000Reynolds number where transition to turbulence starts
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
PressureDifferencedp_nominaldpDamper_nominal + dpFixed_n...Pressure drop at nominal mass flow rate [Pa]
PressureDifferencedpDamper_nominal Pressure drop of fully open damper at nominal mass flow rate [Pa]
PressureDifferencedpFixed_nominal0Pressure drop of duct and resistances other than the damper in series, at nominal mass flow rate [Pa]
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
Booleanfrom_dpfalse= true, use m_flow = f(dp) else dp = f(m_flow)
Booleanlinearizedfalse= true, use linear relation between m_flow and dp for any flow rate
Booleanuse_constant_densitytrueSet to true to use constant density for flow friction
Dynamics
Filtered opening
Booleanuse_inputFiltertrue= true, if opening is filtered with a 2nd order CriticalDamping filter
TimeriseTime120Rise time of the filter (time to reach 99.6 % of an opening step) [s]
InitinitModelica.Blocks.Types.Init.I...Type of initialization (no init/steady state/initial state/initial output)
Realy_start1Initial position of actuator
Damper coefficients
Reala-1.51Coefficient a for damper characteristics [1]
Realb0.105*90Coefficient b for damper characteristics [1]
RealyL15/90Lower value for damper curve [1]
RealyU55/90Upper value for damper curve [1]
Realk10.45Loss coefficient for y=1 (pressure drop divided by dynamic pressure) [1]
Reall0.0001Damper leakage, ratio of flow coefficients k(y=0)/k(y=1)

Connectors

TypeNameDescription
FluidPort_aport_aFluid connector a (positive design flow direction is from port_a to port_b)
FluidPort_bport_bFluid connector b (positive design flow direction is from port_a to port_b)
input RealInputyActuator position (0: closed, 1: open)
output RealOutputy_actualActual actuator position

Modelica definition

partial model PartialDamperExponential "Partial model for air dampers with exponential opening characteristics" extends Buildings.Fluid.BaseClasses.PartialResistance( final dp_nominal=dpDamper_nominal+dpFixed_nominal, final m_flow_turbulent=if use_deltaM then deltaM * m_flow_nominal else eta_default*ReC*sqrt(A)*facRouDuc); extends Buildings.Fluid.Actuators.BaseClasses.ActuatorSignal; parameter Modelica.Units.SI.PressureDifference dpDamper_nominal(displayUnit= "Pa") "Pressure drop of fully open damper at nominal mass flow rate"; parameter Modelica.Units.SI.PressureDifference dpFixed_nominal(displayUnit= "Pa") = 0 "Pressure drop of duct and resistances other than the damper in series, at nominal mass flow rate"; 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"; final parameter Modelica.Units.SI.Velocity v_nominal=(2/rho_default/k1* dpDamper_nominal)^0.5 "Nominal face velocity"; final parameter Modelica.Units.SI.Area A=m_flow_nominal/rho_default/v_nominal "Face area"; parameter Boolean roundDuct = false "Set to true for round duct, false for square cross section"; parameter Real ReC=4000 "Reynolds number where transition to turbulence starts"; parameter Real a(unit="1")=-1.51 "Coefficient a for damper characteristics"; parameter Real b(unit="1")=0.105*90 "Coefficient b for damper characteristics"; parameter Real yL(unit="1") = 15/90 "Lower value for damper curve"; parameter Real yU(unit="1") = 55/90 "Upper value for damper curve"; final parameter Real k0(min=0, unit="1") = 2 * rho_default * (A / kDamMin)^2 "Loss coefficient for y=0 (pressure drop divided by dynamic pressure)"; parameter Real k1(min=0, unit="1") = 0.45 "Loss coefficient for y=1 (pressure drop divided by dynamic pressure)"; parameter Real l(min=1e-10, max=1) = 0.0001 "Damper leakage, ratio of flow coefficients k(y=0)/k(y=1)"; parameter Boolean use_constant_density=true "Set to true to use constant density for flow friction"; Medium.Density rho "Medium density"; final parameter Real kFixed = if dpFixed_nominal > Modelica.Constants.eps then m_flow_nominal / sqrt(dpFixed_nominal) else Modelica.Constants.inf "Flow coefficient of fixed resistance that may be in series with damper, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2)."; Real kDam "Flow coefficient of damper, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2)"; Real k "Flow coefficient of damper plus fixed resistance, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2)"; protected parameter Medium.Density rho_default=Medium.density(sta_default) "Density, used to compute fluid volume"; parameter Real facRouDuc= if roundDuct then sqrt(Modelica.Constants.pi)/2 else 1 "Shape factor used to compute the hydraulic diameter for round ducts"; parameter Real kL = Buildings.Fluid.Actuators.BaseClasses.exponentialDamper( y=yL, a=a, b=b, cL=cL, cU=cU, yL=yL, yU=yU)^2 "Loss coefficient at the lower limit of the exponential characteristics"; parameter Real kU = Buildings.Fluid.Actuators.BaseClasses.exponentialDamper( y=yU, a=a, b=b, cL=cL, cU=cU, yL=yL, yU=yU)^2 "Loss coefficient at the upper limit of the exponential characteristics"; parameter Real[3] cL={ (Modelica.Math.log(k0) - b - a)/yL^2, (-b*yL - 2*Modelica.Math.log(k0) + 2*b + 2*a)/yL, Modelica.Math.log(k0)} "Polynomial coefficients for curve fit for y < yl"; parameter Real[3] cU={ (Modelica.Math.log(k1) - a)/(yU^2 - 2*yU + 1), (-b*yU^2 - 2*Modelica.Math.log(k1)*yU - (-2*b - 2*a)*yU - b)/(yU^2 - 2*yU + 1), (Modelica.Math.log(k1)*yU^2 + b*yU^2 + (-2*b - 2*a)*yU + b + a)/(yU^2 - 2*yU + 1)} "Polynomial coefficients for curve fit for y > yu"; parameter Real kDamMax = (2 * rho_default / k1)^0.5 * A "Flow coefficient of damper fully open, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2)"; parameter Real kTotMax = if dpFixed_nominal > Modelica.Constants.eps then sqrt(1 / (1 / kFixed^2 + 1 / kDamMax^2)) else kDamMax "Flow coefficient of damper fully open plus fixed resistance, with unit=(kg.m)^(1/2)"; parameter Real kDamMin = l * kDamMax "Flow coefficient of damper fully closed, with unit=(kg.m)^(1/2)"; parameter Real kTotMin = if dpFixed_nominal > Modelica.Constants.eps then sqrt(1 / (1 / kFixed^2 + 1 / kDamMin^2)) else kDamMin "Flow coefficient of damper fully closed + fixed resistance, with unit=(kg.m)^(1/2)"; initial equation assert(dpDamper_nominal > Modelica.Constants.eps, "In " + getInstanceName() + ": dpDamper_nominal must be strictly greater than zero."); assert(dpFixed_nominal >= 0, "In " + getInstanceName() + ": dpFixed_nominal must be greater than zero."); assert(yL < yU, "In " + getInstanceName() + ": yL must be strictly lower than yU."); assert(m_flow_turbulent > 0, "In " + getInstanceName() + ": m_flow_turbulent must be strictly greater than zero."); assert(k1 >= 0.2, "In " + getInstanceName() + ": k1 must be greater than 0.2. k1=" + String(k1)); assert(k1 < kU, "In " + getInstanceName() + ": k1 must be strictly lower than exp(a + b * (1 - yU)). k1=" + String(k1) + ", exp(...) = " + String(kU)); assert(k0 <= 1e10, "In " + getInstanceName() + ": k0 must be lower than 1e10. k0=" + String(k0)); assert(k0 > kL, "In " + getInstanceName() + ": k0 must be strictly higher than exp(a + b * (1 - yL)). k0=" + String(k0) + ", exp(...) = " + String(kL)); equation rho = if use_constant_density then rho_default else Medium.density(Medium.setState_phX( port_a.p, inStream(port_a.h_outflow), inStream(port_a.Xi_outflow))); // flow coefficient, k = m_flow/sqrt(dp) kDam=sqrt(2*rho)*A/Buildings.Fluid.Actuators.BaseClasses.exponentialDamper( y=y_actual, a=a, b=b, cL=cL, cU=cU, yL=yL, yU=yU); k = if dpFixed_nominal > Modelica.Constants.eps then sqrt(1/(1/kFixed^2 + 1/kDam^2)) else kDam; 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 for flow from port_1 → port_2 is a parameter. The flow coefficient for the bypass flow from port_3 → port_2 is computed as

         Kv(port_3 → port_2)
  fraK = ----------------------
         Kv(port_1 → port_2)

where 0 < fraK ≤ 1 is a parameter with a default value of fraK=0.7.

Since this model uses two way valves to construct a three way valve, see Buildings.Fluid.Actuators.BaseClasses.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.ActuatorSignal (Partial model that implements the filtered opening for valves and dampers), Buildings.Fluid.Actuators.BaseClasses.ValveParameters (Model with parameters for valves).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
RealfraK0.7Fraction Kv(port_3&rarr;port_2)/Kv(port_1&rarr;port_2)
Reall[2]{0.0001,0.0001}Valve leakage, l=Kv(y=0)/Kv(y=1)
Flow Coefficient
CvTypesCvDataBuildings.Fluid.Types.CvType...Selection of flow coefficient
RealKv Kv (metric) flow coefficient [m3/h/(bar)^(1/2)]
RealCv Cv (US) flow coefficient [USG/min/(psi)^(1/2)]
AreaAv Av (metric) flow coefficient [m2]
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]
PressureDifferencedpValve_nominal Nominal pressure drop of fully open valve, used if CvData=Buildings.Fluid.Types.CvTypes.OpPoint [Pa]
PressureDifferencedpFixed_nominal[2]{0,0}Nominal pressure drop of pipes and other equipment in flow legs at port_1 and port_3 [Pa]
Dynamics
Conservation equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Type of energy balance: dynamic (3 initialization options) or steady state
MassFlowRatemDyn_flow_nominalm_flow_nominalNominal mass flow rate for dynamic momentum and energy balance [kg/s]
Nominal condition
Timetau10Time constant at nominal flow for dynamic energy and momentum balance [s]
Filtered opening
Booleanuse_inputFiltertrue= true, if opening is filtered with a 2nd order CriticalDamping filter
TimeriseTime120Rise time of the filter (time to reach 99.6 % of an opening step) [s]
InitinitModelica.Blocks.Types.Init.I...Type of initialization (no init/steady state/initial state/initial output)
Realy_start1Initial position of actuator
Initialization
AbsolutePressurep_startMedium.p_defaultStart value of pressure [Pa]
TemperatureT_startMedium.T_defaultStart value of temperature [K]
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
ExtraPropertyC_nominal[Medium.nC]fill(1E-2, Medium.nC)Nominal value of trace substances. (Set to typical order of magnitude.)
Advanced
Booleanfrom_dptrue= true, use m_flow = f(dp) else dp = f(m_flow)
PortFlowDirectionportFlowDirection_1Modelica.Fluid.Types.PortFlo...Flow direction for port_1
PortFlowDirectionportFlowDirection_2Modelica.Fluid.Types.PortFlo...Flow direction for port_2
PortFlowDirectionportFlowDirection_3Modelica.Fluid.Types.PortFlo...Flow direction for port_3
BooleanverifyFlowReversalfalse=true, to assert that the flow does not reverse when portFlowDirection_* does not equal Bidirectional
MassFlowRatem_flow_smallm_flow_nominal*1e-4Small mass flow rate for checking flow reversal [kg/s]
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]

Connectors

TypeNameDescription
FluidPort_aport_1First port, typically inlet
FluidPort_bport_2Second port, typically outlet
FluidPort_aport_3Third port, can be either inlet or outlet
input RealInputyActuator position (0: closed, 1: open)
output RealOutputy_actualActual actuator position

Modelica definition

partial model PartialThreeWayValve "Partial three way valve" extends Buildings.Fluid.BaseClasses.PartialThreeWayResistance( m_flow_small = m_flow_nominal*1e-4, final mDyn_flow_nominal = m_flow_nominal, redeclare replaceable Buildings.Fluid.Actuators.BaseClasses.PartialTwoWayValve res1 constrainedby Buildings.Fluid.Actuators.BaseClasses.PartialTwoWayValve( deltaM=deltaM, dp(start=dpValve_nominal/2), from_dp=from_dp, final linearized=linearized[1], final homotopyInitialization=homotopyInitialization, final CvData=Buildings.Fluid.Types.CvTypes.OpPoint, final m_flow_nominal=m_flow_nominal, final dpValve_nominal=dpValve_nominal, final dpFixed_nominal=dpFixed_nominal[1], final use_inputFilter=false, final riseTime=riseTime), redeclare FixedResistances.LosslessPipe res2( m_flow_nominal=m_flow_nominal), redeclare replaceable Buildings.Fluid.Actuators.BaseClasses.PartialTwoWayValve res3 constrainedby Buildings.Fluid.Actuators.BaseClasses.PartialTwoWayValve( deltaM=deltaM, dp(start=dpValve_nominal/2), from_dp=from_dp, final linearized=linearized[2], final homotopyInitialization=homotopyInitialization, final CvData=Buildings.Fluid.Types.CvTypes.OpPoint, final m_flow_nominal=m_flow_nominal, final dpValve_nominal=dpValve_nominal/fraK^2, final dpFixed_nominal=dpFixed_nominal[2], final use_inputFilter=false, final riseTime=riseTime)); extends Buildings.Fluid.Actuators.BaseClasses.ActuatorSignal; extends Buildings.Fluid.Actuators.BaseClasses.ValveParameters( rhoStd=Medium.density_pTX(101325, 273.15+4, Medium.X_default)); constant Boolean homotopyInitialization = true "= true, use homotopy method"; parameter Modelica.Units.SI.PressureDifference dpFixed_nominal[2]( each displayUnit="Pa", each min=0) = {0,0} "Nominal pressure drop of pipes and other equipment in flow legs at port_1 and port_3"; parameter Real fraK(min=0, max=1) = 0.7 "Fraction Kv(port_3&rarr;port_2)/Kv(port_1&rarr;port_2)"; parameter Real[2] l(each min=0, each max=1) = {0.0001, 0.0001} "Valve leakage, l=Kv(y=0)/Kv(y=1)"; parameter Real deltaM = 0.02 "Fraction of nominal flow rate where linearization starts, if y=1"; parameter Boolean[2] linearized = {false, false} "= true, use linear relation between m_flow and dp for any flow rate"; protected Modelica.Blocks.Math.Feedback inv "Inversion of control signal"; Modelica.Blocks.Sources.Constant uni(final k=1) "Outputs one for bypass valve"; initial equation assert(homotopyInitialization, "In " + getInstanceName() + ": The constant homotopyInitialization has been modified from its default value. This constant will be removed in future releases.", level = AssertionLevel.warning); 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, quick opening or pressure-independent.

To prevent the derivative d/dP (m_flow) to be infinite near the origin, this model linearizes the pressure drop versus flow relation ship. The region in which it is linearized is parameterized by

  m_turbulent_flow = deltaM * m_flow_nominal

Because the parameterization contains Kv_SI, the values for deltaM and dp_nominal need not be changed if the valve size changes.

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.

Options

This model allows different parameterization of the flow resistance. The different parameterizations are described in Buildings.Fluid.Actuators.BaseClasses.ValveParameters.

Implementation

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.BaseClasses.PartialResistance (Partial model for a hydraulic resistance), Buildings.Fluid.Actuators.BaseClasses.ValveParameters (Model with parameters for valves), Buildings.Fluid.Actuators.BaseClasses.ActuatorSignal (Partial model that implements the filtered opening for valves and dampers).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
MassFlowRatem_flow_turbulentdeltaM*abs(m_flow_nominal)Turbulent flow if |m_flow| >= m_flow_turbulent [kg/s]
Reall0.0001Valve leakage, l=Kv(y=0)/Kv(y=1)
RealkFixedif dpFixed_nominal > Modelic...Flow coefficient of fixed resistance that may be in series with valve, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2).
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
PressureDifferencedp_nominaldpValve_nominal + dpFixed_no...Pressure drop at nominal mass flow rate [Pa]
PressureDifferencedpValve_nominal Nominal pressure drop of fully open valve, used if CvData=Buildings.Fluid.Types.CvTypes.OpPoint [Pa]
PressureDifferencedpFixed_nominal0Pressure drop of pipe and other resistances that are in series [Pa]
Flow Coefficient
CvTypesCvDataBuildings.Fluid.Types.CvType...Selection of flow coefficient
RealKv Kv (metric) flow coefficient [m3/h/(bar)^(1/2)]
RealCv Cv (US) flow coefficient [USG/min/(psi)^(1/2)]
AreaAv Av (metric) flow coefficient [m2]
Pressure-flow linearization
RealdeltaM0.02Fraction of nominal flow rate where linearization starts, if y=1
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
Booleanfrom_dpfalse= true, use m_flow = f(dp) else dp = f(m_flow)
Booleanlinearizedfalse= 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]
Dynamics
Filtered opening
Booleanuse_inputFiltertrue= true, if opening is filtered with a 2nd order CriticalDamping filter
TimeriseTime120Rise time of the filter (time to reach 99.6 % of an opening step) [s]
InitinitModelica.Blocks.Types.Init.I...Type of initialization (no init/steady state/initial state/initial output)
Realy_start1Initial position of actuator

Connectors

TypeNameDescription
FluidPort_aport_aFluid connector a (positive design flow direction is from port_a to port_b)
FluidPort_bport_bFluid connector b (positive design flow direction is from port_a to port_b)
input RealInputyActuator position (0: closed, 1: open)
output RealOutputy_actualActual actuator position

Modelica definition

partial model PartialTwoWayValve "Partial model for a two way valve" extends Buildings.Fluid.BaseClasses.PartialResistance( final dp_nominal=dpValve_nominal + dpFixed_nominal, dp(nominal=6000), final m_flow_turbulent = deltaM * abs(m_flow_nominal)); extends Buildings.Fluid.Actuators.BaseClasses.ValveParameters( rhoStd=Medium.density_pTX(101325, 273.15+4, Medium.X_default)); extends Buildings.Fluid.Actuators.BaseClasses.ActuatorSignal; parameter Modelica.Units.SI.PressureDifference dpFixed_nominal( displayUnit="Pa", min=0) = 0 "Pressure drop of pipe and other resistances that are in series"; parameter Real l(min=1e-10, max=1) = 0.0001 "Valve leakage, l=Kv(y=0)/Kv(y=1)"; input Real phi "Ratio actual to nominal mass flow rate of valve, phi=Kv(y)/Kv(y=1)"; parameter Real kFixed(unit="", min=0) = if dpFixed_nominal > Modelica.Constants.eps then m_flow_nominal / sqrt(dpFixed_nominal) else 0 "Flow coefficient of fixed resistance that may be in series with valve, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2)."; Real kVal(unit="", min=Modelica.Constants.small) "Flow coefficient of valve, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2)."; Real k(unit="", min=Modelica.Constants.small) "Flow coefficient of valve and pipe in series, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2)."; initial equation assert(dpFixed_nominal > -Modelica.Constants.eps, "In " + getInstanceName() + ": Model requires dpFixed_nominal >= 0 but received dpFixed_nominal = " + String(dpFixed_nominal) + " Pa."); end PartialTwoWayValve;

Buildings.Fluid.Actuators.BaseClasses.PartialTwoWayValveKv Buildings.Fluid.Actuators.BaseClasses.PartialTwoWayValveKv

Partial model for a two way valve using a Kv characteristic

Buildings.Fluid.Actuators.BaseClasses.PartialTwoWayValveKv

Information

Partial model for valves with different opening characteristics, such as linear, equal percentage or quick opening. This partial extends from Buildings.Fluid.Actuators.BaseClasses.PartialTwoWayValve and also contains the governing equations for these three two way valve models.

Implementation

Models that extend this model need to provide a binding equation for the flow function phi. An example of such a code can be found in Buildings.Fluid.Actuators.Valves.TwoWayLinear.

Extends from Buildings.Fluid.Actuators.BaseClasses.PartialTwoWayValve (Partial model for a two way valve).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Reall0.0001Valve leakage, l=Kv(y=0)/Kv(y=1)
RealkFixedif dpFixed_nominal > Modelic...Flow coefficient of fixed resistance that may be in series with valve, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2).
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
PressureDifferencedpValve_nominal Nominal pressure drop of fully open valve, used if CvData=Buildings.Fluid.Types.CvTypes.OpPoint [Pa]
PressureDifferencedpFixed_nominal0Pressure drop of pipe and other resistances that are in series [Pa]
Flow Coefficient
CvTypesCvDataBuildings.Fluid.Types.CvType...Selection of flow coefficient
RealKv Kv (metric) flow coefficient [m3/h/(bar)^(1/2)]
RealCv Cv (US) flow coefficient [USG/min/(psi)^(1/2)]
AreaAv Av (metric) flow coefficient [m2]
Pressure-flow linearization
RealdeltaM0.02Fraction of nominal flow rate where linearization starts, if y=1
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
Booleanfrom_dpfalse= true, use m_flow = f(dp) else dp = f(m_flow)
Booleanlinearizedfalse= 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]
Dynamics
Filtered opening
Booleanuse_inputFiltertrue= true, if opening is filtered with a 2nd order CriticalDamping filter
TimeriseTime120Rise time of the filter (time to reach 99.6 % of an opening step) [s]
InitinitModelica.Blocks.Types.Init.I...Type of initialization (no init/steady state/initial state/initial output)
Realy_start1Initial position of actuator

Connectors

TypeNameDescription
FluidPort_aport_aFluid connector a (positive design flow direction is from port_a to port_b)
FluidPort_bport_bFluid connector b (positive design flow direction is from port_a to port_b)
input RealInputyActuator position (0: closed, 1: open)
output RealOutputy_actualActual actuator position

Modelica definition

partial model PartialTwoWayValveKv "Partial model for a two way valve using a Kv characteristic" extends Buildings.Fluid.Actuators.BaseClasses.PartialTwoWayValve; equation kVal = phi*Kv_SI; if (dpFixed_nominal > Modelica.Constants.eps) then k = sqrt(1/(1/kFixed^2 + 1/kVal^2)); else k = kVal; end if; if linearized then // This implementation yields m_flow_nominal = phi*kv_SI * sqrt(dp_nominal) // if m_flow = m_flow_nominal and dp = dp_nominal m_flow*m_flow_nominal_pos = k^2 * dp; else if homotopyInitialization then if from_dp then m_flow=homotopy(actual=Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_dp(dp=dp, k=k, m_flow_turbulent=m_flow_turbulent), simplified=m_flow_nominal_pos*dp/dp_nominal_pos); else dp=homotopy(actual=Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_m_flow(m_flow=m_flow, k=k, m_flow_turbulent=m_flow_turbulent), simplified=dp_nominal_pos*m_flow/m_flow_nominal_pos); end if; else // do not use homotopy if from_dp then m_flow=Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_dp(dp=dp, k=k, m_flow_turbulent=m_flow_turbulent); else dp=Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_m_flow(m_flow=m_flow, k=k, m_flow_turbulent=m_flow_turbulent); end if; end if; // homotopyInitialization end if; // linearized end PartialTwoWayValveKv;

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 protected 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. The value of Kv_SI is computed based on the parameters Av, Kv, Cv, or, if CvData = Buildings.Fluid.Types.CvTypes.OpPoint, based on m_flow_nominal and dpValve_nominal. Conversely, if CvData <> Buildings.Fluid.Types.CvTypes.OpPoint, then dpValve_nominal is computed based on Av, Kv, or Cv, and the nominal mass flow rate m_flow_nominal. Therefore, if CvData <> Buildings.Fluid.Types.CvTypes.OpPoint, then specifying a value for dpValve_nominal is a syntax error.

Parameters

TypeNameDefaultDescription
Flow Coefficient
CvTypesCvDataBuildings.Fluid.Types.CvType...Selection of flow coefficient
RealKv Kv (metric) flow coefficient [m3/h/(bar)^(1/2)]
RealCv Cv (US) flow coefficient [USG/min/(psi)^(1/2)]
AreaAv Av (metric) flow coefficient [m2]
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]
PressureDifferencedpValve_nominal Nominal pressure drop of fully open valve, used if CvData=Buildings.Fluid.Types.CvTypes.OpPoint [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) "Kv (metric) flow coefficient [m3/h/(bar)^(1/2)]"; parameter Real Cv( fixed= if CvData==Buildings.Fluid.Types.CvTypes.Cv then true else false) "Cv (US) flow coefficient [USG/min/(psi)^(1/2)]"; parameter Modelica.Units.SI.Area Av(fixed=if CvData == Buildings.Fluid.Types.CvTypes.Av then true else false) "Av (metric) flow coefficient"; parameter Real deltaM = 0.02 "Fraction of nominal flow rate where linearization starts, if y=1"; parameter Modelica.Units.SI.MassFlowRate m_flow_nominal "Nominal mass flow rate"; parameter Modelica.Units.SI.PressureDifference dpValve_nominal( displayUnit="Pa", min=0, fixed=if CvData == Buildings.Fluid.Types.CvTypes.OpPoint then true else false) "Nominal pressure drop of fully open valve, used if CvData=Buildings.Fluid.Types.CvTypes.OpPoint"; parameter Modelica.Units.SI.Density rhoStd "Inlet density for which valve coefficients are defined"; protected parameter Real Kv_SI( min=0, fixed= false) "Flow coefficient for fully open valve in SI units, Kv=m_flow/sqrt(dp) [kg/s/(Pa)^(1/2)]"; initial equation if CvData == Buildings.Fluid.Types.CvTypes.OpPoint then Kv_SI = m_flow_nominal/sqrt(dpValve_nominal); Kv = Kv_SI/(rhoStd/3600/sqrt(1E5)); Cv = Kv_SI/(rhoStd*0.0631/1000/sqrt(6895)); Av = Kv_SI/sqrt(rhoStd); elseif CvData == Buildings.Fluid.Types.CvTypes.Kv then Kv_SI = Kv*rhoStd/3600/sqrt(1E5) "Unit conversion m3/(h*sqrt(bar)) to kg/(s*sqrt(Pa))"; Cv = Kv_SI/(rhoStd*0.0631/1000/sqrt(6895)); Av = Kv_SI/sqrt(rhoStd); dpValve_nominal = (m_flow_nominal/Kv_SI)^2; elseif CvData == Buildings.Fluid.Types.CvTypes.Cv then Kv_SI = Cv*rhoStd*0.0631/1000/sqrt(6895) "Unit conversion USG/(min*sqrt(psi)) to kg/(s*sqrt(Pa))"; Kv = Kv_SI/(rhoStd/3600/sqrt(1E5)); Av = Kv_SI/sqrt(rhoStd); dpValve_nominal = (m_flow_nominal/Kv_SI)^2; else assert(CvData == Buildings.Fluid.Types.CvTypes.Av, "Invalid value for CvData. Obtained CvData = " + String(CvData) + "."); Kv_SI = Av*sqrt(rhoStd); Kv = Kv_SI/(rhoStd/3600/sqrt(1E5)); Cv = Kv_SI/(rhoStd*0.0631/1000/sqrt(6895)); dpValve_nominal = (m_flow_nominal/Kv_SI)^2; end if; end ValveParameters;

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

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

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

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" extends Modelica.Icons.Function; 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=Kv(y=0)/Kv(y=1)"; 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 * der_y; else if (y > (3/2 * delta)) then der_phi := R^(y-1)*Modelica.Math.log(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.equalPercentage 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 Buildings.Fluid.Actuators.Valves.TwoWayEqualPercentage.

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.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realy Valve opening signal, y=1 is fully open
RealR Rangeability, R=50...100 typically
Reall Valve leakage, l=Kv(y=0)/Kv(y=1)
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); extends Modelica.Icons.Function; 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=Kv(y=0)/Kv(y=1)"; 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 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 Buildings.Fluid.Actuators.Dampers.Exponential.

For yL < y < yU, the damper characteristics is

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

Outside this range, the damper characteristic is defined by a quadratic polynomial.

Note that this implementation returns sqrt(kd(y)) instead of kd(y). This is done for numerical reason since otherwise kd(y) may be an iteration variable, which may cause a lot of warnings and slower convergence if the solver attempts kd(y) < 0 during the iterative solution procedure.

Extends from Modelica.Icons.Function (Icon for functions).

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
RealkThetaSqRtSquare root of loss coefficient, sqrt(pressure drop/dynamic pressure)

Modelica definition

function exponentialDamper "Damper opening characteristics for an exponential damper" extends Modelica.Icons.Function; input Real y(min=0, max=1, 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 kThetaSqRt(min=0) "Square root of loss coefficient, sqrt(pressure drop/dynamic pressure)"; protected Real yC(min=0, max=1, unit="") "y constrained to 0 <= y <= 1 to avoid numerical problems"; algorithm if y < yL then yC :=max(0, y); kThetaSqRt := sqrt(Modelica.Math.exp(cL[3] + yC * (cL[2] + yC * cL[1]))); else if (y > yU) then yC := min(1, y); kThetaSqRt := sqrt(Modelica.Math.exp(cU[3] + yC * (cU[2] + yC * cU[1]))); else kThetaSqRt := sqrt(Modelica.Math.exp(a+b*(1-y))) "y=0 is closed"; end if; end if; end exponentialDamper;