This tutorial presents the theory of heat transfer in buildings and examplifies it on a toy model.
Objectives:
import numpy as np
import pandas as pd
import matplotlib.pyplot as plt
import dm4bem
Figure 1. Simple ventilated room (5 two-layer walls and 1 glass window) equiped with an HVAC system which acts as a proportional controller.
Let’s consider a cubic building with an HVAC systems acting as a proportional controller.
The dimensions and surface areas of the building are:
l = 3 # m length of the cubic room
Sg = l**2 # m² surface of the glass wall
Sc = Si = 5 * Sg # m² surface of concrete & insulation of the 5 walls
The thermophysical properties of the air (in SI units) are:
air = {'Density': 1.2, # kg/m³
'Specific heat': 1000} # J/(kg·K)
# pd.DataFrame.from_dict(air, orient='index', columns=['air'])
pd.DataFrame(air, index=['Air'])
Density | Specific heat | |
---|---|---|
Air | 1.2 | 1000 |
The thermophysical properties (thermal conductivities, densities and specific heat capacities) and the geometry (widths and surface areas) of the three materials (i.e., concrete, insulation, glass) in SI units are:
concrete = {'Conductivity': 1.400,
'Density': 2300.0,
'Specific heat': 880,
'Width': 0.2,
'Surface': 5 * l**2}
insulation = {'Conductivity': 0.027,
'Density': 55.0,
'Specific heat': 1210,
'Width': 0.08,
'Surface': 5 * l**2}
glass = {'Conductivity': 1.4,
'Density': 2500,
'Specific heat': 1210,
'Width': 0.04,
'Surface': l**2}
wall = pd.DataFrame.from_dict({'Layer_out': concrete,
'Layer_in': insulation,
'Glass': glass},
orient='index')
wall
Conductivity | Density | Specific heat | Width | Surface | |
---|---|---|---|---|---|
Layer_out | 1.400 | 2300.0 | 880 | 0.20 | 45 |
Layer_in | 0.027 | 55.0 | 1210 | 0.08 | 45 |
Glass | 1.400 | 2500.0 | 1210 | 0.04 | 9 |
where Slices
is the number of meshes used in the numerical discretization.
The radiative properties of the surfaces are:
# radiative properties
ε_wLW = 0.85 # long wave emmisivity: wall surface (concrete)
ε_gLW = 0.90 # long wave emmisivity: glass pyrex
α_wSW = 0.25 # short wave absortivity: white smooth surface
α_gSW = 0.38 # short wave absortivity: reflective blue glass
τ_gSW = 0.30 # short wave transmitance: reflective blue glass
The Stefan-Boltzmann constant is:
σ = 5.67e-8 # W/(m²⋅K⁴) Stefan-Bolzmann constant
print(f'σ = {σ} W/(m²⋅K⁴)')
σ = 5.67e-08 W/(m²⋅K⁴)
Conventional values for the convection coeficients for indoor and outdoor convection in W/(m²⋅K) are:
h = pd.DataFrame([{'in': 8., 'out': 25}], index=['h']) # W/(m²⋅K)
h
in | out | |
---|---|---|
h | 8.0 | 25 |
Thermal networks (or circuits) are weighted directed graphs in which:
Figure 2. Basic thermal network.
A thermal network has at least one oriented branch, $q$, and one node, $\theta$.
In a node there are a thermal capacity, $C_i$, (which can be positive or zero) and a heat flow rate source, $\dot Q_i$, (which can be zero).
On a branch, there are a conductane, $G_j > 0$, (which needs to be strictely pozitve) and a temperature source, $T_j$ (which can be zero).
The problem of analysis of thermal circuits (or the simulation problem, or the direct problem) is:
given:
find the temperature vector $\theta$ and the flow rate vector $q.$
For the toy model shown in Figure 1, heat transfert is:
The HVAC system is modelled as a proportional controller. There is long wave radiative exchange between the wall and the glass window. The sources are:
Figure 3. Heat processes for the cubic building shown in Figure 1.
Figure 4. Thermal circuit for the cubic building shown in Figure 1 and the heat processes shown in Figure 2. Note: space discretization of the walls is done for simplicity.
Figure 3 shows the models of:
The sources are:
Note: The known values, i.e. the elements of the circuit (the conductances $G$ and capacities $C$) and the sources (of temperature $T$ and of flow rate $\Phi$ or $\dot{Q}$) are noted in uppercase (majuscule) letters. The unknow variables, i.e. the temperatures in the nodes $\theta$ and the flow rates on the branches $q$, are noted in lowercase (minuscule) letters.
The conductances 1, 2, 3, and 4 of the thermal circuit from Figure 3 model the heat transfer by conduction. Conduction conductances, in W/K, are of the form: $$G_{cd} = \frac{\lambda}{w}S$$ where:
# conduction
G_cd = wall['Conductivity'] / wall['Width'] * wall['Surface']
pd.DataFrame(G_cd, columns=['Conductance'])
Conductance | |
---|---|
Layer_out | 315.0000 |
Layer_in | 15.1875 |
Glass | 315.0000 |
The conductances 0, 6 and 7 model the heat transfer by [convection](https://en.m.wikipedia.org/wiki/Convection_(heat_transfer). Convection conductances, in W/K, are of the form: $$G_{cv} = {h S}$$ where:
Table 1. Surface thermal resistances [Dal Zotto et al. 2014, p. 251]
Type of wall | Indoor surface | Outdoor surface |
---|---|---|
$h_i$,W/(m²·K) | $h_o$,W/(m²·K) | |
Vertical (tilt > 60°) | 7.7 | 25 |
Horizontal (tilt < 60°) | ||
- Upward heat flow | 10 | 25 |
- Downward heat flow | 5.9 | 25 |
# convection
Gw = h * wall['Surface'][0] # wall
Gg = h * wall['Surface'][2] # glass
The majority of methods used for modelling the radiative heat exchange use the view factors between surfaces. The view factor $F_{i,j}$ is defined as the proportion of radiation leaving surface $i$ that is intercepted by surface $j$. The view factors can be estimated by differential areas or for different configurations of surfaces [5].
The view factors need to satisfy the summation rule $$\sum_{j=0}^{n-1} F_{i,j} = 1$$ and the reciprocity theorem: $$F_{i,j} S_i = F_{j,i} S_j$$ where $S_{i}$ and $S_{j}$ are the surface areas.
For a convex surface $i$, the self-viewing factor is zero, $$F_{i,i} = 0$$
Two simplified relations are used to calculate the view factors for buildings.
In the first one, the view factors are defined by: $$\begin{cases} F_{i,j} = \frac{S_i}{S_T}\\ F_{i,i} = 0 \end{cases}$$
where $S_{T} = \sum_{j=0}^{n-1} S_j$, i.e. the surface $S_j$ is included in the total surface $S_T$. In this method, the reciprocity theorem is satisfied, $$F_{i,j} S_i = F_{j,i} S_j = \frac{S_i S_j}{S_T}$$ but summation rule isn't, $$\sum_{j=0}^{n-1} F_{i,j} = \sum_{j=0, j \neq i}^{n-1} \frac{S_j}{S_T} = \frac {S_T - S_i}{S_T} \neq 1$$
In this case, the heat balance for each surface would be wrong.
In the second one, the view factors are defined by: $$\begin{cases} F_{i,j} = \frac{S_j}{S_T -S_i}\\ F_{i,i} = 0 \end{cases}$$
where $S_{T} = \sum_{j=0}^{n-1} S_j$, i.e. the surface $S_i$ is not included in the total surface $S_{T,i} = S_T - S_i$.
In this case, the reciprocty theorem is generally not respected: $$F_{i, j} S_i = \frac{S_j}{S_T - S_i} S_i \neq F_{j, i} S_j = \frac{S_i}{S_T - S_j} S_j$$ but the summation rule is respected: $$ \sum_{i=0}^{n-1} F_{i, j} = \frac{1}{S_T - S_i} \sum_{\substack{j=0\\i\neq j}}^{n-1} S_j = 1 $$
# view factor wall-glass
Fwg = glass['Surface'] / concrete['Surface']
Note: The view factor between two surfaces, $j,k$ that are in the same plan (e.g. a window and a wall) is zero, $$F_{j,k} = F_{k,j}=0$$
Therefore the total surface $S_{T,i}$ should be: $$S_{T,i} = \sum_{j=0}^{n-1} S_j - \sum_k S_k$$ i.e. the surfaces $S_k$ in the same plan with the surface $S_i$ are not included in $S_{T,i}$.
The view factor between the top surface of finite wall $w$ tilted relative to an infinite plane of the ground $g$ is [5, 6, eq. 4.18]:
$$ F_{w,g} = \frac {1 - \cos \beta}{2}$$Therefore, the view factor between the tilted wall $w$ and the sky dome $s$ is [6, eq. 4.17]: $$ F_{w,s} = 1 - F_{w,g} = \frac {1 + \cos \beta}{2}$$
The long-wave heat exchange between surfaces may be modelled by using the concept of radiosity and then linearizing the radiative heat exchange.
Figure 5. Radiative long-wave heat exchange between two surfaces: a) modeled by emmitance (source) and radiosity (nodes); b) modeled by linearization of emmitance (temperature sources) and radiosity (temperature nodes).
For two surfaces, shown by temperature nodes 4 and 5 in Figure 3 and by nodes 1 and 2 in Figure 4, the conductances, in m², for radiative heat exchange expressed by using the emmitance (or the radiant excitance) of the black body, the radiosity, and the reciprocity of view factors are:
$$G_{1}^{r} = \frac{\varepsilon_1}{1 - \varepsilon_1} S_1$$$$G_{1,2}^{r} = F_{1,2} S_1 = F_{2,1} S_2$$$$G_{2}^{r} = \frac{\varepsilon_2}{1 - \varepsilon_2} S_2$$where:
The net flows leaving the surfaces 1 and 2 are:
$$q_{net,1} = \frac{\varepsilon_1}{1 - \varepsilon_1} S_1 (M^o_1 - J_1)= G^r_1 (M_1^o - J_1)$$$$q_{net,2} = \frac{\varepsilon_2}{1 - \varepsilon_2} S_2 (M^o_2 - J_2)= G^r_2 (M_2^o - J_2)$$respectively, where:
The net flow between surfaces 1 and 2 is:
$$q_{1,2} = F_{1,2} S_1 (J_1 - J_2) = F_{2,1} S_2 (J_1 - J_2)= G_{1,2}^r (J_1 - J_2)$$In order to express the long-wave radiative exchange as a function of temperature differences, a linearization of the difference of temperatures $T_1^4 - T_2^4$ may be used:
$$T_1^4 - T_2^4 = (T_1^2 + T_2^2)(T_1^2 - T_2^2) = (T_1^2 + T_2^2)(T_1 + T_2)(T_1 - T_2) = 4 \bar{T}^3 (T_1 - T_2)$$where the mean temperature $\bar{T}$, measured in kelvin, is:
$$\bar{T} =\sqrt[3]{ \frac{(T_1^2 + T_2^2)(T_1 + T_2)}{4}}$$The evaluation of mean temperaure, $\bar{T}$, requires the values of the surface tempetratures, $T_1$ and $T_2$. An initial guess can be used (and then an iterative process, for a more precise evaluation).
After linearization, the conductances, in W/K, for radiative heat exchange are:
$$G_{1} = 4 \sigma \bar{T}^3 \frac{\varepsilon_1}{1 - \varepsilon_1} S_1$$$$G_{1,2} = 4 \sigma \bar{T}^3 F_{1,2} S_1 = 4 \sigma \bar{T}^3 F_{2,1} S_2$$$$G_{2} = 4 \sigma \bar{T}^3 \frac{\varepsilon_2}{1 - \varepsilon_2} S_2$$# long wave radiation
Tm = 20 + 273 # K, mean temp for radiative exchange
GLW1 = 4 * σ * Tm**3 * ε_wLW / (1 - ε_wLW) * wall['Surface']['Layer_in']
GLW12 = 4 * σ * Tm**3 * Fwg * wall['Surface']['Layer_in']
GLW2 = 4 * σ * Tm**3 * ε_gLW / (1 - ε_gLW) * wall['Surface']['Glass']
The equivalent conductance, in W/K, for the radiative long-wave heat exchange between the wall and the glass window is: $$G = \frac{1}{1/G_1 + 1/G_{1,2} + 1/G_2}$$
GLW = 1 / (1 / GLW1 + 1 / GLW12 + 1 / GLW2)
Note: Resistances in series or parallel can be replaced by their equivalent resistance.
The volumetric flow rate of the air, in m³/s, is:
$$\dot{V}_a = \frac{\mathrm{ACH}}{3600} V_a$$where:
# ventilation flow rate
Va = l**3 # m³, volume of air
ACH = 1 # air changes per hour
Va_dot = ACH / 3600 * Va # m³/s, air infiltration
The net flow rate that the building receives by advection, i.e., introducing outdoor air at temperature $T_o$ and extracting indoor air at temperature $\theta_i$ by ventilation and/or air infiltration, is:
$$q_v = \dot{m}_a c_a (T_o - \theta_i) = \rho_a c_a \dot{V}_a (T_o - \theta_i)$$where:
Therefore, the conductance of advection by ventilation and/or infiltration, in W/K, is:
$$G_v = \rho_a c_a \dot{V}_a$$# ventilation & advection
Gv = air['Density'] * air['Specific heat'] * Va_dot
Table 2. Typical values for the ventilation rates (in air changes per hour, ACH) as a function of the position of windows (H. Recknagel, E. Spenger, E_R Schramek (2013), Table 1.12.1-4)
Position of windows | Ventilation rate, ACH [h⁻ⁱ] |
---|---|
Window closed, doors closed | 0 to 0.5 |
Tilted window, venetian blind closed | 0.3 to 1.5 |
Tilted window, whitout venetian blind | 0.8 to 4.0 |
Window half opened | 5 to 10 |
Window fully open | 9 to 15 |
Window and French window fully open (cross ventilation) | about 40 |
In the simplest representation, the HVAC system can be considered as a proportional controller that adjusts the heat flow rate $q_{HVAC}$ in order to control the indoor temperature $\theta_i$ at its setpoint value $T_{i,sp}$. The heat flow-rate, in W, injected by the HVAC system into the controlled space is:
$$ q_{HVAC} = K_p (T_{i, sp} - \theta_i)$$where:
This equation shows that the proportional controller can be modelled by a source of temperature, $T_{i, sp}$, and a conductance, $K_p$. If the controller gain tends towards:
Note: Respecting the sign convention, the flow rate $q_{HVAC}$ is oriented from the lower to the higher potential of the temperature source $T_{i,sp}$.
# P-controler gain
Kp = 1e4 # almost perfect controller Kp -> ∞
Kp = 1e-3 # no controller Kp -> 0
Kp = 0
If conductances are connected to temperature nodes which have no capacity and/or flow rate source, then the conductances can be considered in series or parallel (depending on the connection). Let's consider, for example, the outdoor side of the glass window (Figure 3): the outdoor convection conductance and the conduction conductance (corresponding to half of the width of the glass) are in series:
$$ G_{gs} = \frac{1}{1/G_{g,cv.out } + 1/(2 G_{g,cd})} = \frac{1}{\frac{1}{h_{out} S_g} + \frac{w / 2}{\lambda S_g}} $$# glass: convection outdoor & conduction
Ggs = float(1 / (1 / Gg.loc['h', 'out'] + 1 / (2 * G_cd['Glass'])))
The thermal capacities of the wall, in J/kg, are:
$$C_w= m_w c_w= \rho_w c_w w_w S_w$$where:
C = wall['Density'] * wall['Specific heat'] * wall['Surface'] * wall['Width']
C['Air'] = air['Density'] * air['Specific heat'] * Va
pd.DataFrame(C, columns=['Capacity'])
Capacity | |
---|---|
Layer_out | 18216000.0 |
Layer_in | 239580.0 |
Glass | 1089000.0 |
Air | 32400.0 |
The temperature sources model temperatures which vary independently of what happens in the themal circuit; they are inputs of the physical model. Generally, the temperature sources are:
The hourly values of outdoor temperatures can be obtained from weather data files downloadable from the Repository of free climate data for building performance simulation or from Weather data for EnergyPlus® (see the tutorial on Weather data and solar radiation).
If the adjacent spaces are controlled by a HVAC system, it means that their temperature can be considered independent of the studied thermal zone(s); therefore, they can be modelled by a temperature source.
Setpoint temperature does not depend on the heat transfer processes of the analyzed thermal zone. If the HVAC system can deliver the heat flow rate:
$$ q_{HVAC} = K_p (T_{i, sp} - \theta_i)$$where:
then the setpoint for indoor temperature, $T_{i, sp}$, may be modelled by a source of temperature.
The heat flow rate sources model flow rates which vary idependently of what happens in the themal circuit. They are inputs of the physical model. Generally, the heat flow rate sources are:
The direct, diffuse and reflected components of the solar radiation on a tilted surface can be estimated from weather data by using the function sol_rad_tilt_surf
from the module dm4bem
(see the tutorial on Weather data and solar radiation).
The radiation absorbed by the outdoor surface of the wall is:
$$\Phi_o = \alpha_{w,SW} S_w E_{tot}$$where:
The total shortwave incident irradiance on the wall $i$, $E_i$, may be estimated as a function of the direct solar irradiance incident on the surface of the walls, $E_{i}^{o}$:
$$S_i E_i = S_i E_{i}^{o} + \sum_{j=1}^{n} F_{j,i} S_j \rho_j E_j$$where:
By taking into account the reciprocity of the view factors: $S_i F_{i,j} = S_j F_{j,i}$, the set of previous equation becomes:
$$ \begin{bmatrix} 1 - \rho_1 F_{1,1} & - \rho_2 F_{1,2} & ... & - \rho_n F_{1,n}\\ - \rho_1 F_{2,1} & 1 - \rho_2 F_{2,2} & ... & - \rho_n F_{2,n} \\ ... & ... & ... & ... \\ - \rho_1 F_{n,1} & - \rho_2 F_{n,1} & ... & 1 - \rho_n F_{n,n} \end{bmatrix} \begin{bmatrix} E_1\\ E_2\\ ...\\ E_n \end{bmatrix} = \begin{bmatrix} E_{1}^{o}\\ E_{2}^{o}\\ ...\\ E_{n}^{o} \end{bmatrix} $$or
$$(I - \rho \circ F) E = E^o$$The unknown total irradiances on walls, in W/m², are then
$$ E = (I - \rho \circ F)^{-1} E^o$$where:
$I =\begin{bmatrix} 1 & 0 & ... & 0 \\ 0 & 1 & ... & 0 \\ ... & ... & ... & ...\\ 0 & 0 & ... & 1 \end{bmatrix}, $ is the identity matrix;
$\rho = \begin{bmatrix} \rho_1\\ \rho_2\\ ...\\ \rho_n \end{bmatrix}$ - vector of reflectances, $0 \le \rho_{i,j} \le 1$;
$F = \begin{bmatrix} F_{1,1} & F_{1,2} & ... & F_{1,n}\\ F_{2,1} & F_{2,2} & ... & F_{2,n} \\ ... & ... & ... & ...\\ F_{n,1} & F_{n,2} & ... & F_{n,n} \end{bmatrix}$ - matrix of view factors, $0 \le F_{i,j} \le 1$;
$E^o = \begin{bmatrix} E_{1}^{o}\\ E_{2}^{o}\\ ...\\ E_{n}^{o} \end{bmatrix}$ - vector of direct solar irradiances, W/m²;
$E = \begin{bmatrix} E_1\\ E_2\\ ...\\ E_n \end{bmatrix}$ - vector of unknown total irradiances, W/m².
The radiative short wave (i.e. solar) heat flow rate on each surface is:
$$ \Phi = S E $$where:
$\Phi = \begin{bmatrix} \Phi_1\\ \Phi_2\\ ...\\ \Phi_n \end{bmatrix}$ - vector of total heat flow rates due to solar radiation, W;
$S =\begin{bmatrix} S_1 & 0 & ... & 0 \\ 0 & S_2 & ... & 0 \\ ... & ... & ... & ...\\ 0 & 0 & ... & S_n \end{bmatrix}$ - matrix of surface areas of walls $i$, m².
Internal flow rates are generated by occupants and by the electrical equipment (with values given for offices, commercial spaces, etc.).
The analysis of a thermal circuit, or the direct problem (Ghiaus 2022), means to find the temperatures in the nodes, $\theta$, and the heat flows on the branches, $q$, i.e. to solve for $\theta$ and $q$ the system of Differential-Algebraic Equations (DAE) (Figures 3 and 4):
$$\left\{\begin{array}{ll} C \dot{\theta} = -(A^T G A) \theta + A^T G b + f\\ q = G (-A \theta + b) \end{array}\right.$$where:
$\theta$ is the temperature vector of size $n_\theta$ equal to the number of nodes;
$q$ - heat flow vector of size $n_q$ equal to the number of branches;
$A$ - incidence matrix of size $n_q$ rows and $n_{\theta}$ columns, where $n_q$ is the number of flow branches and $n_{\theta}$ is the number of temperature nodes. It shows how the temperature nodes are connected by oriented branches of heat flows:
$G$ - conductance diagonal matrix of size $n_q \times n_q$, where $n_q$ is the number of flow branches: diagonal matrix containing the conductances. Each branch $k$ needs to contain a conductance $0 < G_{k,k} < \infty $.
$C$ - capacity diagonal matrix of size $n_θ \times n_θ$, where $n_θ$ is the number of temperature nodes: diagonal matrix containing the capacities. If there is no capacity in the node n, then $C_{n, n} = 0$.
$b$ - temperature source vector of size $n_q$: if there is no temperature source on the branch m, then $b_m = 0$.
$f$ - heat flow source vector of size $n_θ$: if there is no heat flow source in the node n, then $f_n = 0$.
The resolution is first done for temperatures, $\theta$, by solving the equation $$C \dot{\theta} = -(A^T G A) \theta + A^T G b + f$$ which, generally, is a system of differential-algebraic equations (DAE). Then, the heat flow rates are found from the equation $$q = G (-A \theta + b)$$
\begin{aligned} &C \dot{\theta} = -(A^T G A) \theta + A^T G b + f \\ &q = G (-A \theta + b) \end{aligned}Figure 6. Matrices of the system of Differential-Algebraic Equations (DAE)
The incidence matrix is:
$A_{kl} = \begin{cases}\phantom{-} 0 & \text{if branch } q_k \text{ is not connected to node } \theta_l \\ +1 & \text{if branch } q_k \text{ enters into node } \theta_l\\ -1 & \text{if branch } q_k \text{ gets out of node } \theta_l \end{cases}$
For the themal circuit shown in Figure 3,
$ A = \begin{cases} A_{0,0} = 1\\ A_{1,0} = -1, A_{1,1} = 1\\ ...\\ A_{11,6} = 1\\ \end{cases}$
A = np.zeros([12, 8]) # n° of branches X n° of nodes
A[0, 0] = 1 # branch 0: -> node 0
A[1, 0], A[1, 1] = -1, 1 # branch 1: node 0 -> node 1
A[2, 1], A[2, 2] = -1, 1 # branch 2: node 1 -> node 2
A[3, 2], A[3, 3] = -1, 1 # branch 3: node 2 -> node 3
A[4, 3], A[4, 4] = -1, 1 # branch 4: node 3 -> node 4
A[5, 4], A[5, 5] = -1, 1 # branch 5: node 4 -> node 5
A[6, 4], A[6, 6] = -1, 1 # branch 6: node 4 -> node 6
A[7, 5], A[7, 6] = -1, 1 # branch 7: node 5 -> node 6
A[8, 7] = 1 # branch 8: -> node 7
A[9, 5], A[9, 7] = 1, -1 # branch 9: node 5 -> node 7
A[10, 6] = 1 # branch 10: -> node 6
A[11, 6] = 1 # branch 11: -> node 6
# np.set_printoptions(suppress=False)
# pd.DataFrame(A)
The conductance matrix of the themal circuit shown in Figure 3 is diagonal:
$ G = \begin{cases} G_{0,0} = G_{w,out} & \text{convection outside surface of the wall}\\ G_{1,1} = G_{2,2} = 2G_{cd,Layer\,out} & \text{conduction in half width of the outer layer}\\ G_{3,3} = G_{4,4} = 2G_{cd,Layer\,in} & \text{conduction in half width of the inner layer}\\ G_{5,5} = G_{LW} & \text{long-wave radiation walls - window}\\ G_{6,6} = G_{w,in} & \text{convection inside surface of the wall}\\ G_{7,7} = G_{g,in} & \text{convection inside surface of the glass}\\ G_{8,8} = G_{gs} & \text{convection outside surface of the glass}\\ & \text{and conduction in half width of the glass}\\ G_{9,9} = 2G_{cd,glass} & \text{conduction in half width of the glass}\\ G_{10,10} = G_v & \text{advection by ventilation}\\ G_{11,11} = K_p & \text{gain of proportional controller} \end{cases}$
G = np.diag(np.hstack(
[Gw['out'],
2 * G_cd['Layer_out'], 2 * G_cd['Layer_out'],
2 * G_cd['Layer_in'], 2 * G_cd['Layer_in'],
GLW,
Gw['in'],
Gg['in'],
Ggs,
2 * G_cd['Glass'],
Gv,
Kp]))
# np.set_printoptions(precision=3, threshold=16, suppress=True)
# pd.set_option("display.precision", 1)
# pd.DataFrame(G)
The capacity matrix of the themal circuit shown in Figure 3 is diagonal:
$ C = \begin{cases} C_{1,1} = C_{Layer\,out} & \text{outer layer of the wall}\\ C_{3,3} = C_{Layer\,in} & \text{inner layer of the wall}\\ C_{6,6} = C_{Air} & \text{air of the room}\\ C_{7,7} = C_{Glass} & \text{glass of the windows}\\ \end{cases}$
The thermal capacities of the air and of the glass can be neglected or not.
neglect_air_glass = False
if neglect_air_glass:
C = np.diag([0, C['Layer_out'], 0, C['Layer_in'], 0, 0,
0, 0])
else:
C = np.diag([0, C['Layer_out'], 0, C['Layer_in'], 0, 0,
C['Air'], C['Glass']])
# pd.set_option("display.precision", 3)
# pd.DataFrame(C)
The vector of temperature sources is $b$, of size $n_q$, the number of branches (in this example 12). An element of the vector $b$ corresponding to a branch without a source is zero. If the flow in a source is from the low potential to the high potential of the source (i.e. from - to +), then the source is positive. If the flow rate in the temperature source is from high potential to low potential (i.e. from + to -), then the source is negative (see passive sign convention).
For the thermal circuit shown in Figure 3,
$$b = [\begin{matrix} T_o &0 &0 &0 &0 &0 &0 &0 &T_o &0 &T_o &T_{i,sp} \end{matrix}]^T$$i.e. $b_0 = b_8 = b_{10} = T_o$ and $b_{11} = T_{i,sp}$ where:
Since the temperature sorces $T_o$ and $T_{i,sp}$ are time series, in vector $b$ the branches which contain temperature sources are designated by $1$ and the branches without any temeprature source by $0$.
b = np.zeros(12) # branches
b[[0, 8, 10, 11]] = 1 # branches with temperature sources
print(f'b = ', b)
b = [1. 0. 0. 0. 0. 0. 0. 0. 1. 0. 1. 1.]
The vector of heat sources is $f$, of size $n_{\theta}$, the number of nodes (in this example 8). An element of the vector $f$ corresponding to a node without a heat source is zero.
For the thermal circuit shown in Figure 3,
$$f = [\begin{matrix} \Phi_o &0 &0 &0 &\Phi_i &0 &\dot{Q_a} &\Phi_a \end{matrix}]^T$$i.e. $f_0 = \Phi_o$, $f_4 = \Phi_i$, $f_6 = \dot{Q_a}$, and $f_7 = \Phi_a$, where:
Since the flow rate sorces $\Phi_o$, $\Phi_i$, $\dot{Q}_a$ and $\Phi_a$ are time series, in vector $f$ the nodes which contain flow rate sources are designated by $1$ and the nodes without any flow rate source by $0$.
f = np.zeros(8) # nodes
f[[0, 4, 6, 7]] = 1 # nodes with heat-flow sources
print(f'f = ', f)
f = [1. 0. 0. 0. 1. 0. 1. 1.]
The vector of outputs is $y$, of size $n_{\theta}$, the number of nodes (in this example 8). The non-zero values of $y$ indicate the nodes which are the outputs of the model.
For the thermal circuit shown in Figure 3, if the output is the indoor air temperature, then the output vector is:
$$y = [\begin{matrix} 0 &0 &0 &0 &0 &0 &\theta_6 &0 \end{matrix}]^T$$In vector $y$, the nodes for which the temperatures are outputs are noted by $1$ and the other nodes by $0$.
y = np.zeros(8) # nodes
y[[6]] = 1 # nodes (temperatures) of interest
print(f'y = ', y)
y = [0. 0. 0. 0. 0. 0. 1. 0.]
The differential-algebraic system of equations (DAE)
$$C \dot{\theta} = -(A^T G A) \theta + A^T G b + f$$is transformed in state-space representation (Ghiaus 2013):
$$\left\{\begin{array}{rr} \dot{\theta}_s=A_s \theta_s + B_s u\\ y = C_s \theta_s + D_s u \end{array}\right.$$where:
$\theta_s$ is the vector of state variables which are the temperatures of nodes containing capacities; the elements are in the same order as in the vector of temperatures, $\theta$; its dimension, $\dim \theta_s$, is equal to the number of capacities from the thermal network; for the circuit presented in Figure 3, $\theta_s = [\theta_1, \theta_3, \theta_6, \theta_7]^T$;
$u = \begin{bmatrix} b_T \\ f_Q\end{bmatrix}$ - vector of inputs of dimension $\dim u$ equal to the number of sources (of temperaure, $b_T$, and heat flows, $f_Q$) of the thermal network, where:
$y$ - vector of outputs, a subset of vector $\theta$ representing temperature nodes which are of interest; for the circuit presented in Figure 3, $y = \theta_6$, the indoor temperature.
$A_s$ - state matrix, of dimension $\dim A_s = \dim {\theta_s} \times \dim {\theta_s}$;
$B_s$ - input matrix, of dimension $\dim B_s = \dim {\theta_s} \times \dim u$;
$C_s$ - output matrix, of dimension $\dim C_s = \dim y \times \dim {\theta_s}$;
$D_s$ - feedthrough (or feedforward) matrix, of dimension $\dim D_s = \dim y \times \dim u$.
Note: The subscript $s$ of the matrices $A_s, B_s, C_s, D_s$ is used to differentiante the matrices $A_s, C_s$ of the state-space represenation of the matrices $A, C$ of the system of DAE.
The state-space representation, i.e., matrices $A_s, B_s, C_s, D_s$ is obtained from the system of DAE, i.e., matrices and vectors $A, G, b, C, f, y$ (Ghiaus 2013).
[As, Bs, Cs, Ds] = dm4bem.tc2ss(A, G, C, b, f, y)
θs = ['θ1', 'θ3', 'θ6', 'θ7'] # state temperature nodes
uT = ['q0', 'q8', 'q10', 'q11'] # temperature sources
uQ = ['θ0', 'θ4', 'θ6', 'θ7'] # flow sources
u = uT + uQ # inputs
y = ['θ6'] # output
pd.DataFrame(As, index=θs, columns=θs)
θ1 | θ3 | θ6 | θ7 | |
---|---|---|---|---|
θ1 | -0.000024 | 0.000002 | 0.000000 | 0.000000 |
θ3 | 0.000121 | -0.000239 | 0.000107 | 0.000011 |
θ6 | 0.000000 | 0.000790 | -0.003925 | 0.002857 |
θ7 | 0.000000 | 0.000002 | 0.000085 | -0.000240 |
pd.DataFrame(Bs, index=θs, columns=u)
q0 | q8 | q10 | q11 | θ0 | θ4 | θ6 | θ7 | |
---|---|---|---|---|---|---|---|---|
θ1 | 0.000022 | 0.000000 | 0.000000 | 0.0 | 1.970654e-08 | 0.000000e+00 | 0.000000 | 0.000000e+00 |
θ3 | 0.000000 | 0.000000 | 0.000000 | 0.0 | 0.000000e+00 | 2.931594e-07 | 0.000000 | 0.000000e+00 |
θ6 | 0.000000 | 0.000000 | 0.000278 | 0.0 | 0.000000e+00 | 2.600003e-05 | 0.000031 | 0.000000e+00 |
θ7 | 0.000000 | 0.000152 | 0.000000 | 0.0 | 0.000000e+00 | 8.022402e-08 | 0.000000 | 9.182736e-07 |
pd.DataFrame(Cs, index=y, columns=θs)
θ1 | θ3 | θ6 | θ7 | |
---|---|---|---|---|
θ6 | 0.0 | 0.0 | 1.0 | 0.0 |
pd.DataFrame(Ds, index=y, columns=u)
q0 | q8 | q10 | q11 | θ0 | θ4 | θ6 | θ7 | |
---|---|---|---|---|---|---|---|---|
θ6 | 0.0 | 0.0 | 0.0 | 0.0 | 0.0 | 0.0 | 0.0 | 0.0 |
Steady-state means that the term $C \dot \theta = 0$ in the system of DAE.
In steady-state, the model can be checked if it is incorrect. Let's consider that:
b = np.zeros(12) # temperature sources
b[[0, 8, 10]] = 10 # outdoor temperature
b[[11]] = 20 # indoor set-point temperature
f = np.zeros(8) # flow-rate sources
Note: Steady-state analysis is a test of falsification (refutability) of the model, not a verification and validation. If the model does not pass the steady-state test, it means that it is wrong. If the model passes the steady-state test, it does not mean that it is correct. For example, the values of the capacities in matrix $C$ or of the conductances in matrix $G$ can be wrong even when the steady-state test is passed.
The value of temperature in steady-state is obtained from the system of DAE by considering that $C \dot{\theta} = 0$:
$$\theta_{ss} = (A^T G A)^{-1}(A^T G b + f)$$For the conditions mentioned above, in steady-state, all temperatures $\theta_0 ... \theta_7$, including the indoor air temperature $\theta_6$, are equal to $T_o = 10 \, \mathrm{^\circ C}$.
θ = np.linalg.inv(A.T @ G @ A) @ (A.T @ G @ b + f)
print(f'θ = {θ} °C')
θ = [10. 10. 10. 10. 10. 10. 10. 10.] °C
The input vector $u$ is obtained by stacking the vectors $b_T$ and $f_Q$:
$$u = \begin{bmatrix} b_T \\ f_Q\end{bmatrix}$$where:
Note: Zero in vectors $b$ and $f$ indicates that there is no source on the branch or in the node, respectively. However, a source can have the value zero.
bT = np.array([10, 10, 10, 20]) # [To, To, To, Tisp]
fQ = np.array([0, 0, 0, 0]) # [Φo, Φi, Qa, Φa]
u = np.hstack([bT, fQ])
print(f'u = {u}')
u = [10 10 10 20 0 0 0 0]
The steady-state value of the output of the state-space representation is obtained when $\dot \theta_{C} = 0$:
$$y_{ss} = (-C_s A_s^{-1} B_s + D_s) u$$yss = (-Cs @ np.linalg.inv(As) @ Bs + Ds) @ u
print(f'yss = {yss} °C')
yss = [10.] °C
The error between the steady-state values obtained from the system of DAE, $\theta_6$, and the output of the state-space representation, $y_{ss}$,
$$\varepsilon = \left | \theta_6 - y_{ss} \right |$$is practically zero; the slight difference is due to numerical errors.
print(f'Max error between DAE and state-space: \
{max(abs(θ[6] - yss)):.2e} °C')
Max error between DAE and state-space: 8.88e-15 °C
The condition for numerical stability of Euler explicit integration method is
$$\left | \lambda_i \Delta t + 1 \right | < 1, \forall \lambda_i, $$i.e. in the complex plane, $\lambda_i \Delta t$ is inside a circle of radius 1 centered in {-1, 0j}, where:
For positive real eigenvalues $\left \{ \lambda \in \Re |\lambda >0 \right \}$, which is the case of thermal networks, the above condition becomes
$$- \lambda_i \Delta t - 1 < 1, \forall \lambda_i, $$or
$$ 0 < \Delta t < -\frac{2}{\min \lambda_i} = 2 \min -\frac{1}{\lambda_i} = 2 \min T_i$$where $T_i$ are the time constants, $T_i = - \frac{1}{\lambda_i} $
λ = np.linalg.eig(As)[0] # eigenvalues of matrix As
λ = np.sort(λ)
print('Time constants:')
print([f'{T:.2f} s' for T in -1 / λ])
print('\n2 x Time constants:')
print([f'{T:.2f} s' for T in -2 / λ])
dtmax = 2 * min(-1. / λ)
print(f'\nMaximum time step: {dtmax:.2f} s = {dtmax / 60:.2f} min')
Time constants: ['249.30 s', '4093.20 s', '6729.11 s', '44033.06 s'] 2 x Time constants: ['498.60 s', '8186.41 s', '13458.22 s', '88066.12 s'] Maximum time step: 498.60 s = 8.31 min
Let's chose a time step smaller than $\Delta t_{max} = \min (-2 / \lambda_i) $.
# time step
dt = np.floor(dtmax / 60) * 60 # s
print(f'dt = {dt} s = {dt / 60:.0f} min')
dt = 480.0 s = 8 min
The settling time is roughly 4 times the larger time constant.
# settling time
time_const = np.array([int(x) for x in sorted(-1 / λ)])
print('4 * Time constants: \n', 4 * time_const, 's \n')
t_settle = 4 * max(-1 / λ)
print(f'Settling time: \
{t_settle:.0f} s = \
{t_settle / 60:.1f} min = \
{t_settle / (3600):.2f} h = \
{t_settle / (3600 * 24):.2f} days')
4 * Time constants: [ 996 16372 26916 176132] s Settling time: 176132 s = 2935.5 min = 48.93 h = 2.04 days
Let's obtain the dynamic response of the system to a step input.
The duration of the simulation needs to be larger than the estimated settling time. This requires a corresponding number of time steps in the time vector.
# Step response
# -------------
# duration: next multiple of 3600 s that is larger than t_settle
duration = np.ceil(t_settle / 3600) * 3600
n = int(np.floor(duration / dt)) # number of time steps
t = np.arange(0, n * dt, dt) # time vector for n time steps
print(f'Duration = {duration} s')
print(f'Number of time steps = {n}')
# pd.DataFrame(t, columns=['time'])
Duration = 176400.0 s Number of time steps = 367
In dynamic simulation, the inputs are time series, e.g., the oudoor temperature will have $n$ values $T_o = [T_{o(0)}, T_{o(1)}, ..., T_{o(n-1)}]$ at discrete time $t = [t_0, t_1, ... , t_{n-1}]$.
The input vector $u$ of the state-space representation is obtained by stacking the vectors $b_T$ and $f_Q$ of the system of Differential Algebraic Equations:
$$u = \begin{bmatrix} b_T \\ f_Q\end{bmatrix}$$where:
and $$b_T = [T_o, T_o, T_o, T_{i,sp}]^T$$ corresponding to branches 0, 8, 10 and 11;
and
$$f_Q = [\Phi_o, \Phi_i, \dot{Q}_a, \Phi_a]^T$$corresponding to nodes 0, 4, 6, and 7.
For the thermal circuit shown in Figure 3, the time series of the input vector, $u = [u_0, u_1, ... , u_{n-1}]^T$, is:
$$u = \begin{bmatrix} T_o\\ T_o\\ T_o\\ T_{i,sp}\\ \Phi_o\\ \Phi_i\\ \dot{Q}_a\\ \Phi_a \end{bmatrix} = \begin{bmatrix} T_{o(0)} & T_{o(1)}& ... & T_{o(n-1)}\\ T_{o(0)} & T_{o(1)}& ... & T_{o(n-1)}\ \\ T_{o(0)} & T_{o(1)}& ... & T_{o(n-1)}\ \\ T_{i,sp(0)} & T_{i,sp(1)}& ... & T_{i,sp(n-1)}\ \\ \Phi_{o(0)} & \Phi_{o(1)} & ... & \Phi_{o(n-1)}\\ \Phi_{i(0)} & \Phi_{i(1)} & ... & \Phi_{i(n-1)}\\ \dot{Q}_{a(0)} & \dot{Q}_{a(1)} & ... & \dot{Q}_{a(n-1)}\\ \Phi_{a(0)} & \Phi_{a(1)} & ... & \Phi_{a(n-1)} \end{bmatrix}$$where:
$T_o = [T_{o(0)}, T_{o(1)}, ..., T_{o(n-1)}]$ is the time series of the oudoor temperature at discrete time $t = [t_0, t_1, ... , t_{n-1}]$.
$T_{i, sp} = [T_{{i, sp}(0)}, T_{{i, sp}(1)}, ..., T_{{i, sp}(n-1)}]$ is the time series of the setpoint indoor temperature at discrete time $t = [t_0, t_1, ... , t_{n-1}]$.
$\Phi_o = [\Phi_{o(0)}, \Phi_{o(1)}, ..., \Phi_{o(n-1)}]$ is the time series of the solar radiation absorbed by the outdoor surface of the wall at discrete time $t = [t_0, t_1, ... , t_{n-1}]$.
$\Phi_i = [\Phi_{i(0)}, \Phi_{i(1)}, ..., \Phi_{i(n-1)}]$ is the time series of the solar radiation absorbed by the indoor surface of the wall at discrete time $t = [t_0, t_1, ... , t_{n-1}]$.
$\dot{Q}_a = [\dot{Q}_{a(0)}, \dot{Q}_{a(1)}, ..., \dot{Q}_{a(n-1)}]$ is the time series of the auxiliary heat gains (i.e., occupants, electrical devices, etc.) at discrete time $t = [t_0, t_1, ... , t_{n-1}]$.
$\Phi_a = [\Phi_{a(0)}, \Phi_{a(1)}, ..., \Phi_{a(n-1)}]$ is the time series of the solar radiation absorbed by the glass at discrete time $t = [t_0, t_1, ... , t_{n-1}]$.
Let's consider a step response in the conditions used for steady-state analysis, i.e. $T_o = 10 \, \mathrm{^\circ C}$, $T_{i,sp} = 20 \, \mathrm{^\circ C}$, and all the flow sources zero (including the HVAC system).
# input vector
#u = np.zeros([8, n]) # u = [To To To Tisp Φo Φi Qa Φa]
#u[0:3, :] = 10 * np.ones([3, n]) # To = 10 for n time steps
#u[3, :] = 20 * np.ones([1, n]) # Tisp = 20 for n time steps
#pd.DataFrame(u)
# input vector
u = np.zeros([n, 8]) # u = [To To To Tisp Φo Φi Qa Φa]
u[:, 0:3] = 10 * np.ones([n, 3]) # To = 10 for n time steps
u[:, 3] = 20 * np.ones([n]) # Tisp = 20 for n time steps
pd.DataFrame(u)
0 | 1 | 2 | 3 | 4 | 5 | 6 | 7 | |
---|---|---|---|---|---|---|---|---|
0 | 10.0 | 10.0 | 10.0 | 20.0 | 0.0 | 0.0 | 0.0 | 0.0 |
1 | 10.0 | 10.0 | 10.0 | 20.0 | 0.0 | 0.0 | 0.0 | 0.0 |
2 | 10.0 | 10.0 | 10.0 | 20.0 | 0.0 | 0.0 | 0.0 | 0.0 |
3 | 10.0 | 10.0 | 10.0 | 20.0 | 0.0 | 0.0 | 0.0 | 0.0 |
4 | 10.0 | 10.0 | 10.0 | 20.0 | 0.0 | 0.0 | 0.0 | 0.0 |
... | ... | ... | ... | ... | ... | ... | ... | ... |
362 | 10.0 | 10.0 | 10.0 | 20.0 | 0.0 | 0.0 | 0.0 | 0.0 |
363 | 10.0 | 10.0 | 10.0 | 20.0 | 0.0 | 0.0 | 0.0 | 0.0 |
364 | 10.0 | 10.0 | 10.0 | 20.0 | 0.0 | 0.0 | 0.0 | 0.0 |
365 | 10.0 | 10.0 | 10.0 | 20.0 | 0.0 | 0.0 | 0.0 | 0.0 |
366 | 10.0 | 10.0 | 10.0 | 20.0 | 0.0 | 0.0 | 0.0 | 0.0 |
367 rows × 8 columns
The state-space model
$$\left\{\begin{array}{rr} \dot{\theta}_C=A_s \theta_C + B_s u\\ y = C_s \theta_C + D_s u \end{array}\right.$$is integrated in time by using Euler forward (or explicit) method for numerical integration:
$$ \theta_{s,k+1} = (I + \Delta t A) \theta_{s,k} + \Delta t B u_k $$and Euler backward (or implicit) method for numerical integration:
$$\theta_{s,k+1} = (I - \Delta t A)^{-1} ( \theta_{s,k} + \Delta t B u_k )$$where $k = 0, ... , n - 1$.
# initial conditions
n_s = As.shape[0] # number of state variables
θ_exp = np.zeros([n_s, t.shape[0]]) # explicit Euler in time t
θ_imp = np.zeros([n_s, t.shape[0]]) # implicit Euler in time t
# time integration
I = np.eye(n_s) # identity matrix
for k in range(n - 1):
θ_exp[:, k + 1] = (I + dt * As) @ θ_exp[:, k]\
+ dt * Bs @ u[k, :]
θ_imp[:, k + 1] = np.linalg.inv(I - dt * As) @ (θ_imp[:, k]\
+ dt * Bs @ u[k, :])
Then, we obtain the outputs
$$ y = C_s \theta_s + D_s u$$for explicit and for implicit Euler methods, respectively.
# outputs
y_exp = Cs @ θ_exp + Ds @ u.T
y_imp = Cs @ θ_imp + Ds @ u.T
The results of explicit and implicit Euler integration are practically identical.
fig, ax = plt.subplots()
ax.plot(t / 3600, y_exp.T, t / 3600, y_imp.T)
ax.set(xlabel='Time, $t$ / h',
ylabel='Temperatue, $θ_i$ / °C',
title='Step input: outdoor temperature $T_o$')
ax.legend(['Explicit', 'Implicit'])
ax.grid()
plt.show()
Figure 7. Step response to outdoor temperature by using Euler
implicit and explicit integration.
The value the indoor temperature obtained after the settling time is almost equal to the value obtained in steady-state.
print('Steady-state indoor temperature obtained with:')
print(f'- DAE model: {float(θ[6]):.4f} °C')
print(f'- state-space model: {float(yss):.4f} °C')
print(f'- steady-state response to step input: {float(y_exp[:, -2]):.4f} °C')
Steady-state indoor temperature obtained with: - DAE model: 10.0000 °C - state-space model: 10.0000 °C - steady-state response to step input: 9.9640 °C
The simulation will be done from start_date
to end_date
indicated in the format MM-DD HH:MM:SS
(month, day, hour:minute:second).
start_date = '02-01 12:00:00'
end_date = '02-07 18:00:00'
The weather data are for a year. The choice of 2000
for the year is arbitrary; it used in order to respect the format YYYY-MM-DD HH:MM:SS
.
start_date = '2000-' + start_date
end_date = '2000-' + end_date
print(f'{start_date} \tstart date')
print(f'{end_date} \tend date')
2000-02-01 12:00:00 start date 2000-02-07 18:00:00 end date
Dynamic simulation needs time series of weather data for air temperature, direct solar radiation on a normal surface and diffuse solar radiation on an horizontal surface (see the tutorial on Weather data and solar radiation).
From the weather data, we select:
filename = './weather_data/FRA_Lyon.074810_IWEC.epw'
[data, meta] = dm4bem.read_epw(filename, coerce_year=None)
weather = data[["temp_air", "dir_n_rad", "dif_h_rad"]]
del data
The year is set to 2000
by convention and the data is selected from start to end.
weather.index = weather.index.map(lambda t: t.replace(year=2000))
weather = weather.loc[start_date:end_date]
For the surface orientation given by slope
, azimuth
and latitude
, and the albedo
of the surface in front of the wall, by using the weather data, we can calculate the:
for hourly solar irradiance on a tilted surface.
surface_orientation = {'slope': 90,
'azimuth': 0,
'latitude': 45}
albedo = 0.2
rad_surf = dm4bem.sol_rad_tilt_surf(
weather, surface_orientation, albedo)
# pd.DataFrame(rad_surf)
The total solar irradiance $E_{tot}$, in W/m², is the sum of direct, diffuse, and reflected components.
rad_surf['Φtot'] = rad_surf.sum(axis=1)
The weather data is at the time-step of 1h. It needs to be resampled at time step $\Delta t$ used for numerical integration.
# resample weather data
data = pd.concat([weather['temp_air'], rad_surf['Φtot']], axis=1)
data = data.resample(str(dt) + 'S').interpolate(method='linear')
data = data.rename(columns={'temp_air': 'To'})
data = data.rename_axis('Time')
# pd.DataFrame(data)
Let's consider the indoor temperature setpoint $T_{i,sp} = 20 \, \mathrm{^\circ C}$ and the auxiliary heat flow $\dot{Q}_a = 0 \, \mathrm{W}$ constant for the whole duration of the simulation.
data['Ti'] = 20 * np.ones(data.shape[0])
data['Qa'] = 0 * np.ones(data.shape[0])
# pd.DataFrame(data)
The input is formed by the vectors of time series of temperature sources $\left [ T_o, T_o ,T_o, T_{i,sp} \right ]^T$ and vectors of time series of the heat flow sources $\left [ \Phi_o, \Phi_i, \dot{Q_a}, \Phi_a \right ]^T$:
$$u = \begin{bmatrix} T_o\\ T_o\\ T_o\\ T_{i,sp}\\ \Phi_o\\ \Phi_i\\ \dot{Q}_a\\ \Phi_a \end{bmatrix} = \begin{bmatrix} T_{o(0)} & T_{o(1)}& ... & T_{o(n-1)}\\ T_{o(0)} & T_{o(1)}& ... & T_{o(n-1)}\ \\ T_{o(0)} & T_{o(1)}& ... & T_{o(n-1)}\ \\ T_{i,sp(0)} & T_{i,sp(1)}& ... & T_{i,sp(n-1)}\ \\ \Phi_{o,(0)} & \Phi_{o,(1)} & ... & \Phi_{o,(n-1)}\\ \Phi_{i,(0)} & \Phi_{i,(1)} & ... & \Phi_{i,(n-1)}\\ \dot{Q}_{a(0)} & \dot{Q}_{a(1)} & ... & \dot{Q}_{a(n-1)}\\ \Phi_{a,(0)} & \Phi_{a,(1)} & ... & \Phi_{a,(n-1)} \end{bmatrix}$$where:
$T_o$: the time series vector of outdoor temperatures (from weather data), °C.
$T_{i,sp}$: time series vector of indoor setpoint temperatures, °C.
$\Phi_o$: time series vector of solar (i.e. short wave) radiation, in W, absorbed by the outdoor surface of the wall:
$$\Phi_o = \alpha_{w,SW} S_w E_{tot}$$where:
$\Phi_i$: time series vector of short wave (i.e. solar) radiation, in W, absorbed by the indoor surfaces of the wall:
$$\Phi_i = \tau_{g,SW} \alpha_{w,SW} S_g E_{tot}$$where:
$\dot{Q}_a$: time vector of auxiliary heat flows (from occupants, electrical devices, etc.), W.
$\Phi_a$: time series vector of short wave (i.e. solar) radiation, in W, absorbed by the window glass:
$$\Phi_a = \alpha_{g,SW} S_g E_{tot}$$where:
# input vector
To = data['To']
Ti = data['Ti']
Φo = α_wSW * wall['Surface']['Layer_out'] * data['Φtot']
Φi = τ_gSW * α_wSW * wall['Surface']['Glass'] * data['Φtot']
Qa = data['Qa']
Φa = α_gSW * wall['Surface']['Glass'] * data['Φtot']
u = pd.concat([To, To, To, Ti, Φo, Φi, Qa, Φa], axis=1)
u.columns.values[[4, 5, 7]] = ['Φo', 'Φi', 'Φa']
pd.DataFrame(u)
To | To | To | Ti | Φo | Φi | Qa | Φa | |
---|---|---|---|---|---|---|---|---|
Time | ||||||||
2000-02-01 12:00:00+01:00 | 10.00 | 10.00 | 10.00 | 20.0 | 803.250000 | 48.195000 | 0.0 | 244.188000 |
2000-02-01 12:08:00+01:00 | 10.20 | 10.20 | 10.20 | 20.0 | 975.898057 | 58.553883 | 0.0 | 296.673009 |
2000-02-01 12:16:00+01:00 | 10.40 | 10.40 | 10.40 | 20.0 | 1148.546115 | 68.912767 | 0.0 | 349.158019 |
2000-02-01 12:24:00+01:00 | 10.60 | 10.60 | 10.60 | 20.0 | 1321.194172 | 79.271650 | 0.0 | 401.643028 |
2000-02-01 12:32:00+01:00 | 10.80 | 10.80 | 10.80 | 20.0 | 1493.842229 | 89.630534 | 0.0 | 454.128038 |
... | ... | ... | ... | ... | ... | ... | ... | ... |
2000-02-07 17:28:00+01:00 | 4.68 | 4.68 | 4.68 | 20.0 | 163.144887 | 9.788693 | 0.0 | 49.596046 |
2000-02-07 17:36:00+01:00 | 4.56 | 4.56 | 4.56 | 20.0 | 122.358665 | 7.341520 | 0.0 | 37.197034 |
2000-02-07 17:44:00+01:00 | 4.44 | 4.44 | 4.44 | 20.0 | 81.572443 | 4.894347 | 0.0 | 24.798023 |
2000-02-07 17:52:00+01:00 | 4.32 | 4.32 | 4.32 | 20.0 | 40.786222 | 2.447173 | 0.0 | 12.399011 |
2000-02-07 18:00:00+01:00 | 4.20 | 4.20 | 4.20 | 20.0 | 0.000000 | 0.000000 | 0.0 | 0.000000 |
1126 rows × 8 columns
The initial value of the state-vector can be zero or different from zero.
θ_exp = 20 * np.ones([As.shape[0], u.shape[0]])
Explicit Euler integration in time,
$$ \theta_{s,k+1} = (I + \Delta t A) \theta_{s,k} + \Delta t B u_k $$where $k = 0, ... , n - 1$,
for k in range(u.shape[0] - 1):
θ_exp[:, k + 1] = (I + dt * As) @ θ_exp[:, k]\
+ dt * Bs @ u.iloc[k, :]
yields the time variation of state variable $\theta$, from which we obtain the variation of the output (i.e. indoor temperature):
$$y = C_s \theta_s + D_s u$$and the variation of the heat flow of the HVAC system:
$$q_{HVAC} = K_p (T_{i,sp} - \theta_i) = K_p (T_{i,sp} - y)$$where $K_p$ is the gain of the P-controller and $T_{i,sp}$ is the HVAC-setpoint for the indoor temperature.
y_exp = Cs @ θ_exp + Ds @ u.to_numpy().T
q_HVAC = Kp * (data['Ti'] - y_exp[0, :])
data['θi_exp'] = y_exp.T
data['q_HVAC'] = q_HVAC.T
fig, axs = plt.subplots(2, 1)
data[['To', 'θi_exp']].plot(ax=axs[0],
xticks=[],
ylabel='Temperature, $θ$ / °C')
axs[0].legend(['$θ_{outdoor}$', '$θ_{indoor}$'],
loc='upper right')
data[['Φtot', 'q_HVAC']].plot(ax=axs[1],
ylabel='Heat rate, $q$ / W')
axs[1].set(xlabel='Time')
axs[1].legend(['$Φ_{total}$', '$q_{HVAC}$'],
loc='upper right')
plt.show()
Figure 6. Simulation in free-running with weather data using Euler explicit method of integration. a) Indoor and outdoor temperatures. b) Solar and HVAC heat flow rates.
t = dt * np.arange(data.shape[0]) # time vector
fig, axs = plt.subplots(2, 1)
# plot outdoor and indoor temperature
axs[0].plot(t / 3600 / 24, data['To'], label='$θ_{outdoor}$')
axs[0].plot(t / 3600 / 24, y_exp[0, :], label='$θ_{indoor}$')
axs[0].set(ylabel='Temperatures, $θ$ / °C',
title='Simulation for weather')
axs[0].legend(loc='upper right')
# plot total solar radiation and HVAC heat flow
axs[1].plot(t / 3600 / 24, data['Φtot'], label='$Φ_{total}$')
axs[1].plot(t / 3600 / 24, q_HVAC, label='$q_{HVAC}$')
axs[1].set(xlabel='Time, $t$ / day',
ylabel='Heat flows, $q$ / W')
axs[1].legend(loc='upper right')
fig.tight_layout()
Figure 7. Simulation in free-running with weather data using Euler explicit method of integration. a) Indoor and outdoor temperatures. b) Solar and HVAC heat flow rates.
Interchange the materials of the layers of the wall. Discuss the step responses and the simuation for weather. Give arguments for the advantages and the disadvanted of indoor and outdoor insulation.
The time step depends on:
Kp
:C['Air']
and of the glass $C_g =$ C['Glass']
are considered, then the time step is small;The controller models an HVAC system able to heat (when $q_{HVAC} > 0$) and to cool (when $q_{HVAC} < 0$).
C. Ghiaus (2013) Causality issue in the heat balance method for calculating the design heating and cooling loads, Energy 50: 292-301, https://doi.org/10.1016/j.energy.2012.10.024, open access preprint: HAL-03605823
C. Ghiaus (2021). Dynamic Models for Energy Control of Smart Homes, in S. Ploix M. Amayri, N. Bouguila (eds.) Towards Energy Smart Homes, Online ISBN: 978-3-030-76477-7, Print ISBN: 978-3-030-76476-0, Springer, pp. 163-198 (ref.)
DOI 10.1007/978-3-030-76477-7_5, open access preprint: HAL 03578578
J.A. Duffie, W. A. Beckman, N. Blair (2020) Solar Engineering of Thermal Processes, 5th ed. John Wiley & Sons, Inc. ISBN 9781119540281
Réglementation Thermique 2005. Méthode de calcul Th-CE. Annexe à l’arrêté du 19 juillet 2006
H. Recknagel, E. Sprenger, E.-R. Schramek (2013) Génie climatique, 5e edition, Dunod, Paris. ISBN 978-2-10-070451-4
J.R. Howell et al. (2021) Thermal Radiation Heat Transfer 7th edition, ISBN 978-0-367-34707-0, A Catalogue of Configuration Factors
J. Widén, J. Munkhammar (2019) Solar Radiation Theory, Uppsala University