Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors

Models for heat transfer outside boreholes

Information

This package contains functions to evaluate temperature response factors used by Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.GroundTemperatureResponse to evaluate borehole wall temperatures.

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

Package Content

Name Description
Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.cylindricalHeatSource cylindricalHeatSource Cylindrical heat source solution from Carslaw and Jaeger
Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.cylindricalHeatSource_Integrand cylindricalHeatSource_Integrand Integrand function for cylindrical heat source evaluation
Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource finiteLineSource Finite line source solution of Claesson and Javed
Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Erfint finiteLineSource_Erfint Integral of the error function
Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Integrand finiteLineSource_Integrand Integrand function for finite line source evaluation
Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.gFunction gFunction Evaluate the g-function of a bore field
Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.infiniteLineSource infiniteLineSource Infinite line source model for borehole heat exchangers
Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.shaGFunction shaGFunction Returns a SHA1 encryption of the formatted arguments for the g-function generation
Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.timeGeometric timeGeometric Geometric expansion of time steps
Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.Validation Validation Example models for ThermalResponseFactors

Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.cylindricalHeatSource Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.cylindricalHeatSource

Cylindrical heat source solution from Carslaw and Jaeger

Information

This function evaluates the cylindrical heat source solution. This solution gives the relation between the constant heat transfer rate (per unit length) injected by a cylindrical heat source of infinite length and the temperature raise in the medium. The cylindrical heat source solution is defined by

image

where ΔT(t,r) is the temperature raise after a time t of constant heat injection and at a distance r from the cylindrical source, Q' is the heat injection rate per unit length, ks is the soil thermal conductivity, Fo is the Fourier number, aSois is the ground thermal diffusivity, rb is the radius of the cylindrical source and G is the cylindrical heat source solution.

The cylindrical heat source solution is given by:

image

The integral is solved numerically, with the integrand defined in Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.cylindricalHeatSource_Integrand.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Timet Time [s]
ThermalDiffusivityaSoi Ground thermal diffusivity [m2/s]
Distancedis Radial distance between borehole axes [m]
RadiusrBor Radius of emitting borehole [m]

Outputs

TypeNameDescription
RealGThermal response factor of borehole 1 on borehole 2

Modelica definition

function cylindricalHeatSource "Cylindrical heat source solution from Carslaw and Jaeger" extends Modelica.Icons.Function; input Modelica.Units.SI.Time t "Time"; input Modelica.Units.SI.ThermalDiffusivity aSoi "Ground thermal diffusivity"; input Modelica.Units.SI.Distance dis "Radial distance between borehole axes"; input Modelica.Units.SI.Radius rBor "Radius of emitting borehole"; output Real G "Thermal response factor of borehole 1 on borehole 2"; protected Real Fo = aSoi*t/rBor^2 "Fourier number"; Real p = dis/rBor "Fourier number"; algorithm G := Modelica.Math.Nonlinear.quadratureLobatto( function Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.cylindricalHeatSource_Integrand( Fo=Fo, p=p), a = 1e-12, b = 100, tolerance = 1e-6); end cylindricalHeatSource;

Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.cylindricalHeatSource_Integrand Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.cylindricalHeatSource_Integrand

Integrand function for cylindrical heat source evaluation

Information

Integrand of the cylindrical heat source solution for use in Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.cylindricalHeatSource.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realu Normalized integration variable
RealFo Fourier number
Realp Ratio of distance over radius

Outputs

TypeNameDescription
RealyValue of integrand

Modelica definition

function cylindricalHeatSource_Integrand "Integrand function for cylindrical heat source evaluation" extends Modelica.Icons.Function; input Real u "Normalized integration variable"; input Real Fo "Fourier number"; input Real p "Ratio of distance over radius"; output Real y "Value of integrand"; algorithm y := 1.0/(u^2*Modelica.Constants.pi^2)*(exp(-u^2*Fo) - 1.0) /(Buildings.Utilities.Math.Functions.besselJ1(u)^2+Buildings.Utilities.Math.Functions.besselY1(u)^2) *(Buildings.Utilities.Math.Functions.besselJ0(p*u)*Buildings.Utilities.Math.Functions.besselY1(u) -Buildings.Utilities.Math.Functions.besselJ1(u)*Buildings.Utilities.Math.Functions.besselY0(p*u)); end cylindricalHeatSource_Integrand;

Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource

