 
This package contains base classes that are used to construct the models in Buildings.Rooms.
Extends from Modelica.Icons.BasesPackage (Icon for packages containing base classes).
| Name | Description | 
|---|---|
|  ExteriorBoundaryConditions | Model for convection and radiation bounary condition of exterior constructions | 
|  ExteriorBoundaryConditionsWithWindow | Model for exterior boundary conditions for constructions with a window | 
|  HeatGain | Model to convert internal heat gain signals | 
|  InfraredRadiationExchange | Infrared radiation heat exchange between the room facing surfaces | 
|  InfraredRadiationGainDistribution | Infrared radiative heat gain distribution between the room facing surfaces | 
|  MixedAir | Model for room air that is completely mixed | 
|  RadiationAdapter | Model to connect between signals and heat port for radiative gains of the room | 
|  RadiationTemperature | Radiative temperature of the room | 
|  SolarRadiationExchange | Solar radiation heat exchange between the room facing surfaces | 
|  PartialSurfaceInterface | Partial model that is used for infrared radiation balance | 
|  SkyRadiationExchange | Radiative heat exchange with the sky and the ambient | 
|  ConstructionNumbers | Data records for construction data | 
|  ConstructionRecords | Data records for construction data | 
|  ParameterConstructionWithWindow | Record for exterior constructions that have a window | 
|  ParameterConstruction | Record for exterior constructions that have no window | 
|  PartialParameterConstruction | Partial record for constructions | 
|  Examples | Collection of models that illustrate model use and test models | 
 Buildings.Rooms.BaseClasses.ExteriorBoundaryConditions
Buildings.Rooms.BaseClasses.ExteriorBoundaryConditions
 
The model computes the infrared, solar, and convective heat exchange between these surfaces and the exterior temperature and the sky temperature. Input into this model are weather data that may be obtained from Buildings.BoundaryConditions.WeatherData.
In this model, the solar radiation data are converted from horizontal irradiation to irradiation on tilted surfaces using models from the package Buildings.BoundaryConditions.SolarIrradiation. The convective heat transfer between the exterior surface of the opaque constructions is computed using Buildings.HeatTransfer.Convection.
The heat transfer of windows are not computed in this model. They are implemented in Buildings.Rooms.BaseClasses.ExteriorBoundaryConditionsWithWindow.
| Type | Name | Default | Description | 
|---|---|---|---|
| Angle | lat | Latitude [rad] | |
| Angle | til[nCon] | Surface tilt (0 if the surface is a roof) [rad] | |
| Angle | azi[nCon] | Surface azimuth [rad] | |
| Area | AOpa[nCon] | Areas of exterior constructions (excluding the window area) [m2] | |
| Boolean | linearizeRadiation | Set to true to linearize emissive power | |
| Emissivity | absIR[nCon] | Infrared absorptivity of building surface [1] | |
| Emissivity | absSol[nCon] | Solar absorptivity of building surface [1] | |
| Exterior constructions | |||
| Integer | nCon | Number of exterior constructions | |
| Convective heat transfer | |||
| ExteriorConvection | conMod | Buildings.HeatTransfer.Types... | Convective heat transfer model for opaque part of the constructions | 
| CoefficientOfHeatTransfer | hFixed | 10.0 | Constant convection coefficient for opaque part of the constructions [W/(m2.K)] | 
| Type | Name | Description | 
|---|---|---|
| HeatPort_a | opa_a[nCon] | Heat port at surface a of opaque construction | 
| Bus | weaBus | 
model ExteriorBoundaryConditions 
  "Model for convection and radiation bounary condition of exterior constructions"
  parameter Integer nCon(min=1) "Number of exterior constructions";
  parameter Modelica.SIunits.Angle lat "Latitude";
  parameter Modelica.SIunits.Angle til[nCon] 
    "Surface tilt (0 if the surface is a roof)";
  parameter Modelica.SIunits.Angle azi[nCon] "Surface azimuth";
  parameter Modelica.SIunits.Area AOpa[nCon] 
    "Areas of exterior constructions (excluding the window area)";
  parameter Boolean linearizeRadiation 
    "Set to true to linearize emissive power";
  parameter Modelica.SIunits.Emissivity absIR[nCon] 
    "Infrared absorptivity of building surface";
  parameter Modelica.SIunits.Emissivity absSol[nCon] 
    "Solar absorptivity of building surface";
  Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a opa_a[nCon] 
    "Heat port at surface a of opaque construction"; 
  parameter Buildings.HeatTransfer.Types.ExteriorConvection conMod=
  Buildings.HeatTransfer.Types.ExteriorConvection.TemperatureWind 
    "Convective heat transfer model for opaque part of the constructions";
  parameter Modelica.SIunits.CoefficientOfHeatTransfer hFixed=10.0 
    "Constant convection coefficient for opaque part of the constructions";
  // The convection coefficients are not final to allow a user to individually
  // assign them.
  // We reassign the tilt since a roof has been declared in the room model as the
  // ceiling (of the room)
  HeatTransfer.Convection.Exterior conOpa[nCon](
    final A=AOpa,
    final til=Modelica.Constants.pi*ones(nCon) - til,
    final azi=azi,
    each conMod=conMod,
    each hFixed=hFixed) "Convection model for opaque part of the wall"; 
  SkyRadiationExchange skyRadExc(
    final n=nCon,
    each final A=AOpa,
    each final absIR=absIR,
    vieFacSky={(Modelica.Constants.pi - til[i])./Modelica.Constants.pi for i in 1:nCon}) 
    "Infrared radiative heat exchange with sky";
  BoundaryConditions.WeatherData.Bus weaBus; 
  BoundaryConditions.SolarIrradiation.DirectTiltedSurface HDirTil[
            nCon](
    each final lat=lat,
    final til=til,
    final azi=azi) "Direct solar irradiation on the surface";
  BoundaryConditions.SolarIrradiation.DiffusePerez HDifTil[nCon](
    each final lat=lat,
    final til=til,
    final azi=azi) "Diffuse solar irradiation";
  Modelica.Blocks.Math.Add HTotConExt[nCon](
    final k1=absSol .* AOpa,
    final k2=absSol .* AOpa) "Total solar irradiation";
  Buildings.HeatTransfer.Sources.PrescribedHeatFlow solHeaGaiConExt[nCon] 
    "Total solar heat gain of the surface"; 
protected 
  Buildings.HeatTransfer.Sources.PrescribedTemperature TAirConExt[
    nCon] "Outside air temperature for exterior constructions";
  Modelica.Blocks.Routing.Replicator repConExt(nout=nCon) "Signal replicator"; 
  Modelica.Blocks.Routing.Replicator repConExt1(
                                               nout=nCon) "Signal replicator";
  Modelica.Blocks.Routing.Replicator repConExt2(
                                               nout=nCon) "Signal replicator"; 
equation 
  connect(conOpa.solid, opa_a);
  connect(skyRadExc.port, opa_a); 
  connect(TAirConExt.port, conOpa.fluid);
  connect(repConExt.y, TAirConExt.T);
  connect(repConExt.u, weaBus.TDryBul);
  connect(skyRadExc.TOut, weaBus.TDryBul);
  connect(skyRadExc.TBlaSky, weaBus.TBlaSky); 
  for i in 1:nCon loop
  connect(weaBus, HDirTil[i].weaBus);
  connect(HDifTil[i].weaBus, weaBus); 
   end for;
  connect(HTotConExt.y, solHeaGaiConExt.Q_flow);
  connect(solHeaGaiConExt.port, opa_a);
  connect(HDirTil.H, HTotConExt.u1);
  connect(HDifTil.H, HTotConExt.u2);
  connect(repConExt2.u, weaBus.winDir);
  connect(repConExt1.u, weaBus.winSpe);
  connect(repConExt1.y, conOpa.v);
  connect(repConExt2.y, conOpa.dir); 
end ExteriorBoundaryConditions;
 
 Buildings.Rooms.BaseClasses.ExteriorBoundaryConditionsWithWindow
Buildings.Rooms.BaseClasses.ExteriorBoundaryConditionsWithWindow
 
The model computes the infrared, solar, and convective heat exchange between these surfaces and the exterior temperature and the sky temperature. Input into this model are weather data that may be obtained from Buildings.BoundaryConditions.WeatherData.
This model extends Buildings.Rooms.BaseClasses.ExteriorBoundaryConditions, which models the boundary conditions for the opaque constructions, and then implements the boundary condition for windows by using the model Buildings.HeatTransfer.Windows.ExteriorHeatTransfer.
Extends from Buildings.Rooms.BaseClasses.ExteriorBoundaryConditions (Model for convection and radiation bounary condition of exterior constructions).
| Type | Name | Default | Description | 
|---|---|---|---|
| Angle | lat | Latitude [rad] | |
| Angle | til[nCon] | Surface tilt (0 if the surface is a roof) [rad] | |
| Angle | azi[nCon] | Surface azimuth [rad] | |
| Area | AOpa[nCon] | Areas of exterior constructions (excluding the window area) [m2] | |
| Boolean | linearizeRadiation | Set to true to linearize emissive power | |
| Emissivity | absIR[nCon] | Infrared absorptivity of building surface [1] | |
| Emissivity | absSol[nCon] | Solar absorptivity of building surface [1] | |
| Area | AWin[nCon] | Areas of window glass and frame [m2] | |
| Exterior constructions | |||
| Integer | nCon | Number of exterior constructions | |
| Convective heat transfer | |||
| ExteriorConvection | conMod | Buildings.HeatTransfer.Types... | Convective heat transfer model for opaque part of the constructions | 
| CoefficientOfHeatTransfer | hFixed | 10.0 | Constant convection coefficient for opaque part of the constructions [W/(m2.K)] | 
| Frame | |||
| Real | fFra[nCon] | Fraction of window frame divided by total window area | |
| Emissivity | absSolFra[nCon] | Solar absorptivity of window frame [1] | |
| Emissivity | absIRFra[nCon] | Infrared absorptivity of window frame [1] | |
| Shading | |||
| Emissivity | absIRSha_air[nCon] | Infrared absorptivity of shade surface that faces air [1] | |
| Emissivity | absIRSha_glass[nCon] | Infrared absorptivity of shade surface that faces glass [1] | |
| TransmissionCoefficient | tauIRSha_air[nCon] | Infrared transmissivity of shade for radiation coming from the exterior or the room [1] | |
| TransmissionCoefficient | tauIRSha_glass[nCon] | Infrared transmissivity of shade for radiation coming from the glass [1] | |
| Boolean | haveExteriorShade[nCon] | Set to true if window has exterior shade (at surface a) | |
| Boolean | haveInteriorShade[nCon] | Set to true if window has interior shade (at surface b) | |
| Type | Name | Description | 
|---|---|---|
| HeatPort_a | opa_a[nCon] | Heat port at surface a of opaque construction | 
| Bus | weaBus | |
| input RealInput | uSha[nCon] | Control signal for the shading device, 0: unshaded; 1: fully shaded | 
| input RealInput | QAbsSolSha_flow[nCon] | Solar radiation absorbed by shade [W] | 
| output RadiosityOutflow | JOutUns[nCon] | Outgoing radiosity that connects to unshaded part of glass at exterior side [W] | 
| input RadiosityInflow | JInUns[nCon] | Incoming radiosity that connects to unshaded part of glass at exterior side [W] | 
| output RadiosityOutflow | JOutSha[nCon] | Outgoing radiosity that connects to shaded part of glass at exterior side [W] | 
| input RadiosityInflow | JInSha[nCon] | Incoming radiosity that connects to shaded part of glass at exterior side [W] | 
| HeatPort_a | glaUns[nCon] | Heat port at unshaded glass of exterior-facing surface | 
| HeatPort_a | glaSha[nCon] | Heat port at shaded glass of exterior-facing surface | 
| HeatPort_a | fra[nCon] | Heat port at frame of exterior-facing surface | 
| output RealOutput | HDir[nCon] | Direct solar irradition on tilted surface [W/m2] | 
| output RealOutput | HDif[nCon] | Diffuse solar irradiation on tilted surface [W/m2] | 
| output RealOutput | inc[nCon] | Incidence angle [rad] | 
model ExteriorBoundaryConditionsWithWindow 
  "Model for exterior boundary conditions for constructions with a window"
  extends Buildings.Rooms.BaseClasses.ExteriorBoundaryConditions;
  parameter Modelica.SIunits.Area AWin[nCon] "Areas of window glass and frame";
  parameter Real fFra[nCon](each min=0, each max=1) 
    "Fraction of window frame divided by total window area";
  parameter Modelica.SIunits.Emissivity absSolFra[nCon] 
    "Solar absorptivity of window frame";
  parameter Modelica.SIunits.Emissivity absIRFra[nCon] 
    "Infrared absorptivity of window frame";
  parameter Modelica.SIunits.Emissivity absIRSha_air[nCon] 
    "Infrared absorptivity of shade surface that faces air";
  parameter Modelica.SIunits.Emissivity absIRSha_glass[nCon] 
    "Infrared absorptivity of shade surface that faces glass";
  parameter Modelica.SIunits.TransmissionCoefficient tauIRSha_air[nCon] 
    "Infrared transmissivity of shade for radiation coming from the exterior or the room";
  parameter Modelica.SIunits.TransmissionCoefficient tauIRSha_glass[nCon] 
    "Infrared transmissivity of shade for radiation coming from the glass";
  parameter Boolean haveExteriorShade[nCon] 
    "Set to true if window has exterior shade (at surface a)";
  parameter Boolean haveInteriorShade[nCon] 
    "Set to true if window has interior shade (at surface b)";
  final parameter Boolean windowHasShade=
    haveExteriorShade[1] or haveInteriorShade[1] 
    "Set to true if window system has a shade";
  Modelica.Blocks.Interfaces.RealInput uSha[nCon](min=0, max=1) if 
       windowHasShade 
    "Control signal for the shading device, 0: unshaded; 1: fully shaded"; 
  Modelica.Blocks.Interfaces.RealInput QAbsSolSha_flow[nCon](
    final unit="W", quantity="Power") "Solar radiation absorbed by shade"; 
  HeatTransfer.Windows.ExteriorHeatTransfer conExtWin[nCon](
    final A=AWin,
    final fFra=fFra,
    each final linearizeRadiation = linearizeRadiation,
    final vieFacSky={(Modelica.Constants.pi - til[i]) ./ Modelica.Constants.pi for i in 1:nCon},
    final absIRSha_air=absIRSha_air,
    final absIRSha_glass=absIRSha_glass,
    final tauIRSha_air=tauIRSha_air,
    final tauIRSha_glass=tauIRSha_glass,
    final haveExteriorShade=haveExteriorShade,
    final haveInteriorShade=haveInteriorShade) 
    "Exterior convection of the window"; 
  SkyRadiationExchange skyRadExcWin(
    final n=nCon,
    each final absIR=absIRFra,
    vieFacSky={(Modelica.Constants.pi - til[i]) ./ Modelica.Constants.pi for i in 
            1:nCon},
    each final A=AWin .* fFra) 
    "Infrared radiative heat exchange between window frame and sky";
  HeatTransfer.Interfaces.RadiosityOutflow JOutUns[nCon] 
    "Outgoing radiosity that connects to unshaded part of glass at exterior side";
  HeatTransfer.Interfaces.RadiosityInflow JInUns[nCon] 
    "Incoming radiosity that connects to unshaded part of glass at exterior side";
  HeatTransfer.Interfaces.RadiosityOutflow JOutSha[nCon] if 
       windowHasShade 
    "Outgoing radiosity that connects to shaded part of glass at exterior side";
  HeatTransfer.Interfaces.RadiosityInflow JInSha[nCon] if 
       windowHasShade 
    "Incoming radiosity that connects to shaded part of glass at exterior side";
  Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a glaUns[nCon] 
    "Heat port at unshaded glass of exterior-facing surface";
  Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a glaSha[nCon] if 
       windowHasShade "Heat port at shaded glass of exterior-facing surface";
  Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a fra[nCon](T(each nominal=
          300, each start=283.15)) 
    "Heat port at frame of exterior-facing surface";
  Modelica.Blocks.Math.Add HTotConExtWinFra[nCon](
     final k1=fFra .* absSolFra .* AWin,
     final k2=fFra .* absSolFra .* AWin) 
    "Total solar irradiation on window frame";
  Buildings.HeatTransfer.Sources.PrescribedHeatFlow solHeaGaiConWin[nCon] 
    "Total solar heat gain of the window frame";
  Modelica.Blocks.Interfaces.RealOutput HDir[nCon](
     each final quantity="RadiantEnergyFluenceRate",
     each final unit="W/m2") "Direct solar irradition on tilted surface";
  Modelica.Blocks.Interfaces.RealOutput HDif[nCon](
     each final quantity="RadiantEnergyFluenceRate",
     each final unit="W/m2") "Diffuse solar irradiation on tilted surface";
  Modelica.Blocks.Interfaces.RealOutput inc[nCon](
    each final quantity="Angle",
    each final unit="rad",
    each displayUnit="deg") "Incidence angle"; 
