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    if (trace_use) call da_trace_entry("da_evapo_adj")
23 
24    QRTH = 1.e-6
25    beta = 0.0486   ! original
26 
27    do K = kts, kte
28  
29       if ( DT(k) <= 0. ) cycle
30  
31       if ( SCR3(K) > QRTH .and. QV_B(k) < SCR5(k) ) then
32          PRE(k) = beta * ( QV_B(K)-SCR5(K) ) * ( SCR6(k)*SCR3(K) )**0.65
33       else if ( SCR3(K) <= QRTH .and. SCR3(k) > 0. &
34                           .and. QV_B(k) < SCR5(k) ) then
35          PRE(k) = beta * ( QV_B(K)-SCR5(K) ) * ( SCR6(k)*QRTH )**0.65
36       else
37          PRE(k) = 0.
38       end if
39 
40    end do
41 
42    do K = kts, kte
43 
44       if (DT(k) <= 0.) cycle
45 
46       SCR39(k) = 0.
47       if ( PRE(k) < -SCR3(k)/DT(k) ) then
48          SCR39(k) = -PRE9(k) / DT(k)
49          PRE9(k)  = 0.
50       end if
51 
52       SCR59(k) = 0.
53       SCR69(k) = 0.
54       if ( SCR3(k) > QRTH .and. QV_B(k) < SCR5(k) ) then
55          TMP  = beta * ( QV_B(k)-SCR5(k) ) * 0.65 * ( SCR6(k)*SCR3(k) )**(-0.35)
56          TMP2 = beta * ( SCR6(k)*SCR3(k) )**0.65
57 
58          QV_A(k) = QV_A(k) + TMP2 * PRE9(k)
59          SCR59(k) = -TMP2 * PRE9(k)
60          SCR39(k) = SCR39(k) + TMP * SCR6(k) * PRE9(k)
61          SCR69(k) = TMP * SCR3(k) * PRE9(k)
62       else if (SCR3(k) <= QRTH .and. SCR3(k) > 0. .and. QV_B(k) < SCR5(k) ) then
63          TMP  = beta * ( QV_B(k)-SCR5(k) ) * 0.65 * ( SCR6(k)*QRTH )**(-0.35)
64          TMP2 = beta * ( SCR6(k)*QRTH )**0.65
65 
66          QV_A(k) = QV_A(k) + TMP2 * PRE9(k)
67          SCR59(k) = -TMP2 * PRE9(k)
68          SCR69(k) = TMP * QRTH * PRE9(k)
69       end if
70 
71    end do
72 
73    if (trace_use) call da_trace_exit("da_evapo_adj")
74 
75 end subroutine da_evapo_adj