da_check_rh.inc

References to this file elsewhere.
1 subroutine da_check_rh(xb, xa, xp)
2 
3    !---------------------------------------------------------------------------
4    ! Purpose: Remove RH over 100% or under 10%
5    !          Made the modification to those levels, which RH are less than 95%
6    !---------------------------------------------------------------------------
7 
8    implicit none
9 
10    type (xb_type), intent(in)           :: xb    ! First guess state
11    type (x_type),  intent(inout)        :: xa    ! Analysis increment
12    type (xpose_type), intent(in)        :: xp    ! domain decomposition vars.
13 
14    integer   :: imod(xp%kts:xp%kte)
15    real      :: rhtol(xp%kts:xp%kte)
16    real      :: x_qs(xp%kts:xp%kte)
17    real      :: dz(xp%kts:xp%kte)
18 
19    integer :: is,ie, js,je, ks,ke, i, j, k
20    real    :: tol_adjust_moist, tol_moist, oldrha, each_moist, es, weight
21    real    :: upper_modify_rh, lower_modify_rh
22 
23    if (trace_use) call da_trace_entry("da_check_rh")
24 
25    is = xp % its
26    ie = xp % ite
27    js = xp % jts
28    je = xp % jte
29    ks = xp % kts
30    ke = xp % kte
31 
32    upper_modify_rh = 95.0
33    lower_modify_rh = 11.0
34 
35    do i=is,ie
36       do j=js,je
37 
38          tol_adjust_moist = 0.0
39          tol_moist        = 0.0
40 
41          dz(1)             = xb%ztop-xb%hf(i,j,1)
42          imod(1)           = 0
43          x_qs(1)           = 0.0
44          rhtol(1)          = xb%rh(i,j,1)+xa%rh(i,j,1)
45 
46          do k=ks+1,ke
47             dz(k)=xb%hf(i,j,k-1)-xb%hf(i,j,k)
48 
49             imod(k)           = 0
50             x_qs(k)           = 0.0
51             rhtol(k)          = xb%rh(i,j,k)+xa%rh(i,j,k)
52          end do
53 
54          do k=ks,ke
55             if (rhtol(k) .gt. maximum_rh) then
56                oldrha       = xa%rh(i,j,k)
57                ! modify xa%rh
58                xa%rh(i,j,k) = maximum_rh - xb%rh(i,j,k)
59 
60                call da_tp_to_qs(xb%t(i,j,k)+xa%t(i,j,k), &
61                   xb%p(i,j,k)+xa%p(i,j,k), es, x_qs(k))
62 
63                ! calculate xa%q
64                call da_tprh_to_q_lin1(xb%t(i,j,k), xb%p(i,j,k), &
65                   xb%es(i,j,k), xb%q(i,j,k), xb%rh(i,j,k),  xa%t(i,j,k), &
66                   xa%p(i,j,k), xa%rh(i,j,k), xa%q(i,j,k))
67 
68                tol_adjust_moist = tol_adjust_moist + x_qs(k)*(oldrha - &
69                   xa%rh(i,j,k))* dz(k)*(xb%rho(i,j,k)+xa%rho(i,j,k))
70                imod(k)=-1
71             end if
72 
73             if (rhtol(k) .lt. minimum_rh) then
74                oldrha           = xa%rh(i,j,k)
75                xa%rh(i,j,k)     = minimum_rh - xb%rh(i,j,k)
76                call da_tp_to_qs(xb%t(i,j,k)+xa%t(i,j,k), &
77                   xb%p(i,j,k)+xa%p(i,j,k), es, x_qs(k))
78 
79                call da_tprh_to_q_lin1(xb%t(i,j,k), xb%p(i,j,k), &
80                   xb%es(i,j,k), xb%q(i,j,k), xb%rh(i,j,k),  xa%t(i,j,k), &
81                   xa%p(i,j,k), xa%rh(i,j,k), xa%q(i,j,k))
82 
83 
84                tol_adjust_moist = tol_adjust_moist + x_qs(k)*(oldrha - &
85                   xa%rh(i,j,k))* dz(k)*(xb%rho(i,j,k)+xa%rho(i,j,k))
86                imod(k)=-1
87             end if
88          end do
89 
90          if (tol_adjust_moist .gt. 0.0) then
91             do k=ks,ke
92                if (rhtol(k) .lt. upper_modify_rh .and. imod(k) .eq. 0) then
93                   call da_tp_to_qs(xb%t(i,j,k)+xa%t(i,j,k), &
94                                     xb%p(i,j,k)+xa%p(i,j,k),es,x_qs(k))
95 
96                   each_moist   = rhtol(k)*x_qs(k)* &
97                                  dz(k)*(xb%rho(i,j,k)+xa%rho(i,j,k))
98                   tol_moist    = tol_moist + each_moist
99                   imod(k)      = 1
100                end if
101             end do
102          end if
103 
104          if (tol_adjust_moist .lt. 0.0) then
105             do k=ks,ke
106                if (rhtol(k) .gt. lower_modify_rh .and. imod(k) .eq. 0) then
107                   call da_tp_to_qs(xb%t(i,j,k)+xa%t(i,j,k), &
108                                     xb%p(i,j,k)+xa%p(i,j,k), es, x_qs(k))
109 
110                   each_moist   = rhtol(k)*x_qs(k)* &
111                                  dz(k)*(xb%rho(i,j,k)+xa%rho(i,j,k))
112                   tol_moist    = tol_moist + each_moist
113                   imod(k)      = 1
114                end if
115             end do
116          end if
117 
118          if (tol_moist > 0) then
119            weight       = tol_adjust_moist/tol_moist
120            do k=ks,ke
121              if (imod(k) .eq. 1) then
122                xa%rh(i,j,k) = xa%rh(i,j,k)+(xb%rh(i,j,k)+xa%rh(i,j,k))*weight
123                call da_tprh_to_q_lin1(xb%t(i,j,k), xb%p(i,j,k), xb%es(i,j,k), &
124                                       xb%q(i,j,k), xb%rh(i,j,k),  xa%t(i,j,k), &
125                                       xa%p(i,j,k), xa%rh(i,j,k), xa%q(i,j,k))
126 
127              end if
128            end do
129          end if
130       end do
131    end do
132 
133    if (trace_use) call da_trace_exit("da_check_rh")
134 
135 end subroutine da_check_rh
136 
137