1 subroutine da_check_rh(grid)
3 !---------------------------------------------------------------------------
4 ! Purpose: Remove RH over 100% or under 10%
5 ! Made the modification to those levels, which RH are less than 95%
6 !---------------------------------------------------------------------------
10 type (domain), intent(inout) :: grid
12 integer :: imod(kts:kte)
13 real :: rhtol(kts:kte)
17 integer :: i, j, k, ij
18 real :: tol_adjust_moist, tol_moist, oldrha, each_moist, es, weight
19 real :: upper_modify_rh, lower_modify_rh
21 if (trace_use) call da_trace_entry("da_check_rh")
24 upper_modify_rh = 95.0
25 lower_modify_rh = 11.0
27 ! To get the grid%xa%rh for the moisture adjustments
29 if (cv_options_hum == cv_options_hum_specific_humidity) then
30 call da_tpq_to_rh_lin (grid)
33 !$OMP PARALLEL DO SCHEDULE (DYNAMIC, 1) &
34 !$OMP PRIVATE ( i, j, k, tol_adjust_moist, tol_moist) &
35 !$OMP PRIVATE ( weight, oldrha, each_moist, imod, dz, x_qs, rhtol, es)
36 ! do ij = 1 , grid%num_tiles
39 !do j=grid%j_start(ij), grid%j_end(ij)
42 tol_adjust_moist = 0.0
46 dz(k)=grid%xb%hf(i,j,k+1)-grid%xb%hf(i,j,k)
50 rhtol(k) = grid%xb%rh(i,j,k)+grid%xa%rh(i,j,k)
54 if (rhtol(k) .gt. maximum_rh) then
55 oldrha = grid%xa%rh(i,j,k)
57 grid%xa%rh(i,j,k) = maximum_rh - grid%xb%rh(i,j,k)
59 call da_tp_to_qs(grid%xb%t(i,j,k)+grid%xa%t(i,j,k), &
60 grid%xb%p(i,j,k)+grid%xa%p(i,j,k), es, x_qs(k))
63 call da_tprh_to_q_lin1(grid%xb%t(i,j,k), grid%xb%p(i,j,k), &
64 grid%xb%es(i,j,k), grid%xb%q(i,j,k), grid%xb%rh(i,j,k), grid%xa%t(i,j,k), &
65 grid%xa%p(i,j,k), grid%xa%rh(i,j,k), grid%xa%q(i,j,k))
67 tol_adjust_moist = tol_adjust_moist + x_qs(k)*(oldrha - &
68 grid%xa%rh(i,j,k))* dz(k)*(grid%xb%rho(i,j,k)+grid%xa%rho(i,j,k))
72 if (rhtol(k) .lt. minimum_rh) then
73 oldrha = grid%xa%rh(i,j,k)
74 grid%xa%rh(i,j,k) = minimum_rh - grid%xb%rh(i,j,k)
75 call da_tp_to_qs(grid%xb%t(i,j,k)+grid%xa%t(i,j,k), &
76 grid%xb%p(i,j,k)+grid%xa%p(i,j,k), es, x_qs(k))
78 call da_tprh_to_q_lin1(grid%xb%t(i,j,k), grid%xb%p(i,j,k), &
79 grid%xb%es(i,j,k), grid%xb%q(i,j,k), grid%xb%rh(i,j,k), grid%xa%t(i,j,k), &
80 grid%xa%p(i,j,k), grid%xa%rh(i,j,k), grid%xa%q(i,j,k))
83 tol_adjust_moist = tol_adjust_moist + x_qs(k)*(oldrha - &
84 grid%xa%rh(i,j,k))* dz(k)*(grid%xb%rho(i,j,k)+grid%xa%rho(i,j,k))
89 if (tol_adjust_moist .gt. 0.0) then
91 if (rhtol(k) .lt. upper_modify_rh .and. imod(k) .eq. 0) then
92 call da_tp_to_qs(grid%xb%t(i,j,k)+grid%xa%t(i,j,k), &
93 grid%xb%p(i,j,k)+grid%xa%p(i,j,k),es,x_qs(k))
95 each_moist = rhtol(k)*x_qs(k)* &
96 dz(k)*(grid%xb%rho(i,j,k)+grid%xa%rho(i,j,k))
97 tol_moist = tol_moist + each_moist
103 if (tol_adjust_moist .lt. 0.0) then
105 if (rhtol(k) .gt. lower_modify_rh .and. imod(k) .eq. 0) then
106 call da_tp_to_qs(grid%xb%t(i,j,k)+grid%xa%t(i,j,k), &
107 grid%xb%p(i,j,k)+grid%xa%p(i,j,k), es, x_qs(k))
109 each_moist = rhtol(k)*x_qs(k)* &
110 dz(k)*(grid%xb%rho(i,j,k)+grid%xa%rho(i,j,k))
111 tol_moist = tol_moist + each_moist
117 if (tol_moist > 0) then
118 weight = tol_adjust_moist/tol_moist
120 if (imod(k) .eq. 1) then
121 grid%xa%rh(i,j,k) = grid%xa%rh(i,j,k)+(grid%xb%rh(i,j,k)+grid%xa%rh(i,j,k))*weight
123 ! To guarantee the adjusted relative humidity will not be out of the range (YRG, 10/21/2008):
124 oldrha = grid%xa%rh(i,j,k)+grid%xb%rh(i,j,k)
125 if ( (oldrha > maximum_rh) ) then
126 grid%xa%rh(i,j,k) = maximum_rh - grid%xb%rh(i,j,k)
127 if (print_detail_xa ) &
128 write(unit=stdout, fmt='(3I5," Warning=== Adjusted RH > maximum_rh=",f10.2,&
129 & " total_rh, xb%rh, xa%rh:",3f10.2)') &
130 i, j, k, maximum_rh, oldrha, grid%xb%rh(i,j,k), grid%xa%rh(i,j,k)
131 else if ( oldrha < minimum_rh ) then
132 grid%xa%rh(i,j,k) = minimum_rh - grid%xb%rh(i,j,k)
133 if (print_detail_xa ) &
134 write(unit=stdout, fmt='(3I5," Warning=== Adjusted RH < minimum_rh=",f10.2,&
135 & " total_rh, xb%rh, xa%rh:",3f10.2)') &
136 i, j, k, minimum_rh, oldrha, grid%xb%rh(i,j,k), grid%xa%rh(i,j,k)
138 ! ...........................................................................................
140 call da_tprh_to_q_lin1(grid%xb%t(i,j,k), grid%xb%p(i,j,k), grid%xb%es(i,j,k), &
141 grid%xb%q(i,j,k), grid%xb%rh(i,j,k), grid%xa%t(i,j,k), &
142 grid%xa%p(i,j,k), grid%xa%rh(i,j,k), grid%xa%q(i,j,k))
151 !$OMP END PARALLEL DO
153 if (trace_use) call da_trace_exit("da_check_rh")
155 end subroutine da_check_rh