chgres_cube  1.4.0
 All Data Structures Files Functions Variables
grib2_util.F90
Go to the documentation of this file.
1 
4 
12 module grib2_util
13 
14 use esmf
15 
16 use model_grid, only : i_input, j_input
17 
18 implicit none
19 
20 contains
21 
29  subroutine rh2spfh(rh_sphum,p,t)
30 
31  implicit none
32  real,parameter :: alpha=-9.477e-4 , & !K^-1,
33  tnot=273.15, & !K
34  lnot=2.5008e6, & !JKg^-1
35  rv=461.51, & !JKg^-1K^-1
36  esnot=611.21 !Pa
37 
38  real(esmf_kind_r4), intent(inout), dimension(i_input,j_input) ::rh_sphum
39  real(esmf_kind_r8), intent(in) :: p, t(i_input,j_input)
40 
41  real, dimension(i_input,j_input) :: es, e, rh
42 
43  print*,"- CONVERT RH TO SPFH AT LEVEL ", p
44 
45  rh = rh_sphum
46  !print *, 'T = ', T, ' RH = ', RH, ' P = ', P
47  es = esnot * exp( lnot/rv * ((t-tnot)/(t*tnot) + alpha * log(t/tnot) - alpha * (t-tnot)/ t))
48  !print *, 'es = ', es
49  e = rh * es / 100.0
50  !print *, 'e = ', e
51  rh_sphum = 0.622 * e / p
52  !print *, 'q = ', sphum
53 
54  !if (P .eq. 100000.0) THEN
55  ! print *, 'T = ', T, ' RH = ', RH, ' P = ', P, ' es = ', es, ' e = ', e, ' q = ', sphum
56  !end if
57 
58 end subroutine rh2spfh
59 
70 subroutine convert_omega(omega,p,t,q,clb,cub)
71 
72  implicit none
73  real(esmf_kind_r8), pointer :: omega(:,:,:), p(:,:,:), t(:,:,:), q(:,:,:),omtmp,ptmp
74 
75  integer :: clb(3), cub(3), i ,j, k
76 
77  real, parameter :: rd = 287.15_esmf_kind_r8, & !JKg^-1K^-1
78  rv=461.51_esmf_kind_r8, & !JKg^-1K^-1
79  g = 9.81_esmf_kind_r8 ! ms^-2
80 
81  real(esmf_kind_r8) :: tv, w
82 
83  do k = clb(3),cub(3)
84  do j = clb(2),cub(2)
85  do i = clb(1),cub(1)
86  tv = t(i,j,k)*(1+rd/rv*q(i,j,k))
87  omtmp=>omega(i,j,k)
88  ptmp=>p(i,j,k)
89 
90  w = -1 * omtmp * rd * tv / (ptmp * g)
91  omega(i,j,k)=w
92  enddo
93  enddo
94  enddo
95 
96 end subroutine convert_omega
97 
98  end module grib2_util
Sets up the ESMF grid objects for the input data grid and target FV3 grid.
Definition: model_grid.F90:9
subroutine rh2spfh(rh_sphum, p, t)
Convert relative humidity to specific humidity.
Definition: grib2_util.F90:29
Utilities for use when reading grib2 data.
Definition: grib2_util.F90:12
subroutine convert_omega(omega, p, t, q, clb, cub)
Convert omega to vertical velocity.
Definition: grib2_util.F90:70