g-Function Module

pygfunction.gfunction.equal_inlet_temperature(boreholes, UTubes, m_flow_borehole, cp_f, time, alpha, kind='linear', nSegments=8, segment_ratios=<function segment_ratios>, use_similarities=True, disTol=0.01, tol=1e-06, dtype=<class 'numpy.float64'>, disp=False, **kwargs)

Evaluate the g-function with equal inlet fluid temperatures.

This function superimposes the finite line source (FLS) solution to estimate the g-function of a geothermal bore field. Each borehole is modeled as a series of finite line source segments, as proposed in 1.

Parameters
boreholeslist of Borehole objects

List of boreholes included in the bore field.

UTubeslist of pipe objects

Model for pipes inside each borehole.

m_flow_boreholefloat or array

Fluid mass flow rate per borehole (in kg/s).

cp_ffloat

Fluid specific isobaric heat capacity (in J/kg.K).

timefloat or array

Values of time (in seconds) for which the g-function is evaluated.

alphafloat

Soil thermal diffusivity (in m2/s).

nSegmentsint or list, optional

Number of line segments used per borehole, or list of number of line segments used for each borehole. Default is 8.

segment_ratiosarray, list of arrays, or callable, optional

Ratio of the borehole length represented by each segment. The sum of ratios must be equal to 1. The shape of the array is of (nSegments,) or list of (nSegments[i],). If segment_ratios==None, segments of equal lengths are considered. If a callable is provided, it must return an array of size (nSegments,) when provided with nSegments (of type int) as an argument, or an array of size (nSegments[i],) when provided with an element of nSegments (of type list). Default is utilities.segment_ratios().

kindstring, optional

Interpolation method used for segment-to-segment thermal response factors. See documentation for scipy.interpolate.interp1d. Default is ‘linear’.

use_similaritiesbool, optional

True if similarities are used to limit the number of FLS evaluations. Default is True.

disTolfloat, optional

Relative tolerance on radial distance. Two distances (d1, d2) between two pairs of boreholes are considered equal if the difference between the two distances (abs(d1-d2)) is below tolerance. Default is 0.01.

tolfloat, optional

Relative tolerance on length and depth. Two lengths H1, H2 (or depths D1, D2) are considered equal if abs(H1 - H2)/H2 < tol. Default is 1.0e-6.

dtypenumpy dtype, optional

numpy data type used for matrices and vectors. Should be one of numpy.single or numpy.double. Default is numpy.double.

dispbool, optional

Set to true to print progression messages. Default is False.

Returns
gFunctionfloat or array

Values of the g-function

References

1

Cimmino, M. (2015). The effects of borehole thermal resistances and fluid flow rate on the g-functions of geothermal bore fields. International Journal of Heat and Mass Transfer, 91, 1119-1127.

class pygfunction.gfunction.gFunction(boreholes_or_network, alpha, time=None, method='equivalent', boundary_condition=None, options={})

Bases: object

Class for the calculation and visualization of the g-functions of geothermal bore fields.

This class superimposes the finite line source (FLS) solution to estimate the g-function of a geothermal bore field. Each borehole is modeled as a series of finite line source segments, as proposed in 2. Different boundary conditions and solution methods are implemented.

Notes

  • The ‘equivalent’ solver does not support the ‘MIFT’ boundary condition when boreholes are connected in series.

  • The ‘equivalent’ solver does not support inclined boreholes.

  • The g-function is linearized for times t < r_b**2 / (25 * self.alpha). The g-function value is then interpolated between 0 and its value at the threshold.

References

2(1,2,3)

Cimmino, M., & Bernier, M. (2014). A semi-analytical method to generate g-functions for geothermal bore fields. International Journal of Heat and Mass Transfer, 70, 641-650.

3

Cimmino, M. (2015). The effects of borehole thermal resistances and fluid flow rate on the g-functions of geothermal bore fields. International Journal of Heat and Mass Transfer, 91, 1119-1127.

4

Cimmino, M. (2018). Fast calculation of the g-functions of geothermal borehole fields using similarities in the evaluation of the finite line source solution. Journal of Building Performance Simulation, 11 (6), 655-668.

5

Cimmino, M. (2019). Semi-analytical method for g-function calculation of bore fields with series- and parallel-connected boreholes. Science and Technology for the Built Environment, 25 (8), 1007-1022.

