Name | Description |
---|---|
der_equalPercentage | Derivative of valve opening characteristics for equal percentage valve |
equalPercentage | Valve opening characteristics for equal percentage valve |
Examples | Collection of models that illustrate model use and test models |
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 |
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 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;
Type | Name | Default | Description |
---|---|---|---|
replaceable package Medium | PartialMedium | Medium in the component | |
MassFlowRate | m_small_flow | Mass flow rate where transition to laminar occurs [kg/s] | |
Initialization | |||
MassFlowRate | m_flow | Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction) [kg/s] | |
Pressure | dp | Pressure difference between port_a and port_b [Pa] | |
Advanced | |||
Temp | flowDirection | Modelica_Fluid.Types.FlowDir... | Unidirectional (port_a -> port_b) or bidirectional flow component |
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 |
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 "Damper position (0: closed, 1: open)"; end PartialActuator;
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.
Type | Name | Default | Description |
---|---|---|---|
replaceable package Medium | PartialMedium | Medium in the component | |
MassFlowRate | m_small_flow | eta0*ReC*sqrt(A)*facRouDuc | Mass flow rate where transition to laminar occurs [kg/s] |
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 |
Initialization | |||
MassFlowRate | m_flow | Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction) [kg/s] | |
Pressure | dp | Pressure difference between port_a and port_b [Pa] | |
Advanced | |||
Temp | flowDirection | Modelica_Fluid.Types.FlowDir... | Unidirectional (port_a -> port_b) or bidirectional flow component |
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 |
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(m_small_flow=eta0*ReC*sqrt(A)*facRouDuc); 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"; Real kDam(unit="(kg*m)^(1/2)", start=1) "Flow coefficient for damper, k=m_flow/sqrt(dp)"; protected Real kTheta(min=0) "Flow coefficient, kTheta = pressure drop divided by dynamic pressure"; protected parameter Real facRouDuc= if roundDuct then sqrt(Modelica.Constants.pi)/2 else 1; equation assert(y >= (15/90) and y <= (55/90), "Damper characteristics not implemented for angles outside 15...55 degree.\n" + "Received y = " + realString(y) + ". Corresponds to " + realString(y*90) + " degrees."); kTheta = exp(a+b*(1-y)) "y=0 is closed, but theta=1 is closed in ASHRAE-825"; kDam = sqrt(kTheta/2/medium_a.d) / A "flow coefficient for resistance base model"; end PartialDamperExponential;
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 k_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
k_SI(port_1 -> port_2) fraK = ---------------------- k_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.
Type | Name | Default | Description |
---|---|---|---|
replaceable package Medium | PartialMedium | 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 | k_SI | Flow coefficient for fully open valve in SI units, k=m_flow/sqrt(dp) [(kg*m)^(1/2)] | |
Real | fraK | 0.7 | Fraction k_SI(port_1->port_2)/k_SI(port_3->port_2) |
Real | l[2] | {0,0} | Valve leakage, l=Cv(y=0)/Cvs |
Pressure-flow linearization | |||
Real | deltaM | 0.02 | Fraction of nominal flow rate where linearization starts, if y=1 |
Pressure | dp0 | 6000 | Nominal pressure drop [Pa] |
Advanced | |||
Boolean | from_dp | true | = true, use m_flow = f(dp) else dp = f(m_flow) |
Temp | flowDirection | Modelica_Fluid.Types.FlowDir... | Unidirectional (port_1 -> port_2) or bidirectional flow component |
Boolean | linearized[2] | {false,false} | = true, use linear relation between m_flow and dp for any flow rate |
Type | Name | Description |
---|---|---|
FluidPort_b | port_1 | |
FluidPort_b | port_2 | |
input RealInput | y | Damper position (0: closed, 1: open) |
partial model PartialThreeWayValve "Partial three way valve" extends Buildings.Fluids.BaseClasses.PartialThreeWayResistance( port_3( m_flow(start=0,min=if allowFlowReversal then -Modelica.Constants.inf else 0)), redeclare FixedResistances.LosslessPipe res2( redeclare package Medium = Medium, flowDirection=Modelica_Fluid.Types.FlowDirection.Bidirectional)); parameter Real k_SI(min=0, unit="(kg*m)^(1/2)") "Flow coefficient for fully open valve in SI units, k=m_flow/sqrt(dp)"; parameter Real fraK(min=0, max=1) = 0.7 "Fraction k_SI(port_1->port_2)/k_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 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;
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.
The parameter k_SI 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 assume infinity 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 * k_SI * sqrt(dp0)Because the parameterization contains k_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.
Type | Name | Default | Description |
---|---|---|---|
replaceable package Medium | PartialMedium | Medium in the component | |
Real | k_SI | Flow coefficient for fully open valve in SI units, k=m_flow/sqrt(dp) [(kg*m)^(1/2)] | |
Real | l | 0 | Valve leakage, l=Cv(y=0)/Cvs |
Initialization | |||
MassFlowRate | m_flow | Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction) [kg/s] | |
Pressure | dp | Pressure difference between port_a and port_b [Pa] | |
Pressure-flow linearization | |||
Real | deltaM | 0.02 | Fraction of nominal flow rate where linearization starts, if y=1 |
Pressure | dp0 | 6000 | Nominal pressure drop [Pa] |
Advanced | |||
Temp | flowDirection | Modelica_Fluid.Types.FlowDir... | Unidirectional (port_a -> port_b) or bidirectional flow component |
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 |
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 PartialActuator(final m_small_flow = deltaM * k_SI * sqrt(dp0)); parameter Real k_SI(min=0, unit="(kg*m)^(1/2)") "Flow coefficient for fully open valve in SI units, k=m_flow/sqrt(dp)"; parameter Real l(min=0, max=1) = 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 Modelica.SIunits.Pressure dp0 = 6000 "Nominal pressure drop"; Real phi "Ratio actual to nominal mass flow rate, phi=Cv(y)/Cv(y=1)"; equation k = phi * k_SI; end PartialTwoWayValve;