Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

add use_rh and calrh for grib2 humidity processing #566

Merged
merged 14 commits into from
Sep 1, 2021
Merged
Show file tree
Hide file tree
Changes from 12 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
173 changes: 173 additions & 0 deletions sorc/chgres_cube.fd/grib2_util.F90
Original file line number Diff line number Diff line change
Expand Up @@ -17,6 +17,11 @@ module grib2_util

implicit none

public :: rh2spfh_gfs
public :: rh2spfh_nam
public :: fpvsnew


contains

!> Convert relative humidity to specific humidity.
Expand Down Expand Up @@ -57,6 +62,174 @@ subroutine rh2spfh(rh_sphum,p,t)

end subroutine RH2SPFH



!> Convert relative humidity to specific humidity (NAM formula)
!!
!! @param[inout] rh_sphum rel humidity on input. spec hum on output.
!! @param[in] p pressure in Pa
!! @param[in] t temperature
!! @author Jili Dong NCEP/EMC
subroutine rh2spfh_nam(rh_sphum,p,t)

implicit none

real, parameter :: PQ0=379.90516
real, parameter :: A2=17.2693882
real, parameter :: A3=273.16
real, parameter :: A4=35.86


real(esmf_kind_r4), intent(inout), dimension(i_input,j_input) ::rh_sphum
real(esmf_kind_r8), intent(in) :: p, t(i_input,j_input)

real, dimension(i_input,j_input) :: QC, rh

print*,"- CONVERT RH TO SPFH AT LEVEL ", p

rh = rh_sphum

QC = PQ0/P*EXP(A2*(T-A3)/(T-A4))

!print *, 'T = ', T, ' RH = ', RH, ' P = ', P
rh_sphum = rh*QC/100.0
!print *, 'q = ', sphum

end subroutine RH2SPFH_NAM

!> Convert relative humidity to specific humidity (GFS formula)
!!
!! @param[inout] rh_sphum rel humidity on input. spec hum on output.
!! @param[in] p pressure in Pa
!! @param[in] t temperature
!! @author Jili Dong NCEP/EMC
subroutine rh2spfh_gfs(rh_sphum,p,t)

implicit none

integer kind_phys

parameter (kind_phys = selected_real_kind(13,60)) ! the '60' maps to 64-bit real


real(kind=kind_phys),parameter:: con_rd =2.8705e+2 ! gas constant air (J/kg/K)
real(kind=kind_phys),parameter:: con_rv =4.6150e+2 ! gas constant H2O (J/kg/K)

real(kind=kind_phys),parameter:: con_eps =con_rd/con_rv
real(kind=kind_phys),parameter:: con_epsm1 =con_rd/con_rv-1.




real(esmf_kind_r4), intent(inout), dimension(i_input,j_input) ::rh_sphum
real(esmf_kind_r8), intent(in) :: p, t(i_input,j_input)

real, dimension(i_input,j_input) :: QC, rh
real :: ES
integer :: i,j

print*,"- CONVERT RH TO SPFH AT LEVEL ", p

rh = rh_sphum

do j=1,j_input
do i=1,i_input
ES = MIN(FPVSNEW(T(I,J)),P)
QC(i,j) = CON_EPS*ES/(P+CON_EPSM1*ES)
rh_sphum(i,j) = rh(i,j)*QC(i,j)/100.0
end do
end do


!print *, 'T = ', T, ' RH = ', RH, ' P = ', P
!print *, 'q = ', sphum


end subroutine RH2SPFH_GFS


!> Compute saturation vapor pressure
!!
!! @param[in] t temperature in Kelvin
!! @return fpvsnew Saturation vapor pressure
!! @author N Phillips w/NMC2X2
!! @date 30 dec 82