protected 
  Buildings.HeatTransfer.Sources.PrescribedTemperature TAirConExtWin[
    nCon] "Outside air temperature for window constructions";
  Modelica.Blocks.Routing.Replicator repConExtWin(final nout=nCon) 
    "Signal replicator";
  Modelica.Blocks.Routing.Replicator repConExtWinVWin(final nout=nCon) 
    "Signal replicator";
  Modelica.Blocks.Routing.Replicator repConExtWinTSkyBla(final nout=nCon) 
    "Signal replicator"; 
equation 
  connect(uSha, conExtWin.uSha);
  connect(JInUns,conExtWin. JInUns);
  connect(conExtWin.JOutUns,JOutUns);
  connect(conExtWin.glaUns,glaUns);
  connect(conExtWin.glaSha,glaSha);
  connect(conExtWin.JOutSha,JOutSha);
  connect(conExtWin.JInSha,JInSha);
  connect(conExtWin.frame,fra);
  connect(TAirConExtWin.port,conExtWin. air);
  connect(TAirConExtWin.T,repConExtWin. y);
  connect(repConExtWin.u, weaBus.TDryBul);
  connect(repConExtWinVWin.y,conExtWin. vWin);
  connect(repConExtWinVWin.u, weaBus.winSpe);
  connect(HTotConExtWinFra.y, solHeaGaiConWin.Q_flow);
  connect(solHeaGaiConWin.port, fra);
  connect(HDirTil.H, HDir);
  connect(HDifTil.H, HDif);
  connect(HDirTil.inc, inc);
  connect(HTotConExtWinFra.u2, HDifTil.H);
  connect(HTotConExtWinFra.u1, HDirTil.H);
  connect(conExtWin.QAbs_flow, QAbsSolSha_flow);
  connect(skyRadExcWin.TOut, weaBus.TDryBul);
  connect(skyRadExcWin.TBlaSky, weaBus.TBlaSky);
  connect(skyRadExcWin.port, fra);
  connect(repConExtWin.y, conExtWin.TOut);
  connect(repConExtWinTSkyBla.y, conExtWin.TBlaSky);
  connect(repConExtWinTSkyBla.u, weaBus.TBlaSky); 
end ExteriorBoundaryConditionsWithWindow;
 
 Buildings.Rooms.BaseClasses.HeatGain
Buildings.Rooms.BaseClasses.HeatGain
 
Extends from Buildings.BaseClasses.BaseIcon (Base icon).
| Type | Name | Default | Description | 
|---|---|---|---|
| Area | AFlo | Floor area [m2] | |
| Temperature | TWat | 273.15 + 34 | Temperature at which water vapor is released [K] | 
| Type | Name | Description | 
|---|---|---|
| HeatPort_a | QCon_flow | Convective heat gain | 
| input RealInput | qGai_flow[3] | Radiant, convective and latent heat input into room (positive if heat gain) | 
| output RealOutput | QRad_flow | Radiative heat gain | 
| FluidPort_b | QLat_flow | Latent heat gain | 
model HeatGain "Model to convert internal heat gain signals"
  extends Buildings.BaseClasses.BaseIcon;
  replaceable package Medium =
     Modelica.Media.Interfaces.PartialMedium "Medium in the component";
  parameter Modelica.SIunits.Area AFlo "Floor area";
  parameter Modelica.SIunits.Temperature TWat=273.15+34 
    "Temperature at which water vapor is released";
  Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a QCon_flow 
    "Convective heat gain"; 
public 
  Modelica.Blocks.Interfaces.RealInput qGai_flow[3] 
    "Radiant, convective and latent heat input into room (positive if heat gain)";
  Modelica.Blocks.Interfaces.RealOutput QRad_flow "Radiative heat gain";
  Modelica.Fluid.Interfaces.FluidPort_b QLat_flow(redeclare final package
      Medium =
        Medium) "Latent heat gain"; 
protected 
  Modelica.SIunits.MassFlowRate mWat_flow 
    "Water vapor flow rate released by latent gain";
  parameter Modelica.SIunits.SpecificEnergy h_fg = Medium.enthalpyOfCondensingGas(TWat) 
    "Latent heat of water vapor";
 // constant Medium.MassFraction Xi[Medium.nXi] = {1}
 //   "Species concentration (water vapor only)";
protected 
  parameter Real s[Medium.nX](fixed=false) 
    "Vector with zero everywhere except where water vapor is";
initial algorithm 
  for i in 1:Medium.nX loop
    if ( Modelica.Utilities.Strings.isEqual(string1=Medium.substanceNames[i],
                                            string2="water",
                                            caseSensitive=false)) then
      s[i] :=1;
    else
      s[i] :=0;
    end if;
   end for;
   assert(Medium.nX == 1 or abs(1-sum(s)) < 1E-4, "Substance 'water' is not present in medium '"
                  + Medium.mediumName + "'.\n"
                  + "Change medium model to one that has 'water' as a substance.");
equation 
  QRad_flow = qGai_flow[1]*AFlo;
  QCon_flow.Q_flow = -qGai_flow[2]*AFlo;
  // Interface to fluid port
  // If a medium does not contain water vapor, then h_fg is equal to zero.
  if Medium.nXi == 0 or (h_fg == 0) then
    mWat_flow = 0;
  else
    mWat_flow = qGai_flow[3]*AFlo/h_fg;
  end if;
  QLat_flow.C_outflow  = fill(0, Medium.nC);
  QLat_flow.h_outflow  = h_fg;
  QLat_flow.Xi_outflow = s[1:Medium.nXi];
  QLat_flow.m_flow     = if (h_fg > 0) then 
                           -qGai_flow[3]*AFlo/h_fg else 
                            0;
end HeatGain;
 
 Buildings.Rooms.BaseClasses.InfraredRadiationExchange
Buildings.Rooms.BaseClasses.InfraredRadiationExchange
 
This model computes the infrared radiative heat transfer between the interior surfaces of a room. Each opaque surface emits radiation according to
Ei = σ Ai εi (Ti)4,
where
σ
is the Stefan-Boltzmann constant,
Ai 
is the surface area,
εi 
is the absorptivity in the infrared spectrum, and
Ti
is the surface temperature.
If the parameter linearizeRadidation is set to true,
then the term (Ti)4 is replaced with
T03 Ti,
where T0 = 20°C is a parameter.
The incoming radiation at surface i is
Gi = -∑j Fj,i Jj
where Fj,i is the view factor from surface j to surface i, Jj is the radiosity leaving surface j and the sum is over all surfaces. For opaque surfaces, it follows from the first law that the radiosity Ji is
Ji = -Ei - (1-εi) Gi.
For windows, the outgoing radiosity is an input into this model because the window model computes this quantity directly.
For each surface i, the heat balance is
0 = Qi + Ji + Gi.
For opaque surfaces, the heat flow rate Qi is set to be equal to the heat flow rate at the heat port. For the glass of the windows, the radiosity outflow at the connector is set to the radiosity Gi that is leaving the surface.
The view factor from surface i to j is approximated as
Fi,j = Aj ⁄ ∑k Ak.
Extends from Buildings.Rooms.BaseClasses.PartialSurfaceInterface (Partial model that is used for infrared radiation balance).
| Type | Name | Default | Description | 
|---|---|---|---|
| ParameterConstruction | datConExt[NConExt] | Data for exterior construction | |
| ParameterConstructionWithWindow | datConExtWin[NConExtWin] | Data for exterior construction with window | |
| ParameterConstruction | datConPar[NConPar] | Data for partition construction | |
| ParameterConstruction | datConBou[NConBou] | Data for construction boundary | |
| Generic | surBou[NSurBou] | Record for data of surfaces whose heat conduction is modeled outside of this room | |
| Boolean | linearizeRadiation | Set to true to linearize emissive power | |
| Exterior constructions | |||
| Integer | nConExt | Number of exterior constructions | |
| Integer | nConExtWin | Number of window constructions | |
| Partition constructions | |||
| Integer | nConPar | Number of partition constructions | |
| Boundary constructions | |||
| Integer | nConBou | Number of constructions that have their outside surface exposed to the boundary of this room | |
| Integer | nSurBou | Number of surface heat transfer models that connect to constructions that are modeled outside of this room | |
| Advanced | |||
| Boolean | homotopyInitialization | true | = true, use homotopy method | 
| Type | Name | Description | 
|---|---|---|
| HeatPort_a | conExt[NConExt] | Heat port that connects to room-side surface of exterior constructions | 
| HeatPort_a | conExtWin[NConExtWin] | Heat port that connects to room-side surface of exterior constructions that contain a window | 
| HeatPort_a | conExtWinFra[NConExtWin] | Heat port that connects to room-side surface of window frame | 
| HeatPort_a | conPar_a[NConPar] | Heat port that connects to room-side surface a of partition constructions | 
| HeatPort_a | conPar_b[NConPar] | Heat port that connects to room-side surface b of partition constructions | 
| HeatPort_a | conBou[NConBou] | Heat port that connects to room-side surface of constructions that expose their other surface to the outside | 
| HeatPort_a | conSurBou[NSurBou] | Heat port to surfaces of models that compute the heat conduction outside of this room | 
| input RadiosityInflow | JInConExtWin[NConExtWin] | Incoming radiosity that connects to non-frame part of the window [W] | 
| output RadiosityOutflow | JOutConExtWin[NConExtWin] | Outgoing radiosity that connects to non-frame part of the window [W] | 
model InfraredRadiationExchange 
  "Infrared radiation heat exchange between the room facing surfaces"
  extends Buildings.Rooms.BaseClasses.PartialSurfaceInterface;
  parameter Boolean linearizeRadiation 
    "Set to true to linearize emissive power";
  parameter Boolean homotopyInitialization=true "= true, use homotopy method";
  HeatTransfer.Interfaces.RadiosityInflow JInConExtWin[NConExtWin] 
    "Incoming radiosity that connects to non-frame part of the window";
  HeatTransfer.Interfaces.RadiosityOutflow JOutConExtWin[NConExtWin] 
    "Outgoing radiosity that connects to non-frame part of the window"; 
protected 
  final parameter Integer NOpa=NConExt + 2*NConExtWin + 2*NConPar + NConBou +
      NSurBou "Number of opaque surfaces, including the window frame";
  final parameter Integer nOpa=nConExt + 2*nConExtWin + 2*nConPar + nConBou +
      nSurBou "Number of opaque surfaces, including the window frame";
  final parameter Integer NWin=NConExtWin "Number of window surfaces";
  final parameter Integer nWin=nConExtWin "Number of window surfaces";
  final parameter Integer NTot=NOpa + NWin "Total number of surfaces";
  final parameter Integer nTot=nOpa + nWin "Total number of surfaces";
  final parameter Real epsOpa[nOpa](
    min=0,
    max=1,
    fixed=false) "Absorptivity of opaque surfaces";
  final parameter Real rhoOpa[nOpa](
    min=0,
    max=1,
    fixed=false) "Reflectivity of opaque surfaces";
  final parameter Modelica.SIunits.Area AOpa[nOpa](fixed=false) 
    "Surface area of opaque surfaces";
  final parameter Modelica.SIunits.Area A[nTot](fixed=false) "Surface areas";
  final parameter Real kOpa[nOpa](unit="W/K4", fixed=false) 
    "Product sigma*epsilon*A for opaque surfaces";
  final parameter Real F[nTot, nTot](
    min=0,
    max=1,
    fixed=false) "View factor from surface i to j";
  Modelica.SIunits.HeatFlowRate J[nTot](
    max=0,
    start=A .* 0.8*Modelica.Constants.sigma*293.15^4,
    each nominal=10*0.8*Modelica.Constants.sigma*293.15^4) 
    "Radiosity leaving the surface";
  Modelica.SIunits.HeatFlowRate G[nTot](
    min=0,
    start=A .* 0.8*Modelica.Constants.sigma*293.15^4,
    each nominal=10*0.8*Modelica.Constants.sigma*293.15^4) 
    "Radiosity entering the surface";
  constant Real T30(unit="K3") = 293.15^3 "Nominal temperature";
  constant Real T40(unit="K4") = 293.15^4 "Nominal temperature";
  Modelica.SIunits.Temperature TOpa[nOpa](each start=293.15, each nominal=
        293.15) "Temperature of opaque surfaces";
  Real T4Opa[nOpa](
    each unit="K4",
    each start=T40,
    each nominal=293.15^4) "Forth power of temperature of opaque surfaces";
  Modelica.SIunits.HeatFlowRate Q_flow[nTot] "Heat flow rate at surfaces";
  parameter Modelica.SIunits.Temperature T0=293.15 
    "Temperature used to linearize radiative heat transfer";
  final parameter Real T03(
    min=0,
    unit="K3") = T0^3 "3rd power of temperature T0";
  Modelica.SIunits.HeatFlowRate sumEBal "Sum of energy balance, should be zero";
