LBL logo

Buildings.Fluid.Movers.BaseClasses.Characteristics

Functions for fan or pump characteristics

Information

This package implements performance curves for fans and pumps, and records for parameter that can be used with these performance curves.

The following performance curves are implemented:
Independent variable Dependent variable Record for performance data Function
Volume flow rate Pressure flowParameters pressure
Relative volumetric flow rate Efficiency efficiencyParameters efficiency
Volume flow rate Power powerParameters power

Package Content

NameDescription
Buildings.Fluid.Movers.BaseClasses.Characteristics.flowParameters flowParameters Record for flow parameters
Buildings.Fluid.Movers.BaseClasses.Characteristics.flowParametersInternal flowParametersInternal Record for flow parameters with prescribed size
Buildings.Fluid.Movers.BaseClasses.Characteristics.efficiencyParameters efficiencyParameters Record for efficiency parameters
Buildings.Fluid.Movers.BaseClasses.Characteristics.powerParameters powerParameters Record for electrical power parameters
Buildings.Fluid.Movers.BaseClasses.Characteristics.pressure pressure Flow vs. head characteristics for fan or pump pressure raise
Buildings.Fluid.Movers.BaseClasses.Characteristics.flowApproximationAtOrigin flowApproximationAtOrigin Approximation for fan or pump pressure raise at origin
Buildings.Fluid.Movers.BaseClasses.Characteristics.power power Flow vs. electrical power characteristics for fan or pump
Buildings.Fluid.Movers.BaseClasses.Characteristics.efficiency efficiency Flow vs. efficiency characteristics for fan or pump

Buildings.Fluid.Movers.BaseClasses.Characteristics.flowParameters Buildings.Fluid.Movers.BaseClasses.Characteristics.flowParameters

Record for flow parameters

Information

Data record for performance data that describe volume flow rate versus pressure rise. The volume flow rate V_flow must be increasing, i.e., V_flow[i] < V_flow[i+1]. Both vectors, V_flow and dp must have the same size.

Extends from Modelica.Icons.Record (Icon for records).

Parameters

TypeNameDefaultDescription
VolumeFlowRateV_flow[:] Volume flow rate at user-selected operating points [m3/s]
Pressuredp[size(V_flow, 1)] Fan or pump total pressure at these flow rates [Pa]

Modelica definition

record flowParameters "Record for flow parameters"
  extends Modelica.Icons.Record;
  parameter Modelica.SIunits.VolumeFlowRate V_flow[:](each min=0) 
    "Volume flow rate at user-selected operating points";
  parameter Modelica.SIunits.Pressure dp[size(V_flow,1)](
     each min=0, each displayUnit="Pa") 
    "Fan or pump total pressure at these flow rates";
end flowParameters;

Buildings.Fluid.Movers.BaseClasses.Characteristics.flowParametersInternal Buildings.Fluid.Movers.BaseClasses.Characteristics.flowParametersInternal

Record for flow parameters with prescribed size

Information

Data record for performance data that describe volume flow rate versus pressure rise. The volume flow rate V_flow must be increasing, i.e., V_flow[i] < V_flow[i+1]. Both vectors, V_flow and dp must have the same size.

This record is identical to Buildings.Fluid.Movers.BaseClasses.Characteristic.flowParameters, except that it takes the size of the array as a parameter. This is required in Dymola 2014. Otherwise, the array size would need to be computed in Buildings.Fluid.Movers.BaseClasses.FlowMachineInterface in the initial algorithm section, which is not supported.

Extends from Modelica.Icons.Record (Icon for records).

Parameters

TypeNameDefaultDescription
Integern Number of elements in each array
VolumeFlowRateV_flow[n] Volume flow rate at user-selected operating points [m3/s]
Pressuredp[n] Fan or pump total pressure at these flow rates [Pa]

Modelica definition

record flowParametersInternal 
  "Record for flow parameters with prescribed size"
  extends Modelica.Icons.Record;
  parameter Integer n "Number of elements in each array";
  parameter Modelica.SIunits.VolumeFlowRate V_flow[n](each min=0) 
    "Volume flow rate at user-selected operating points";
  parameter Modelica.SIunits.Pressure dp[n](
     each min=0, each displayUnit="Pa") 
    "Fan or pump total pressure at these flow rates";
end flowParametersInternal;

Buildings.Fluid.Movers.BaseClasses.Characteristics.efficiencyParameters Buildings.Fluid.Movers.BaseClasses.Characteristics.efficiencyParameters

Record for efficiency parameters

Information

Data record for performance data that describe volume flow rate versus efficiency. The volume flow rate r_V must be increasing, i.e., r_V[i] < r_V[i+1]. Both vectors, r_V and eta must have the same size.

Extends from Modelica.Icons.Record (Icon for records).

