Buildings.Fluid.Sensors

Package with sensor models

Information


Package Sensors consists of idealized sensor components that provide variables of a medium model and/or fluid ports as output signals. These signals can be, e.g., further processed with components of the Modelica.Blocks library. Also more realistic sensor models can be built, by further processing (e.g., by attaching block Modelica.Blocks.FirstOrder to model the time constant of the sensor).

For the thermodynamic state variables temperature, specific entalpy, specific entropy and density the fluid library provides two different types of sensors: regular one port and two port sensors.

Package Content

NameDescription
Buildings.Fluid.Sensors.Conversions Conversions Package with conversions for sensor models
Buildings.Fluid.Sensors.Density Density Ideal one port density sensor
Buildings.Fluid.Sensors.DensityTwoPort DensityTwoPort Ideal two port density sensor
Buildings.Fluid.Sensors.EnthalpyFlowRate EnthalpyFlowRate Ideal enthalphy flow rate sensor
Buildings.Fluid.Sensors.LatentEnthalpyFlowRate LatentEnthalpyFlowRate Ideal enthalphy flow rate sensor that outputs the latent enthalpy flow rate only
Buildings.Fluid.Sensors.MassFlowRate MassFlowRate Ideal sensor for mass flow rate
Buildings.Fluid.Sensors.MassFraction MassFraction Ideal one port mass fraction sensor
Buildings.Fluid.Sensors.Pressure Pressure Ideal pressure sensor
Buildings.Fluid.Sensors.RelativeHumidity RelativeHumidity Ideal one port relative humidity sensor
Buildings.Fluid.Sensors.RelativePressure RelativePressure Ideal relative pressure sensor
Buildings.Fluid.Sensors.RelativeTemperature RelativeTemperature Ideal relative temperature sensor
Buildings.Fluid.Sensors.SensibleEnthalpyFlowRate SensibleEnthalpyFlowRate Ideal enthalphy flow rate sensor that outputs the sensible enthalpy flow rate only
Buildings.Fluid.Sensors.SpecificEnthalpy SpecificEnthalpy Ideal one port specific enthalpy sensor
Buildings.Fluid.Sensors.SpecificEnthalpyTwoPort SpecificEnthalpyTwoPort Ideal two port sensor for the specific enthalpy
Buildings.Fluid.Sensors.SpecificEntropy SpecificEntropy Ideal one port specific entropy sensor
Buildings.Fluid.Sensors.SpecificEntropyTwoPort SpecificEntropyTwoPort Ideal two port sensor for the specific entropy
Buildings.Fluid.Sensors.Temperature Temperature Ideal one port temperature sensor
Buildings.Fluid.Sensors.TemperatureDryBulbDynamic TemperatureDryBulbDynamic Ideal temperature sensor
Buildings.Fluid.Sensors.TemperatureTwoPort TemperatureTwoPort Ideal two port temperature sensor
Buildings.Fluid.Sensors.TemperatureWetBulb TemperatureWetBulb Ideal wet bulb temperature sensor
Buildings.Fluid.Sensors.TraceSubstances TraceSubstances Ideal one port trace substances sensor
Buildings.Fluid.Sensors.TraceSubstancesTwoPort TraceSubstancesTwoPort Ideal two port sensor for trace substance
Buildings.Fluid.Sensors.VolumeFlowRate VolumeFlowRate Ideal sensor for volume flow rate
Buildings.Fluid.Sensors.Examples Examples Collection of models that illustrate model use and test models


Buildings.Fluid.Sensors.Density Buildings.Fluid.Sensors.Density

Ideal one port density sensor

Buildings.Fluid.Sensors.Density

Information


This component monitors the density of the fluid passing its port. The sensor is ideal, i.e. it does not influence the fluid.

If using the one port sensor please read the Information first.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor (Partial component to model a sensor that measures a potential variable), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the sensor

Connectors

TypeNameDescription
FluidPort_aport 
output RealOutputdDensity in port medium [kg/m3]

Modelica definition

model Density "Ideal one port density sensor"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor;
  extends Modelica.Icons.RotationalSensor;
  Modelica.Blocks.Interfaces.RealOutput d(final quantity="Density",
                                          final unit="kg/m3",
                                          min=0) "Density in port medium";

equation 
  d = Medium.density(Medium.setState_phX(port.p, inStream(port.h_outflow), inStream(port.Xi_outflow)));
end Density;

Buildings.Fluid.Sensors.DensityTwoPort Buildings.Fluid.Sensors.DensityTwoPort

Ideal two port density sensor

Buildings.Fluid.Sensors.DensityTwoPort

Information


This component monitors the density of the fluid flowing from port_a to port_b. The sensor is ideal, i.e. it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor (Partial component to model sensors that measure flow properties), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_smallsystem.m_flow_smallFor bi-directional flow, density is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required) [kg/s]

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)
output RealOutputdDensity of the passing fluid [kg/m3]

Modelica definition

model DensityTwoPort "Ideal two port density sensor"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor;
  extends Modelica.Icons.RotationalSensor;
  Modelica.Blocks.Interfaces.RealOutput d(final quantity="Density",
                                          final unit="kg/m3",
                                          min=0) "Density of the passing fluid";
  parameter Medium.MassFlowRate m_flow_small(min=0) = system.m_flow_small 
    "For bi-directional flow, density is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required)";