Finite line source solution of Claesson and Javed

Information

This function evaluates the finite line source solution. This solution gives the relation between the constant heat transfer rate (per unit length) injected by a line source of finite length H1 buried at a distance D1 from a constant temperature surface (T=0) and the average temperature raise over a line of finite length H2 buried at a distance D2 from the constant temperature surface. The finite line source solution is defined by:

image

where ΔT1-2(t,r,H1,D1,H2,D2) is the temperature raise after a time t of constant heat injection and at a distance r from the line heat source, Q' is the heat injection rate per unit length, ks is the soil thermal conductivity and hFLS is the finite line source solution.

The finite line source solution is given by:

image

where αs is the ground thermal diffusivity and erfint is the integral of the error function, defined in Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_erfint. The integral is solved numerically, with the integrand defined in Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Integrand.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Timet Time [s]
ThermalDiffusivityaSoi Ground thermal diffusivity [m2/s]
Distancedis Radial distance between borehole axes [m]
Heightlen1 Length of emitting borehole [m]
HeightburDep1 Buried depth of emitting borehole [m]
Heightlen2 Length of receiving borehole [m]
HeightburDep2 Buried depth of receiving borehole [m]
BooleanincludeRealSourcetrueTrue if contribution of real source is included
BooleanincludeMirrorSourcetrueTrue if contribution of mirror source is included

Outputs

TypeNameDescription
Realh_21Thermal response factor of borehole 1 on borehole 2

Modelica definition

function finiteLineSource "Finite line source solution of Claesson and Javed" extends Modelica.Icons.Function; input Modelica.Units.SI.Time t "Time"; input Modelica.Units.SI.ThermalDiffusivity aSoi "Ground thermal diffusivity"; input Modelica.Units.SI.Distance dis "Radial distance between borehole axes"; input Modelica.Units.SI.Height len1 "Length of emitting borehole"; input Modelica.Units.SI.Height burDep1 "Buried depth of emitting borehole"; input Modelica.Units.SI.Height len2 "Length of receiving borehole"; input Modelica.Units.SI.Height burDep2 "Buried depth of receiving borehole"; input Boolean includeRealSource = true "True if contribution of real source is included"; input Boolean includeMirrorSource = true "True if contribution of mirror source is included"; output Real h_21 "Thermal response factor of borehole 1 on borehole 2"; protected Real lowBou(unit="m-1") "Lower bound of integration"; // Upper bound is infinite Real uppBou(unit="m-1") = max(100.0, 10.0/dis) "Upper bound of integration"; Modelica.Units.SI.Distance disMin "Minimum distance between sources and receiving line"; Modelica.Units.SI.Time timTre "Time treshold for evaluation of the solution"; algorithm h_21 := 0; if t > 0 and (includeRealSource or includeMirrorSource) then // Find the minimum distance between the line source and the line where the // finite line source solution is evaluated. if includeRealSource then if (burDep2 + len2) < burDep1 then disMin := sqrt(dis^2 + (burDep1 - burDep2 - len2)^2); elseif burDep2 > (burDep1 + len1) then disMin := sqrt(dis^2 + (burDep1 - burDep2 + len1)^2); else disMin := dis; end if; else disMin := sqrt(dis^2 + (burDep1 + burDep2)^2); end if; // The traveled distance of the temperature front is assumed to be: // d = 5*sqrt(aSoi*t). // The solution is only evaluated at times when the traveled distance is // greater than the minimum distance. timTre := disMin^2/(25*aSoi); if t >= timTre then lowBou := 1.0/sqrt(4*aSoi*t); h_21 := Modelica.Math.Nonlinear.quadratureLobatto( function Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Integrand( dis=dis, len1=len1, burDep1=burDep1, len2=len2, burDep2=burDep2, includeRealSource=includeRealSource, includeMirrorSource=includeMirrorSource), lowBou, uppBou, 1.0e-6); else // Linearize the solution at times below the time treshold. lowBou := 1.0/sqrt(4*aSoi*timTre); h_21 := t/timTre*Modelica.Math.Nonlinear.quadratureLobatto( function Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Integrand( dis=dis, len1=len1, burDep1=burDep1, len2=len2, burDep2=burDep2, includeRealSource=includeRealSource, includeMirrorSource=includeMirrorSource), lowBou, uppBou, 1.0e-6); end if; end if; end finiteLineSource;

Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Erfint Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Erfint

