Update version info for release v4.6.1 (#2122)
[WRF.git] / var / da / da_gen_be / da_get_trh.inc
blob1f313c40c5e7ad648985839b9d1160e2048cdec9
1 subroutine da_get_trh( input_file, dim1, dim2, dim3, k, temp, rh )
3    !---------------------------------------------------------------------------
4    ! Purpose: Calculates T, RH from input WRF file.
5    !---------------------------------------------------------------------------
7    implicit none
9    character(len=200), intent(in)  :: input_file       ! NETCDF file name.
10    integer,            intent(in)  :: dim1, dim2, dim3          ! Dimensions.
11    integer,            intent(in)  :: k                         ! Model level.
12    real,               intent(out) :: temp(1:dim1,1:dim2)       ! Temperature.
13    real,               intent(out) :: rh(1:dim1,1:dim2)         ! Relative humidity.
15    character(len=10) :: var                       ! Variable to search for. var = "T"
16    integer           :: i, j                      ! Loop counters.
18    real              :: thetap(1:dim1,1:dim2)     ! Perturbation potential temperature.
19    real              :: pb(1:dim1,1:dim2)         ! Base state pressure.
20    real              :: pp(1:dim1,1:dim2)         ! Pressure perturbation.
21    real              :: x(1:dim1,1:dim2)          ! Vapor mixing ratio.
23    real              :: theta                     ! Potential temperature.
24    real              :: p                         ! Pressure.
25    real              :: q                         ! Specific humidity.
26    real              :: t_c                       ! Temp(Celsius).
27    real              :: es                        ! Saturation vapor pressure.
28    real              :: qs                        ! Saturation specific humidity.
30    if (trace_use) call da_trace_entry("da_get_trh")
32    var = "T" ! Perturbation potential temperature in WRF.
33    call da_get_field( input_file, var, 3, dim1, dim2, dim3, k, thetap)
35    var = "PB"  ! Base state pressure in WRF.
36    call da_get_field( input_file, var, 3, dim1, dim2, dim3, k, pb)
38    var = "P" ! Perturbation pressure in WRF.
39    call da_get_field( input_file, var, 3, dim1, dim2, dim3, k, pp)
41    var = "QVAPOR"  ! Water vapor mixing ratio.
42    call da_get_field( input_file, var, 3, dim1, dim2, dim3, k, x)
44    do j = 1, dim2
45       do i = 1, dim1
46          ! Convert p', theta' to T:
47          theta = t0 + thetap(i,j)                 ! Theta = Theta0 + Thetap
48          p = pb(i,j) + pp(i,j)                     ! p = base p + pert p.
49          temp(i,j) = theta *( p/base_pres)**kappa ! Theta to T.
51          ! Convert to specific humidity.
52          q = x(i,j) /( 1.0 + x(i,j))
54          ! Calculate relative humidity:
55          t_c = temp(i,j) - t_kelvin
56          es = es_alpha * exp( es_beta * t_c /( t_c + es_gamma))
57          qs = rd_over_rv * es /( p - rd_over_rv1 * es)
58          rh(i,j) = 100.0 * q / qs
59       end do
60    end do
62    if (trace_use) call da_trace_exit("da_get_trh")
64 end subroutine da_get_trh