Buildings.Fluid.Geothermal.Boreholes.BaseClasses

Base classes for Borehole

Information

This package contains base classes that are used to construct the models in Buildings.Fluid.Geothermal.Boreholes.

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

Package Content

Name Description
Buildings.Fluid.Geothermal.Boreholes.BaseClasses.BoreholeSegment BoreholeSegment Vertical segment of a borehole
Buildings.Fluid.Geothermal.Boreholes.BaseClasses.ExtendableArray ExtendableArray class used to create the external object: ExtendableArray
Buildings.Fluid.Geothermal.Boreholes.BaseClasses.HexInternalElement HexInternalElement Internal part of a borehole
Buildings.Fluid.Geothermal.Boreholes.BaseClasses.SingleUTubeBoundaryCondition SingleUTubeBoundaryCondition Prescribed temperature at the outer boundary of a single U tube borehole
Buildings.Fluid.Geothermal.Boreholes.BaseClasses.convectionResistance convectionResistance Thermal resistance between the fluid and the tube
Buildings.Fluid.Geothermal.Boreholes.BaseClasses.exchangeValues exchangeValues Store the value u at the element n in the external object, and return the value of the element i
Buildings.Fluid.Geothermal.Boreholes.BaseClasses.powerSeries powerSeries Power series used to compute far-field temperature
Buildings.Fluid.Geothermal.Boreholes.BaseClasses.singleUTubeResistances singleUTubeResistances Thermal resistances for single U-tube
Buildings.Fluid.Geothermal.Boreholes.BaseClasses.temperatureDrop temperatureDrop Calculate the temperature drop of the soil at the external boundary of the cylinder
Buildings.Fluid.Geothermal.Boreholes.BaseClasses.Examples Examples Example models to test base classes

Buildings.Fluid.Geothermal.Boreholes.BaseClasses.BoreholeSegment Buildings.Fluid.Geothermal.Boreholes.BaseClasses.BoreholeSegment

Vertical segment of a borehole

Buildings.Fluid.Geothermal.Boreholes.BaseClasses.BoreholeSegment

Information

Horizontal layer that is used to model a U-tube borehole heat exchanger. This model combines three models, each simulating a different aspect of a borehole heat exchanger.

The instance pipFil computes the heat transfer in the pipes and the filling material. This computation is done using the model Buildings.Fluid.Geothermal.Boreholes.BaseClasses.HexInternalElement.

The instance soi computes transient and steady state heat transfer in the soil using a vertical cylinder. The computation is done using the model Buildings.HeatTransfer.Conduction.SingleLayerCylinder.

The model TBouCon computes the far-field temperature boundary condition, i.e., the temperature at the outer surface of the above cylindrical heat transfer computation. The computation is done using the model Buildings.Fluid.Geothermal.Boreholes.BaseClasses.SingleUTubeBoundaryCondition.

Extends from Buildings.Fluid.Interfaces.PartialFourPortInterface (Partial model transporting fluid between two ports without storing mass or energy), Buildings.Fluid.Interfaces.TwoPortFlowResistanceParameters (Parameters for flow resistance for models with two ports), Buildings.Fluid.Interfaces.LumpedVolumeDeclarations (Declarations for lumped volumes).

Parameters