Integral of the error function

Information

This function evaluates the integral of the error function, given by:

image

Extends from Modelica.Math.Nonlinear.Interfaces.partialScalarFunction (Interface for a function with one input and one output Real signal).

Inputs

TypeNameDefaultDescription
Realu Independent variable

Outputs

TypeNameDescription
RealyDependent variable y=f(u)

Modelica definition

function finiteLineSource_Erfint "Integral of the error function" extends Modelica.Math.Nonlinear.Interfaces.partialScalarFunction; algorithm y := u*Modelica.Math.Special.erf(u) - 1/sqrt(Modelica.Constants.pi)*(1 - exp(-u^2)); end finiteLineSource_Erfint;

Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Integrand Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Integrand

Integrand function for finite line source evaluation

Information

Integrand of the cylindrical heat source solution for use in Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realu Integration variable [1/m]
Distancedis Radial distance between borehole axes [m]
Heightlen1 Length of emitting borehole [m]
HeightburDep1 Buried depth of emitting borehole [m]
Heightlen2 Length of receiving borehole [m]
HeightburDep2 Buried depth of receiving borehole [m]
BooleanincludeRealSourcetruetrue if contribution of real source is included
BooleanincludeMirrorSourcetruetrue if contribution of mirror source is included

Outputs

TypeNameDescription
RealyValue of integrand [m]

Modelica definition

function finiteLineSource_Integrand "Integrand function for finite line source evaluation" extends Modelica.Icons.Function; input Real u(unit="1/m") "Integration variable"; input Modelica.Units.SI.Distance dis "Radial distance between borehole axes"; input Modelica.Units.SI.Height len1 "Length of emitting borehole"; input Modelica.Units.SI.Height burDep1 "Buried depth of emitting borehole"; input Modelica.Units.SI.Height len2 "Length of receiving borehole"; input Modelica.Units.SI.Height burDep2 "Buried depth of receiving borehole"; input Boolean includeRealSource = true "true if contribution of real source is included"; input Boolean includeMirrorSource = true "true if contribution of mirror source is included"; output Real y(unit="m") "Value of integrand"; protected Real f "Intermediate variable"; algorithm if includeRealSource then f := sum({ +Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Erfint( (burDep2 - burDep1 + len2)*u), -Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Erfint( (burDep2 - burDep1)*u), +Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Erfint( (burDep2 - burDep1 - len1)*u), -Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Erfint( (burDep2 - burDep1 + len2 - len1)*u)}); else f := 0; end if; if includeMirrorSource then f := f + sum({ +Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Erfint( (burDep2 + burDep1 + len2)*u), -Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Erfint( (burDep2 + burDep1)*u), +Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Erfint( (burDep2 + burDep1 + len1)*u), -Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource_Erfint( (burDep2 + burDep1 + len2 + len1)*u)}); end if; y := 0.5/(len2*u^2)*f*exp(-dis^2*u^2); end finiteLineSource_Integrand;

Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.gFunction Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.gFunction

Evaluate the g-function of a bore field

Information

This function implements the g-function evaluation method introduced by Cimmino and Bernier (see: Cimmino and Bernier (2014), and Cimmino (2018)) based on the g-function function concept first introduced by Eskilson (1987). The g-function gives the relation between the variation of the borehole wall temperature at a time t and the heat extraction and injection rates at all times preceding time t as

image

where Tb is the borehole wall temperature, Tg is the undisturbed ground temperature, Q is the heat injection rate into the ground through the borehole wall per unit borehole length, ks is the soil thermal conductivity and g is the g-function.

