Solvers Module

class pygfunction.solvers.Detailed(boreholes, network, time, boundary_condition, m_flow_borehole=None, m_flow_network=None, cp_f=None, nSegments=8, segment_ratios=<function segment_ratios>, approximate_FLS=False, mQuad=11, nFLS=10, linear_threshold=None, disp=False, profiles=False, kind='linear', dtype=<class 'numpy.float64'>, **other_options)

Bases: _BaseSolver

Detailed solver for the evaluation of the g-function.

This solver 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.

networknetwork object

Model of the network.

timefloat or array

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

boundary_conditionstr

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

  • ‘UHTR’ :

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

  • ‘UBWT’ :

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

  • ‘MIFT’ :

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

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().

m_flow_borehole(nInlets,) array or (nMassFlow, nInlets,) array, optional

Fluid mass flow rate into each circuit of the network. If a (nMassFlow, nInlets,) array is supplied, the (nMassFlow, nMassFlow,) variable mass flow rate g-functions will be evaluated using the method of Cimmino (2024) [5]. Only required for the ‘MIFT’ boundary condition. Only one of ‘m_flow_borehole’ and ‘m_flow_network’ can be provided. Default is None.

m_flow_networkfloat or (nMassFlow,) array, optional

Fluid mass flow rate into the network of boreholes. If an array is supplied, the (nMassFlow, nMassFlow,) variable mass flow rate g-functions will be evaluated using the method of Cimmino (2024) [5]. Only required for the ‘MIFT’ boundary condition. Only one of ‘m_flow_borehole’ and ‘m_flow_network’ can be provided. Default is None.

cp_ffloat, optional

Fluid specific isobaric heat capacity (in J/kg.degC). Only required for the ‘MIFT’ boundary condition. Default is None.

approximate_FLSbool, optional

Set to true to use the approximation of the FLS solution of Cimmino (2021) [4]. This approximation does not require the numerical evaluation of any integral. 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.

kindstring, 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.

References

initialize(**kwargs)

Split boreholes into segments.

Returns:
nSourcesint

Number of finite line heat sources in the borefield used to initialize the matrix of segment-to-segment thermal response factors (of size: nSources x nSources).

thermal_response_factors(time, alpha, kind='linear')

Evaluate the segment-to-segment thermal response factors for all pairs of segments in the borefield at all time steps using the finite line source solution.

This method returns a scipy.interpolate.interp1d object of the matrix of thermal response factors, containing a copy of the matrix accessible by h_ij.y[:nSources,:nSources,:nt+1]. The first index along the third axis corresponds to time t=0. The interp1d object can be used to obtain thermal response factors at any intermediate time by h_ij(t)[:nSources,:nSources].

Attributes:
timefloat or array

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

alphafloat

Soil thermal diffusivity (in m2/s).

kindstring, optional

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

Returns:
h_ijinterp1d

interp1d object (scipy.interpolate) of the matrix of segment-to-segment thermal response factors.

class pygfunction.solvers.Equivalent(boreholes, network, time, boundary_condition, m_flow_borehole=None, m_flow_network=None, cp_f=None, nSegments=8, segment_ratios=<function segment_ratios>, approximate_FLS=False, mQuad=11, nFLS=10, linear_threshold=None, disp=False, profiles=False, kind='linear', dtype=<class 'numpy.float64'>, **other_options)

Bases: _BaseSolver

Equivalent solver for the evaluation of the g-function.

This solver uses hierarchical agglomerative clustering to identify groups of boreholes that are expected to have similar borehole wall temperatures and heat extraction rates, as proposed by Prieto and Cimmino (2021) [9]. Each group of boreholes is represented by a single equivalent borehole. The FLS solution is adapted to evaluate thermal interactions between groups of boreholes. This greatly reduces the number of evaluations of the FLS solution and the size of the system of equations to evaluate the g-function.

Parameters:
boreholeslist of Borehole objects

List of boreholes included in the bore field.

networknetwork object

Model of the network.

timefloat or array

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

boundary_conditionstr

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

  • ‘UHTR’ :

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

  • ‘UBWT’ :

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

  • ‘MIFT’ :

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

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().

m_flow_borehole(nInlets,) array or (nMassFlow, nInlets,) array, optional

Fluid mass flow rate into each circuit of the network. If a (nMassFlow, nInlets,) array is supplied, the (nMassFlow, nMassFlow,) variable mass flow rate g-functions will be evaluated using the method of Cimmino (2024) [11]. Only required for the ‘MIFT’ boundary condition. Only one of ‘m_flow_borehole’ and ‘m_flow_network’ can be provided. Default is None.

m_flow_networkfloat or (nMassFlow,) array, optional

Fluid mass flow rate into the network of boreholes. If an array is supplied, the (nMassFlow, nMassFlow,) variable mass flow rate g-functions will be evaluated using the method of Cimmino (2024) [11]. Only required for the ‘MIFT’ boundary condition. Only one of ‘m_flow_borehole’ and ‘m_flow_network’ can be provided. Default is None.

cp_ffloat, optional

Fluid specific isobaric heat capacity (in J/kg.degC). Only required for the ‘MIFT’ boundary condition. Default is None.

approximate_FLSbool, optional

Set to true to use the approximation of the FLS solution of Cimmino (2021) [10]. 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.

kindstring, 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.

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.

References

find_groups(tol=1e-06)