TypeNameDefaultDescription
replaceable package Medium1PartialMediumMedium 1 in the component
replaceable package Medium2PartialMediumMedium 2 in the component
replaceable package MediumPartialMediumMedium in the component
RadiusrBor0.1Radius of the borehole [m]
HeighthSeg Height of the element [m]
LengthxC0.05Shank spacing, defined as the distance between the center of a pipe and the center of the borehole [m]
Nominal condition
MassFlowRatem1_flow_nominalm_flow_nominalNominal mass flow rate [kg/s]
MassFlowRatem2_flow_nominalm_flow_nominalNominal mass flow rate [kg/s]
PressureDifferencedp_nominal Pressure difference [Pa]
MassFlowRatem_flow_nominal Nominal mass flow rate [kg/s]
Soil
GenericmatSoiredeclare parameter Building...Thermal properties of soil
RadiusrExt3Radius of the soil used for the external boundary condition [m]
TemperatureTExt_start283.15Initial far field temperature [K]
IntegernSta10Number of state variables in the soil
TimesamplePeriod604800Sample period for the external boundary condition [s]
Filling material
GenericmatFilredeclare parameter Building...Thermal properties of the filling material
TemperatureTFil_start283.15Initial temperature of the filling material [K]
Tubes
RadiusrTub0.02Radius of the tubes [m]
ThermalConductivitykTub0.5Thermal conductivity of the tubes [W/(m.K)]
LengtheTub0.002Thickness of the tubes [m]
Assumptions
BooleanallowFlowReversal1allowFlowReversal= false to simplify equations, assuming, but not enforcing, no flow reversal for medium 1
BooleanallowFlowReversal2allowFlowReversal= false to simplify equations, assuming, but not enforcing, no flow reversal for medium 2
BooleanallowFlowReversaltrue= true to allow flow reversal, false restricts to design direction (port_a -> port_b)
Advanced
MassFlowRatem1_flow_smallm_flow_smallSmall mass flow rate for regularization of zero flow [kg/s]
MassFlowRatem2_flow_smallm_flow_smallSmall mass flow rate for regularization of zero flow [kg/s]
MassFlowRatem_flow_small1E-4*abs(m_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
Flow resistance
BooleancomputeFlowResistancetrue=true, compute flow resistance. Set to false to assume no friction
Booleanfrom_dpfalse= true, use m_flow = f(dp) else dp = f(m_flow)
BooleanlinearizeFlowResistancefalse= true, use linear relation between m_flow and dp for any flow rate
RealdeltaM0.1Fraction of nominal flow rate where flow transitions to laminar
Dynamics
Equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Type of energy balance: dynamic (3 initialization options) or steady state
DynamicsmassDynamicsenergyDynamicsType of mass balance: dynamic (3 initialization options) or steady state
RealmSenFac1Factor for scaling the sensible thermal mass of the volume
Initialization
AbsolutePressurep_startMedium.p_defaultStart value of pressure [Pa]
TemperatureT_startTFil_startStart value of temperature [K]
MassFractionX_start[Medium.nX]Medium.X_defaultStart value of mass fractions m_i/m [kg/kg]
ExtraPropertyC_start[Medium.nC]fill(0, Medium.nC)Start value of trace substances
ExtraPropertyC_nominal[Medium.nC]fill(1E-2, Medium.nC)Nominal value of trace substances. (Set to typical order of magnitude.)

Connectors

TypeNameDescription
replaceable package Medium1Medium 1 in the component
replaceable package Medium2Medium 2 in the component
FluidPort_aport_a1Fluid connector a1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_bport_b1Fluid connector b1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_aport_a2Fluid connector a2 (positive design flow direction is from port_a2 to port_b2)
FluidPort_bport_b2Fluid connector b2 (positive design flow direction is from port_a2 to port_b2)
replaceable package MediumMedium in the component

Modelica definition

model BoreholeSegment "Vertical segment of a borehole" extends Buildings.Fluid.Interfaces.PartialFourPortInterface( redeclare final package Medium1 = Medium, redeclare final package Medium2 = Medium, final m1_flow_nominal = m_flow_nominal, final m2_flow_nominal = m_flow_nominal, final m1_flow_small = m_flow_small, final m2_flow_small = m_flow_small, final allowFlowReversal1=allowFlowReversal, final allowFlowReversal2=allowFlowReversal); extends Buildings.Fluid.Interfaces.TwoPortFlowResistanceParameters; extends Buildings.Fluid.Interfaces.LumpedVolumeDeclarations(T_start=TFil_start); replaceable package Medium = Modelica.Media.Interfaces.PartialMedium "Medium in the component"; constant Boolean homotopyInitialization = true "= true, use homotopy method"; replaceable parameter Buildings.HeatTransfer.Data.Soil.Generic matSoi "Thermal properties of soil"; replaceable parameter Buildings.HeatTransfer.Data.BoreholeFillings.Generic matFil "Thermal properties of the filling material"; parameter Modelica.SIunits.MassFlowRate m_flow_nominal "Nominal mass flow rate"; parameter Modelica.SIunits.MassFlowRate m_flow_small(min=0) = 1E-4*abs(m_flow_nominal) "Small mass flow rate for regularization of zero flow"; parameter Modelica.SIunits.Radius rTub=0.02 "Radius of the tubes"; parameter Modelica.SIunits.ThermalConductivity kTub=0.5 "Thermal conductivity of the tubes"; parameter Modelica.SIunits.Length eTub=0.002 "Thickness of the tubes"; parameter Modelica.SIunits.Temperature TFil_start=283.15 "Initial temperature of the filling material"; parameter Modelica.SIunits.Radius rExt=3 "Radius of the soil used for the external boundary condition"; parameter Modelica.SIunits.Temperature TExt_start=283.15 "Initial far field temperature"; parameter Integer nSta(min=1) = 10 "Number of state variables in the soil"; parameter Modelica.SIunits.Time samplePeriod=604800 "Sample period for the external boundary condition"; parameter Modelica.SIunits.Radius rBor=0.1 "Radius of the borehole"; parameter Modelica.SIunits.Height hSeg "Height of the element"; parameter Modelica.SIunits.Length xC=0.05 "Shank spacing, defined as the distance between the center of a pipe and the center of the borehole"; parameter Boolean allowFlowReversal = true "= true to allow flow reversal, false restricts to design direction (port_a -> port_b)"; Buildings.Fluid.Geothermal.Boreholes.BaseClasses.HexInternalElement pipFil( redeclare final package Medium = Medium, final matFil=matFil, final matSoi=matSoi, final hSeg=hSeg, final rTub=rTub, final eTub=eTub, final kTub=kTub, final kSoi=matSoi.k, final xC=xC, final rBor=rBor, final TFil_start=TFil_start, final m1_flow_nominal=m_flow_nominal, final m2_flow_nominal=m_flow_nominal, final dp1_nominal=dp_nominal, final dp2_nominal=0, final from_dp1=from_dp, final from_dp2=from_dp, final linearizeFlowResistance1=linearizeFlowResistance, final linearizeFlowResistance2=linearizeFlowResistance, final deltaM1=deltaM, final deltaM2=deltaM, final m1_flow_small=m_flow_small, final m2_flow_small=m_flow_small, final allowFlowReversal1=allowFlowReversal, final allowFlowReversal2=allowFlowReversal, final homotopyInitialization=homotopyInitialization, final energyDynamics=energyDynamics, final massDynamics=massDynamics, final p1_start=p_start, T1_start=T_start, X1_start=X_start, C1_start=C_start, C1_nominal=C_nominal, final p2_start=p_start, T2_start=T_start, X2_start=X_start, C2_start=C_start, C2_nominal=C_nominal) "Internal part of the borehole including the pipes and the filling material"; Buildings.HeatTransfer.Conduction.SingleLayerCylinder soi( final material=matSoi, final h=hSeg, final nSta=nSta, final r_a=rBor, final r_b=rExt, final steadyStateInitial=false, final TInt_start=TFil_start, final TExt_start=TExt_start) "Heat conduction in the soil"; Buildings.Fluid.Geothermal.Boreholes.BaseClasses.SingleUTubeBoundaryCondition TBouCon( final matSoi=matSoi, final rExt=rExt, final hSeg=hSeg, final TExt_start=TExt_start, final samplePeriod=samplePeriod) "Thermal boundary condition for the far-field"; protected Modelica.Thermal.HeatTransfer.Sensors.HeatFlowSensor heaFlo; initial equation assert(homotopyInitialization, "In " + getInstanceName() + ": The constant homotopyInitialization has been modified from its default value. This constant will be removed in future releases.", level = AssertionLevel.warning); equation connect(pipFil.port_b1, port_b1); connect(pipFil.port_a2, port_a2); connect(pipFil.port_b2, port_b2); connect(pipFil.port, heaFlo.port_a); connect(heaFlo.port_b, soi.port_a); connect(soi.port_b, TBouCon.port); connect(port_a1, pipFil.port_a1); connect(heaFlo.Q_flow, TBouCon.Q_flow); end BoreholeSegment;

Buildings.Fluid.Geothermal.Boreholes.BaseClasses.ExtendableArray

class used to create the external object: ExtendableArray

Information

Class derived from ExternalObject having two local external function definition, named destructor and constructor respectively.

These functions create and release an external object that allows the storage of real parameters in an array of extendable dimension.

Extends from ExternalObject.

Modelica definition

class ExtendableArray "class used to create the external object: ExtendableArray" extends ExternalObject; function constructor "Construct an extendable array that can be used to store double values" output ExtendableArray table; external "C" table = initArray(); end constructor; function destructor "Release storage of table and close the external object" input ExtendableArray table; external "C" freeArray(table); end destructor; end ExtendableArray;

Buildings.Fluid.Geothermal.Boreholes.BaseClasses.HexInternalElement Buildings.Fluid.Geothermal.Boreholes.BaseClasses.HexInternalElement

Internal part of a borehole

Buildings.Fluid.Geothermal.Boreholes.BaseClasses.HexInternalElement

Information

Model for the heat transfer between the fluid and within the borehole filling. This model computes the dynamic response of the fluid in the tubes, the heat transfer between the fluid and the borehole filling, and the heat storage within the fluid and the borehole filling.

This model computes the different thermal resistances present in a single-U-tube borehole using the method of Bauer et al. (2011) and computing explicitly the fluid-to-ground thermal resistance Rb and the grout-to-grout resistance Ra as defined by Hellstroem (1991) using the multipole method. The multipole method is implemented in Buildings.Fluid.Geothermal.Boreholes.BaseClasses.singleUTubeResistances. The convection resistance is calculated using the Dittus-Boelter correlation as implemented in Buildings.Fluid.Geothermal.Boreholes.BaseClasses.convectionResistance.

The figure below shows the thermal network set up by Bauer et al. (2010).

image

References

G. Hellström. Ground heat storage: thermal analyses of duct storage systems (Theory). Dept. of Mathematical Physics, University of Lund, Sweden, 1991.

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.Interfaces.FourPortHeatMassExchanger (Model transporting two fluid streams between four ports with storing mass or energy).

Parameters

TypeNameDefaultDescription
replaceable package Medium1PartialMediumMedium 1 in the component
replaceable package Medium2PartialMediumMedium 2 in the component
replaceable package MediumModelica.Media.Interfaces.Pa...Medium in the component
ThermalConductivitykSoi Thermal conductivity of the soil used for the calculation of the internal interference resistance [W/(m.K)]
HeighthSeg Height of the element [m]
RadiusrBor Radius of the borehole [m]
LengthxC0.05Shank spacing, defined as half the center-to-center distance between the two pipes [m]
Nominal condition
MassFlowRatem1_flow_nominal Nominal mass flow rate [kg/s]
MassFlowRatem2_flow_nominal Nominal mass flow rate [kg/s]
PressureDifferencedp1_nominal Pressure difference [Pa]
PressureDifferencedp2_nominal Pressure difference [Pa]
Filling material
GenericmatFilredeclare parameter Building...Thermal properties of the filling material
TemperatureTFil_start283.15Initial temperature of the filling material [K]
Soil
GenericmatSoiredeclare parameter Building...Thermal properties of soil
Pipes
RadiusrTub0.02Radius of the tubes [m]
ThermalConductivitykTub0.5Thermal conductivity of the tubes [W/(m.K)]
LengtheTub0.002Thickness of the tubes [m]
Assumptions
BooleanallowFlowReversal1true= false to simplify equations, assuming, but not enforcing, no flow reversal for medium 1
BooleanallowFlowReversal2true= false to simplify equations, assuming, but not enforcing, no flow reversal for medium 2
Advanced
MassFlowRatem1_flow_small1E-4*abs(m1_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
MassFlowRatem2_flow_small1E-4*abs(m2_flow_nominal)Small mass flow rate for regularization of zero flow [kg/s]
Diagnostics
Booleanshow_Tfalse= true, if actual temperature at port is computed
Flow resistance
Medium 1
Booleanfrom_dp1false= true, use m_flow = f(dp) else dp = f(m_flow)
BooleanlinearizeFlowResistance1false= true, use linear relation between m_flow and dp for any flow rate
RealdeltaM10.1Fraction of nominal flow rate where flow transitions to laminar
Medium 2
Booleanfrom_dp2false= true, use m_flow = f(dp) else dp = f(m_flow)
BooleanlinearizeFlowResistance2false= true, use linear relation between m_flow and dp for any flow rate
RealdeltaM20.1Fraction of nominal flow rate where flow transitions to laminar
Dynamics
Nominal condition
Timetau1Modelica.Constants.pi*rTub^2...Time constant at nominal flow [s]
Timetau2Modelica.Constants.pi*rTub^2...Time constant at nominal flow [s]
Equations
DynamicsenergyDynamicsModelica.Fluid.Types.Dynamic...Type of energy balance: dynamic (3 initialization options) or steady state
DynamicsmassDynamicsenergyDynamicsType of mass balance: dynamic (3 initialization options) or steady state
Initialization
Medium 1
AbsolutePressurep1_startMedium1.p_defaultStart value of pressure [Pa]
TemperatureT1_startTFil_startStart value of temperature [K]
MassFractionX1_start[Medium1.nX]Medium1.X_defaultStart value of mass fractions m_i/m [kg/kg]
ExtraPropertyC1_start[Medium1.nC]fill(0, Medium1.nC)Start value of trace substances
ExtraPropertyC1_nominal[Medium1.nC]fill(1E-2, Medium1.nC)Nominal value of trace substances. (Set to typical order of magnitude.)
Medium 2
AbsolutePressurep2_startMedium2.p_defaultStart value of pressure [Pa]
TemperatureT2_startTFil_startStart value of temperature [K]
MassFractionX2_start[Medium2.nX]Medium2.X_defaultStart value of mass fractions m_i/m [kg/kg]
ExtraPropertyC2_start[Medium2.nC]fill(0, Medium2.nC)Start value of trace substances
ExtraPropertyC2_nominal[Medium2.nC]fill(1E-2, Medium2.nC)Nominal value of trace substances. (Set to typical order of magnitude.)

Connectors

TypeNameDescription
replaceable package Medium1Medium 1 in the component
replaceable package Medium2Medium 2 in the component
FluidPort_aport_a1Fluid connector a1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_bport_b1Fluid connector b1 (positive design flow direction is from port_a1 to port_b1)
FluidPort_aport_a2Fluid connector a2 (positive design flow direction is from port_a2 to port_b2)
FluidPort_bport_b2Fluid connector b2 (positive design flow direction is from port_a2 to port_b2)
replaceable package MediumMedium in the component
HeatPort_aportHeat port that connects to filling material

Modelica definition

model HexInternalElement "Internal part of a borehole" extends Buildings.Fluid.Interfaces.FourPortHeatMassExchanger( redeclare final package Medium1 = Medium, redeclare final package Medium2 = Medium, T1_start=TFil_start, T2_start=TFil_start, final tau1=Modelica.Constants.pi*rTub^2*hSeg*rho1_nominal/m1_flow_nominal, final tau2=Modelica.Constants.pi*rTub^2*hSeg*rho2_nominal/m2_flow_nominal, vol1(final energyDynamics=energyDynamics, final massDynamics=massDynamics, final prescribedHeatFlowRate=false, final V=m2_flow_nominal*tau2/rho2_nominal, final m_flow_small=m1_flow_small), final vol2(final energyDynamics=energyDynamics, final massDynamics=massDynamics, final prescribedHeatFlowRate=false, final V=m1_flow_nominal*tau1/rho1_nominal, final m_flow_small=m2_flow_small)); replaceable package Medium = Modelica.Media.Interfaces.PartialMedium "Medium in the component"; replaceable parameter Buildings.HeatTransfer.Data.BoreholeFillings.Generic matFil "Thermal properties of the filling material"; replaceable parameter Buildings.HeatTransfer.Data.Soil.Generic matSoi "Thermal properties of soil"; parameter Modelica.SIunits.Radius rTub=0.02 "Radius of the tubes"; parameter Modelica.SIunits.ThermalConductivity kTub=0.5 "Thermal conductivity of the tubes"; parameter Modelica.SIunits.Length eTub=0.002 "Thickness of the tubes"; parameter Modelica.SIunits.ThermalConductivity kSoi "Thermal conductivity of the soil used for the calculation of the internal interference resistance"; parameter Modelica.SIunits.Temperature TFil_start=283.15 "Initial temperature of the filling material"; parameter Modelica.SIunits.Height hSeg "Height of the element"; parameter Modelica.SIunits.Radius rBor "Radius of the borehole"; parameter Modelica.SIunits.Length xC=0.05 "Shank spacing, defined as half the center-to-center distance between the two pipes"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_a port "Heat port that connects to filling material"; Modelica.Thermal.HeatTransfer.Components.HeatCapacitor capFil1( final C=Co_fil/2, T(final start=TFil_start, fixed=(energyDynamics == Modelica.Fluid.Types.Dynamics.FixedInitial)), der_T(fixed=(energyDynamics == Modelica.Fluid.Types.Dynamics.SteadyStateInitial))) "Heat capacity of the filling material"; Modelica.Thermal.HeatTransfer.Components.HeatCapacitor capFil2( final C=Co_fil/2, T(final start=TFil_start, fixed=(energyDynamics == Modelica.Fluid.Types.Dynamics.FixedInitial)), der_T(fixed=(energyDynamics == Modelica.Fluid.Types.Dynamics.SteadyStateInitial))) "Heat capacity of the filling material"; protected final parameter Modelica.SIunits.SpecificHeatCapacity cpFil=matFil.c "Specific heat capacity of the filling material"; final parameter Modelica.SIunits.ThermalConductivity kFil=matFil.k "Thermal conductivity of the filling material"; final parameter Modelica.SIunits.Density dFil=matFil.d "Density of the filling material"; parameter Modelica.SIunits.HeatCapacity Co_fil=dFil*cpFil*hSeg*Modelica.Constants.pi*(rBor^2 - 2*(rTub + eTub)^2) "Heat capacity of the filling material"; parameter Modelica.SIunits.SpecificHeatCapacity cpMed= Medium.specificHeatCapacityCp(Medium.setState_pTX( Medium.p_default, Medium.T_default, Medium.X_default)) "Specific heat capacity of the fluid"; parameter Modelica.SIunits.ThermalConductivity kMed= Medium.thermalConductivity(Medium.setState_pTX( Medium.p_default, Medium.T_default, Medium.X_default)) "Thermal conductivity of the fluid"; parameter Modelica.SIunits.DynamicViscosity mueMed=Medium.dynamicViscosity( Medium.setState_pTX( Medium.p_default, Medium.T_default, Medium.X_default)) "Dynamic viscosity of the fluid"; parameter Modelica.SIunits.ThermalResistance Rgb_val(fixed=false) "Thermal resistance between grout zone and borehole wall"; parameter Modelica.SIunits.ThermalResistance Rgg_val(fixed=false) "Thermal resistance between the two grout zones"; parameter Modelica.SIunits.ThermalResistance RCondGro_val(fixed=false) "Thermal resistance of the pipe wall"; parameter Real x(fixed=false) "Capacity location"; Modelica.Thermal.HeatTransfer.Components.ConvectiveResistor RConv1 "Pipe convective resistance"; Modelica.Thermal.HeatTransfer.Components.ConvectiveResistor RConv2 "Pipe convective resistance"; Modelica.Thermal.HeatTransfer.Components.ThermalResistor Rpg1( final R=RCondGro_val) "Grout thermal resistance"; Modelica.Thermal.HeatTransfer.Components.ThermalResistor Rpg2( final R=RCondGro_val) "Grout thermal resistance"; Modelica.Thermal.HeatTransfer.Components.ThermalResistor Rgb1( final R=Rgb_val) "Grout thermal resistance"; Modelica.Thermal.HeatTransfer.Components.ThermalResistor Rgb2( final R=Rgb_val) "Grout thermal resistance"; Modelica.Thermal.HeatTransfer.Components.ThermalResistor Rgg( final R=Rgg_val) "Grout thermal resistance"; Modelica.Blocks.Sources.RealExpression RVol1(y= convectionResistance( hSeg=hSeg, rTub=rTub, kMed=kMed, mueMed=mueMed, cpMed=cpMed, m_flow=m1_flow, m_flow_nominal=m1_flow_nominal)) "Convective and thermal resistance at fluid 1"; Modelica.Blocks.Sources.RealExpression RVol2(y= convectionResistance( hSeg=hSeg, rTub=rTub, kMed=kMed, mueMed=mueMed, cpMed=cpMed, m_flow=m2_flow, m_flow_nominal=m2_flow_nominal)) "Convective and thermal resistance at fluid 2"; initial equation (Rgb_val, Rgg_val, RCondGro_val, x) = singleUTubeResistances( hSeg=hSeg, rBor=rBor, rTub=rTub, eTub=eTub, xC=xC, kSoi=matSoi.k, kFil=matFil.k, kTub=kTub); equation connect(vol1.heatPort, RConv1.fluid); connect(RConv1.solid, Rpg1.port_a); connect(Rpg1.port_b, capFil1.port); connect(capFil1.port, Rgb1.port_a); connect(capFil1.port, Rgg.port_a); connect(Rgb1.port_b, port); connect(RConv2.solid, Rpg2.port_a); connect(Rpg2.port_b, capFil2.port); connect(RConv2.fluid, vol2.heatPort); connect(capFil2.port, Rgb2.port_a); connect(Rgg.port_b, capFil2.port); connect(Rgb2.port_b, port); connect(RVol1.y, RConv1.Rc); connect(RVol2.y, RConv2.Rc); end HexInternalElement;

Buildings.Fluid.Geothermal.Boreholes.BaseClasses.SingleUTubeBoundaryCondition Buildings.Fluid.Geothermal.Boreholes.BaseClasses.SingleUTubeBoundaryCondition

Prescribed temperature at the outer boundary of a single U tube borehole

Buildings.Fluid.Geothermal.Boreholes.BaseClasses.SingleUTubeBoundaryCondition

Information

This model computes the temperature boundary condition at the outer boundary of the borehole. It takes as an input the heat flow rate at the center of the borehole. This heat flow rate is averaged over the sample period. At each sampling interval, typically every one week, a new temperature boundary conditions is computed using the analytical solution to a line source heat transfer problem.

Implementation

The computation of the temperature change of the boundary is computed using the function Buildings.Fluid.Geothermal.Boreholes.BaseClasses.temperatureDrop.

Parameters

TypeNameDefaultDescription
GenericmatSoiredeclare parameter Building...Thermal properties of the soil
RadiusrExt3Distance from the brine where the calculation is performed [m]
HeighthSeg10Height of the segment [m]
TemperatureTExt_start283.15Initial external temperature [K]
TimesamplePeriod604800Period between two samples [s]

Connectors

TypeNameDescription
input RealInputQ_flowHeat flow rate at the center of the borehole, positive if heat is added to soil [W]
HeatPort_bportHeat port

Modelica definition

model SingleUTubeBoundaryCondition "Prescribed temperature at the outer boundary of a single U tube borehole" replaceable parameter Buildings.HeatTransfer.Data.Soil.Generic matSoi "Thermal properties of the soil"; parameter Modelica.SIunits.Radius rExt=3 "Distance from the brine where the calculation is performed"; parameter Modelica.SIunits.Height hSeg=10 "Height of the segment"; parameter Modelica.SIunits.Temperature TExt_start=283.15 "Initial external temperature"; parameter Modelica.SIunits.Time samplePeriod=604800 "Period between two samples"; ExtendableArray table=ExtendableArray() "Extentable array, used to store history of rate of heat flows"; Modelica.SIunits.HeatFlowRate QAve_flow "Average heat flux over a time period"; Modelica.Blocks.Interfaces.RealInput Q_flow(unit="W") "Heat flow rate at the center of the borehole, positive if heat is added to soil"; Modelica.Thermal.HeatTransfer.Interfaces.HeatPort_b port "Heat port"; protected final parameter Modelica.SIunits.SpecificHeatCapacity c= matSoi.c "Specific heat capacity of the soil"; final parameter Modelica.SIunits.ThermalConductivity k= matSoi.k "Thermal conductivity of the soil"; final parameter Modelica.SIunits.Density d = matSoi.d "Density of the soil"; Modelica.SIunits.Energy UOld "Internal energy at the previous period"; Modelica.SIunits.Energy U "Current internal energy, defined as U=0 for t=tStart"; final parameter Modelica.SIunits.Time startTime(fixed=false) "Start time of the simulation"; Integer iSam(min=1) "Counter for how many time the model was sampled. Defined as iSam=1 when called at t=0"; initial equation U = 0; UOld = 0; startTime = time; iSam = 1; port.T = TExt_start; QAve_flow = 0; equation der(U) = Q_flow; when sample(startTime, samplePeriod) then QAve_flow = (U-pre(UOld))/samplePeriod; UOld = U; port.T = TExt_start + Buildings.Fluid.Geothermal.Boreholes.BaseClasses.temperatureDrop( table=table, iSam=pre(iSam), Q_flow=QAve_flow, samplePeriod=samplePeriod, rExt=rExt, hSeg=hSeg, k=k, d=d, c=c); iSam = pre(iSam)+1; end when; end SingleUTubeBoundaryCondition;

Buildings.Fluid.Geothermal.Boreholes.BaseClasses.convectionResistance

Thermal resistance between the fluid and the tube

Information

This model computes the convection resistance in the pipes of a borehole segment with heigth hSeg.

The correlation of Dittus-Boelter (1930) is used to find the convection heat transfer coefficient as

Nu = 0.023   Re0.8   Prn,

where Nu is the Nusselt number, Re is the Reynolds number and Pr is the Prandlt number. We selected n=0.35, as the reference uses n=0.4 for heating and n=0.3 for cooling. Dittus-Boelter's correlation is valid for turbulent flow in cylindrical smooth pipe.

+

References

Dittus P.W. and L.M.K Boelter, (1930). Heat transfer in automobile radiators of the tubular type. Univ Calif Pub Eng, 2(13):443-461. (Reprinted in Int. J. Comm. Heat Mass Transf. 12 (1985), 3:22). DOI:10.1016/0735-1933(85)90003-X.

Inputs

TypeNameDefaultDescription
HeighthSeg Height of the element [m]
RadiusrTub Tube radius [m]
ThermalConductivitykMed Thermal conductivity of the fluid [W/(m.K)]
DynamicViscositymueMed 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
ThermalResistanceRThermal resistance between the fluid and the tube [K/W]

Modelica definition

function convectionResistance "Thermal resistance between the fluid and the tube" // Geometry of the borehole input Modelica.SIunits.Height hSeg "Height of the element"; input Modelica.SIunits.Radius rTub "Tube radius"; // Thermal properties input Modelica.SIunits.ThermalConductivity kMed "Thermal conductivity of the fluid"; input Modelica.SIunits.DynamicViscosity mueMed "Dynamic viscosity of the fluid"; input Modelica.SIunits.SpecificHeatCapacity cpMed "Specific heat capacity of the fluid"; // Mass flow rates input Modelica.SIunits.MassFlowRate m_flow "Mass flow rate"; input Modelica.SIunits.MassFlowRate m_flow_nominal "Nominal mass flow rate"; // Outputs output Modelica.SIunits.ThermalResistance R "Thermal resistance between the fluid and the tube"; protected Modelica.SIunits.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"; algorithm // ********** Convection resistance ********** // Dittus-Boelter: 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) k := 2/(mueMed*Modelica.Constants.pi*rTub); // Convection h := 0.023*kMed*(cpMed*mueMed/kMed)^(0.35)/(2*rTub)* Buildings.Utilities.Math.Functions.regNonZeroPower( x=m_flow*k, n=0.8, delta=0.01*m_flow_nominal*k); R := 1/(2*Modelica.Constants.pi*rTub*hSeg*h); end convectionResistance;

Buildings.Fluid.Geothermal.Boreholes.BaseClasses.exchangeValues

Store the value u at the element n in the external object, and return the value of the element i

Information

External function that stores the input value x as the element a[iX] in the array a = [a[1], a[2], ...], and that returns the element a[iY]. The size of the array a is automatically enlarged as needed.

Inputs

TypeNameDefaultDescription
ExtendableArraytable External object
IntegeriX One-based index where u needs to be stored in the array of the external object
Realx Value to store in the external object
IntegeriY One-based index of the element that needs to be returned from C to Modelica

Outputs

TypeNameDescription
RealyValue of the i-th element

Modelica definition

function exchangeValues "Store the value u at the element n in the external object, and return the value of the element i" input ExtendableArray table "External object"; input Integer iX "One-based index where u needs to be stored in the array of the external object"; input Real x "Value to store in the external object"; input Integer iY "One-based index of the element that needs to be returned from C to Modelica"; output Real y "Value of the i-th element"; external"C" y = exchangeValues(table, iX, x, iY); end exchangeValues;

Buildings.Fluid.Geothermal.Boreholes.BaseClasses.powerSeries

Power series used to compute far-field temperature

Information

This function computes the power series that is used to compute the far-field temperature.

Inputs

TypeNameDefaultDescription
Realu u
IntegerN Number of coefficients

Outputs

TypeNameDescription
RealWPower series

Modelica definition

function powerSeries "Power series used to compute far-field temperature" input Real u "u"; input Integer N "Number of coefficients"; output Real W "Power series"; algorithm W := -0.5772 - Modelica.Math.log(u) + sum((-1)^(j + 1)*u^j/(j*Buildings.Utilities.Math.Functions.factorial(j)) for j in 1:N); end powerSeries;

Buildings.Fluid.Geothermal.Boreholes.BaseClasses.singleUTubeResistances

Thermal resistances for single U-tube

Information

This model computes the thermal resistances of 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 Hellstroem (1991) using the multipole method.

The figure below shows the thermal network set up by Bauer et al. (2011).

image

The different resistances are calculated as follows. The grout zone and bore hole wall thermal resistance are related as

Rgb1U = ( 1 - x1U ) Rg1U.

The thermal resistance between the two grout zones are

Rgg1U = 2 Rgb1U ( Rar1U - 2 x1U Rg1U ) / ( 2 Rgb1U - Rar1U + 2 x1U Rg1U ).

Thermal resistance between the pipe wall to the capacity in the grout is

RCondGro = x1U Rg1U + log ( ( rTub + eTub ) /rTub ) / ( 2 π hSeg kTub ).

The capacities are located at

x1U =log ( (rBor2 + 2 ( rTub + eTub ) 2)1⁄2/ ( 2 ( rTub + eTub ) ) ) / log ( rBor/ ( √2 ( rTub + eTub ))).

The thermal resistance between the outer borehole wall and one tube is

Rg1U =2 Rb ⁄ hSeg.

The thermal resistance between the two pipe outer walls is

Rar1U =Ra ⁄ hSeg.

The fluid to ground thermal resistance Rb and the grout to grout thermal resistance Ra are calculated with the multipole method (Hellstroem (1991)) as

Rb =1/ ( 4 π kFil ) ( log ( rBor/ ( rTub + eTub ) ) + log ( rBor/ ( 2 xC ) ) + σ log ( rBor4/ ( rBor4 - xC4 ) ) ) - 1/ ( 4 π kFil ) ( ( rTub + eTub ) 2/ ( 4 xC2 ) ( 1 - σ 4 xC4/ ( rBor4 - xC4 ) ) 2 ) / ( ( 1 + β ) / ( 1 - β ) + ( rTub + eTub ) 2/ ( 4 xC2 ) ( 1 + σ 16 xC4 rBor4/ ( rBor4 - xC4 ) 2)).

Ra = 1/ ( π kFil ) ( log ( 2 xC/rTub ) + σ log (( rBor2 + xC2 ) / ( rBor2 - xC2 ) ) ) - 1/ ( π kFil ) ( rTub2/ ( 4 xC2 ) ( 1 + σ 4 rBor4 xC2/ ( rBor4 - xC4 ) ) / ( ( 1 + β ) / ( 1 - β ) - rTub2/ ( 4 xC2 ) + σ 2 rTub2 rBor2 ( rBor4 + xC4 ) / ( rBor4 - xC4 ) 2)),

with σ = ( kFil - kSoi ) / ( kFil + kSoi ) and β = 2 π kFil RCondPipe, where kFil and kSoi are the conductivity of the filling material and of the ground, rTub+eTub and rBor are the pipe and the borehole outside radius and xC is the shank spacing, which is equal to the distance between the center of borehole and the center of the pipe.

Note: The value of Rgg1U may be negative as long as

1/Rgg1U + 1/(2   Rgb1U) > 0,

in which case the laws of thermodynamics are not violated. See Bauer et al. (2011) for details.

References

G. Hellström. Ground heat storage: thermal analyses of duct storage systems (Theory). Dept. of Mathematical Physics, University of Lund, Sweden, 1991.

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.

Inputs

TypeNameDefaultDescription
HeighthSeg Height of the element [m]
RadiusrBor Radius of the borehole [m]
RadiusrTub Radius of the tube [m]
LengtheTub Thickness of the tubes [m]
LengthxC Shank spacing, defined as the distance between the center of a pipe and the center of the borehole [m]
ThermalConductivitykSoi Thermal conductivity of the soi [W/(m.K)]
ThermalConductivitykFil Thermal conductivity of the grout [W/(m.K)]
ThermalConductivitykTub Thermal conductivity of the tube [W/(m.K)]

Outputs

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

Modelica definition

function singleUTubeResistances "Thermal resistances for single U-tube" // Geometry of the borehole input Modelica.SIunits.Height hSeg "Height of the element"; input Modelica.SIunits.Radius rBor "Radius of the borehole"; // Geometry of the pipe input Modelica.SIunits.Radius rTub "Radius of the tube"; input Modelica.SIunits.Length eTub "Thickness of the tubes"; input Modelica.SIunits.Length xC "Shank spacing, defined as the distance between the center of a pipe and the center of the borehole"; // Thermal properties input Modelica.SIunits.ThermalConductivity kSoi "Thermal conductivity of the soi"; input Modelica.SIunits.ThermalConductivity kFil "Thermal conductivity of the grout"; input Modelica.SIunits.ThermalConductivity kTub "Thermal conductivity of the tube"; // Outputs output Modelica.SIunits.ThermalResistance Rgb "Thermal resistance between the grout zone and the borehole wall"; output Modelica.SIunits.ThermalResistance Rgg "Thermal resistance between the two grout zones"; output Modelica.SIunits.ThermalResistance RCondGro "Thermal resistance between the pipe wall ant the capacity in the grout"; output Real x "Capacity location"; protected Boolean test=false "Thermodynamic test for R and x value"; Modelica.SIunits.ThermalResistance Rg "Thermal resistance between outer borehole wall and one tube"; Modelica.SIunits.ThermalResistance Rar "Thermal resistance between the two pipe outer walls"; Modelica.SIunits.ThermalResistance RCondPipe "Thermal resistance of the pipe wall"; Real Rb "Fluid-to-grout resistance, as defined by Hellstroem. Resistance from the fluid in the pipe to the borehole wall"; Real Ra "Grout-to-grout resistance (2D) as defined by Hellstroem. Interaction between the different grout part"; // Help variables Real sigma "Help variable as defined by Hellstroem"; Real beta "Help variable as defined by Hellstroem"; Real R_1delta_LS "One leg of the triangle resistance network, corresponding to the line source solution"; Real R_1delta_MP "One leg of the triangle resistance network, corresponding to the multipole solution"; Real Ra_LS "Grout-to-grout resistance calculated with the line-source approximation"; Integer i=1 "Loop counter"; algorithm // ********** Rb and Ra from multipole ********** // Help variables RCondPipe :=Modelica.Math.log((rTub + eTub)/rTub)/(2*Modelica.Constants.pi*hSeg*kTub); beta :=2*Modelica.Constants.pi*kFil*RCondPipe; sigma :=(kFil - kSoi)/(kFil + kSoi); R_1delta_LS :=1/(2*Modelica.Constants.pi*kFil)*(log(rBor/(rTub + eTub)) + log(rBor/(2*xC)) + sigma*log(rBor^4/(rBor^4 - xC^4))); R_1delta_MP :=R_1delta_LS - 1/(2*Modelica.Constants.pi*kFil)*((rTub + eTub)^2/ (4*xC^2)*(1 - sigma*4*xC^4/(rBor^4 - xC^4))^2)/((1 + beta)/(1 - beta) + ( rTub + eTub)^2/(4*xC^2)*(1 + sigma*16*xC^4*rBor^4/(rBor^4 - xC^4)^2)); Ra_LS :=1/(Modelica.Constants.pi*kFil)*(log(2*xC/rTub) + sigma*log(( rBor^2 + xC^2)/(rBor^2 - xC^2))); //Rb and Ra Rb :=R_1delta_MP/2; Ra :=Ra_LS - 1/(Modelica.Constants.pi*kFil)*(rTub^2/(4*xC^2)*(1 + sigma* 4*rBor^4*xC^2/(rBor^4 - xC^4))/((1 + beta)/(1 - beta) - rTub^2/(4*xC^2) + sigma*2*rTub^2*rBor^2*(rBor^4 + xC^4)/(rBor^4 - xC^4)^2)); //Conversion of Rb (resp. Ra) to Rg (resp. Rar) of Bauer: Rg :=2*Rb/hSeg; Rar :=Ra/hSeg; /* **************** Simplification of Bauer for single U-tube ************************ //Thermal resistance between: Outer wall and one tube Rg := Modelica.Math.acosh((rBor^2 + (rTub + eTub)^2 - xC^2)/(2*rBor*(rTub + eTub)))/(2*Modelica.Constants.pi*hSeg*kFil)*(1.601 - 0.888*xC/rBor); //Thermal resistance between: The two pipe outer walls Rar := Modelica.Math.acosh((2*xC^2 - (rTub + eTub)^2)/(rTub + eTub)^2)/(2* Modelica.Constants.pi*hSeg*kFil); *************************************************************************************** */ // ********** Resistances and capacity location according to Bauer ********** while test == false and i <= 15 loop // Capacity location (with correction factor in case that the test is negative) x := Modelica.Math.log(sqrt(rBor^2 + 2*(rTub + eTub)^2)/(2*(rTub + eTub)))/ Modelica.Math.log(rBor/(sqrt(2)*(rTub + eTub)))*((15 - i + 1)/15); //Thermal resistance between the grout zone and bore hole wall Rgb := (1 - x)*Rg; //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 > 0); i := i + 1; end while; // Check for errors. assert(test, "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.Geothermal.Boreholes.BaseClasses.singleUTubeResistances is hSeg = " + String(hSeg) + " m rBor = " + String(rBor) + " m rTub = " + String(rTub) + " m eTub = " + String(eTub) + " m xC = " + String(xC) + " m kSoi = " + String(kSoi) + " W/m/K kFil = " + String(kFil) + " W/m/K kTub = " + String(kTub) + " W/m/K Computed x = " + String(x) + " K/W Rgb = " + String(Rgb) + " K/W Rgg = " + String(Rgg) + " K/W"); //Conduction resistance in grout from pipe wall to capacity in grout RCondGro := x*Rg + RCondPipe; end singleUTubeResistances;

Buildings.Fluid.Geothermal.Boreholes.BaseClasses.temperatureDrop

Calculate the temperature drop of the soil at the external boundary of the cylinder

Information

This function calculates the temperature drop of the soil at the outer boundary of the cylinder. The analytical formula of Hart and Couvillion (1986) for constant heat extraction is adapted to a non-constant heat flux. To adapt the formula for a variable rate of heat extraction, different constant heat extraction rates, starting at different time instances, are super-imposed. To obtain the temperature drop at the time t=n*Δt, the effects of constant rate of heat extractions are super-imposed as

ΔT ( r , t=n Δt )= 1 ⁄ ( 4 π k ) ∑ W(u(r, t= i Δt)) (qn-i+1-qn-i),

where r is the radius for which the temperature is computed, k is the thermal conductivity of the material, W is a solution of the heat conduction in polar coordinates and qi=Qi/h is the specific rate of heat extraction per unit length at time t=i Δt. The value of W is obtained using

W(u)=[-0.5772 - ln(u) + u - u2/(2   2!) +u3/(3   3!) - u4/(4   4!) + ....].

where u(r,t)= c ρ r2 ⁄ (4 t k) , ρ is the mass density and c is the specific heat capacity per unit mass.

Implementation

The rate of heat flow Qi is obtained from the function Buildings.Fluid.Geothermal.Boreholes.BaseClasses.exchangeValues.

References

Hart and Couvillion, (1986). Earth Coupled Heat Transfer. Publication of the National Water Well Association.

Inputs

TypeNameDefaultDescription
ExtendableArraytable External object that contains the history terms of the heat flux
IntegeriSam Counter for how many time the model was sampled. Defined as iSam=1 when called at t=0
HeatFlowRateQ_flow Heat flow rate to be stored in the external object [W]
TimesamplePeriod Period between two samples [s]
RadiusrExt External radius of the cylinder [m]
HeighthSeg Height of the cylinder [m]
ThermalConductivityk Thermal conductivity of the soil [W/(m.K)]
Densityd Density of the soil [kg/m3]
SpecificHeatCapacityc Specific heat capacity of the soil [J/(kg.K)]

Outputs

TypeNameDescription
TemperatureDifferencedTTemperature drop of the soil [K]

Modelica definition

function temperatureDrop "Calculate the temperature drop of the soil at the external boundary of the cylinder" input ExtendableArray table "External object that contains the history terms of the heat flux"; input Integer iSam(min=1) "Counter for how many time the model was sampled. Defined as iSam=1 when called at t=0"; input Modelica.SIunits.HeatFlowRate Q_flow "Heat flow rate to be stored in the external object"; input Modelica.SIunits.Time samplePeriod "Period between two samples"; input Modelica.SIunits.Radius rExt "External radius of the cylinder"; input Modelica.SIunits.Height hSeg "Height of the cylinder"; input Modelica.SIunits.ThermalConductivity k "Thermal conductivity of the soil"; input Modelica.SIunits.Density d "Density of the soil"; input Modelica.SIunits.SpecificHeatCapacity c "Specific heat capacity of the soil"; output Modelica.SIunits.TemperatureDifference dT "Temperature drop of the soil"; protected Modelica.SIunits.Time minSamplePeriod= rExt^2/(4*(k/c/d)*3.8) "Minimal length of the sampling period"; Modelica.SIunits.HeatFlowRate QL_flow "Intermediate variable for heat flow rate at the lower bound of the time interval"; Modelica.SIunits.HeatFlowRate QU_flow "Intermediate variable for heat flow rate at the upper bound of the time interval"; algorithm assert(rExt*rExt/(4*(k/c/d)*samplePeriod)<=3.8, "The samplePeriod has to be bigger than " + String(minSamplePeriod) + " for convergence purpose. samplePeriod = "+ String(samplePeriod)); if iSam == 1 then // First call, at t=0 dT := 0; QL_flow := Buildings.Fluid.Geothermal.Boreholes.BaseClasses.exchangeValues( table=table, iX=iSam, x=Q_flow, iY=iSam); else dT := 0; // The first evaluation is at iSam=2, in which case we have one term of the sum, // and t=samplePeriod=(iSam-1)*samplePeriod for i in 1:(iSam-1) loop QL_flow := Buildings.Fluid.Geothermal.Boreholes.BaseClasses.exchangeValues( table=table, iX=iSam, x=Q_flow, iY=iSam+1-i); QU_flow := Buildings.Fluid.Geothermal.Boreholes.BaseClasses.exchangeValues( table=table, iX=iSam, x=Q_flow, iY=iSam-i); // The division by hSeg is because QU_flow and QL_flow are in [W], but the equation // requires [W/m], i.e., heat flow rate per unit length of the line source. dT := dT + 1/(4*Modelica.Constants.pi*k)* Buildings.Fluid.Geothermal.Boreholes.BaseClasses.powerSeries( u=c*d/(4*k*i*samplePeriod)*rExt^2, N=10)* (QL_flow-QU_flow)/hSeg; end for; end if; end temperatureDrop;