Buildings.Fluids.HeatExchangers.BaseClasses

Package with base classes for heat exchanger models

Information


This package contains base classes that are used to construct the models in 
Buildings.Fluids.HeatExchangers.

Extends from Modelica_Fluid.Icons.BaseClassLibrary (Icon for library).

Package Content

NameDescription
Buildings.Fluids.HeatExchangers.BaseClasses.CoilHeader CoilHeader Header for a heat exchanger register
Buildings.Fluids.HeatExchangers.BaseClasses.CoilRegister CoilRegister Register for a heat exchanger
Buildings.Fluids.HeatExchangers.BaseClasses.DuctManifoldFixedResistance DuctManifoldFixedResistance Manifold for a heat exchanger air duct connection
Buildings.Fluids.HeatExchangers.BaseClasses.DuctManifoldNoResistance DuctManifoldNoResistance Duct manifold without resistance
Buildings.Fluids.HeatExchangers.BaseClasses.Examples Examples Collection of models that illustrate model use and test models
Buildings.Fluids.HeatExchangers.BaseClasses.HADryCoil HADryCoil Sensible convective heat transfer model for air to water coil
Buildings.Fluids.HeatExchangers.BaseClasses.HexElement HexElement Element of a heat exchanger
Buildings.Fluids.HeatExchangers.BaseClasses.MassExchange MassExchange Block to compute the latent heat transfer based on the Lewis number
Buildings.Fluids.HeatExchangers.BaseClasses.PartialDuctManifold PartialDuctManifold Partial manifold for heat exchanger duct connection
Buildings.Fluids.HeatExchangers.BaseClasses.PartialDuctPipeManifold PartialDuctPipeManifold Partial heat exchanger duct and pipe manifold
Buildings.Fluids.HeatExchangers.BaseClasses.PartialPipeManifold PartialPipeManifold Partial pipe manifold for a heat exchanger
Buildings.Fluids.HeatExchangers.BaseClasses.PipeManifoldFixedResistance PipeManifoldFixedResistance Pipe manifold for a heat exchanger connection
Buildings.Fluids.HeatExchangers.BaseClasses.PipeManifoldNoResistance PipeManifoldNoResistance Manifold for heat exchanger register


Buildings.Fluids.HeatExchangers.BaseClasses.CoilHeader Buildings.Fluids.HeatExchangers.BaseClasses.CoilHeader

Header for a heat exchanger register

Buildings.Fluids.HeatExchangers.BaseClasses.CoilHeader

Information


Header for a heat exchanger coil.

This model connects the flow between its ports without modeling flow friction. Currently, the ports are connected without redistributing the flow. In latter versions, the model may be changed to define different flow reroutings in the coil header.


Extends from Buildings.BaseClasses.BaseIcon (Base icon).

Parameters

TypeNameDefaultDescription
IntegernPipPar Number of parallel pipes in each register
Initialization
MassFlowRatemStart_flow_a Guess value for mass flow rate at port_a [kg/s]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)

Connectors

TypeNameDescription
FluidPort_aport_a[nPipPar]Fluid connector a for medium (positive design flow direction is from port_a to port_b)
FluidPort_bport_b[nPipPar]Fluid connector b for medium (positive design flow direction is from port_a to port_b)

Modelica definition

model CoilHeader "Header for a heat exchanger register"
  extends Buildings.BaseClasses.BaseIcon;

  outer Modelica_Fluid.System system "System wide properties";

  replaceable package Medium = 
      Modelica.Media.Interfaces.PartialMedium "Medium in the component";
  parameter Boolean allowFlowReversal = system.allowFlowReversal 
    "= true to allow flow reversal, false restricts to design direction (port_a -> port_b)";

  parameter Integer nPipPar(min=1) "Number of parallel pipes in each register";
  parameter Modelica.SIunits.MassFlowRate mStart_flow_a 
    "Guess value for mass flow rate at port_a";

  Modelica_Fluid.Interfaces.FluidPort_a port_a[nPipPar](
        redeclare each final package Medium = Medium,
        each m_flow(start=mStart_flow_a/nPipPar, min=if allowFlowReversal then -Modelica.Constants.inf else 0)) 
    "Fluid connector a for medium (positive design flow direction is from port_a to port_b)";

  Modelica_Fluid.Interfaces.FluidPort_b port_b[nPipPar](
        redeclare each final package Medium = Medium,
        each m_flow(start=-mStart_flow_a/nPipPar, max=if allowFlowReversal then +Modelica.Constants.inf else 0)) 
    "Fluid connector b for medium (positive design flow direction is from port_a to port_b)";

equation 
  connect(port_a, port_b);
end CoilHeader;

Buildings.Fluids.HeatExchangers.BaseClasses.CoilRegister Buildings.Fluids.HeatExchangers.BaseClasses.CoilRegister

Register for a heat exchanger

Buildings.Fluids.HeatExchangers.BaseClasses.CoilRegister

Information


Register of a heat exchanger with dynamics on the fluids and the solid. The register represents one array of pipes that are perpendicular to the air stream. The hA value for both fluids is an input. The driving force for the heat transfer is the temperature difference between the fluid volumes and the solid in each heat exchanger element.


Extends from Buildings.Fluids.Interfaces.FourPortFlowResistanceParameters (Parameters for flow resistance for models with four ports).

Parameters

