Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions

Package with functions for evaluation of borehole thermal resistances

Information

This package contains functions to evaluate borehole internal resistances used by models in Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.

Extends from Modelica.Icons.VariantsPackage (Icon for package containing variants).

Package Content

Name Description
Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.convectionResistanceCircularPipe convectionResistanceCircularPipe Thermal resistance from the fluid in pipes and the grout zones (Bauer et al. 2011)
Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.internalResistancesOneUTube internalResistancesOneUTube Thermal resistances for single U-tube, according to Bauer et al. (2011)
Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.internalResistancesTwoUTube internalResistancesTwoUTube Thermal resistances for double U-tube, according to Bauer et al (2011)
Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.multipoleFluidTemperature multipoleFluidTemperature Fluid temperatures from multipole solution
Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.multipoleFmk multipoleFmk Complex matrix F_mk from Claesson and Hellstrom (2011)
Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.multipoleThermalResistances multipoleThermalResistances Thermal resistances from multipole solution
Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.partialInternalResistances partialInternalResistances Partial model for borehole resistance calculation
Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.Validation Validation Models to validate borehole thermal resistances functions

Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.convectionResistanceCircularPipe Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.convectionResistanceCircularPipe

Thermal resistance from the fluid in pipes and the grout zones (Bauer et al. 2011)

Information

This model computes the convection resistance in the pipes of a borehole segment with heigth hSeg using correlations suggested by Bergman et al. (2011).

If the flow is laminar (Re ≤ 2300, with Re being the Reynolds number of the flow), the Nusselt number of the flow is assumed to be constant at 3.66. If the flow is turbulent (Re > 2300), the correlation of Dittus-Boelter is used to find the convection heat transfer coefficient as

Nu = 0.023   Re0.8   Prn,

where Nu is the Nusselt number and Pr is the Prandlt number. A value of n=0.35 is used, as the reference uses n=0.4 for heating and n=0.3 for cooling. To ensure that the function is continuously differentiable, a smooth transition between the laminar and turbulent values is created for the range 2300 < Re < 2400.

References

Bergman, T. L., Incropera, F. P., DeWitt, D. P., & Lavine, A. S. (2011). Fundamentals of heat and mass transfer (7th ed.). New York: John Wiley & Sons.

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

Inputs

TypeNameDefaultDescription
HeighthSeg Height of the element [m]
RadiusrTub Tube radius [m]
LengtheTub Tube thickness [m]
ThermalConductivitykMed Thermal conductivity of the fluid [W/(m.K)]
DynamicViscositymuMed Dynamic viscosity of the fluid [Pa.s]
SpecificHeatCapacitycpMed Specific heat capacity of the fluid [J/(kg.K)]
MassFlowRatem_flow Mass flow rate [kg/s]
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]

Outputs

TypeNameDescription
ThermalResistanceRFluPipConvection resistance (or conduction in fluid if no mass flow) [K/W]

Modelica definition

function convectionResistanceCircularPipe "Thermal resistance from the fluid in pipes and the grout zones (Bauer et al. 2011)" extends Modelica.Icons.Function; // Geometry of the borehole input Modelica.Units.SI.Height hSeg "Height of the element"; input Modelica.Units.SI.Radius rTub "Tube radius"; input Modelica.Units.SI.Length eTub "Tube thickness"; // thermal properties input Modelica.Units.SI.ThermalConductivity kMed "Thermal conductivity of the fluid"; input Modelica.Units.SI.DynamicViscosity muMed "Dynamic viscosity of the fluid"; input Modelica.Units.SI.SpecificHeatCapacity cpMed "Specific heat capacity of the fluid"; input Modelica.Units.SI.MassFlowRate m_flow "Mass flow rate"; input Modelica.Units.SI.MassFlowRate m_flow_nominal "Nominal mass flow rate"; // Outputs output Modelica.Units.SI.ThermalResistance RFluPip "Convection resistance (or conduction in fluid if no mass flow)"; protected parameter Modelica.Units.SI.Radius rTub_in=rTub - eTub "Pipe inner radius"; Modelica.Units.SI.CoefficientOfHeatTransfer h "Convective heat transfer coefficient of the fluid"; Real k(unit="s/kg") "Coefficient used in the computation of the convective heat transfer coefficient"; Modelica.Units.SI.MassFlowRate m_flow_abs= Buildings.Utilities.Math.Functions.spliceFunction( m_flow, -m_flow, m_flow, m_flow_nominal/30); Real Re "Reynolds number"; Real NuTurb "Nusselt at Re=2400"; Real Nu "Nusselt"; algorithm // Convection resistance and Reynolds number k := 2/(muMed*Modelica.Constants.pi*rTub_in); Re := m_flow_abs*k; if Re>=2400 then // Turbulent, fully-developped flow in a smooth circular pipe with the // Dittus-Boelter correlation: h = 0.023*k_f*Re*Pr/(2*rTub) // Re = rho*v*DTub / mue_f // = m_flow/(pi r^2) * DTub/mue_f = 2*m_flow / ( mue*pi*rTub) Nu := 0.023*(cpMed*muMed/kMed)^(0.35)* Buildings.Utilities.Math.Functions.regNonZeroPower( x=Re, n=0.8, delta=0.01*m_flow_nominal*k); else // Laminar, fully-developped flow in a smooth circular pipe with uniform // imposed temperature: Nu=3.66 for Re<=2300. For 2300<Re<2400, a smooth // transition is created with the splice function. NuTurb := 0.023*(cpMed*muMed/kMed)^(0.35)*(2400)^(0.8); Nu := Buildings.Utilities.Math.Functions.spliceFunction(NuTurb,3.66,Re-(2300+2400)/2,((2300+2400)/2)-2300); end if; h := Nu*kMed/(2*rTub_in); RFluPip := 1/(2*Modelica.Constants.pi*rTub_in*hSeg*h); end convectionResistanceCircularPipe;

Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.internalResistancesOneUTube Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.internalResistancesOneUTube

