LBL logo

Buildings.Fluid.SolarCollectors.BaseClasses

Package with base classes for Buildings.Fluid.SolarCollectors

Information

This package contains base classes that are used to construct the models in Buildings.Fluid.SolarCollectors.

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

Package Content

NameDescription
Buildings.Fluid.SolarCollectors.BaseClasses.ASHRAEHeatLoss ASHRAEHeatLoss Calculate the heat loss of a solar collector per ASHRAE standard 93
Buildings.Fluid.SolarCollectors.BaseClasses.ASHRAESolarGain ASHRAESolarGain Calculate the solar heat gain of a solar collector per ASHRAE Standard 93
Buildings.Fluid.SolarCollectors.BaseClasses.EN12975HeatLoss EN12975HeatLoss Calculate the heat loss of a solar collector per EN12975
Buildings.Fluid.SolarCollectors.BaseClasses.EN12975SolarGain EN12975SolarGain Model calculating solar gains per the EN12975 standard
Buildings.Fluid.SolarCollectors.BaseClasses.PartialHeatLoss PartialHeatLoss Partial heat loss model on which ASHRAEHeatLoss and EN12975HeatLoss are based
PartialParameters Partial model for parameters
Buildings.Fluid.SolarCollectors.BaseClasses.PartialSolarCollector PartialSolarCollector Partial model for solar collectors
Buildings.Fluid.SolarCollectors.BaseClasses.IAM IAM Function for incident angle modifer
Buildings.Fluid.SolarCollectors.BaseClasses.Examples Examples Collection of models that illustrate model use and test models

Buildings.Fluid.SolarCollectors.BaseClasses.ASHRAEHeatLoss Buildings.Fluid.SolarCollectors.BaseClasses.ASHRAEHeatLoss

Calculate the heat loss of a solar collector per ASHRAE standard 93

Buildings.Fluid.SolarCollectors.BaseClasses.ASHRAEHeatLoss

Information

This component computes the heat loss from the solar thermal collector to the environment. It is designed for use with ratings data collected in accordance with ASHRAE Standard 93. A negative QLos[i] indicates that heat is being lost to the environment.

This model calculates the heat lost from a multiple-segment model using ratings data based solely on the inlet temperature. As a result, the slope from the ratings data is converted to a UA value which, for a given number of segments, yields the same heat loss as the ratings data would at nominal conditions. The first three equations, which perform calculations at nominal conditions without a UA value, are based on equations 6.17.1 through 6.17.3 in Duffie and Beckman (2006). The UA value is identified using the system of equations below:

QUse,nom = Gnom Ac FR(τα) + FRUL Ac (TIn,nom - TEnv,nom )
TFluid,nom[nSeg]=TIn,nom+QUse,nom/(m Flow,nomCp)
QLos,nom=-FRUL Ac (TIn,nom -TEnv,nom)
TFluid,nom[i] = TFluid,nom[i-1] + (Gnom F R(τα) Ac/nSeg - UA/nSeg (TFluid,nom[i-1] -TEnv,nom))/(mFlow,nom cp)
QLoss,UA=UA/nSeg (TFluid,nom[i]-TEnv,nom)
sum(QLoss,UA[1:nSeg])=QLoss,nom

where QUse,nom is the useful heat gain at nominal conditions, Gnom is the nominal solar irradiance, Ac is the area of the collector, FR(τα) is the collector maximum efficiency, FR UL is the collector heat loss coefficient,TIn,nom is the nominal inlet temperature, T Env,nom is the ambient temperature at nominal conditions, T Fluid,nom[i] is the temperature of fluid in a given segment of the collector, mFlow,nom is the fluid flow at nominal conditions, Cp is the specific heat of the heated fluid, Q Loss,nom is the heat loss identified using the default value UA is the identified heat loss coefficient for a multiple-segment equivalent solar collector, nSeg is the number of segments in the simulation, and Q Loss,UA is the heat loss identified using the UA value.

The effective UA value is calculated at the beginning of the simulation and used as a constant through the rest of the simulation. The actual heat loss from the collector is calculated using

-QLoss[i] = UA/nSeg (TFluid[i] - TEnv)

where QLoss[i] is the heat loss from a given segment, UA is the heat loss coefficient for a multiple segments model, nSeg is the number of segments in the simulation, TFluid[i] is the temperature of the fluid in a given segment, and TEnv is the temperature of the surrounding air.

This model reduces the heat loss rate to 0 W when the fluid temperature is within 1 degree C of the minimum temperature of the medium model. The calucation is performed using the Buildings.Utilities.Math.Functions.smoothHeaviside function.

References

J.A. Duffie and W.A. Beckman 2006, Solar Engineering of Thermal Processes (3rd Edition), John Wiley & Sons, Inc.
ASHRAE 93-2010 -- Methods of Testing to Determine the Thermal Performance of Solar Collectors (ANSI approved)

Extends from Buildings.Fluid.SolarCollectors.BaseClasses.PartialHeatLoss (Partial heat loss model on which ASHRAEHeatLoss and EN12975HeatLoss are based).

Parameters

TypeNameDefaultDescription
AreaA_c Area of the collector [m2]
IntegernSeg3Number of segments
Realy_intercept Y intercept (Maximum efficiency)
replaceable package MediumPartialMediumMedium in the component
CoefficientOfHeatTransferslope Slope from ratings data [W/(m2.K)]
Nominal condition
IrradianceG_nominal Irradiance at nominal conditions [W/m2]
TemperatureDifferencedT_nominal Ambient temperature at nomincal conditions [K]
MassFlowRatem_flow_nominal Fluid flow rate at nominal conditions [kg/s]

Connectors

TypeNameDescription
input RealInputTEnvTemperature of surrounding environment [K]
input RealInputTFlu[nSeg]Temperature of the heat transfer fluid [K]
output RealOutputQLos[nSeg]Limited heat loss rate at current conditions [W]

Modelica definition

block ASHRAEHeatLoss 
  "Calculate the heat loss of a solar collector per ASHRAE standard 93"

  extends Buildings.Fluid.SolarCollectors.BaseClasses.PartialHeatLoss;

  parameter Modelica.SIunits.CoefficientOfHeatTransfer slope 
    "Slope from ratings data";