Identify groups of boreholes that can be represented by a single equivalent borehole for the calculation of the g-function.

Hierarchical agglomerative clustering is applied to the superposed steady-state finite line source solution (i.e. the steady-state dimensionless borehole wall temperature due to a uniform heat extraction equal for all boreholes). The number of clusters is evaluated by cutting the dendrogram at the half-height of the longest branch and incrementing the number of intercepted branches by the value of the kClusters parameter.

Parameters:
tolfloat

Tolerance on the temperature to identify the maxiumum number of equivalent boreholes. Default is 1e-6.

Returns:
nSourcesint

Number of heat sources in the bore field.

initialize(disTol=0.01, tol=1e-06, kClusters=1, **kwargs)

Initialize paramteters. Identify groups for equivalent boreholes.

Returns:
nSourcesint

Number of finite line heat sources in the borefield used to initialize the matrix of segment-to-segment thermal response factors (of size: nSources x nSources).

segment_lengths()

Return the length of all segments in the bore field.

The segments lengths are used for the energy balance in the calculation of the g-function. For equivalent boreholes, the length of segments is multiplied by the number of boreholes in the group.

Returns:
Harray

Array of segment lengths (in m).

thermal_response_factors(time, alpha, kind='linear')

Evaluate the segment-to-segment thermal response factors for all pairs of segments in the borefield at all time steps using the finite line source solution.

This method returns a scipy.interpolate.interp1d object of the matrix of thermal response factors, containing a copy of the matrix accessible by h_ij.y[:nSources,:nSources,:nt+1]. The first index along the third axis corresponds to time t=0. The interp1d object can be used to obtain thermal response factors at any intermediate time by h_ij(t)[:nSources,:nSources].

Parameters:
timefloat or array

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

alphafloat

Soil thermal diffusivity (in m2/s).

kindstring, optional

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

Returns:
h_ijinterp1d

interp1d object (scipy.interpolate) of the matrix of segment-to-segment thermal response factors.

class pygfunction.solvers.Similarities(boreholes, network, time, boundary_condition, m_flow_borehole=None, m_flow_network=None, cp_f=None, nSegments=8, segment_ratios=<function segment_ratios>, approximate_FLS=False, mQuad=11, nFLS=10, linear_threshold=None, disp=False, profiles=False, kind='linear', dtype=<class 'numpy.float64'>, **other_options)

Bases: _BaseSolver

Similarities solver for the evaluation of the g-function.

This solver 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 [12]. The number of evaluations of the FLS solution is decreased by identifying similar pairs of boreholes, for which the same FLS value can be applied [14].

Parameters:
boreholeslist of Borehole objects

List of boreholes included in the bore field.

networknetwork object

Model of the network.

timefloat or array

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

boundary_conditionstr

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

  • ‘UHTR’ :

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

  • ‘UBWT’ :

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

  • ‘MIFT’ :

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

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().

m_flow_borehole(nInlets,) array or (nMassFlow, nInlets,) array, optional

Fluid mass flow rate into each circuit of the network. If a (nMassFlow, nInlets,) array is supplied, the (nMassFlow, nMassFlow,) variable mass flow rate g-functions will be evaluated using the method of Cimmino (2024) [17]. Only required for the ‘MIFT’ boundary condition. Only one of ‘m_flow_borehole’ and ‘m_flow_network’ can be provided. Default is None.

m_flow_networkfloat or (nMassFlow,) array, optional

Fluid mass flow rate into the network of boreholes. If an array is supplied, the (nMassFlow, nMassFlow,) variable mass flow rate g-functions will be evaluated using the method of Cimmino (2024) [17]. Only required for the ‘MIFT’ boundary condition. Only one of ‘m_flow_borehole’ and ‘m_flow_network’ can be provided. Default is None.

cp_ffloat, optional

Fluid specific isobaric heat capacity (in J/kg.degC). Only required for the ‘MIFT’ boundary condition. Default is None.

approximate_FLSbool, optional

Set to true to use the approximation of the FLS solution of Cimmino (2021) [16]. This approximation does not require the numerical evaluation of any integral. 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.

kindstring, 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.

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.

References

find_similarities()

Find similarities in the FLS solution for groups of boreholes.

This function identifies pairs of boreholes for which the evaluation of the Finite Line Source (FLS) solution is equivalent.

initialize(disTol=0.01, tol=1e-06, **kwargs)

Split boreholes into segments and identify similarities in the borefield.

Returns:
nSourcesint

Number of finite line heat sources in the borefield used to initialize the matrix of segment-to-segment thermal response factors (of size: nSources x nSources).

thermal_response_factors(time, alpha, kind='linear')

Evaluate the segment-to-segment thermal response factors for all pairs of segments in the borefield at all time steps using the finite line source solution.

This method returns a scipy.interpolate.interp1d object of the matrix of thermal response factors, containing a copy of the matrix accessible by h_ij.y[:nSources,:nSources,:nt+1]. The first index along the third axis corresponds to time t=0. The interp1d object can be used to obtain thermal response factors at any intermediate time by h_ij(t)[:nSources,:nSources].

Attributes:
timefloat or array

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

alphafloat

Soil thermal diffusivity (in m2/s).

kindstring, optional

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

Returns:
h_ijinterp1d

interp1d object (scipy.interpolate) of the matrix of segment-to-segment thermal response factors.