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 ! To get the xa%rh for the moisture adjustments (YRG, 04/02/2007):
36 
37       if ( cv_options_hum == 1  ) then
38          call DA_TPQ_To_RH_Lin( xb, xp, xa )
39       endif
40 
41    do i=is,ie
42       do j=js,je
43 
44          tol_adjust_moist = 0.0
45          tol_moist        = 0.0
46 
47          DO k=ks,ke
48             dz(k)=xb%hf(i,j,k+1)-xb%hf(i,j,k)
49 
50             imod(k)           = 0
51             x_qs(k)           = 0.
52             rhtol(k)          = xb%rh(i,j,k)+xa%rh(i,j,k)
53          ENDDO
54 
55          do k=ks,ke
56             if (rhtol(k) .gt. maximum_rh) then
57                oldrha       = xa%rh(i,j,k)
58                ! modify xa%rh
59                xa%rh(i,j,k) = maximum_rh - xb%rh(i,j,k)
60 
61                call da_tp_to_qs(xb%t(i,j,k)+xa%t(i,j,k), &
62                   xb%p(i,j,k)+xa%p(i,j,k), es, x_qs(k))
63 
64                ! calculate xa%q
65                call da_tprh_to_q_lin1(xb%t(i,j,k), xb%p(i,j,k), &
66                   xb%es(i,j,k), xb%q(i,j,k), xb%rh(i,j,k),  xa%t(i,j,k), &
67                   xa%p(i,j,k), xa%rh(i,j,k), xa%q(i,j,k))
68 
69                tol_adjust_moist = tol_adjust_moist + x_qs(k)*(oldrha - &
70                   xa%rh(i,j,k))* dz(k)*(xb%rho(i,j,k)+xa%rho(i,j,k))
71                imod(k)=-1
72             end if
73 
74             if (rhtol(k) .lt. minimum_rh) then
75                oldrha           = xa%rh(i,j,k)
76                xa%rh(i,j,k)     = minimum_rh - xb%rh(i,j,k)
77                call da_tp_to_qs(xb%t(i,j,k)+xa%t(i,j,k), &
78                   xb%p(i,j,k)+xa%p(i,j,k), es, x_qs(k))
79 
80                call da_tprh_to_q_lin1(xb%t(i,j,k), xb%p(i,j,k), &
81                   xb%es(i,j,k), xb%q(i,j,k), xb%rh(i,j,k),  xa%t(i,j,k), &
82                   xa%p(i,j,k), xa%rh(i,j,k), xa%q(i,j,k))
83 
84 
85                tol_adjust_moist = tol_adjust_moist + x_qs(k)*(oldrha - &
86                   xa%rh(i,j,k))* dz(k)*(xb%rho(i,j,k)+xa%rho(i,j,k))
87                imod(k)=-1
88             end if
89          end do
90 
91          if (tol_adjust_moist .gt. 0.0) then
92             do k=ks,ke
93                if (rhtol(k) .lt. upper_modify_rh .and. imod(k) .eq. 0) then
94                   call da_tp_to_qs(xb%t(i,j,k)+xa%t(i,j,k), &
95                                     xb%p(i,j,k)+xa%p(i,j,k),es,x_qs(k))
96 
97                   each_moist   = rhtol(k)*x_qs(k)* &
98                                  dz(k)*(xb%rho(i,j,k)+xa%rho(i,j,k))
99                   tol_moist    = tol_moist + each_moist
100                   imod(k)      = 1
101                end if
102             end do
103          end if
104 
105          if (tol_adjust_moist .lt. 0.0) then
106             do k=ks,ke
107                if (rhtol(k) .gt. lower_modify_rh .and. imod(k) .eq. 0) then
108                   call da_tp_to_qs(xb%t(i,j,k)+xa%t(i,j,k), &
109                                     xb%p(i,j,k)+xa%p(i,j,k), es, x_qs(k))
110 
111                   each_moist   = rhtol(k)*x_qs(k)* &
112                                  dz(k)*(xb%rho(i,j,k)+xa%rho(i,j,k))
113                   tol_moist    = tol_moist + each_moist
114                   imod(k)      = 1
115                end if
116             end do
117          end if
118 
119          if (tol_moist > 0) then
120            weight       = tol_adjust_moist/tol_moist
121            do k=ks,ke
122              if (imod(k) .eq. 1) then
123                xa%rh(i,j,k) = xa%rh(i,j,k)+(xb%rh(i,j,k)+xa%rh(i,j,k))*weight
124                call da_tprh_to_q_lin1(xb%t(i,j,k), xb%p(i,j,k), xb%es(i,j,k), &
125                                       xb%q(i,j,k), xb%rh(i,j,k),  xa%t(i,j,k), &
126                                       xa%p(i,j,k), xa%rh(i,j,k), xa%q(i,j,k))
127 
128              end if
129            end do
130          end if
131       end do
132    end do
133 
134    if (trace_use) call da_trace_exit("da_check_rh")
135 
136 end subroutine da_check_rh
137 
138