protected 
  final parameter Modelica.SIunits.ThermalConductance UA(start = -slope*A_c, fixed = false) 
    "Coefficient describing heat loss to ambient conditions";
initial equation 
   //Identifies useful heat gain at nominal conditions
   QUse_nominal = G_nominal * A_c * y_intercept + slope * A_c * (dT_nominal);
   //Identifies TFlu[nSeg] at nominal conditions
   m_flow_nominal * Cp_default * (dT_nominal_fluid[nSeg]) = QUse_nominal;
   //Identifies heat lost to environment at nominal conditions
   QLos_nominal = -slope * A_c * (dT_nominal);
   //Governing equation for the first segment (i=1)
   G_nominal * y_intercept * A_c/nSeg - UA/nSeg * (dT_nominal) = m_flow_nominal * Cp_default
   * (dT_nominal_fluid[1]);
   //Loop with the governing equations for segments 2:nSeg-1
   for i in 2:nSeg-1 loop
     G_nominal * y_intercept * A_c/nSeg - UA/nSeg * (dT_nominal_fluid[i-1]+dT_nominal) =
     m_flow_nominal * Cp_default * (dT_nominal_fluid[i]-dT_nominal_fluid[i-1]);
   end for;
   for i in 1:nSeg loop
     nSeg * QLosUA[i] = UA * (dT_nominal_fluid[i]+dT_nominal);
   end for;
   sum(QLosUA[1:nSeg]) = QLos_nominal;
equation 
   for i in 1:nSeg loop
     -QLosInt[i] * nSeg = UA * (TFlu[i] - TEnv);
   end for;
end ASHRAEHeatLoss;

Buildings.Fluid.SolarCollectors.BaseClasses.ASHRAESolarGain Buildings.Fluid.SolarCollectors.BaseClasses.ASHRAESolarGain

Calculate the solar heat gain of a solar collector per ASHRAE Standard 93

Buildings.Fluid.SolarCollectors.BaseClasses.ASHRAESolarGain

Information

This component computes the solar heat gain of the solar thermal collector. It only calculates the solar heat gain without considering the heat loss to the environment. This model uses ratings data according to ASHRAE93. The solar heat gain is calculated using Equations 555 - 559 in the referenced EnergyPlus documentation.

The solar radiation absorbed by the panel is identified using Eq 559 from the EnergyPlus documentation. It is

QFlow[i]=Ac/nSeg (FR(τα) K(τα)net (GDir (1-shaCoe)+GDif,Sky+GDif,Gnd))

where QFlow[i] is the heat gain in each segment, A c is the area of the collector, nSeg is the user-specified number of segments in the simulation, FR(τα) is the maximum collector efficiency, K(τα)net> is the incidence angle modifier, GDir is the direct solar radiation, shaCoe is the user-specified shading coefficient, GDif,Sky is the diffuse solar radiation from the sky, and GDif,Gnd is the diffuse radiation from the ground.

The solar radiation equation indicates that the collector is divided into multiple segments. The number of segments used in the simulation is specified by the user (parameter: nSeg). The area of an individual segment is identified by dividing the collector area by the total number of segments. The term shaCoe is used to define the percentage of the collector that is shaded.

The incidence angle modifier used in the solar radiation equation is found using Eq 556 from the EnergyPlus documentation. It is

K(τα),net=(GBeam K(τα), Beam+GDif,Sky K(τα),Sky+GDif,Gnd K(τα),Gnd) / (Gbeam+GDif,Sky+G Dif,Gnd)

where K(τα),net is the net incidence angle modified, GBeam is the beam radiation, K(τα),Beam is the incidence angle modifier for beam radiation, GDif,Sky is the diffuse radiation from the sky, K(τα),Sky is the incidence angle modifier for radiation from the sky, GDif, Gnd is the diffuse radiation from the ground, and K(τα),Gnd is the incidence angle modifier for diffuse radiation from the ground.

Each incidence angle modifier is calculated using Eq 555 from the EnergyPlus documentation. It is

K(τα),x=1+b0 (1/cos(θ)-1)+b1 (1/cos(θ)-1)2

where x can refer to beam, sky or ground. θ is the incidence angle. For beam radiation θ is found via standard geometry. The incidence angle for sky and ground diffuse radiation are found using, respectively, Eq 557 and 558 from the EnergyPlus documentation. They are

θsky=59.68-0.1388 til+0.001497 til2
θgnd=90.0-0.5788 til+0.002693 til2

where θsky is the incidence angle for diffuse radiation from the sky, til is the tilt of the solar thermal collector, and θgnd is the incidence angle for diffuse radiation from the ground.

These two equations must be evaluated in degrees. The necessary unit conversions are made internally.

This model reduces the heat gain rate to 0 W when the fluid temperature is within 1 degree C of the maximum temperature of the medium model. The calucation is performed using the Buildings.Utilities.Math.Functions.smoothHeaviside function.

References

EnergyPlus 7.0.0 Engineering Reference, October 13, 2011.
ASHRAE 93-2010 -- Methods of Testing to Determine the Thermal Performance of Solar Collectors (ANSI approved)

Extends from Modelica.Blocks.Interfaces.BlockIcon (Basic graphical layout of input/output block), SolarCollectors.BaseClasses.PartialParameters (Partial model for parameters).

Parameters

TypeNameDefaultDescription
AreaA_c Area of the collector [m2]
IntegernSeg3Number of segments
Realy_intercept Y intercept (Maximum efficiency)
RealB0 1st incident angle modifer coefficient
RealB1 2nd incident angle modifer coefficient
Angletil Surface tilt [rad]
replaceable package MediumModelica.Media.Interfaces.Pa...Medium in the system
Shading
Booleanuse_shaCoe_infalseEnable input connector for shaCoe
RealshaCoe0Shading coefficient 0.0: no shading, 1.0: full shading

Connectors