initial equation 
  // The next loops build the array epsOpa, AOpa and kOpa that simplify
  // the model equations.
  // These arrays store the values of the constructios in the following order
  // [x[1:NConExt] x[1:NConPar] x[1: NConPar] x[1: NConBou] x[1: NSurBou] x[1: NConExtWin] x[1: NConExtWin]]
  // where x is epsOpa, AOpa or kOpa.
  // The last two entries are for the opaque wall that contains a window, and for the window frame.
  for i in 1:nConExt loop
    epsOpa[i] = epsConExt[i];
    AOpa[i] = AConExt[i];
    kOpa[i] = Modelica.Constants.sigma*epsConExt[i]*AOpa[i];
  end for;
  for i in 1:nConPar loop
    epsOpa[i + nConExt] = epsConPar_a[i];
    AOpa[i + nConExt] = AConPar[i];
    kOpa[i + nConExt] = Modelica.Constants.sigma*epsConPar_a[i]*AOpa[i +
      nConExt];
    epsOpa[i + nConExt + nConPar] = epsConPar_b[i];
    AOpa[i + nConExt + nConPar] = AConPar[i];
    kOpa[i + nConExt + nConPar] = Modelica.Constants.sigma*epsConPar_b[i]*AOpa[
      i + nConExt + nConPar];
  end for;
  for i in 1:nConBou loop
    epsOpa[i + nConExt + 2*nConPar] = epsConBou[i];
    AOpa[i + nConExt + 2*nConPar] = AConBou[i];
    kOpa[i + nConExt + 2*nConPar] = Modelica.Constants.sigma*epsConBou[i]*AOpa[
      i + nConExt + 2*nConPar];
  end for;
  for i in 1:nSurBou loop
    epsOpa[i + nConExt + 2*nConPar + nConBou] = epsSurBou[i];
    AOpa[i + nConExt + 2*nConPar + nConBou] = ASurBou[i];
    kOpa[i + nConExt + 2*nConPar + nConBou] = Modelica.Constants.sigma*
      epsSurBou[i]*AOpa[i + nConExt + 2*nConPar + nConBou];
  end for;
  for i in 1:nConExtWin loop
    // Opaque part of construction that has a window embedded
    epsOpa[i + nConExt + 2*nConPar + nConBou + nSurBou] = epsConExtWinOpa[i];
    AOpa[i + nConExt + 2*nConPar + nConBou + nSurBou] = AConExtWinOpa[i];
    kOpa[i + nConExt + 2*nConPar + nConBou + nSurBou] = Modelica.Constants.sigma
      *epsConExtWinOpa[i]*AOpa[i + nConExt + 2*nConPar + nConBou + nSurBou];
    // Window frame
    epsOpa[i + nConExt + 2*nConPar + nConBou + nSurBou + nConExtWin] =
      epsConExtWinFra[i];
    AOpa[i + nConExt + 2*nConPar + nConBou + nSurBou + nConExtWin] =
      AConExtWinFra[i];
    kOpa[i + nConExt + 2*nConPar + nConBou + nSurBou + nConExtWin] = Modelica.Constants.sigma
      *epsConExtWinFra[i]*AOpa[i + nConExt + 2*nConPar + nConBou + nSurBou +
      nConExtWin];
  end for;
  // Vector with all surface areas.
  // The next loops build the array A that simplifies
  // the model equations.
  // These array stores the values of the constructios in the following order
  // [AOpa[1:nConExt] AOpa[1:nConPar] AOpa[1: nConPar] AOpa[1: nConBou] AOpa[1: nSurBou]
  //  AOpa[1: nConExtWin] AOpa[1: nConExtWin] AGla[1: nConExtWin]]
  // since nWin=nConExtWin.
  for i in 1:nOpa loop
    A[i] = AOpa[i];
  end for;
  for i in 1:nWin loop
    A[i + nOpa] = AConExtWinGla[i];
  end for;
  // Reflectivity for opaque surfaces
  rhoOpa = 1 .- epsOpa;
  // View factors from surface i to surface j
  for i in 1:nTot loop
    for j in 1:nTot loop
      F[i, j] = A[j]/sum((A[k]) for k in 1:nTot);
    end for;
  end for;
  // Test whether the view factors add up to one, or the sum is zero in case there
  // is only one construction
  for i in 1:nTot loop
    assert((abs(1 - sum(F[i, j] for j in 1:nTot))) < 1E-10,
      "Program error: Sum 1 of view factors is " + String(sum(F[i, j] for j in 
      1:nTot)));
  end for;
  ////////////////////////////////////////////////////////////////////
equation 
  // If the room has no window, then the incoming radiosity from the window
  // is not connected. In this situation, we set it to zero.
  // This approach is easier than using a conditional connector, since
  // this port carries a flow variable, and hence the sign of the radiosity
  // would change in a connect statement.
  if (cardinality(JInConExtWin) == 0) then
    JInConExtWin = zeros(NConExtWin);
  end if;
  // Assign temperature of opaque surfaces
  for i in 1:nConExt loop
    TOpa[i] = conExt[i].T;
  end for;
  for i in 1:nConPar loop
    TOpa[i + nConExt] = conPar_a[i].T;
    TOpa[i + nConExt + nConPar] = conPar_b[i].T;
  end for;
  for i in 1:nConBou loop
    TOpa[i + nConExt + 2*nConPar] = conBou[i].T;
  end for;
  for i in 1:nSurBou loop
    TOpa[i + nConExt + 2*nConPar + nConBou] = conSurBou[i].T;
  end for;
  for i in 1:nConExtWin loop
    TOpa[i + nConExt + 2*nConPar + nConBou + nSurBou] = conExtWin[i].T;
    TOpa[i + nConExt + 2*nConPar + nConBou + nConExtWin + nSurBou] =
      conExtWinFra[i].T;
  end for;
  // Incoming radiosity at each surface
  // is equal to the negative of the outgoing radiosity of
  // all other surfaces times the view factor
  G = -transpose(F)*J;
  // Outgoing radiosity
  // Opaque surfaces.
  // If kOpa[j]=absIR[j]*A[j] < 1E-28, then A < 1E-20 and the surface is
  // from a dummy construction. In this situation, we set T40=293.15^4 to
  // avoid a singularity.
  for j in 1:nOpa loop
    //   T4Opa[j] = if (kOpa[j] > 1E-28) then (Q_flow[j]-epsOpa[j] * G[j])/kOpa[j] else T40;
    T4Opa[j] = (-J[j] - rhoOpa[j]*G[j])/kOpa[j];
  end for;
  // 4th power of temperature
  if linearizeRadiation then
    TOpa = (T4Opa .+ 3*T40)/(4*T30);
    // Based on T4 = 4*T30*T-3*T40
  else
    if homotopyInitialization then
      TOpa = homotopy(actual=Buildings.Utilities.Math.Functions.powerLinearized(
        x=T4Opa,
        x0=243.15^4,
        n=0.25), simplified=(T4Opa .+ 3*T40)/(4*T30));
    else
      TOpa = Buildings.Utilities.Math.Functions.powerLinearized(
        x=T4Opa,
        x0=243.15^4,
        n=0.25);
    end if;
  end if;
  // Assign radiosity that comes from window
  // and that leaves window.
  // J < 0 because it leaves the surface
  // G > 0 because it strikes the surface
  // JIn > 0 because it enters the model
  // JOut < 0 because it leaves the model
  for j in 1:nWin loop
    J[j + nOpa] = -JInConExtWin[j];
    G[j + nOpa] = -JOutConExtWin[j];
  end for;
  // Net heat exchange
  Q_flow = -J - G;
  // Assign heat exchange to connectors
  for i in 1:nConExt loop
    Q_flow[i] = conExt[i].Q_flow;
  end for;
  if nConExt == 0 then
    conExt[1].T = T0;
  end if;
  for i in 1:nConPar loop
    Q_flow[i + nConExt] = conPar_a[i].Q_flow;
    Q_flow[i + nConExt + nConPar] = conPar_b[i].Q_flow;
  end for;
  if nConPar == 0 then
    conPar_a[1].T = T0;
    conPar_b[1].T = T0;
  end if;
  for i in 1:nConBou loop
    Q_flow[i + nConExt + 2*nConPar] = conBou[i].Q_flow;
  end for;
  if nConBou == 0 then
    conBou[1].T = T0;
  end if;
  for i in 1:nSurBou loop
    Q_flow[i + nConExt + 2*nConPar + nConBou] = conSurBou[i].Q_flow;
  end for;
  if nSurBou == 0 then
    conSurBou[1].T = T0;
  end if;
  for i in 1:nConExtWin loop
    Q_flow[i + nConExt + 2*nConPar + nConBou + nSurBou] = conExtWin[i].Q_flow;
    Q_flow[i + nConExt + 2*nConPar + nConBou + nSurBou + nConExtWin] =
      conExtWinFra[i].Q_flow;
  end for;
  if nConExtWin == 0 then
    conExtWin[1].T = T0;
    conExtWinFra[1].T = T0;
    JOutConExtWin[1] = 0;
  end if;
  // Sum of energy balance
  // Remove sumEBal and assert statement for final release
  sumEBal = sum(conExt.Q_flow) + sum(conPar_a.Q_flow) + sum(conPar_b.Q_flow) +
    sum(conBou.Q_flow) + sum(conSurBou.Q_flow) + sum(conExtWin.Q_flow) + sum(
    conExtWinFra.Q_flow) + (sum(JInConExtWin) + sum(JOutConExtWin));
  assert(abs(sumEBal) < 1E-1,
    "Program error: Energy is not conserved in InfraredRadiationExchange." +
    "\n  Sum of all energy is " + String(sumEBal));
end InfraredRadiationExchange;
 
 Buildings.Rooms.BaseClasses.InfraredRadiationGainDistribution
Buildings.Rooms.BaseClasses.InfraredRadiationGainDistribution
 
Qi = Q Ai εi ⁄ ∑k Ak εk.
For opaque surfaces, the heat flow rate Qi is set to be equal to the heat flow rate at the heat port. For the glass of the windows, the heat flow rate Qi is set to the radiosity Ji that will strike the glass or the window shade.
Extends from Buildings.Rooms.BaseClasses.PartialSurfaceInterface (Partial model that is used for infrared radiation balance).
| Type | Name | Default | Description | 
|---|---|---|---|
| ParameterConstruction | datConExt[NConExt] | Data for exterior construction | |
| ParameterConstructionWithWindow | datConExtWin[NConExtWin] | Data for exterior construction with window | |
| ParameterConstruction | datConPar[NConPar] | Data for partition construction | |
| ParameterConstruction | datConBou[NConBou] | Data for construction boundary | |
| Generic | surBou[NSurBou] | Record for data of surfaces whose heat conduction is modeled outside of this room | |
| Boolean | haveShade | Set to true if a shade is present | |
| Exterior constructions | |||
| Integer | nConExt | Number of exterior constructions | |
| Integer | nConExtWin | Number of window constructions | |
| Partition constructions | |||
| Integer | nConPar | Number of partition constructions | |
| Boundary constructions | |||
| Integer | nConBou | Number of constructions that have their outside surface exposed to the boundary of this room | |
| Integer | nSurBou | Number of surface heat transfer models that connect to constructions that are modeled outside of this room | |
| Type | Name | Description | 
|---|---|---|
| HeatPort_a | conExt[NConExt] | Heat port that connects to room-side surface of exterior constructions | 
| HeatPort_a | conExtWin[NConExtWin] | Heat port that connects to room-side surface of exterior constructions that contain a window | 
| HeatPort_a | conExtWinFra[NConExtWin] | Heat port that connects to room-side surface of window frame | 
| HeatPort_a | conPar_a[NConPar] | Heat port that connects to room-side surface a of partition constructions | 
| HeatPort_a | conPar_b[NConPar] | Heat port that connects to room-side surface b of partition constructions | 
| HeatPort_a | conBou[NConBou] | Heat port that connects to room-side surface of constructions that expose their other surface to the outside | 
| HeatPort_a | conSurBou[NSurBou] | Heat port to surfaces of models that compute the heat conduction outside of this room | 
| input RealInput | uSha[NConExtWin] | Control signal for the shading device (removed if no shade is present) | 
| input RealInput | Q_flow | Radiative heat input into room (positive if heat gain) | 
| output RadiosityOutflow | JOutConExtWin[NConExtWin] | Outgoing radiosity that connects to shaded and unshaded part of glass [W] | 
model InfraredRadiationGainDistribution "Infrared radiative heat gain distribution between the room facing surfaces" extends Buildings.Rooms.BaseClasses.PartialSurfaceInterface; parameter Boolean haveShade "Set to true if a shade is present";Modelica.Blocks.Interfaces.RealInput uSha[NConExtWin](each min=0, each max=1) if haveShade "Control signal for the shading device (removed if no shade is present)"; Modelica.Blocks.Interfaces.RealInput Q_flow "Radiative heat input into room (positive if heat gain)"; HeatTransfer.Interfaces.RadiosityOutflow[NConExtWin] JOutConExtWin "Outgoing radiosity that connects to shaded and unshaded part of glass"; protected Real fraConExt[NConExt] = AEpsConExt/sumAEps "Fraction of infrared radiant heat gain absorbed by exterior constructions"; Real fraConExtWinOpa[NConExtWin] = AEpsConExtWinOpa/sumAEps "Fraction of infrared radiant heat gain absorbed by opaque part of exterior constructions that have a window"; Real fraConExtWinGla[NConExtWin] = (AEpsConExtWinSha + AEpsConExtWinUns)/sumAEps "Fraction of infrared radiant heat gain absorbed by opaque part of glass constructions that have a window"; Real fraConExtWinFra[NConExtWin] = AEpsConExtWinFra/sumAEps "Fraction of infrared radiant heat gain absorbed by window frame of exterior constructions that have a window"; Real fraConPar_a[NConPar] = AEpsConPar_a/sumAEps "Fraction of infrared radiant heat gain absorbed by partition constructions surface a"; Real fraConPar_b[NConPar] = AEpsConPar_b/sumAEps "Fraction of infrared radiant heat gain absorbed by partition constructions surface b"; Real fraConBou[NConBou] = AEpsConBou/sumAEps "Fraction of infrared radiant heat gain absorbed by constructions with exterior boundary conditions exposed to outside of room model"; Real fraSurBou[NSurBou] = AEpsSurBou/sumAEps "Fraction of infrared radiant heat gain absorbed by surface models of constructions that are modeled outside of this room"; parameter Real AEpsConExt[NConExt] = {AConExt[i]*epsConExt[i] for i in 1:NConExt} "Absorptivity times area of exterior constructions"; parameter Real AEpsConExtWinOpa[NConExtWin] = {AConExtWinOpa[i]*epsConExtWinOpa[i] for i in 1:NConExtWin} "Absorptivity times area of opaque part of exterior constructions that contain a window"; Real AEpsConExtWinUns[NConExtWin] = {shaSig[i].yCom * AConExtWinGla[i]*epsConExtWinUns[i] for i in 1:NConExtWin} "Absorptivity times area of unshaded window constructions"; Real AEpsConExtWinSha[NConExtWin] = {shaSig[i].y * AConExtWinGla[i]*epsConExtWinSha[i] for i in 1:NConExtWin} "Absorptivity times area of shaded window constructions"; parameter Real AEpsConExtWinFra[NConExtWin] = {AConExtWinFra[i]*epsConExtWinFra[i] for i in 1:NConExtWin} "Absorptivity times area of window frame"; parameter Real AEpsConPar_a[NConPar] = {AConPar[i]*epsConPar_a[i] for i in 1:NConPar} "Absorptivity times area of partition constructions surface a"; parameter Real AEpsConPar_b[NConPar] = {AConPar[i]*epsConPar_b[i] for i in 1:NConPar} "Absorptivity times area of partition constructions surface b"; parameter Real AEpsConBou[NConBou] = {AConBou[i]*epsConBou[i] for i in 1:NConBou} "Absorptivity times area of constructions with exterior boundary conditions exposed to outside of room model"; parameter Real AEpsSurBou[NSurBou] = {ASurBou[i]*epsSurBou[i] for i in 1:NSurBou} "Absorptivity times area of surface models of constructions that are modeled outside of this room"; parameter Real sumAEpsNoWin(fixed=false) "Sum of absorptivity times area of all constructions except for windows"; Real sumAEps "Sum of absorptivity times area of all constructions including windows"; Buildings.HeatTransfer.Windows.BaseClasses.ShadingSignal shaSig[NConExtWin]( each final haveShade=haveShade) "Block to constrain the shading control signal to be strictly within (0, 1) if a shade is present"; initial equation sumAEpsNoWin = sum(AEpsConExt)+sum(AEpsConExtWinOpa)+sum(AEpsConExtWinFra) +sum(AEpsConPar_a)+sum(AEpsConPar_b)+sum(AEpsConBou)+sum(AEpsSurBou); equation connect(uSha, shaSig.u); sumAEps = sumAEpsNoWin + sum(AEpsConExtWinUns) + sum(AEpsConExtWinSha); // Infrared radiative heat flow conExt.Q_flow = -fraConExt*Q_flow; conExtWin.Q_flow = -fraConExtWinOpa*Q_flow; conPar_a.Q_flow = -fraConPar_a*Q_flow; conPar_b.Q_flow = -fraConPar_b*Q_flow; conBou.Q_flow = -fraConBou*Q_flow; conSurBou.Q_flow = -fraSurBou*Q_flow; // This model makes the simplification that the shade, the glass and the frame have // the same absorptivity in the infrared region JOutConExtWin = -fraConExtWinGla*Q_flow; conExtWinFra.Q_flow = -fraConExtWinFra*Q_flow; // Check for conservation of energy assert(abs(1 - sum(fraConExt) - sum(fraConExtWinOpa)- sum(fraConExtWinGla) - sum(fraConExtWinFra) - sum(fraConPar_a) - sum(fraConPar_b) - sum(fraConBou) - sum(fraSurBou)) < 1E-5, "Programming error: Radiation balance is wrong. Check equations.");end InfraredRadiationGainDistribution; 
 Buildings.Rooms.BaseClasses.MixedAir
