Calculation of effective bore field thermal resistanceΒΆ
This example demonstrates the use of the
pipes.field_thermal_resistance()
function to evaluate the bore field
thermal resistance. The concept effective bore field thermal is detailed by
Cimmino 1
The following script evaluates the effective bore field thermal resistance for fields of 1 to 5 series-connected boreholes with fluid flow rates ranging from 0.01 kg/s ti 1.00 kg/s.
The script is located in: pygfunction/examples/bore_field_thermal_resistance.py
1# -*- coding: utf-8 -*-
2""" Example of calculation of effective bore field thermal resistance.
3
4 The effective bore field thermal resistance of fields of up to 5 boreholes
5 of equal lengths connected in series is calculated for various fluid flow
6 rates.
7
8"""
9import matplotlib.pyplot as plt
10import numpy as np
11from matplotlib.ticker import AutoMinorLocator
12from scipy import pi
13
14import pygfunction as gt
15
16
17def main():
18 # -------------------------------------------------------------------------
19 # Simulation parameters
20 # -------------------------------------------------------------------------
21
22 # Number of boreholes
23 nBoreholes = 5
24 # Borehole dimensions
25 D = 4.0 # Borehole buried depth (m)
26 # Borehole length (m)
27 H = 150.
28 r_b = 0.075 # Borehole radius (m)
29 B = 7.5 # Borehole spacing (m)
30
31 # Pipe dimensions
32 r_out = 0.02 # Pipe outer radius (m)
33 r_in = 0.015 # Pipe inner radius (m)
34 D_s = 0.05 # Shank spacing (m)
35 epsilon = 1.0e-6 # Pipe roughness (m)
36
37 # Pipe positions
38 # Single U-tube [(x_in, y_in), (x_out, y_out)]
39 pos_pipes = [(-D_s, 0.), (D_s, 0.)]
40
41 # Ground properties
42 k_s = 2.0 # Ground thermal conductivity (W/m.K)
43
44 # Grout properties
45 k_g = 1.0 # Grout thermal conductivity (W/m.K)
46
47 # Pipe properties
48 k_p = 0.4 # Pipe thermal conductivity (W/m.K)
49
50 # Fluid properties
51 # Total fluid mass flow rate per borehole (kg/s), from 0.01 kg/s to 1 kg/s
52 m_flow_network = 10**np.arange(-2, 0.001, 0.05)
53 # The fluid is propylene-glycol (20 %) at 20 degC
54 fluid = gt.media.Fluid('MPG', 20.)
55 cp_f = fluid.cp # Fluid specific isobaric heat capacity (J/kg.K)
56 rho_f = fluid.rho # Fluid density (kg/m3)
57 mu_f = fluid.mu # Fluid dynamic viscosity (kg/m.s)
58 k_f = fluid.k # Fluid thermal conductivity (W/m.K)
59
60 # -------------------------------------------------------------------------
61 # Borehole field
62 # -------------------------------------------------------------------------
63
64 boreField = []
65 bore_connectivity = []
66 for i in range(nBoreholes):
67 x = i*B
68 borehole = gt.boreholes.Borehole(H, D, r_b, x, 0.)
69 boreField.append(borehole)
70 # Boreholes are connected in series: The index of the upstream
71 # borehole is that of the previous borehole
72 bore_connectivity.append(i - 1)
73
74 # -------------------------------------------------------------------------
75 # Evaluate the effective bore field thermal resistance
76 # -------------------------------------------------------------------------
77
78 # Initialize result array
79 R = np.zeros((nBoreholes, len(m_flow_network)))
80 for i in range(nBoreholes):
81 for j, m_flow_network_j in enumerate(m_flow_network):
82 nBoreholes = i + 1
83 # Boreholes are connected in series
84 m_flow_borehole = m_flow_network_j
85 # Boreholes are single U-tube
86 m_flow_pipe = m_flow_borehole
87
88 # Pipe thermal resistance
89 R_p = gt.pipes.conduction_thermal_resistance_circular_pipe(
90 r_in, r_out, k_p)
91 # Fluid to inner pipe wall thermal resistance (Single U-tube)
92 h_f = gt.pipes.convective_heat_transfer_coefficient_circular_pipe(
93 m_flow_pipe, r_in, mu_f, rho_f, k_f, cp_f, epsilon)
94 R_f = 1.0/(h_f*2*pi*r_in)
95
96 # Single U-tube, same for all boreholes in the bore field
97 UTubes = []
98 for borehole in boreField:
99 SingleUTube = gt.pipes.SingleUTube(
100 pos_pipes, r_in, r_out, borehole, k_s, k_g, R_f + R_p)
101 UTubes.append(SingleUTube)
102 network = gt.networks.Network(
103 boreField[:nBoreholes],
104 UTubes[:nBoreholes],
105 bore_connectivity=bore_connectivity[:nBoreholes])
106
107 # Effective bore field thermal resistance
108 R_field = gt.networks.network_thermal_resistance(
109 network, m_flow_network_j, cp_f)
110 # Add to result array
111 R[i,j] = R_field
112
113 # -------------------------------------------------------------------------
114 # Plot bore field thermal resistances
115 # -------------------------------------------------------------------------
116
117 # Configure figure and axes
118 fig = gt.utilities._initialize_figure()
119
120 ax1 = fig.add_subplot(111)
121 # Axis labels
122 ax1.set_xlabel(r'$\dot{m}$ [kg/s]')
123 ax1.set_ylabel(r'$R^*_{field}$ [m.K/W]')
124 # Axis limits
125 ax1.set_xlim([0., 1.])
126 ax1.set_ylim([0., 1.])
127
128 gt.utilities._format_axes(ax1)
129
130 # Bore field thermal resistances
131 ax1.plot(m_flow_network, R[0,:], '-', label='1 borehole')
132 ax1.plot(m_flow_network, R[2,:], '--', label='3 boreholes')
133 ax1.plot(m_flow_network, R[4,:], '-.', label='5 boreholes')
134 ax1.legend()
135 # Adjust to plot window
136 plt.tight_layout()
137
138 return
139
140
141# Main function
142if __name__ == '__main__':
143 main()
References
- 1
Cimmino, M. (2018). g-Functions for bore fields with mixed parallel and series connections considering the axial fluid temperature variations. Proceedings of the IGSHPA Sweden Research Track 2018. Stockholm, Sweden. pp. 262-270.