TypeNameDescription
replaceable package MediumMedium in the system
input RealInputshaCoe_inShading coefficient
input RealInputTFlu[nSeg][K]
input RealInputHSkyDifTilDiffuse solar irradiation on a tilted surfce from the sky [W/m2]
input RealInputHGroDifTilDiffuse solar irradiation on a tilted surfce from the ground [W/m2]
input RealInputincAngIncidence angle of the sun beam on a tilted surface [rad]
input RealInputHDirTilDirect solar irradiation on a tilted surfce [W/m2]
output RealOutputQSol_flow[nSeg]Solar heat gain [W]

Modelica definition

block ASHRAESolarGain 
  "Calculate the solar heat gain of a solar collector per ASHRAE Standard 93"
  extends Modelica.Blocks.Interfaces.BlockIcon;
  extends SolarCollectors.BaseClasses.PartialParameters;

  parameter Real B0 "1st incident angle modifer coefficient";
  parameter Real B1 "2nd incident angle modifer coefficient";
  parameter Boolean use_shaCoe_in = false "Enable input connector for shaCoe";

  parameter Real shaCoe(
    min=0.0,
    max=1.0) = 0 "Shading coefficient 0.0: no shading, 1.0: full shading";

  parameter Modelica.SIunits.Angle til "Surface tilt";

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

  Modelica.Blocks.Interfaces.RealInput shaCoe_in if use_shaCoe_in 
    "Shading coefficient";
   Modelica.Blocks.Interfaces.RealInput TFlu[nSeg](
   unit = "K",
   displayUnit="degC",
   quantity="Temperature");
   Modelica.Blocks.Interfaces.RealInput HSkyDifTil(
     unit="W/m2", quantity="RadiantEnergyFluenceRate") 
    "Diffuse solar irradiation on a tilted surfce from the sky";
  Modelica.Blocks.Interfaces.RealInput HGroDifTil(
    unit="W/m2", quantity="RadiantEnergyFluenceRate") 
    "Diffuse solar irradiation on a tilted surfce from the ground";
  Modelica.Blocks.Interfaces.RealInput incAng(
    quantity="Angle",
    unit="rad",
    displayUnit="degree") "Incidence angle of the sun beam on a tilted surface";
  Modelica.Blocks.Interfaces.RealInput HDirTil(
    unit="W/m2", quantity="RadiantEnergyFluenceRate") 
    "Direct solar irradiation on a tilted surfce";
  Modelica.Blocks.Interfaces.RealOutput QSol_flow[nSeg](final unit="W") 
    "Solar heat gain";

protected 
  final parameter Real iamSky( fixed = false) 
    "Incident angle modifier for diffuse solar radiation from the sky";
  final parameter Real iamGro( fixed = false) 
    "Incident angle modifier for diffuse solar radiation from the ground";
  final parameter Modelica.SIunits.Angle incAngSky( fixed = false) 
    "Incident angle of diffuse radiation from the sky";
  final parameter Modelica.SIunits.Angle incAngGro( fixed = false) 
    "Incident angle of diffuse radiation from the ground";
  final parameter Real tilDeg(
    unit = "deg") = Modelica.SIunits.Conversions.to_deg(til) 
    "Surface tilt angle in degrees";
  final parameter Modelica.SIunits.HeatFlux HTotMin = 1 
    "Minimum HTot to avoid div/0";
  final parameter Real HMinDel = 0.001 
    "Delta of the smoothing function for HTot";

  Real iamBea "Incident angle modifier for director solar radiation";
  Real iam "Weighted incident angle modifier";

  Modelica.Blocks.Interfaces.RealInput shaCoe_internal 
    "Internally used shading coefficient";

initial equation 
  // E+ Equ (557)
  incAngSky = Modelica.SIunits.Conversions.from_deg(59.68 - 0.1388*(tilDeg) +
  0.001497*(tilDeg)^2);
  // Diffuse radiation from the sky
  // E+ Equ (555)
  iamSky = SolarCollectors.BaseClasses.IAM(incAngSky, B0, B1);
  // E+ Equ (558)
  incAngGro = Modelica.SIunits.Conversions.from_deg(90 - 0.5788*(tilDeg)+
  0.002693*(tilDeg)^2);
  // Diffuse radiation from the ground
  // E+ Equ (555)
  iamGro = SolarCollectors.BaseClasses.IAM(
    incAngGro,
    B0,
    B1);

equation 

  connect(shaCoe_internal, shaCoe_in);

  if not use_shaCoe_in then
    shaCoe_internal = shaCoe;
  end if;

  // E+ Equ (555)
  iamBea = Buildings.Utilities.Math.Functions.smoothHeaviside(1/3*
  Modelica.Constants.pi- incAng, Modelica.Constants.pi/60)*
  SolarCollectors.BaseClasses.IAM(
    incAng,
    B0,
    B1);
  // E+ Equ (556)
  iam = Buildings.Utilities.Math.Functions.smoothHeaviside(
      1/3*Modelica.Constants.pi-incAng,Modelica.Constants.pi/60)*
      ((HDirTil*iamBea + HSkyDifTil*iamSky + HGroDifTil*iamGro)/
      Buildings.Utilities.Math.Functions.smoothMax((HDirTil +
      HSkyDifTil + HGroDifTil), HTotMin, HMinDel));
  // Modified from EnergyPlus Equ (559) by applying shade effect for
  //direct solar radiation
  // Only solar heat gain is considered here
  for i in 1 : nSeg loop
    QSol_flow[i] = A_c/nSeg*(y_intercept*iam*(HDirTil*(1.0 -
    shaCoe_internal) + HSkyDifTil + HGroDifTil))*
    Buildings.Utilities.Math.Functions.smoothHeaviside(
     (Medium.T_max+1)-TFlu[i],1);
  end for;

end ASHRAESolarGain;

Buildings.Fluid.SolarCollectors.BaseClasses.EN12975HeatLoss Buildings.Fluid.SolarCollectors.BaseClasses.EN12975HeatLoss

Calculate the heat loss of a solar collector per EN12975

Buildings.Fluid.SolarCollectors.BaseClasses.EN12975HeatLoss

Information