TypeNameDefaultDescription
IntegernPipPar2Number of parallel pipes in each register
IntegernPipSeg3Number of pipe segments per register used for discretization
BooleanallowCondensationtrueSet to false to compute sensible heat transfer only
Nominal condition
Pressuredp1_nominal Pressure [Pa]
Pressuredp2_nominal Pressure [Pa]
ThermalConductanceUA_nominal Thermal conductance at nominal flow, used to compute time constant [W/K]
MassFlowRatem1_flow_nominal Mass flow rate medim 1 [kg/s]
MassFlowRatem2_flow_nominal Mass flow rate medium 2 [kg/s]
Timetau120Time constant at nominal flow for medium 1 [s]
Timetau21Time constant at nominal flow for medium 2 [s]
Timetau_m60Time constant of metal at nominal UA value [s]
Fluid 1
BooleansteadyState_1falseSet to true for steady state model for fluid 1
Fluid 2
BooleansteadyState_2falseSet to true for steady state model for fluid 2
Flow resistance
Medium 1
BooleancomputeFlowResistance1true=true, compute flow resistance. Set to false to assume no friction
Booleanfrom_dp1true= true, use m_flow = f(dp) else dp = f(m_flow)
BooleanlinearizeFlowResistance1false= true, use linear relation between m_flow and dp for any flow rate
RealdeltaM10.1Fraction of nominal flow rate where flow transitions to laminar
Medium 2
BooleancomputeFlowResistance2true=true, compute flow resistance. Set to false to assume no friction
Booleanfrom_dp2true= true, use m_flow = f(dp) else dp = f(m_flow)
BooleanlinearizeFlowResistance2false= true, use linear relation between m_flow and dp for any flow rate
RealdeltaM20.1Fraction of nominal flow rate where flow transitions to laminar
Assumptions
BooleanallowFlowReversal1system.allowFlowReversal= true to allow flow reversal in medium 1, false restricts to design direction (port_a -> port_b)
BooleanallowFlowReversal2system.allowFlowReversal= true to allow flow reversal in medium 2, false restricts to design direction (port_a -> port_b)
Dynamics
DynamicsenergyDynamics1Modelica_Fluid.Types.Dynamic...Default formulation of energy balances for volume 1
DynamicsenergyDynamics2Modelica_Fluid.Types.Dynamic...Default formulation of energy balances for volume 2

Connectors

TypeNameDescription
FluidPort_aport_a1[nPipPar]Fluid connector a for medium 1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_bport_b1[nPipPar]Fluid connector b for medium 1 (positive design flow direction is from port_a to port_b)
FluidPort_aport_a2[nPipPar, nPipSeg]Fluid connector a for medium 2 (positive design flow direction is from port_a2 to port_b2)
FluidPort_bport_b2[nPipPar, nPipSeg]Fluid connector b for medium 2 (positive design flow direction is from port_a to port_b)
input RealInputGc_2Signal representing the convective thermal conductance medium 2 in [W/K]
input RealInputGc_1Signal representing the convective thermal conductance medium 1 in [W/K]

Modelica definition

model CoilRegister "Register for a heat exchanger"
  import Modelica.Constants;
  extends Buildings.Fluids.Interfaces.FourPortFlowResistanceParameters(
     final computeFlowResistance1=true, final computeFlowResistance2=true);
  replaceable package Medium1 = 
      Modelica.Media.Interfaces.PartialMedium "Medium 1 in the component";
  replaceable package Medium2 = 
      Modelica.Media.Interfaces.PartialMedium "Medium 2 in the component";
  outer Modelica_Fluid.System system "System wide properties";

  parameter Boolean allowFlowReversal1 = system.allowFlowReversal 
    "= true to allow flow reversal in medium 1, false restricts to design direction (port_a -> port_b)";
  parameter Boolean allowFlowReversal2 = system.allowFlowReversal 
    "= true to allow flow reversal in medium 2, false restricts to design direction (port_a -> port_b)";

  parameter Integer nPipPar(min=1)=2 
    "Number of parallel pipes in each register";
  parameter Integer nPipSeg(min=1)=3 
    "Number of pipe segments per register used for discretization";
  final parameter Integer nEle = nPipPar * nPipSeg 
    "Number of heat exchanger elements";

  Buildings.Fluids.HeatExchangers.BaseClasses.HexElement[
                      nPipPar, nPipSeg] ele(
    redeclare each package Medium1 = Medium1,
    redeclare each package Medium2 = Medium2,
    each allowFlowReversal1=allowFlowReversal1,
    each allowFlowReversal2=allowFlowReversal2,
    each tau1=tau1/nPipSeg,
    each m1_flow_nominal=m1_flow_nominal/nPipPar,
    each tau2=tau2,
    each m2_flow_nominal=m2_flow_nominal/nPipPar/nPipSeg,
    each tau_m=tau_m,
    each UA_nominal=UA_nominal/nPipPar/nPipSeg,
    each energyDynamics1=energyDynamics1,
    each energyDynamics2=energyDynamics2,
    each allowCondensation=allowCondensation,
    each from_dp1=from_dp1,
    each linearizeFlowResistance1=linearizeFlowResistance1,
    each deltaM1=deltaM1,
    each from_dp2=from_dp2,
    each linearizeFlowResistance2=linearizeFlowResistance2,
    each deltaM2=deltaM2,
    each dp1_nominal=dp1_nominal,
    each dp2_nominal=dp2_nominal) "Element of a heat exchanger";

  Modelica_Fluid.Interfaces.FluidPort_a[nPipPar] port_a1(
        redeclare each package Medium = Medium1,
        each m_flow(start=0, min=if allowFlowReversal1 then -Constants.inf else 0)) 
    "Fluid connector a for medium 1 (positive design flow direction is from port_a1 to port_b1)";
  Modelica_Fluid.Interfaces.FluidPort_b[nPipPar] port_b1(
        redeclare each package Medium = Medium1,
        each m_flow(start=0, max=if allowFlowReversal1 then +Constants.inf else 0)) 
    "Fluid connector b for medium 1 (positive design flow direction is from port_a to port_b)";
  Modelica_Fluid.Interfaces.FluidPort_a[nPipPar,nPipSeg] port_a2(
        redeclare each package Medium = Medium2,
        each m_flow(start=0, min=if allowFlowReversal2 then -Constants.inf else 0)) 
    "Fluid connector a for medium 2 (positive design flow direction is from port_a2 to port_b2)";
  Modelica_Fluid.Interfaces.FluidPort_b[nPipPar,nPipSeg] port_b2(
        redeclare each package Medium = Medium2,
        each m_flow(start=0, max=if allowFlowReversal2 then +Constants.inf else 0)) 
    "Fluid connector b for medium 2 (positive design flow direction is from port_a to port_b)";

  parameter Modelica.SIunits.ThermalConductance UA_nominal 
    "Thermal conductance at nominal flow, used to compute time constant";

  parameter Modelica.SIunits.MassFlowRate m1_flow_nominal 
    "Mass flow rate medim 1";
  parameter Modelica.SIunits.MassFlowRate m2_flow_nominal 
    "Mass flow rate medium 2";

  parameter Modelica.SIunits.Time tau1=20 
    "Time constant at nominal flow for medium 1";
  parameter Modelica.SIunits.Time tau2=1 
    "Time constant at nominal flow for medium 2";
  parameter Boolean steadyState_1=false 
    "Set to true for steady state model for fluid 1";
  parameter Boolean steadyState_2=false 
    "Set to true for steady state model for fluid 2";
  Modelica.SIunits.HeatFlowRate Q1_flow 
    "Heat transfered from solid into medium 1";
  Modelica.SIunits.HeatFlowRate Q2_flow 
    "Heat transfered from solid into medium 2";
  parameter Modelica.SIunits.Time tau_m=60 
    "Time constant of metal at nominal UA value";
  parameter Boolean allowCondensation = true 
    "Set to false to compute sensible heat transfer only";
  parameter Modelica_Fluid.Types.Dynamics energyDynamics1=
    Modelica_Fluid.Types.Dynamics.DynamicFreeInitial 
    "Default formulation of energy balances for volume 1";
  parameter Modelica_Fluid.Types.Dynamics energyDynamics2=
    Modelica_Fluid.Types.Dynamics.DynamicFreeInitial 
    "Default formulation of energy balances for volume 2";
  Modelica.Blocks.Interfaces.RealInput Gc_2 
    "Signal representing the convective thermal conductance medium 2 in [W/K]";
  Modelica.Blocks.Interfaces.RealInput Gc_1 
    "Signal representing the convective thermal conductance medium 1 in [W/K]";
protected 
  Modelica.Blocks.Math.Gain gai_1(k=1/nEle) 
    "Gain medium-side 1 to take discretization into account";
  Modelica.Blocks.Math.Gain gai_2(k=1/nEle) 
    "Gain medium-side 2 to take discretization into account";
equation 
  Q1_flow = sum(ele[i,j].Q1_flow for i in 1:nPipPar, j in 1:nPipSeg);
  Q2_flow = sum(ele[i,j].Q2_flow for i in 1:nPipPar, j in 1:nPipSeg);
  for i in 1:nPipPar loop
    // liquid side (pipes)
    connect(ele[i,1].port_a1,       port_a1[i]);
    connect(ele[i,nPipSeg].port_b1, port_b1[i]);
    for j in 1:nPipSeg-1 loop
      connect(ele[i,j].port_b1, ele[i,j+1].port_a1);
    end for;
    // gas side (duct)                                                                                      //water connections
    for j in 1:nPipSeg loop
      connect(ele[i,j].port_a2, port_a2[i,j]);
      connect(ele[i,j].port_b2, port_b2[i,j]);
    end for;
  end for;

  connect(Gc_1, gai_1.u);
  connect(Gc_2, gai_2.u);
  for i in 1:nPipPar loop

     for j in 1:nPipSeg loop
      connect(gai_1.y, ele[i,j].Gc_1);
      connect(gai_2.y, ele[i,j].Gc_2);
     end for;
  end for;

end CoilRegister;

Buildings.Fluids.HeatExchangers.BaseClasses.DuctManifoldFixedResistance Buildings.Fluids.HeatExchangers.BaseClasses.DuctManifoldFixedResistance

Manifold for a heat exchanger air duct connection

Buildings.Fluids.HeatExchangers.BaseClasses.DuctManifoldFixedResistance

Information


Duct manifold with a fixed flow resistance.

This model causes the flow to be distributed equally into each flow path by using a fixed flow resistance for each flow path.


Extends from PartialDuctManifold (Partial manifold for heat exchanger duct connection).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
IntegernPipPar Number of parallel pipes in each register
IntegernPipSeg Number of pipe segments per register used for discretization
Booleanuse_dhfalseSet to true to specify hydraulic diameter
Lengthdh1Hydraulic diameter of duct [m]
RealReC4000Reynolds number where transition to turbulent starts
RealdeltaM0.3Fraction of nominal mass flow rate where transition to turbulent occurs
Lengthdl0.3Length of mixing volume [m]
Initialization
MassFlowRatemStart_flow_a Guess value for mass flow rate at port_a [kg/s]
Nominal Condition
MassFlowRatem_flow_nominal Mass flow rate at port_a [kg/s]
Pressuredp_nominal Pressure [Pa]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Dynamics
DynamicsenergyDynamicsModelica_Fluid.Types.Dynamic...Default formulation of energy balances for volume
Advanced
Booleanlinearizedfalse= true, use linear relation between m_flow and dp for any flow rate
Booleanfrom_dpfalse= true, use m_flow = f(dp) else dp = f(m_flow)

Connectors

TypeNameDescription
FluidPort_aport_aFluid connector a for medium (positive design flow direction is from port_a to port_b)
FluidPort_bport_b[nPipPar, nPipSeg]Fluid connector b for medium (positive design flow direction is from port_a to port_b)

Modelica definition

model DuctManifoldFixedResistance 
  "Manifold for a heat exchanger air duct connection"
  extends PartialDuctManifold;


  parameter Boolean use_dh = false "Set to true to specify hydraulic diameter";

  parameter Modelica.SIunits.MassFlowRate m_flow_nominal 
    "Mass flow rate at port_a";
  parameter Modelica.SIunits.Pressure dp_nominal(min=0) "Pressure";
  parameter Modelica.SIunits.Length dh=1 "Hydraulic diameter of duct";
  parameter Real ReC=4000 
    "Reynolds number where transition to turbulent starts";
  parameter Boolean linearized = false 
    "= true, use linear relation between m_flow and dp for any flow rate";
  parameter Real deltaM(min=0) = 0.3 
    "Fraction of nominal mass flow rate where transition to turbulent occurs";
  parameter Boolean from_dp = false 
    "= true, use m_flow = f(dp) else dp = f(m_flow)";

  Fluids.FixedResistances.FixedResistanceDpM[nPipPar,nPipSeg] fixRes(
    redeclare each package Medium = Medium,
    each m_flow_nominal=m_flow_nominal/nPipPar/nPipSeg,
    each m_flow(start=mStart_flow_a/nPipPar/nPipSeg),
    each dp_nominal=dp_nominal,
    each dh=dh/sqrt(nPipPar*nPipSeg),
    each from_dp=from_dp,
    each deltaM=deltaM,
    each ReC=ReC,
    each use_dh=use_dh,
    each linearized=linearized) "Fixed resistance for each duct";
  parameter Modelica.SIunits.Length dl = 0.3 "Length of mixing volume";
  Fluids.MixingVolumes.MixingVolume vol(redeclare package Medium = Medium,
    final V=dh*dh*dl,
    final nPorts=1+nPipPar*nPipSeg,
    final energyDynamics=energyDynamics,
    final massDynamics=energyDynamics);
  parameter Modelica_Fluid.Types.Dynamics energyDynamics=
    Modelica_Fluid.Types.Dynamics.DynamicFreeInitial 
    "Default formulation of energy balances for volume";
equation 
  for i in 1:nPipPar loop
    for j in 1:nPipSeg loop
      connect(vol.ports[1+(i-1)*nPipSeg+j], fixRes[i, j].port_a);
    end for;
  end for;

  connect(port_a, vol.ports[1]);
  connect(fixRes.port_b, port_b);
end DuctManifoldFixedResistance;

Buildings.Fluids.HeatExchangers.BaseClasses.DuctManifoldNoResistance Buildings.Fluids.HeatExchangers.BaseClasses.DuctManifoldNoResistance

Duct manifold without resistance

Buildings.Fluids.HeatExchangers.BaseClasses.DuctManifoldNoResistance

Information


Duct manifold without flow resistance.

This model connects the flows between the ports without modeling flow friction. The model is used in conjunction with a manifold which contains pressure drop elements and that is added to the other side of the heat exchanger registers.


Extends from PartialDuctManifold (Partial manifold for heat exchanger duct connection).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
IntegernPipPar Number of parallel pipes in each register
IntegernPipSeg Number of pipe segments per register used for discretization
Initialization
MassFlowRatemStart_flow_a Guess value for mass flow rate at port_a [kg/s]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)

Connectors

TypeNameDescription
FluidPort_aport_aFluid connector a for medium (positive design flow direction is from port_a to port_b)
FluidPort_bport_b[nPipPar, nPipSeg]Fluid connector b for medium (positive design flow direction is from port_a to port_b)

Modelica definition

model DuctManifoldNoResistance "Duct manifold without resistance"
  extends PartialDuctManifold;

equation 
  for i in 1:nPipPar loop
    for j in 1:nPipSeg loop
    connect(port_a, port_b[i, j]);
    end for;
  end for;
end DuctManifoldNoResistance;

Buildings.Fluids.HeatExchangers.BaseClasses.HADryCoil Buildings.Fluids.HeatExchangers.BaseClasses.HADryCoil

Sensible convective heat transfer model for air to water coil

Buildings.Fluids.HeatExchangers.BaseClasses.HADryCoil

Information


Model for sensible convective heat transfer coefficients for an air to water coil.

This model computes the convective heat transfer coefficient for an air to water coil. The parameters allow a user to enable or disable, individually for each medium, the mass flow and/or the temperature dependence of the convective heat transfer coefficients. For a detailed explanation of the equation, see the references below.

References


Extends from Buildings.BaseClasses.BaseIcon (Base icon).

Parameters

TypeNameDefaultDescription
Realr0.5Ratio between air-side and water-side convective heat transfer coefficient
Realn_w0.85Water-side exponent for convective heat transfer coefficient, h~m_flow^n
Realn_a0.8Air-side exponent for convective heat transfer coefficient, h~m_flow^n
Nominal condition
ThermalConductanceUA_nominal Thermal conductance at nominal flow [W/K]
MassFlowRatem_flow_nominal_w Water mass flow rate [kg/s]
MassFlowRatem_flow_nominal_a Air mass flow rate [kg/s]
ThermalConductancehA_nominal_wUA_nominal*(r + 1)/rWater side convective heat transfer coefficient [W/K]
ThermalConductancehA_nominal_ar*hA_nominal_wAir side convective heat transfer coefficient, including fin resistance [W/K]
TemperatureT0_wModelica.SIunits.Conversions...Water temperature [K]
TemperatureT0_aModelica.SIunits.Conversions...Air temperature [K]
Advanced
Modeling detail
BooleanwaterSideFlowDependenttrueSet to false to make water-side hA independent of mass flow rate
BooleanairSideFlowDependenttrueSet to false to make air-side hA independent of mass flow rate
BooleanwaterSideTemperatureDependenttrueSet to false to make water-side hA independent of temperature
BooleanairSideTemperatureDependenttrueSet to false to make air-side hA independent of temperature

Connectors

TypeNameDescription
input RealInputm1_flowMass flow rate medium 1
input RealInputm2_flowMass flow rate medium 2
input RealInputT_1Temperature medium 1
input RealInputT_2Temperature medium 2
output RealOutputhA_1Convective heat transfer medium 1
output RealOutputhA_2Convective heat transfer medium 2

Modelica definition

model HADryCoil 
  "Sensible convective heat transfer model for air to water coil"
  extends Buildings.BaseClasses.BaseIcon;
  parameter Modelica.SIunits.ThermalConductance UA_nominal(min=0) 
    "Thermal conductance at nominal flow";

  parameter Modelica.SIunits.MassFlowRate m_flow_nominal_w 
    "Water mass flow rate";
  parameter Modelica.SIunits.MassFlowRate m_flow_nominal_a "Air mass flow rate";

  Modelica.Blocks.Interfaces.RealInput m1_flow "Mass flow rate medium 1";
  Modelica.Blocks.Interfaces.RealInput m2_flow "Mass flow rate medium 2";
  Modelica.Blocks.Interfaces.RealInput T_1 "Temperature medium 1";
  Modelica.Blocks.Interfaces.RealInput T_2 "Temperature medium 2";

  Modelica.Blocks.Interfaces.RealOutput hA_1 
    "Convective heat transfer medium 1";
  Modelica.Blocks.Interfaces.RealOutput hA_2 
    "Convective heat transfer medium 2";


  parameter Real r(min=0, max=1)=0.5 
    "Ratio between air-side and water-side convective heat transfer coefficient";
  parameter Modelica.SIunits.ThermalConductance hA_nominal_w(min=0)=UA_nominal * (r+1)/r 
    "Water side convective heat transfer coefficient";
  parameter Modelica.SIunits.ThermalConductance hA_nominal_a(min=0)=r * hA_nominal_w 
    "Air side convective heat transfer coefficient, including fin resistance";
  parameter Real n_w(min=0, max=1)=0.85 
    "Water-side exponent for convective heat transfer coefficient, h~m_flow^n";
  parameter Real n_a(min=0, max=1)=0.8 
    "Air-side exponent for convective heat transfer coefficient, h~m_flow^n";
  parameter Modelica.SIunits.Temperature T0_w=
          Modelica.SIunits.Conversions.from_degC(20) "Water temperature";
  parameter Modelica.SIunits.Temperature T0_a=
          Modelica.SIunits.Conversions.from_degC(20) "Air temperature";
  parameter Boolean waterSideFlowDependent = true 
    "Set to false to make water-side hA independent of mass flow rate";
  parameter Boolean airSideFlowDependent = true 
    "Set to false to make air-side hA independent of mass flow rate";
  parameter Boolean waterSideTemperatureDependent = true 
    "Set to false to make water-side hA independent of temperature";
  parameter Boolean airSideTemperatureDependent = true 
    "Set to false to make air-side hA independent of temperature";
protected 
  Real x_a(min=0) 
    "Factor for air side temperature dependent variation of heat transfer coefficient";
  Real x_w(min=0) 
    "Factor for water side temperature dependent variation of heat transfer coefficient";
  Real s_w(min=0, nominal=0.01) 
    "Coefficient for temperature dependence of water side heat transfer coefficient";
  Real fm_w "Fraction of actual to nominal mass flow rate";
  Real fm_a "Fraction of actual to nominal mass flow rate";
equation 
  fm_w = if waterSideFlowDependent then 
              m1_flow / m_flow_nominal_w else 1;
  fm_a = if airSideFlowDependent then 
              m2_flow / m_flow_nominal_a else 1;
  s_w =  if waterSideTemperatureDependent then 
            0.014/(1+0.014*Modelica.SIunits.Conversions.to_degC(T_1)) else 
              1;
  x_w = if waterSideTemperatureDependent then 
         1 + s_w * (T_1-T0_w) else 
              1;
  x_a = if airSideTemperatureDependent then 
         1 + 4.769E-3 * (T_2-T0_a) else 
              1;
  if ( waterSideFlowDependent == true) then
    hA_1 = x_w * hA_nominal_w
               * Buildings.Utilities.Math.Functions.regNonZeroPower(fm_w, n_w, 0.1);
  else
    hA_1 = x_w * hA_nominal_w;
  end if;

  if ( airSideFlowDependent == true) then
    hA_2 = x_a * hA_nominal_a
               * Buildings.Utilities.Math.Functions.regNonZeroPower(fm_a, n_a, 0.1);
  else
    hA_2 = x_a * hA_nominal_a;
  end if;
end HADryCoil;

Buildings.Fluids.HeatExchangers.BaseClasses.HexElement Buildings.Fluids.HeatExchangers.BaseClasses.HexElement

Element of a heat exchanger

Buildings.Fluids.HeatExchangers.BaseClasses.HexElement

Information


Element of a heat exchanger with dynamics on the fluids and the solid. The hA value for both fluids is an input. The driving force for the heat transfer is the temperature difference between the fluid volumes and the solid.

The heat capacity C of the metal is assigned as follows. Suppose the metal temperature is governed by

     dT
  C ---- = hA_1 (T_1 - T) + hA_2 (T_2 - T)
     dt
where hA are the convective heat transfer coefficients that also take into account heat conduction in the heat exchanger fins and T_1 and T_2 are the medium temperatures. Assuming hA_1=hA_2, this equation can be rewritten as
     dT
  C ---- = 2 UA_nominal ( (T_1 - T) + (T_2 - T) )
     dt
where UA_nominal is the UA value at nominal condition. Hence we set the heat capacity of the metal to C = 2 * UA_nominal * tau_m.


Extends from Buildings.Fluids.Interfaces.PartialDynamicFourPortTransformer (Partial model transporting two fluid streams between four ports with storing mass or energy).

Parameters

TypeNameDefaultDescription
replaceable package Medium1PartialMediumMedium 1 in the component
replaceable package Medium2PartialMediumMedium 2 in the component
HeatCapacityC2*UA_nominal*tau_mHeat capacity of metal (= cp*m) [J/K]
BooleanallowCondensationtrueSet to false to compute sensible heat transfer only
Nominal condition
MassFlowRatem1_flow_nominal Nominal mass flow rate [kg/s]
MassFlowRatem2_flow_nominalm1_flow_nominalNominal mass flow rate [kg/s]
Pressuredp1_nominal Pressure [Pa]
Pressuredp2_nominal Pressure [Pa]
Timetau160Time constant at nominal flow [s]
Timetau260Time constant at nominal flow [s]
ThermalConductanceUA_nominal Thermal conductance at nominal flow, used to compute time constant [W/K]
Timetau_m60Time constant of metal at nominal UA value [s]
Initialization
MassFlowRatem1_flow.start0Mass flow rate from port_a1 to port_b1 (m1_flow > 0 is design flow direction) [kg/s]
Pressuredp1.start0Pressure difference between port_a1 and port_b1 [Pa]
MassFlowRatem2_flow.start0Mass flow rate from port_a2 to port_b2 (m2_flow > 0 is design flow direction) [kg/s]
Pressuredp2.start0Pressure difference between port_a2 and port_b2 [Pa]
Assumptions
BooleanallowFlowReversal1system.allowFlowReversal= true to allow flow reversal in medium 1, false restricts to design direction (port_a -> port_b)
BooleanallowFlowReversal2system.allowFlowReversal= true to allow flow reversal in medium 2, false restricts to design direction (port_a -> port_b)
Dynamics
DynamicsenergyDynamics1Modelica_Fluid.Types.Dynamic...Default formulation of energy balances for volume 1
DynamicsenergyDynamics2Modelica_Fluid.Types.Dynamic...Default formulation of energy balances for volume 2
Advanced
MassFlowRatem1_flow_small1E-4*m1_flow_nominalSmall mass flow rate for regularization of zero flow [kg/s]
MassFlowRatem2_flow_small1E-4*m2_flow_nominalSmall mass flow rate for regularization of zero flow [kg/s]
Diagnostics
Booleanshow_V_flowtrue= true, if volume flow rate at inflowing port is computed
Initialization
AbsolutePressurep_a1_startsystem.p_startGuess value for inlet pressure [Pa]
AbsolutePressurep_b1_startp_a1_startGuess value for outlet pressure [Pa]
AbsolutePressurep_a2_startsystem.p_startGuess value for inlet pressure [Pa]
AbsolutePressurep_b2_startp_a2_startGuess value for outlet pressure [Pa]
Flow resistance
Medium 1
Booleanfrom_dp1true= true, use m_flow = f(dp) else dp = f(m_flow)
BooleanlinearizeFlowResistance1false= true, use linear relation between m_flow and dp for any flow rate
RealdeltaM10.1Fraction of nominal flow rate where flow transitions to laminar
Medium 2
Booleanfrom_dp2true= true, use m_flow = f(dp) else dp = f(m_flow)
BooleanlinearizeFlowResistance2false= true, use linear relation between m_flow and dp for any flow rate
RealdeltaM20.1Fraction of nominal flow rate where flow transitions to laminar

Connectors

TypeNameDescription
FluidPort_aport_a1Fluid connector a1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_bport_b1Fluid connector b1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_aport_a2Fluid connector a2 (positive design flow direction is from port_a2 to port_b2)
FluidPort_bport_b2Fluid connector b2 (positive design flow direction is from port_a2 to port_b2)
input RealInputGc_1Signal representing the convective thermal conductance medium 1 in [W/K]
input RealInputGc_2Signal representing the convective thermal conductance medium 2 in [W/K]

Modelica definition

model HexElement "Element of a heat exchanger"
  extends Buildings.Fluids.Interfaces.PartialDynamicFourPortTransformer(
    vol1(redeclare package Medium = Medium1,
          V=m1_flow_nominal*tau1/rho1_nominal,
          nPorts=2,
          final energyDynamics=energyDynamics1,
          final massDynamics=energyDynamics1),
    vol2(redeclare package Medium = Medium2,
          nPorts = 2,
          V=m2_flow_nominal*tau2/rho2_nominal,
          final energyDynamics=energyDynamics2,
          final massDynamics=energyDynamics2));
  // Note that we MUST declare the value of vol2.V here.
  // Otherwise, if the class of vol2 is redeclared at a higher level,
  // it will overwrite the assignment of V in the base class
  // PartialDynamicFourPortTransformer, which will cause V=0.
  // Assigning the values for vol1 here is optional, but we added
  // it to be consistent in the implementation of vol1 and vol2.

  parameter Modelica.SIunits.ThermalConductance UA_nominal 
    "Thermal conductance at nominal flow, used to compute time constant";
  parameter Modelica.SIunits.Time tau_m(min=0) = 60 
    "Time constant of metal at nominal UA value";
  parameter Modelica.SIunits.HeatCapacity C=2*UA_nominal*tau_m 
    "Heat capacity of metal (= cp*m)";

  Modelica.Blocks.Interfaces.RealInput Gc_1 
    "Signal representing the convective thermal conductance medium 1 in [W/K]";
  Modelica.Blocks.Interfaces.RealInput Gc_2 
    "Signal representing the convective thermal conductance medium 2 in [W/K]";

  parameter Boolean allowCondensation = true 
    "Set to false to compute sensible heat transfer only";

  MassExchange masExc(
       redeclare package Medium = Medium2) if allowCondensation 
    "Model for mass exchange";

  parameter Modelica_Fluid.Types.Dynamics energyDynamics1=
    Modelica_Fluid.Types.Dynamics.DynamicFreeInitial 
    "Default formulation of energy balances for volume 1";
  parameter Modelica_Fluid.Types.Dynamics energyDynamics2=
    Modelica_Fluid.Types.Dynamics.DynamicFreeInitial 
    "Default formulation of energy balances for volume 2";

  Modelica.Thermal.HeatTransfer.Components.HeatCapacitor mas(
                                                  C=C, T(stateSelect=StateSelect.always)) 
    "Mass of metal";
  Modelica.Thermal.HeatTransfer.Components.Convection con1(dT(min=-200)) 
    "Convection (and conduction) on fluid side 1";
  Modelica.Thermal.HeatTransfer.Components.Convection con2(dT(min=-200)) 
    "Convection (and conduction) on fluid side 2";
  Modelica.Thermal.HeatTransfer.Sensors.TemperatureSensor temSen(
    T(final quantity="ThermodynamicTemperature",
      final unit = "K", displayUnit = "degC", min=0)) 
    "Temperature sensor of metal";
  Modelica.Thermal.HeatTransfer.Sensors.HeatFlowSensor heaFloSen_1 
    "Heat input into fluid 1";
  Modelica.Thermal.HeatTransfer.Sensors.HeatFlowSensor heaFloSen_2 
    "Heat input into fluid 1";
equation 
  connect(Gc_1, con1.Gc);
  connect(Gc_2, con2.Gc);
  connect(temSen.T, masExc.TSur);
  connect(Gc_2, masExc.Gc);
  connect(masExc.mWat_flow, vol2.mWat_flow);
  connect(masExc.TLiq, vol2.TWat);
  connect(vol2.XWat, masExc.XInf);
  connect(con1.solid,mas. port);
  connect(con2.solid,mas. port);
  connect(mas.port,temSen. port);
  connect(con1.fluid,heaFloSen_1. port_a);
  connect(con2.fluid,heaFloSen_2. port_a);
  connect(heaFloSen_1.port_b, vol1.heatPort);
  connect(heaFloSen_2.port_b, vol2.heatPort);
end HexElement;

Buildings.Fluids.HeatExchangers.BaseClasses.MassExchange Buildings.Fluids.HeatExchangers.BaseClasses.MassExchange

Block to compute the latent heat transfer based on the Lewis number

Buildings.Fluids.HeatExchangers.BaseClasses.MassExchange

Information


This model computes the mass transfer based on similarity laws between the convective sensible heat transfer coefficient and the mass transfer coefficient.

Using the Lewis number which is defined as the ratio between the heat and mass diffusion coefficients, one can obtain the ratio between convection heat transfer coefficient h in (W/(m^2*K)) and mass transfer coefficient h_m in (m/s) as follows:

  h
 --- = rho * c_p * Le^(1-n),
 h_m
where rho is the mass density, c_p is the specific heat capacity of the bulk medium and n is a coefficient from the boundary layer analysis, which is typically n=1/3. From this equation, we can compute the water vapor mass flow rate n_A in (kg/s) as
  n_A = (Gc) / c_p / Le^(1-n) * (X_s - X_inf),
where Gc is the sensible heat conductivity in (W/K) and X_s and X_inf are the water vapor mass per unit volume in the boundary layer and in the bulk of the medium. In this model, X_s is the saturation water vapor pressure corresponding to the temperature T_sur which is an input.


Extends from Buildings.BaseClasses.BaseIcon (Base icon).

Parameters

TypeNameDefaultDescription
RealLe1Lewis number (around 1 for water vapor in air)
Realn1/3Exponent in bondary layer ratio, delta/delta_t = Pr^n

Connectors

TypeNameDescription
input RealInputXInfWater mass fraction of medium
input RealInputTSurSurface temperature [K]
output RealOutputmWat_flowWater flow rate [kg/s]
output RealOutputTLiqTemperature at which condensate drains from system [K]
input RealInputGcSignal representing the convective (sensible) thermal conductance in [W/K]

Modelica definition

model MassExchange 
  "Block to compute the latent heat transfer based on the Lewis number"
  import Buildings;
  extends Buildings.BaseClasses.BaseIcon;
  replaceable package Medium = Modelica.Media.Interfaces.PartialMedium 
    "Fluid medium model";

  Modelica.Blocks.Interfaces.RealInput XInf "Water mass fraction of medium";
  Modelica.Blocks.Interfaces.RealInput TSur(final quantity="ThermodynamicTemperature",
                                            final unit = "K", displayUnit = "degC", min=0) 
    "Surface temperature";
  Modelica.Blocks.Interfaces.RealOutput mWat_flow(final unit = "kg/s") 
    "Water flow rate";
  Modelica.Blocks.Interfaces.RealOutput TLiq(final quantity="ThermodynamicTemperature",
                                             final unit = "K", displayUnit = "degC", min=0) 
    "Temperature at which condensate drains from system";
  Modelica.Blocks.Interfaces.RealInput Gc 
    "Signal representing the convective (sensible) thermal conductance in [W/K]";
  parameter Real Le = 1 "Lewis number (around 1 for water vapor in air)";
  parameter Real n = 1/3 
    "Exponent in bondary layer ratio, delta/delta_t = Pr^n";
public 
  Buildings.Utilities.Psychrometrics.HumidityRatio_pWat humRatPre(use_p_in=
        false) "Model to convert water vapor pressure into humidity ratio";
  Buildings.Utilities.Psychrometrics.DewPointTemperature_pWat TDewPoi 
    "Model to compute the water vapor pressure at the dew point";
  Modelica.Blocks.Math.Gain gain(k=1/cpLe) 
    "Constant to convert from heat transfer to mass transfer";
  Modelica.Blocks.Math.Product mWat "Water flow rate";
  Modelica.Blocks.Math.Min min;
  Modelica.Blocks.Sources.Constant zero(k=0) "Constant for zero";
  Modelica.Blocks.Math.Add delX(k2=-1, k1=1) "Humidity difference";
protected 
 parameter Medium.ThermodynamicState sta0 = Medium.setState_phX(h=Medium.h_default,
       p=Medium.p_default, X=Medium.X_default);
 parameter Modelica.SIunits.SpecificHeatCapacity cp=Medium.specificHeatCapacityCp(sta0) 
    "Density, used to compute fluid volume";
 parameter Real cpLe(unit="J/(kg.K)") = cp * Le^(1-n);
equation 
  connect(TSur, TDewPoi.T);
  connect(TDewPoi.p_w, humRatPre.p_w);
  connect(TSur, TLiq);
  connect(Gc, gain.u);
  connect(gain.y, mWat.u2);
  connect(mWat.y, mWat_flow);
  connect(zero.y,min. u1);
  connect(delX.u2,XInf);
  connect(humRatPre.XWat, delX.u1);
  connect(delX.y, min.u2);
  connect(min.y, mWat.u1);
end MassExchange;

Buildings.Fluids.HeatExchangers.BaseClasses.PartialDuctManifold Buildings.Fluids.HeatExchangers.BaseClasses.PartialDuctManifold

Partial manifold for heat exchanger duct connection

Buildings.Fluids.HeatExchangers.BaseClasses.PartialDuctManifold

Information


Partial duct manifold for a heat exchanger.

This model defines the duct connection to a heat exchanger. It is extended by other models that model the flow connection between the ports with and without flow friction.


Extends from PartialDuctPipeManifold (Partial heat exchanger duct and pipe manifold).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
IntegernPipPar Number of parallel pipes in each register
IntegernPipSeg Number of pipe segments per register used for discretization
Initialization
MassFlowRatemStart_flow_a Guess value for mass flow rate at port_a [kg/s]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)

Connectors

TypeNameDescription
FluidPort_aport_aFluid connector a for medium (positive design flow direction is from port_a to port_b)
FluidPort_bport_b[nPipPar, nPipSeg]Fluid connector b for medium (positive design flow direction is from port_a to port_b)

Modelica definition

partial model PartialDuctManifold 
  "Partial manifold for heat exchanger duct connection"
  extends PartialDuctPipeManifold;
  parameter Integer nPipSeg(min=1) 
    "Number of pipe segments per register used for discretization";

  Modelica_Fluid.Interfaces.FluidPort_b[nPipPar,nPipSeg] port_b(
        redeclare each package Medium = Medium,
        each m_flow(start=-mStart_flow_a/nPipSeg/nPipPar,
             max=if allowFlowReversal then +Modelica.Constants.inf else 0)) 
    "Fluid connector b for medium (positive design flow direction is from port_a to port_b)";
end PartialDuctManifold;

Buildings.Fluids.HeatExchangers.BaseClasses.PartialDuctPipeManifold Buildings.Fluids.HeatExchangers.BaseClasses.PartialDuctPipeManifold

Partial heat exchanger duct and pipe manifold

Buildings.Fluids.HeatExchangers.BaseClasses.PartialDuctPipeManifold

Information


Partial heat exchanger manifold. This partial model defines ports and parameters that are used for air-side and water-side heat exchanger manifolds.


Extends from Buildings.BaseClasses.BaseIcon (Base icon).

Parameters

TypeNameDefaultDescription
IntegernPipPar Number of parallel pipes in each register
Initialization
MassFlowRatemStart_flow_a Guess value for mass flow rate at port_a [kg/s]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)

Connectors

TypeNameDescription
FluidPort_aport_aFluid connector a for medium (positive design flow direction is from port_a to port_b)

Modelica definition

partial model PartialDuctPipeManifold 
  "Partial heat exchanger duct and pipe manifold"
  extends Buildings.BaseClasses.BaseIcon;

  outer Modelica_Fluid.System system "System wide properties";

  replaceable package Medium = 
      Modelica.Media.Interfaces.PartialMedium "Medium in the component";

  parameter Boolean allowFlowReversal = system.allowFlowReversal 
    "= true to allow flow reversal, false restricts to design direction (port_a -> port_b)";

  parameter Integer nPipPar(min=1) "Number of parallel pipes in each register";

  parameter Modelica.SIunits.MassFlowRate mStart_flow_a 
    "Guess value for mass flow rate at port_a";

  Modelica_Fluid.Interfaces.FluidPort_a port_a(
        redeclare package Medium = Medium,
        m_flow(start=mStart_flow_a, min=if allowFlowReversal then -Modelica.Constants.inf else 0)) 
    "Fluid connector a for medium (positive design flow direction is from port_a to port_b)";
end PartialDuctPipeManifold;

Buildings.Fluids.HeatExchangers.BaseClasses.PartialPipeManifold Buildings.Fluids.HeatExchangers.BaseClasses.PartialPipeManifold

Partial pipe manifold for a heat exchanger

Buildings.Fluids.HeatExchangers.BaseClasses.PartialPipeManifold

Information


Partial pipe manifold for a heat exchanger.

This model defines the pipe connection to a heat exchanger. It is extended by other models that model the flow connection between the ports with and without flow friction.


Extends from PartialDuctPipeManifold (Partial heat exchanger duct and pipe manifold).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
IntegernPipPar Number of parallel pipes in each register
Initialization
MassFlowRatemStart_flow_a Guess value for mass flow rate at port_a [kg/s]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)

Connectors

TypeNameDescription
FluidPort_aport_aFluid connector a for medium (positive design flow direction is from port_a to port_b)
FluidPort_bport_b[nPipPar]Fluid connector b for medium (positive design flow direction is from port_a to port_b)

Modelica definition

partial model PartialPipeManifold 
  "Partial pipe manifold for a heat exchanger"
  extends PartialDuctPipeManifold;
  Modelica_Fluid.Interfaces.FluidPort_b[nPipPar] port_b(
        redeclare each package Medium = Medium,
        each m_flow(start=-mStart_flow_a/nPipPar, max=if allowFlowReversal then +Modelica.Constants.inf else 0)) 
    "Fluid connector b for medium (positive design flow direction is from port_a to port_b)";
end PartialPipeManifold;

Buildings.Fluids.HeatExchangers.BaseClasses.PipeManifoldFixedResistance Buildings.Fluids.HeatExchangers.BaseClasses.PipeManifoldFixedResistance

Pipe manifold for a heat exchanger connection

Buildings.Fluids.HeatExchangers.BaseClasses.PipeManifoldFixedResistance

Information


Pipe manifold with a fixed flow resistance.

This model causes the flow to be distributed equally into each flow path by using a fixed flow resistance for each flow path.


Extends from PartialPipeManifold (Partial pipe manifold for a heat exchanger).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
IntegernPipPar Number of parallel pipes in each register
Booleanuse_dhfalseSet to true to specify hydraulic diameter
Lengthdh0.025Hydraulic diameter for each pipe [m]
RealReC4000Reynolds number where transition to turbulent starts
RealdeltaM0.3Fraction of nominal mass flow rate where transition to turbulent occurs
Initialization
MassFlowRatemStart_flow_a Guess value for mass flow rate at port_a [kg/s]
Nominal Condition
MassFlowRatem_flow_nominal Mass flow rate at port_a [kg/s]
Pressuredp_nominal Pressure [Pa]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
Booleanlinearizedfalse= true, use linear relation between m_flow and dp for any flow rate
Booleanfrom_dpfalse= true, use m_flow = f(dp) else dp = f(m_flow)

Connectors

TypeNameDescription
FluidPort_aport_aFluid connector a for medium (positive design flow direction is from port_a to port_b)
FluidPort_bport_b[nPipPar]Fluid connector b for medium (positive design flow direction is from port_a to port_b)

Modelica definition

model PipeManifoldFixedResistance 
  "Pipe manifold for a heat exchanger connection"
  extends PartialPipeManifold;

  parameter Modelica.SIunits.MassFlowRate m_flow_nominal 
    "Mass flow rate at port_a";
  parameter Modelica.SIunits.Pressure dp_nominal(min=0) "Pressure";

  parameter Boolean use_dh = false "Set to true to specify hydraulic diameter";
  parameter Modelica.SIunits.Length dh=0.025 "Hydraulic diameter for each pipe";
  parameter Real ReC=4000 
    "Reynolds number where transition to turbulent starts";
  parameter Boolean linearized = false 
    "= true, use linear relation between m_flow and dp for any flow rate";
  parameter Real deltaM(min=0) = 0.3 
    "Fraction of nominal mass flow rate where transition to turbulent occurs";
  parameter Boolean from_dp = false 
    "= true, use m_flow = f(dp) else dp = f(m_flow)";

  Fluids.FixedResistances.FixedResistanceDpM[nPipPar] fixRes(
    redeclare each package Medium = Medium,
    each m_flow_nominal=m_flow_nominal/nPipPar,
    each m_flow(start=mStart_flow_a),
    each dp_nominal=dp_nominal,
    each dh=dh,
    each from_dp=from_dp,
    each deltaM=deltaM,
    each ReC=ReC,
    each use_dh=use_dh,
    each linearized=linearized) "Fixed resistance for each duct";
equation 
  for i in 1:nPipPar loop
    connect(port_a, fixRes[i].port_a);
    connect(fixRes[i].port_b, port_b[i]);
  end for;
end PipeManifoldFixedResistance;

Buildings.Fluids.HeatExchangers.BaseClasses.PipeManifoldNoResistance Buildings.Fluids.HeatExchangers.BaseClasses.PipeManifoldNoResistance

Manifold for heat exchanger register

Buildings.Fluids.HeatExchangers.BaseClasses.PipeManifoldNoResistance

Information


Pipe manifold without flow resistance.

This model connects the flows between the ports without modeling flow friction. The model is used in conjunction with a manifold which contains pressure drop elements and that is added to the other side of the heat exchanger registers.


Extends from PartialPipeManifold (Partial pipe manifold for a heat exchanger).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
IntegernPipPar Number of parallel pipes in each register
BooleanconnectAllPressurestrue 
Initialization
MassFlowRatemStart_flow_a Guess value for mass flow rate at port_a [kg/s]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)

Connectors

TypeNameDescription
FluidPort_aport_aFluid connector a for medium (positive design flow direction is from port_a to port_b)
FluidPort_bport_b[nPipPar]Fluid connector b for medium (positive design flow direction is from port_a to port_b)

Modelica definition

model PipeManifoldNoResistance "Manifold for heat exchanger register"
  extends PartialPipeManifold;
 parameter Boolean connectAllPressures=true;
  Modelica_Fluid.Fittings.MultiPort mulPor(
      redeclare package Medium = Medium,
      final nPorts_b=nPipPar);
equation 
  connect(port_a, mulPor.port_a);
  connect(mulPor.ports_b, port_b);
end PipeManifoldNoResistance;

HTML-documentation generated by Dymola Fri May 15 10:14:36 2009.