LBL logo

Buildings.Fluid.Movers.BaseClasses

Package with base classes for Buildings.Fluid.Movers

Information

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

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

Package Content

NameDescription
Buildings.Fluid.Movers.BaseClasses.Characteristics Characteristics Functions for fan or pump characteristics
Buildings.Fluid.Movers.BaseClasses.ControlledFlowMachine ControlledFlowMachine Partial model for fan or pump with ideally controlled mass flow rate or head as input signal
Buildings.Fluid.Movers.BaseClasses.PrescribedFlowMachine PrescribedFlowMachine Partial model for fan or pump with speed (y or Nrpm) as input signal
Buildings.Fluid.Movers.BaseClasses.PartialFlowMachine PartialFlowMachine Partial model to interface fan or pump models with the medium
Buildings.Fluid.Movers.BaseClasses.IdealSource IdealSource Base class for pressure and mass flow source with optional power input
PowerInterface Partial model to compute power draw and heat dissipation of fans and pumps
Buildings.Fluid.Movers.BaseClasses.FlowMachineInterface FlowMachineInterface Partial model with performance curves for fans or pumps


Buildings.Fluid.Movers.BaseClasses.ControlledFlowMachine Buildings.Fluid.Movers.BaseClasses.ControlledFlowMachine

Partial model for fan or pump with ideally controlled mass flow rate or head as input signal

Buildings.Fluid.Movers.BaseClasses.ControlledFlowMachine

Information

This model describes a fan or pump that takes as an input the head or the mass flow rate.