Buildings.Rooms.BaseClasses.MixedAir
 
The main components that are used in this model are as follows:
Extends from Buildings.Fluid.Interfaces.LumpedVolumeDeclarations (Declarations for lumped volumes), Buildings.Rooms.BaseClasses.PartialSurfaceInterface (Partial model that is used for infrared radiation balance).
| Type | Name | Default | Description | 
|---|---|---|---|
| replaceable package Medium | PartialMedium | Medium in the component | |
| ParameterConstruction | datConExt[NConExt] | Data for exterior construction | |
| ParameterConstructionWithWindow | datConExtWin[NConExtWin] | Data for exterior construction with window | |
| ParameterConstruction | datConPar[NConPar] | Data for partition construction | |
| ParameterConstruction | datConBou[NConBou] | Data for construction boundary | |
| Generic | surBou[NSurBou] | Record for data of surfaces whose heat conduction is modeled outside of this room | |
| Volume | V | Volume [m3] | |
| Area | AFlo | Floor area [m2] | |
| Length | hRoo | Average room height [m] | |
| Boolean | isFloorSurBou[NSurBou] | surBou.isFloor | Flag to indicate if floor for constructions that are modeled outside of this room | 
| Emissivity | tauGlaSol[NConExtWin] | Transmissivity of window [1] | |
| Boolean | linearizeRadiation | Set to true to linearize emissive power | |
| Exterior constructions | |||
| Integer | nConExt | Number of exterior constructions | |
| Integer | nConExtWin | Number of window constructions | |
| Partition constructions | |||
| Integer | nConPar | Number of partition constructions | |
| Boundary constructions | |||
| Integer | nConBou | Number of constructions that have their outside surface exposed to the boundary of this room | |
| Integer | nSurBou | Number of surface heat transfer models that connect to constructions that are modeled outside of this room | |
| Convective heat transfer | |||
| InteriorConvection | conMod | Buildings.HeatTransfer.Types... | Convective heat transfer model for opaque constructions | 
| CoefficientOfHeatTransfer | hFixed | 3.0 | Constant convection coefficient for opaque constructions [W/(m2.K)] | 
| Nominal condition | |||
| MassFlowRate | m_flow_nominal | Nominal mass flow rate [kg/s] | |
| Dynamics | |||
| Equations | |||
| Dynamics | energyDynamics | Modelica.Fluid.Types.Dynamic... | Formulation of energy balance | 
| Dynamics | massDynamics | energyDynamics | Formulation of mass balance | 
| Initialization | |||
| AbsolutePressure | p_start | Medium.p_default | Start value of pressure [Pa] | 
| Temperature | T_start | Medium.T_default | Start value of temperature [K] | 
| MassFraction | X_start[Medium.nX] | Medium.X_default | Start value of mass fractions m_i/m [kg/kg] | 
| ExtraProperty | C_start[Medium.nC] | fill(0, Medium.nC) | Start value of trace substances | 
| ExtraProperty | C_nominal[Medium.nC] | fill(1E-2, Medium.nC) | Nominal value of trace substances. (Set to typical order of magnitude.) | 
| Type | Name | Description | 
|---|---|---|
| HeatPort_a | conExt[NConExt] | Heat port that connects to room-side surface of exterior constructions | 
| HeatPort_a | conExtWin[NConExtWin] | Heat port that connects to room-side surface of exterior constructions that contain a window | 
| HeatPort_a | conExtWinFra[NConExtWin] | Heat port that connects to room-side surface of window frame | 
| HeatPort_a | conPar_a[NConPar] | Heat port that connects to room-side surface a of partition constructions | 
| HeatPort_a | conPar_b[NConPar] | Heat port that connects to room-side surface b of partition constructions | 
| HeatPort_a | conBou[NConBou] | Heat port that connects to room-side surface of constructions that expose their other surface to the outside | 
| HeatPort_a | conSurBou[NSurBou] | Heat port to surfaces of models that compute the heat conduction outside of this room | 
| HeatPort_a | glaUns[NConExtWin] | Heat port that connects to room-side surface of unshaded glass | 
| HeatPort_a | glaSha[NConExtWin] | Heat port that connects to room-side surface of shaded glass | 
| VesselFluidPorts_b | ports[nPorts] | Fluid inlets and outlets | 
| HeatPort_a | heaPorAir | Heat port to air volume | 
| input RadiosityInflow | JInSha[NConExtWin] | Incoming radiosity that connects to shaded part of glass [W] | 
| output RadiosityOutflow | JOutSha[NConExtWin] | Outgoing radiosity that connects to shaded part of glass [W] | 
| input RadiosityInflow | JInUns[NConExtWin] | Incoming radiosity that connects to unshaded part of glass [W] | 
| output RadiosityOutflow | JOutUns[NConExtWin] | Outgoing radiosity that connects to unshaded part of glass [W] | 
| input RealInput | qGai_flow[3] | Radiant, convective and latent heat input into room (positive if heat gain) | 
| input RealInput | uSha[NConExtWin] | Control signal for the shading device (removed if no shade is present) | 
| input RealInput | QAbsSolSha_flow[NConExtWin] | Solar radiation absorbed by shade [W] | 
| input RealInput | JInConExtWin[NConExtWin] | Solar radiation transmitted from the outside through the glazing system [W] | 
| output RealOutput | HOutConExtWin[NConExtWin] | Outgoing solar radiation that strikes window [W/m2] | 
| HeatPort_a | heaPorRad | Heat port for radiative heat gain and radiative temperature | 
model MixedAir "Model for room air that is completely mixed"
  extends Buildings.Fluid.Interfaces.LumpedVolumeDeclarations;
  extends Buildings.Rooms.BaseClasses.PartialSurfaceInterface;
  parameter Modelica.SIunits.Volume V "Volume";
  // Port definitions
  parameter Integer nPorts=0 "Number of fluid ports of this model";
  parameter Modelica.SIunits.Area AFlo "Floor area";
  parameter Modelica.SIunits.Length hRoo "Average room height";
  parameter Buildings.HeatTransfer.Types.InteriorConvection conMod=
  Buildings.HeatTransfer.Types.InteriorConvection.Temperature 
    "Convective heat transfer model for opaque constructions";
  parameter Modelica.SIunits.CoefficientOfHeatTransfer hFixed=3.0 
    "Constant convection coefficient for opaque constructions";
  final parameter Boolean isFloorConExt[NConExt]=
    datConExt.isFloor "Flag to indicate if floor for exterior constructions";
  final parameter Boolean isFloorConExtWin[NConExtWin]=
    datConExtWin.isFloor "Flag to indicate if floor for constructions";
  final parameter Boolean isFloorConPar_a[NConPar]=
    datConPar.isFloor "Flag to indicate if floor for constructions";
  final parameter Boolean isFloorConPar_b[NConPar]=
    datConPar.isCeiling "Flag to indicate if floor for constructions";
  final parameter Boolean isFloorConBou[NConBou]=
    datConBou.isFloor 
    "Flag to indicate if floor for constructions with exterior boundary conditions exposed to outside of room model";
  parameter Boolean isFloorSurBou[NSurBou]=
    surBou.isFloor 
    "Flag to indicate if floor for constructions that are modeled outside of this room";
  parameter Modelica.SIunits.Emissivity tauGlaSol[NConExtWin] 
    "Transmissivity of window";
  Buildings.Fluid.MixingVolumes.MixingVolume vol(V=AFlo*hRoo,
    redeclare package Medium = Medium,
    final energyDynamics=energyDynamics,
    final massDynamics=massDynamics,
    final p_start=p_start,
    final T_start=T_start,
    final X_start=X_start,
    final C_start=C_start,
    final C_nominal=C_nominal,
    final m_flow_nominal = m_flow_nominal,
    final prescribedHeatFlowRate = true,
    final nPorts=nPorts+1) "Room air volume"; 
  // Heat ports that are needed to connect to the window glass
  Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a glaUns[NConExtWin] if 
     haveConExtWin 
    "Heat port that connects to room-side surface of unshaded glass";
  Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a glaSha[NConExtWin] if 
       haveShade "Heat port that connects to room-side surface of shaded glass";
  Modelica.Fluid.Vessels.BaseClasses.VesselFluidPorts_b ports[nPorts](
      redeclare each final package Medium = Medium) "Fluid inlets and outlets"; 
  Buildings.HeatTransfer.Convection.Interior convConExt[
                                     NConExt](
    final A=AConExt,
    final til =  datConExt.til,
    each conMod=conMod,
    each hFixed=hFixed) if 
       haveConExt "Convective heat transfer";
  Buildings.HeatTransfer.Convection.Interior convConExtWin[
                                        NConExtWin](
    final A=AConExtWinOpa,
    final til =  datConExtWin.til,
    each conMod=conMod,
    each hFixed=hFixed) if 
       haveConExtWin "Convective heat transfer";
  HeatTransfer.Windows.InteriorHeatTransfer convConWin[NConExtWin](
    final fFra=fFra,
    each final linearizeRadiation = linearizeRadiation,
    final haveExteriorShade=haveExteriorShade,
    final haveInteriorShade=haveInteriorShade,
    final A=AConExtWinGla + AConExtWinFra,
    final absIRSha_air=epsConExtWinSha,
    final absIRSha_glass=epsConExtWinUns,
    final tauIRSha_air=tauIRSha_air,
    final tauIRSha_glass=tauIRSha_glass) if 
       haveConExtWin "Model for convective heat transfer at window"; 
  // For conPar_a, we use for the tilt pi-tilt since it is the
  // surface that is on the other side of the construction
  Buildings.HeatTransfer.Convection.Interior convConPar_a[
                                       nConPar](
    final A=AConPar,
    final til=Modelica.Constants.pi .- datConPar.til,
    each conMod=conMod,
    each hFixed=hFixed) if 
       haveConPar "Convective heat transfer";
  Buildings.HeatTransfer.Convection.Interior convConPar_b[
                                       nConPar](
    final A=AConPar,
    final til =  datConPar.til,
    each conMod=conMod,
    each hFixed=hFixed) if 
       haveConPar "Convective heat transfer";
  Buildings.HeatTransfer.Convection.Interior convConBou[
                                     nConBou](
    final A=AConBou,
    final til =  datConBou.til,
    each conMod=conMod,
    each hFixed=hFixed) if 
       haveConBou "Convective heat transfer";
  Buildings.HeatTransfer.Convection.Interior convSurBou[
                                     nSurBou](
    final A=ASurBou,
    final til =  surBou.til,
    each conMod=conMod,
    each hFixed=hFixed) if 
       haveSurBou "Convective heat transfer";
