Buildings.Airflow.Multizone.BaseClasses

Package with base classes for Buildings.Airflow.Multizone

Information


This package contains base classes that are used to construct the models in Buildings.Airflow.Multizone.

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

Package Content

NameDescription
Buildings.Airflow.Multizone.BaseClasses.DoorDiscretized DoorDiscretized Door model using discretization along height coordinate
Buildings.Airflow.Multizone.BaseClasses.powerLaw powerLaw Power law used in orifice equations
Buildings.Airflow.Multizone.BaseClasses.PowerLawResistance PowerLawResistance Flow resistance that uses the power law
Buildings.Airflow.Multizone.BaseClasses.TwoWayFlowElement TwoWayFlowElement Flow resistance that uses the power law
Buildings.Airflow.Multizone.BaseClasses.TwoWayFlowElementBuoyancy TwoWayFlowElementBuoyancy Flow resistance that uses the power law
Buildings.Airflow.Multizone.BaseClasses.ZonalFlow ZonalFlow Flow across zonal boundaries of a room
ErrorControl Interface that defines parameters for error control


Buildings.Airflow.Multizone.BaseClasses.DoorDiscretized Buildings.Airflow.Multizone.BaseClasses.DoorDiscretized

Door model using discretization along height coordinate

Buildings.Airflow.Multizone.BaseClasses.DoorDiscretized

Information


This is a partial model for the bi-directional air flow through a door.

To compute the bi-directional flow, the door is discretize along the height coordinate, and uses an orifice equation to compute the flow for each compartment.

The compartment area dA is a variable, which allows using the model for a door that can be open or closed.

Extends from Buildings.Airflow.Multizone.BaseClasses.TwoWayFlowElementBuoyancy (Flow resistance that uses the power law).

Parameters

