wrf svn trunk commit r4103
[wrffire.git] / wrfv2_fire / var / da / da_physics / da_trh_to_td.inc
bloba86a732d2fe1a98683a51653ee0805487bc9e54e
1 subroutine da_trh_to_td (grid)
3    !---------------------------------------------------------------------
4    !
5    !                       function f_td_from_rh
6    !                     **************************
7    !
8    !  purpose:
9    !  -------
10    !     compute dew point from temperature and relative humidity
11    !
12    !   method:
13    !   ------
14    !     invert the relation
15    !
16    !     rh = 100.0 * exp (l_over_rv * (1.0/t - 1.0/td))
17    !
18    !   input:
19    !   -----
20    !      t_k:   temperature       in k
21    !      rh:    relative humidity in %
22    !
23    !   output:
24    !   ------
25    !      td:    dew point in k
26    !
27    !   references:
28    !   -----------
29    !    R. R. Rogers and M. K. Yau, 1989: a short course in cloud physics,
30    !                                   3nd edition, pergamon press, page 14-19.
31    !
32    !   verification set:
33    !   -----------------
34    !    t_k  = 268.15 k,  
35    !    td_k = 262.55 k
36    !    rh   = 65 %, 
37    !    p_pa = 80000  pa, 
38    !    qv   = 2.11e-03 kg/kg,
39    !
40    !  modifications:
41    !   ------------
42    !    parallel implementation. -al bourgeoits
43    ! 
44    !-------------------------------------------------------------------------
46    implicit none
48    type (domain), intent(inout) :: grid
50    integer :: i, j, k, ij
52    real    :: invdifftd, invtd
54    if (trace_use_dull) call da_trace_entry("da_trh_to_td")
56    !$OMP PARALLEL DO &
57    !$OMP PRIVATE( ij, i, j, k, invdifftd, invtd )
58    do ij = 1 , grid%num_tiles
60    do k=kts,kte
61       do j=grid%j_start(ij), grid%j_end(ij)
62          do i=its,ite
63             if (grid%xb%rh(i,j,k) < 10.0) then
64                grid%xb%rh(i,j,k) = 10.0
65             else if (grid%xb%rh(i,j,k) > 105.0) then
66                grid%xb%rh(i,j,k) = 105.0
67             end if
69             invdifftd = log (grid%xb%rh(i,j,k)/100.0) / l_over_rv
71             invtd = 1/grid%xb%t(i,j,k)  - invdifftd
73             grid%xb%td(i,j,k)  = 1.0 / invtd
75             if (grid%xb%td(i,j,k) > grid%xb%t(i,j,k)) &
76                grid%xb%td(i,j,k) = grid%xb%t(i,j,k)
77          end do
78       end do
79    end do
81    end do
82    !$OMP END PARALLEL DO
84    if (trace_use_dull) call da_trace_exit("da_trh_to_td")
86 end subroutine da_trh_to_td