6

Prieto, C., & Cimmino, M. (2021). Thermal interactions in large irregular fields of geothermal boreholes: the method of equivalent borehole. Journal of Building Performance Simulation, 14 (4), 446-460.

7

Cimmino, M. (2021). An approximation of the finite line source solution to model thermal interactions between geothermal boreholes. International Communications in Heat and Mass Transfer, 127, 105496.

Attributes
boreholes_or_networklist of Borehole objects or Network object

List of boreholes included in the bore field, or network of boreholes and pipes.

alphafloat

Soil thermal diffusivity (in m2/s).

timefloat or array, optional

Values of time (in seconds) for which the g-function is evaluated. The g-function is only evaluated at initialization if a value is provided. Default is None.

methodstr, optional

Method for the evaluation of the g-function. Should be one of

  • ‘similarities’ :

    The accelerated method of Cimmino (2018) 4, using similarities in the bore field to decrease the number of evaluations of the FLS solution.

  • ‘detailed’ :

    The classical superposition of the FLS solution. The FLS solution is evaluated for all pairs of segments in the bore field.

  • ‘equivalent’ :

    The equivalent borehole method of Prieto and Cimmino (2021) 6. Boreholes are assembled into groups of boreholes that share similar borehole wall temperatures and heat extraction rates. Each group is represented by an equivalent borehole and the group-to-group thermal interactions are calculated by the FLS solution. This is an approximation of the ‘similarities’ method.

Default is ‘equivalent’.

boundary_conditionstr, optional

Boundary condition for the evaluation of the g-function. Should be one of

  • ‘UHTR’ :

    Uniform heat transfer rate. This is corresponds to boundary condition BC-I as defined by Cimmino and Bernier (2014) 2.

  • ‘UBWT’ :

    Uniform borehole wall temperature. This is corresponds to boundary condition BC-III as defined by Cimmino and Bernier (2014) 2.

  • ‘MIFT’ :

    Mixed inlet fluid temperatures. This boundary condition was introduced by Cimmino (2015) 3 for parallel-connected boreholes and extended to mixed configurations by Cimmino (2019) 5.

If not given, chosen to be ‘UBWT’ if a list of boreholes is provided or ‘MIFT’ if a Network object is provided.

optionsdict, optional

A dictionary of solver options. All methods accept the following generic options:

nSegmentsint or list, optional

Number of line segments used per borehole, or list of number of line segments used for each borehole. Default is 8.

segment_ratiosarray, list of arrays, or callable, optional

Ratio of the borehole length represented by each segment. The sum of ratios must be equal to 1. The shape of the array is of (nSegments,) or list of (nSegments[i],). If segment_ratios==None, segments of equal lengths are considered. If a callable is provided, it must return an array of size (nSegments,) when provided with nSegments (of type int) as an argument, or an array of size (nSegments[i],) when provided with an element of nSegments (of type list). Default is utilities.segment_ratios().

approximate_FLSbool, optional

Set to true to use the approximation of the FLS solution of Cimmino (2021) 7. This approximation does not require the numerical evaluation of any integral. When using the ‘equivalent’ solver, the approximation is only applied to the thermal response at the borehole radius. Thermal interaction between boreholes is evaluated using the FLS solution. Default is False.

nFLSint, optional

Number of terms in the approximation of the FLS solution. This parameter is unused if approximate_FLS is set to False. Default is 10. Maximum is 25.

mQuadint, optional

Number of Gauss-Legendre sample points for the integral over \(u\) in the inclined FLS solution. Default is 11.

linear_thresholdfloat, optional

Threshold time (in seconds) under which the g-function is linearized. The g-function value is then interpolated between 0 and its value at the threshold. If linear_threshold==None, the g-function is linearized for times t < r_b**2 / (25 * self.alpha). Default is None.

dispbool, optional

Set to true to print progression messages. Default is False.

profilesbool, optional

Set to true to keep in memory the temperatures and heat extraction rates. Default is False.

kindstr, optional

Interpolation method used for segment-to-segment thermal response factors. See documentation for scipy.interpolate.interp1d. Default is ‘linear’.

dtypenumpy dtype, optional