This component computes the heat loss from the solar thermal collector to the environment. It is designed anticipating ratings data collected in accordance with EN12975. A negative QLos[i] indicates that heat is being lost to the environment.

This model calculates the heat lost from a multiple-segment model using ratings data based on the mean collector temperature. As a result, the slope from the ratings data must be converted to a UA value which, for a given number of segments, returns the same heat loss as the ratings data would at nominal conditions. The UA value is identified using the system of equations below:

QUse,nom = Gnom Ac FR (τα) - C1 Ac (TMean,nom- TEnv,nom)-C2 Ac (TMean,nom -TEnv,nom)2
TFluid,nom[nSeg]=TMean,nom+QUse,nom/ (mFlow,nom Cp)
QLos,nom=-C1 Ac (TMean,nom-T Env,nom)-C2 Ac (TMean,nom- TEnv,nom)2
TFluid,nom[i] = TFluid,nom[i-1] + (Gnom FR(τα) Ac/nSeg - UA/nSeg (TFluid, nom[i-1]-TEnv,nom))/(mFlow,nom cp )
QLoss,UA=UA/nSeg (TFluid,nom[i]-TEnv,nom )
sum(QLoss,UA[1:nSeg])=QLoss,nom

where QUse,nom is the useful heat gain at nominal conditions,Gnom is the nominal solar irradiance, Ac is the area of the collector, FR (τα) is the collector maximum efficiency, C1 is the collector heat loss coefficient, C2 is the temperature dependance of the heat loss coefficient, TMean, nom is the nominal mean temperature of the solar collector, T Env,nom is the ambient temperature at nominal conditions, TFluid,nom[i] is the temperature of fluid in a given segment of the collector, mFlow,nom is the fluid flow at nominal conditions, Cp is the specific heat of the heated fluid, QLoss,nom is the heat loss identified using the default value UA is the identified heat loss coefficient for a multiple-segment equivalent solar collector, nSeg is the number of segments in the simulation, and QLoss,UA is the heat loss identified using the UA value.

The effective UA value is calculated at the beginning of the simulation and used as a constant through the rest of the simulation. The actual heat loss from the collector is calculated using

QLoss[i] = UA/nSeg (TFluid[i] - TEnv).

where QLoss[i] is the heat loss from a given segment, UA is the heat loss coefficient for a multiple segments model, nSeg is the number of segments in the simulation,TFluid [i] is the temperature of the fluid in a given segment, and TEnv is the temperature of the surrounding air.

This model reduces the heat loss rate to 0 W when the fluid temperature is within 1 degree C of the minimum temperature of the medium model. The calucation is performed using the Buildings.Utilities.Math.Functions.smoothHeaviside function.

References

CEN 2006, European Standard 12975-1:2006, European Committee for Standardization

Extends from Buildings.Fluid.SolarCollectors.BaseClasses.PartialHeatLoss (Partial heat loss model on which ASHRAEHeatLoss and EN12975HeatLoss are based).

Parameters

TypeNameDefaultDescription
AreaA_c Area of the collector [m2]
IntegernSeg3Number of segments
Realy_intercept Y intercept (Maximum efficiency)
replaceable package MediumPartialMediumMedium in the component
CoefficientOfHeatTransferC1 C1 from ratings data [W/(m2.K)]
RealC2 C2 from ratings data [W/(m2.K2)]
Nominal condition
IrradianceG_nominal Irradiance at nominal conditions [W/m2]
TemperatureDifferencedT_nominal Ambient temperature at nomincal conditions [K]
MassFlowRatem_flow_nominal Fluid flow rate at nominal conditions [kg/s]

Connectors

TypeNameDescription
input RealInputTEnvTemperature of surrounding environment [K]
input RealInputTFlu[nSeg]Temperature of the heat transfer fluid [K]
output RealOutputQLos[nSeg]Limited heat loss rate at current conditions [W]

Modelica definition

block EN12975HeatLoss 
  "Calculate the heat loss of a solar collector per EN12975"

  extends Buildings.Fluid.SolarCollectors.BaseClasses.PartialHeatLoss;

  parameter Modelica.SIunits.CoefficientOfHeatTransfer C1 
    "C1 from ratings data";
  parameter Real C2(
  final unit = "W/(m2.K2)") "C2 from ratings data";

protected 
  final parameter Modelica.SIunits.ThermalConductance UA(
     fixed = false,
     start=QLos_nominal/(dT_nominal)) 
    "Coefficient describing heat loss to ambient conditions";
initial equation 
   //Identifies useful heat gain at nominal conditions
   QUse_nominal = G_nominal * A_c * y_intercept -C1 * A_c *
      (dT_nominal) - C2 * A_c * (dT_nominal)^2;
   //Identifies TFlu[nSeg] at nominal conditions
   m_flow_nominal * Cp_default * (dT_nominal_fluid[nSeg]) = QUse_nominal;
   //Identifies heat lost to environment at nominal conditions
   QLos_nominal = -C1 * A_c * (dT_nominal)-C2 * A_c * (dT_nominal)^2;
   //Governing equation for the first segment (i=1)
   G_nominal * y_intercept * A_c/nSeg - UA/nSeg * (dT_nominal)
     = m_flow_nominal * Cp_default * (dT_nominal_fluid[1]);
   //Loop with the governing equations for segments 2:nSeg-1
   for i in 2:nSeg-1 loop
     G_nominal * y_intercept * A_c/nSeg - UA/nSeg * (dT_nominal_fluid[i-1]
     +dT_nominal)= m_flow_nominal * Cp_default * (dT_nominal_fluid[i]-
     dT_nominal_fluid[i-1]);
   end for;
   for i in 1:nSeg loop
     nSeg * QLosUA[i] = UA * (dT_nominal_fluid[i]+dT_nominal);
   end for;
   sum(QLosUA) = QLos_nominal;
equation 
   for i in 1:nSeg loop
     QLosInt[i] * nSeg = UA * (TFlu[i] - TEnv);
   end for;
end EN12975HeatLoss;

Buildings.Fluid.SolarCollectors.BaseClasses.EN12975SolarGain Buildings.Fluid.SolarCollectors.BaseClasses.EN12975SolarGain

