Nodes of different colours represent the following:
Solid arrows point from a parent (sub)module to the submodule which is descended from it. Dashed arrows point from a module being used to the module or program unit using it. Where possible, edges connecting nodes are given different colours to make them easier to distinguish in large graphs.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | ZoneNum |
Nodes of different colours represent the following:
Solid arrows point from a procedure to one which it calls. Dashed arrows point from an interface to procedures which implement that interface. This could include the module procedures in a generic interface or the implementation in a submodule of an interface in a parent module. Where possible, edges connecting nodes are given different colours to make them easier to distinguish in large graphs.
Nodes of different colours represent the following:
Solid arrows point from a procedure to one which it calls. Dashed arrows point from an interface to procedures which implement that interface. This could include the module procedures in a generic interface or the implementation in a submodule of an interface in a parent module. Where possible, edges connecting nodes are given different colours to make them easier to distinguish in large graphs.
SUBROUTINE CalcUCSDUI(ZoneNum)
! SUBROUTINE INFORMATION:
! AUTHOR Fred Buhl
! DATE WRITTEN August 2005
! MODIFIED Brent Griffith June 2008 for new interpolation and time history
! RE-ENGINEERED na
! PURPOSE OF THIS SUBROUTINE:
! Using the UCSD UFAD interior zone model, this subroutine calculates the occupied subzone height,
! surface heat transfer coefficients, the occupied subzone temperature, and the upper subzone temperature.
! METHODOLOGY EMPLOYED:
! The zone is divided into 2 subzones with a variable transition height.
! REFERENCES:
! The model is described in the EnergyPlus Engineering Reference in Anna Liu's UCSD PhD thesis.
! USE STATEMENTS:
USE ScheduleManager, ONLY: GetCurrentScheduleValue
USE DataZoneEquipment, ONLY: ZoneEquipConfig
USE Psychrometrics, ONLY: PsyRhoAirFnPbTdbW, PsyCpAirFnWTdb
USE DataHeatBalFanSys
USE DataHVACGlobals, ONLY: TimestepSys, UseZoneTimeStepHistory
USE InternalHeatGains, ONLY: SumInternalConvectionGainsByTypes, SumReturnAirConvectionGainsByTypes
IMPLICIT NONE ! Enforce explicit typing of all variables in this routine
! SUBROUTINE ARGUMENT DEFINITIONS:
INTEGER, INTENT(IN) :: ZoneNum ! index number for the specified zone
! SUBROUTINE PARAMETER DEFINITIONS:
! na
! INTERFACE BLOCK SPECIFICATIONS:
! na
! DERIVED TYPE DEFINITIONS:
! na
! SUBROUTINE LOCAL VARIABLE DECLARATIONS:
LOGICAL :: MIXFLAG = .FALSE. ! if true treat as a mixed zone
REAL(r64) :: CeilingHeight ! zone ceiling height above floor [m]
INTEGER :: UINum ! index to underfloor interior zone model data
REAL(r64) :: GainsFrac ! fraction of occupied subzone heat gains that remain in the subzone;
! that is, don't go into the plumes
! REAL(r64) :: NumPLPP ! number of plumes per person
REAL(r64) :: HeightThermostat ! height of the thermostat above the floor [m]
REAL(r64) :: HeightComfort ! height at which comfort temperature is calculated
REAL(r64) :: TempDiffCritRep ! Minimum temperature difference between upper and occupied subzones for reporting
REAL(r64) :: ConvGainsOccSubzone ! convective heat gains into the lower (occupied) subzone [W]
REAL(r64) :: ConvGainsUpSubzone ! convective heat gains into the upper subzone [W]
REAL(r64) :: ConvGains ! total zone convective gains (exclusing surfaces) [W]
INTEGER :: ZoneEquipConfigNum ! ZoneEquipConfig index for this UFAD zone
REAL(r64) :: SumSysMCp ! Sum of system mass flow rate * specific heat for this zone [W/K]
REAL(r64) :: SumSysMCpT ! Sum of system mass flow rate * specific heat * temperature for this zone [W]
REAL(r64) :: SumSysM ! Sum of systems mass flow rate [kg/s]
REAL(r64) :: NodeTemp ! inlet node temperature [K]
REAL(r64) :: MassFlowRate ! system mass flow rate [kg/s]
REAL(r64) :: CpAir ! specific heat of air [J/kgK]
INTEGER :: InNodeIndex ! inlet node index in ZoneEquipConfig
REAL(r64) :: SumMCp ! mass flow rate * specific heat for this zone for infiltration, ventilation, mixing [W/K]
REAL(r64) :: SumMCpT ! mass flow rate * specific heat* temp for this zone for infiltration, ventilation, mixing [W]
REAL(r64) :: MCP_Total ! total mass flow rate * specific heat for this zone [W/K]
REAL(r64) :: MCpT_Total ! total mass flow rate * specific heat* temp for this zone [W]
REAL(r64) :: NumberOfPlumes
REAL(r64) :: PowerInPlumes ! [W]
REAL(r64) :: PowerPerPlume=0.0d0 ! power generating each plume [W]
REAL(r64) :: HeightFrac ! Fractional height of transition between occupied and upper subzones
REAL(r64) :: TotSysFlow ! [m3/s]
REAL(r64) :: NumDiffusersPerPlume
REAL(r64) :: NumDiffusers
REAL(r64) :: TSupK ! supply yemperature [K]
REAL(r64) :: Gamma ! dimensionless height parameter; higher gamma means interface height will be
! higher, smaller gamma means interface height will be lower.
REAL(r64) :: DiffArea ! diffuser effective area [m2]
REAL(r64) :: ThrowAngle ! diffuser slot angle relative to vertical [radians]
REAL(r64) :: SourceHeight ! height of plume sources above the floor [m]
INTEGER :: Ctd
REAL(r64) :: AirCap
REAL(r64) :: TempHistTerm
REAL(r64) :: ZTAveraged
REAL(r64) :: HeightUpSubzoneAve ! Height of center of upper air subzone
REAL(r64) :: HeightOccupiedSubzoneAve ! Height of center of occupied air subzone
REAL(r64) :: ZoneMult ! total zone multiplier
INTEGER :: ZoneNodeNum ! node number of the HVAC zone node
REAL(r64) :: TempDepCoef=0.0d0 ! Formerly CoefSumha, coef in zone temp equation with dimensions of h*A
REAL(r64) :: TempIndCoef=0.0d0 ! Formerly CoefSumhat, coef in zone temp equation with dimensions of h*A(T1
INTEGER, DIMENSION(28) :: IntGainTypesOccupied = (/IntGainTypeOf_People, &
IntGainTypeOf_WaterHeaterMixed, &
IntGainTypeOf_WaterHeaterStratified, &
IntGainTypeOf_ThermalStorageChilledWaterMixed, &
IntGainTypeOf_ThermalStorageChilledWaterStratified, &
IntGainTypeOf_ElectricEquipment, &
IntGainTypeOf_GasEquipment, &
IntGainTypeOf_HotWaterEquipment, &
IntGainTypeOf_SteamEquipment, &
IntGainTypeOf_OtherEquipment, &
IntGainTypeOf_ZoneBaseboardOutdoorTemperatureControlled, &
IntGainTypeOf_GeneratorFuelCell, &
IntGainTypeOf_WaterUseEquipment, &
IntGainTypeOf_GeneratorMicroCHP, &
IntGainTypeOf_ElectricLoadCenterTransformer, &
IntGainTypeOf_ElectricLoadCenterInverterSimple, &
IntGainTypeOf_ElectricLoadCenterInverterFunctionOfPower, &
IntGainTypeOf_ElectricLoadCenterInverterLookUpTable, &
IntGainTypeOf_ElectricLoadCenterStorageBattery, &
IntGainTypeOf_ElectricLoadCenterStorageSimple, &
IntGainTypeOf_PipeIndoor, &
IntGainTypeOf_RefrigerationCase, &
IntGainTypeOf_RefrigerationCompressorRack, &
IntGainTypeOf_RefrigerationSystemAirCooledCondenser ,&
IntGainTypeOf_RefrigerationSystemSuctionPipe, &
IntGainTypeOf_RefrigerationSecondaryReceiver, &
IntGainTypeOf_RefrigerationSecondaryPipe, &
IntGainTypeOf_RefrigerationWalkIn/)
INTEGER, DIMENSION(2) :: IntGainTypesUpSubzone = (/IntGainTypeOf_DaylightingDeviceTubular , &
IntGainTypeOf_Lights/)
REAL(r64) :: RetAirGains
! Exact solution or Euler method
If (ZoneAirSolutionAlgo .NE. Use3rdOrder) Then
If (ShortenTimeStepSysRoomAir .and. TimeStepSys .LT. TimeStepZone) Then
If (PreviousTimeStep < TimeStepZone) Then
Zone1OC(ZoneNum) = ZoneM2OC(ZoneNum)
Zone1MX(ZoneNum) = ZoneM2MX(ZoneNum)
Else
Zone1OC(ZoneNum) = ZoneMXOC(ZoneNum)
Zone1MX(ZoneNum) = ZoneMXMX(ZoneNum)
End If
Else
Zone1OC(ZoneNum) = ZTOC(ZoneNum)
Zone1MX(ZoneNum) = ZTMX(ZoneNum)
End If
End If
MIXFLAG = .FALSE.
UFHcIn = HConvIn
SumSysMCp = 0.0d0
SumSysMCpT = 0.0d0
TotSysFlow = 0.0d0
TSupK = 0.0d0
SumSysM = 0.0d0
ZoneMult = Zone(ZoneNum)%Multiplier * Zone(ZoneNum)%ListMultiplier
CeilingHeight = ZoneCeilingHeight((ZoneNum-1)*2 + 2) - ZoneCeilingHeight((ZoneNum-1)*2 + 1)
UINum = ZoneUFPtr(ZoneNum)
HeightThermostat = ZoneUCSDUI(UINum)%ThermostatHeight
HeightComfort = ZoneUCSDUI(UINum)%ComfortHeight
TempDiffCritRep = ZoneUCSDUI(UINum)%TempTrigger
DiffArea = ZoneUCSDUI(UINum)%DiffArea
ThrowAngle = DegToRadians*ZoneUCSDUI(UINum)%DiffAngle
SourceHeight = 0.0d0
NumDiffusers = ZoneUCSDUI(UINum)%DiffusersPerZone
PowerPerPlume = ZoneUCSDUI(UINum)%PowerPerPlume
! gains from occupants, task lighting, elec equip, gas equip, other equip, hot water equip, steam equip,
! baseboards (nonthermostatic), water heater skin loss
CALL SumInternalConvectionGainsByTypes(ZoneNum, IntGainTypesOccupied, ConvGainsOccSubzone)
! Add heat to return air if zonal system (no return air) or cycling system (return air frequently very
! low or zero)
IF (Zone(ZoneNum)%NoHeatToReturnAir) THEN
CALL SumReturnAirConvectionGainsByTypes(ZoneNum, IntGainTypesOccupied, RetAirGains)
ConvGainsOccSubzone = ConvGainsOccSubzone + RetAirGains
END IF
! gains from lights (ceiling), tubular daylighting devices, high temp radiant heaters
CALL SumInternalConvectionGainsByTypes(ZoneNum, IntGainTypesUpSubzone, ConvGainsUpSubzone)
ConvGainsUpSubzone = ConvGainsUpSubzone + SumConvHTRadSys(ZoneNum)
IF (Zone(ZoneNum)%NoHeatToReturnAir) THEN
CALL SumReturnAirConvectionGainsByTypes(ZoneNum, IntGainTypesUpSubzone, RetAirGains)
ConvGainsUpSubzone = ConvGainsUpSubzone + RetAirGains
END IF
ConvGains = ConvGainsOccSubzone + ConvGainsUpSubzone + SysDepZoneLoadsLagged(ZoneNum)
ZoneEquipConfigNum = ZoneUCSDUI(UINum)%ZoneEquipPtr
IF (ZoneEquipConfigNum > 0) THEN
DO InNodeIndex = 1,ZoneEquipConfig(ZoneEquipConfigNum)%NumInletNodes
NodeTemp = Node(ZoneEquipConfig(ZoneEquipConfigNum)%InletNode(InNodeIndex))%Temp
MassFlowRate = Node(ZoneEquipConfig(ZoneEquipConfigNum)%InletNode(InNodeIndex))%MassFlowRate
CpAir = PsyCpAirFnWTdb(ZoneAirHumRat(ZoneNum), NodeTemp)
SumSysMCp = SumSysMCp + MassFlowRate * CpAir
SumSysMCpT = SumSysMCpT + MassFlowRate * CpAir * NodeTemp
TotSysFlow = TotSysFlow + MassFlowRate / PsyRhoAirFnPbTdbW(OutBaroPress,NodeTemp,ZoneAirHumRat(ZoneNum))
TSupK = TSupK + MassFlowRate * NodeTemp
SumSysM = SumSysM + MassFlowRate
END DO
IF (TotSysFlow > 0.0d0) THEN
TSupK = TSupK/SumSysM + KelvinConv
ELSE
TSupK = 0.0d0
END IF
END IF
! mass flow times specific heat for infiltration, ventilation, mixing, earth tube
SumMCp = MCPI(ZoneNum) + MCPV(ZoneNum) + MCPM(ZoneNum) + MCPE(ZoneNum) + MCPC(ZoneNum) + MdotCPOA(ZoneNum)
! mass flow times specific heat times temperature for infiltration, ventilation, mixing, earth tube
SumMCpT = MCPTI(ZoneNum) + MCPTV(ZoneNum) + MCPTM(ZoneNum) + MCPTE(ZoneNum) + MCPTC(ZoneNum) + &
MdotCPOA(ZoneNum)*Zone(ZoneNum)%OutDryBulbTemp
MCP_Total = SumMCp + SumSysMCp
MCpT_Total = SumMCpT + SumSysMCpT
! For the York MIT diffusers (variable area) the area varies with the flow rate. Assume 400 ft/min velocity
! at the diffuser, and a design flow rate of 150 cfm (.0708 m3/s). Then the design area for each diffuser is
! 150 ft3/min / 400 ft/min = .375 ft2 = .035 m2. This is adjusted each time step by
! (TotSysFlow/(NumDiffusers*.0708))*.035
IF (ZoneUCSDUI(UINum)%DiffuserType .EQ. VarArea) THEN
DiffArea = .035d0*TotSysFlow/(.0708d0*NumDiffusers)
END IF
! initial estimate of convective transfer from surfaces; assume HeightFrac is 0.5.
CALL HcUCSDUF(ZoneNum,0.5d0)
PowerInPlumes = ConvGains + HAT_OC - HA_OC*ZTOC(ZoneNum) + HAT_MX - HA_MX*ZTMX(ZoneNum)
IF (PowerPerPlume > 0.0d0 .AND. PowerInPlumes > 0.0d0) THEN
NumberOfPlumes = PowerInPlumes / PowerPerPlume
NumDiffusersPerPlume = NumDiffusers / NumberOfPlumes
ELSE
NumberOfPlumes = 1.0d0
NumDiffusersPerPlume = 1.0d0
END IF
IF ((PowerInPlumes <= 0.0d0) .OR. (TotSysFlow .EQ. 0.0d0) .OR. (TsupK-KelvinConv) > MAT(ZoneNum)) THEN
! The system will mix
HeightFrac = 0.0d0
ELSE
Gamma = (TotSysFlow*COS(ThrowAngle))**1.5d0 / (NumberOfPlumes*(NumDiffusersPerPlume*DiffArea)**1.25d0 &
* (0.0281d0*0.001d0*PowerInPlumes)**0.5d0)
IF (ZoneUCSDUI(UINum)%CalcTransHeight) THEN
HeightFrac = ((NumDiffusersPerPlume*DiffArea)**0.5d0 * (7.43d0*log(Gamma) - 1.35d0) + 0.5d0*SourceHeight) / CeilingHeight
ELSE
HeightFrac = ZoneUCSDUI(UINum)%TransHeight / CeilingHeight
END IF
HeightFrac = MAX(0.0d0,MIN(1.0d0,HeightFrac))
DO Ctd = 1,4
CALL HcUCSDUF(ZoneNum,HeightFrac)
PowerInPlumes = ConvGains + HAT_OC - HA_OC*ZTOC(ZoneNum) + HAT_MX - HA_MX*ZTMX(ZoneNum)
IF (PowerPerPlume > 0.0d0 .AND. PowerInPlumes > 0.0d0) THEN
NumberOfPlumes = PowerInPlumes / PowerPerPlume
NumDiffusersPerPlume = NumDiffusers / NumberOfPlumes
ELSE
NumberOfPlumes = 1.0d0
NumDiffusersPerPlume = 1.0d0
END IF
IF (PowerInPlumes .LE. 0.0d0) EXIT
Gamma = (TotSysFlow*COS(ThrowAngle))**1.5d0 / (NumberOfPlumes*(NumDiffusersPerPlume*DiffArea)**1.25d0 &
* (0.0281d0*0.001d0*PowerInPlumes)**0.5d0)
IF (ZoneUCSDUI(UINum)%CalcTransHeight) THEN
HeightFrac = ((NumDiffusersPerPlume*DiffArea)**0.5d0 * (7.43d0*log(Gamma) - 1.35d0) + 0.5d0*SourceHeight) / CeilingHeight
ELSE
HeightFrac = ZoneUCSDUI(UINum)%TransHeight / CeilingHeight
END IF
HeightFrac = MAX(0.0d0,MIN(1.0d0,HeightFrac))
HeightTransition(ZoneNum) = HeightFrac * CeilingHeight
GainsFrac = ZoneUCSDUI(UINum)%A_Kc * Gamma**ZoneUCSDUI(UINum)%B_Kc + ZoneUCSDUI(UINum)%C_Kc + &
ZoneUCSDUI(UINum)%D_Kc * Gamma + ZoneUCSDUI(UINum)%E_Kc * Gamma**2
GainsFrac = MAX(0.6d0,MIN(GainsFrac,1.0d0))
AIRRATOC(ZoneNum) = Zone(ZoneNum)%Volume*(HeightTransition(ZoneNum)-MIN(HeightTransition(ZoneNum),0.2d0)) &
/CeilingHeight*ZoneVolCapMultpSens &
*PsyRhoAirFnPbTdbW(OutBaroPress,MATOC(ZoneNum),ZoneAirHumRat(ZoneNum)) &
*PsyCpAirFnWTdb(ZoneAirHumRat(ZoneNum),MATOC(ZoneNum))/(TimeStepSys*SecInHour)
AIRRATMX(ZoneNum) = Zone(ZoneNum)%Volume*(CeilingHeight-HeightTransition(ZoneNum)) &
/CeilingHeight*ZoneVolCapMultpSens &
*PsyRhoAirFnPbTdbW(OutBaroPress,MATMX(ZoneNum),ZoneAirHumRat(ZoneNum)) &
*PsyCpAirFnWTdb(ZoneAirHumRat(ZoneNum),MATMX(ZoneNum))/(TimeStepSys*SecInHour)
IF (UseZoneTimeStepHistory) THEN
ZTM3OC(ZoneNum) = XM3TOC(ZoneNum)
ZTM2OC(ZoneNum) = XM2TOC(ZoneNum)
ZTM1OC(ZoneNum) = XMATOC(ZoneNum)
ZTM3MX(ZoneNum) = XM3TMX(ZoneNum)
ZTM2MX(ZoneNum) = XM2TMX(ZoneNum)
ZTM1MX(ZoneNum) = XMATMX(ZoneNum)
ELSE
ZTM3OC(ZoneNum) = DSXM3TOC(ZoneNum)
ZTM2OC(ZoneNum) = DSXM2TOC(ZoneNum)
ZTM1OC(ZoneNum) = DSXMATOC(ZoneNum)
ZTM3MX(ZoneNum) = DSXM3TMX(ZoneNum)
ZTM2MX(ZoneNum) = DSXM2TMX(ZoneNum)
ZTM1MX(ZoneNum) = DSXMATMX(ZoneNum)
ENDIF
AirCap = AIRRATOC(ZoneNum)
TempHistTerm = AirCap*(3.0d0*ZTM1OC(ZoneNum)-(3.0d0/2.0d0)*ZTM2OC(ZoneNum)+(1.0d0/3.0d0)*ZTM3OC(ZoneNum))
TempDepCoef = GainsFrac*HA_OC + MCP_Total
TempIndCoef = GainsFrac*(ConvGains+HAT_OC+HAT_MX-HA_MX*ZTMX(ZoneNum))+MCPT_Total+NonAirSystemResponse(ZoneNum)/ZoneMult
SELECT CASE (ZoneAirSolutionAlgo)
CASE (Use3rdOrder)
ZTOC(ZoneNum) = (TempHistTerm + GainsFrac*(ConvGains + HAT_OC + HAT_MX - HA_MX*ZTMX(ZoneNum)) &
+ MCPT_Total + NonAirSystemResponse(ZoneNum)/ZoneMult) &
/ ((11.0d0/6.0d0)*AirCap + GainsFrac*HA_OC + MCP_Total)
CASE (UseAnalyticalSolution)
If (TempDepCoef .eq. 0.0d0) Then ! B=0
ZTOC(ZoneNum) = Zone1OC(ZoneNum) + TempIndCoef/AirCap
Else
ZTOC(ZoneNum) = (Zone1OC(ZoneNum)-TempIndCoef/TempDepCoef)*exp(MIN(700.d0,-TempDepCoef/AirCap))+TempIndCoef/TempDepCoef
End If
CASE (UseEulerMethod)
ZTOC(ZoneNum) = (AirCap*Zone1OC(ZoneNum)+TempIndCoef)/(AirCap+TempDepCoef)
END SELECT
AirCap = AIRRATMX(ZoneNum)
TempHistTerm = AirCap*(3.0d0*ZTM1MX(ZoneNum)-(3.0d0/2.0d0)*ZTM2MX(ZoneNum)+(1.0d0/3.0d0)*ZTM3MX(ZoneNum))
TempDepCoef = (1.0d0-GainsFrac)*HA_MX + MCP_Total
TempIndCoef = (1.0d0-GainsFrac)*(ConvGains + HAT_OC + HAT_MX - HA_OC*ZTOC(ZoneNum)) + ZTOC(ZoneNum)*MCP_Total
SELECT CASE (ZoneAirSolutionAlgo)
CASE (Use3rdOrder)
ZTMX(ZoneNum) = (TempHistTerm + (1.0d0-GainsFrac)*(ConvGains + HAT_OC + HAT_MX - HA_OC*ZTOC(ZoneNum)) + &
ZTOC(ZoneNum)*MCP_Total) / ((11.0d0/6.0d0)*AirCap + (1.0d0-GainsFrac)*HA_MX + MCP_Total)
CASE (UseAnalyticalSolution)
If (TempDepCoef .eq. 0.0d0) Then ! B=0
ZTMX(ZoneNum) = Zone1MX(ZoneNum) + TempIndCoef/AirCap
Else
ZTMX(ZoneNum) = (Zone1MX(ZoneNum)-TempIndCoef/TempDepCoef)*exp(MIN(700.d0,-TempDepCoef/AirCap))+TempIndCoef/TempDepCoef
End If
CASE (UseEulerMethod)
ZTMX(ZoneNum) = (AirCap*Zone1MX(ZoneNum)+TempIndCoef)/(AirCap+TempDepCoef)
END SELECT
ZTFLOOR(ZoneNum) = ZTOC(ZoneNum)
END DO
IF (PowerInPlumes .LE. 0.0d0) THEN
HeightFrac = 0.0d0
AirModel(ZoneNum)%SimAirModel = .FALSE.
ZoneUFGamma(ZoneNum) = 0.0d0
ZoneUFPowInPlumes(ZoneNum) = 0.0d0
ELSE
AirModel(ZoneNum)%SimAirModel = .TRUE.
ZoneUFGamma(ZoneNum) = Gamma
ZoneUFPowInPlumes(ZoneNum) = PowerInPlumes
END IF
END IF
!=============================== M I X E D Calculation ==============================================
IF(ZTMX(ZoneNum) < ZTOC(ZoneNum) .or. MCP_Total <= 0.0d0 .or. &
HeightFrac*CeilingHeight < ThickOccupiedSubzoneMin) THEN
MIXFLAG = .TRUE.
HeightFrac = 0.0d0
AvgTempGrad(ZoneNum) = 0.0d0
MaxTempGrad(ZoneNum) = 0.0d0
AirModel(ZoneNum)%SimAirModel = .FALSE.
AirCap = AIRRAT(ZoneNum)
TempHistTerm = AirCap*(3.0d0*ZTM1(ZoneNum)-(3.0d0/2.0d0)*ZTM2(ZoneNum)+(1.0d0/3.0d0)*ZTM3(ZoneNum))
DO Ctd = 1,3
TempDepCoef = HA_MX + HA_OC + MCP_Total
TempIndCoef = ConvGains + HAT_MX + HAT_OC + MCpT_Total
SELECT CASE (ZoneAirSolutionAlgo)
CASE (Use3rdOrder)
ZTAveraged = (TempHistTerm + ConvGains + HAT_MX + HAT_OC + MCpT_Total)/ &
((11.0d0/6.0d0)*AirCap + HA_MX + HA_OC + MCP_Total)
CASE (UseAnalyticalSolution)
If (TempDepCoef .eq. 0.0d0) Then ! B=0
ZTAveraged = ZoneT1(ZoneNum) + TempIndCoef/AirCap
Else
ZTAveraged = (ZoneT1(ZoneNum)-TempIndCoef/TempDepCoef)*exp(MIN(700.d0,-TempDepCoef/AirCap))+TempIndCoef/TempDepCoef
End If
CASE (UseEulerMethod)
ZTAveraged = (AirCap*ZoneT1(ZoneNum)+TempIndCoef)/(AirCap+TempDepCoef)
END SELECT
ZTOC(ZoneNum) = ZTAveraged
ZTMX(ZoneNum) = ZTAveraged
ZTFLOOR(ZoneNum) = ZTAveraged
CALL HcUCSDUF(ZoneNum,HeightFrac)
TempDepCoef = HA_MX + HA_OC + MCP_Total
TempIndCoef = ConvGains + HAT_MX + HAT_OC + MCpT_Total
SELECT CASE (ZoneAirSolutionAlgo)
CASE (Use3rdOrder)
ZTAveraged = (TempHistTerm + ConvGains + HAT_MX + HAT_OC + MCpT_Total)/ &
((11.0d0/6.0d0)*AirCap + HA_MX + HA_OC + MCP_Total)
CASE (UseAnalyticalSolution)
If (TempDepCoef .eq. 0.0d0) Then ! B=0
ZTAveraged = ZoneT1(ZoneNum) + TempIndCoef/AirCap
Else
ZTAveraged = (ZoneT1(ZoneNum)-TempIndCoef/TempDepCoef)*exp(MIN(700.d0,-TempDepCoef/AirCap))+TempIndCoef/TempDepCoef
End If
CASE (UseEulerMethod)
ZTAveraged = (AirCap*ZoneT1(ZoneNum)+TempIndCoef)/(AirCap+TempDepCoef)
END SELECT
ZTOC(ZoneNum) = ZTAveraged
ZTMX(ZoneNum) = ZTAveraged
ZTFLOOR(ZoneNum) = ZTAveraged
END DO
END IF
!=========================================================================================
! Comfort temperature and temperature at the thermostat/temperature control sensor
HeightTransition(ZoneNum) = HeightFrac * CeilingHeight
HeightUpSubzoneAve = (CeilingHeight + HeightTransition(ZoneNum)) / 2.d0
HeightOccupiedSubzoneAve = HeightTransition(ZoneNum) / 2.d0
! Comfort temperature
IF (MIXFLAG) THEN
TCMF(ZoneNum) = ZTAveraged
ELSE
IF (HeightComfort < HeightOccupiedSubzoneAve) THEN
TCMF(ZoneNum) = ZTOC(ZoneNum)
ELSEIF (HeightComfort >= HeightOccupiedSubzoneAve .AND. HeightComfort < HeightUpSubzoneAve) THEN
TCMF(ZoneNum) = (ZTOC(ZoneNum) * (HeightUpSubzoneAve - HeightComfort) &
+ ZTMX(ZoneNum) * (HeightComfort - HeightOccupiedSubzoneAve)) &
/ (HeightUpSubzoneAve - HeightOccupiedSubzoneAve)
ELSEIF (HeightComfort >= HeightUpSubzoneAve .AND. HeightComfort <= CeilingHeight) THEN
TCMF(ZoneNum) = ZTMX(ZoneNum)
ELSE
CALL ShowFatalError ('UFAD comfort height is above ceiling or below floor in Zone: '// &
TRIM(Zone(ZoneNum)%Name))
ENDIF
ENDIF
! Temperature at the thermostat/temperature control sensor
IF (MIXFLAG) THEN
TempTstatAir(ZoneNum) = ZTAveraged
ELSE
IF (HeightThermostat < HeightOccupiedSubzoneAve) THEN
TempTstatAir(ZoneNum) = ZTOC(ZoneNum)
ELSEIF (HeightThermostat >= HeightOccupiedSubzoneAve .AND. HeightThermostat < HeightUpSubzoneAve) THEN
TempTstatAir(ZoneNum) = (ZTOC(ZoneNum) * (HeightUpSubzoneAve - HeightThermostat) &
+ ZTMX(ZoneNum) * (HeightThermostat - HeightOccupiedSubzoneAve)) &
/ (HeightUpSubzoneAve - HeightOccupiedSubzoneAve)
ELSEIF (HeightThermostat >= HeightUpSubzoneAve .AND. HeightThermostat <= CeilingHeight) THEN
TempTstatAir(ZoneNum) = ZTMX(ZoneNum)
ELSE
CALL ShowFatalError ('Underfloor air distribution thermostat height is above ceiling or below floor in Zone: '// &
TRIM(Zone(ZoneNum)%Name))
ENDIF
ENDIF
! Temperature gradients
IF ((HeightUpSubzoneAve - HeightOccupiedSubzoneAve) > 0.1d0) THEN
AvgTempGrad(ZoneNum) = (ZTMX(ZoneNum)-ZTOC(ZoneNum))/(HeightUpSubzoneAve - HeightOccupiedSubzoneAve)
ELSE
AvgTempGrad(ZoneNum) = 0.0d0
ENDIF
IF (MIXFLAG) THEN
ZoneUFMixedFlag(ZoneNum) = 1
AirModel(ZoneNum)%SimAirModel = .FALSE.
ELSE
ZoneUFMixedFlag(ZoneNum) = 0
AirModel(ZoneNum)%SimAirModel = .TRUE.
ENDIF
IF (ZoneEquipConfigNum > 0) THEN
ZoneNodeNum = Zone(ZoneNum)%SystemZoneNodeNumber
Node(ZoneNodeNum)%Temp = ZTMX(ZoneNum)
ENDIF
IF (MIXFLAG) THEN
Phi(ZoneNum) = 1.0d0
ELSE
Phi(ZoneNum) = (ZTOC(ZoneNum) - (TSupK-KelvinConv)) / (ZTMX(ZoneNum) - (TSupK-KelvinConv))
END IF
! Mixed for reporting purposes
IF ((MIXFLAG) .OR. ((ZTMX(ZoneNum)-ZTOC(ZoneNum)).LT.TempDiffCritRep)) THEN
ZoneUFMixedFlagRep(ZoneNum) = 1.d0
HeightTransition(ZoneNum) = 0.0d0
AvgTempGrad(ZoneNum) = 0.0d0
ELSE
ZoneUFMixedFlagRep(ZoneNum) = 0.0d0
ENDIF
RETURN
END SUBROUTINE CalcUCSDUI