protected 
  Medium.Density rho_a_inflow "Density of inflowing fluid at port_a";
  Medium.Density rho_b_inflow 
    "Density of inflowing fluid at port_b or rho_a_inflow, if uni-directional flow";
equation 
  if allowFlowReversal then
     rho_a_inflow = Medium.density(Medium.setState_phX(port_b.p, port_b.h_outflow, port_b.Xi_outflow));
     rho_b_inflow = Medium.density(Medium.setState_phX(port_a.p, port_a.h_outflow, port_a.Xi_outflow));
     d = Modelica.Fluid.Utilities.regStep(port_a.m_flow, rho_a_inflow, rho_b_inflow, m_flow_small);
  else
     d = Medium.density(Medium.setState_phX(port_b.p, port_b.h_outflow, port_b.Xi_outflow));
     rho_a_inflow = d;
     rho_b_inflow = d;
  end if;
end DensityTwoPort;

Buildings.Fluid.Sensors.EnthalpyFlowRate Buildings.Fluid.Sensors.EnthalpyFlowRate

Ideal enthalphy flow rate sensor

Buildings.Fluid.Sensors.EnthalpyFlowRate

Information


This component monitors the enthalphy flow rate of the medium in the flow between fluid ports. The sensor is ideal, i.e., it does not influence the fluid.

For a sensor that measures the latent enthalpy flow rate, use Buildings.Fluid.Sensors.LatentEnthalpyFlowRate.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor (Partial component to model sensors that measure flow properties), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)

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)
output RealOutputH_flowEnthalpy flow rate, positive if from port_a to port_b [W]

Modelica definition

model EnthalpyFlowRate "Ideal enthalphy flow rate sensor"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor;
  extends Modelica.Icons.RotationalSensor;
  Modelica.Blocks.Interfaces.RealOutput H_flow(unit="W") 
    "Enthalpy flow rate, positive if from port_a to port_b";

equation 
  if allowFlowReversal then
     H_flow = port_a.m_flow * actualStream(port_a.h_outflow);
  else
     H_flow = port_a.m_flow * port_b.h_outflow;
  end if;
end EnthalpyFlowRate;

Buildings.Fluid.Sensors.LatentEnthalpyFlowRate Buildings.Fluid.Sensors.LatentEnthalpyFlowRate

Ideal enthalphy flow rate sensor that outputs the latent enthalpy flow rate only

Buildings.Fluid.Sensors.LatentEnthalpyFlowRate

Information


This component monitors the latent enthalphy flow rate of the medium in the flow between fluid ports. In particular, if the enthalpy flow rate is

  HTotal_flow = HSensible_flow + HLatent_flow,
where HSensible_flow = m_flow * (1-X[water]) * cp_air, then this sensor outputs H_flow = HLatent_flow. This sensor can only be used with medium models that implement the function enthalpyOfCondensingGas(state).

For a sensor that measures HTotal_flow, use Buildings.Fluid.Sensors.EnthalpyFlowRate.
For a sensor that measures HSensible_flow, use Buildings.Fluid.Sensors.SensibleEnthalpyFlowRate.

The sensor is ideal, i.e., it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor (Partial component to model sensors that measure flow properties), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)

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)
output RealOutputH_flowLatent enthalpy flow rate, positive if from port_a to port_b [W]

Modelica definition

model LatentEnthalpyFlowRate 
  "Ideal enthalphy flow rate sensor that outputs the latent enthalpy flow rate only"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor;
  extends Modelica.Icons.RotationalSensor;
  // redeclare Medium with a more restricting base class. This improves the error
  // message if a user selects a medium that does not contain the function
  // enthalpyOfLiquid(.)
  replaceable package Medium =
      Modelica.Media.Interfaces.PartialCondensingGases;

  Modelica.Blocks.Interfaces.RealOutput H_flow(unit="W") 
    "Latent enthalpy flow rate, positive if from port_a to port_b";
  Medium.MassFraction XiActual[Medium.nXi] "Medium mass fraction in port_a";
  Medium.SpecificEnthalpy hActual "Medium enthalpy in port_a";
  Medium.ThermodynamicState sta "Medium properties in port_a";

protected 
  parameter Integer i_w(min=1, fixed=false) "Index for water substance";
initial algorithm 
  i_w :=1;
    for i in 1:Medium.nXi loop
      if Modelica.Utilities.Strings.isEqual(Medium.substanceNames[i], "Water") then
        i_w :=i;
      end if;
    end for;
equation 
  if allowFlowReversal then
     XiActual = actualStream(port_a.Xi_outflow);
     hActual = actualStream(port_a.h_outflow);
     sta = Medium.setState_phX(port_a.p, hActual, XiActual);
  else
     XiActual = port_b.Xi_outflow;
     hActual = port_b.h_outflow;
     sta = Medium.setState_phX(port_b.p, port_b.h_outflow, XiActual);
  end if;
     // Compute H_flow as difference between total enthalpy and enthalpy on non-condensing gas.
     // This is needed to compute the liquid vs. gas fraction of water, using the equations
     // provided by the medium model
     H_flow = port_a.m_flow * (hActual -
       (1-XiActual[i_w]) * Medium.enthalpyOfNonCondensingGas(Medium.temperature(sta)));