protected 
  Modelica.Thermal.HeatTransfer.Components.ThermalCollector theConConExt(final m=
        nConExt) if 
       haveConExt 
    "Thermal collector to convert from vector to scalar connector";
  Modelica.Thermal.HeatTransfer.Components.ThermalCollector theConConExtWin(
      final m=nConExtWin) if 
       haveConExtWin 
    "Thermal collector to convert from vector to scalar connector";
  Modelica.Thermal.HeatTransfer.Components.ThermalCollector theConConWin(final m=
        nConExtWin) if 
       haveConExtWin 
    "Thermal collector to convert from vector to scalar connector";
  Modelica.Thermal.HeatTransfer.Components.ThermalCollector theConConPar_a(final m=
        nConPar) if 
       haveConPar 
    "Thermal collector to convert from vector to scalar connector";
  Modelica.Thermal.HeatTransfer.Components.ThermalCollector theConConPar_b(final m=
        nConPar) if 
       haveConPar 
    "Thermal collector to convert from vector to scalar connector";
  Modelica.Thermal.HeatTransfer.Components.ThermalCollector theConConBou(final m=
        nConBou) if 
       haveConBou 
    "Thermal collector to convert from vector to scalar connector";
  Modelica.Thermal.HeatTransfer.Components.ThermalCollector theConSurBou(final m=
        nSurBou) if 
       haveSurBou 
    "Thermal collector to convert from vector to scalar connector";
public 
  Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a heaPorAir 
    "Heat port to air volume"; 
  final parameter Modelica.SIunits.TransmissionCoefficient tauIRSha_air[NConExtWin]=
    datConExtWin.glaSys.shade.tauIR_a 
    "Infrared transmissivity of shade for radiation coming from the exterior or the room";
        final parameter Modelica.SIunits.TransmissionCoefficient tauIRSha_glass[
                                                                          NConExtWin]=
    datConExtWin.glaSys.shade.tauIR_b 
    "Infrared transmissivity of shade for radiation coming from the glass";
  final parameter Boolean haveExteriorShade[NConExtWin]=
    {datConExtWin[i].glaSys.haveExteriorShade for i in 1:NConExtWin} 
    "Set to true if window has exterior shade (at surface a)";
  final parameter Boolean haveInteriorShade[NConExtWin]=
    {datConExtWin[i].glaSys.haveInteriorShade for i in 1:NConExtWin} 
    "Set to true if window has interior shade (at surface b)";
  final parameter Boolean haveShade = haveExteriorShade[1] or haveInteriorShade[1] 
    "Set to true if the windows have a shade";
  final parameter Real fFra[NConExtWin](each min=0, each max=1) = datConExtWin.fFra 
    "Fraction of window frame divided by total window area";
  parameter Boolean linearizeRadiation 
    "Set to true to linearize emissive power";
  parameter Modelica.SIunits.MassFlowRate m_flow_nominal(min=0) 
    "Nominal mass flow rate";
  HeatTransfer.Interfaces.RadiosityInflow JInSha[NConExtWin] if 
                                       haveShade 
    "Incoming radiosity that connects to shaded part of glass"; 
  HeatTransfer.Interfaces.RadiosityOutflow JOutSha[NConExtWin] if 
                                         haveShade 
    "Outgoing radiosity that connects to shaded part of glass";
  HeatTransfer.Interfaces.RadiosityInflow JInUns[NConExtWin] if 
     haveConExtWin "Incoming radiosity that connects to unshaded part of glass";
  HeatTransfer.Interfaces.RadiosityOutflow JOutUns[NConExtWin] if 
     haveConExtWin "Outgoing radiosity that connects to unshaded part of glass"; 
  SolarRadiationExchange solRadExc(
    final nConExt=nConExt,
    final nConExtWin=nConExtWin,
    final nConPar=nConPar,
    final nConBou=nConBou,
    final nSurBou=nSurBou,
    final datConExt = datConExt,
    final datConExtWin = datConExtWin,
    final datConPar = datConPar,
    final datConBou = datConBou,
    final surBou = surBou,
    final isFloorConExt=isFloorConExt,
    final isFloorConExtWin=isFloorConExtWin,
    final isFloorConPar_a=isFloorConPar_a,
    final isFloorConPar_b=isFloorConPar_b,
    final isFloorConBou=isFloorConBou,
    final isFloorSurBou=isFloorSurBou,
    final tauGla=tauGlaSol) if 
       haveConExtWin "Solar radiative heat exchange"; 
  InfraredRadiationGainDistribution irRadGai(
    final nConExt=nConExt,
    final nConExtWin=nConExtWin,
    final nConPar=nConPar,
    final nConBou=nConBou,
    final nSurBou=nSurBou,
    final datConExt = datConExt,
    final datConExtWin = datConExtWin,
    final datConPar = datConPar,
    final datConBou = datConBou,
    final surBou = surBou,
    final haveShade=haveShade) 
    "Distribution for infrared radiative heat gains (e.g., due to equipment and people)";
  InfraredRadiationExchange irRadExc(
    final nConExt=nConExt,
    final nConExtWin=nConExtWin,
    final nConPar=nConPar,
    final nConBou=nConBou,
    final nSurBou=nSurBou,
    final datConExt = datConExt,
    final datConExtWin = datConExtWin,
    final datConPar = datConPar,
    final datConBou = datConBou,
    final surBou = surBou,
    final linearizeRadiation = linearizeRadiation) 
    "Infrared radiative heat exchange"; 
  Modelica.Blocks.Interfaces.RealInput qGai_flow[3] 
    "Radiant, convective and latent heat input into room (positive if heat gain)"; 
  Modelica.Blocks.Interfaces.RealInput uSha[NConExtWin](each min=0, each max=1) if 
       haveShade 
    "Control signal for the shading device (removed if no shade is present)"; 
  Modelica.Blocks.Interfaces.RealInput QAbsSolSha_flow[NConExtWin](
    final unit="W", quantity="Power") if 
       haveConExtWin "Solar radiation absorbed by shade";
  Modelica.Blocks.Interfaces.RealInput JInConExtWin[NConExtWin](final unit="W",
      quantity="Power") if haveConExtWin 
    "Solar radiation transmitted from the outside through the glazing system"; 
  Modelica.Blocks.Interfaces.RealOutput HOutConExtWin[NConExtWin](unit="W/m2") if 
       haveConExtWin "Outgoing solar radiation that strikes window"; 
  HeatGain heaGai(redeclare package Medium = Medium, final AFlo=AFlo) 
    "Model to convert internal heat gains";
protected 
  RadiationAdapter radiationAdapter;
  Modelica.Blocks.Math.Add add;
public 
  RadiationTemperature radTem(
    final nConExt=nConExt,
    final nConExtWin=nConExtWin,
    final nConPar=nConPar,
    final nConBou=nConBou,
    final nSurBou=nSurBou,
    final datConExt=datConExt,
    final datConExtWin=datConExtWin,
    final datConPar=datConPar,
    final datConBou=datConBou,
    final surBou=surBou,
    final haveShade=haveShade) "Radiative temperature of the room";
  Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a heaPorRad 
    "Heat port for radiative heat gain and radiative temperature"; 
equation 
  connect(convConExt.solid, conExt);
  connect(convConPar_a.solid, conPar_a);
  connect(convConPar_b.solid, conPar_b);
  connect(convConBou.solid, conBou);
  connect(convSurBou.solid, conSurBou);
  connect(convConExt.fluid, theConConExt.port_a);
  connect(convConPar_a.fluid, theConConPar_a.port_a);
  connect(convConPar_b.fluid, theConConPar_b.port_a);
  connect(convConBou.fluid, theConConBou.port_a);
  connect(convSurBou.fluid, theConSurBou.port_a);
  connect(theConConExt.port_b,vol. heatPort);
  connect(theConConExtWin.port_b,vol. heatPort);
  connect(theConConPar_a.port_b,vol. heatPort);
  connect(theConConPar_b.port_b,vol. heatPort);
  connect(theConConBou.port_b,vol. heatPort);
  connect(theConSurBou.port_b,vol. heatPort);
  connect(vol.heatPort, heaPorAir);
  connect(convConWin.JOutUns, JOutUns);
  connect(JInUns, convConWin.JInUns);
  connect(convConWin.JOutSha, JOutSha);
  connect(convConWin.JInSha, JInSha);
  connect(uSha, convConWin.uSha);
  connect(conExt, irRadExc.conExt);
  connect(conExtWinFra, irRadExc.conExtWinFra);
  connect(conPar_a, irRadExc.conPar_a);
  connect(conPar_b, irRadExc.conPar_b); 
  connect(conBou, irRadExc.conBou); 
  connect(conSurBou, irRadExc.conSurBou);
  connect(irRadExc.JOutConExtWin, convConWin.JInRoo);
  connect(convConWin.JOutRoo, irRadExc.JInConExtWin);
  connect(irRadGai.JOutConExtWin, convConWin.JInRoo);
  connect(irRadGai.conExt, conExt);
  connect(irRadGai.conExtWinFra, conExtWinFra);
  connect(irRadGai.conPar_a, conPar_a);
  connect(irRadGai.conPar_b, conPar_b);
  connect(irRadGai.conBou, conBou);
  connect(irRadGai.conSurBou, conSurBou); 
  connect(irRadGai.uSha, uSha);
  connect(heaGai.QCon_flow,vol. heatPort);
  connect(heaGai.qGai_flow, qGai_flow);
  connect(convConExtWin.fluid, theConConExtWin.port_a);
  connect(convConExtWin.solid, conExtWin);
  connect(conExtWin, irRadExc.conExtWin);
  connect(conExtWin, irRadGai.conExtWin);
  connect(heaGai.QLat_flow,vol. ports[nPorts+1]); 
  for i in 1:nPorts loop
    connect(ports[i],vol. ports[i]); 
  end for;
  connect(convConWin.air, theConConWin.port_a);
  connect(theConConWin.port_b,vol. heatPort);
  connect(glaUns, convConWin.glaUns);
  connect(convConWin.glaSha, glaSha); 
  connect(conExt, solRadExc.conExt);
  connect(conExtWinFra, solRadExc.conExtWinFra);
  connect(conPar_a, solRadExc.conPar_a);
  connect(conPar_b, solRadExc.conPar_b);
  connect(conBou, solRadExc.conBou);
  connect(conSurBou, solRadExc.conSurBou);
  connect(conExtWin, solRadExc.conExtWin);
  connect(QAbsSolSha_flow, convConWin.QAbs_flow);
  connect(solRadExc.JInConExtWin, JInConExtWin);
  connect(solRadExc.HOutConExtWin,HOutConExtWin);
  connect(uSha, radTem.uSha);
  connect(conExt, radTem.conExt);
  connect(conExtWin, radTem.conExtWin);
  connect(conExtWinFra, radTem.conExtWinFra);
  connect(conPar_a, radTem.conPar_a);
  connect(conPar_b, radTem.conPar_b); 
  connect(conBou, radTem.conBou); 
  connect(conSurBou, radTem.conSurBou); 
  connect(radTem.glaUns, glaUns);
  connect(radTem.glaSha, glaSha);
  connect(convConWin.sha, radTem.sha); 
  connect(radTem.TRad, radiationAdapter.TRad);
  connect(radiationAdapter.rad, heaPorRad); 
  connect(radiationAdapter.QRad_flow, add.u1);
  connect(heaGai.QRad_flow, add.u2);
  connect(add.y, irRadGai.Q_flow);
  connect(conExtWinFra, convConWin.frame); 
end MixedAir;
 
 Buildings.Rooms.BaseClasses.RadiationAdapter
Buildings.Rooms.BaseClasses.RadiationAdapter
 
Extends from Buildings.BaseClasses.BaseIcon (Base icon).
| Type | Name | Description | 
|---|---|---|
| HeatPort_a | rad | Port for radiative heat gain and radiation temperature | 
| input RealInput | TRad | Radiation temperature of room | 
| output RealOutput | QRad_flow | Radiative heat gain | 
model RadiationAdapter "Model to connect between signals and heat port for radiative gains of the room" extends Buildings.BaseClasses.BaseIcon;Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a rad "Port for radiative heat gain and radiation temperature"; public Modelica.Blocks.Interfaces.RealInput TRad "Radiation temperature of room"; Modelica.Blocks.Interfaces.RealOutput QRad_flow "Radiative heat gain"; equation QRad_flow = rad.Q_flow; rad.T = TRad;end RadiationAdapter; 
 Buildings.Rooms.BaseClasses.RadiationTemperature
Buildings.Rooms.BaseClasses.RadiationTemperature
 