Thermal resistances for single U-tube, according to Bauer et al. (2011)

Information

This model computes the different thermal resistances present in a single-U-tube borehole using the method of Bauer et al. (2011). It also computes the fluid-to-ground thermal resistance Rb and the grout-to-grout thermal resistance Ra as defined by Claesson and Hellstrom (2011) using the multipole method.

References

J. Claesson and G. Hellstrom. Multipole method to calculate borehole thermal resistances in a borehole heat exchanger. HVAC&R Research, 17(6): 895-911, 2011.

D. Bauer, W. Heidemann, H. Müller-Steinhagen, and H.-J. G. Diersch. Thermal resistance and capacity models for borehole heat exchangers. International Journal of Energy Research, 35:312–320, 2011.

Extends from Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.partialInternalResistances (Partial model for borehole resistance calculation).

Inputs

TypeNameDefaultDescription
Booleanuse_RbfalseTrue if the value Rb should be used instead of calculated
RealRb Borehole thermal resistance [(m.K)/W]
HeighthSeg Height of the element [m]
RadiusrBor Radius of the borehole [m]
RadiusrTub Radius of the tube [m]
LengtheTub Thickness of the tubes [m]
Lengthsha Shank spacing, defined as the distance between the center of a pipe and the center of the borehole [m]
ThermalConductivitykFil Thermal conductivity of the grout [W/(m.K)]
ThermalConductivitykSoi Thermal conductivity of the soi [W/(m.K)]
ThermalConductivitykTub Thermal conductivity of the tube [W/(m.K)]
ThermalConductivitykMed Thermal conductivity of the fluid [W/(m.K)]
DynamicViscositymuMed Dynamic viscosity of the fluid [Pa.s]
SpecificHeatCapacitycpMed Specific heat capacity of the fluid [J/(kg.K)]
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]

Outputs

TypeNameDescription
RealxCapacity location
ThermalResistanceRgbThermal resistance between grout zone and borehole wall [K/W]
ThermalResistanceRggThermal resistance between the two grout zones [K/W]
ThermalResistanceRCondGroThermal resistance between: pipe wall to capacity in grout [K/W]

Modelica definition