Extends from Buildings.Fluid.Movers.BaseClasses.PartialFlowMachine (Partial model to interface fan or pump models with the medium), Buildings.Fluid.Movers.BaseClasses.PowerInterface (Partial model to compute power draw and heat dissipation of fans and pumps).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
BooleanaddPowerToMediumtrueSet to false to avoid any power (=heat and flow work) being added to medium (may give simpler equations)
Densityrho_nominalMedium.density(sta_nominal)Nominal fluid density [kg/m3]
Booleancontrol_m_flow = false to control head instead of m_flow
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
Initialization
MassFlowRatem_flow.start0Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction) [kg/s]
Pressuredp.start0Pressure difference between port_a and port_b [Pa]
Characteristics
Booleanuse_powerCharacteristicfalseUse powerCharacteristic (vs. efficiencyCharacteristic)
BooleanmotorCooledByFluidtrueIf true, then motor heat is added to fluid stream
efficiencyParametersmotorEfficiency Normalized volume flow rate vs. efficiency
efficiencyParametershydraulicEfficiency Normalized volume flow rate vs. efficiency
Dynamics
Equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Formulation of energy balance
DynamicsmassDynamicsenergyDynamicsFormulation of mass balance
BooleandynamicBalancetrueSet to true to use a dynamic balance, which often leads to smaller systems of equations
Nominal condition
Timetau1Time constant of fluid volume for nominal flow, used if dynamicBalance=true [s]
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.)
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_small1E-4*abs(m_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
BooleanhomotopyInitializationtrue= true, use homotopy method
Diagnostics
Booleanshow_V_flowfalse= true, if volume flow rate at inflowing port is computed
Booleanshow_Tfalse= true, if actual temperature at port is computed (may lead to events)

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

Modelica definition

model ControlledFlowMachine 
  "Partial model for fan or pump with ideally controlled mass flow rate or head as input signal"
  extends Buildings.Fluid.Movers.BaseClasses.PartialFlowMachine(
   final show_V_flow = false,
   preSou(final control_m_flow=control_m_flow));

  extends Buildings.Fluid.Movers.BaseClasses.PowerInterface(
     final use_powerCharacteristic = false,
     final rho_nominal = Medium.density(sta_nominal));

  import cha = Buildings.Fluid.Movers.BaseClasses.Characteristics;
//  parameter Medium.MassFlowRate m_flow_nominal
//    "Nominal mass flow rate, used as flow rate if control_m_flow";
//  parameter Modelica.SIunits.MassFlowRate m_flow_max = m_flow_nominal
//    "Maximum mass flow rate (at zero head)";
  // what to control
  constant Boolean control_m_flow "= false to control head instead of m_flow";

  Real r_V(start=1) 
    "Ratio V_flow/V_flow_max = V_flow/V_flow(dp=0, N=N_nominal)";

protected 
  final parameter Medium.AbsolutePressure p_a_nominal(displayUnit="Pa") = Medium.p_default 
    "Nominal inlet pressure for predefined fan or pump characteristics";
  parameter Medium.ThermodynamicState sta_nominal = Medium.setState_pTX(T=T_start,
     p=p_a_nominal, X=X_start[1:Medium.nXi]);

  Modelica.Blocks.Sources.RealExpression PToMedium_flow(y=Q_flow + WFlo) if  addPowerToMedium 
    "Heat and work input into medium";
initial equation 
  V_flow_max=m_flow_nominal/rho_nominal;
equation 
  r_V = VMachine_flow/V_flow_max;
  etaHyd = cha.efficiency(data=hydraulicEfficiency, r_V=r_V, d=hydDer);
  etaMot = cha.efficiency(data=motorEfficiency,     r_V=r_V, d=motDer);
  dpMachine = -dp;
  VMachine_flow = -port_b.m_flow/rho_in;

  connect(PToMedium_flow.y, prePow.Q_flow);
end ControlledFlowMachine;

Buildings.Fluid.Movers.BaseClasses.PrescribedFlowMachine Buildings.Fluid.Movers.BaseClasses.PrescribedFlowMachine

Partial model for fan or pump with speed (y or Nrpm) as input signal

Buildings.Fluid.Movers.BaseClasses.PrescribedFlowMachine

Information

This is the base model for fans and pumps that take as input a control signal in the form of the pump speed Nrpm or the normalized pump speed y=Nrpm/N_nominal.

Extends from Buildings.Fluid.Movers.BaseClasses.FlowMachineInterface (Partial model with performance curves for fans or pumps), Buildings.Fluid.Movers.BaseClasses.PartialFlowMachine (Partial model to interface fan or pump models with the medium).

Parameters

TypeNameDefaultDescription
Densityrho_nominalMedium.density_pTX(Medium.p_...Nominal fluid density [kg/m3]
replaceable package MediumPartialMediumMedium in the component
BooleanaddPowerToMediumtrueSet to false to avoid any power (=heat and flow work) being added to medium (may give simpler equations)
Characteristics
Booleanuse_powerCharacteristicfalseUse powerCharacteristic (vs. efficiencyCharacteristic)
BooleanmotorCooledByFluidtrueIf true, then motor heat is added to fluid stream
efficiencyParametersmotorEfficiency Normalized volume flow rate vs. efficiency
efficiencyParametershydraulicEfficiency Normalized volume flow rate vs. efficiency
AngularVelocity_rpmN_nominal1500Nominal rotational speed for flow characteristic [1/min]
flowParameterspressure Volume flow rate vs. total pressure rise
powerParameterspower Volume flow rate vs. electrical power consumption
Initialization
Realr_V.start1Ratio V_flow/V_flow_max [1]
MassFlowRatem_flow.start0Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction) [kg/s]
Pressuredp.start0Pressure difference between port_a and port_b [Pa]
Nominal condition
MassFlowRatem_flow_nominalmax(pressure.V_flow)*rho_nom...Nominal mass flow rate [kg/s]
Advanced
BooleanhomotopyInitializationtrue= true, use homotopy method
MassFlowRatem_flow_small1E-4*abs(m_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
Diagnostics
Booleanshow_V_flowfalse= true, if volume flow rate at inflowing port is computed
Booleanshow_Tfalse= true, if actual temperature at port is computed (may lead to events)
Dynamics
Filtered speed
BooleanfilteredSpeedtrue= true, if speed is filtered with a 2nd order CriticalDamping filter
TimeriseTime30Rise time of the filter (time to reach 99.6 % of the speed) [s]
InitinitModelica.Blocks.Types.Init.I...Type of initialization (no init/steady state/initial state/initial output)
RealN_start0Initial value of speed
Equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Formulation of energy balance
DynamicsmassDynamicsenergyDynamicsFormulation of mass balance
BooleandynamicBalancetrueSet to true to use a dynamic balance, which often leads to smaller systems of equations
Nominal condition
Timetau1Time constant of fluid volume for nominal flow, used if dynamicBalance=true [s]
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.)
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)