end LatentEnthalpyFlowRate;

Buildings.Fluid.Sensors.MassFlowRate Buildings.Fluid.Sensors.MassFlowRate

Ideal sensor for mass flow rate

Buildings.Fluid.Sensors.MassFlowRate

Information


This component monitors the mass flow rate flowing from port_a to port_b. The sensor is ideal, i.e., it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor (Partial component to model sensors that measure flow properties), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)

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)
output RealOutputm_flowMass flow rate from port_a to port_b [kg/s]

Modelica definition

model MassFlowRate "Ideal sensor for mass flow rate"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor;
  extends Modelica.Icons.RotationalSensor;
  Modelica.Blocks.Interfaces.RealOutput m_flow(quantity="MassFlowRate",
                                               final unit="kg/s") 
    "Mass flow rate from port_a to port_b";

equation 
  m_flow = port_a.m_flow;
end MassFlowRate;

Buildings.Fluid.Sensors.MassFraction Buildings.Fluid.Sensors.MassFraction

Ideal one port mass fraction sensor

Buildings.Fluid.Sensors.MassFraction

Information


This component monitors the mass fraction contained in the fluid passing its port. The sensor is ideal, i.e. it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor (Partial component to model a sensor that measures a potential variable), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the sensor
StringsubstanceName"water"Name of species substance

Connectors

TypeNameDescription
FluidPort_aport 
output RealOutputXMass fraction in port medium

Modelica definition

model MassFraction "Ideal one port mass fraction sensor"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor;
  extends Modelica.Icons.RotationalSensor;
  parameter String substanceName = "water" "Name of species substance";

  Modelica.Blocks.Interfaces.RealOutput X "Mass fraction in port medium";

protected 
  parameter Integer ind(fixed=false) 
    "Index of species in vector of auxiliary substances";
  Medium.MassFraction XiVec[Medium.nXi](
      quantity=Medium.extraPropertiesNames) 
    "Trace substances vector, needed because indexed argument for the operator inStream is not supported";
initial algorithm 
  ind:= -1;
  for i in 1:Medium.nX loop
    if ( Modelica.Utilities.Strings.isEqual(Medium.substanceNames[i], substanceName)) then
      ind := i;
    end if;
  end for;
  assert(ind > 0, "Species with name '" + substanceName + "' is not present in medium '"
         + Medium.mediumName + "'.\n"
         + "Check sensor parameter and medium model.");
equation 
  XiVec = inStream(port.Xi_outflow);
  X = if ind > Medium.nXi then (1-sum(XiVec)) else XiVec[ind];
end MassFraction;

Buildings.Fluid.Sensors.Pressure Buildings.Fluid.Sensors.Pressure

Ideal pressure sensor

Buildings.Fluid.Sensors.Pressure

Information


This component monitors the absolute pressure at its fluid port. The sensor is ideal, i.e., it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor (Partial component to model a sensor that measures a potential variable), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the sensor

Connectors

TypeNameDescription
FluidPort_aport 
output RealOutputpPressure at port [Pa]

Modelica definition

model Pressure "Ideal pressure sensor"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor;
  extends Modelica.Icons.RotationalSensor;
  Modelica.Blocks.Interfaces.RealOutput p(final quantity="Pressure",
                                          final unit="Pa",
                                          min=0) "Pressure at port";
equation 
  p = port.p;
end Pressure;

Buildings.Fluid.Sensors.RelativeHumidity Buildings.Fluid.Sensors.RelativeHumidity

Ideal one port relative humidity sensor

Buildings.Fluid.Sensors.RelativeHumidity

Information


This component monitors the relative humidity contained in the fluid passing its port. The sensor is ideal, i.e. it does not influence the fluid.

Note that this sensor can only be used with media that contain the variable phi, which is typically the case for moist air models.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor (Partial component to model a sensor that measures a potential variable), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the sensor

Connectors

TypeNameDescription
FluidPort_aport 
output RealOutputphiRelative humidity in port medium

Modelica definition

model RelativeHumidity "Ideal one port relative humidity sensor"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor;
  extends Modelica.Icons.RotationalSensor;

  Modelica.Blocks.Interfaces.RealOutput phi "Relative humidity in port medium";

  Medium.BaseProperties med "Medium state at dry bulb temperature";

equation 
  med.p = port.p;
  med.h = inStream(port.h_outflow);
  med.Xi = inStream(port.Xi_outflow);
  phi = med.phi;

end RelativeHumidity;

Buildings.Fluid.Sensors.RelativePressure Buildings.Fluid.Sensors.RelativePressure

Ideal relative pressure sensor

Buildings.Fluid.Sensors.RelativePressure

Information


The relative pressure "port_a.p - port_b.p" is determined between the two ports of this component and is provided as output signal. The sensor should be connected in parallel with other equipment, no flow through the sensor is allowed.

Extends from Modelica.Icons.TranslationalSensor (Icon representing translational measurement device).

Connectors

TypeNameDescription
FluidPort_aport_a 
FluidPort_bport_b 
output RealOutputp_relRelative pressure signal [Pa]

Modelica definition