The g-function is constructed from the combination of the combination of the finite line source (FLS) solution (see Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource), the cylindrical heat source (CHS) solution (see Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.cylindricalHeatSource), and the infinite line source (ILS) solution (see Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.infiniteLineSource). To obtain the g-function of a bore field, each borehole is divided into a series of nSeg segments of equal length, each modeled as a line source of finite length. The finite line source solution is superimposed in space to obtain a system of equations that gives the relation between the heat injection rate at each of the segments and the borehole wall temperature at each of the segments. The system is solved to obtain the uniform borehole wall temperature required at any time to maintain a constant total heat injection rate (Qtot = 2πksHtot) into the bore field. The uniform borehole wall temperature is then equal to the finite line source based g-function.

Since this g-function is based on line sources of heat, rather than cylinders, the g-function is corrected to consider the cylindrical geometry. The correction factor is then the difference between the cylindrical heat source solution and the infinite line source solution, as proposed by Li et al. (2014) as

g(t) = gFLS + (gCHS - gILS)

Implementation

The calculation of the g-function is separated into two regions: the short-time region and the long-time region. In the short-time region, corresponding to times t < 1 hour, heat interaction between boreholes and axial variations of heat injection rate are not considered. The g-function is calculated using only one borehole and one segment. In the long-time region, corresponding to times t > 1 hour, all boreholes are represented as series of nSeg line segments and the g-function is evaluated as described above.

References

Cimmino, M. and 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.

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. DOI: 10.1080/19401493.2017.1423390.

Eskilson, P. 1987. Thermal analysis of heat extraction boreholes. Ph.D. Thesis. Department of Mathematical Physics. University of Lund. Sweden.

Li, M., Li, P., Chan, V. and Lai, A.C.K. 2014. Full-scale temperature response function (G-function) for heat transfer by borehole heat exchangers (GHEs) from sub-hour to decades. Applied Energy 136: 197-205.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
IntegernBor Number of boreholes
PositioncooBor[nBor, 2] Coordinates of boreholes [m]
HeighthBor Borehole length [m]
HeightdBor Borehole buried depth [m]
RadiusrBor Borehole radius [m]
ThermalDiffusivityaSoi Ground thermal diffusivity used in g-function evaluation [m2/s]
IntegernSeg Number of line source segments per borehole
IntegernTimSho Number of time steps in short time region
IntegernTimLon Number of time steps in long time region
RealttsMax Maximum adimensional time for gfunc calculation
RealrelTol0.02Relative tolerance on distance between boreholes

Outputs

TypeNameDescription
TimetGFun[nTimSho + nTimLon]Time of g-function evaluation [s]
Realg[nTimSho + nTimLon]g-function

Modelica definition