function internalResistancesOneUTube "Thermal resistances for single U-tube, according to Bauer et al. (2011)" extends Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.partialInternalResistances; // Outputs output Modelica.Units.SI.ThermalResistance Rgb "Thermal resistance between grout zone and borehole wall"; output Modelica.Units.SI.ThermalResistance Rgg "Thermal resistance between the two grout zones"; output Modelica.Units.SI.ThermalResistance RCondGro "Thermal resistance between: pipe wall to capacity in grout"; protected Real[2,2] RDelta(each unit="(m.K)/W") "Delta-circuit thermal resistances"; Real[2,2] R(each unit="(m.K)/W") "Internal thermal resistances"; Modelica.Units.SI.Position[2] xPip={-sha,sha} "x-Coordinates of pipes"; Modelica.Units.SI.Position[2] yPip={0.,0.} "y-Coordinates of pipes"; Modelica.Units.SI.Radius[2] rPip={rTub,rTub} "Outer radius of pipes"; Real[2] RFluPip(each unit="(m.K)/W") = {RCondPipe+RConv, RCondPipe+RConv} "Fluid to pipe wall thermal resistances"; Modelica.Units.SI.ThermalResistance Rg "Thermal resistance between outer borehole wall and one tube"; Modelica.Units.SI.ThermalResistance Rar "Thermal resistance between the two pipe outer walls"; Real Ra(unit="(m.K)/W") "Grout-to-grout resistance (2D) as defined by Hellstrom. Interaction between the different grout parts"; algorithm // Internal thermal resistances (RDelta, R) := Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.multipoleThermalResistances( 2, 3, xPip, yPip, rBor, rPip, kFil, kSoi, RFluPip); // Rb and Ra Rb_multipole := 1./(1./RDelta[1,1] + 1./RDelta[2,2]); Rb_internal := if use_Rb then Rb else Rb_multipole; // The short-circuit resistance in weigthed by the ratio between the used // value of Rb and the theoretical value Ra := (R[1,1] + R[2,2] - 2*R[1,2])*Rb_internal/Rb_multipole; // Conversion of Rb (resp. Ra) to Rg (resp. Rar) of Bauer: Rg :=(2*Rb_internal-RCondPipe-RConv)/hSeg; Rar :=(Ra-2*(RCondPipe + RConv))/hSeg; // If any of the internal delta-circuit resistances is negative, then // the location (x) of the thermal capacity is set to zero to limit // instabilities in the calculations. Otherwise, calculations follow the // method of Bauer et al. (2011). if (RDelta[1,2] < 0) then //Thermal resistance between the grout zone and borehole wall Rgb := Rg; //Conduction resistance in grout from pipe wall to capacity in grout RCondGro := RCondPipe/hSeg; //Thermal resistance between the two grout zones Rgg := 2*Rgb*Rar/(2*Rgb - Rar); test := true; else // ********** Resistances and capacity location according to Bauer ********** while test == false and i <= 10 loop // Capacity location (with correction factor in case that the test is // negative) x := Modelica.Math.log(sqrt(rBor^2 + 2*rTub^2)/(2*rTub))/ Modelica.Math.log(rBor/(sqrt(2)*rTub))*((15 - i + 1)/15); //Thermal resistance between the grout zone and borehole wall Rgb := (1 - x)*Rg; //Conduction resistance in grout from pipe wall to capacity in grout RCondGro := x*Rg + RCondPipe/hSeg; //Thermal resistance between the two grout zones Rgg := 2*Rgb*(Rar - 2*x*Rg)/(2*Rgb - Rar + 2*x*Rg); // Thermodynamic test to check if negative R values make sense. If not, // decrease x-value. test := ((1/Rgg + 1/2/Rgb)^(-1) > 0); i := i + 1; end while; end if; assert(test, "In " + getInstanceName() + ":\n" + "Maximum number of iterations exceeded. Check the borehole geometry. The tubes may be too close to the borehole wall. Input to the function Buildings.Fluid.HeatExchangers.Boreholes.BaseClasses.singleUTubeResistances is hSeg = " + String(hSeg) + " m rBor = " + String(rBor) + " m rTub = " + String(rTub) + " m eTub = " + String(eTub) + " m kSoi = " + String(kSoi) + " W/m/K kFil = " + String(kFil) + " W/m/K kTub = " + String(kTub) + " W/m/K i = " + String(i) + " Computed x = " + String(x) + " K/W Rgb = " + String(Rgb) + " K/W Rgg = " + String(Rgg) + " K/W"); end internalResistancesOneUTube;

Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.internalResistancesTwoUTube Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.internalResistancesTwoUTube

Thermal resistances for double U-tube, according to Bauer et al (2011)

Information

This model computes the different thermal resistances present in a double U-tube borehole using the method of Bauer et al. (2011). It also computes the fluid-to-ground thermal resistance Rb and the grout-to-grout thermal resistance Ra as defined by Claesson and Hellstrom (2011) using the multipole method.

References

J. Claesson and G. Hellstrom. Multipole method to calculate borehole thermal resistances in a borehole heat exchanger. HVAC&R Research, 17(6): 895-911, 2011.

D. Bauer, W. Heidemann, H. Müller-Steinhagen, and H.-J. G. Diersch. Thermal resistance and capacity models for borehole heat exchangers. International Journal of Energy Research, 35:312–320, 2011.

Extends from Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.partialInternalResistances (Partial model for borehole resistance calculation).

Inputs

TypeNameDefaultDescription
Booleanuse_RbfalseTrue if the value Rb should be used instead of calculated
RealRb Borehole thermal resistance [(m.K)/W]
HeighthSeg Height of the element [m]
RadiusrBor Radius of the borehole [m]
RadiusrTub Radius of the tube [m]
LengtheTub Thickness of the tubes [m]
Lengthsha Shank spacing, defined as the distance between the center of a pipe and the center of the borehole [m]
ThermalConductivitykFil Thermal conductivity of the grout [W/(m.K)]
ThermalConductivitykSoi Thermal conductivity of the soi [W/(m.K)]
ThermalConductivitykTub Thermal conductivity of the tube [W/(m.K)]
ThermalConductivitykMed Thermal conductivity of the fluid [W/(m.K)]
DynamicViscositymuMed Dynamic viscosity of the fluid [Pa.s]
SpecificHeatCapacitycpMed Specific heat capacity of the fluid [J/(kg.K)]
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]

Outputs

TypeNameDescription
RealxCapacity location
ThermalResistanceRgbThermal resistance between a grout capacity and the borehole wall, as defined by Bauer et al (2010) [K/W]
ThermalResistanceRgg1Thermal resistance between two neightbouring grout capacities, as defined by Bauer et al (2010) [K/W]
ThermalResistanceRgg2Thermal resistance between two grout capacities opposite to each other, as defined by Bauer et al (2010) [K/W]
ThermalResistanceRCondGroThermal resistance between a pipe wall and the grout capacity, as defined by Bauer et al (2010) [K/W]

Modelica definition