model RelativePressure "Ideal relative pressure sensor"
  extends Modelica.Icons.TranslationalSensor;
  replaceable package Medium =
    Modelica.Media.Interfaces.PartialMedium "Medium in the sensor";

  Modelica.Fluid.Interfaces.FluidPort_a port_a(m_flow(min=0),
                                p(start=Medium.p_default),
                                redeclare package Medium = Medium);
  Modelica.Fluid.Interfaces.FluidPort_b port_b(m_flow(min=0),
                                p(start=Medium.p_default),
                                redeclare package Medium = Medium);

  Modelica.Blocks.Interfaces.RealOutput p_rel(final quantity="Pressure",
                                              final unit="Pa",
                                              displayUnit="bar") 
    "Relative pressure signal";
equation 
  // Zero flow equations for connectors
  port_a.m_flow = 0;
  port_b.m_flow = 0;

  // No contribution of specific quantities
  port_a.h_outflow = 0;
  port_b.h_outflow = 0;
  port_a.Xi_outflow = zeros(Medium.nXi);
  port_b.Xi_outflow = zeros(Medium.nXi);
  port_a.C_outflow  = zeros(Medium.nC);
  port_b.C_outflow  = zeros(Medium.nC);

  // Relative pressure
  p_rel = port_a.p - port_b.p;
end RelativePressure;

Buildings.Fluid.Sensors.RelativeTemperature Buildings.Fluid.Sensors.RelativeTemperature

Ideal relative temperature sensor

Buildings.Fluid.Sensors.RelativeTemperature

Information


The relative temperature "T(port_a) - T(port_b)" is determined between the two ports of this component and is provided as output signal. The sensor should be connected in parallel with other equipment, no flow through the sensor is allowed.

Extends from Modelica.Icons.TranslationalSensor (Icon representing translational measurement device).

Connectors

TypeNameDescription
FluidPort_aport_a 
FluidPort_bport_b 
output RealOutputT_relRelative temperature signal [K]

Modelica definition

model RelativeTemperature "Ideal relative temperature sensor"
  extends Modelica.Icons.TranslationalSensor;
  replaceable package Medium =
    Modelica.Media.Interfaces.PartialMedium "Medium in the sensor";
  Modelica.Fluid.Interfaces.FluidPort_a port_a(m_flow(min=0),
                                p(start=Medium.p_default),
                                redeclare package Medium = Medium);
  Modelica.Fluid.Interfaces.FluidPort_b port_b(m_flow(min=0),
                                p(start=Medium.p_default),
                                redeclare package Medium = Medium);

  Modelica.Blocks.Interfaces.RealOutput T_rel(final quantity="ThermodynamicTemperature",
                                              final unit = "K", displayUnit = "degC", min=0) 
    "Relative temperature signal";
equation 
  // Zero flow equations for connectors
  port_a.m_flow = 0;
  port_b.m_flow = 0;

  // No contribution of specific quantities
  port_a.h_outflow = 0;
  port_b.h_outflow = 0;
  port_a.Xi_outflow = zeros(Medium.nXi);
  port_b.Xi_outflow = zeros(Medium.nXi);
  port_a.C_outflow  = zeros(Medium.nC);
  port_b.C_outflow  = zeros(Medium.nC);

  // Relative temperature
  T_rel = Medium.temperature(Medium.setState_phX(port_a.p, inStream(port_a.h_outflow), inStream(port_a.Xi_outflow))) -
          Medium.temperature(Medium.setState_phX(port_b.p, inStream(port_b.h_outflow), inStream(port_b.Xi_outflow)));
end RelativeTemperature;

Buildings.Fluid.Sensors.SensibleEnthalpyFlowRate Buildings.Fluid.Sensors.SensibleEnthalpyFlowRate

Ideal enthalphy flow rate sensor that outputs the sensible enthalpy flow rate only

Buildings.Fluid.Sensors.SensibleEnthalpyFlowRate

Information


This component monitors the sensible enthalphy flow rate of the medium in the flow between fluid ports. In particular, if the enthalpy flow rate is

  HTotal_flow = HSensible_flow + HLatent_flow,
where HSensible_flow = m_flow * (1-X[water]) * cp_air, then this sensor outputs H_flow = HSensible_flow. This sensor can only be used with medium models that implement the function enthalpyOfNonCondensingGas(state).

For a sensor that measures HTotal_flow, use Buildings.Fluid.Sensors.EnthalpyFlowRate.
For a sensor that measures HLatent_flow, use Buildings.Fluid.Sensors.LatentEnthalpyFlowRate.

The sensor is ideal, i.e., it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor (Partial component to model sensors that measure flow properties), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)

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)
output RealOutputH_flowSensible enthalpy flow rate, positive if from port_a to port_b [W]

Modelica definition

model SensibleEnthalpyFlowRate 
  "Ideal enthalphy flow rate sensor that outputs the sensible enthalpy flow rate only"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor;
  extends Modelica.Icons.RotationalSensor;
  // redeclare Medium with a more restricting base class. This improves the error
  // message if a user selects a medium that does not contain the function
  // enthalpyOfLiquid(.)
  replaceable package Medium =
      Modelica.Media.Interfaces.PartialCondensingGases;

  Modelica.Blocks.Interfaces.RealOutput H_flow(unit="W") 
    "Sensible enthalpy flow rate, positive if from port_a to port_b";
  Medium.MassFraction XiActual[Medium.nXi] "Medium mass fraction in port_a";
  Medium.ThermodynamicState sta "Medium properties in port_a";