function gFunction "Evaluate the g-function of a bore field" extends Modelica.Icons.Function; input Integer nBor "Number of boreholes"; input Modelica.Units.SI.Position cooBor[nBor,2] "Coordinates of boreholes"; input Modelica.Units.SI.Height hBor "Borehole length"; input Modelica.Units.SI.Height dBor "Borehole buried depth"; input Modelica.Units.SI.Radius rBor "Borehole radius"; input Modelica.Units.SI.ThermalDiffusivity aSoi "Ground thermal diffusivity used in g-function evaluation"; input Integer nSeg "Number of line source segments per borehole"; input Integer nTimSho "Number of time steps in short time region"; input Integer nTimLon "Number of time steps in long time region"; input Real ttsMax "Maximum adimensional time for gfunc calculation"; input Real relTol = 0.02 "Relative tolerance on distance between boreholes"; output Modelica.Units.SI.Time tGFun[nTimSho + nTimLon] "Time of g-function evaluation"; output Real g[nTimSho+nTimLon] "g-function"; protected Modelica.Units.SI.Time ts=hBor^2/(9*aSoi) "Characteristic time"; Modelica.Units.SI.Time tSho_min=1 "Minimum time for short time calculations"; Modelica.Units.SI.Time tSho_max=3600 "Maximum time for short time calculations"; Modelica.Units.SI.Time tLon_min=tSho_max "Minimum time for long time calculations"; Modelica.Units.SI.Time tLon_max=ts*ttsMax "Maximum time for long time calculations"; Modelica.Units.SI.Time tSho[nTimSho] "Time vector for short time calculations"; Modelica.Units.SI.Time tLon[nTimLon] "Time vector for long time calculations"; Modelica.Units.SI.Distance dis "Separation distance between boreholes"; Modelica.Units.SI.Distance dis_mn "Separation distance for comparison"; Modelica.Units.SI.Radius rLin=0.0005*hBor "Radius for evaluation of same-borehole line source solutions"; Real hSegRea[nSeg] "Real part of the FLS solution"; Real hSegMir[2*nSeg-1] "Mirror part of the FLS solution"; Modelica.Units.SI.Height dSeg "Buried depth of borehole segment"; Integer Done[nBor, nBor] = fill(0,nBor,nBor) "Matrix for tracking of FLS evaluations"; Real A[nSeg*nBor+1, nSeg*nBor+1] "Coefficient matrix for system of equations"; Real B[nSeg*nBor+1] "Coefficient vector for system of equations"; Real X[nSeg*nBor+1] "Solution vector for system of equations"; Real FLS "Finite line source solution"; Real ILS "Infinite line source solution"; Real CHS "Cylindrical heat source solution"; algorithm // Generate geometrically expanding time vectors tSho := Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.timeGeometric( tSho_min, tSho_max, nTimSho); tLon := Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.timeGeometric( tLon_min, tLon_max, nTimLon); // Concatenate the short- and long-term parts tGFun := cat(1, {0}, tSho[1:nTimSho - 1], tLon); // ----------------------- // Short time calculations // ----------------------- g[1] := 0.; for k in 1:nTimSho loop // Finite line source solution FLS := Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource( tSho[k], aSoi, rLin, hBor, dBor, hBor, dBor); // Infinite line source solution ILS := 0.5* Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.infiniteLineSource( tSho[k], aSoi, rLin); // Cylindrical heat source solution CHS := 2*Modelica.Constants.pi* Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.cylindricalHeatSource( tSho[k], aSoi, rBor, rBor); // Correct finite line source solution for cylindrical geometry g[k+1] := FLS + (CHS - ILS); end for; // ---------------------- // Long time calculations // ---------------------- // Initialize coefficient matrix A for m in 1:nBor loop for u in 1:nSeg loop // Tb coefficient in spatial superposition equations A[(m-1)*nSeg+u,nBor*nSeg+1] := -1; // Q coefficient in heat balance equation A[nBor*nSeg+1,(m-1)*nSeg+u] := 1; end for; end for; // Initialize coefficient vector B // The total heat extraction rate is constant B[nBor*nSeg+1] := nBor*nSeg; // Evaluate thermal response matrix at all times for k in 1:nTimLon-1 loop for i in 1:nBor loop for j in i:nBor loop // Distance between boreholes if i == j then // If same borehole, distance is the radius dis := rLin; else dis := sqrt((cooBor[i,1] - cooBor[j,1])^2 + (cooBor[i,2] - cooBor[j,2])^2); end if; // Only evaluate the thermal response factors if not already evaluated if Done[i,j] < k then // Evaluate Real and Mirror parts of FLS solution // Real part for m in 1:nSeg loop hSegRea[m] := Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource( tLon[k + 1], aSoi, dis, hBor/nSeg, dBor, hBor/nSeg, dBor + (m - 1)*hBor/nSeg, includeMirrorSource=false); end for; // Mirror part for m in 1:(2*nSeg-1) loop hSegMir[m] := Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.finiteLineSource( tLon[k + 1], aSoi, dis, hBor/nSeg, dBor, hBor/nSeg, dBor + (m - 1)*hBor/nSeg, includeRealSource=false); end for; // Apply to all pairs that have the same separation distance for m in 1:nBor loop for n in m:nBor loop if m == n then dis_mn := rLin; else dis_mn := sqrt((cooBor[m,1] - cooBor[n,1])^2 + (cooBor[m,2] - cooBor[n,2])^2); end if; if abs(dis_mn - dis) < relTol*dis then // Add thermal response factor to coefficient matrix A for u in 1:nSeg loop for v in 1:nSeg loop A[(m-1)*nSeg+u,(n-1)*nSeg+v] := hSegRea[abs(u-v)+1] + hSegMir[u+v-1]; A[(n-1)*nSeg+v,(m-1)*nSeg+u] := hSegRea[abs(u-v)+1] + hSegMir[u+v-1]; end for; end for; // Mark current pair as evaluated Done[m,n] := k; Done[n,m] := k; end if; end for; end for; end if; end for; end for; // Solve the system of equations X := Modelica.Math.Matrices.solve(A,B); // The g-function is equal to the borehole wall temperature g[nTimSho+k+1] := X[nBor*nSeg+1]; end for; // Correct finite line source solution for cylindrical geometry for k in 2:nTimLon loop // Infinite line source ILS := 0.5* Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.infiniteLineSource( tLon[k], aSoi, rLin); // Cylindrical heat source CHS := 2*Modelica.Constants.pi* Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.cylindricalHeatSource( tLon[k], aSoi, rBor, rBor); g[nTimSho+k] := g[nTimSho+k] + (CHS - ILS); end for; end gFunction;

Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.infiniteLineSource Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.infiniteLineSource

