Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_physics / da_evapo_lin.inc
blob2424a26cdb2aa73738f7dc8d387623c9baca8b34
1 subroutine da_evapo_lin(dt,scr3,scr5,qv_b,pre,scr6, scr39,scr59,qv_a,pre9,scr69, kts,kte,kms,kme)
3    !-----------------------------------------------------------------------
4    ! Purpose: Rainwater evaporation
5    !-----------------------------------------------------------------------
7    implicit none
9    integer, intent(in)  :: kts, kte, kms, kme
10    real,    intent(in)  :: dt(kms:kme)
11    real,    intent(in)  :: scr3(kms:kme)
12    real,    intent(in)  :: scr5(kms:kme)
13    real,    intent(in)  :: scr6(kms:kme)
14    real,    intent(in)  :: qv_b(kms:kme)
15    real,    intent(out) :: pre(kms:kme)
16    real,    intent(out) :: pre9(kms:kme)
17    real,    intent(in)  :: scr39(kms:kme)
18    real,    intent(in)  :: scr59(kms:kme)
19    real,    intent(in)  :: scr69(kms:kme)
20    real,    intent(in)  :: qv_a(kms:kme)
22    integer :: k
23    real    :: beta, qrth
24    real    :: tmp, tmp2
26    if (trace_use) call da_trace_entry("da_evapo_lin")
28    qrth = 1.0e-6
29    beta = 0.0486   ! original
31    do k = kts, kte
32       if (dt(k) <= 0.0) cycle
34       if ( scr3(k) > qrth .and. qv_b(k) < scr5(k) ) then
35          tmp  = beta * ( qv_b(k)-scr5(k) )* 0.65 * ( scr6(k)*scr3(k) )**(-0.35)
36          tmp2 = beta * ( scr6(k)*scr3(k) )**0.65
37          pre9(k) = tmp * ( scr69(k)*scr3(k)+scr6(k)*scr39(k) ) + &
38                    tmp2 * ( qv_a(k)-scr59(k) )
39          pre(k)  = beta * ( qv_b(k)-scr5(k) ) * ( scr6(k)*scr3(k) )**0.65
40       else if ( scr3(k) <= qrth .and. scr3(k) > 0.0 .and. qv_b(k) < scr5(k) ) then
41          tmp  = beta * ( qv_b(k)-scr5(k) ) * 0.65 * ( scr6(k)*qrth )**(-0.35)
42          tmp2 = beta * ( scr6(k)*qrth )**0.65
43          pre9(k) = tmp * ( scr69(k)*qrth ) + tmp2 * ( qv_a(k)-scr59(k) )
44          pre(k)  = beta * ( qv_b(k)-scr5(k) ) * ( scr6(k)*qrth )**0.65
45       else
46          pre9(k) = 0.0
47          pre(k) = 0.0
48       end if
50       if ( pre(k) < -scr3(k)/dt(k) ) then
51          pre9(k) = -scr39(k) / dt(k)
52          pre(k)  = -scr3(k) / dt(k)
53       end if
54    end do
56    if (trace_use) call da_trace_exit("da_evapo_lin")
58 end subroutine da_evapo_lin