protected 
  parameter Integer i_w(min=1, fixed=false) "Index for water substance";
initial algorithm 
  i_w :=1;
    for i in 1:Medium.nXi loop
      if Modelica.Utilities.Strings.isEqual(Medium.substanceNames[i], "Water") then
        i_w :=i;
      end if;
    end for;
equation 
  if allowFlowReversal then
     XiActual = actualStream(port_a.Xi_outflow);
     sta = Medium.setState_phX(port_a.p, actualStream(port_a.h_outflow), XiActual);
  else
     XiActual = port_b.Xi_outflow;
     sta = Medium.setState_phX(port_b.p, port_b.h_outflow, XiActual);
  end if;
     H_flow = port_a.m_flow * (1-XiActual[i_w]) * Medium.enthalpyOfNonCondensingGas(Medium.temperature(sta));
end SensibleEnthalpyFlowRate;

Buildings.Fluid.Sensors.SpecificEnthalpy Buildings.Fluid.Sensors.SpecificEnthalpy

Ideal one port specific enthalpy sensor

Buildings.Fluid.Sensors.SpecificEnthalpy

Information


This component monitors the specific enthalpy of the fluid passing its port. The sensor is ideal, i.e. it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor (Partial component to model a sensor that measures a potential variable), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the sensor

Connectors

TypeNameDescription
FluidPort_aport 
output RealOutputh_outSpecific enthalpy in port medium [J/kg]

Modelica definition

model SpecificEnthalpy "Ideal one port specific enthalpy sensor"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor;
  extends Modelica.Icons.RotationalSensor;
  Modelica.Blocks.Interfaces.RealOutput h_out(final quantity="SpecificEnergy",
                                              final unit="J/kg") 
    "Specific enthalpy in port medium";

equation 
  h_out = inStream(port.h_outflow);
end SpecificEnthalpy;

Buildings.Fluid.Sensors.SpecificEnthalpyTwoPort Buildings.Fluid.Sensors.SpecificEnthalpyTwoPort

Ideal two port sensor for the specific enthalpy

Buildings.Fluid.Sensors.SpecificEnthalpyTwoPort

Information


This component monitors the specific enthalpy of a passing fluid. The sensor is ideal, i.e. it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor (Partial component to model sensors that measure flow properties), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_smallsystem.m_flow_smallFor bi-directional flow, specific enthalpy is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required) [kg/s]

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)
output RealOutputh_outSpecific enthalpy of the passing fluid [J/kg]

Modelica definition

model SpecificEnthalpyTwoPort 
  "Ideal two port sensor for the specific enthalpy"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor;
  extends Modelica.Icons.RotationalSensor;
  Modelica.Blocks.Interfaces.RealOutput h_out(final quantity="SpecificEnergy",
                                              final unit="J/kg") 
    "Specific enthalpy of the passing fluid";

  parameter Medium.MassFlowRate m_flow_small(min=0) = system.m_flow_small 
    "For bi-directional flow, specific enthalpy is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required)";

equation 
  if allowFlowReversal then
     h_out = Modelica.Fluid.Utilities.regStep(port_a.m_flow, port_b.h_outflow, port_a.h_outflow, m_flow_small);
  else
     h_out = port_b.h_outflow;
  end if;
end SpecificEnthalpyTwoPort;

Buildings.Fluid.Sensors.SpecificEntropy Buildings.Fluid.Sensors.SpecificEntropy

Ideal one port specific entropy sensor

Buildings.Fluid.Sensors.SpecificEntropy

Information


This component monitors the specific entropy of the fluid passing its port. The sensor is ideal, i.e. it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor (Partial component to model a sensor that measures a potential variable), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the sensor

Connectors

TypeNameDescription
FluidPort_aport 
output RealOutputsSpecific entropy in port medium [J/(kg.K)]

Modelica definition

model SpecificEntropy "Ideal one port specific entropy sensor"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor;
  extends Modelica.Icons.RotationalSensor;
  Modelica.Blocks.Interfaces.RealOutput s(final quantity="SpecificEntropy",
                                          final unit="J/(kg.K)") 
    "Specific entropy in port medium";

equation 
  s = Medium.specificEntropy(Medium.setState_phX(port.p, inStream(port.h_outflow), inStream(port.Xi_outflow)));
end SpecificEntropy;

Buildings.Fluid.Sensors.SpecificEntropyTwoPort Buildings.Fluid.Sensors.SpecificEntropyTwoPort

Ideal two port sensor for the specific entropy

Buildings.Fluid.Sensors.SpecificEntropyTwoPort

Information


This component monitors the specific entropy of the passing fluid. The sensor is ideal, i.e. it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor (Partial component to model sensors that measure flow properties), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_smallsystem.m_flow_smallFor bi-directional flow, specific entropy is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required) [kg/s]

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)
output RealOutputsSpecific entropy of the passing fluid [J/(kg.K)]

Modelica definition