Infinite line source model for borehole heat exchangers

Information

This function evaluates the infinite line source solution. This solution gives the relation between the constant heat transfer rate (per unit length) injected by a line heat source of infinite length and the temperature raise in the medium. The infinite line source solution is defined by

image

where ΔT(t,r) is the temperature raise after a time t of constant heat injection and at a distance r from the line source, Q' is the heat injection rate per unit length, ks is the soil thermal conductivity and hILS is the infinite line source solution.

The infinite line source solution is given by the exponential integral

image

where αs is the ground thermal diffusivity. The exponential integral is implemented in Buildings.Utilities.Math.Functions.exponentialIntegralE1.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Realt Time
RealaSoi Ground thermal diffusivity
Realdis Radial distance between borehole axes

Outputs

TypeNameDescription
Realh_ilsThermal response factor of borehole 1 on borehole 2

Modelica definition

function infiniteLineSource "Infinite line source model for borehole heat exchangers" extends Modelica.Icons.Function; input Real t "Time"; input Real aSoi "Ground thermal diffusivity"; input Real dis "Radial distance between borehole axes"; output Real h_ils "Thermal response factor of borehole 1 on borehole 2"; algorithm h_ils := if t > 0.0 then Buildings.Utilities.Math.Functions.exponentialIntegralE1(dis^2/(4*aSoi*t)) else 0.0; end infiniteLineSource;

Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.shaGFunction Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.shaGFunction

Returns a SHA1 encryption of the formatted arguments for the g-function generation

Information

This function returns the SHA1 encryption of its arguments.

Implementation

Each argument is formatted in exponential notation with four significant digits, for example 1.234e+001, with no spaces or other separating characters between each argument value. To prevent too long strings that can cause buffer overflows, the sha encoding of each argument is computed and added to the next string that is parsed.

The SHA1 encryption is computed using Buildings.Utilities.Cryptographics.sha.

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
IntegernBor Number of boreholes
PositioncooBor[nBor, 2] Coordinates of boreholes [m]
HeighthBor Borehole length [m]
HeightdBor Borehole buried depth [m]
RadiusrBor Borehole radius [m]
ThermalDiffusivityaSoi Ground thermal diffusivity used in g-function evaluation [m2/s]
IntegernSeg Number of line source segments per borehole
IntegernTimSho Number of time steps in short time region
IntegernTimLon Number of time steps in long time region
RealttsMax Maximum adimensional time for gfunc calculation

Outputs

TypeNameDescription
StringshaSHA1 encryption of the g-function arguments

Modelica definition