This model computes the radiative temperature in the room. For a room with windows but no shade, the radiative temperature is computed as
Trad = ∑i (Ai εi Ti) ⁄ ∑i (Ai εi)
where Trad is the radiative temperature of the room, Ai are the surface areas of the room, εi are the infrared emissivities of the surfaces, and Ti are the surface temperatures.If a the windows have a shade, then the equation is modified to take the actual shaded and non-shaded surface area into account. In this situation, the shaded part of a window has a infrared radiative power of
E = A ( u εs Ts + (1-u) εg τs Tgs)
where A is the surface area of the glass, u is the control signal of the shade, εs is the infrared absorptivity of the shade, Ts is the temperature of the shade, εg is the infrared absorptivity of the glass, τs is the infrared transmittance of the shade, and Tgs is the glass temperature behind the shade.E = A (1-u) εg Tgn
where Tgn is the glass temperature of the non-shaded part of the window.Extends from Buildings.Rooms.BaseClasses.PartialSurfaceInterface (Partial model that is used for infrared radiation balance).
| Type | Name | Default | Description | 
|---|---|---|---|
| ParameterConstruction | datConExt[NConExt] | Data for exterior construction | |
| ParameterConstructionWithWindow | datConExtWin[NConExtWin] | Data for exterior construction with window | |
| ParameterConstruction | datConPar[NConPar] | Data for partition construction | |
| ParameterConstruction | datConBou[NConBou] | Data for construction boundary | |
| Generic | surBou[NSurBou] | Record for data of surfaces whose heat conduction is modeled outside of this room | |
| Boolean | haveShade | Set to true if the windows have a shade | |
| Exterior constructions | |||
| Integer | nConExt | Number of exterior constructions | |
| Integer | nConExtWin | Number of window constructions | |
| Partition constructions | |||
| Integer | nConPar | Number of partition constructions | |
| Boundary constructions | |||
| Integer | nConBou | Number of constructions that have their outside surface exposed to the boundary of this room | |
| Integer | nSurBou | Number of surface heat transfer models that connect to constructions that are modeled outside of this room | |
| Type | Name | Description | 
|---|---|---|
| HeatPort_a | conExt[NConExt] | Heat port that connects to room-side surface of exterior constructions | 
| HeatPort_a | conExtWin[NConExtWin] | Heat port that connects to room-side surface of exterior constructions that contain a window | 
| HeatPort_a | conExtWinFra[NConExtWin] | Heat port that connects to room-side surface of window frame | 
| HeatPort_a | conPar_a[NConPar] | Heat port that connects to room-side surface a of partition constructions | 
| HeatPort_a | conPar_b[NConPar] | Heat port that connects to room-side surface b of partition constructions | 
| HeatPort_a | conBou[NConBou] | Heat port that connects to room-side surface of constructions that expose their other surface to the outside | 
| HeatPort_a | conSurBou[NSurBou] | Heat port to surfaces of models that compute the heat conduction outside of this room | 
| HeatPort_a | glaUns[NConExtWin] | Heat port that connects to room-side surface of unshaded glass | 
| HeatPort_a | glaSha[NConExtWin] | Heat port that connects to room-side surface of shaded glass | 
| HeatPort_a | sha[NConExtWin] | Heat port that connects to shade | 
| input RealInput | uSha[NConExtWin] | Control signal for the shading device (removed if no shade is present) | 
| output RealOutput | TRad | Radiative temperature [K] | 
model RadiationTemperature "Radiative temperature of the room" extends Buildings.Rooms.BaseClasses.PartialSurfaceInterface;Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a glaUns[NConExtWin] if haveConExtWin "Heat port that connects to room-side surface of unshaded glass"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a glaSha[NConExtWin] if haveShade "Heat port that connects to room-side surface of shaded glass"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a sha[NConExtWin] if haveShade "Heat port that connects to shade"; parameter Boolean haveShade "Set to true if the windows have a shade";Modelica.Blocks.Interfaces.RealInput uSha[NConExtWin](each min=0, each max=1) if haveShade "Control signal for the shading device (removed if no shade is present)"; Modelica.Blocks.Interfaces.RealOutput TRad(min=0, unit="K", displayUnit="degC") "Radiative temperature"; protected final parameter Integer NOpa = NConExt+2*NConExtWin+2*NConPar+NConBou+NSurBou "Number of opaque surfaces, including the window frame"; final parameter Integer NWin = NConExtWin "Number of window surfaces"; final parameter Integer NTot = NOpa + NWin "Total number of surfaces"; final parameter Modelica.SIunits.Area AGla[NWin] = datConExtWin.AGla "Surface area of opaque surfaces"; final parameter Real epsGla[NWin](min=0, max=1)= {datConExtWin[i].glaSys.glass[datConExtWin[i].glaSys.nLay].absIR_b for i in 1:NWin} "Absorptivity of glass"; final parameter Real epsSha[NWin](min=0, max=1)= {datConExtWin[i].glaSys.shade.absIR_a for i in 1:NWin} "Absorptivity of shade"; final parameter Real tauSha[NWin](min=0, max=1)= {(if datConExtWin[i].glaSys.haveInteriorShade then datConExtWin[i].glaSys.shade.tauIR_a else 1) for i in 1:NWin} "Transmissivity of shade"; final parameter Modelica.SIunits.Area epsAOpa[NOpa](fixed=false) "Product of area times absorptivity of opaque surfaces"; final parameter Modelica.SIunits.Area epsAGla[NWin](fixed=false) "Product of area times absorptivity of window surfaces"; final parameter Modelica.SIunits.Area epsASha[NWin](fixed=false) "Product of area times absorptivity of window shade"; final parameter Modelica.SIunits.Area epsTauASha[NWin](fixed=false) "Product of area times glass absorptivity times shade transmittance"; Modelica.SIunits.Temperature TOpa[NOpa](each start=293.15, each nominal=293.15) "Temperature of opaque surfaces"; Modelica.SIunits.Temperature TGlaUns[NWin](each start=293.15, each nominal=293.15) "Temperature of unshaded part of glass"; Modelica.SIunits.Temperature TGlaSha[NWin](each start=293.15, each nominal=293.15) "Temperature of shaded part of glass"; Modelica.SIunits.Temperature TSha[NWin](each start=293.15, each nominal=293.15) "Temperature of shade"; // Internal connectors, used because of the conditionally removed connectors Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a glaUns_internal[NConExtWin] "Heat port that connects to room-side surface of unshaded glass"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a glaSha_internal[NConExtWin] "Heat port that connects to room-side surface of shaded glass"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a sha_internal[NConExtWin] "Heat port that connects to shade"; Modelica.Blocks.Interfaces.RealInput uSha_internal[NConExtWin](each min=0, each max=1) "Control signal for the shading device"; initial equation // The next loops build the array epsAOpa that simplifies // the model equations. // The arrays stores the values of the constructios in the following order // [x[1:NConExt] x[1:NConPar] x[1: NConPar] x[1: NConBou] x[1: NSurBou] x[1: NConExtWin] x[1: NConExtWin]] // where x is epsOpa, AOpa or kOpa. // The last two entries are for the opaque wall that contains a window, and for the window frame. for i in 1:NConExt loop epsAOpa[i] = epsConExt[i] * AConExt[i]; end for; for i in 1:NConPar loop epsAOpa[i+NConExt] = epsConPar_a[i] * AConPar[i]; epsAOpa[i+NConExt+NConPar] = epsConPar_b[i] * AConPar[i]; end for; for i in 1:NConBou loop epsAOpa[i+NConExt+2*NConPar] = epsConBou[i] * AConBou[i]; end for; for i in 1:NSurBou loop epsAOpa[i+NConExt+2*NConPar+NConBou] = epsSurBou[i] * ASurBou[i]; end for; for i in 1:NConExtWin loop // Opaque part of construction that has a window embedded epsAOpa[i+NConExt+2*NConPar+NConBou+NSurBou] = epsConExtWinOpa[i] * AConExtWinOpa[i]; // Window frame epsAOpa[i+NConExt+2*NConPar+NConBou+NSurBou+NConExtWin] = epsConExtWinFra[i] * AConExtWinFra[i]; end for; // Window glass for i in 1:NConExtWin loop // Window glass epsAGla[i] = AGla[i] * epsGla[i]; // Window shade epsASha[i] = AGla[i] * epsSha[i]; // Emitted from glas and transmitted through window shade epsTauASha[i] = AGla[i] * epsGla[i] * tauSha[i]; end for; //////////////////////////////////////////////////////////////////// equation // Conditional connnector connect(glaUns, glaUns_internal); connect(glaSha, glaSha_internal); connect(sha, sha_internal); connect(uSha, uSha_internal); if not haveConExtWin then glaUns_internal.T = fill(293.15, NConExtWin); end if; if not haveShade then glaSha_internal.T = fill(293.15, NConExtWin); sha_internal.T = fill(293.15, NConExtWin); uSha_internal = fill(0, NConExtWin); end if; // Assign temperature of opaque surfaces for i in 1:NConExt loop TOpa[i] = conExt[i].T; end for; for i in 1:NConPar loop TOpa[i+NConExt] = conPar_a[i].T; TOpa[i+NConExt+NConPar] = conPar_b[i].T; end for; for i in 1:NConBou loop TOpa[i+NConExt+2*NConPar] = conBou[i].T; end for; for i in 1:NSurBou loop TOpa[i+NConExt+2*NConPar+NConBou] = conSurBou[i].T; end for; for i in 1:NConExtWin loop TOpa[i+NConExt+2*NConPar+NConBou+NSurBou] = conExtWin[i].T; TOpa[i+NConExt+2*NConPar+NConBou+NConExtWin+NSurBou] = conExtWinFra[i].T; end for; // Assign temperature of glass and shade for i in 1:NConExtWin loop TGlaUns[i] = glaUns_internal[i].T; TGlaSha[i] = glaSha_internal[i].T; TSha[i] = sha_internal[i].T; end for; // Compute radiative temperature if haveShade then TRad = (sum(epsAOpa[i] * TOpa[i] for i in 1:NOpa) + sum( ( uSha_internal[i] * (epsASha[i] * TSha[i] + epsTauASha[i] * TGlaSha[i]) + (1-uSha_internal[i]) * epsAGla[i] * TGlaUns[i]) for i in 1:NConExtWin)) / (sum(epsAOpa) + sum( ( uSha_internal[i] * (epsASha[i] + epsTauASha[i]) + (1-uSha_internal[i]) * epsAGla[i]) for i in 1:NConExtWin)); else TRad = (sum(epsAOpa[i] * TOpa[i] for i in 1:NOpa) + sum(epsAGla .* TGlaUns)) / (sum(epsAOpa) + sum(epsAGla)); end if; // Assign heat exchange to connectors for i in 1:NConExt loop 0 = conExt[i].Q_flow; end for; for i in 1:NConPar loop 0 = conPar_a[i].Q_flow; 0 = conPar_b[i].Q_flow; end for; for i in 1:NConBou loop 0 = conBou[i].Q_flow; end for; for i in 1:NSurBou loop 0 = conSurBou[i].Q_flow; end for; for i in 1:NConExtWin loop 0 = conExtWin[i].Q_flow; 0 = conExtWinFra[i].Q_flow; end for; /* for i in 1:NConExtWin loop 0 = glaUns_internal[i].Q_flow; 0 = glaSha_internal[i].Q_flow; 0 = sha_internal[i].Q_flow; end for; */end RadiationTemperature; 
 Buildings.Rooms.BaseClasses.SolarRadiationExchange
Buildings.Rooms.BaseClasses.SolarRadiationExchange
 
This model computes the distribution of the solar radiation gain to the room surfaces. Let Nw denote the number of windows, Nf the number of floor elements and Nn the number of non-floor elements such as ceiling, wall and window elements. Input to the model are the solar radiosities Ji, i ∈ {1, … , Nw}, that were transmitted through the window. The total incoming solar radiation is therefore
H = ∑i=1Nw Jini
It is assumed that H first hits the floor where some of it is absorbed, and some of it is diffusely reflected to all other surfaces. Only the first reflection is taken into account and the location of the floor patch relative to the window is neglected.
Hence, the radiation that is absorbed by each floor patch i ∈ {1, …, Nf}, and may be partially transmitted in the unusual case that the floor contains a window, is
Qi = H (εi+τi) Ai ⁄ ∑j=1Nf Aj.
The sum of the radiation that is reflected by the floor is thereforeJf = H ∑i=1Nf (1-εi-τi) Ai ⁄ ∑j=1Nf Aj.
This reflected radiosity is then distributed to all non-floor areas i ∈ {1, …, Nn} using
Qi = Jf Ai (εi+τi) ⁄ ∑k=1Nn Ak (εk+τk)
For opaque surfaces, the heat flow rate Qi is set to be equal to the heat flow rate at the heat port. For the glass of the windows, the heat flow rate Qi is set to the radiosity Jouti that will strike the glass or the window shade as diffuse solar radiation.
Extends from Buildings.Rooms.BaseClasses.PartialSurfaceInterface (Partial model that is used for infrared radiation balance).
| Type | Name | Default | Description | 
|---|---|---|---|
| ParameterConstruction | datConExt[NConExt] | Data for exterior construction | |
| ParameterConstructionWithWindow | datConExtWin[NConExtWin] | Data for exterior construction with window | |
| ParameterConstruction | datConPar[NConPar] | Data for partition construction | |
| ParameterConstruction | datConBou[NConBou] | Data for construction boundary | |
| Generic | surBou[NSurBou] | Record for data of surfaces whose heat conduction is modeled outside of this room | |
| Boolean | isFloorConExt[NConExt] | Flag to indicate if floor for exterior constructions | |
| Boolean | isFloorConExtWin[NConExtWin] | Flag to indicate if floor for constructions | |
| Boolean | isFloorConPar_a[NConPar] | Flag to indicate if floor for constructions | |
| Boolean | isFloorConPar_b[NConPar] | Flag to indicate if floor for constructions | |
| Boolean | isFloorConBou[NConBou] | Flag to indicate if floor for constructions with exterior boundary conditions exposed to outside of room model | |
| Boolean | isFloorSurBou[NSurBou] | Flag to indicate if floor for constructions that are modeled outside of this room | |
| Emissivity | tauGla[NConExtWin] | Transmissivity of window [1] | |
| Exterior constructions | |||
| Integer | nConExt | Number of exterior constructions | |
| Integer | nConExtWin | Number of window constructions | |
| Partition constructions | |||
| Integer | nConPar | Number of partition constructions | |
| Boundary constructions | |||
| Integer | nConBou | Number of constructions that have their outside surface exposed to the boundary of this room | |
| Integer | nSurBou | Number of surface heat transfer models that connect to constructions that are modeled outside of this room | |
| Type | Name | Description | 
|---|---|---|
| HeatPort_a | conExt[NConExt] | Heat port that connects to room-side surface of exterior constructions | 
| HeatPort_a | conExtWin[NConExtWin] | Heat port that connects to room-side surface of exterior constructions that contain a window | 
| HeatPort_a | conExtWinFra[NConExtWin] | Heat port that connects to room-side surface of window frame | 
| HeatPort_a | conPar_a[NConPar] | Heat port that connects to room-side surface a of partition constructions | 
| HeatPort_a | conPar_b[NConPar] | Heat port that connects to room-side surface b of partition constructions | 
| HeatPort_a | conBou[NConBou] | Heat port that connects to room-side surface of constructions that expose their other surface to the outside | 
| HeatPort_a | conSurBou[NSurBou] | Heat port to surfaces of models that compute the heat conduction outside of this room | 
| output RealOutput | HOutConExtWin[NConExtWin] | Outgoing solar radiation that strikes window per unit area [W/m2] | 
| input RealInput | JInConExtWin[NConExtWin] | Solar radiation transmitted by window per unit area [W] | 
model SolarRadiationExchange 
  "Solar radiation heat exchange between the room facing surfaces"
  extends Buildings.Rooms.BaseClasses.PartialSurfaceInterface(
  final epsConExt = datConExt.layers.absSol_b,
  final epsConExtWinOpa = datConExtWin.layers.absSol_b,
  final epsConExtWinUns={(1-datConExtWin[i].glaSys.glass[datConExtWin[i].glaSys.nLay].tauSol
                     -datConExtWin[i].glaSys.glass[datConExtWin[i].glaSys.nLay].rhoSol_b) for i in 1:NConExtWin},
  final epsConExtWinSha = {(1-datConExtWin[i].glaSys.glass[datConExtWin[i].glaSys.nLay].tauSol
                       -datConExtWin[i].glaSys.glass[datConExtWin[i].glaSys.nLay].rhoSol_b) for i in 1:NConExtWin},
  final epsConExtWinFra = datConExtWin.glaSys.absSolFra,
  final epsConPar_a = datConPar.layers.absSol_a,
  final epsConPar_b = datConPar.layers.absSol_b,
  final epsConBou = datConBou.layers.absSol_b,
  final epsSurBou = surBou.absSol);
  // In the above declaration, we simplified the assignment of epsConExtWinSha.
  // An exact formulation would need to take into account the transmission and reflection
  // of the shade for the solar radiation that strikes the window from the room-side.
  // The simplification leads to too low a value of epsConExtWinSha. Since epsConExtWinSha
  // is used as a weight for how much solar radiation hits the window from the room-side,
  // underestimating epsConExtWinSha does not seem to cause concerns. The reason is that
  // the model assumes diffuse reflection, whereas in reality, reflection of the solar
  // radiation at the floor is likely specular, and therefore less radiation would hit
  // the window from the room-side.
  parameter Boolean isFloorConExt[NConExt] 
    "Flag to indicate if floor for exterior constructions";
  parameter Boolean isFloorConExtWin[NConExtWin] 
    "Flag to indicate if floor for constructions";
  parameter Boolean isFloorConPar_a[NConPar] 
    "Flag to indicate if floor for constructions";
  parameter Boolean isFloorConPar_b[NConPar] 
    "Flag to indicate if floor for constructions";
  parameter Boolean isFloorConBou[NConBou] 
    "Flag to indicate if floor for constructions with exterior boundary conditions exposed to outside of room model";
  parameter Boolean isFloorSurBou[NSurBou] 
    "Flag to indicate if floor for constructions that are modeled outside of this room";
  parameter Modelica.SIunits.Emissivity tauGla[NConExtWin] 
    "Transmissivity of window";
  Modelica.Blocks.Interfaces.RealOutput HOutConExtWin[NConExtWin](unit="W/m2") 
    "Outgoing solar radiation that strikes window per unit area";
  Modelica.Blocks.Interfaces.RealInput JInConExtWin[NConExtWin](unit="W") 
    "Solar radiation transmitted by window per unit area"; 
  Modelica.SIunits.HeatFlowRate JOutConExtWin[NConExtWin] 
    "Outgoing solar radiation that strikes the window";
  Modelica.SIunits.HeatFlowRate HTot 
    "Total solar radiation that enters the room";