Model calculating solar gains per the EN12975 standard

Buildings.Fluid.SolarCollectors.BaseClasses.EN12975SolarGain

Information

This component computes the solar heat gain of the solar thermal collector. It only calculates the solar heat gain without considering the heat loss to the environment. This model performs calculations using ratings data from EN12975. The solar heat gain is calculated using Equation 559 in the referenced EnergyPlus documentation. The calculation is modified somewhat to use coefficients from EN12975.

The equation used to calculate solar gain is a modified version of Eq 559 from the EnergyPlus documentation. It is

QFlow[i] = Ac/nSeg FR(τα) (K (τα),Beam GBeam (1-shaCoe)+KDiff G Diff),

where QFlow[i] is the heat gained in each segment, A c is the area of the collector, nSeg is the number of segments in the collector, FR(τα) is the maximum efficiency of the collector, K(τα),Beam is the incidence angle modifier for beam radiation, GBeam is the current beam radiation on the collector, shaCoe is the shading coefficient, KDiff is the incidence angle modifier for diffuse radiation and GDiff is the diffuse radiation striking the surface.

The solar radiation equation indicates that the collector is divided into multiple segments. The number of segments used in the simulation is specified using the parameter nSeg. The area of an individual segment is identified by dividing the collector area by the total number of segments. The parameter shaCoe is used to define the percentage of the collector which is shaded. The main difference between this model and the ASHRAE model is the handling of diffuse radiation. The ASHRAE model contains calculated incidence angle modifiers for both sky and ground diffuse radiation while this model uses a coefficient from test data for diffuse radiation.

The incidence angle modifier for beam radiation is calculated using Eq 555 from the EnergyPlus documentation, as

K(τα),Beam=1+b0 (1/cos(θ)-1)+b1 (1/cos(θ)-1)2,

where K(τα),Beam is the incidence angle modifier for beam radiation, b0 is the first incidence angle modifier coefficient, θ is the incidence angle and b1 is the second incidence angle modifier coefficient.

This model reduces the heat gain rate to 0 W when the fluid temperature is within 1 degree C of the maximum temperature of the medium model. The calucation is performed using the Buildings.Utilities.Math.Functions.smoothHeaviside function.

References

EnergyPlus 7.0.0 Engineering Reference, October 13, 2011.
CEN 2006, European Standard 12975-1:2006, European Committee for Standardization

Extends from Modelica.Blocks.Interfaces.BlockIcon (Basic graphical layout of input/output block), SolarCollectors.BaseClasses.PartialParameters (Partial model for parameters).

Parameters

TypeNameDefaultDescription
AreaA_c Area of the collector [m2]
IntegernSeg3Number of segments
Realy_intercept Y intercept (Maximum efficiency)
RealB0 1st incident angle modifer coefficient
RealB1 2nd incident angle modifer coefficient
RealiamDiff Incidence angle modifier for diffuse radiation
replaceable package MediumModelica.Media.Interfaces.Pa...Medium in the system
Shading
Booleanuse_shaCoe_infalseEnables an input connector for shaCoe
RealshaCoe0Shading coefficient 0.0: no shading, 1.0: full shading

Connectors

TypeNameDescription
replaceable package MediumMedium in the system
input RealInputshaCoe_inTime varying input for the shading coefficient
input RealInputHSkyDifTilDiffuse solar irradiation on a tilted surfce from the sky [W/m2]
input RealInputincAngIncidence angle of the sun beam on a tilted surface [rad]
input RealInputHDirTilDirect solar irradiation on a tilted surfce [W/m2]
output RealOutputQSol_flow[nSeg]Solar heat gain [W]
input RealInputTFlu[nSeg][K]

Modelica definition

model EN12975SolarGain 
  "Model calculating solar gains per the EN12975 standard"
  extends Modelica.Blocks.Interfaces.BlockIcon;
  extends SolarCollectors.BaseClasses.PartialParameters;

  parameter Real B0 "1st incident angle modifer coefficient";
  parameter Real B1 "2nd incident angle modifer coefficient";
  parameter Boolean use_shaCoe_in = false 
    "Enables an input connector for shaCoe";
  parameter Real shaCoe(
    min=0.0,
    max=1.0) = 0 "Shading coefficient 0.0: no shading, 1.0: full shading";
  parameter Real iamDiff "Incidence angle modifier for diffuse radiation";

  replaceable package Medium = Modelica.Media.Interfaces.PartialMedium 
    "Medium in the system";
  Modelica.Blocks.Interfaces.RealInput shaCoe_in if use_shaCoe_in 
    "Time varying input for the shading coefficient";

  Modelica.Blocks.Interfaces.RealInput HSkyDifTil(
    unit="W/m2", quantity="RadiantEnergyFluenceRate") 
    "Diffuse solar irradiation on a tilted surfce from the sky";
  Modelica.Blocks.Interfaces.RealInput incAng(
    quantity="Angle",
    unit="rad",
    displayUnit="degree") "Incidence angle of the sun beam on a tilted surface";
  Modelica.Blocks.Interfaces.RealInput HDirTil(
    unit="W/m2", quantity="RadiantEnergyFluenceRate") 
    "Direct solar irradiation on a tilted surfce";
  Modelica.Blocks.Interfaces.RealOutput QSol_flow[nSeg](final unit="W") 
    "Solar heat gain";
  Modelica.Blocks.Interfaces.RealInput TFlu[nSeg](
     unit="K",
     displayUnit="degC",
     quantity="Temperature");

protected 
  Real iamBea "Incidence angle modifier for director solar radiation";
  Modelica.Blocks.Interfaces.RealInput shaCoe_internal "Internally used shaCoe";

