da_cloud_model.inc
References to this file elsewhere.
1 subroutine da_cloud_model (TB, PB, QB, QCWB, QRNB, ZB, ZFB, DT, kts, kte)
2
3 !-----------------------------------------------------------------
4 ! Purpose: Calculate DT (=dz/w) using cumulus parameterization
5 ! of a one-dimensional cloud model.
6 !-----------------------------------------------------------------
7
8 ! Calculate DT
9
10 implicit none
11
12 integer, intent(in) :: kts, kte
13 real, intent(in), dimension(kts:kte) :: TB, PB, QB, QCWB, QRNB, ZB
14 real, intent(in), dimension(kts:kte+1) :: ZFB
15 real, intent(out), dimension(kts:kte) :: DT
16
17 integer :: k
18 real :: P0, Z0, T0, Q0
19 real :: PLCL, ZLCL, TLCL, QLCL
20 integer :: KCB, KCT
21 real :: PCT, ZCT
22 real, dimension(kts:kte) :: ZC, TC, QC, PP, QT
23 real, dimension(kts:kte) :: TCV, TBV, B
24 real :: G, ALPHA, RC, MU, XX, YY
25 real, dimension(kts:kte+1) :: W0, W
26
27 if (trace_use) call da_trace_entry("da_cloud_model")
28
29 G=9.81
30 ALPHA=0.5
31 RC=100.0
32 MU=0.183/RC
33
34 do k = kts, kte+1
35 W0(k)=0.0
36 W(k)=0.0
37 end do
38
39 do k = kts, kte
40 PP(k)=PB(k)/100.0
41 DT(k)=0.0
42 end do
43
44 P0 = PP(kts)
45 Z0 = ZB(kts)
46 T0 = MAX(TB(kts),303.)
47
48 call da_qfrmrh (P0, T0, 95., Q0)
49
50 call da_lcl (P0, Z0, T0, Q0, PLCL, ZLCL, TLCL, QLCL)
51
52 call da_qfrmrh (PLCL, TLCL, 95., QLCL)
53
54 call da_cumulus (ZLCL, TLCL, QLCL, PLCL, PP, TB, &
55 ZC, TC, QC, KCB, KCT, PCT, ZCT, kts, kte)
56
57 do k = KCB, KCT
58 TCV(k) = TC(k) * (1. + 0.608 * QC(k))
59 TBV(k) = TB(k) * (1. + 0.608 * QB(k))
60
61 B(k) = (TCV(k)-TBV(k)) / TBV(k)
62
63 QT(k) = QC(k) + QCWB(k) + QRNB(k)
64 end do
65
66 W0(KCB) = 0.0
67 do k = KCB+1, KCT+1
68 XX = 1.0+2.0*MU*(ZFB(k)-ZFB(k-1))
69 YY = 2.0*G*(B(k-1)/(1.+ALPHA) - QT(k-1)) * (ZFB(k)-ZFB(k-1))
70 W0(k) = (W0(k-1)+YY) / XX
71 end do
72
73 do k = KCB, KCT+1
74 if (W0(k) >= 0.) then
75 W(k) = sqrt(W0(k))
76 end if
77 end do
78
79
80 do k = KCT, KCB+1, -1
81 if (W(k) >= 0.01) then
82 DT(k) = (ZB(k)-ZB(k-1))/W(k)
83 else
84 DT(k) = 0.0
85 end if
86 end do
87
88 if (trace_use) call da_trace_exit("da_cloud_model")
89
90 end subroutine da_cloud_model
91
92