function internalResistancesTwoUTube "Thermal resistances for double U-tube, according to Bauer et al (2011)" extends Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.partialInternalResistances; // Outputs output Modelica.Units.SI.ThermalResistance Rgb "Thermal resistance between a grout capacity and the borehole wall, as defined by Bauer et al (2010)"; output Modelica.Units.SI.ThermalResistance Rgg1 "Thermal resistance between two neightbouring grout capacities, as defined by Bauer et al (2010)"; output Modelica.Units.SI.ThermalResistance Rgg2 "Thermal resistance between two grout capacities opposite to each other, as defined by Bauer et al (2010)"; output Modelica.Units.SI.ThermalResistance RCondGro "Thermal resistance between a pipe wall and the grout capacity, as defined by Bauer et al (2010)"; protected Real[4,4] RDelta(each unit="(m.K)/W") "Delta-circuit thermal resistances"; Real[4,4] R(each unit="(m.K)/W") "Internal thermal resistances"; Modelica.Units.SI.Position[4] xPip={-sha,sha,0.,0.} "x-Coordinates of pipes"; Modelica.Units.SI.Position[4] yPip={0.,0.,-sha,sha} "y-Coordinates of pipes"; Modelica.Units.SI.Radius[4] rPip={rTub,rTub,rTub,rTub} "Outer radius of pipes"; Real[4] RFluPip(each unit="(m.K)/W") = {RCondPipe+RConv, RCondPipe+RConv, RCondPipe+RConv, RCondPipe+RConv} "Fluid to pipe wall thermal resistances"; Real Ra( unit="(m.K)/W") "Grout-to-grout resistance (2D) as defined by Hellstrom. Interaction between the different grout parts"; Modelica.Units.SI.ThermalResistance Rg "Thermal resistance between outer borehole wall and one tube"; Modelica.Units.SI.ThermalResistance Rar1 "Thermal resistance between the two closest pipe outer walls"; Modelica.Units.SI.ThermalResistance Rar2 "Thermal resistance between the two farthest pipe outer walls"; algorithm // Internal thermal resistances (RDelta, R) := Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.multipoleThermalResistances( 4, 3, xPip, yPip, rBor, rPip, kFil, kSoi, RFluPip); // Rb and Ra Rb_multipole := 1./(1./RDelta[1,1] + 1./RDelta[2,2] + 1./RDelta[3,3] + 1./RDelta[4,4]); Rb_internal := if use_Rb then Rb else Rb_multipole; // The short-circuit resistance in weigthed by the ratio between the used // value of Rb and the theoretical value Ra := (R[1,1] + R[3,3] - 2*R[1,3])*Rb_internal/Rb_multipole; // ------ Calculation according to Bauer et al. (2010) Rg := (4*Rb_internal - RCondPipe - RConv)/hSeg; Rar1 := ((2 + sqrt(2))*Rg*hSeg*(Ra - RCondPipe)/(Rg*hSeg + Ra - RCondPipe))/hSeg; Rar2 := sqrt(2)*Rar1; // If any of the internal delta-circuit resistances is negative, then // the location (x) of the thermal capacity is set to zero to limit // instabilities in the calculations. Otherwise, calculations follow the // method of Bauer et al. (2011). if (RDelta[1,2] < 0) or (RDelta[1,3] < 0) then //Thermal resistance between the grout zone and borehole wall Rgb := Rg; //Conduction resistance in grout from pipe wall to capacity in grout RCondGro := RCondPipe/hSeg; //Thermal resistance between the two grout zones Rgg1 := 2*Rgb*(Rar1)/(2*Rgb - Rar1); Rgg2 := 2*Rgb*(Rar2)/(2*Rgb - Rar2); test := true; else // ********** Resistances and capacity location according to Bauer ********** while test == false and i <= 16 loop // Capacity location (with correction factor in case that the test is negative) x := Modelica.Math.log(sqrt(rBor^2 + 4*rTub^2)/(2*sqrt(2)*rTub))/ Modelica.Math.log(rBor/(2*rTub))*((15 - i + 1)/15); //Thermal resistance between the grout zone and borehole wall Rgb := (1 - x)*Rg; //Conduction resistance in grout from pipe wall to capacity in grout RCondGro := x*Rg + RCondPipe/hSeg; //Thermal resistance between the two grout zones Rgg1 := 2*Rgb*(Rar1 - 2*x*Rg)/(2*Rgb - Rar1 + 2*x*Rg); Rgg2 := 2*Rgb*(Rar2 - 2*x*Rg)/(2*Rgb - Rar2 + 2*x*Rg); // Thermodynamic test to check if negative R values make sense. If not, decrease x-value. test := (((1/Rgg1 + 1/2/Rgb)^(-1) > 0) and ((1/Rgg2 + 1/2/Rgb)^(-1) > 0)); i := i + 1; end while; end if; assert(test, "In " + getInstanceName() + ":\n" + "Maximum number of iterations exceeded. Check the borehole geometry. The tubes may be too close to the borehole wall. Input to the function Buildings.Fluid.HeatExchangers.Boreholes.BaseClasses.doubleUTubeResistances is hSeg = " + String(hSeg) + " m rBor = " + String(rBor) + " m rTub = " + String(rTub) + " m eTub = " + String(eTub) + " m sha = " + String(sha) + " m kSoi = " + String(kSoi) + " W/m/K kFil = " + String(kFil) + " W/m/K kTub = " + String(kTub) + " W/m/K Computed x = " + String(x) + " m Rgb = " + String(Rgb) + " K/W Rgg1 = " + String(Rgg1) + " K/W Rgg2 = " + String(Rgg2) + " K/W"); end internalResistancesTwoUTube;

Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.multipoleFluidTemperature Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.multipoleFluidTemperature