!
!-------------------------------------------------------------------------------------
!
elemental function fpvsnew(t)
!
implicit none
integer,parameter:: nxpvs=7501
real,parameter:: con_ttp =2.7316e+2 ! temp at H2O 3pt
real,parameter:: con_psat =6.1078e+2 ! pres at H2O 3pt
real,parameter:: con_cvap =1.8460e+3 ! spec heat H2O gas (J/kg/K)
real,parameter:: con_cliq =4.1855e+3 ! spec heat H2O liq
real,parameter:: con_hvap =2.5000e+6 ! lat heat H2O cond
real,parameter:: con_rv =4.6150e+2 ! gas constant H2O
real,parameter:: con_csol =2.1060e+3 ! spec heat H2O ice
real,parameter:: con_hfus =3.3358e+5 ! lat heat H2O fusion
real,parameter:: tliq=con_ttp
real,parameter:: tice=con_ttp-20.0
real,parameter:: dldtl=con_cvap-con_cliq
real,parameter:: heatl=con_hvap
real,parameter:: xponal=-dldtl/con_rv
real,parameter:: xponbl=-dldtl/con_rv+heatl/(con_rv*con_ttp)
real,parameter:: dldti=con_cvap-con_csol
real,parameter:: heati=con_hvap+con_hfus
real,parameter:: xponai=-dldti/con_rv
real,parameter:: xponbi=-dldti/con_rv+heati/(con_rv*con_ttp)
real tr,w,pvl,pvi
real fpvsnew
real(esmf_kind_r8),intent(in):: t
integer jx
real xj,x,tbpvs(nxpvs),xp1
real xmin,xmax,xinc,c2xpvs,c1xpvs
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
xmin=180.0
xmax=330.0
xinc=(xmax-xmin)/(nxpvs-1)
! c1xpvs=1.-xmin/xinc
c2xpvs=1./xinc
c1xpvs=1.-xmin*c2xpvs
! xj=min(max(c1xpvs+c2xpvs*t,1.0),real(nxpvs,krealfp))
xj=min(max(c1xpvs+c2xpvs*t,1.0),float(nxpvs))
jx=min(xj,float(nxpvs)-1.0)
x=xmin+(jx-1)*xinc

tr=con_ttp/x
if(x>=tliq) then
tbpvs(jx)=con_psat*(tr**xponal)*exp(xponbl*(1.-tr))
elseif(x<tice) then
tbpvs(jx)=con_psat*(tr**xponai)*exp(xponbi*(1.-tr))
else
w=(t-tice)/(tliq-tice)
pvl=con_psat*(tr**xponal)*exp(xponbl*(1.-tr))
pvi=con_psat*(tr**xponai)*exp(xponbi*(1.-tr))
tbpvs(jx)=w*pvl+(1.-w)*pvi
endif

xp1=xmin+(jx-1+1)*xinc

tr=con_ttp/xp1
if(xp1>=tliq) then
tbpvs(jx+1)=con_psat*(tr**xponal)*exp(xponbl*(1.-tr))
elseif(xp1<tice) then
tbpvs(jx+1)=con_psat*(tr**xponai)*exp(xponbi*(1.-tr))
else
w=(t-tice)/(tliq-tice)
pvl=con_psat*(tr**xponal)*exp(xponbl*(1.-tr))
pvi=con_psat*(tr**xponai)*exp(xponbi*(1.-tr))
tbpvs(jx+1)=w*pvl+(1.-w)*pvi
endif

fpvsnew=tbpvs(jx)+(xj-jx)*(tbpvs(jx+1)-tbpvs(jx))
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
end function fpvsnew



!> Convert omega to vertical velocity.
!!
!! @param[inout] omega on input, vertical velocity on output
Expand Down
27 changes: 23 additions & 4 deletions sorc/chgres_cube.fd/input_data.F90
Original file line number Diff line number Diff line change
Expand Up @@ -2459,7 +2459,7 @@ subroutine read_input_atm_grib2_file(localpet)

use wgrib2api

use grib2_util, only : rh2spfh, convert_omega
use grib2_util, only : rh2spfh, rh2spfh_gfs, convert_omega

implicit none

Expand Down Expand Up @@ -2490,7 +2490,9 @@ subroutine read_input_atm_grib2_file(localpet)
logical :: lret
logical :: conv_omega=.false., &
hasspfh=.true., &
isnative=.false.
isnative=.false., &
use_rh=.false.


real(esmf_kind_r8), allocatable :: rlevs(:)
real(esmf_kind_r4), allocatable :: dummy2d(:,:)
Expand Down Expand Up @@ -2603,11 +2605,22 @@ subroutine read_input_atm_grib2_file(localpet)
if (localpet==0) print*, "- LEVEL AFTER SORT = ",slevs(i)
enddo
endif

! Is SPFH on full levels Jili Dong
do vlev = 1, lev_input
iret = grb2_inq(the_file,inv_file,':SPFH:',slevs(vlev))
if (iret <= 0) then
use_rh = .TRUE.
if (localpet == 0) print*,':SPFH:'//slevs(vlev)//' not existing. Use RH instead.'
JiliDong-NOAA marked this conversation as resolved.
Show resolved Hide resolved
exit
end if
end do


if (localpet == 0) print*,"- FIND SPFH OR RH IN FILE"
iret = grb2_inq(the_file,inv_file,trim(trac_names_grib_1(1)),trac_names_grib_2(1),lvl_str_space)

