Buildings.Fluid.FixedResistances

Package with models for fixed flow resistances

Information

This package contains component models for fixed flow resistances. By fixed flow resistance, we mean resistances that do not change the flow coefficient

k = m ⁄ √Δp.

For models of valves and air dampers, see Buildings.Fluid.Actuators. For models of flow resistances as part of the building constructions, see Buildings.Airflow.Multizone.

The model Buildings.Fluid.FixedResistances.PressureDrop is a fixed flow resistance that takes as parameter a nominal flow rate and a nominal pressure drop. The actual resistance is scaled using the above equation.

The model Buildings.Fluid.FixedResistances.HydraulicDiameter is a fixed flow resistance that takes as parameter a nominal flow rate and a hydraulic diameter. The actual resistance is scaled using the above equation.

The model Buildings.Fluid.FixedResistances.LosslessPipe is an ideal pipe segment with no pressure drop. It is primarily used in models in which the above pressure drop model need to be replaced by a model with no pressure drop.

The model Buildings.Fluid.FixedResistances.Junction can be used to model flow splitters or flow merges.

Extends from Modelica.Icons.VariantsPackage (Icon for package containing variants).

Package Content

Name Description
Buildings.Fluid.FixedResistances.CheckValve CheckValve Check valve that avoids flow reversal
Buildings.Fluid.FixedResistances.HydraulicDiameter HydraulicDiameter Fixed flow resistance with hydraulic diameter and m_flow as parameter
Buildings.Fluid.FixedResistances.Junction Junction Flow splitter with fixed resistance at each port
Buildings.Fluid.FixedResistances.LosslessPipe LosslessPipe Pipe with no flow friction and no heat transfer
Buildings.Fluid.FixedResistances.Pipe Pipe Pipe with finite volume discretization along flow path
Buildings.Fluid.FixedResistances.PlugFlowPipe PlugFlowPipe Pipe model using spatialDistribution for temperature delay
Buildings.Fluid.FixedResistances.PlugFlowPipeDiscretized PlugFlowPipeDiscretized Discretized pipe model using spatialDistribution for temperature delay
Buildings.Fluid.FixedResistances.PressureDrop PressureDrop Fixed flow resistance with dp and m_flow as parameter
Buildings.Fluid.FixedResistances.BuriedPipes BuriedPipes Package with models for buried pipes
Buildings.Fluid.FixedResistances.Examples Examples Collection of models that illustrate model use and test models
Buildings.Fluid.FixedResistances.Validation Validation Collection of validation models
Buildings.Fluid.FixedResistances.BaseClasses BaseClasses Package with base classes for Buildings.Fluid.FixedResistances

Buildings.Fluid.FixedResistances.CheckValve Buildings.Fluid.FixedResistances.CheckValve

Check valve that avoids flow reversal

Buildings.Fluid.FixedResistances.CheckValve

Information

Implementation of a hydraulic check valve. Note that the small reverse flows can still occur with this model.

Main equations

The basic flow function

ṁ = sign(Δp) k √ Δp  ,

with regularization near the origin, is used to compute the pressure drop. The flow coefficient

k = ṁ ⁄ √ Δp  

is increased from l*KV_Si to KV_Si, where KV_Si is equal to Kv but in SI units. Therefore, the flow coefficient k is set to a value close to zero for negative pressure differences, thereby restricting reverse flow to a small value. The flow coefficient k saturates to its maximum value at the pressure dpValve_closing. For larger pressure drops, the pressure drop is a quadratic function of the flow rate.

Typical use and important parameters

The parameters m_flow_nominal and dpValve_nominal determine the flow coefficient of the check valve when it is fully opened. A typical value for a nominal flow rate of 1 m/s is dpValve_nominal = 3400 Pa. The leakage ratio l determines the minimum flow coefficient, for negative pressure differences. The parameter dpFixed_nominal allows to include a series pressure drop with a fixed flow coefficient into the model. The parameter dpValve_closing determines when the flow coefficient starts to increase, which is typically in the order of dpValve_nominal.

Implementation

The check valve implementation approximates the physics where a forward pressure difference opens the valve such that the valve opening increases, causing a growing orifice area and thus increasing the flow coefficient. Near dp=dpValve_closing, the valve is fully open and the flow coefficient saturates to the flow coefficient value determined by dpValve_nominal and m_flow_nominal. For typical valve diameters, the check valve is only fully open near nominal mass flow rate. Therefore, the model sets dpValve_closing=dpValve_nominal/2 by default.

Extends from Buildings.Fluid.BaseClasses.PartialResistance (Partial model for a hydraulic resistance), Buildings.Fluid.Actuators.BaseClasses.ValveParameters (Model with parameters for valves).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
MassFlowRatem_flow_turbulentdeltaM*abs(m_flow_nominal)Turbulent flow if |m_flow| >= m_flow_turbulent [kg/s]
Reall0.001Valve leakage, l=Kv(y=0)/Kv(y=1)
RealkFixedif dpFixed_nominal > Modelic...Flow coefficient of fixed resistance that may be in series with valve, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2).
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
PressureDifferencedp_nominaldpValve_nominal + dpFixed_no...Pressure drop at nominal mass flow rate [Pa]
PressureDifferencedpValve_nominal Nominal pressure drop of fully open valve, used if CvData=Buildings.Fluid.Types.CvTypes.OpPoint [Pa]
PressureDifferencedpFixed_nominal0Pressure drop of pipe and other resistances that are in series [Pa]
PressureDifferencedpValve_closingdpValve_nominal/2Pressure drop when the check valve starts to close [Pa]
Flow Coefficient
CvTypesCvDataBuildings.Fluid.Types.CvType...Selection of flow coefficient
RealKv Kv (metric) flow coefficient [m3/h/(bar)^(1/2)]
RealCv Cv (US) flow coefficient [USG/min/(psi)^(1/2)]
AreaAv Av (metric) flow coefficient [m2]
Pressure-flow linearization
RealdeltaM0.02Fraction of nominal flow rate where linearization starts, if y=1
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
Booleanfrom_dptrue= true, use m_flow = f(dp) else dp = f(m_flow)
Booleanlinearizedfalse= true, use linear relation between m_flow and dp for any flow rate
Nominal condition
DensityrhoStdMedium.density_pTX(101325, 2...Inlet density for which valve coefficients are defined [kg/m3]

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)

Modelica definition

model CheckValve "Check valve that avoids flow reversal" extends Buildings.Fluid.BaseClasses.PartialResistance( dp(nominal=2000), final dp_nominal=dpValve_nominal + dpFixed_nominal, final m_flow_turbulent=deltaM*abs(m_flow_nominal), final from_dp=true, final linearized=false, allowFlowReversal=true); extends Buildings.Fluid.Actuators.BaseClasses.ValveParameters( rhoStd=Medium.density_pTX(101325, 273.15 + 4, Medium.X_default)); parameter Modelica.Units.SI.PressureDifference dpFixed_nominal( displayUnit="Pa", min=0) = 0 "Pressure drop of pipe and other resistances that are in series"; parameter Modelica.Units.SI.PressureDifference dpValve_closing= dpValve_nominal/2 "Pressure drop when the check valve starts to close"; parameter Real l(min=1e-10, max=1)=0.001 "Valve leakage, l=Kv(y=0)/Kv(y=1)"; parameter Real kFixed(unit="", min=0)= if dpFixed_nominal > Modelica.Constants.eps then m_flow_nominal/sqrt(dpFixed_nominal) else 0 "Flow coefficient of fixed resistance that may be in series with valve, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2)."; Real k(min=Modelica.Constants.small) "Flow coefficient of valve and pipe in series in allowed/forward direction, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2)."; protected Real a "Scaled pressure variable"; Real cv "Twice differentiable Heaviside check valve characteristic"; Real kCv "Smoothed restriction characteristic"; initial equation assert(dpFixed_nominal > -Modelica.Constants.eps, "In " + getInstanceName() + ": We require dpFixed_nominal >= 0. Received dpFixed_nominal = " + String(dpFixed_nominal) + " Pa."); assert(l > -Modelica.Constants.eps, "In " + getInstanceName() + ": We require l >= 0. Received l = " + String(l)); equation a = dp/dpValve_closing; cv = smooth(2, max(0, min(1, a^3*(10+a*(-15+6*a))))); kCv = Kv_SI*(cv*(1-l) + l); if (dpFixed_nominal > Modelica.Constants.eps) then k = sqrt(1/(1/kFixed^2 + 1/kCv^2)); else k = kCv; end if; if homotopyInitialization then m_flow = homotopy( actual=Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_dp( dp = dp, k = k, m_flow_turbulent = m_flow_turbulent), simplified = m_flow_nominal_pos*dp/dp_nominal_pos); else m_flow = Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_dp( dp = dp, k = k, m_flow_turbulent = m_flow_turbulent); end if; end CheckValve;

Buildings.Fluid.FixedResistances.HydraulicDiameter Buildings.Fluid.FixedResistances.HydraulicDiameter

Fixed flow resistance with hydraulic diameter and m_flow as parameter

Buildings.Fluid.FixedResistances.HydraulicDiameter

Information

This is a model of a flow resistance with a fixed flow coefficient. The mass flow rate is computed as

ṁ = k √Δp,

where k is a constant and Δp is the pressure drop. The constant k is equal to k=m_flow_nominal/sqrt(dp_nominal), where m_flow_nominal is a parameter.

Assumptions

In the region abs(m_flow) < m_flow_turbulent, the square root is replaced by a differentiable function with finite slope. The value of m_flow_turbulent is computed as m_flow_turbulent = eta_nominal*dh/4*π*ReC, where eta_nominal is the dynamic viscosity, obtained from the medium model. The parameter dh is the hydraulic diameter and ReC=4000 is the critical Reynolds number, which both can be set by the user.

Important parameters

By default, the pressure drop at nominal flow rate is computed as

dp_nominal = fac * dpStraightPipe_nominal,

where dpStraightPipe_nominal is a parameter that is automatically computed based on the nominal mass flow rate, hydraulic diameter, pipe roughness and medium properties. The hydraulic diameter dh is by default computed based on the flow velocity v_nominal and the nominal mass flow rate m_flow_nominal. Hence, users should change the default values of dh or v_nominal if they are not applicable for their model.

The factor fac takes into account additional resistances such as for bends. The default value of 2 can be changed by the user.

The parameter from_dp is used to determine whether the mass flow rate is computed as a function of the pressure drop (if from_dp=true), or vice versa. This setting can affect the size of the nonlinear system of equations.

If the parameter linearized is set to true, then the pressure drop is computed as a linear function of the mass flow rate.

Setting allowFlowReversal=false can lead to simpler equations. However, this should only be set to false if one can guarantee that the flow never reverses its direction. This can be difficult to guarantee, as pressure imbalance after the initialization, or due to medium expansion and contraction, can lead to reverse flow.

If the parameter show_T is set to true, then the model will compute the temperature at its ports. Note that this can lead to state events when the mass flow rate approaches zero, which can increase computing time.

Notes

For more detailed models that compute the actual flow friction, models from the package Modelica.Fluid can be used and combined with models from the Buildings library.

For a model that uses dp_nominal as a parameter rather than geoemetric data, use Buildings.Fluid.FixedResistances.PressureDrop.

Implementation

The pressure drop is computed by calling a function in the package Buildings.Fluid.BaseClasses.FlowModels, This package contains regularized implementations of the equation

ṁ = sign(Δp) k √ Δp  

and its inverse function.

To decouple the energy equation from the mass equations, the pressure drop is a function of the mass flow rate, and not the volume flow rate. This leads to simpler equations.

Extends from Buildings.Fluid.FixedResistances.PressureDrop (Fixed flow resistance with dp and m_flow as parameter).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Lengthdhsqrt(4*m_flow_nominal/rho_de...Hydraulic diameter (assuming a round cross section area) [m]
Lengthlength Length of the pipe [m]
RealReC4000Reynolds number where transition to turbulence starts
Lengthroughness2.5e-5Absolute roughness of pipe, with a default for a smooth steel pipe (dummy if use_roughness = false) [m]
Realfac2Factor to take into account resistance of bends etc., fac=dp_nominal/dpStraightPipe_nominal
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
PressureDifferencedp_nominalfac*dpStraightPipe_nominalPressure drop at nominal mass flow rate [Pa]
Velocityv_nominalif rho_default < 500 then 1....Velocity at m_flow_nominal (used to compute default value for hydraulic diameter dh) [m/s]
Transition to laminar
RealdeltaMeta_default*dh/4*Modelica.Co...Fraction of nominal mass flow rate where transition to turbulent occurs
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
Booleanfrom_dpfalse= true, use m_flow = f(dp) else dp = f(m_flow)
Booleanlinearizedfalse= true, use linear relation between m_flow and dp for any flow rate

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)

Modelica definition

model HydraulicDiameter "Fixed flow resistance with hydraulic diameter and m_flow as parameter" extends Buildings.Fluid.FixedResistances.PressureDrop( final deltaM = eta_default*dh/4*Modelica.Constants.pi*ReC/m_flow_nominal_pos, final dp_nominal=fac*dpStraightPipe_nominal); parameter Modelica.Units.SI.Length dh=sqrt(4*m_flow_nominal/rho_default/ v_nominal/Modelica.Constants.pi) "Hydraulic diameter (assuming a round cross section area)"; parameter Modelica.Units.SI.Length length "Length of the pipe"; parameter Real ReC(min=0)=4000 "Reynolds number where transition to turbulence starts"; parameter Modelica.Units.SI.Velocity v_nominal=if rho_default < 500 then 1.5 else 0.15 "Velocity at m_flow_nominal (used to compute default value for hydraulic diameter dh)"; parameter Modelica.Units.SI.Length roughness(min=0) = 2.5e-5 "Absolute roughness of pipe, with a default for a smooth steel pipe (dummy if use_roughness = false)"; parameter Real fac(min=1) = 2 "Factor to take into account resistance of bends etc., fac=dp_nominal/dpStraightPipe_nominal"; final parameter Modelica.Units.SI.PressureDifference dpStraightPipe_nominal( displayUnit="Pa") = Modelica.Fluid.Pipes.BaseClasses.WallFriction.Detailed.pressureLoss_m_flow( m_flow=m_flow_nominal, rho_a=rho_default, rho_b=rho_default, mu_a=mu_default, mu_b=mu_default, length=length, diameter=dh, roughness=roughness, m_flow_small=m_flow_small) "Pressure loss of a straight pipe at m_flow_nominal"; Modelica.Units.SI.Velocity v=m_flow/(rho_default*ARound) "Flow velocity (assuming a round cross section area)"; protected parameter Modelica.Units.SI.Area ARound=dh^2*Modelica.Constants.pi/4 "Cross sectional area (assuming a round cross section area)"; parameter Medium.ThermodynamicState state_default= Medium.setState_pTX( T=Medium.T_default, p=Medium.p_default, X=Medium.X_default[1:Medium.nXi]) "Default state"; parameter Modelica.Units.SI.Density rho_default=Medium.density(state_default) "Density at nominal condition"; parameter Modelica.Units.SI.DynamicViscosity mu_default= Medium.dynamicViscosity(state_default) "Dynamic viscosity at nominal condition"; end HydraulicDiameter;

Buildings.Fluid.FixedResistances.Junction Buildings.Fluid.FixedResistances.Junction

Flow splitter with fixed resistance at each port

Buildings.Fluid.FixedResistances.Junction

Information

Model of a flow junction with an optional fixed resistance in each flow leg and an optional mixing volume at the junction.

The pressure drop is implemented using the model Buildings.Fluid.FixedResistances.PressureDrop. If its nominal pressure drop is set to zero, then the pressure drop model will be removed. For example, the pressure drop declaration

  m_flow_nominal={ 0.1, 0.1,  -0.2},
  dp_nominal =   {500,    0, -6000}

would model a flow mixer that has the nominal flow rates and associated pressure drops as shown in the figure below. Note that port_3 is set to negative values. The negative values indicate that at the nominal conditions, fluid is leaving the component.

image

If energyDynamics <> Modelica.Fluid.Types.Dynamics.SteadyState, then at the flow junction, a fluid volume is modeled. The fluid volume is implemented using the model Buildings.Fluid.Delays.DelayFirstOrder. The fluid volume has the size

  V = sum(abs(m_flow_nominal[:])/3)*tau/rho_nominal

where tau is a parameter and rho_nominal is the density of the medium in the volume at nominal condition. Setting energyDynamics=Modelica.Fluid.Types.Dynamics.FixedInitial can help reducing the size of the nonlinear system of equations.

Extends from Buildings.Fluid.BaseClasses.PartialThreeWayResistance (Flow splitter with partial resistance model at each port).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Nominal condition
MassFlowRatem_flow_nominal[3] Mass flow rate. Set negative at outflowing ports. [kg/s]
Pressuredp_nominal[3] Pressure drop at nominal mass flow rate, set to zero or negative number at outflowing ports. [Pa]
Transition to laminar
RealdeltaM0.3Fraction of nominal mass flow rate where transition to turbulent occurs
Dynamics
Conservation equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Type of energy balance: dynamic (3 initialization options) or steady state
MassFlowRatemDyn_flow_nominalsum(abs(m_flow_nominal[:])/3)Nominal mass flow rate for dynamic momentum and energy balance [kg/s]
Nominal condition
Timetau10Time constant at nominal flow for dynamic energy and momentum balance [s]
Initialization
AbsolutePressurep_startMedium.p_defaultStart value of pressure [Pa]
TemperatureT_startMedium.T_defaultStart value of temperature [K]
MassFractionX_start[Medium.nX]Medium.X_defaultStart value of mass fractions m_i/m [kg/kg]
ExtraPropertyC_start[Medium.nC]fill(0, Medium.nC)Start value of trace substances
ExtraPropertyC_nominal[Medium.nC]fill(1E-2, Medium.nC)Nominal value of trace substances. (Set to typical order of magnitude.)
Advanced
Booleanfrom_dptrue= true, use m_flow = f(dp) else dp = f(m_flow)
PortFlowDirectionportFlowDirection_1Modelica.Fluid.Types.PortFlo...Flow direction for port_1
PortFlowDirectionportFlowDirection_2Modelica.Fluid.Types.PortFlo...Flow direction for port_2
PortFlowDirectionportFlowDirection_3Modelica.Fluid.Types.PortFlo...Flow direction for port_3
BooleanverifyFlowReversalfalse=true, to assert that the flow does not reverse when portFlowDirection_* does not equal Bidirectional
MassFlowRatem_flow_smallmDyn_flow_nominal*1e-4Small mass flow rate for checking flow reversal [kg/s]
Booleanlinearizedfalse= true, use linear relation between m_flow and dp for any flow rate

Connectors

TypeNameDescription
FluidPort_aport_1First port, typically inlet
FluidPort_bport_2Second port, typically outlet
FluidPort_aport_3Third port, can be either inlet or outlet

Modelica definition

model Junction "Flow splitter with fixed resistance at each port" extends Buildings.Fluid.BaseClasses.PartialThreeWayResistance( m_flow_small=mDyn_flow_nominal*1e-4, mDyn_flow_nominal = sum(abs(m_flow_nominal[:])/3), redeclare Buildings.Fluid.FixedResistances.PressureDrop res1( from_dp=from_dp, final m_flow_nominal=m_flow_nominal[1], final dp_nominal=dp_nominal[1], linearized=linearized, homotopyInitialization=homotopyInitialization, deltaM=deltaM), redeclare Buildings.Fluid.FixedResistances.PressureDrop res2( from_dp=from_dp, final m_flow_nominal=m_flow_nominal[2], final dp_nominal=dp_nominal[2], linearized=linearized, homotopyInitialization=homotopyInitialization, deltaM=deltaM), redeclare Buildings.Fluid.FixedResistances.PressureDrop res3( from_dp=from_dp, final m_flow_nominal=m_flow_nominal[3], final dp_nominal=dp_nominal[3], linearized=linearized, homotopyInitialization=homotopyInitialization, deltaM=deltaM)); constant Boolean homotopyInitialization = true "= true, use homotopy method"; parameter Modelica.Units.SI.MassFlowRate[3] m_flow_nominal "Mass flow rate. Set negative at outflowing ports."; parameter Modelica.Units.SI.Pressure[3] dp_nominal(each displayUnit="Pa") "Pressure drop at nominal mass flow rate, set to zero or negative number at outflowing ports."; parameter Real deltaM(min=0) = 0.3 "Fraction of nominal mass flow rate where transition to turbulent occurs"; parameter Boolean linearized = false "= true, use linear relation between m_flow and dp for any flow rate"; initial equation assert(homotopyInitialization, "In " + getInstanceName() + ": The constant homotopyInitialization has been modified from its default value. This constant will be removed in future releases.", level = AssertionLevel.warning); end Junction;

Buildings.Fluid.FixedResistances.LosslessPipe Buildings.Fluid.FixedResistances.LosslessPipe

Pipe with no flow friction and no heat transfer

Buildings.Fluid.FixedResistances.LosslessPipe

Information

Model of a pipe with no flow resistance, no heat loss and no transport delay. This model can be used to replace a replaceable pipe model in flow legs in which no friction should be modeled. This is for example done in the outlet port of the base class for three way valves, Buildings.Fluid.Actuators.BaseClasses.PartialThreeWayValve.