Fluid temperatures from multipole solution

Information

This model evaluates the fluid temperatures using the multipole method of Claesson and Hellstrom (2011).

References

J. Claesson and G. Hellstrom. Multipole method to calculate borehole thermal resistances in a borehole heat exchanger. HVAC&R Research, 17(6): 895-911, 2011.

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

Inputs

TypeNameDefaultDescription
IntegernPip Number of pipes
IntegerJ Number of multipoles
PositionxPip[nPip] x-Coordinates of pipes [m]
PositionyPip[nPip] y-Coordinates of pipes [m]
RealQPip_flow[nPip] Heat flow in pipes [W/m]
TemperatureTBor Average borehole wall temperature [K]
RadiusrBor Borehole radius [m]
RadiusrPip[nPip] Outter radius of pipes [m]
ThermalConductivitykFil Thermal conductivity of grouting material [W/(m.K)]
ThermalConductivitykSoi Thermal conductivity of soil material [W/(m.K)]
RealRFluPip[nPip] Fluid to pipe wall thermal resistances [(m.K)/W]
Realeps1.0e-5Iteration relative accuracy
Integerit_max100Maximum number of iterations

Outputs

TypeNameDescription
TemperatureTFlu[nPip]Fluid temperature in pipes [K]

Modelica definition