model SpecificEntropyTwoPort 
  "Ideal two port sensor for the specific entropy"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor;
  extends Modelica.Icons.RotationalSensor;
  Modelica.Blocks.Interfaces.RealOutput s(final quantity="SpecificEntropy",
                                          final unit="J/(kg.K)") 
    "Specific entropy of the passing fluid";
  parameter Medium.MassFlowRate m_flow_small(min=0) = system.m_flow_small 
    "For bi-directional flow, specific entropy is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required)";

protected 
  Medium.SpecificEntropy s_a_inflow 
    "Specific entropy of inflowing fluid at port_a";
  Medium.SpecificEntropy s_b_inflow 
    "Specific entropy of inflowing fluid at port_b or s_a_inflow, if uni-directional flow";
equation 
  if allowFlowReversal then
     s_a_inflow = Medium.specificEntropy(Medium.setState_phX(port_b.p, port_b.h_outflow, port_b.Xi_outflow));
     s_b_inflow = Medium.specificEntropy(Medium.setState_phX(port_a.p, port_a.h_outflow, port_a.Xi_outflow));
     s = Modelica.Fluid.Utilities.regStep(port_a.m_flow, s_a_inflow, s_b_inflow, m_flow_small);
  else
     s = Medium.specificEntropy(Medium.setState_phX(port_b.p, port_b.h_outflow, port_b.Xi_outflow));
     s_a_inflow = s;
     s_b_inflow = s;
  end if;
end SpecificEntropyTwoPort;

Buildings.Fluid.Sensors.Temperature Buildings.Fluid.Sensors.Temperature

Ideal one port temperature sensor

Buildings.Fluid.Sensors.Temperature

Information


This component monitors the temperature of the fluid passing its port. The sensor is ideal, i.e. it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor (Partial component to model a sensor that measures a potential variable).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the sensor

Connectors

TypeNameDescription
FluidPort_aport 
output RealOutputTTemperature in port medium [K]

Modelica definition

model Temperature "Ideal one port temperature sensor"
    extends Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor;

  Modelica.Blocks.Interfaces.RealOutput T(final quantity="ThermodynamicTemperature",
                                          final unit = "K", min=0) 
    "Temperature in port medium";

equation 
  T = Medium.temperature(Medium.setState_phX(port.p, inStream(port.h_outflow), inStream(port.Xi_outflow)));
end Temperature;

Buildings.Fluid.Sensors.TemperatureDryBulbDynamic Buildings.Fluid.Sensors.TemperatureDryBulbDynamic

Ideal temperature sensor

Buildings.Fluid.Sensors.TemperatureDryBulbDynamic

Information


This component monitors the temperature of the medium in the flow between fluid ports. The sensor does not influence the fluid. Its output is computed using a first order differential equation.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor (Partial component to model sensors that measure flow properties).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Timetau10Time constant [s]
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
Initialization
InitinitTypeModelica.Blocks.Types.Init.N...Type of initialization (InitialState and InitialOutput are identical)
TemperatureT_startMedium.T_defaultInitial or guess value of output (= state) [K]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_small1E-4*m_flow_nominalFor bi-directional flow, temperature is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required) [kg/s]

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)
output RealOutputTTemperature of the passing fluid [K]

Modelica definition

model TemperatureDryBulbDynamic "Ideal temperature sensor"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor;

  parameter Modelica.SIunits.Time tau(min=0) = 10 "Time constant";
  Modelica.Blocks.Interfaces.RealOutput T( final quantity="ThermodynamicTemperature",
                                           final unit="K",
                                           min = 0,
                                           start=T_start) 
    "Temperature of the passing fluid";

  parameter Medium.MassFlowRate m_flow_nominal(min=0) "Nominal mass flow rate";
  parameter Medium.MassFlowRate m_flow_small(min=0) = 1E-4*m_flow_nominal 
    "For bi-directional flow, temperature is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required)";

  parameter Modelica.Blocks.Types.Init initType = Modelica.Blocks.Types.Init.NoInit 
    "Type of initialization (InitialState and InitialOutput are identical)";
  parameter Modelica.SIunits.Temperature T_start=Medium.T_default 
    "Initial or guess value of output (= state)";
// temporarily made public for testing of large model protected
  Medium.Temperature T_a_inflow "Temperature of inflowing fluid at port_a";
  Medium.Temperature T_b_inflow 
    "Temperature of inflowing fluid at port_b or T_a_inflow, if uni-directional flow";
  Medium.Temperature TMed "Medium to which sensor is exposed";
initial equation 
  if initType == Modelica.Blocks.Types.Init.SteadyState then
    der(T) = 0;
  elseif initType == Modelica.Blocks.Types.Init.InitialState or 
         initType == Modelica.Blocks.Types.Init.InitialOutput then
    T = T_start;
  end if;
equation 
  if allowFlowReversal then
     T_a_inflow = Medium.temperature(Medium.setState_phX(port_b.p, port_b.h_outflow, port_b.Xi_outflow));
     T_b_inflow = Medium.temperature(Medium.setState_phX(port_a.p, port_a.h_outflow, port_a.Xi_outflow));
     TMed = Modelica.Fluid.Utilities.regStep(port_a.m_flow, T_a_inflow, T_b_inflow, m_flow_small);
  else
     TMed = Medium.temperature(Medium.setState_phX(port_b.p, port_b.h_outflow, port_b.Xi_outflow));
     T_a_inflow = TMed;
     T_b_inflow = TMed;
  end if;
  der(T)  = (TMed-T)/tau;