Parameters

TypeNameDefaultDescription
Realr_V[:] Volumetric flow rate divided by nominal flow rate at user-selected operating points
Realeta[size(r_V, 1)] Fan or pump efficiency at these flow rates

Modelica definition

record efficiencyParameters "Record for efficiency parameters"
  extends Modelica.Icons.Record;
  parameter Real  r_V[:](each min=0, each max=1, each displayUnit="1") 
    "Volumetric flow rate divided by nominal flow rate at user-selected operating points";
  parameter Real eta[size(r_V,1)](
     each min=0, each max=1, each displayUnit="1") 
    "Fan or pump efficiency at these flow rates";
end efficiencyParameters;

Buildings.Fluid.Movers.BaseClasses.Characteristics.powerParameters Buildings.Fluid.Movers.BaseClasses.Characteristics.powerParameters

Record for electrical power parameters

Information

Data record for performance data that describe volume flow rate versus electrical power. The volume flow rate V_flow must be increasing, i.e., V_flow[i] < V_flow[i+1]. Both vectors, V_flow and P must have the same size.

Extends from Modelica.Icons.Record (Icon for records).

Parameters

TypeNameDefaultDescription
VolumeFlowRateV_flow[:]{0}Volume flow rate at user-selected operating points [m3/s]
PowerP[size(V_flow, 1)]{0}Fan or pump electrical power at these flow rates [W]

Modelica definition

record powerParameters "Record for electrical power parameters"
  extends Modelica.Icons.Record;
  parameter Modelica.SIunits.VolumeFlowRate V_flow[:](each min=0)= {0} 
    "Volume flow rate at user-selected operating points";
  parameter Modelica.SIunits.Power P[size(V_flow,1)](
     each min=0) = {0} "Fan or pump electrical power at these flow rates";
end powerParameters;

Buildings.Fluid.Movers.BaseClasses.Characteristics.pressure Buildings.Fluid.Movers.BaseClasses.Characteristics.pressure

Flow vs. head characteristics for fan or pump pressure raise

Information

This function computes the fan static pressure raise as a function of volume flow rate and revolution in the form

Δp = rN2   s(V/rN, d) - Δpr ,

where Δp is the pressure rise, rN is the normalized fan speed, V is the volume flow rate and d are performance data for fan or pump power consumption at rN=1. The term

Δpr = V   Δpmax ⁄ Vmax   δ

models the flow resistance of the fan, approximated using a linear equation. This is done for numerical reasons to avoid a singularity at rN=0. Since δ is small, the contribution of this term is small. The fan and pump models in Buildings.Fluid.Movers modify the user-supplied performance data to add the term Δpr prior to computing the performance curve. Thus, at full speed, the fan or pump can operate exactly at the user-supplied performance data.

Implementation

The function s(·, ·) is a cubic hermite spline. If the data d define a monotone decreasing sequence, then s(·, d) is a monotone decreasing function.

For rN < δ, the polynomial is replaced with an other model to avoid a singularity at the origin. The composite model is once continuously differentiable in all input variables.

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

Inputs

TypeNameDefaultDescription
flowParametersInternaldata Pressure performance data
VolumeFlowRateV_flow Volumetric flow rate [m3/s]
Realr_N Relative revolution, r_N=N/N_nominal [1]
VolumeFlowRateVDelta_flow Small volume flow rate [m3/s]
PressuredpDelta Small pressure [Pa]
VolumeFlowRateV_flow_max Maximum volume flow rate at r_N=1 and dp=0 [m3/s]
PressuredpMax Maximum pressure at r_N=1 and V_flow=0 [Pa]
Reald[:] Derivatives at support points for spline interpolation
Realdelta Small value used to transition to other fan curve
RealcBar[2] Coefficients for linear approximation of pressure vs. flow rate
RealkRes Linear coefficient for fan-internal pressure drop [kg/(s.m4)]

Outputs

TypeNameDescription
PressuredpPressure raise [Pa]

Modelica definition

function pressure 
  "Flow vs. head characteristics for fan or pump pressure raise"
  extends Modelica.Icons.Function;
  input Buildings.Fluid.Movers.BaseClasses.Characteristics.flowParametersInternal
                                                                                  data 
    "Pressure performance data";
  input Modelica.SIunits.VolumeFlowRate V_flow "Volumetric flow rate";
  input Real r_N(unit="1") "Relative revolution, r_N=N/N_nominal";
  input Modelica.SIunits.VolumeFlowRate VDelta_flow "Small volume flow rate";
  input Modelica.SIunits.Pressure dpDelta "Small pressure";

  input Modelica.SIunits.VolumeFlowRate V_flow_max 
    "Maximum volume flow rate at r_N=1 and dp=0";
  input Modelica.SIunits.Pressure dpMax(min=0) 
    "Maximum pressure at r_N=1 and V_flow=0";

  input Real d[:] "Derivatives at support points for spline interpolation";
  input Real delta "Small value used to transition to other fan curve";
  input Real cBar[2] 
    "Coefficients for linear approximation of pressure vs. flow rate";
  input Real kRes(unit="kg/(s.m4)") 
    "Linear coefficient for fan-internal pressure drop";
  output Modelica.SIunits.Pressure dp "Pressure raise";