function multipoleFluidTemperature "Fluid temperatures from multipole solution" extends Modelica.Icons.Function; input Integer nPip "Number of pipes"; input Integer J "Number of multipoles"; input Modelica.Units.SI.Position xPip[nPip] "x-Coordinates of pipes"; input Modelica.Units.SI.Position yPip[nPip] "y-Coordinates of pipes"; input Real QPip_flow[nPip](each unit="W/m") "Heat flow in pipes"; input Modelica.Units.SI.Temperature TBor "Average borehole wall temperature"; input Modelica.Units.SI.Radius rBor "Borehole radius"; input Modelica.Units.SI.Radius rPip[nPip] "Outter radius of pipes"; input Modelica.Units.SI.ThermalConductivity kFil "Thermal conductivity of grouting material"; input Modelica.Units.SI.ThermalConductivity kSoi "Thermal conductivity of soil material"; input Real RFluPip[nPip](each unit="(m.K)/W") "Fluid to pipe wall thermal resistances"; input Real eps=1.0e-5 "Iteration relative accuracy"; input Integer it_max=100 "Maximum number of iterations"; output Modelica.Units.SI.Temperature TFlu[nPip] "Fluid temperature in pipes"; protected Real pikFil(unit="(m.K)/W")=1/(2*Modelica.Constants.pi*kFil) "Coefficient based on grout thermal conductivity"; Real sigma=(kFil - kSoi)/(kFil + kSoi) "Thermal conductivity ratio"; Real betaPip[nPip]=2*Modelica.Constants.pi*kFil*RFluPip "Dimensionless fluid to outter pipe wall thermal resistance"; Complex zPip_i "Position of pipe i"; Complex zPip_j "Position of pipe j"; Complex P_nj "Multipole of order j for pipe n"; Real PRea[nPip,J] "Matrix of real part of multipoles"; Real PIma[nPip,J] "Matrix of imaginary part of multipole"; Complex P_nj_new "New value of multipole of order j for pipe n"; Real PRea_new[nPip,J] "New value of real part of multipoles"; Real PIma_new[nPip,J] "New value of imaginary part of multipoles"; Complex F_mk "Coefficient F of order k of pipe m"; Real FRea[nPip,J] "Real part of matrix F_mk"; Real FIma[nPip,J] "Imaginary part of matrix F_mk"; Real R0[nPip,nPip](each unit="(m.K)/W") "Line source approximation of thermal resistances"; Complex deltaTFlu "Fluid temperature difference with line source approximation"; Real rbm "Intermediate coefficient"; Modelica.Units.SI.Distance dz "Pipe to pipe distance"; Real coeff[nPip,J] "Coefficient for multiplication with matrix F_mk"; Real diff "Difference in subsequent multipole evaluations"; Real diff_max "Maximum difference in subsequent multipole evaluations"; Real diff_min "Minimum difference in subsequent multipole evaluations"; Real diff0 "Difference in subsequent multipole evaluations"; Integer it "Iteration counter"; Real eps_max "Convergence variable"; algorithm // Thermal resistance matrix from 0th order multipole for i in 1:nPip loop zPip_i := Complex(xPip[i], yPip[i]); rbm :=rBor^2/(rBor^2 - Modelica.ComplexMath.abs(zPip_i)^2); R0[i, i] := pikFil*(log(rBor/rPip[i]) + betaPip[i] + sigma*log(rbm)); for j in 1:nPip loop zPip_j := Complex(xPip[j], yPip[j]); if i <> j then dz :=Modelica.ComplexMath.abs(zPip_i - zPip_j); rbm :=rBor^2/Modelica.ComplexMath.abs(rBor^2 - zPip_j* Modelica.ComplexMath.conj(zPip_i)); R0[i, j] := pikFil*(log(rBor/dz) + sigma*log(rbm)); end if; end for; end for; // Initialize maximum error and iteration counter eps_max := 1.0e99; it := 0; // Multipoles if J > 0 then for m in 1:nPip loop for k in 1:J loop coeff[m, k] := -(1 - k*betaPip[m])/(1 + k*betaPip[m]); PRea[m, k] := 0; PIma[m, k] := 0; end for; end for; while eps_max > eps and it < it_max loop it := it + 1; (FRea, FIma) := Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.multipoleFmk( nPip, J, QPip_flow, PRea, PIma, rBor, rPip, xPip, yPip, kFil, kSoi); for m in 1:nPip loop for k in 1:J loop F_mk := Complex(FRea[m, k], FIma[m, k]); P_nj_new := coeff[m, k]*Modelica.ComplexMath.conj(F_mk); PRea_new[m, k] := Modelica.ComplexMath.real(P_nj_new); PIma_new[m, k] := Modelica.ComplexMath.imag(P_nj_new); end for; end for; diff_max := 0; diff_min := 1e99; for m in 1:nPip loop for k in 1:J loop P_nj := Complex(PRea[m, k], PIma[m, k]); P_nj_new := Complex(PRea_new[m, k], PIma_new[m, k]); diff_max :=max(diff_max, Modelica.ComplexMath.abs(P_nj_new - P_nj)); diff_min :=min(diff_min, Modelica.ComplexMath.abs(P_nj_new - P_nj)); end for; end for; diff := diff_max - diff_min; if it == 1 then diff0 :=diff; end if; eps_max := diff/diff0; PRea := PRea_new; PIma := PIma_new; end while; end if; // Fluid Temperatures TFlu := TBor .+ R0*QPip_flow; if J > 0 then for m in 1:nPip loop zPip_i :=Complex(xPip[m], yPip[m]); deltaTFlu := Complex(0, 0); for n in 1:nPip loop zPip_j :=Complex(xPip[n], yPip[n]); for j in 1:J loop P_nj := Complex(PRea[n, j], PIma[n, j]); if n <> m then // Second term deltaTFlu := deltaTFlu + P_nj*(rPip[n]/(zPip_i - zPip_j))^j; end if; // Third term deltaTFlu := deltaTFlu + sigma*P_nj*(rPip[n]* Modelica.ComplexMath.conj(zPip_i)/(rBor^2 - zPip_j* Modelica.ComplexMath.conj(zPip_i)))^j; end for; end for; TFlu[m] := TFlu[m] + Modelica.ComplexMath.real(deltaTFlu); end for; end if; end multipoleFluidTemperature;

Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.multipoleFmk Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.multipoleFmk

Complex matrix F_mk from Claesson and Hellstrom (2011)

Information

This model evaluates the complex coefficient matrix F_mk from Claesson and Hellstrom (2011).

References

J. Claesson and G. Hellstrom. Multipole method to calculate borehole thermal resistances in a borehole heat exchanger. HVAC&R Research, 17(6): 895-911, 2011.

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

Inputs

TypeNameDefaultDescription
IntegernPip Number of pipes
IntegerJ Number of multipoles
RealQPip_flow[nPip] Heat flow in pipes [W/m]
RealPRea[nPip, J] Multipoles (Real part)
RealPIma[nPip, J] Multipoles (Imaginary part)
RadiusrBor Borehole radius [m]
RadiusrPip[nPip] Outter radius of pipes [m]
PositionxPip[nPip] x-Coordinates of pipes [m]
PositionyPip[nPip] y-Coordinates of pipes [m]
ThermalConductivitykFil Thermal conductivity of grouting material [W/(m.K)]
ThermalConductivitykSoi Thermal conductivity of soil material [W/(m.K)]

Outputs

TypeNameDescription
RealFRea[nPip, J]Multipole coefficients
RealFIma[nPip, J]Multipole coefficients

Modelica definition

