da_evapo_adj.inc

References to this file elsewhere.
1 subroutine da_evapo_adj(DT,SCR3,SCR5,QV_B,PRE,SCR6,  &
2                         SCR39,SCR59,QV_A,PRE9,SCR69, &
3                         kts,kte,kms,kme)
4 
5 
6    !-----------------------------------------------------------------------
7    ! Purpose: Rainwater evaporation
8    !-----------------------------------------------------------------------
9 
10    implicit none
11 
12    integer, intent(in)                     :: kts, kte, kms, kme
13    real, dimension(kms:kme), intent(in)    :: DT, SCR3, SCR5, SCR6, QV_B
14    real, dimension(kms:kme), intent(inout) :: PRE, PRE9
15    real, dimension(kms:kme), intent(inout) :: QV_A
16    real, dimension(kms:kme), intent(out)   :: SCR39, SCR59, SCR69
17 
18    integer             :: k
19    real                :: beta, QRTH
20    real                :: TMP, TMP2
21 
22    QRTH = 1.e-6
23    beta = 0.0486   ! original
24 
25    do K = kts, kte
26  
27       if ( DT(k) <= 0. ) cycle
28  
29       if ( SCR3(K) > QRTH .and. QV_B(k) < SCR5(k) ) then
30          PRE(k) = beta * ( QV_B(K)-SCR5(K) ) * ( SCR6(k)*SCR3(K) )**0.65
31       else if ( SCR3(K) <= QRTH .and. SCR3(k) > 0. &
32                           .and. QV_B(k) < SCR5(k) ) then
33          PRE(k) = beta * ( QV_B(K)-SCR5(K) ) * ( SCR6(k)*QRTH )**0.65
34       else
35          PRE(k) = 0.
36       end if
37 
38    end do
39 
40    do K = kts, kte
41 
42       if (DT(k) <= 0.) cycle
43 
44       SCR39(k) = 0.
45       if ( PRE(k) < -SCR3(k)/DT(k) ) then
46          SCR39(k) = -PRE9(k) / DT(k)
47          PRE9(k)  = 0.
48       end if
49 
50       SCR59(k) = 0.
51       SCR69(k) = 0.
52       if ( SCR3(k) > QRTH .and. QV_B(k) < SCR5(k) ) then
53          TMP  = beta * ( QV_B(k)-SCR5(k) ) * 0.65 * ( SCR6(k)*SCR3(k) )**(-0.35)
54          TMP2 = beta * ( SCR6(k)*SCR3(k) )**0.65
55 
56          QV_A(k) = QV_A(k) + TMP2 * PRE9(k)
57          SCR59(k) = -TMP2 * PRE9(k)
58          SCR39(k) = SCR39(k) + TMP * SCR6(k) * PRE9(k)
59          SCR69(k) = TMP * SCR3(k) * PRE9(k)
60       else if (SCR3(k) <= QRTH .and. SCR3(k) > 0. .and. QV_B(k) < SCR5(k) ) then
61          TMP  = beta * ( QV_B(k)-SCR5(k) ) * 0.65 * ( SCR6(k)*QRTH )**(-0.35)
62          TMP2 = beta * ( SCR6(k)*QRTH )**0.65
63 
64          QV_A(k) = QV_A(k) + TMP2 * PRE9(k)
65          SCR59(k) = -TMP2 * PRE9(k)
66          SCR69(k) = TMP * QRTH * PRE9(k)
67       end if
68 
69    end do
70 
71 end subroutine da_evapo_adj