protected 
   Integer dimD(min=2)=size(data.V_flow, 1) "Dimension of data vector";

  function performanceCurve "Performance curve away from the origin"
    input Modelica.SIunits.VolumeFlowRate V_flow "Volumetric flow rate";
    input Real r_N(unit="1") "Relative revolution, r_N=N/N_nominal";
    input Real d[dimD] "Coefficients for polynomial of pressure vs. flow rate";
    input Buildings.Fluid.Movers.BaseClasses.Characteristics.flowParametersInternal
                                                                                    data 
      "Pressure performance data";
    input Integer dimD "Dimension of data vector";

    output Modelica.SIunits.Pressure dp "Pressure raise";

  protected 
    Modelica.SIunits.VolumeFlowRate rat "Ratio of V_flow/r_N";
    Integer i "Integer to select data interval";
  algorithm 
    rat := V_flow/r_N;
    i :=1;
    // Since the coefficients for the spline were evaluated for
    // rat_nominal = V_flow_nominal/r_N_nominal = V_flow_nominal/1, we use
    // V_flow_nominal below
    for j in 1:dimD-1 loop
       if rat > data.V_flow[j] then
         i := j;
       end if;
    end for;
    // Extrapolate or interpolate the data
    dp:=r_N^2*Buildings.Utilities.Math.Functions.cubicHermiteLinearExtrapolation(
                x=rat,
                x1=data.V_flow[i],
                x2=data.V_flow[i + 1],
                y1=data.dp[i],
                y2=data.dp[i + 1],
                y1d=d[i],
                y2d=d[i+1]);
  end performanceCurve;

algorithm 
  if r_N >= delta then
     dp := performanceCurve(V_flow=V_flow, r_N=r_N, d=d,
                            data=data, dimD=dimD);
  elseif r_N <= delta/2 then
    dp := flowApproximationAtOrigin(r_N=r_N, V_flow=V_flow,
                                    VDelta_flow=  VDelta_flow, dpDelta=dpDelta,
                                    delta=delta, cBar=cBar);
  else
    dp := Modelica.Fluid.Utilities.regStep(x=r_N-0.75*delta,
                                           y1=performanceCurve(V_flow=V_flow, r_N=r_N, d=d,
                                                               data=data, dimD=dimD),
                                           y2=flowApproximationAtOrigin(r_N=r_N, V_flow=V_flow,
                                                   VDelta_flow=VDelta_flow, dpDelta=dpDelta,
                                                   delta=delta, cBar=cBar),
                                           x_small=delta/4);
  end if;
  dp := dp - V_flow*kRes;
end pressure;

Buildings.Fluid.Movers.BaseClasses.Characteristics.flowApproximationAtOrigin Buildings.Fluid.Movers.BaseClasses.Characteristics.flowApproximationAtOrigin

Approximation for fan or pump pressure raise at origin

Information

This function computes the fan static pressure raise as a function of volume flow rate and revolution near the origin. It is used to avoid a singularity in the pump or fan curve if the revolution approaches zero.

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

Inputs

TypeNameDefaultDescription
VolumeFlowRateV_flow Volumetric flow rate [m3/s]
Realr_N Relative revolution, r_N=N/N_nominal [1]
VolumeFlowRateVDelta_flow Small volume flow rate [m3/s]
PressuredpDelta Small pressure [Pa]
Realdelta Small value used to transition to other fan curve
RealcBar[2] Coefficients for linear approximation of pressure vs. flow rate

Outputs

TypeNameDescription
PressuredpPressure raise [Pa]

Modelica definition

function flowApproximationAtOrigin 
  "Approximation for fan or pump pressure raise at origin"
  extends Modelica.Icons.Function;
  input Modelica.SIunits.VolumeFlowRate V_flow "Volumetric flow rate";
  input Real r_N(unit="1") "Relative revolution, r_N=N/N_nominal";
  input Modelica.SIunits.VolumeFlowRate VDelta_flow "Small volume flow rate";
  input Modelica.SIunits.Pressure dpDelta "Small pressure";
  input Real delta "Small value used to transition to other fan curve";
  input Real cBar[2] 
    "Coefficients for linear approximation of pressure vs. flow rate";
  output Modelica.SIunits.Pressure dp "Pressure raise";