function multipoleFmk "Complex matrix F_mk from Claesson and Hellstrom (2011)" extends Modelica.Icons.Function; input Integer nPip "Number of pipes"; input Integer J "Number of multipoles"; input Real QPip_flow[nPip](each unit="W/m") "Heat flow in pipes"; input Real PRea[nPip,J] "Multipoles (Real part)"; input Real PIma[nPip,J] "Multipoles (Imaginary part)"; input Modelica.Units.SI.Radius rBor "Borehole radius"; input Modelica.Units.SI.Radius rPip[nPip] "Outter radius of pipes"; input Modelica.Units.SI.Position xPip[nPip] "x-Coordinates of pipes"; input Modelica.Units.SI.Position yPip[nPip] "y-Coordinates of pipes"; input Modelica.Units.SI.ThermalConductivity kFil "Thermal conductivity of grouting material"; input Modelica.Units.SI.ThermalConductivity kSoi "Thermal conductivity of soil material"; output Real FRea[nPip,J] "Multipole coefficients"; output Real FIma[nPip,J] "Multipole coefficients"; protected Complex zPip_i "Position of pipe i"; Complex zPip_j "Position of pipe j"; Complex P_nj "Multipole of order j for pipe n"; Real pikFil(unit="(m.K)/W")=1/(2*Modelica.Constants.pi*kFil) "Coefficient based on grout thermal conductivity"; Real sigma=(kFil - kSoi)/(kFil + kSoi) "Thermal conductivity ratio"; Complex f "Intermedia value of multipole coefficient"; Integer j_pend "Maximum loop index in fourth term"; algorithm for m in 1:nPip loop zPip_i := Complex(xPip[m], yPip[m]); for k in 1:J loop f := Complex(0, 0); for n in 1:nPip loop zPip_j := Complex(xPip[n], yPip[n]); // First term if m <> n then f := f + QPip_flow[n]*pikFil/k*(rPip[m]/(zPip_j - zPip_i))^k; end if; // Second term f := f + sigma*QPip_flow[n]*pikFil/k*(rPip[m]*Modelica.ComplexMath.conj( zPip_j)/(rBor^2 - zPip_i*Modelica.ComplexMath.conj(zPip_j)))^k; for j in 1:J loop P_nj := Complex(PRea[n, j], PIma[n, j]); // Third term if m <> n then f := f + P_nj*Buildings.Utilities.Math.Functions.binomial(j + k - 1, j - 1)*rPip[n]^j*(-rPip[m])^k/(zPip_i - zPip_j)^(j + k); end if; //Fourth term j_pend := min(k, j); for jp in 0:j_pend loop f := f + sigma*Modelica.ComplexMath.conj(P_nj)* Buildings.Utilities.Math.Functions.binomial(j, jp)* Buildings.Utilities.Math.Functions.binomial(j + k - jp - 1, j - 1)* rPip[n]^j*rPip[m]^k*zPip_i^(j - jp)*Modelica.ComplexMath.conj( zPip_j)^(k - jp)/(rBor^2 - zPip_i*Modelica.ComplexMath.conj( zPip_j))^(k + j - jp); end for; end for; end for; FRea[m, k] := Modelica.ComplexMath.real(f); FIma[m, k] := Modelica.ComplexMath.imag(f); end for; end for; end multipoleFmk;

Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.multipoleThermalResistances Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.multipoleThermalResistances

Thermal resistances from multipole solution

Information

This model evaluates the delta-circuit borehole thermal resistances using the multipole method of Claesson and Hellstrom (2011).

References

J. Claesson and G. Hellstrom. Multipole method to calculate borehole thermal resistances in a borehole heat exchanger. HVAC&R Research, 17(6): 895-911, 2011.

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

Inputs

TypeNameDefaultDescription
IntegernPip Number of pipes
IntegerJ Number of multipoles
PositionxPip[nPip] x-Coordinates of pipes [m]
PositionyPip[nPip] y-Coordinates of pipes [m]
RadiusrBor Borehole radius [m]
RadiusrPip[nPip] Outter radius of pipes [m]
ThermalConductivitykFil Thermal conductivity of grouting material [W/(m.K)]
ThermalConductivitykSoi Thermal conductivity of soil material [W/(m.K)]
RealRFluPip[nPip] Fluid to pipe wall thermal resistances [(m.K)/W]
TemperatureTBor0Average borehole wall temperature [K]

Outputs

TypeNameDescription
RealRDelta[nPip, nPip]Delta-circuit thermal resistances [(m.K)/W]
RealR[nPip, nPip]Internal thermal resistances [(m.K)/W]

Modelica definition