end TemperatureDryBulbDynamic;

Buildings.Fluid.Sensors.TemperatureTwoPort Buildings.Fluid.Sensors.TemperatureTwoPort

Ideal two port temperature sensor

Buildings.Fluid.Sensors.TemperatureTwoPort

Information


This component monitors the temperature of the passing fluid. The sensor is ideal, i.e. it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor (Partial component to model sensors that measure flow properties).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_smallsystem.m_flow_smallFor bi-directional flow, temperature is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required) [kg/s]

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)
output RealOutputTTemperature of the passing fluid [K]

Modelica definition

model TemperatureTwoPort "Ideal two port temperature sensor"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor;

  Modelica.Blocks.Interfaces.RealOutput T( final quantity="ThermodynamicTemperature",
                                           final unit="K",
                                           min = 0) 
    "Temperature of the passing fluid";
  parameter Medium.MassFlowRate m_flow_small(min=0) = system.m_flow_small 
    "For bi-directional flow, temperature is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required)";

protected 
  Medium.Temperature T_a_inflow "Temperature of inflowing fluid at port_a";
  Medium.Temperature T_b_inflow 
    "Temperature of inflowing fluid at port_b or T_a_inflow, if uni-directional flow";
equation 
  if allowFlowReversal then
     T_a_inflow = Medium.temperature(Medium.setState_phX(port_b.p, port_b.h_outflow, port_b.Xi_outflow));
     T_b_inflow = Medium.temperature(Medium.setState_phX(port_a.p, port_a.h_outflow, port_a.Xi_outflow));
     T = Modelica.Fluid.Utilities.regStep(port_a.m_flow, T_a_inflow, T_b_inflow, m_flow_small);
  else
     T = Medium.temperature(Medium.setState_phX(port_b.p, port_b.h_outflow, port_b.Xi_outflow));
     T_a_inflow = T;
     T_b_inflow = T;
  end if;
end TemperatureTwoPort;

Buildings.Fluid.Sensors.TemperatureWetBulb Buildings.Fluid.Sensors.TemperatureWetBulb

Ideal wet bulb temperature sensor

Buildings.Fluid.Sensors.TemperatureWetBulb

Information


This component monitors the wet bulb temperature of the medium in the flow between fluid ports. The sensor is ideal, i.e., it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor (Partial component to model sensors that measure flow properties).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_small1e-4For bi-directional flow, specific enthalpy is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required) [kg/s]

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)
output RealOutputTWet bulb temperature in port medium [K]

Modelica definition

model TemperatureWetBulb "Ideal wet bulb temperature sensor"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor;

  Modelica.Blocks.Interfaces.RealOutput T(
    start=Medium.T_default,
    final quantity="Temperature",
    final unit="K") "Wet bulb temperature in port medium";
  parameter Medium.MassFlowRate m_flow_small(min=0) = 1e-4 
    "For bi-directional flow, specific enthalpy is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required)";
protected 
  Buildings.Utilities.Psychrometrics.Twb_TdbXi wetBulMod(redeclare package
      Medium = Medium) "Block for wet bulb temperature";
  Modelica.SIunits.SpecificEnthalpy h "Specific enthalpy";
  Medium.MassFraction Xi[Medium.nXi] 
    "Species vector, needed because indexed argument for the operator inStream is not supported";
  //Medium.MassFraction X[Medium.nX] "Species vector";
equation 

  if allowFlowReversal then
    h  = Modelica.Fluid.Utilities.regStep(port_a.m_flow, port_b.h_outflow, port_a.h_outflow, m_flow_small);
    Xi = Modelica.Fluid.Utilities.regStep(port_a.m_flow, port_b.Xi_outflow, port_a.Xi_outflow, m_flow_small);
  else
    h = port_b.h_outflow;
    Xi = port_b.Xi_outflow;
  end if;
  // Compute wet bulb temperature
  wetBulMod.Tdb = Medium.T_phX(port_a.p, h, cat(1,Xi,{1-sum(Xi)}));
  wetBulMod.Xi = Xi;
  wetBulMod.p  = port_a.p;
  T = wetBulMod.Twb;
end TemperatureWetBulb;

Buildings.Fluid.Sensors.TraceSubstances Buildings.Fluid.Sensors.TraceSubstances

Ideal one port trace substances sensor

Buildings.Fluid.Sensors.TraceSubstances

Information


This component monitors the trace substances contained in the fluid passing its port. The sensor is ideal, i.e. it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor (Partial component to model a sensor that measures a potential variable), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the sensor
StringsubstanceName"CO2"Name of trace substance

Connectors

TypeNameDescription
FluidPort_aport 
output RealOutputCTrace substance in port medium

Modelica definition

model TraceSubstances "Ideal one port trace substances sensor"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialAbsoluteSensor;
  extends Modelica.Icons.RotationalSensor;
  parameter String substanceName = "CO2" "Name of trace substance";

  Modelica.Blocks.Interfaces.RealOutput C "Trace substance in port medium";

protected 
  parameter Integer ind(fixed=false) 
    "Index of species in vector of auxiliary substances";
  parameter Real s[Medium.nC](fixed=false) 
    "Vector with zero everywhere except where species is";
  Medium.ExtraProperty CVec[Medium.nC](
      quantity=Medium.extraPropertiesNames) 
    "Trace substances vector, needed because indexed argument for the operator inStream is not supported";