Connectors

TypeNameDescription
output RealOutputN_actual[1/min]
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)
HeatPort_aheatPort 

Modelica definition

partial model PrescribedFlowMachine 
  "Partial model for fan or pump with speed (y or Nrpm) as input signal"
  extends Buildings.Fluid.Movers.BaseClasses.FlowMachineInterface(
    V_flow_max(start=V_flow_nominal),
    final rho_nominal = Medium.density_pTX(Medium.p_default, Medium.T_default, Medium.X_default));

  extends Buildings.Fluid.Movers.BaseClasses.PartialFlowMachine(
      final show_V_flow = false,
      final m_flow_nominal = max(pressure.V_flow)*rho_nominal,
      preSou(final control_m_flow=false));

  // Models
protected 
  Modelica.Blocks.Sources.RealExpression dpMac(y=-dpMachine) 
    "Pressure drop of the pump or fan";

  Modelica.Blocks.Sources.RealExpression PToMedium_flow(y=Q_flow + WFlo) if  addPowerToMedium 
    "Heat and work input into medium";
equation 
 VMachine_flow = -port_b.m_flow/rho;
 rho = rho_in;

  connect(preSou.dp_in, dpMac.y);
  connect(PToMedium_flow.y, prePow.Q_flow);
end PrescribedFlowMachine;

Buildings.Fluid.Movers.BaseClasses.PartialFlowMachine Buildings.Fluid.Movers.BaseClasses.PartialFlowMachine

Partial model to interface fan or pump models with the medium

Buildings.Fluid.Movers.BaseClasses.PartialFlowMachine

Information

This is the base model for fans and pumps. It provides an interface between the equations that compute head and power consumption, and the implementation of the energy and pressure balance of the fluid.

Depending on the value of the parameter dynamicBalance, the fluid volume is computed using a dynamic balance or a steady-state balance.

The parameter addPowerToMedium determines whether any power is added to the fluid. The default is addPowerToMedium=true, and hence the outlet enthalpy is higher than the inlet enthalpy if the flow device is operating. The setting addPowerToMedium=false is physically incorrect (since the flow work, the flow friction and the fan heat do not increase the enthalpy of the medium), but this setting does in some cases lead to simpler equations and more robust simulation, in particular if the mass flow is equal to zero.