equation 
  connect(shaCoe_internal, shaCoe_in);

  // E+ Equ (555)
  iamBea = Buildings.Utilities.Math.Functions.smoothHeaviside(x=Modelica.Constants.pi
    /3 - incAng, delta=Modelica.Constants.pi/60)*
    SolarCollectors.BaseClasses.IAM(
    incAng,
    B0,
    B1);
  // Modified from EnergyPlus Equ (559) by applying shade effect for direct solar radiation
  // Only solar heat gain is considered here

  if not use_shaCoe_in then
    shaCoe_internal = shaCoe;
  end if;

  for i in 1 : nSeg loop
  QSol_flow[i] = A_c/nSeg*(y_intercept*(iamBea*HDirTil*(1.0 - shaCoe_internal) + iamDiff *
  HSkyDifTil))*Buildings.Utilities.Math.Functions.smoothHeaviside((Medium.T_max+1)-TFlu[i],1);
  end for;
end EN12975SolarGain;

Buildings.Fluid.SolarCollectors.BaseClasses.PartialHeatLoss Buildings.Fluid.SolarCollectors.BaseClasses.PartialHeatLoss

Partial heat loss model on which ASHRAEHeatLoss and EN12975HeatLoss are based

Buildings.Fluid.SolarCollectors.BaseClasses.PartialHeatLoss

Information

This component is a partial model used as the base for Buildings.Fluid.SolarCollectors.BaseClasses.ASHRAEHeatLoss and Buildings.Fluid.SolarCollectors.BaseClasses.EN12975HeatLoss. It contains the input, output and parameter declarations which are common to both models. More detailed information is available in the documentation for the extending classes.

Extends from Modelica.Blocks.Interfaces.BlockIcon (Basic graphical layout of input/output block), SolarCollectors.BaseClasses.PartialParameters (Partial model for parameters).

Parameters

TypeNameDefaultDescription
AreaA_c Area of the collector [m2]
IntegernSeg3Number of segments
Realy_intercept Y intercept (Maximum efficiency)
replaceable package MediumModelica.Media.Interfaces.Pa...Medium in the component
Nominal condition
IrradianceG_nominal Irradiance at nominal conditions [W/m2]
TemperatureDifferencedT_nominal Ambient temperature at nomincal conditions [K]
MassFlowRatem_flow_nominal Fluid flow rate at nominal conditions [kg/s]

Connectors

TypeNameDescription
replaceable package MediumMedium in the component
input RealInputTEnvTemperature of surrounding environment [K]
input RealInputTFlu[nSeg]Temperature of the heat transfer fluid [K]
output RealOutputQLos[nSeg]Limited heat loss rate at current conditions [W]

Modelica definition

block PartialHeatLoss 
  "Partial heat loss model on which ASHRAEHeatLoss and EN12975HeatLoss are based"
  extends Modelica.Blocks.Interfaces.BlockIcon;
  extends SolarCollectors.BaseClasses.PartialParameters;

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

  parameter Modelica.SIunits.Irradiance G_nominal 
    "Irradiance at nominal conditions";
  parameter Modelica.SIunits.TemperatureDifference dT_nominal 
    "Ambient temperature at nomincal conditions";
  parameter Modelica.SIunits.MassFlowRate m_flow_nominal 
    "Fluid flow rate at nominal conditions";

  Modelica.Blocks.Interfaces.RealInput TEnv(
    quantity="ThermodynamicTemperature",
    unit="K",
    displayUnit="degC") "Temperature of surrounding environment";

public 
  Modelica.Blocks.Interfaces.RealInput TFlu[nSeg](
    quantity="ThermodynamicTemperature",
    unit = "K",
    displayUnit="degC") "Temperature of the heat transfer fluid";
  Modelica.Blocks.Interfaces.RealOutput QLos[nSeg](
    quantity="HeatFlowRate",
    unit="W",
    displayUnit="W") "Limited heat loss rate at current conditions";

protected 
  final parameter Modelica.SIunits.HeatFlowRate QUse_nominal(fixed = false) 
    "Useful heat gain at nominal conditions";
  final parameter Modelica.SIunits.HeatFlowRate QLos_nominal(fixed = false) 
    "Heat loss at nominal conditions";
  final parameter Modelica.SIunits.HeatFlowRate QLosUA[nSeg](fixed = false) 
    "Heat loss at current conditions";
  final parameter Modelica.SIunits.Temperature dT_nominal_fluid[nSeg](each 
  start = 293.15, fixed = false) 
    "Temperature of each semgent in the collector at nominal conditions";
  Medium.ThermodynamicState sta_default=Medium.setState_pTX(
      T=Medium.T_default,
      p=Medium.p_default,
      X=Medium.X_default);
  Modelica.SIunits.SpecificHeatCapacity Cp_default = Medium.specificHeatCapacityCp(sta_default) 
    "Specific heat capacity of the fluid";
  Modelica.SIunits.HeatFlowRate QLosInt[nSeg] 
    "Heat loss rate at current conditions";

equation 
  for i in 1:nSeg loop
    QLos[i] = QLosInt[i] * Buildings.Utilities.Math.Functions.smoothHeaviside(
     TFlu[i]-(Medium.T_min+1), 1);
  end for;

end PartialHeatLoss;

Buildings.Fluid.SolarCollectors.BaseClasses.PartialParameters

Partial model for parameters

Information

Partial parameters used in all solar collector models

Parameters

TypeNameDefaultDescription
AreaA_c Area of the collector [m2]
IntegernSeg3Number of segments
Realy_intercept Y intercept (Maximum efficiency)

Modelica definition

block PartialParameters "Partial model for parameters"

  parameter Modelica.SIunits.Area A_c "Area of the collector";
  parameter Integer nSeg(min=3)=3 "Number of segments";
  parameter Real y_intercept "Y intercept (Maximum efficiency)";

end PartialParameters;

Buildings.Fluid.SolarCollectors.BaseClasses.PartialSolarCollector Buildings.Fluid.SolarCollectors.BaseClasses.PartialSolarCollector

Partial model for solar collectors

Buildings.Fluid.SolarCollectors.BaseClasses.PartialSolarCollector

Information

This component is a partial model of a solar thermal collector. It can be expanded to create solar collector models based on either ASHRAE93 or EN12975 ratings data.

Notice