Extends from Buildings.Fluid.Interfaces.PartialTwoPortInterface (Partial model with two ports and declaration of quantities that are used by many models).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
MassFlowRatem_flow_small1E-4*abs(m_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed

Connectors

TypeNameDescription
FluidPort_aport_aFluid connector a (positive design flow direction is from port_a to port_b)
FluidPort_bport_bFluid connector b (positive design flow direction is from port_a to port_b)

Modelica definition

model LosslessPipe "Pipe with no flow friction and no heat transfer" extends Buildings.Fluid.Interfaces.PartialTwoPortInterface; final parameter Boolean from_dp=true "Used to satisfy replaceable models"; equation dp=0; // Isenthalpic state transformation (no storage and no loss of energy) port_a.h_outflow = if allowFlowReversal then inStream(port_b.h_outflow) else Medium.h_default; port_b.h_outflow = inStream(port_a.h_outflow); // Mass balance (no storage) port_a.m_flow + port_b.m_flow = 0; // Transport of substances port_a.Xi_outflow = if allowFlowReversal then inStream(port_b.Xi_outflow) else Medium.X_default[1:Medium.nXi]; port_b.Xi_outflow = inStream(port_a.Xi_outflow); port_a.C_outflow = if allowFlowReversal then inStream(port_b.C_outflow) else zeros(Medium.nC); port_b.C_outflow = inStream(port_a.C_outflow); end LosslessPipe;

Buildings.Fluid.FixedResistances.Pipe Buildings.Fluid.FixedResistances.Pipe

Pipe with finite volume discretization along flow path

Buildings.Fluid.FixedResistances.Pipe

Information

Model of a pipe with flow resistance and optional heat exchange with environment.

Heat loss calculation

There are two possible configurations:

  1. If useMultipleHeatPorts=false (default option), the pipe uses a single heat port for the heat exchange with the environment. Note that if the heat port is unconnected, then all volumes are still connected through the heat conduction elements conPipWal. Therefore, they exchange a small amount of heat, which is not physical. To avoid this, set useMultipleHeatPorts=true.
  2. If useMultipleHeatPorts=true, then one heat port for each segment of the pipe is used for the heat exchange with the environment. If the heat port is unconnected, then the pipe has no heat loss.

Pressure drop calculation

The default value for the parameter diameter is computed such that the flow velocity is equal to v_nominal=0.15 for a mass flow rate of m_flow_nominal. Both parameters, diameter and v_nominal, can be overwritten by the user. The default value for dp_nominal is two times the pressure drop that the pipe would have if it were straight with no fittings. The factor of two that takes into account the pressure loss of fittings can be overwritten. These fittings could also be explicitly modeled outside of this component using models from the package Modelica.Fluid.Fittings. For mass flow rates other than m_flow_nominal, the model Buildings.Fluid.FixedResistances.PressureDrop is used to compute the pressure drop.

For a steady-state model of a flow resistance, use Buildings.Fluid.FixedResistances.PressureDrop instead of this model.

Extends from Buildings.Fluid.FixedResistances.BaseClasses.Pipe (Model of a pipe with finite volume discretization along the flow path).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
IntegernSeg10Number of volume segments
LengththicknessIns Thickness of insulation [m]
ThermalConductivitylambdaIns Heat conductivity of insulation [W/(m.K)]
Lengthdiametersqrt(4*m_flow_nominal/rho_de...Pipe diameter (without insulation) [m]
Lengthlength Length of the pipe [m]
Velocityv_nominal0.15Velocity at m_flow_nominal (used to compute default diameter) [m/s]
Lengthroughness2.5e-5Absolute roughness of pipe, with a default for a smooth steel pipe (dummy if use_roughness = false) [m]
BooleanuseMultipleHeatPortsfalse= true to use one heat port for each segment of the pipe, false to use a single heat port for the entire pipe
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
PressureDifferencedp_nominal2*dpStraightPipe_nominalPressure difference [Pa]
Dynamics
Conservation equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Type of energy balance: dynamic (3 initialization options) or steady state
RealmSenFac1Factor for scaling the sensible thermal mass of the volume
Advanced
Dynamics
DynamicsmassDynamicsenergyDynamicsType of mass balance: dynamic (3 initialization options) or steady state, must be steady state if energyDynamics is steady state
MassFlowRatem_flow_small1E-4*abs(m_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
Initialization
AbsolutePressurep_startMedium.p_defaultStart value of pressure [Pa]
TemperatureT_startMedium.T_defaultStart value of temperature [K]
MassFractionX_start[Medium.nX]Medium.X_defaultStart value of mass fractions m_i/m [kg/kg]
ExtraPropertyC_start[Medium.nC]fill(0, Medium.nC)Start value of trace substances
ExtraPropertyC_nominal[Medium.nC]fill(1E-2, Medium.nC)Nominal value of trace substances. (Set to typical order of magnitude.)
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Flow resistance
Booleanfrom_dpfalse= true, use m_flow = f(dp) else dp = f(m_flow)
BooleanlinearizeFlowResistancefalse= true, use linear relation between m_flow and dp for any flow rate
RealdeltaM0.1Fraction of nominal flow rate where flow transitions to laminar
RealReC4000Reynolds number where transition to turbulence starts

Connectors

TypeNameDescription
FluidPort_aport_aFluid connector a (positive design flow direction is from port_a to port_b)
FluidPort_bport_bFluid connector b (positive design flow direction is from port_a to port_b)
HeatPort_aheatPortSingle heat port that connects to outside of pipe wall (default, enabled when useMultipleHeatPorts=false)
HeatPorts_aheatPorts[nSeg]Multiple heat ports that connect to outside of pipe wall (enabled if useMultipleHeatPorts=true)

Modelica definition

model Pipe "Pipe with finite volume discretization along flow path" extends Buildings.Fluid.FixedResistances.BaseClasses.Pipe( diameter=sqrt(4*m_flow_nominal/rho_default/v_nominal/Modelica.Constants.pi), dp_nominal=2*dpStraightPipe_nominal, preDro(dp(nominal=length*10))); // Because dp_nominal is a non-literal value, we set // dp.nominal=100 instead of the default dp.nominal=dp_nominal, // because the latter is ignored by Dymola 2012 FD 01. parameter Modelica.Units.SI.Velocity v_nominal=0.15 "Velocity at m_flow_nominal (used to compute default diameter)"; parameter Modelica.Units.SI.Length roughness(min=0) = 2.5e-5 "Absolute roughness of pipe, with a default for a smooth steel pipe (dummy if use_roughness = false)"; final parameter Modelica.Units.SI.PressureDifference dpStraightPipe_nominal( displayUnit="Pa") = Modelica.Fluid.Pipes.BaseClasses.WallFriction.Detailed.pressureLoss_m_flow( m_flow=m_flow_nominal, rho_a=rho_default, rho_b=rho_default, mu_a=mu_default, mu_b=mu_default, length=length, diameter=diameter, roughness=roughness, m_flow_small=m_flow_small) "Pressure loss of a straight pipe at m_flow_nominal"; parameter Boolean useMultipleHeatPorts=false "= true to use one heat port for each segment of the pipe, false to use a single heat port for the entire pipe"; Modelica.Thermal.HeatTransfer.Components.ThermalConductor conPipWal[nSeg]( each G=2*Modelica.Constants.pi*lambdaIns*length/nSeg/Modelica.Math.log(( diameter/2.0 + thicknessIns)/(diameter/2.0))) "Thermal conductance through pipe wall"; Modelica.Thermal.HeatTransfer.Components.ThermalCollector colAllToOne(m=nSeg) if not useMultipleHeatPorts "Connector to assign multiple heat ports to one heat port"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a heatPort if not useMultipleHeatPorts "Single heat port that connects to outside of pipe wall (default, enabled when useMultipleHeatPorts=false)"; Modelica.Fluid.Interfaces.HeatPorts_a heatPorts[nSeg] if useMultipleHeatPorts "Multiple heat ports that connect to outside of pipe wall (enabled if useMultipleHeatPorts=true)"; equation connect(conPipWal.port_b, vol.heatPort); if useMultipleHeatPorts then connect(heatPorts, conPipWal.port_a); else connect(colAllToOne.port_a, conPipWal.port_a); connect(colAllToOne.port_b, heatPort); end if; end Pipe;

Buildings.Fluid.FixedResistances.PlugFlowPipe Buildings.Fluid.FixedResistances.PlugFlowPipe

Pipe model using spatialDistribution for temperature delay

Buildings.Fluid.FixedResistances.PlugFlowPipe

Information

Pipe with heat loss using the time delay based heat losses and transport of the fluid using a plug flow model, applicable for simulation of long pipes such as in district heating and cooling systems.

This model takes into account transport delay along the pipe length idealized as a plug flow. The model also includes thermal inertia of the pipe wall.

Implementation

The spatialDistribution operator is used for the temperature wave propagation through the length of the pipe. This operator is contained in Buildings.Fluid.FixedResistances.BaseClasses.PlugFlow.

The model Buildings.Fluid.FixedResistances.BaseClasses.PlugFlowHeatLoss implements a heat loss in design direction, but leaves the enthalpy unchanged in opposite flow direction. Therefore it is used in front of and behind the time delay.

The pressure drop is implemented using Buildings.Fluid.FixedResistances.HydraulicDiameter.

The thermal capacity of the pipe wall is implemented as a mixing volume of the fluid in the pipe, of which the thermal capacity is equal to that of the pipe wall material. In addition, this mixing volume allows the hydraulic separation of subsequent pipes.
The mixing volume is either split between the inlet and outlet ports (port_a and port_b) or lumped in at the outlet (port_b) if have_symmetry is set to false. This mixing volume can be removed from this model with the Boolean parameter have_pipCap, in cases where the pipe wall heat capacity is negligible and a state is not needed at the pipe outlet (see the note below about numerical Jacobians).

Note that in order to model a branched network it is recommended to use Buildings.Fluid.FixedResistances.Junction at each junction and to configure that junction model with a state (energyDynamics <> Modelica.Fluid.Types.Dynamics.SteadyState), see for instance Buildings.Fluid.FixedResistances.Validation.PlugFlowPipes.PlugFlowAIT. This will avoid the numerical Jacobian that is otherwise created when the inlet ports of two instances of the plug flow model are connected together.

Assumptions

References

Full details on the model implementation and experimental validation can be found in:

van der Heijde, B., Fuchs, M., Ribas Tugores, C., Schweiger, G., Sartor, K., Basciotti, D., Müller, D., Nytsch-Geusen, C., Wetter, M. and Helsen, L. (2017).
Dynamic equation-based thermo-hydraulic pipe model for district heating and cooling systems.
Energy Conversion and Management, vol. 151, p. 158-169. doi: 10.1016/j.enconman.2017.08.072.

Extends from Buildings.Fluid.FixedResistances.BaseClasses.PlugFlowPipe (Pipe model using spatialDistribution for temperature delay).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
RealReC4000Reynolds number where transition to turbulence starts
Realfac1Factor to take into account flow resistance of bends etc., fac=dp_nominal/dpStraightPipe_nominal
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
Velocityv_nominal1.5Velocity at m_flow_nominal (used to compute default value for hydraulic diameter dh) [m/s]
Material
Lengthdhsqrt(4*m_flow_nominal/rho_de...Hydraulic diameter (assuming a round cross section area) [m]
Heightroughness2.5e-5Average height of surface asperities (default: smooth steel pipe) [m]
Lengthlength Pipe length [m]
SpecificHeatCapacitycPip2300Specific heat of pipe wall material. 2300 for PE, 500 for steel [J/(kg.K)]
DensityrhoPip930Density of pipe wall material. 930 for PE, 8000 for steel [kg/m3]
Lengththickness0.0035Pipe wall thickness [m]
Thermal resistance
LengthdIns Thickness of pipe insulation, used to compute R [m]
ThermalConductivitykIns Heat conductivity of pipe insulation, used to compute R [W/(m.K)]
RealR1/(kIns*2*Modelica.Constants...Thermal resistance per unit length from fluid to boundary temperature [(m.K)/W]
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
MassFlowRatem_flow_small1E-4*abs(m_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
Booleanfrom_dpfalse= true, use m_flow = f(dp) else dp = f(m_flow)
Booleanhave_pipCaptrue= true, a mixing volume is added that corresponds to the heat capacity of the pipe wall
Booleanhave_symmetrytrue= false, the mixing volume is only on port_b, which improve performances, but reduces dynamic accuracy.
Booleanlinearizedfalse= true, use linear relation between m_flow and dp for any flow rate
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
Initialization
TemperatureT_start_inMedium.T_defaultInitialization temperature at pipe inlet [K]
TemperatureT_start_outT_start_inInitialization temperature at pipe outlet [K]
BooleaninitDelayfalseInitialize delay for a constant mass flow rate if true, otherwise start from 0
MassFlowRatem_flow_start0Initial value of mass flow rate through pipe [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)
HeatPort_aheatPortHeat transfer to or from surroundings (positive if pipe is colder than surrounding)

Modelica definition

model PlugFlowPipe "Pipe model using spatialDistribution for temperature delay" extends Buildings.Fluid.FixedResistances.BaseClasses.PlugFlowPipe( redeclare final Buildings.Fluid.FixedResistances.HydraulicDiameter res( final dh=dh, final from_dp=from_dp, final length=length, final roughness=roughness, final fac=fac, final ReC=ReC, final v_nominal=v_nominal, final homotopyInitialization=homotopyInitialization, final linearized=linearized, dp(nominal=fac*200*length))); end PlugFlowPipe;

Buildings.Fluid.FixedResistances.PlugFlowPipeDiscretized Buildings.Fluid.FixedResistances.PlugFlowPipeDiscretized

Discretized pipe model using spatialDistribution for temperature delay

Buildings.Fluid.FixedResistances.PlugFlowPipeDiscretized

Information

Wrapper around Buildings.Fluid.FixedResistances.PlugFlowPipe which allows to specify nSeg successive segments of pipes (connected in series).

This wrapper simplifies use-cases where different segments of the same pipe might have different boundary conditions. This would be the case, for instance, for sufficiently long stretches of buried pipes.

To reduce coupled nonlinear equations, the pipe flow resistance is aggregated to a single instance of Buildings.Fluid.FixedResistances.HydraulicDiameter rather than being instantiated separately for each segment.

Extends from Buildings.Fluid.Interfaces.PartialTwoPort (Partial component with two ports).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
IntegernSeg1Number of axial segment
RealReC4000Reynolds number where transition to turbulence starts
Realfac1Factor to take into account flow resistance of bends etc., fac=dp_nominal/dpStraightPipe_nominal
Material
Lengthdh Hydraulic diameter [m]
Heightroughness2.5e-5Average height of surface asperities (default: smooth steel pipe) [m]
LengthtotLensum(segLen)Total pipe length (used to compute segment length) [m]
LengthsegLen[nSeg]fill(totLen/nSeg, nSeg)Pipe segment length [m]
Lengththickness0.0035Pipe wall thickness [m]
SpecificHeatCapacitycPip2300Specific heat of pipe wall material. 2300 for PE, 500 for steel [J/(kg.K)]
DensityrhoPip930Density of pipe wall material. 930 for PE, 8000 for steel [kg/m3]
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
Thermal resistance
LengthdIns Thickness of pipe insulation, used to compute R [m]
ThermalConductivitykIns Heat conductivity of pipe insulation, used to compute R [W/(m.K)]
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
Booleanfrom_dpfalse= true, use m_flow = f(dp) else dp = f(m_flow)
Booleanhave_pipCaptrue= true, a mixing volume is added to each segment that corresponds to the heat capacity of the pipe segment wall
Booleanhave_symmetrytrue= false, the mixing volume is only on port_b of each segment, which improve performances, but reduces dynamic accuracy
MassFlowRatem_flow_small1E-4*abs(m_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
Booleanlinearizedfalse= true, use linear relation between m_flow and dp for any flow rate
Initialization
TemperatureT_start_in[nSeg]fill(Medium.T_default, nSeg)Initialization temperature at pipe inlet [K]
TemperatureT_start_out[nSeg]T_start_inInitialization temperature at pipe outlet [K]
BooleaninitDelayfalseInitialize delay for a constant mass flow rate if true, otherwise start from 0
MassFlowRatem_flow_start0Initial value of mass flow rate through pipe [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)
HeatPort_aheatPorts[nSeg]Heat transfer to or from surrounding for each pipe segment (positive if pipe is colder than surrounding)

Modelica definition

model PlugFlowPipeDiscretized "Discretized pipe model using spatialDistribution for temperature delay" extends Buildings.Fluid.Interfaces.PartialTwoPort; constant Boolean homotopyInitialization = true "= true, use homotopy method"; parameter Integer nSeg(min=1) = 1 "Number of axial segment"; parameter Boolean from_dp=false "= true, use m_flow = f(dp) else dp = f(m_flow)"; parameter Boolean have_pipCap=true "= true, a mixing volume is added to each segment that corresponds to the heat capacity of the pipe segment wall"; parameter Boolean have_symmetry=true "= false, the mixing volume is only on port_b of each segment, which improve performances, but reduces dynamic accuracy"; parameter Modelica.Units.SI.Length dh "Hydraulic diameter"; parameter Real ReC=4000 "Reynolds number where transition to turbulence starts"; parameter Modelica.Units.SI.Height roughness=2.5e-5 "Average height of surface asperities (default: smooth steel pipe)"; parameter Modelica.Units.SI.Length totLen=sum(segLen) "Total pipe length (used to compute segment length)"; parameter Modelica.Units.SI.Length segLen[nSeg]=fill(totLen/nSeg, nSeg) "Pipe segment length"; parameter Modelica.Units.SI.MassFlowRate m_flow_nominal "Nominal mass flow rate"; parameter Modelica.Units.SI.Length dIns "Thickness of pipe insulation, used to compute R"; parameter Modelica.Units.SI.ThermalConductivity kIns "Heat conductivity of pipe insulation, used to compute R"; parameter Modelica.Units.SI.Length thickness=0.0035 "Pipe wall thickness"; parameter Modelica.Units.SI.SpecificHeatCapacity cPip=2300 "Specific heat of pipe wall material. 2300 for PE, 500 for steel"; parameter Modelica.Units.SI.Density rhoPip(displayUnit="kg/m3") = 930 "Density of pipe wall material. 930 for PE, 8000 for steel"; parameter Modelica.Units.SI.Temperature T_start_in[nSeg]=fill(Medium.T_default, nSeg) "Initialization temperature at pipe inlet"; parameter Modelica.Units.SI.Temperature T_start_out[nSeg]=T_start_in "Initialization temperature at pipe outlet"; parameter Boolean initDelay = false "Initialize delay for a constant mass flow rate if true, otherwise start from 0"; parameter Modelica.Units.SI.MassFlowRate m_flow_start=0 "Initial value of mass flow rate through pipe"; parameter Real fac=1 "Factor to take into account flow resistance of bends etc., fac=dp_nominal/dpStraightPipe_nominal"; parameter Modelica.Units.SI.MassFlowRate m_flow_small=1E-4*abs(m_flow_nominal) "Small mass flow rate for regularization of zero flow"; parameter Boolean linearized = false "= true, use linear relation between m_flow and dp for any flow rate"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a heatPorts[nSeg] "Heat transfer to or from surrounding for each pipe segment (positive if pipe is colder than surrounding)"; Modelica.Units.SI.HeatFlowRate QEnv_flow=sum(heatPorts.Q_flow) "Heat transfer to or from surroundings of total pipe length (positive if pipe is colder than surrounding)"; Modelica.Units.SI.PressureDifference dp(displayUnit="Pa") = res.dp "Pressure difference between port_a and port_b"; Modelica.Units.SI.MassFlowRate m_flow=port_a.m_flow "Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction)"; final parameter Modelica.Units.SI.Velocity v_nominal=m_flow_nominal/(APip* rho_default) "Velocity at m_flow_nominal"; Modelica.Units.SI.Velocity v=pipSeg[1].v "Flow velocity of medium in pipe"; FixedResistances.HydraulicDiameter res( redeclare final package Medium = Medium, final m_flow_nominal=m_flow_nominal, final dh=dh, final from_dp=from_dp, final length=totLen, final roughness=roughness, final fac=fac, final ReC=ReC, final v_nominal=v_nominal, final allowFlowReversal=allowFlowReversal, final show_T=false, final homotopyInitialization=homotopyInitialization, final linearized=linearized, dp(nominal= if rho_default > 500 then totLen * fac * 200 else totLen * fac * 2)) "Pressure drop calculation for this pipe"; protected parameter Modelica.Units.SI.Length rInt=dh/2 "Pipe interior radius"; parameter Modelica.Units.SI.Area APip=Modelica.Constants.pi*rInt^2 "Pipe hydraulic cross-sectional area"; parameter Modelica.Units.SI.Density rho_default=Medium.density_pTX( p=Medium.p_default, T=Medium.T_default, X=Medium.X_default) "Default density (e.g., rho_liquidWater = 995, rho_air = 1.2)"; Buildings.Fluid.FixedResistances.BaseClasses.PlugFlowPipe pipSeg[nSeg]( redeclare each final package Medium = Medium, redeclare each final FixedResistances.LosslessPipe res, final length=segLen, final T_start_in=T_start_in, final T_start_out=T_start_out, each final dIns=dIns, each final kIns=kIns, each final cPip=cPip, each final rhoPip=rhoPip, each final dh=dh, each final v_nominal=v_nominal, each final allowFlowReversal=allowFlowReversal, each final m_flow_nominal=m_flow_nominal, each final thickness=thickness, each final m_flow_small=m_flow_small, each final initDelay=initDelay, each final have_pipCap=have_pipCap, each final have_symmetry=have_symmetry) "Pipe segments"; equation connect(res.port_b, port_b); for i in 2:nSeg loop connect(pipSeg[i-1].port_b, pipSeg[i].port_a); end for; connect(pipSeg.heatPort, heatPorts); connect(pipSeg[nSeg].port_b, res.port_a); connect(port_a, pipSeg[1].port_a); end PlugFlowPipeDiscretized;

Buildings.Fluid.FixedResistances.PressureDrop Buildings.Fluid.FixedResistances.PressureDrop

Fixed flow resistance with dp and m_flow as parameter

Buildings.Fluid.FixedResistances.PressureDrop

Information

Model of a flow resistance with a fixed flow coefficient. The mass flow rate is

ṁ = k √Δp,

where k is a constant and Δp is the pressure drop. The constant k is equal to k=m_flow_nominal/sqrt(dp_nominal), where m_flow_nominal and dp_nominal are parameters.

Assumptions

In the region abs(m_flow) < m_flow_turbulent, the square root is replaced by a differentiable function with finite slope. The value of m_flow_turbulent is computed as m_flow_turbulent = deltaM * abs(m_flow_nominal), where deltaM=0.3 and m_flow_nominal are parameters that can be set by the user.

The figure below shows the pressure drop for the parameters m_flow_nominal=5 kg/s, dp_nominal=10 Pa and deltaM=0.3.

image

Important parameters

The parameter from_dp is used to determine whether the mass flow rate is computed as a function of the pressure drop (if from_dp=true), or vice versa. This setting can affect the size of the nonlinear system of equations.

If the parameter linearized is set to true, then the pressure drop is computed as a linear function of the mass flow rate.

Setting allowFlowReversal=false can lead to simpler equations. However, this should only be set to false if one can guarantee that the flow never reverses its direction. This can be difficult to guarantee, as pressure imbalance after the initialization, or due to medium expansion and contraction, can lead to reverse flow.

If the parameter show_T is set to true, then the model will compute the temperature at its ports. Note that this can lead to state events when the mass flow rate approaches zero, which can increase computing time.

Notes

For more detailed models that compute the actual flow friction, models from the package Modelica.Fluid can be used and combined with models from the Buildings library.

For a model that uses the hydraulic parameter and flow velocity at nominal conditions as a parameter, use Buildings.Fluid.FixedResistances.HydraulicDiameter.

Implementation

The pressure drop is computed by calling a function in the package Buildings.Fluid.BaseClasses.FlowModels, This package contains regularized implementations of the equation

ṁ = sign(Δp) k √ Δp  

and its inverse function.

To decouple the energy equation from the mass equations, the pressure drop is a function of the mass flow rate, and not the volume flow rate. This leads to simpler equations.

Extends from Buildings.Fluid.BaseClasses.PartialResistance (Partial model for a hydraulic resistance).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
MassFlowRatem_flow_turbulentif computeFlowResistance the...Turbulent flow if |m_flow| >= m_flow_turbulent [kg/s]
Nominal condition
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
PressureDifferencedp_nominal Pressure drop at nominal mass flow rate [Pa]
Transition to laminar
RealdeltaM0.3Fraction of nominal mass flow rate where transition to turbulent occurs
Assumptions
BooleanallowFlowReversaltrue= false to simplify equations, assuming, but not enforcing, no flow reversal
Advanced
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
Booleanfrom_dpfalse= true, use m_flow = f(dp) else dp = f(m_flow)
Booleanlinearizedfalse= true, use linear relation between m_flow and dp for any flow rate

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)

Modelica definition

model PressureDrop "Fixed flow resistance with dp and m_flow as parameter" extends Buildings.Fluid.BaseClasses.PartialResistance( final m_flow_turbulent = if computeFlowResistance then deltaM * m_flow_nominal_pos else 0); parameter Real deltaM(min=1E-6) = 0.3 "Fraction of nominal mass flow rate where transition to turbulent occurs"; final parameter Real k = if computeFlowResistance then m_flow_nominal_pos / sqrt(dp_nominal_pos) else 0 "Flow coefficient, k=m_flow/sqrt(dp), with unit=(kg.m)^(1/2)"; protected final parameter Boolean computeFlowResistance=(dp_nominal_pos > Modelica.Constants.eps) "Flag to enable/disable computation of flow resistance"; final parameter Real coeff= if linearized and computeFlowResistance then if from_dp then k^2/m_flow_nominal_pos else m_flow_nominal_pos/k^2 else 0 "Precomputed coefficient to avoid division by parameter"; initial equation if computeFlowResistance then assert(m_flow_turbulent > 0, "m_flow_turbulent must be bigger than zero."); end if; assert(m_flow_nominal_pos > 0, "m_flow_nominal_pos must be non-zero. Check parameters."); equation // Pressure drop calculation if computeFlowResistance then if linearized then if from_dp then m_flow = dp*coeff; else dp = m_flow*coeff; end if; else if homotopyInitialization then if from_dp then m_flow=homotopy( actual=Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_dp( dp=dp, k=k, m_flow_turbulent=m_flow_turbulent), simplified=m_flow_nominal_pos*dp/dp_nominal_pos); else dp=homotopy( actual=Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_m_flow( m_flow=m_flow, k=k, m_flow_turbulent=m_flow_turbulent), simplified=dp_nominal_pos*m_flow/m_flow_nominal_pos); end if; // from_dp else // do not use homotopy if from_dp then m_flow=Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_dp( dp=dp, k=k, m_flow_turbulent=m_flow_turbulent); else dp=Buildings.Fluid.BaseClasses.FlowModels.basicFlowFunction_m_flow( m_flow=m_flow, k=k, m_flow_turbulent=m_flow_turbulent); end if; // from_dp end if; // homotopyInitialization end if; // linearized else // do not compute flow resistance dp = 0; end if; // computeFlowResistance end PressureDrop;