Extends from Buildings.Fluid.Interfaces.LumpedVolumeDeclarations (Declarations for lumped volumes), Buildings.Fluid.Interfaces.PartialTwoPortInterface (Partial model transporting fluid between two ports without storing mass or energy).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
BooleanaddPowerToMediumtrueSet to false to avoid any power (=heat and flow work) being added to medium (may give simpler equations)
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
Initialization
MassFlowRatem_flow.start0Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction) [kg/s]
Pressuredp.start0Pressure difference between port_a and port_b [Pa]
Dynamics
Equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Formulation of energy balance
DynamicsmassDynamicsenergyDynamicsFormulation of mass balance
BooleandynamicBalancetrueSet to true to use a dynamic balance, which often leads to smaller systems of equations
Nominal condition
Timetau1Time constant of fluid volume for nominal flow, used if dynamicBalance=true [s]
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.)
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_small1E-4*abs(m_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
BooleanhomotopyInitializationtrue= true, use homotopy method
Diagnostics
Booleanshow_V_flowfalse= true, if volume flow rate at inflowing port is computed
Booleanshow_Tfalse= true, if actual temperature at port is computed (may lead to events)

Connectors

TypeNameDescription
HeatPort_aheatPort 

Modelica definition

partial model PartialFlowMachine 
  "Partial model to interface fan or pump models with the medium"
  extends Buildings.Fluid.Interfaces.LumpedVolumeDeclarations;
  import Modelica.Constants;

  extends Buildings.Fluid.Interfaces.PartialTwoPortInterface(show_T=false,
    port_a(
      h_outflow(start=h_outflow_start),
      final m_flow(min = if allowFlowReversal then -Constants.inf else 0)),
    port_b(
      h_outflow(start=h_outflow_start),
      p(start=p_start),
      final m_flow(max = if allowFlowReversal then +Constants.inf else 0)),
      final showDesignFlowDirection=false);

  Delays.DelayFirstOrder vol(
    redeclare package Medium = Medium,
    tau=tau,
    energyDynamics=if dynamicBalance then energyDynamics else Modelica.Fluid.Types.Dynamics.SteadyState,
    massDynamics=if dynamicBalance then massDynamics else Modelica.Fluid.Types.Dynamics.SteadyState,
    T_start=T_start,
    X_start=X_start,
    C_start=C_start,
    m_flow_nominal=m_flow_nominal,
    p_start=p_start,
    prescribedHeatFlowRate=true,
    allowFlowReversal=allowFlowReversal,
    nPorts=2) "Fluid volume for dynamic model";
   parameter Boolean dynamicBalance = true 
    "Set to true to use a dynamic balance, which often leads to smaller systems of equations";

  parameter Boolean addPowerToMedium=true 
    "Set to false to avoid any power (=heat and flow work) being added to medium (may give simpler equations)";

  parameter Modelica.SIunits.Time tau=1 
    "Time constant of fluid volume for nominal flow, used if dynamicBalance=true";

  // Models
  Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a heatPort;

protected 
  Modelica.SIunits.Density rho_in "Density of inflowing fluid";

  Buildings.Fluid.Movers.BaseClasses.IdealSource preSou(
  redeclare package Medium = Medium,
    allowFlowReversal=allowFlowReversal) "Pressure source";

  Buildings.HeatTransfer.Sources.PrescribedHeatFlow prePow if addPowerToMedium 
    "Prescribed power (=heat and flow work) flow for dynamic model";

  parameter Medium.ThermodynamicState sta_start=Medium.setState_pTX(
      T=T_start, p=p_start, X=X_start);
  parameter Modelica.SIunits.SpecificEnthalpy h_outflow_start = Medium.specificEnthalpy(sta_start) 
    "Start value for outflowing enthalpy";

equation 
  // For computing the density, we assume that the fan operates in the design flow direction.
  rho_in = Medium.density(
       Medium.setState_phX(port_a.p, inStream(port_a.h_outflow), inStream(port_a.Xi_outflow)));
  connect(prePow.port, vol.heatPort);

  connect(vol.heatPort, heatPort);
  connect(port_a, vol.ports[1]);
  connect(vol.ports[2], preSou.port_a);
  connect(preSou.port_b, port_b);
end PartialFlowMachine;

Buildings.Fluid.Movers.BaseClasses.IdealSource Buildings.Fluid.Movers.BaseClasses.IdealSource

Base class for pressure and mass flow source with optional power input

Buildings.Fluid.Movers.BaseClasses.IdealSource

Information

Model of a fictious pipe that is used as a base class for a pressure source or to prescribe a mass flow rate.

Note that for fans and pumps with dynamic balance, both the heat and the flow work are added to the volume of air or water. This simplifies the equations compared to adding heat to the volume, and flow work to this model.

Extends from Modelica.Fluid.Interfaces.PartialTwoPortTransport (Partial element transporting fluid between two ports without storage of mass or energy).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Booleancontrol_m_flow = false to control dp instead of m_flow
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_smallsystem.m_flow_smallSmall mass flow rate for regularization of zero flow [kg/s]
Diagnostics
Booleanshow_Tfalse= true, if temperatures at port_a and port_b are computed
Booleanshow_V_flowfalse= 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 RealInputm_flow_inPrescribed mass flow rate
input RealInputdp_inPrescribed outlet pressure

Modelica definition

model IdealSource 
  "Base class for pressure and mass flow source with optional power input"
  extends Modelica.Fluid.Interfaces.PartialTwoPortTransport(show_V_flow=false,
                                                            show_T=false);

  // what to control
  parameter Boolean control_m_flow "= false to control dp instead of m_flow";
  Modelica.Blocks.Interfaces.RealInput m_flow_in if control_m_flow 
    "Prescribed mass flow rate";
  Modelica.Blocks.Interfaces.RealInput dp_in if not control_m_flow 
    "Prescribed outlet pressure";

protected 
  Modelica.Blocks.Interfaces.RealInput m_flow_internal 
    "Needed to connect to conditional connector";
  Modelica.Blocks.Interfaces.RealInput dp_internal 
    "Needed to connect to conditional connector";
equation 

  // Ideal control
  if control_m_flow then
    m_flow = m_flow_internal;
    dp_internal = 0;
  else
    dp = dp_internal;
    m_flow_internal = 0;
  end if;

  connect(dp_internal, dp_in);
  connect(m_flow_internal, m_flow_in);

  // Energy balance (no storage)
  port_a.h_outflow = inStream(port_b.h_outflow);
  port_b.h_outflow = inStream(port_a.h_outflow);

end IdealSource;

Buildings.Fluid.Movers.BaseClasses.PowerInterface

Partial model to compute power draw and heat dissipation of fans and pumps

Buildings.Fluid.Movers.BaseClasses.PowerInterface

Information

This is an interface that implements the functions to compute the power draw and the heat dissipation of fans and pumps. It is used by the model Buildings.Fluid.Movers.BaseClasses.FlowMachineInterface.

Parameters

TypeNameDefaultDescription
Densityrho_nominal Nominal fluid density [kg/m3]
Characteristics
Booleanuse_powerCharacteristicfalseUse powerCharacteristic (vs. efficiencyCharacteristic)
BooleanmotorCooledByFluidtrueIf true, then motor heat is added to fluid stream
efficiencyParametersmotorEfficiency Normalized volume flow rate vs. efficiency
efficiencyParametershydraulicEfficiency Normalized volume flow rate vs. efficiency
Advanced
BooleanhomotopyInitializationtrue= true, use homotopy method

Modelica definition

partial model PowerInterface 
  "Partial model to compute power draw and heat dissipation of fans and pumps"

  import Modelica.Constants;

  parameter Boolean use_powerCharacteristic = false 
    "Use powerCharacteristic (vs. efficiencyCharacteristic)";

  parameter Boolean motorCooledByFluid = true 
    "If true, then motor heat is added to fluid stream";
  parameter Boolean homotopyInitialization = true "= true, use homotopy method";

  parameter Buildings.Fluid.Movers.BaseClasses.Characteristics.efficiencyParameters
      motorEfficiency(r_V={1}, eta={0.7}) 
    "Normalized volume flow rate vs. efficiency";
  parameter Buildings.Fluid.Movers.BaseClasses.Characteristics.efficiencyParameters
      hydraulicEfficiency(r_V={1}, eta={0.7}) 
    "Normalized volume flow rate vs. efficiency";

  parameter Modelica.SIunits.Density rho_nominal "Nominal fluid density";

  Modelica.SIunits.Power PEle "Electrical power input";
  Modelica.SIunits.Power WHyd 
    "Hydraulic power input (converted to flow work and heat)";
  Modelica.SIunits.Power WFlo "Flow work";
  Modelica.SIunits.HeatFlowRate Q_flow "Heat input from fan or pump to medium";
  Real eta(min=0, max=1) "Global efficiency";
  Real etaHyd(min=0, max=1) "Hydraulic efficiency";
  Real etaMot(min=0, max=1) "Motor efficiency";

  Modelica.SIunits.Pressure dpMachine(displayUnit="Pa") "Pressure increase";
  Modelica.SIunits.VolumeFlowRate VMachine_flow "Volume flow rate";
protected 
  parameter Modelica.SIunits.VolumeFlowRate V_flow_max(fixed=false) 
    "Maximum volume flow rate, used for smoothing";
  //Modelica.SIunits.HeatFlowRate QThe_flow "Heat input into the medium";
  parameter Modelica.SIunits.VolumeFlowRate delta_V_flow = 1E-3*V_flow_max 
    "Factor used for setting heat input into medium to zero at very small flows";
  final parameter Real motDer[size(motorEfficiency.r_V, 1)](fixed=false) 
    "Coefficients for polynomial of pressure vs. flow rate";
  final parameter Real hydDer[size(hydraulicEfficiency.r_V,1)](fixed=false) 
    "Coefficients for polynomial of pressure vs. flow rate";

  Modelica.SIunits.HeatFlowRate QThe_flow 
    "Heat input from fan or pump to medium";

initial algorithm 
 // Compute derivatives for cubic spline
 motDer :=
   if use_powerCharacteristic then 
     zeros(size(motorEfficiency.r_V, 1))
   elseif ( size(motorEfficiency.r_V, 1) == 1)  then 
       {0}
   else 
      Buildings.Utilities.Math.Functions.splineDerivatives(
      x=motorEfficiency.r_V,
      y=motorEfficiency.eta);
  hydDer :=
     if use_powerCharacteristic then 
       zeros(size(hydraulicEfficiency.r_V, 1))
     elseif ( size(hydraulicEfficiency.r_V, 1) == 1)  then 
       {0}
     else 
       Buildings.Utilities.Math.Functions.splineDerivatives(
                   x=hydraulicEfficiency.r_V,
                   y=hydraulicEfficiency.eta);
equation 
  eta = etaHyd * etaMot;
  WFlo = eta * PEle;
  // Flow work
  WFlo = dpMachine*VMachine_flow;
  // Hydraulic power (transmitted by shaft), etaHyd = WFlo/WHyd
  etaHyd * WHyd   = WFlo;
  // Heat input into medium
  QThe_flow +  WFlo = if motorCooledByFluid then PEle else WHyd;
  // At m_flow = 0, the solver may still obtain positive values for QThe_flow.
  // The next statement sets the heat input into the medium to zero for very small flow rates.
  if homotopyInitialization then
    Q_flow = homotopy(actual=Buildings.Utilities.Math.Functions.spliceFunction(pos=QThe_flow, neg=0,
                       x=noEvent(abs(VMachine_flow))-2*delta_V_flow, deltax=delta_V_flow),
                     simplified=0);
  else
    Q_flow = Buildings.Utilities.Math.Functions.spliceFunction(pos=QThe_flow, neg=0,
                       x=noEvent(abs(VMachine_flow))-2*delta_V_flow, deltax=delta_V_flow);
  end if;
end PowerInterface;

Buildings.Fluid.Movers.BaseClasses.FlowMachineInterface Buildings.Fluid.Movers.BaseClasses.FlowMachineInterface

Partial model with performance curves for fans or pumps

Buildings.Fluid.Movers.BaseClasses.FlowMachineInterface

Information

This is an interface that implements the functions to compute the head, power draw and efficiency of fans and pumps. It is used by the model PrescribedFlowMachine.

The nominal hydraulic characteristic (volume flow rate versus total pressure) is given by a set of data points. A cubic hermite spline with linear extrapolation is used to compute the performance at other operating points.

The fan or pump energy balance can be specified in two alternative ways: