updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_setup_structures / da_cloud_model.inc
blobc10ce2bcf4fead142dbadb050c47bca9a59df5f4
1 subroutine da_cloud_model (TB, PB, QB, QCWB, QRNB, ZB, ZFB, DT, kts, kte)
3    !-----------------------------------------------------------------
4    ! Purpose: Calculate DT (=dz/w) using cumulus parameterization 
5    !          of a one-dimensional cloud model.
6    !-----------------------------------------------------------------
8    ! Calculate DT
10    implicit none
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
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                       :: ALPHA, RC, MU, XX, YY
25    real, dimension(kts:kte+1) :: W0, W
27    if (trace_use) call da_trace_entry("da_cloud_model")
29    ALPHA=0.5
30    RC=100.0
31    MU=0.183/RC
33    do k = kts, kte+1
34       W0(k)=0.0  
35       W(k)=0.0  
36    end do
38    do k = kts, kte
39       PP(k)=PB(k)/100.0
40       DT(k)=0.0
41    end do
43    P0 = PP(kts)
44    Z0 = ZB(kts)
45    T0 = MAX(TB(kts),303.0)
47    call da_qfrmrh (P0, T0, 95.0, Q0)
49    call da_lcl (P0, Z0, T0, Q0, PLCL, ZLCL, TLCL, QLCL)
51    call da_qfrmrh (PLCL, TLCL, 95.0, QLCL)
53    call da_cumulus (ZLCL, TLCL, QLCL, PLCL, PP, TB,            &
54                   ZC, TC, QC, KCB, KCT, PCT, ZCT, kts, kte)
56    do k = KCB, KCT
57       TCV(k) = TC(k) * (1.0 + 0.608 * QC(k))
58       TBV(k) = TB(k) * (1.0 + 0.608 * QB(k))
59    
60       B(k) = (TCV(k)-TBV(k)) / TBV(k)
62       QT(k) = QC(k) + QCWB(k) + QRNB(k)
63    end do
65    W0(KCB) = 0.0
66    do k = KCB+1, KCT+1
67       XX = 1.0+2.0*MU*(ZFB(k)-ZFB(k-1))
68       YY = 2.0*gravity*(B(k-1)/(1.0+ALPHA) - QT(k-1)) * (ZFB(k)-ZFB(k-1))
69       W0(k) =  (W0(k-1)+YY) / XX
70    end do
71      
72    do k = KCB, KCT+1
73       if (W0(k) >= 0.0) then
74          W(k) = sqrt(W0(k))
75       end if
76    end do
79    do k = KCT, KCB+1, -1
80       if (W(k) >= 0.01) then
81          DT(k) = (ZB(k)-ZB(k-1))/W(k)
82       else
83          DT(k) = 0.0
84       end if
85    end do
87    if (trace_use) call da_trace_exit("da_cloud_model")
89 end subroutine da_cloud_model