1. As mentioned in the reference, the SRCC incident angle modifier equation coefficients are only valid for incident angles of 60 degrees or less. Because these curves behave poorly for angles greater than 60 degrees the model does not calculatue either direct or diffuse solar radiation gains when the incidence angle is greater than 60 degrees.
2. By default, the estimated heat capacity of the collector without fluid is calculated based on the dry mass and the specific heat capacity of copper.

References

EnergyPlus 7.0.0 Engineering Reference, October 13, 2011.
CEN 2006, European Standard 12975-1:2006, European Committee for Standardization

Extends from Buildings.Fluid.Interfaces.LumpedVolumeDeclarations (Declarations for lumped volumes), Buildings.Fluid.Interfaces.TwoPortFlowResistanceParameters (Parameters for flow resistance for models with two ports), Buildings.Fluid.Interfaces.PartialTwoPortInterface (Partial model transporting fluid between two ports without storing mass or energy).

Parameters

TypeNameDefaultDescription
replaceable package MediumPartialMediumMedium in the component
IntegernSeg3Number of segments used to discretize the collector model
Anglelat Latitude [rad]
Angleazi Surface azimuth [rad]
Angletil Surface tilt [rad]
Realrho Ground reflectance
HeatCapacityC385*perPar.mDryHeat capacity of solar collector without fluid (default: cp_copper*mDry*nPanels) [J/K]
Nominal condition
Pressuredp_nominaldp_nominal_finalPressure [Pa]
MassFlowRatem_flow_nominalperPar.mperA_flow_nominal*pe...Nominal mass flow rate [kg/s]
Initialization
MassFlowRatem_flow.start0Mass flow rate from port_a to port_b (m_flow > 0 is design flow direction) [kg/s]
Pressuredp.start0Pressure difference between port_a and port_b [Pa]
Shading
Booleanuse_shaCoe_infalseEnables an input connector for shaCoe
RealshaCoe0Shading coefficient. 0.0: no shading, 1.0: full shading
Area declarations
NumberSelectionnColTypeBuildings.Fluid.SolarCollect...Selection of area specification format
IntegernPanels0Desired number of panels in the simulation
AreatotalArea0Total area of panels in the simulation [m2]
Configuration declarations
SystemConfigurationsysConfigBuildings.Fluid.SolarCollect...Selection of system configuration
Dynamics
Equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Formulation of energy balance
DynamicsmassDynamicsenergyDynamicsFormulation of mass balance
Initialization
AbsolutePressurep_startMedium.p_defaultStart value of pressure [Pa]
TemperatureT_startMedium.T_defaultStart value of temperature [K]
MassFractionX_start[Medium.nX]Medium.X_defaultStart value of mass fractions m_i/m [kg/kg]
ExtraPropertyC_start[Medium.nC]fill(0, Medium.nC)Start value of trace substances
ExtraPropertyC_nominal[Medium.nC]fill(1E-2, Medium.nC)Nominal value of trace substances. (Set to typical order of magnitude.)
Flow resistance
BooleancomputeFlowResistancetrue=true, compute flow resistance. Set to false to assume no friction
Booleanfrom_dpfalse= true, use m_flow = f(dp) else dp = f(m_flow)
BooleanlinearizeFlowResistancefalse= true, use linear relation between m_flow and dp for any flow rate
RealdeltaM0.1Fraction of nominal flow rate where flow transitions to laminar
Assumptions
BooleanallowFlowReversalsystem.allowFlowReversal= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem_flow_small1E-4*abs(m_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
BooleanhomotopyInitializationtrue= true, use homotopy method
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed

Connectors

TypeNameDescription
FluidPort_aport_aFluid connector a (positive design flow direction is from port_a to port_b)
FluidPort_bport_bFluid connector b (positive design flow direction is from port_a to port_b)
input RealInputshaCoe_inShading coefficient
BusweaBusWeather data bus

Modelica definition

model PartialSolarCollector "Partial model for solar collectors"
 extends Buildings.Fluid.Interfaces.LumpedVolumeDeclarations;
  extends Buildings.Fluid.Interfaces.TwoPortFlowResistanceParameters(final dp_nominal = dp_nominal_final);
  extends Buildings.Fluid.Interfaces.PartialTwoPortInterface(
    showDesignFlowDirection=false,
    final m_flow_nominal=perPar.mperA_flow_nominal*perPar.A);
  parameter Integer nSeg(min=3) = 3 
    "Number of segments used to discretize the collector model";

  parameter Modelica.SIunits.Angle lat "Latitude";
  parameter Modelica.SIunits.Angle azi "Surface azimuth";
  parameter Modelica.SIunits.Angle til "Surface tilt";
  parameter Real rho "Ground reflectance";
  parameter Modelica.SIunits.HeatCapacity C=385*perPar.mDry 
    "Heat capacity of solar collector without fluid (default: cp_copper*mDry*nPanels)";

  parameter Boolean use_shaCoe_in = false 
    "Enables an input connector for shaCoe";
  parameter Real shaCoe(
    min=0.0,
    max=1.0) = 0 "Shading coefficient. 0.0: no shading, 1.0: full shading";

  parameter Buildings.Fluid.SolarCollectors.Types.NumberSelection nColType=
  Buildings.Fluid.SolarCollectors.Types.NumberSelection.Number 
    "Selection of area specification format";
  parameter Integer nPanels= 0 "Desired number of panels in the simulation";

  parameter Modelica.SIunits.Area totalArea=0 
    "Total area of panels in the simulation";

  parameter Buildings.Fluid.SolarCollectors.Types.SystemConfiguration sysConfig=
  Buildings.Fluid.SolarCollectors.Types.SystemConfiguration.Series 
    "Selection of system configuration";

  Modelica.Blocks.Interfaces.RealInput shaCoe_in if use_shaCoe_in 
    "Shading coefficient";

  Modelica.Thermal.HeatTransfer.Components.HeatCapacitor heaCap[nSeg](
    T(each start =   T_start), each C=(C*nPanels_internal)/nSeg) if 
       not (energyDynamics == Modelica.Fluid.Types.Dynamics.SteadyState) 
    "Heat capacity for one segment of the the solar collector";

  Buildings.BoundaryConditions.WeatherData.Bus weaBus "Weather data bus";
  Buildings.BoundaryConditions.SolarIrradiation.DiffusePerez HDifTilIso(
    final outSkyCon=true,
    final outGroCon=true,
    final til=til,
    final lat=lat,
    final azi=azi,
    final rho=rho) "Diffuse solar irradiation on a tilted surface";

  Buildings.BoundaryConditions.SolarIrradiation.DirectTiltedSurface HDirTil(
    final til=til,
    final lat=lat,
    final azi=azi) "Direct solar irradiation on a tilted surface";

  Buildings.Fluid.Sensors.MassFlowRate senMasFlo(redeclare package Medium =
    Medium, allowFlowReversal=allowFlowReversal) "Mass flow rate sensor";
  Buildings.Fluid.FixedResistances.FixedResistanceDpM res(
    redeclare final package Medium = Medium,
    final from_dp=from_dp,
    final show_T=show_T,
    final m_flow_nominal=m_flow_nominal,
    final allowFlowReversal=allowFlowReversal,
    final linearized=linearizeFlowResistance,
    final homotopyInitialization=homotopyInitialization,
    use_dh=false,
    deltaM=deltaM,
    final dp_nominal=dp_nominal_final) "Flow resistance";
  Buildings.Fluid.MixingVolumes.MixingVolume vol[nSeg](
    each nPorts=2,
    redeclare package Medium = Medium,
    each final m_flow_nominal=m_flow_nominal,
    each final energyDynamics=energyDynamics,
    each final p_start=p_start,
    each final T_start=T_start,
    each m_flow_small=m_flow_small,
    each homotopyInitialization=homotopyInitialization,
    each final V=(perPar.V*nPanels_internal)/nSeg) 
    "Volume of fluid in one segment of the solar collector";

  Modelica.Thermal.HeatTransfer.Sensors.TemperatureSensor temSen[nSeg] 
    "Temperature sensor";

  Buildings.HeatTransfer.Sources.PrescribedHeatFlow heaGai[nSeg] 
    "Rate of solar heat gain";
  Buildings.HeatTransfer.Sources.PrescribedHeatFlow QLos[nSeg] 
    "Rate of solar heat gain";

protected 
  parameter Buildings.Fluid.SolarCollectors.Data.GenericSolarCollector perPar 
    "Partial performance data";

    Modelica.Blocks.Interfaces.RealInput shaCoe_internal 
    "Internally used shading coefficient";

    final parameter Modelica.SIunits.Pressure dp_nominal_final=
    if sysConfig == Buildings.Fluid.SolarCollectors.Types.SystemConfiguration.Series then 
       nPanels_internal*perPar.dp_nominal
    else 
      perPar.dp_nominal "Nominal pressure loss across the system of collectors";

  parameter Modelica.SIunits.Area TotalArea_internal=
      nPanels_internal * perPar.A "Area used in the simulation";

  parameter Real nPanels_internal=
    if nColType == Buildings.Fluid.SolarCollectors.Types.NumberSelection.Number then 
      nPanels
    else 
      totalArea/perPar.A "Number of panels used in the simulation";

equation 
  connect(shaCoe_internal,shaCoe_in);

  if not use_shaCoe_in then
    shaCoe_internal=shaCoe;
  end if;

  connect(weaBus, HDifTilIso.weaBus);
  connect(weaBus, HDirTil.weaBus);
  connect(port_a, senMasFlo.port_a);
  connect(senMasFlo.port_b, res.port_a);
  connect(heaCap.port, vol.heatPort);
      connect(vol[nSeg].ports[2], port_b);
  connect(vol[1].ports[1], res.port_b);
      for i in 1:(nSeg - 1) loop
    connect(vol[i].ports[2], vol[i + 1].ports[1]);
  end for;
  connect(vol.heatPort, temSen.port);
  connect(heaGai.port, vol.heatPort);
  connect(QLos.port, vol.heatPort);
end PartialSolarCollector;

Buildings.Fluid.SolarCollectors.BaseClasses.IAM

Function for incident angle modifer

Information

Overview

This function computes the incidence angle modifier for solar insolation striking the surface of the solar thermal collector. It is calculated using Eq 555 in the EnergyPlus 7.0.0 Engineering Reference.

Notice

As stated in EnergyPlus7.0.0 the incidence angle equation performs poorly at angles greater than 60 degrees. This model outputs 0 whenever the incidence angle is greater than 60 degrees.

References

EnergyPlus 7.0.0 Engineering Reference, October 13, 2011.

Inputs

TypeNameDefaultDescription
AngleincAng Incident angle [rad]
RealB0 1st incident angle modifer coefficient
RealB1 2nd incident angle modifer coefficient

Outputs

TypeNameDescription
RealincAngModIncident angle modifier coefficient

Modelica definition

function IAM "Function for incident angle modifer"

  input Modelica.SIunits.Angle incAng "Incident angle";
  input Real B0 "1st incident angle modifer coefficient";
  input Real B1 "2nd incident angle modifer coefficient";
  output Real incAngMod "Incident angle modifier coefficient";
  parameter Modelica.SIunits.Angle incAngMin = Modelica.Constants.pi / 2 -0.1 
    "Minimum incidence angle to avoid divison by zero";
  parameter Real cosIncAngMin(min=Modelica.Constants.eps)=
  Modelica.Math.cos(incAngMin) "Cosine of minimum incidence angle";
  parameter Real delta = 0.0001 "Width of the smoothing function";

algorithm 
  // E+ Equ (555)

  incAngMod :=
  Buildings.Utilities.Math.Functions.smoothHeaviside(
  Modelica.Constants.pi/3-incAng, delta)*
  (1 + B0*(1/Buildings.Utilities.Math.Functions.smoothMax(
  Modelica.Math.cos(incAng), Modelica.Math.cos(incAngMin), delta) - 1) +
  B1*(1/Buildings.Utilities.Math.Functions.smoothMax(Modelica.Math.cos(incAng),
  cosIncAngMin, delta) - 1)^2);

end IAM;

Automatically generated Thu Oct 24 15:09:54 2013.