function shaGFunction "Returns a SHA1 encryption of the formatted arguments for the g-function generation" extends Modelica.Icons.Function; input Integer nBor "Number of boreholes"; input Modelica.Units.SI.Position cooBor[nBor,2] "Coordinates of boreholes"; input Modelica.Units.SI.Height hBor "Borehole length"; input Modelica.Units.SI.Height dBor "Borehole buried depth"; input Modelica.Units.SI.Radius rBor "Borehole radius"; input Modelica.Units.SI.ThermalDiffusivity aSoi "Ground thermal diffusivity used in g-function evaluation"; input Integer nSeg "Number of line source segments per borehole"; input Integer nTimSho "Number of time steps in short time region"; input Integer nTimLon "Number of time steps in long time region"; input Real ttsMax "Maximum adimensional time for gfunc calculation"; output String sha "SHA1 encryption of the g-function arguments"; protected constant String formatStrGen = "1.3e" "String format for general parameters"; constant String formatStrCoo = ".2f" "String format for coordinate"; algorithm sha := Buildings.Utilities.Cryptographics.sha(String(nBor, format=formatStrGen)); sha := Buildings.Utilities.Cryptographics.sha(sha + String(hBor, format=formatStrGen)); sha := Buildings.Utilities.Cryptographics.sha(sha + String(dBor, format=formatStrGen)); sha := Buildings.Utilities.Cryptographics.sha(sha + String(rBor, format=formatStrGen)); sha := Buildings.Utilities.Cryptographics.sha(sha + String(aSoi, format=formatStrGen)); sha := Buildings.Utilities.Cryptographics.sha(sha + String(nSeg, format=formatStrGen)); sha := Buildings.Utilities.Cryptographics.sha(sha + String(nTimSho, format=formatStrGen)); sha := Buildings.Utilities.Cryptographics.sha(sha + String(nTimLon, format=formatStrGen)); sha := Buildings.Utilities.Cryptographics.sha(sha + String(ttsMax, format=formatStrGen)); for i in 1:nBor loop sha := Buildings.Utilities.Cryptographics.sha(sha + String(cooBor[i, 1], format=formatStrCoo)); sha := Buildings.Utilities.Cryptographics.sha(sha + String(cooBor[i, 2], format=formatStrCoo)); end for; end shaGFunction;

Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.timeGeometric Buildings.Fluid.Geothermal.Borefields.BaseClasses.HeatTransfer.ThermalResponseFactors.timeGeometric

Geometric expansion of time steps

Information

This function attemps to build a vector of length nTim with a geometric expansion of the time variable between dt and t_max.

If t_max > nTim*dt, then a geometrically expanding vector is built as

t = [dt, dt*(1-r2)/(1-r), ... , dt*(1-rn)/(1-r), ... , tmax],

where r is the geometric expansion factor.

If t_max < nTim*dt, then a linearly expanding vector is built as

t = [dt, 2*dt, ... , n*dt, ... , nTim*dt]

Extends from Modelica.Icons.Function (Icon for functions).

Inputs

TypeNameDefaultDescription
Durationdt Minimum time step [s]
Timet_max Maximum value of time [s]
IntegernTim Number of time values

Outputs

TypeNameDescription
Realt[nTim]Time vector

Modelica definition

function timeGeometric "Geometric expansion of time steps" extends Modelica.Icons.Function; input Modelica.Units.SI.Duration dt "Minimum time step"; input Modelica.Units.SI.Time t_max "Maximum value of time"; input Integer nTim "Number of time values"; output Real t[nTim] "Time vector"; protected Real r(min=1.) "Expansion rate of time values"; Real dr "Error on expansion rate evaluation"; algorithm if t_max > nTim*dt then // Determine expansion rate (r) dr := 1e99; r := 2; while abs(dr) > 1e-10 loop dr := (1+t_max/dt*(r-1))^(1/nTim) - r; r := r + dr; end while; // Assign time values for i in 1:nTim-1 loop t[i] := dt*(1-r^i)/(1-r); end for; t[nTim] := t_max; else // Number of time values too large for chosen parameters: // Use a constant time step for i in 1:nTim loop t[i] := i*dt; end for; end if; end timeGeometric;