numpy data type used for matrices and vectors. Should be one of numpy.single or numpy.double. Default is numpy.double.

The ‘similarities’ solver accepts the following method-specific options:

disTolfloat, optional

Relative tolerance on radial distance. Two distances (d1, d2) between two pairs of boreholes are considered equal if the difference between the two distances (abs(d1-d2)) is below tolerance. Default is 0.01.

tolfloat, optional

Relative tolerance on length and depth. Two lengths H1, H2 (or depths D1, D2) are considered equal if abs(H1 - H2)/H2 < tol. Default is 1.0e-6.

The ‘equivalent’ solver accepts the following method-specific options:

disTolfloat, optional

Relative tolerance on radial distance. Two distances (d1, d2) between two pairs of boreholes are considered equal if the difference between the two distances (abs(d1-d2)) is below tolerance. Default is 0.01.

tolfloat, optional

Relative tolerance on length and depth. Two lengths H1, H2 (or depths D1, D2) are considered equal if abs(H1 - H2)/H2 < tol. Default is 1.0e-6.

kClustersint, optional

Increment on the minimum number of equivalent boreholes determined by cutting the dendrogram of the bore field given by the hierarchical agglomerative clustering method. Increasing the value of this parameter increases the accuracy of the method. Default is 1.

evaluate_g_function(time)

Evaluate the g-function.

Parameters
timefloat or array

Values of time (in seconds) for which the g-function is evaluated.

Returns
gFunctionfloat or array

Values of the g-function

visualize_g_function()

Plot the g-function of the borefield.

Returns
figfigure

Figure object (matplotlib).

visualize_heat_extraction_rate_profiles(time=None, iBoreholes=None, showTilt=True)

Plot the heat extraction rate profiles at chosen time.

Parameters
timefloat

Values of time (in seconds) to plot heat extraction rate profiles. If time is None, heat extraction rates are plotted at the last time step. Default is None.

iBoreholeslist of int

Borehole indices to plot heat extraction rate profiles. If iBoreholes is None, heat extraction rates are plotted for all boreholes. Default is None.

showTiltbool

Set to True to show borehole inclination. Default is True

Returns
figfigure

Figure object (matplotlib).

visualize_heat_extraction_rates(iBoreholes=None, showTilt=True)

Plot the time-variation of the average heat extraction rates.

Parameters
iBoreholeslist of int

Borehole indices to plot heat extraction rates. If iBoreholes is None, heat extraction rates are plotted for all boreholes. Default is None.

showTiltbool

Set to True to show borehole inclination. Default is True

Returns
figfigure

Figure object (matplotlib).

visualize_temperature_profiles(time=None, iBoreholes=None, showTilt=True)

Plot the borehole wall temperature profiles at chosen time.

Parameters
timefloat

Values of time (in seconds) to plot temperature profiles. If time is None, temperatures are plotted at the last time step. Default is None.

iBoreholeslist of int

Borehole indices to plot temperature profiles. If iBoreholes is None, temperatures are plotted for all boreholes. Default is None.

showTiltbool

Set to True to show borehole inclination. Default is True

Returns
figfigure

Figure object (matplotlib).

visualize_temperatures(iBoreholes=None, showTilt=True)

Plot the time-variation of the average borehole wall temperatures.

Parameters
iBoreholeslist of int

Borehole indices to plot temperatures. If iBoreholes is None, temperatures are plotted for all boreholes. Default is None.

showTiltbool

Set to True to show borehole inclination. Default is True

Returns
figfigure

Figure object (matplotlib).

pygfunction.gfunction.mixed_inlet_temperature(network, m_flow_network, cp_f, time, alpha, kind='linear', nSegments=8, segment_ratios=<function segment_ratios>, use_similarities=True, disTol=0.01, tol=1e-06, dtype=<class 'numpy.float64'>, disp=False, **kwargs)

Evaluate the g-function with mixed inlet fluid temperatures.

This function superimposes the finite line source (FLS) solution to estimate the g-function of a geothermal bore field. Each borehole is modeled as a series of finite line source segments, as proposed in 8. The piping configurations between boreholes can be any combination of series and parallel connections.

Parameters
networkNetwork objects

List of boreholes included in the bore field.

m_flow_networkfloat or array