function multipoleThermalResistances "Thermal resistances from multipole solution" extends Modelica.Icons.Function; input Integer nPip "Number of pipes"; input Integer J "Number of multipoles"; input Modelica.Units.SI.Position xPip[nPip] "x-Coordinates of pipes"; input Modelica.Units.SI.Position yPip[nPip] "y-Coordinates of pipes"; input Modelica.Units.SI.Radius rBor "Borehole radius"; input Modelica.Units.SI.Radius rPip[nPip] "Outter radius of pipes"; input Modelica.Units.SI.ThermalConductivity kFil "Thermal conductivity of grouting material"; input Modelica.Units.SI.ThermalConductivity kSoi "Thermal conductivity of soil material"; input Real RFluPip[nPip](each unit="(m.K)/W") "Fluid to pipe wall thermal resistances"; input Modelica.Units.SI.Temperature TBor=0 "Average borehole wall temperature"; output Real RDelta[nPip,nPip](each unit="(m.K)/W") "Delta-circuit thermal resistances"; output Real R[nPip,nPip](each unit="(m.K)/W") "Internal thermal resistances"; protected Real QPip_flow[nPip](each unit="W/m") "Pipe heat transfer rates"; Modelica.Units.SI.Temperature TFlu[nPip] "Fluid temperatures"; Real K[nPip,nPip](each unit="W/(m.K)") "Internal thermal conductances"; algorithm for m in 1:nPip loop for n in 1:nPip loop if n == m then QPip_flow[n] := 1; else QPip_flow[n] := 0; end if; end for; TFlu := Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.multipoleFluidTemperature( nPip, J, xPip, yPip, QPip_flow, TBor, rBor, rPip, kFil, kSoi, RFluPip); for n in 1:nPip loop R[n, m] := TFlu[n]; end for; end for; K := -Modelica.Math.Matrices.inv(R); for m in 1:nPip loop K[m, m] := -K[m, m]; for n in 1:nPip loop if m <> n then K[m, m] := K[m, m] - K[m, n]; end if; end for; end for; for m in 1:nPip loop for n in 1:nPip loop RDelta[m, n] := 1./K[m, n]; end for; end for; end multipoleThermalResistances;

Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.partialInternalResistances Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.partialInternalResistances

Partial model for borehole resistance calculation

Information

This partial function defines the common inputs to functions that calculate the borehole internal resistances.

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

Inputs

TypeNameDefaultDescription
Booleanuse_RbfalseTrue if the value Rb should be used instead of calculated
RealRb Borehole thermal resistance [(m.K)/W]
HeighthSeg Height of the element [m]
RadiusrBor Radius of the borehole [m]
RadiusrTub Radius of the tube [m]
LengtheTub Thickness of the tubes [m]
Lengthsha Shank spacing, defined as the distance between the center of a pipe and the center of the borehole [m]
ThermalConductivitykFil Thermal conductivity of the grout [W/(m.K)]
ThermalConductivitykSoi Thermal conductivity of the soi [W/(m.K)]
ThermalConductivitykTub Thermal conductivity of the tube [W/(m.K)]
ThermalConductivitykMed Thermal conductivity of the fluid [W/(m.K)]
DynamicViscositymuMed Dynamic viscosity of the fluid [Pa.s]
SpecificHeatCapacitycpMed Specific heat capacity of the fluid [J/(kg.K)]
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]

Outputs

TypeNameDescription
RealxCapacity location

Modelica definition

partial function partialInternalResistances "Partial model for borehole resistance calculation" extends Modelica.Icons.Function; // Geometry of the borehole input Boolean use_Rb = false "True if the value Rb should be used instead of calculated"; input Real Rb(unit="(m.K)/W") "Borehole thermal resistance"; input Modelica.Units.SI.Height hSeg "Height of the element"; input Modelica.Units.SI.Radius rBor "Radius of the borehole"; // Geometry of the pipe input Modelica.Units.SI.Radius rTub "Radius of the tube"; input Modelica.Units.SI.Length eTub "Thickness of the tubes"; input Modelica.Units.SI.Length sha "Shank spacing, defined as the distance between the center of a pipe and the center of the borehole"; // Thermal properties input Modelica.Units.SI.ThermalConductivity kFil "Thermal conductivity of the grout"; input Modelica.Units.SI.ThermalConductivity kSoi "Thermal conductivity of the soi"; input Modelica.Units.SI.ThermalConductivity kTub "Thermal conductivity of the tube"; input Modelica.Units.SI.ThermalConductivity kMed "Thermal conductivity of the fluid"; input Modelica.Units.SI.DynamicViscosity muMed "Dynamic viscosity of the fluid"; input Modelica.Units.SI.SpecificHeatCapacity cpMed "Specific heat capacity of the fluid"; input Modelica.Units.SI.MassFlowRate m_flow_nominal "Nominal mass flow rate"; // Outputs output Real x "Capacity location"; protected parameter Real pi = 3.141592653589793 "pi"; parameter Real rTub_in = rTub-eTub "Inner radius of tube"; Real RConv(unit="(m.K)/W")= Buildings.Fluid.Geothermal.Borefields.BaseClasses.Boreholes.BaseClasses.Functions.convectionResistanceCircularPipe( hSeg=hSeg, rTub=rTub, eTub=eTub, kMed=kMed, muMed=muMed, cpMed=cpMed, m_flow=m_flow_nominal, m_flow_nominal=m_flow_nominal)*hSeg; Boolean test=false "thermodynamic test for R and x value"; Real RCondPipe(unit="(m.K)/W") = Modelica.Math.log((rTub)/rTub_in)/(2*Modelica.Constants.pi*kTub) "Thermal resistance of the pipe wall"; Real Rb_internal(unit="(m.K)/W") "Resistance from the fluid in the pipe to the borehole wall"; Real Rb_multipole(unit="(m.K)/W") "Theoretical Fluid-to-borehole-wall resistance evaluated from the multipole method"; Integer i=1 "Loop counter"; end partialInternalResistances;