algorithm 
  dp := r_N * dpDelta + r_N^2 * (cBar[1] + cBar[2]*V_flow);
end flowApproximationAtOrigin;

Buildings.Fluid.Movers.BaseClasses.Characteristics.power Buildings.Fluid.Movers.BaseClasses.Characteristics.power

Flow vs. electrical power characteristics for fan or pump

Information

This function computes the fan power consumption for given volume flow rate, speed and performance data. The power consumption is

P = rN3   s(V, d),

where P is the power consumption, rN is the normalized fan speed, V is the volume flow rate and d are performance data for fan or pump power consumption at rN=1.

Implementation

The function s(·, ·) is a cubic hermite spline. If the data d define a monotone decreasing sequence, then s(·, d) is a monotone decreasing function.

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

Inputs

TypeNameDefaultDescription
powerParametersdata Pressure performance data
VolumeFlowRateV_flow Volumetric flow rate [m3/s]
Realr_N Relative revolution, r_N=N/N_nominal [1]
Reald[:] Derivatives at support points for spline interpolation

Outputs

TypeNameDescription
PowerPPower consumption [W]

Modelica definition

function power 
  "Flow vs. electrical power characteristics for fan or pump"
  extends Modelica.Icons.Function;
  input Buildings.Fluid.Movers.BaseClasses.Characteristics.powerParameters data 
    "Pressure performance data";
  input Modelica.SIunits.VolumeFlowRate V_flow "Volumetric flow rate";
  input Real r_N(unit="1") "Relative revolution, r_N=N/N_nominal";
  input Real d[:] "Derivatives at support points for spline interpolation";
  output Modelica.SIunits.Power P "Power consumption";

protected 
   Integer n=size(data.V_flow, 1) "Dimension of data vector";

   Modelica.SIunits.VolumeFlowRate rat "Ratio of V_flow/r_N";
   Integer i "Integer to select data interval";

algorithm 
  if n == 1 then
    P := r_N^3*data.P[1];
  else
    i :=1;
    // Since the coefficients for the spline were evaluated for
    // rat_nominal = V_flow_nominal/r_N_nominal = V_flow_nominal/1, we use
    // V_flow_nominal below
    for j in 1:n-1 loop
       if V_flow > data.V_flow[j] then
         i := j;
       end if;
    end for;
    // Extrapolate or interpolate the data
    P:=r_N^3*Buildings.Utilities.Math.Functions.cubicHermiteLinearExtrapolation(
                x=V_flow,
                x1=data.V_flow[i],
                x2=data.V_flow[i + 1],
                y1=data.P[i],
                y2=data.P[i + 1],
                y1d=d[i],
                y2d=d[i+1]);
  end if;
end power;

Buildings.Fluid.Movers.BaseClasses.Characteristics.efficiency Buildings.Fluid.Movers.BaseClasses.Characteristics.efficiency

Flow vs. efficiency characteristics for fan or pump

Information

This function computes the fan or pump efficiency for given normalized volume flow rate and performance data. The efficiency is

η = s(rV, d),

where η is the efficiency, rV is the normalized volume flow rate, and d are performance data for fan or pump efficiency.

Implementation

The function s(·, ·) is a cubic hermite spline. If the data d define a monotone decreasing sequence, then s(·, d) is a monotone decreasing function.

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

Inputs

TypeNameDefaultDescription
efficiencyParametersdata Efficiency performance data
Realr_V Volumetric flow rate divided by nominal flow rate [1]
Reald[:] Derivatives at support points for spline interpolation

Outputs

TypeNameDescription
RealetaEfficiency [1]

Modelica definition

function efficiency 
  "Flow vs. efficiency characteristics for fan or pump"
  extends Modelica.Icons.Function;
  input Buildings.Fluid.Movers.BaseClasses.Characteristics.efficiencyParameters
    data "Efficiency performance data";
  input Real r_V(unit="1") "Volumetric flow rate divided by nominal flow rate";
  input Real d[:] "Derivatives at support points for spline interpolation";
  output Real eta(min=0, unit="1") "Efficiency";

protected 
  Integer n = size(data.r_V, 1) "Number of data points";
  Integer i "Integer to select data interval";
algorithm 
  if n == 1 then
    eta := data.eta[1];
  else
    i :=1;
    for j in 1:n-1 loop
       if r_V > data.r_V[j] then
         i := j;
       end if;
    end for;
    // Extrapolate or interpolate the data
    eta:=Buildings.Utilities.Math.Functions.cubicHermiteLinearExtrapolation(
                x=r_V,
                x1=data.r_V[i],
                x2=data.r_V[i + 1],
                y1=data.eta[i],
                y2=data.eta[i + 1],
                y1d=d[i],
                y2d=d[i+1]);
  end if;

end efficiency;

Automatically generated Thu Oct 24 15:09:47 2013.