Total mass flow rate into the network or inlet mass flow rates into each circuit of the network (in kg/s). If a float is supplied, the total mass flow rate is split equally into all circuits.

cp_ffloat or array

Fluid specific isobaric heat capacity (in J/kg.degC). Must be the same for all circuits (a single float can be supplied).

timefloat or array

Values of time (in seconds) for which the g-function is evaluated.

alphafloat

Soil thermal diffusivity (in m2/s).

nSegmentsint or list, optional

Number of line segments used per borehole, or list of number of line segments used for each borehole. Default is 8.

segment_ratiosarray, list of arrays, or callable, optional

Ratio of the borehole length represented by each segment. The sum of ratios must be equal to 1. The shape of the array is of (nSegments,) or list of (nSegments[i],). If segment_ratios==None, segments of equal lengths are considered. Default is None.

kindstring, optional

Interpolation method used for segment-to-segment thermal response factors. See documentation for scipy.interpolate.interp1d. Default is ‘linear’.

use_similaritiesbool, optional

True if similarities are used to limit the number of FLS evaluations. Default is True.

disTolfloat, optional

Relative tolerance on radial distance. Two distances (d1, d2) between two pairs of boreholes are considered equal if the difference between the two distances (abs(d1-d2)) is below tolerance. Default is 0.01.

tolfloat, optional

Relative tolerance on length and depth. Two lengths H1, H2 (or depths D1, D2) are considered equal if abs(H1 - H2)/H2 < tol. Default is 1.0e-6.

dtypenumpy dtype, optional

numpy data type used for matrices and vectors. Should be one of numpy.single or numpy.double. Default is numpy.double.

dispbool, optional

Set to true to print progression messages. Default is False.

Returns
gFunctionfloat or array

Values of the g-function

References

8

Cimmino, M. (2019). Semi-analytical method for g-function calculation of bore fields with series- and parallel-connected boreholes. Science and Technology for the Built Environment, 25 (8), 1007-1022

Examples

>>> b1 = gt.boreholes.Borehole(H=150., D=4., r_b=0.075, x=0., y=0.)
>>> b2 = gt.boreholes.Borehole(H=150., D=4., r_b=0.075, x=5., y=0.)
>>> Utube1 = gt.pipes.SingleUTube(pos=[(-0.05, 0), (0, -0.05)],
                                  r_in=0.015, r_out=0.02,
                                  borehole=b1,k_s=2, k_g=1, R_fp=0.1)
>>> Utube2 = gt.pipes.SingleUTube(pos=[(-0.05, 0), (0, -0.05)],
                                  r_in=0.015, r_out=0.02,
                                  borehole=b1,k_s=2, k_g=1, R_fp=0.1)
>>> bore_connectivity = [-1, 0]
>>> network = gt.networks.Network([b1, b2], [Utube1, Utube2], bore_connectivity)
>>> time = np.array([1.0*10**i for i in range(4, 12)])
>>> m_flow_network = 0.25
>>> cp_f = 4000.
>>> alpha = 1.0e-6
>>> gt.gfunction.mixed_inlet_temperature(network, m_flow_network, cp_f, time, alpha)
array([0.63782415, 1.63304116, 2.72191316, 4.04091713, 5.98240458,
   7.77216202, 8.66195828, 8.77567215])
pygfunction.gfunction.uniform_heat_extraction(boreholes, time, alpha, use_similarities=True, disTol=0.01, tol=1e-06, dtype=<class 'numpy.float64'>, disp=False, **kwargs)

Evaluate the g-function with uniform heat extraction along boreholes.

This function superimposes the finite line source (FLS) solution to estimate the g-function of a geothermal bore field. This boundary condition correponds to BC-I, as defined by 9.

Parameters
boreholeslist of Borehole objects

List of boreholes included in the bore field.

timefloat or array

Values of time (in seconds) for which the g-function is evaluated.

alphafloat

Soil thermal diffusivity (in m2/s).

use_similaritiesbool, optional

True if similarities are used to limit the number of FLS evaluations. Default is True.

disTolfloat, optional

Relative tolerance on radial distance. Two distances (d1, d2) between two pairs of boreholes are considered equal if the difference between the two distances (abs(d1-d2)) is below tolerance. Default is 0.01.

tolfloat, optional