initial algorithm 
  ind:= -1;
  for i in 1:Medium.nC loop
    if ( Modelica.Utilities.Strings.isEqual(Medium.extraPropertiesNames[i], substanceName)) then
      ind := i;
      s[i] :=1;
    else
      s[i] :=0;
    end if;
  end for;
  assert(ind > 0, "Trace substance '" + substanceName + "' is not present in medium '"
         + Medium.mediumName + "'.\n"
         + "Check sensor parameter and medium model.");
equation 
  CVec = inStream(port.C_outflow);
  // We obtain the species concentration with a vector multiplication
  // because Dymola 7.3 cannot find the derivative in the model
  // Buildings.Examples.VAVSystemCTControl.mo
  // if we set C = CVec[ind];
  C = s*CVec;
end TraceSubstances;

Buildings.Fluid.Sensors.TraceSubstancesTwoPort Buildings.Fluid.Sensors.TraceSubstancesTwoPort

Ideal two port sensor for trace substance

Buildings.Fluid.Sensors.TraceSubstancesTwoPort

Information


This component monitors the trace substance of the passing fluid. The sensor is ideal, i.e. it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor (Partial component to model sensors that measure flow properties), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
StringsubstanceName"CO2"Name of trace substance
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_smallsystem.m_flow_smallFor bi-directional flow, trace substance is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required) [kg/s]

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)
output RealOutputCTrace substance of the passing fluid

Modelica definition

model TraceSubstancesTwoPort 
  "Ideal two port sensor for trace substance"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor;
  extends Modelica.Icons.RotationalSensor;
  Modelica.Blocks.Interfaces.RealOutput C 
    "Trace substance of the passing fluid";
  parameter String substanceName = "CO2" "Name of trace substance";
  parameter Medium.MassFlowRate m_flow_small(min=0) = system.m_flow_small 
    "For bi-directional flow, trace substance is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required)";

protected 
  parameter Integer ind(fixed=false) 
    "Index of species in vector of auxiliary substances";
initial algorithm 
  ind:= -1;
  for i in 1:Medium.nC loop
    if ( Modelica.Utilities.Strings.isEqual(Medium.extraPropertiesNames[i], substanceName)) then
      ind := i;
    end if;
  end for;
  assert(ind > 0, "Trace substance '" + substanceName + "' is not present in medium '"
         + Medium.mediumName + "'.\n"
         + "Check sensor parameter and medium model.");
equation 
  if allowFlowReversal then
     C = Modelica.Fluid.Utilities.regStep(port_a.m_flow, port_b.C_outflow[ind], port_a.C_outflow[ind], m_flow_small);
  else
     C = port_b.C_outflow[ind];
  end if;
end TraceSubstancesTwoPort;

Buildings.Fluid.Sensors.VolumeFlowRate Buildings.Fluid.Sensors.VolumeFlowRate

Ideal sensor for volume flow rate

Buildings.Fluid.Sensors.VolumeFlowRate

Information


This component monitors the volume flow rate flowing from port_a to port_b. The sensor is ideal, i.e. it does not influence the fluid.

Extends from Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor (Partial component to model sensors that measure flow properties), Modelica.Icons.RotationalSensor (Icon representing rotational measurement device).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_smallsystem.m_flow_smallFor bi-directional flow, density is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required) [kg/s]

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)
output RealOutputV_flowVolume flow rate from port_a to port_b [m3/s]

Modelica definition

model VolumeFlowRate "Ideal sensor for volume flow rate"
  extends Modelica.Fluid.Sensors.BaseClasses.PartialFlowSensor;
  extends Modelica.Icons.RotationalSensor;
  Modelica.Blocks.Interfaces.RealOutput V_flow(final quantity="VolumeFlowRate",
                                               final unit="m3/s") 
    "Volume flow rate from port_a to port_b";
  parameter Medium.MassFlowRate m_flow_small(min=0) = system.m_flow_small 
    "For bi-directional flow, density is regularized in the region |m_flow| < m_flow_small (m_flow_small > 0 required)";

protected 
  Medium.Density rho_a_inflow "Density of inflowing fluid at port_a";
  Medium.Density rho_b_inflow 
    "Density of inflowing fluid at port_b or rho_a_inflow, if uni-directional flow";
  Medium.Density d "Density of the passing fluid";
equation 
  if allowFlowReversal then
     rho_a_inflow = Medium.density(Medium.setState_phX(port_b.p, port_b.h_outflow, port_b.Xi_outflow));
     rho_b_inflow = Medium.density(Medium.setState_phX(port_a.p, port_a.h_outflow, port_a.Xi_outflow));
     d = Modelica.Fluid.Utilities.regStep(port_a.m_flow, rho_a_inflow, rho_b_inflow, m_flow_small);
  else
     d = Medium.density(Medium.setState_phX(port_b.p, port_b.h_outflow, port_b.Xi_outflow));
     rho_a_inflow = d;
     rho_b_inflow = d;
  end if;
  V_flow = port_a.m_flow/d;
end VolumeFlowRate;

HTML-documentation generated by Dymola Mon Jun 14 14:28:42 2010.