protected 
  final parameter Real k1(unit="1", fixed=false) 
    "Intermediate variable for gain for solar radiation distribution";
  final parameter Real k2(fixed=false) 
    "Intermediate variable for gain for solar radiation distribution";
  Modelica.SIunits.HeatFlowRate Q_flow[NTot] 
    "Total solar radiation that is absorbed by the surfaces (or transmitted back through the glass)";
  final parameter Integer NOpa = NConExt+2*NConExtWin+2*NConPar+NConBou+NSurBou 
    "Number of opaque surfaces, including the window frame";
  final parameter Integer NWin = NConExtWin "Number of window surfaces";
  final parameter Integer NTot = NOpa + NWin "Total number of surfaces";
  final parameter Boolean isFlo[NTot](fixed=false) 
    "Flag, true if a surface is a floor";
  final parameter Real eps[NTot](min=0, max=1, fixed=false) 
    "Solar absorptivity";
  final parameter Real tau[NTot](min=0, max=1, fixed=false) 
    "Solar transmissivity";
  final parameter Modelica.SIunits.Area AFlo(fixed=false) "Total floor area";
  final parameter Modelica.SIunits.Area A[NTot](fixed=false) "Surface areas";
  final parameter Real k[NTot](unit="1", fixed=false) 
    "Gain for solar radiation distribution";
initial equation 
  // The next loops builds arrays that simplify
  // the model equations.
  // These arrays store the values of the constructios in the following order
  // [x[1:NConExt] x[1:NConPar] x[1: NConPar] x[1: NConBou] x[1: NSurBou] x[1: NConExtWin] x[1: NConExtWin]]
  // where x is epsOpa, AOpa or kOpa.
  // The last two entries are for the opaque wall that contains a window, and for the window frame.
  for i in 1:NConExt loop
    eps[i] = epsConExt[i];
    A[i]      = AConExt[i];
    isFlo[i]  = isFloorConExt[i];
  end for;
  for i in 1:NConPar loop
    eps[i+NConExt]           = epsConPar_a[i];
    A[i+NConExt]             = AConPar[i];
    isFlo[i+NConExt]         = isFloorConPar_a[i];
    eps[i+NConExt+NConPar]   = epsConPar_b[i];
    A[i+NConExt+NConPar]     = AConPar[i];
    isFlo[i+NConExt+NConPar] = isFloorConPar_b[i];
  end for;
  for i in 1:NConBou loop
    eps[i+NConExt+2*NConPar]   = epsConBou[i];
    A[i+NConExt+2*NConPar]     = AConBou[i];
    isFlo[i+NConExt+2*NConPar] = isFloorConBou[i];
  end for;
  for i in 1:NSurBou loop
    eps[i+NConExt+2*NConPar+NConBou]   = epsSurBou[i];
    A[i+NConExt+2*NConPar+NConBou]     = ASurBou[i];
    isFlo[i+NConExt+2*NConPar+NConBou] = isFloorSurBou[i];
  end for;
  for i in 1:NConExtWin loop
    // Opaque part of construction that has a window embedded
    eps[i+NConExt+2*NConPar+NConBou+NSurBou]   = epsConExtWinOpa[i];
    A[i+NConExt+2*NConPar+NConBou+NSurBou]     = AConExtWinOpa[i];
    isFlo[i+NConExt+2*NConPar+NConBou+NSurBou] = isFloorConExtWin[i];
    // Window frame
    eps[i+NConExt+2*NConPar+NConBou+NSurBou+NConExtWin]   = epsConExtWinFra[i];
    A[i+NConExt+2*NConPar+NConBou+NSurBou+NConExtWin]     = AConExtWinFra[i];
    isFlo[i+NConExt+2*NConPar+NConBou+NSurBou+NConExtWin] = isFloorConExtWin[i];
  end for;
  // Window glass
  for i in 1:NConExtWin loop
    // We simplify and assume that the shaded and unshaded part of the window
    // have the same solar absorbtance.
    // This simplification allows lumping the solar distribution into
    // a parameter.
    eps[i+NConExt+2*NConPar+NConBou+NSurBou+2*NConExtWin] = epsConExtWinUns[i];
    isFlo[i+NConExt+2*NConPar+NConBou+NSurBou+2*NConExtWin] = isFloorConExtWin[i];
    A[i+NConExt+2*NConPar+NConBou+NSurBou+2*NConExtWin] = AConExtWinGla[i];
  end for;
  // Vector with all surface areas.
  // The next loops build the array A that simplifies
  // the model equations.
  // These array stores the values of the constructios in the following order
  // [AOpa[1:NConExt] AOpa[1:NConPar] AOpa[1: NConPar] AOpa[1: NConBou] AOpa[1: NSurBou]
  //  AOpa[1: NConExtWin] AOpa[1: NConExtWin] AGla[1: NConExtWin]]
  // since NWin=NConExtWin.
  // Solar transmissivity
  for i in 1:NOpa loop
    tau[i] = 0;
  end for;
  for i in 1:NWin loop
    tau[NOpa+i] = tauGla[i];
  end for;
  // Sum of surface areas
  AFlo = sum( (if isFlo[i] then A[i] else 0) for i in 1:NTot);
  // Coefficient that is used for non-floor areas.
  // The expression  max(1E-20, AFlo) is used to prevent a division by zero in case AFlo=0.
  // The situation for AFlo=0 is caught by the assert statement.
  k1 = sum( ( if isFlo[i] then (A[i] * (1-eps[i]-tau[i])) else 0)  for i in 1:NTot) / max(1E-20, AFlo);
  k2 = sum( ( if isFlo[i] then 0 else (A[i] * (eps[i]+tau[i])))  for i in 1:NTot);
  if ( k2 > 1E-10) then
    for i in 1:NTot loop
      if isFlo[i] then
         k[i] = (eps[i]+tau[i]) * A[i] / AFlo;
      else
         k[i] = k1/k2*(eps[i]+tau[i]) * A[i];
      end if;
     end for;
  else
        // This branch only happens if k2=0, i.e., there is no surface other than floors
    for i in 1:NTot loop
      if isFlo[i] then
        k[i] = A[i] / AFlo;
      else
        k[i] = 0;
      end if;
    end for;
  end if;
  // Test whether there is a floor inside this room
  assert( AFlo > 1E-10,
     "Error in parameters of the room model: The geometry is incorrect:\n" +
     "    The room model must have a construction that is a floor,\n" +
     "    and this construction must not have a window.\n" +
     "    The parameters for the room model are such that there is no such construction.\n" +
     "    Revise the model parameters.");
  // Test whether the distribution factors add up to one
  assert(abs(1-sum(k)) < 1E-5,
     "Program error: Sum of solar distribution factors in room is not equal to one. k=" + String(sum(k)));
////////////////////////////////////////////////////////////////////
equation 
  // Incoming radiation
  HTot = sum(JInConExtWin);
  // Radiation that is absorbed by the surfaces
  Q_flow = -k .* HTot;
  // Assign heat exchange to connectors
  for i in 1:NConExt loop
    Q_flow[i] = conExt[i].Q_flow;
  end for;
  for i in 1:NConPar loop
    Q_flow[i+NConExt]         = conPar_a[i].Q_flow;
    Q_flow[i+NConExt+NConPar] = conPar_b[i].Q_flow;
  end for;
  for i in 1:NConBou loop
    Q_flow[i+NConExt+2*NConPar] = conBou[i].Q_flow;
  end for;
  for i in 1:NSurBou loop
    Q_flow[i+NConExt+2*NConPar+NConBou] = conSurBou[i].Q_flow;
  end for;
  for i in 1:NConExtWin loop
    Q_flow[i+NConExt+2*NConPar+NConBou+NSurBou]            = conExtWin[i].Q_flow;
    Q_flow[i+NConExt+2*NConPar+NConBou+NSurBou+NConExtWin] = conExtWinFra[i].Q_flow;
  end for;
  // Windows
  for j in 1:NWin loop
    Q_flow[j+NOpa] = JOutConExtWin[j];
    HOutConExtWin[j] = if (AConExtWinGla[j] > 1E-10) then JOutConExtWin[j] / AConExtWinGla[j] else 0;
  end for;
end SolarRadiationExchange;
 
 Buildings.Rooms.BaseClasses.PartialSurfaceInterface
Buildings.Rooms.BaseClasses.PartialSurfaceInterface
 
There are also parameters that contain the number of constructions,
such as the number of exterior constructions nConExt. 
This parameter may take on the value 0. 
If this parameter were to be used to declare the size of vectors of
component models, then there may be vectors with zero components.
This can cause problems in Dymola 7.4. 
We therefore also introduced the parameter
NConExt = max(1, nConExt)which can be used to set the size of the vector of component models.
There are also parameters that can be used to conditionally remove components,
such as haveConExt, which is set to 
haveConExt = nConExt > 0;
Extends from Buildings.Rooms.BaseClasses.ConstructionRecords (Data records for construction data).
| Type | Name | Default | Description | 
|---|---|---|---|
| ParameterConstruction | datConExt[NConExt] | Data for exterior construction | |
| ParameterConstructionWithWindow | datConExtWin[NConExtWin] | Data for exterior construction with window | |
| ParameterConstruction | datConPar[NConPar] | Data for partition construction | |
| ParameterConstruction | datConBou[NConBou] | Data for construction boundary | |
| Generic | surBou[NSurBou] | Record for data of surfaces whose heat conduction is modeled outside of this room | |
| Exterior constructions | |||
| Integer | nConExt | Number of exterior constructions | |
| Integer | nConExtWin | Number of window constructions | |
| Partition constructions | |||
| Integer | nConPar | Number of partition constructions | |
| Boundary constructions | |||
| Integer | nConBou | Number of constructions that have their outside surface exposed to the boundary of this room | |
| Integer | nSurBou | Number of surface heat transfer models that connect to constructions that are modeled outside of this room | |
| Type | Name | Description | 
|---|---|---|
| HeatPort_a | conExt[NConExt] | Heat port that connects to room-side surface of exterior constructions | 
| HeatPort_a | conExtWin[NConExtWin] | Heat port that connects to room-side surface of exterior constructions that contain a window | 
| HeatPort_a | conExtWinFra[NConExtWin] | Heat port that connects to room-side surface of window frame | 
| HeatPort_a | conPar_a[NConPar] | Heat port that connects to room-side surface a of partition constructions | 
| HeatPort_a | conPar_b[NConPar] | Heat port that connects to room-side surface b of partition constructions | 
| HeatPort_a | conBou[NConBou] | Heat port that connects to room-side surface of constructions that expose their other surface to the outside | 
| HeatPort_a | conSurBou[NSurBou] | Heat port to surfaces of models that compute the heat conduction outside of this room | 
partial model PartialSurfaceInterface "Partial model that is used for infrared radiation balance" import Buildings; extends Buildings.Rooms.BaseClasses.ConstructionRecords;Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a conExt[NConExt] "Heat port that connects to room-side surface of exterior constructions"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a conExtWin[NConExtWin] "Heat port that connects to room-side surface of exterior constructions that contain a window"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a conExtWinFra[NConExtWin] "Heat port that connects to room-side surface of window frame"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a conPar_a[NConPar] "Heat port that connects to room-side surface a of partition constructions"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a conPar_b[NConPar] "Heat port that connects to room-side surface b of partition constructions"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a conBou[NConBou] "Heat port that connects to room-side surface of constructions that expose their other surface to the outside"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a conSurBou[NSurBou] "Heat port to surfaces of models that compute the heat conduction outside of this room"; protected final parameter Modelica.SIunits.Area AConExt[NConExt] = datConExt.A "Areas of exterior constructions"; final parameter Modelica.SIunits.Area AConExtWinOpa[NConExtWin] = datConExtWin.AOpa "Opaque areas of exterior construction that have a window"; final parameter Modelica.SIunits.Area AConExtWinGla[NConExtWin] = (1 .- datConExtWin.fFra) .* datConExtWin.AWin "Glass areas of exterior construction that have a window"; final parameter Modelica.SIunits.Area AConExtWinFra[NConExtWin] = datConExtWin.fFra .* datConExtWin.AWin "Frame areas of exterior construction that have a window"; final parameter Modelica.SIunits.Area AConPar[NConPar] = datConPar.A "Areas of partition constructions"; final parameter Modelica.SIunits.Area AConBou[NConBou] = datConBou.A "Areas of constructions with exterior boundary conditions exposed to outside of room model"; final parameter Modelica.SIunits.Area ASurBou[NSurBou] = surBou.A "Area of surface models of constructions that are modeled outside of this room"; parameter Modelica.SIunits.Emissivity epsConExt[NConExt] = datConExt.layers.absIR_b "Absorptivity of exterior constructions"; parameter Modelica.SIunits.Emissivity epsConExtWinOpa[NConExtWin] = datConExtWin.layers.absIR_b "Absorptivity of opaque part of exterior constructions that contain a window"; parameter Modelica.SIunits.Emissivity epsConExtWinUns[NConExtWin]= {(datConExtWin[i].glaSys.glass[datConExtWin[i].glaSys.nLay].absIR_b) for i in 1:NConExtWin} "Absorptivity of unshaded part of window constructions"; parameter Modelica.SIunits.Emissivity epsConExtWinSha[NConExtWin] = datConExtWin.glaSys.shade.absIR_a "Absorptivity of shaded part of window constructions"; parameter Modelica.SIunits.Emissivity epsConExtWinFra[NConExtWin] = datConExtWin.glaSys.absIRFra "Absorptivity of window frame"; parameter Modelica.SIunits.Emissivity epsConPar_a[NConPar] = datConPar.layers.absIR_a "Absorptivity of partition constructions surface a"; parameter Modelica.SIunits.Emissivity epsConPar_b[NConPar] = datConPar.layers.absIR_b "Absorptivity of partition constructions surface b"; parameter Modelica.SIunits.Emissivity epsConBou[NConBou] = datConBou.layers.absIR_b "Absorptivity of constructions with exterior boundary conditions exposed to outside of room model"; parameter Modelica.SIunits.Emissivity epsSurBou[NSurBou] = surBou.absIR "Absorptivity of surface models of constructions that are modeled outside of this room";end PartialSurfaceInterface; 
 Buildings.Rooms.BaseClasses.SkyRadiationExchange