if (iret <= 0) then
if (iret <= 0 .or. use_rh) then
iret = grb2_inq(the_file,inv_file, ':var0_2','_1_1:',lvl_str_space)
if (iret <= 0) call error_handler("READING ATMOSPHERIC WATER VAPOR VARIABLE.", iret)
hasspfh = .false.
Expand Down Expand Up @@ -2782,7 +2795,13 @@ subroutine read_input_atm_grib2_file(localpet)
endif !iret<=0

if (n==1 .and. .not. hasspfh) then
call rh2spfh(dummy2d,rlevs(vlev),dummy3d(:,:,vlev))
if (trim(external_model) .eq. 'GFS') then
print *,'CALRH GFS'
call rh2spfh_gfs(dummy2d,rlevs(vlev),dummy3d(:,:,vlev))
else
print *,'CALRH non-GFS'
call rh2spfh(dummy2d,rlevs(vlev),dummy3d(:,:,vlev))
end if
endif

print*,'tracer ',vlev, maxval(dummy2d),minval(dummy2d)
Expand Down
7 changes: 7 additions & 0 deletions tests/chgres_cube/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -83,6 +83,13 @@ target_link_libraries(ftst_dint2p

add_test(NAME chgres_cube-ftst_dint2p COMMAND ftst_dint2p)

add_executable(ftst_rh2spfh_gfs ftst_rh2spfh_gfs.F90)
target_link_libraries(ftst_rh2spfh_gfs
chgres_cube_lib)

add_test(NAME chgres_cube-ftst_rh2spfh_gfs COMMAND ftst_rh2spfh_gfs)


add_executable(ftst_quicksort ftst_quicksort.F90)
target_link_libraries(ftst_quicksort
chgres_cube_lib)
Expand Down
59 changes: 59 additions & 0 deletions tests/chgres_cube/ftst_rh2spfh_gfs.F90
Original file line number Diff line number Diff line change
@@ -0,0 +1,59 @@
program ftst_rh2spfh_gfs

use grib2_util

implicit none

real, parameter :: eps=1E-4

real*4, dimension(3,2) :: rh_500hpa, rh_900hpa, rh_tmp
real*4, dimension(3,2) :: spfh_gfs_900hpa_expected, spfh_gfs_500hpa_expected
real*4, dimension(3,2) :: spfh_nam_900hpa_expected, spfh_nam_500hpa_expected
real*8, dimension(3,2) :: temp_500hpa, temp_900hpa
real*8 :: pres

data temp_900hpa /280.9, 282.6, 281.4, 289.9, 286.5, 288.7/
data rh_900hpa /17.6, 75.2, 99.9, 95.4, 86.6, 89.1/

data temp_500hpa /254.1, 259.7, 259.7, 267.2, 262.0, 269.5/
data rh_500hpa /10.7, 9.1, 0.5, 10.1, 4.2, 0.9/

data spfh_gfs_900hpa_expected /1.286259, 6.169457, 7.554718, 12.64650, 9.210890, 10.93505/
data spfh_gfs_500hpa_expected /0.1519615, 0.2257435, 1.2403490E-02, 0.4856212, 0.1288972, 5.2061662E-02/

data spfh_nam_900hpa_expected /1.281873, 6.145997, 7.528146, 12.56592, 9.164917, 10.87111/
data spfh_nam_500hpa_expected /0.1799187, 0.2447683, 1.3448806E-02, 0.4918184, 0.1360911, 5.2174598E-02/

i_input = 3
j_input = 2

! rh to spfh at 500 hPa with GFSv15/16 conversion
rh_tmp=rh_500hpa
pres=50000.0
call rh2spfh_gfs(rh_tmp,pres,temp_500hpa)
if (any((abs(rh_tmp*1000.0-spfh_gfs_500hpa_expected)) .gt. eps)) stop 500

! rh to spfh at 900 hPa with GFSv15/16 conversion
rh_tmp=rh_900hpa
pres=90000.0
call rh2spfh_gfs(rh_tmp,pres,temp_900hpa)
if (any((abs(rh_tmp*1000.0-spfh_gfs_900hpa_expected)) .gt. eps)) stop 900

! rh to spfh at 500 hPa with GFSv17 conversion
rh_tmp=rh_500hpa
pres=50000.0
call rh2spfh_nam(rh_tmp,pres,temp_500hpa)
if (any((abs(rh_tmp*1000.0-spfh_nam_500hpa_expected)) .gt. eps)) stop 501

! rh to spfh at 900 hPa with GFSv17 conversion
rh_tmp=rh_900hpa
pres=90000.0
call rh2spfh_nam(rh_tmp,pres,temp_900hpa)
if (any((abs(rh_tmp*1000.0-spfh_nam_900hpa_expected)) .gt. eps)) stop 901


print*, "OK"
print*, "SUCCESS!"

end program ftst_rh2spfh_gfs