Relative tolerance on length and depth. Two lengths H1, H2 (or depths D1, D2) are considered equal if abs(H1 - H2)/H2 < tol. Default is 1.0e-6.

dtypenumpy dtype, optional

numpy data type used for matrices and vectors. Should be one of numpy.single or numpy.double. Default is numpy.double.

dispbool, optional

Set to true to print progression messages. Default is False.

Returns
gFunctionfloat or array

Values of the g-function

References

9

Cimmino, M., & Bernier, M. (2014). A semi-analytical method to generate g-functions for geothermal bore fields. International Journal of Heat and Mass Transfer, 70, 641-650.

Examples

>>> b1 = gt.boreholes.Borehole(H=150., D=4., r_b=0.075, x=0., y=0.)
>>> b2 = gt.boreholes.Borehole(H=150., D=4., r_b=0.075, x=5., y=0.)
>>> alpha = 1.0e-6
>>> time = np.array([1.0*10**i for i in range(4, 12)])
>>> gt.gfunction.uniform_heat_extraction([b1, b2], time, alpha)
array([ 0.75978163,  1.84860837,  2.98861057,  4.33496051,  6.29199383,
    8.13636888,  9.08401497,  9.20736188])
pygfunction.gfunction.uniform_temperature(boreholes, time, alpha, nSegments=8, segment_ratios=<function segment_ratios>, kind='linear', use_similarities=True, disTol=0.01, tol=1e-06, dtype=<class 'numpy.float64'>, disp=False, **kwargs)

Evaluate the g-function with uniform borehole wall temperature.

This function superimposes the finite line source (FLS) solution to estimate the g-function of a geothermal bore field. Each borehole is modeled as a series of finite line source segments. This boundary condition correponds to BC-III, as defined by 10.

Parameters
boreholeslist of Borehole objects

List of boreholes included in the bore field.

timefloat or array

Values of time (in seconds) for which the g-function is evaluated.

alphafloat

Soil thermal diffusivity (in m2/s).

nSegmentsint or list, optional

Number of line segments used per borehole, or list of number of line segments used for each borehole. Default is 8.

segment_ratiosarray, list of arrays, or callable, optional

Ratio of the borehole length represented by each segment. The sum of ratios must be equal to 1. The shape of the array is of (nSegments,) or list of (nSegments[i],). If segment_ratios==None, segments of equal lengths are considered. If a callable is provided, it must return an array of size (nSegments,) when provided with nSegments (of type int) as an argument, or an array of size (nSegments[i],) when provided with an element of nSegments (of type list). Default is utilities.segment_ratios().

kindstring, optional

Interpolation method used for segment-to-segment thermal response factors. See documentation for scipy.interpolate.interp1d. Default is ‘linear’.

use_similaritiesbool, optional

True if similarities are used to limit the number of FLS evaluations. Default is True.

disTolfloat, optional

Relative tolerance on radial distance. Two distances (d1, d2) between two pairs of boreholes are considered equal if the difference between the two distances (abs(d1-d2)) is below tolerance. Default is 0.01.

tolfloat, optional

Relative tolerance on length and depth. Two lengths H1, H2 (or depths D1, D2) are considered equal if abs(H1 - H2)/H2 < tol. Default is 1.0e-6.

dtypenumpy dtype, optional

numpy data type used for matrices and vectors. Should be one of numpy.single or numpy.double. Default is numpy.double.

dispbool, optional

Set to true to print progression messages. Default is False.

Returns
gFunctionfloat or array

Values of the g-function

References

10

Cimmino, M., & Bernier, M. (2014). A semi-analytical method to generate g-functions for geothermal bore fields. International Journal of Heat and Mass Transfer, 70, 641-650.

Examples

>>> b1 = gt.boreholes.Borehole(H=150., D=4., r_b=0.075, x=0., y=0.)
>>> b2 = gt.boreholes.Borehole(H=150., D=4., r_b=0.075, x=5., y=0.)
>>> alpha = 1.0e-6
>>> time = np.array([1.0*10**i for i in range(4, 12)])
>>> gt.gfunction.uniform_temperature([b1, b2], time, alpha)
array([ 0.75978079,  1.84859851,  2.98852756,  4.33406497,  6.27830732,
    8.05746656,  8.93697282,  9.04925079])