wrf svn trunk commit r4103
[wrffire.git] / wrfv2_fire / var / da / da_physics / da_check_rh.inc
blob8d70db91f15800b4be03cd06e8929399d2ecd928
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    !---------------------------------------------------------------------------
8    implicit none
10    type (domain), intent(inout) :: grid
12    integer   :: imod(kts:kte)
13    real      :: rhtol(kts:kte)
14    real      :: x_qs(kts:kte)
15    real      :: dz(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)
31    end if
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
38    do i=its,ite
39       !do j=grid%j_start(ij), grid%j_end(ij)
40       do j=jts, jte
42          tol_adjust_moist = 0.0
43          tol_moist        = 0.0
45          do k=kts,kte
46             dz(k)=grid%xb%hf(i,j,k+1)-grid%xb%hf(i,j,k)
48             imod(k)           = 0
49             x_qs(k)           = 0.0
50             rhtol(k)          = grid%xb%rh(i,j,k)+grid%xa%rh(i,j,k)
51          end do
53          do k=kts,kte
54             if (rhtol(k) .gt. maximum_rh) then
55                oldrha       = grid%xa%rh(i,j,k)
56                ! modify grid%xa%rh
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))
62                ! calculate grid%xa%q
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))
69                imod(k)=-1
70             end if
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))
85                imod(k)=-1
86             end if
87          end do
89          if (tol_adjust_moist .gt. 0.0) then
90             do k=kts,kte
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
98                   imod(k)      = 1
99                end if
100             end do
101          end if
103          if (tol_adjust_moist .lt. 0.0) then
104             do k=kts,kte
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
112                   imod(k)      = 1
113                end if
114             end do
115          end if
117          if (tol_moist > 0) then
118            weight       = tol_adjust_moist/tol_moist
119            do k=kts,kte
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)
137                endif
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))
144              end if
145            end do
146          end if
147       end do
148    end do
150 !  end do
151    !$OMP END PARALLEL DO
153    if (trace_use) call da_trace_exit("da_check_rh")
155 end subroutine da_check_rh