Buildings.Rooms.BaseClasses.SkyRadiationExchange
 
Extends from Buildings.BaseClasses.BaseIcon (Base icon).
| Type | Name | Default | Description | 
|---|---|---|---|
| Integer | n | Number of constructions | |
| Area | A[n] | Area of exterior constructions [m2] | |
| Real | vieFacSky[n] | View factor to sky (=1 for roofs) | |
| Emissivity | absIR[n] | Infrared absorptivity of building surface [1] | 
| Type | Name | Description | 
|---|---|---|
| HeatPort_a | port[n] | Heat port | 
| input RealInput | TOut | Outside air temperature [K] | 
| input RealInput | TBlaSky | Black body sky temperature [K] | 
model SkyRadiationExchange 
  "Radiative heat exchange with the sky and the ambient"
  import Buildings;
  extends Buildings.BaseClasses.BaseIcon;
  parameter Integer n(min=1) "Number of constructions";
   parameter Modelica.SIunits.Area A[n] "Area of exterior constructions";
  parameter Real vieFacSky[n](min=0, max=1) "View factor to sky (=1 for roofs)";
  parameter Modelica.SIunits.Emissivity absIR[n] 
    "Infrared absorptivity of building surface";
  Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a port[n] "Heat port";
  Modelica.Blocks.Interfaces.RealInput TOut(final quantity="ThermodynamicTemperature",
                                            final unit = "K", min=0) 
    "Outside air temperature";
  Modelica.Blocks.Interfaces.RealInput TBlaSky(
    final quantity="ThermodynamicTemperature",
    final unit="K",
    min=0) "Black body sky temperature"; 
protected 
  parameter Real k[n](unit="W/K4") = {4*A[i]*Modelica.Constants.sigma*absIR[i] for i in 1:n} 
    "Constant for radiative heat exchange";
  Modelica.SIunits.Temperature TEnv[n] "Environment temperature";
  Real TBlaSky4(unit="K4") "Auxiliary variable for radiative heat exchange";
  Real TOut4(unit="K4") "Auxiliary variable for radiative heat exchange";
  Modelica.SIunits.CoefficientOfHeatTransfer h[n] 
    "Radiative heat transfer coefficient";
equation 
  TBlaSky4 = TBlaSky^4;
  TOut4 = TOut^4;
  for i in 1:n loop
    TEnv[i] = (vieFacSky[i] * TBlaSky4 + (1-vieFacSky[i]) * TOut4)^(0.25);
    // h[i] uses TEnv[i] instead of (port[i].T+TEnv[i])/2 to avoid
    // a nonlinear equation system
    h[i]  = k[i] * TEnv[i]^3;
    port[i].Q_flow = h[i] * (port[i].T-TEnv[i]);
  end for;
end SkyRadiationExchange;
 
Record that defines the number of constructions that are used in the room model.
| Type | Name | Default | Description | 
|---|---|---|---|
| Exterior constructions | |||
| Integer | nConExt | Number of exterior constructions | |
| Integer | nConExtWin | Number of window constructions | |
| Partition constructions | |||
| Integer | nConPar | Number of partition constructions | |
| Boundary constructions | |||
| Integer | nConBou | Number of constructions that have their outside surface exposed to the boundary of this room | |
| Integer | nSurBou | Number of surface heat transfer models that connect to constructions that are modeled outside of this room | |
record ConstructionNumbers "Data records for construction data"
  ////////////////////////////////////////////////////////////////////////
  // Number of constructions and surface areas
  parameter Integer nConExt(min=0) "Number of exterior constructions";
  parameter Integer nConExtWin(min=0) "Number of window constructions";
  parameter Integer nConPar(min=0) "Number of partition constructions";
  parameter Integer nConBou(min=0) 
    "Number of constructions that have their outside surface exposed to the boundary of this room";
  parameter Integer nSurBou(min=0) 
    "Number of surface heat transfer models that connect to constructions that are modeled outside of this room";
  // Dimensions of components and connectors
protected 
  parameter Integer NConExt(min=1)=max(1, nConExt) 
    "Number of elements for exterior constructions";
  parameter Integer NConExtWin(min=1)=max(1, nConExtWin) 
    "Number of elements for exterior constructions with windows";
  parameter Integer NConPar(min=1)=max(1, nConPar) 
    "Number of elements for partition constructions";
  parameter Integer NConBou(min=1)=max(1, nConBou) 
    "Number of elements for constructions that have their outside surface exposed to the boundary of this room";
  parameter Integer NSurBou(min=1)=max(1, nSurBou) 
    "Number of elements for surface heat transfer models that connect to constructions that are modeled outside of this room";
  // Flags to conditionally remove components
  final parameter Boolean haveConExt = nConExt > 0 
    "Flag to conditionally remove components";
  final parameter Boolean haveConExtWin = nConExtWin > 0 
    "Flag to conditionally remove components";
  final parameter Boolean haveConPar = nConPar > 0 
    "Flag to conditionally remove components";
  final parameter Boolean haveConBou = nConBou > 0 
    "Flag to conditionally remove components";
  final parameter Boolean haveSurBou = nSurBou > 0 
    "Flag to conditionally remove components";
end ConstructionNumbers;
 
 
Record that defines the number of constructions that are used in the room model.
Extends from Buildings.Rooms.BaseClasses.ConstructionNumbers (Data records for construction data).
| Type | Name | Default | Description | 
|---|---|---|---|
| ParameterConstruction | datConExt[NConExt] | Data for exterior construction | |
| ParameterConstructionWithWindow | datConExtWin[NConExtWin] | Data for exterior construction with window | |
| ParameterConstruction | datConPar[NConPar] | Data for partition construction | |
| ParameterConstruction | datConBou[NConBou] | Data for construction boundary | |
| Generic | surBou[NSurBou] | Record for data of surfaces whose heat conduction is modeled outside of this room | |
| Exterior constructions | |||
| Integer | nConExt | Number of exterior constructions | |
| Integer | nConExtWin | Number of window constructions | |
| Partition constructions | |||
| Integer | nConPar | Number of partition constructions | |
| Boundary constructions | |||
| Integer | nConBou | Number of constructions that have their outside surface exposed to the boundary of this room | |
| Integer | nSurBou | Number of surface heat transfer models that connect to constructions that are modeled outside of this room | |
record ConstructionRecords "Data records for construction data" extends Buildings.Rooms.BaseClasses.ConstructionNumbers;parameter ParameterConstruction datConExt[NConExt]( each A=0, redeclare Buildings.HeatTransfer.Data.OpaqueConstructions.Brick120 layers, each til=0, each azi=0) "Data for exterior construction"; parameter Buildings.Rooms.BaseClasses.ParameterConstructionWithWindow datConExtWin[NConExtWin]( each A=0, redeclare Buildings.HeatTransfer.Data.OpaqueConstructions.Brick120 layers, each til=0, each azi=0, each AWin=0, redeclare Buildings.HeatTransfer.Data.GlazingSystems.SingleClear3 glaSys) "Data for exterior construction with window"; parameter Buildings.Rooms.BaseClasses.ParameterConstruction datConPar[ NConPar]( each A=0, redeclare Buildings.HeatTransfer.Data.OpaqueConstructions.Brick120 layers, each til=0, each azi=0) "Data for partition construction"; parameter Buildings.Rooms.BaseClasses.ParameterConstruction datConBou[NConBou]( each A=0, redeclare Buildings.HeatTransfer.Data.OpaqueConstructions.Brick120 layers, each til=0, each azi=0) "Data for construction boundary"; parameter Buildings.HeatTransfer.Data.OpaqueSurfaces.Generic surBou[ NSurBou](each A=0, each til=0) "Record for data of surfaces whose heat conduction is modeled outside of this room"; end ConstructionRecords; 
 Buildings.Rooms.BaseClasses.ParameterConstructionWithWindow
Buildings.Rooms.BaseClasses.ParameterConstructionWithWindow
 
This data record is used to set the parameters of constructions that do have a window.
The surface azimuth is defined in Buildings.HeatTransfer.Types.Azimuth and the surface tilt is defined in Buildings.HeatTransfer.Types.Tilt
Extends from Buildings.Rooms.BaseClasses.PartialParameterConstruction (Partial record for constructions).
| Type | Name | Default | Description | 
|---|---|---|---|
| Angle | til | Surface tilt [rad] | |
| Angle | azi | Surface azimuth [rad] | |
| Area | A | Heat transfer area of opaque construction and window combined [m2] | |
| Opaque construction | |||
| Generic | layers | redeclare parameter Building... | Material properties of opaque construction | 
| Initialization | |||
| Boolean | steadyStateInitial | false | =true initializes dT(0)/dt=0, false initializes T(0) at fixed temperature using T_a_start and T_b_start | 
| Temperature | T_a_start | 293.15 | Initial temperature at port_a, used if steadyStateInitial = false [K] | 
| Temperature | T_b_start | 293.15 | Initial temperature at port_b, used if steadyStateInitial = false [K] | 
| Glazing system | |||
| Area | AWin | Heat transfer area of window [m2] | |
| Real | fFra | 0.1 | Fraction of window frame divided by total window area | 
| Boolean | linearizeRadiation | true | Set to true to linearize emissive power | 
| Generic | glaSys | redeclare parameter HeatTran... | Material properties of glazing system | 
record ParameterConstructionWithWindow 
  "Record for exterior constructions that have a window"
  extends Buildings.Rooms.BaseClasses.PartialParameterConstruction;
  parameter Modelica.SIunits.Area A 
    "Heat transfer area of opaque construction and window combined";
  parameter Modelica.SIunits.Area AWin "Heat transfer area of window";
  final parameter Modelica.SIunits.Area AOpa = A-AWin 
    "Heat transfer area of opaque construction";
  parameter Real fFra(
    min=0,
    max=1) = 0.1 "Fraction of window frame divided by total window area";
  final parameter Modelica.SIunits.Area AFra = fFra*AWin "Frame area";
  final parameter Modelica.SIunits.Area AGla=AWin - AFra "Glass area";
  parameter Boolean linearizeRadiation = true 
    "Set to true to linearize emissive power";
 replaceable parameter HeatTransfer.Data.GlazingSystems.Generic glaSys 
    "Material properties of glazing system";
end ParameterConstructionWithWindow;
 
 Buildings.Rooms.BaseClasses.ParameterConstruction
Buildings.Rooms.BaseClasses.ParameterConstruction
This data record is used to set the parameters of constructions that do not have a window.
The surface azimuth is defined in Buildings.HeatTransfer.Types.Azimuth and the surface tilt is defined in Buildings.HeatTransfer.Types.Tilt
Extends from Buildings.Rooms.BaseClasses.PartialParameterConstruction (Partial record for constructions).
| Type | Name | Default | Description | 
|---|---|---|---|
| Angle | til | Surface tilt [rad] | |
| Angle | azi | Surface azimuth [rad] | |
| Area | A | Heat transfer area [m2] | |
| Opaque construction | |||
| Generic | layers | redeclare parameter Building... | Material properties of opaque construction | 
| Initialization | |||
| Boolean | steadyStateInitial | false | =true initializes dT(0)/dt=0, false initializes T(0) at fixed temperature using T_a_start and T_b_start | 
| Temperature | T_a_start | 293.15 | Initial temperature at port_a, used if steadyStateInitial = false [K] | 
| Temperature | T_b_start | 293.15 | Initial temperature at port_b, used if steadyStateInitial = false [K] | 
record ParameterConstruction "Record for exterior constructions that have no window" extends Buildings.Rooms.BaseClasses.PartialParameterConstruction; parameter Modelica.SIunits.Area A "Heat transfer area";end ParameterConstruction; 
 Buildings.Rooms.BaseClasses.PartialParameterConstruction
Buildings.Rooms.BaseClasses.PartialParameterConstruction
This data record is used to set the parameters of constructions that do not have a window.
The surface azimuth is defined in Buildings.HeatTransfer.Types.Azimuth and the surface tilt is defined in Buildings.HeatTransfer.Types.Tilt
Extends from Modelica.Icons.Record (Icon for records).
| Type | Name | Default | Description | 
|---|---|---|---|
| Angle | til | Surface tilt [rad] | |
| Angle | azi | Surface azimuth [rad] | |
| Opaque construction | |||
| Generic | layers | redeclare parameter Building... | Material properties of opaque construction | 
| Initialization | |||
| Boolean | steadyStateInitial | false | =true initializes dT(0)/dt=0, false initializes T(0) at fixed temperature using T_a_start and T_b_start | 
| Temperature | T_a_start | 293.15 | Initial temperature at port_a, used if steadyStateInitial = false [K] | 
| Temperature | T_b_start | 293.15 | Initial temperature at port_b, used if steadyStateInitial = false [K] | 
record PartialParameterConstruction "Partial record for constructions" extends Modelica.Icons.Record;replaceable parameter Buildings.HeatTransfer.Data.OpaqueConstructions.Generic layers "Material properties of opaque construction"; parameter Modelica.SIunits.Angle til "Surface tilt"; parameter Modelica.SIunits.Angle azi "Surface azimuth"; final parameter Boolean isFloor=til > 2.74889125 and til < 3.53428875 "Flag, true if construction is a floor"; final parameter Boolean isCeiling=til > -0.392699 and til < 0.392699 "Flag, true if construction is a floor"; final parameter Integer nLay(min=1, fixed=true) = layers.nLay "Number of layers"; final parameter Integer nSta[nLay](min=1)={layers.material[i].nSta for i in 1:nLay} "Number of states"; parameter Boolean steadyStateInitial=false "=true initializes dT(0)/dt=0, false initializes T(0) at fixed temperature using T_a_start and T_b_start"; parameter Modelica.SIunits.Temperature T_a_start=293.15 "Initial temperature at port_a, used if steadyStateInitial = false"; parameter Modelica.SIunits.Temperature T_b_start=293.15 "Initial temperature at port_b, used if steadyStateInitial = false";end PartialParameterConstruction;