TypeNameDefaultDescription
BooleanforceErrorControlOnFlowtrueFlag to force error control on m_flow. Set to true if interested in flow rate
replaceable package MediumPartialMedium 
VelocityvZer0.001Minimum velocity to prevent zero flow. Recommended: 0.001 [m/s]
IntegernCom10Number of compartments for the discretization
Pressuredp_turbulent0.01Pressure difference where laminar and turbulent flow relation coincide. Recommended: 0.01 [Pa]
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]
Densityrho_a1_inflow.start1.2Density of air flowing in from port_a1 [kg/m3]
Densityrho_a2_inflow.start1.2Density of air flowing in from port_a2 [kg/m3]
Geometry
LengthwOpe0.9Width of opening [m]
LengthhOpe2.1Height of opening [m]
LengthhA2.7/2Height of reference pressure zone A [m]
LengthhB2.7/2Height of reference pressure zone B [m]
Orifice characteristics
RealCD0.65Discharge coefficient
Initialization
SpecificEnthalpyh_outflow_a1_startMedium1.h_defaultStart value for enthalpy flowing out of port a1 [J/kg]
SpecificEnthalpyh_outflow_b1_startMedium1.h_defaultStart value for enthalpy flowing out of port b1 [J/kg]
SpecificEnthalpyh_outflow_a2_startMedium2.h_defaultStart value for enthalpy flowing out of port a2 [J/kg]
SpecificEnthalpyh_outflow_b2_startMedium2.h_defaultStart value for enthalpy flowing out of port b2 [J/kg]
Advanced
MassFlowRatem1_flow_small1E-4*abs(m1_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
MassFlowRatem2_flow_small1E-4*abs(m2_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
BooleanhomotopyInitializationtrue= true, use homotopy method
Diagnostics
Booleanshow_V_flowfalse= true, if volume flow rate at inflowing port is computed
Booleanshow_Tfalse= true, if actual temperature at port is computed (may lead to events)

Connectors

TypeNameDescription
FluidPort_aport_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)

Modelica definition

partial model DoorDiscretized 
  "Door model using discretization along height coordinate"
  extends Buildings.Airflow.Multizone.BaseClasses.TwoWayFlowElementBuoyancy;

  parameter Integer nCom=10 "Number of compartments for the discretization";

  parameter Modelica.SIunits.Pressure dp_turbulent(min=0) = 0.01 
    "Pressure difference where laminar and turbulent flow relation coincide. Recommended: 0.01";
  parameter Real CD=0.65 "|Orifice characteristics|Discharge coefficient";

  Modelica.SIunits.Pressure dpAB[nCom](nominal=1) 
    "Pressure difference between compartments";
  Modelica.SIunits.Velocity v[nCom](nominal=0.01) 
    "Velocity in compartment from A to B";
  Modelica.SIunits.Velocity vTop "Velocity at top of opening from A to B";
  Modelica.SIunits.Velocity vBot "Velocity at bottom of opening from A to B";

protected 
  parameter Modelica.SIunits.Length dh=hOpe/nCom "Height of each compartment";
  Modelica.SIunits.AbsolutePressure pA[nCom](nominal=101325) 
    "Pressure in compartments of room A";
  Modelica.SIunits.AbsolutePressure pB[nCom](nominal=101325) 
    "Pressure in compartments of room B";

  Modelica.SIunits.VolumeFlowRate dV_flow[nCom] 
    "Volume flow rate through compartment from A to B";
  Modelica.SIunits.VolumeFlowRate dVAB_flow[nCom] 
    "Volume flow rate through compartment from A to B if positive";
  Modelica.SIunits.VolumeFlowRate dVBA_flow[nCom] 
    "Volume flow rate through compartment from B to A if positive";

  Real m(min=0.5, max=1) "Flow exponent, m=0.5 for turbulent, m=1 for laminar";
  Real kVal "Flow coefficient for each compartment, k = V_flow/ dp^m";
  Modelica.SIunits.Area dA "Compartment area";

  Modelica.SIunits.Density rhoAve "Average air density";

equation 
  dA = A/nCom;
  rhoAve = (rho_a1_inflow + rho_a2_inflow)/2;

  for i in 1:nCom loop
    // pressure drop in each compartment
    pA[i] = port_a1.p + rho_a1_inflow*Modelica.Constants.g_n*(hA - (i - 0.5)*dh);
    pB[i] = port_a2.p + rho_a2_inflow*Modelica.Constants.g_n*(hB - (i - 0.5)*dh);
    dpAB[i] = pA[i] - pB[i];
    // orifice equation
    dV_flow[i] = Buildings.Airflow.Multizone.BaseClasses.powerLaw(
      k=kVal,
      dp=dpAB[i],
      m=m,
      dp_turbulent=dp_turbulent);
    v[i] = dV_flow[i]/dA;
    // assignment of net volume flows
    dVAB_flow[i] = dV_flow[i]*
      Buildings.Utilities.Math.Functions.smoothHeaviside(x=dV_flow[i], delta=
      VZer_flow/nCom) + VZer_flow/nCom;
    dVBA_flow[i] = -dV_flow[i] + dVAB_flow[i] + 2*VZer_flow/nCom;
  end for;
  // add positive and negative flows
  VAB_flow = ones(nCom)*dVAB_flow;
  VBA_flow = ones(nCom)*dVBA_flow;
  vTop = v[nCom];
  vBot = v[1];
end DoorDiscretized;

Buildings.Airflow.Multizone.BaseClasses.powerLaw

Power law used in orifice equations

Information


This model describes the mass flow rate and pressure difference relation of an orifice in the form

    m_flow = k * dp^m,
where k is a variable and m a parameter. For turbulent flow, set m=1/2 and for laminar flow, set m=1.

The model is used as a base for the interzonal air flow models.

Inputs

TypeNameDefaultDescription
Realk Flow coefficient, k = V_flow/ dp^m
Pressuredp Pressure difference [Pa]
Realm Flow exponent, m=0.5 for turbulent, m=1 for laminar
Pressuredp_turbulent0.001Pressure difference where laminar and turbulent flow relation coincide [Pa]

Outputs

TypeNameDescription
VolumeFlowRateV_flowVolume flow rate [m3/s]

Modelica definition

function powerLaw "Power law used in orifice equations"
  input Real k "Flow coefficient, k = V_flow/ dp^m";
  input Modelica.SIunits.Pressure dp "Pressure difference";
  input Real m(min=0.5, max=1) 
    "Flow exponent, m=0.5 for turbulent, m=1 for laminar";
  input Modelica.SIunits.Pressure dp_turbulent(min=0)=0.001 
    "Pressure difference where laminar and turbulent flow relation coincide";

  output Modelica.SIunits.VolumeFlowRate V_flow "Volume flow rate";
protected 
  Modelica.SIunits.Pressure delP "Half-width of transition interval";
  Modelica.SIunits.Pressure pTilFor 
    "Pressure at which transition occurs at forward flow";
  Modelica.SIunits.Pressure pTilRev 
    "Pressure at which transition occurs at reverse flow";

 Real s "Part of the linear slope that is independent of k";

algorithm 
  s :=dp_turbulent^(m - 1);
                    // m can be time dependent, for example for the operable door

  delP :=dp_turbulent/2.0;
  pTilFor :=dp - dp_turbulent;
  pTilRev :=dp + dp_turbulent;

  /*V_flow :=k*spliceFunction(
      spliceFunction(
        abs(dp)^m,
        s*dp,
        pTilFor,
        delP),
      -(abs(dp)^m),
      pTilRev,
      delP);
*/
 if (dp >= dp_turbulent) then
    V_flow :=k*dp^m;
  elseif (dp <= -dp_turbulent) then
    V_flow :=-k*(-dp)^m;
 else
   // The test below avoids computing 0^(-0.5) in
   // Buildings.Airflow.Multizone.BaseClasses.powerLaw:derf,
   // which causes the following error:
   // Model error - power: (abs(dp)) ** (m-1) = (0) ** (-0.5)
   if (abs(dp)<Modelica.Constants.small) then
      V_flow := 0;
   else
     V_flow :=k*spliceFunction(
      spliceFunction(
        abs(dp)^m,
        s*dp,
        pTilFor,
        delP),
      -(abs(dp)^m),
      pTilRev,
      delP);
   end if;

  end if;

end powerLaw;

Buildings.Airflow.Multizone.BaseClasses.PowerLawResistance Buildings.Airflow.Multizone.BaseClasses.PowerLawResistance

Flow resistance that uses the power law

Buildings.Airflow.Multizone.BaseClasses.PowerLawResistance

Information


This model describes the mass flow rate and pressure difference relation of an orifice in the form

    V_flow = k * dp^m,
where k is a variable and m a parameter. For turbulent flow, set m=1/2 and for laminar flow, set m=1.

The model is used as a base for the interzonal air flow models.

Extends from Buildings.Fluid.Interfaces.PartialStaticTwoPortInterface (Partial model transporting fluid between two ports without storing mass or energy), Buildings.Airflow.Multizone.BaseClasses.ErrorControl (Interface that defines parameters for error control).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
BooleanforceErrorControlOnFlowtrueFlag to force error control on m_flow. Set to true if interested in flow rate
Realm Flow exponent, m=0.5 for turbulent, m=1 for laminar
BooleanuseConstantDensitytrueSet to false to use density based on state (as implemented by the Medium model)
Pressuredp_turbulent0.1Pressure difference where laminar and turbulent flow relation coincide. Recommended = 0.1 [Pa]
LengthlWetsqrt(A)Wetted perimeter used for Reynolds number calculation [m]
Nominal condition
MassFlowRatem_flow_nominalrho_nominal*k*dp_turbulentNominal mass flow rate [kg/s]
Initialization
MassFlowRatem_flow.start0Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction) [kg/s]
Pressuredp.start0Pressure difference between port_a and port_b [Pa]
Orifice characteristics
AreaA Area of orifice [m2]
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_small1E-4*abs(m_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
BooleanhomotopyInitializationtrue= true, use homotopy method
Diagnostics
Booleanshow_V_flowtrue= true, if volume flow rate at inflowing port is computed
Booleanshow_Tfalse= true, if actual temperature at port is computed (may lead to events)

Connectors

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

Modelica definition

partial model PowerLawResistance 
  "Flow resistance that uses the power law"
  extends Buildings.Fluid.Interfaces.PartialStaticTwoPortInterface(final 
      m_flow_nominal=rho_nominal*k*dp_turbulent, final show_V_flow=true);
  extends Buildings.Airflow.Multizone.BaseClasses.ErrorControl;

  parameter Modelica.SIunits.Area A "|Orifice characteristics|Area of orifice";

  parameter Real m(min=0.5, max=1) 
    "Flow exponent, m=0.5 for turbulent, m=1 for laminar";
  parameter Boolean useConstantDensity=true 
    "Set to false to use density based on state (as implemented by the Medium model)";
  Modelica.SIunits.Density rho "Fluid density at port_a";
  parameter Modelica.SIunits.Pressure dp_turbulent(min=0) = 0.1 
    "Pressure difference where laminar and turbulent flow relation coincide. Recommended = 0.1";
  Modelica.SIunits.Velocity v(nominal=1) "Average velocity";
  parameter Modelica.SIunits.Length lWet=sqrt(A) 
    "Wetted perimeter used for Reynolds number calculation";
  Real Re "Reynolds number";

protected 
  parameter Real k(fixed=false) "Flow coefficient, k = V_flow/ dp^m";

  parameter Medium.ThermodynamicState sta0=Medium.setState_pTX(
      T=Medium.T_default,
      p=Medium.p_default,
      X=Medium.X_default);
  parameter Modelica.SIunits.Density rho_nominal=Medium.density(sta0) 
    "Density, used to compute fluid volume";
  Modelica.SIunits.DynamicViscosity dynVis "Dynamic viscosity";

  Modelica.SIunits.Mass mExc(start=0) 
    "Air mass exchanged (for purpose of error control only)";

equation 
  if forceErrorControlOnFlow then
    der(mExc) = port_a.m_flow;
  else
    der(mExc) = 0;
  end if;

  rho = if useConstantDensity then rho_nominal else Medium.density(sta_a);
  dynVis = Medium.dynamicViscosity(sta_a);
  port_a.m_flow = rho*Buildings.Airflow.Multizone.BaseClasses.powerLaw(
    k=k,
    dp=dp,
    m=m,
    dp_turbulent=dp_turbulent);
  v = V_flow/A;
  Re = v*lWet*rho/dynVis;

  // Isenthalpic state transformation (no storage and no loss of energy)
  port_a.h_outflow = inStream(port_b.h_outflow);
  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 = inStream(port_b.Xi_outflow);
  port_b.Xi_outflow = inStream(port_a.Xi_outflow);

  port_a.C_outflow = inStream(port_b.C_outflow);
  port_b.C_outflow = inStream(port_a.C_outflow);
end PowerLawResistance;

Buildings.Airflow.Multizone.BaseClasses.TwoWayFlowElement Buildings.Airflow.Multizone.BaseClasses.TwoWayFlowElement

Flow resistance that uses the power law

Buildings.Airflow.Multizone.BaseClasses.TwoWayFlowElement

Information


This is a partial model for models that describe the bi-directional air flow through large openings.

Models that extend this model need to compute mAB_flow and mBA_flow, or alternatively VAB_flow and VBA_flow, and the face area area. The face area is a variable to allow this partial model to be used for doors that can be open or closed as a function of an input signal.

Extends from Buildings.Fluid.Interfaces.PartialStaticFourPortInterface (Partial model transporting fluid between two ports without storing mass or energy), Buildings.Airflow.Multizone.BaseClasses.ErrorControl (Interface that defines parameters for error control).

Parameters

TypeNameDefaultDescription
replaceable package Medium1PartialMediumMedium 1 in the component
replaceable package Medium2PartialMediumMedium 2 in the component
BooleanforceErrorControlOnFlowtrueFlag to force error control on m_flow. Set to true if interested in flow rate
VelocityvZer0.001Minimum velocity to prevent zero flow. Recommended: 0.001 [m/s]
Nominal condition
MassFlowRatem1_flow_nominal10/3600*1.2Nominal mass flow rate [kg/s]
MassFlowRatem2_flow_nominalm1_flow_nominalNominal mass flow rate [kg/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
BooleanallowFlowReversal1false= true to allow flow reversal in medium 1, false restricts to design direction (port_a -> port_b)
BooleanallowFlowReversal2false= true to allow flow reversal in medium 2, false restricts to design direction (port_a -> port_b)
Initialization
SpecificEnthalpyh_outflow_a1_startMedium1.h_defaultStart value for enthalpy flowing out of port a1 [J/kg]
SpecificEnthalpyh_outflow_b1_startMedium1.h_defaultStart value for enthalpy flowing out of port b1 [J/kg]
SpecificEnthalpyh_outflow_a2_startMedium2.h_defaultStart value for enthalpy flowing out of port a2 [J/kg]
SpecificEnthalpyh_outflow_b2_startMedium2.h_defaultStart value for enthalpy flowing out of port b2 [J/kg]
Advanced
MassFlowRatem1_flow_small1E-4*abs(m1_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
MassFlowRatem2_flow_small1E-4*abs(m2_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
BooleanhomotopyInitializationtrue= true, use homotopy method
Diagnostics
Booleanshow_V_flowfalse= true, if volume flow rate at inflowing port is computed
Booleanshow_Tfalse= true, if actual temperature at port is computed (may lead to events)

Connectors

TypeNameDescription
replaceable package Medium1Medium 1 in the component
replaceable package Medium2Medium 2 in the component
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)

Modelica definition

partial model TwoWayFlowElement 
  "Flow resistance that uses the power law"
  extends Buildings.Fluid.Interfaces.PartialStaticFourPortInterface(
    redeclare final package Medium1 = Medium,
    redeclare final package Medium2 = Medium,
    final allowFlowReversal1=false,
    final allowFlowReversal2=false,
    final m1_flow_nominal=10/3600*1.2,
    final m2_flow_nominal=m1_flow_nominal);
  extends Buildings.Airflow.Multizone.BaseClasses.ErrorControl;

  replaceable package Medium = Modelica.Media.Interfaces.PartialMedium;

  Modelica.SIunits.VolumeFlowRate VAB_flow(nominal=0.001) 
    "Volume flow rate from A to B if positive";
  Modelica.SIunits.VolumeFlowRate VBA_flow(nominal=0.001) 
    "Volume flow rate from B to A if positive";
  Modelica.SIunits.MassFlowRate mAB_flow(nominal=0.001) 
    "Mass flow rate from A to B if positive";
  Modelica.SIunits.MassFlowRate mBA_flow(nominal=0.001) 
    "Mass flow rate from B to A if positive";

  Modelica.SIunits.Velocity vAB(nominal=0.01) "Average velocity from A to B";
  Modelica.SIunits.Velocity vBA(nominal=0.01) "Average velocity from B to A";

  Modelica.SIunits.Density rho_a1_inflow(start=1.2) 
    "Density of air flowing in from port_a1";
  Modelica.SIunits.Density rho_a2_inflow(start=1.2) 
    "Density of air flowing in from port_a2";

  Modelica.SIunits.Area A "Face area";

  parameter Modelica.SIunits.Velocity vZer=0.001 
    "Minimum velocity to prevent zero flow. Recommended: 0.001";
protected 
  Modelica.SIunits.VolumeFlowRate VZer_flow(fixed=false) 
    "Minimum net volume flow rate to prevent zero flow";
protected 
  Modelica.SIunits.Mass mExcAB(start=0) 
    "Air mass exchanged (for purpose of error control only)";
  Modelica.SIunits.Mass mExcBA(start=0) 
    "Air mass exchanged (for purpose of error control only)";

equation 
  // enforcing error control on both direction rather than on the sum only
  // gives higher robustness. The reason may be that for bi-directional flow,
  // (VAB_flow - VBA_flow) may be close to zero.
  if forceErrorControlOnFlow then
    der(mExcAB) = mAB_flow;
    der(mExcBA) = mBA_flow;
  else
    der(mExcAB) = 0;
    der(mExcBA) = 0;
  end if;
  rho_a1_inflow = Medium.density(state_a1_inflow);
  rho_a2_inflow = Medium.density(state_a2_inflow);
  VZer_flow = vZer*A;

  mAB_flow = rho_a1_inflow*VAB_flow;
  mBA_flow = rho_a2_inflow*VBA_flow;
  // Average velocity (using the whole orifice area)
  vAB = VAB_flow/A;
  vBA = VBA_flow/A;

  port_a1.m_flow = mAB_flow;
  port_a2.m_flow = mBA_flow;

  // Energy balance (no storage, no heat loss/gain)
  port_a1.h_outflow = inStream(port_b1.h_outflow);
  port_b1.h_outflow = inStream(port_a1.h_outflow);
  port_a2.h_outflow = inStream(port_b2.h_outflow);
  port_b2.h_outflow = inStream(port_a2.h_outflow);

  // Mass balance (no storage)
  port_a1.m_flow = -port_b1.m_flow;
  port_a2.m_flow = -port_b2.m_flow;

  port_a1.Xi_outflow = inStream(port_b1.Xi_outflow);
  port_b1.Xi_outflow = inStream(port_a1.Xi_outflow);
  port_a2.Xi_outflow = inStream(port_b2.Xi_outflow);
  port_b2.Xi_outflow = inStream(port_a2.Xi_outflow);

  // Transport of trace substances
  port_a1.C_outflow = inStream(port_b1.C_outflow);
  port_b1.C_outflow = inStream(port_a1.C_outflow);
  port_a2.C_outflow = inStream(port_b2.C_outflow);
  port_b2.C_outflow = inStream(port_a2.C_outflow);

end TwoWayFlowElement;

Buildings.Airflow.Multizone.BaseClasses.TwoWayFlowElementBuoyancy Buildings.Airflow.Multizone.BaseClasses.TwoWayFlowElementBuoyancy

Flow resistance that uses the power law

Buildings.Airflow.Multizone.BaseClasses.TwoWayFlowElementBuoyancy

Information


This is a partial model for models that describe the bi-directional air flow through large openings.

Models that extend this model need to compute mAB_flow and mBA_flow, or alternatively VAB_flow and VBA_flow, and the face area area. The face area is a variable to allow this partial model to be used for doors that can be open or closed as a function of an input signal.

Extends from Buildings.Airflow.Multizone.BaseClasses.TwoWayFlowElement (Flow resistance that uses the power law).

Parameters

TypeNameDefaultDescription
BooleanforceErrorControlOnFlowtrueFlag to force error control on m_flow. Set to true if interested in flow rate
replaceable package MediumPartialMedium 
VelocityvZer0.001Minimum velocity to prevent zero flow. Recommended: 0.001 [m/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]
Densityrho_a1_inflow.start1.2Density of air flowing in from port_a1 [kg/m3]
Densityrho_a2_inflow.start1.2Density of air flowing in from port_a2 [kg/m3]
Geometry
LengthwOpe0.9Width of opening [m]
LengthhOpe2.1Height of opening [m]
LengthhA2.7/2Height of reference pressure zone A [m]
LengthhB2.7/2Height of reference pressure zone B [m]
Initialization
SpecificEnthalpyh_outflow_a1_startMedium1.h_defaultStart value for enthalpy flowing out of port a1 [J/kg]
SpecificEnthalpyh_outflow_b1_startMedium1.h_defaultStart value for enthalpy flowing out of port b1 [J/kg]
SpecificEnthalpyh_outflow_a2_startMedium2.h_defaultStart value for enthalpy flowing out of port a2 [J/kg]
SpecificEnthalpyh_outflow_b2_startMedium2.h_defaultStart value for enthalpy flowing out of port b2 [J/kg]
Advanced
MassFlowRatem1_flow_small1E-4*abs(m1_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
MassFlowRatem2_flow_small1E-4*abs(m2_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
BooleanhomotopyInitializationtrue= true, use homotopy method
Diagnostics
Booleanshow_V_flowfalse= true, if volume flow rate at inflowing port is computed
Booleanshow_Tfalse= true, if actual temperature at port is computed (may lead to events)

Connectors

TypeNameDescription
FluidPort_aport_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)

Modelica definition

partial model TwoWayFlowElementBuoyancy 
  "Flow resistance that uses the power law"
  extends Buildings.Airflow.Multizone.BaseClasses.TwoWayFlowElement;

  parameter Modelica.SIunits.Length wOpe=0.9 "|Geometry|Width of opening";
  parameter Modelica.SIunits.Length hOpe=2.1 "|Geometry|Height of opening";

  parameter Modelica.SIunits.Length hA=2.7/2 
    "|Geometry|Height of reference pressure zone A";
  parameter Modelica.SIunits.Length hB=2.7/2 
    "|Geometry|Height of reference pressure zone B";

end TwoWayFlowElementBuoyancy;

Buildings.Airflow.Multizone.BaseClasses.ZonalFlow Buildings.Airflow.Multizone.BaseClasses.ZonalFlow

Flow across zonal boundaries of a room

Buildings.Airflow.Multizone.BaseClasses.ZonalFlow

Information


This is a partial model for computing the air exchange between volumes. Models that extend this model need to provide an equation for port_a1.m_flow and port_a2.m_flow.

Extends from Buildings.Fluid.Interfaces.PartialStaticFourPortInterface (Partial model transporting fluid between two ports without storing mass or energy), Buildings.Airflow.Multizone.BaseClasses.ErrorControl (Interface that defines parameters for error control).

Parameters

TypeNameDefaultDescription
replaceable package Medium1PartialMediumMedium 1 in the component
replaceable package Medium2PartialMediumMedium 2 in the component
BooleanforceErrorControlOnFlowtrueFlag to force error control on m_flow. Set to true if interested in flow rate
Nominal condition
MassFlowRatem1_flow_nominal10/3600*1.2Nominal mass flow rate [kg/s]
MassFlowRatem2_flow_nominalm1_flow_nominalNominal mass flow rate [kg/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
BooleanallowFlowReversal1false= true to allow flow reversal in medium 1, false restricts to design direction (port_a -> port_b)
BooleanallowFlowReversal2false= true to allow flow reversal in medium 2, false restricts to design direction (port_a -> port_b)
Initialization
SpecificEnthalpyh_outflow_a1_startMedium1.h_defaultStart value for enthalpy flowing out of port a1 [J/kg]
SpecificEnthalpyh_outflow_b1_startMedium1.h_defaultStart value for enthalpy flowing out of port b1 [J/kg]
SpecificEnthalpyh_outflow_a2_startMedium2.h_defaultStart value for enthalpy flowing out of port a2 [J/kg]
SpecificEnthalpyh_outflow_b2_startMedium2.h_defaultStart value for enthalpy flowing out of port b2 [J/kg]
Advanced
MassFlowRatem1_flow_small1E-4*abs(m1_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
MassFlowRatem2_flow_small1E-4*abs(m2_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
BooleanhomotopyInitializationtrue= true, use homotopy method
Diagnostics
Booleanshow_V_flowfalse= true, if volume flow rate at inflowing port is computed
Booleanshow_Tfalse= true, if actual temperature at port is computed (may lead to events)

Connectors

TypeNameDescription
replaceable package Medium1Medium 1 in the component
replaceable package Medium2Medium 2 in the component
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)

Modelica definition

partial model ZonalFlow "Flow across zonal boundaries of a room"
  extends Buildings.Fluid.Interfaces.PartialStaticFourPortInterface(
     redeclare final package Medium1 = Medium,
     redeclare final package Medium2 = Medium,
     final allowFlowReversal1 = false,
     final allowFlowReversal2 = false,
     final m1_flow_nominal = 10/3600*1.2,
     final m2_flow_nominal = m1_flow_nominal);
  extends Buildings.Airflow.Multizone.BaseClasses.ErrorControl;

   replaceable package Medium = Modelica.Media.Interfaces.PartialMedium;

equation 
  // Energy balance (no storage, no heat loss/gain)
  port_a1.h_outflow = inStream(port_b1.h_outflow);
  port_b1.h_outflow = inStream(port_a1.h_outflow);
  port_a2.h_outflow = inStream(port_b2.h_outflow);
  port_b2.h_outflow = inStream(port_a2.h_outflow);

  // Mass balance (no storage)
  port_a1.m_flow = -port_b1.m_flow;
  port_a2.m_flow = -port_b2.m_flow;

  port_a1.Xi_outflow = inStream(port_b1.Xi_outflow);
  port_b1.Xi_outflow = inStream(port_a1.Xi_outflow);
  port_a2.Xi_outflow = inStream(port_b2.Xi_outflow);
  port_b2.Xi_outflow = inStream(port_a2.Xi_outflow);

  // Transport of trace substances
  port_a1.C_outflow = inStream(port_b1.C_outflow);
  port_b1.C_outflow = inStream(port_a1.C_outflow);
  port_a2.C_outflow = inStream(port_b2.C_outflow);
  port_b2.C_outflow = inStream(port_a2.C_outflow);

end ZonalFlow;

Buildings.Airflow.Multizone.BaseClasses.ErrorControl

Interface that defines parameters for error control

Information


This is an interface that defines parameters used for error control.

Dymola does error control on state variables, such as temperature, pressure and species concentration. Flow variables such as m_flow are typically not checked during the error control. This can give large errors in flow variables, as long as the error on the volume's state variables that are coupled to the flow variables is small. Obtaining accurate flow variables can be achieved by imposing an error control on the exchanged mass, which can be defined as

  dm/dt = m_flow.
By setting enforceErrorControlOnFlow = true, such an equation is imposed by models that extend this class.

Parameters

TypeNameDefaultDescription
BooleanforceErrorControlOnFlowtrueFlag to force error control on m_flow. Set to true if interested in flow rate

Modelica definition

model ErrorControl 
  "Interface that defines parameters for error control"
  parameter Boolean forceErrorControlOnFlow = true 
    "Flag to force error control on m_flow. Set to true if interested in flow rate";
end ErrorControl;

Automatically generated Fri May 06 14:10:58 2011.