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 7 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
196 changes: 196 additions & 0 deletions sorc/chgres_cube.fd/grib2_util.F90
Original file line number Diff line number Diff line change
Expand Up @@ -57,6 +57,202 @@ 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

!
!-------------------------------------------------------------------------------------
!
elemental function fpvsnew(t)
!$$$ Subprogram Documentation Block
!
! Subprogram: fpvsnew Compute saturation vapor pressure
! Author: N Phillips w/NMC2X2 Date: 30 dec 82
!
JiliDong-NOAA marked this conversation as resolved.
Show resolved Hide resolved
! Abstract: Compute saturation vapor pressure from the temperature.
! A linear interpolation is done between values in a lookup table
! computed in gpvs. See documentation for fpvsx for details.
! Input values outside table range are reset to table extrema.
! The interpolation accuracy is almost 6 decimal places.
! On the Cray, fpvs is about 4 times faster than exact calculation.
! This function should be expanded inline in the calling routine.
!
! Program History Log:
! 91-05-07 Iredell made into inlinable function
! 94-12-30 Iredell expand table
! 1999-03-01 Iredell f90 module
! 2001-02-26 Iredell ice phase
!
! Usage: pvs=fpvsnew(t)
!
! Input argument list:
! t Real(krealfp) temperature in Kelvin
!
! Output argument list:
! fpvsnew Real(krealfp) saturation vapor pressure in Pascals
!
! Attributes:
! Language: Fortran 90.
!$$$
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,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
13 changes: 10 additions & 3 deletions sorc/chgres_cube.fd/input_data.F90
Original file line number Diff line number Diff line change
Expand Up @@ -2459,7 +2459,8 @@ 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
use program_setup, only : use_rh, calrh

implicit none

Expand Down Expand Up @@ -2607,7 +2608,7 @@ subroutine read_input_atm_grib2_file(localpet)
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 +2783,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 (calrh .eq. 0) then
call rh2spfh(dummy2d,rlevs(vlev),dummy3d(:,:,vlev))
else if (calrh .eq. 1) then
call rh2spfh_gfs(dummy2d,rlevs(vlev),dummy3d(:,:,vlev))
else
call error_handler("calrh option not supported",1)
end if
endif

print*,'tracer ',vlev, maxval(dummy2d),minval(dummy2d)
Expand Down
3 changes: 3 additions & 0 deletions sorc/chgres_cube.fd/program_setup.F90
Original file line number Diff line number Diff line change
Expand Up @@ -88,13 +88,15 @@ module program_setup
integer, public :: cycle_hour = -999 !< Cycle hour.
integer, public :: regional = 0 !< For regional target grids. When '1' remove boundary halo region from atmospheric/surface data and
!! output atmospheric boundary file. When '2' output boundary file only. Default is '0' (global grids).
integer, public :: calrh = 0 !< rh to spfh calculation. 0:original; 1: GFSv15/16
JiliDong-NOAA marked this conversation as resolved.
Show resolved Hide resolved
integer, public :: halo_bndy = 0 !< Number of row/cols of lateral halo, where pure lateral bndy conditions are applied (regional target grids).
integer, public :: halo_blend = 0 !< Number of row/cols of blending halo, where model tendencies and lateral boundary tendencies are applied. Regional target grids only.
integer, public :: nsoill_out = 4 !< Number of soil levels desired in the output data. chgres_cube can interpolate from 9 input to 4 output levels. DEFAULT: 4.

logical, public :: convert_atm = .false. !< Convert atmospheric data when true.
logical, public :: convert_nst = .false. !< Convert nst data when true.
logical, public :: convert_sfc = .false. !< Convert sfc data when true.
logical, public :: use_rh = .false. !< for grib2, use rh (if true) or spfh (if false)
logical, public :: wam_cold_start = .false. !< When true, cold start for whole atmosphere model.

! Options for replacing vegetation/soil type, veg fraction, and lai with data from the grib2 file
Expand Down Expand Up @@ -181,6 +183,7 @@ subroutine read_setup_namelist(filename)
cycle_year, cycle_mon, cycle_day, &
cycle_hour, convert_atm, &
convert_nst, convert_sfc, &
use_rh, calrh, &
wam_cold_start, &
vgtyp_from_climo, &
sotyp_from_climo, &
Expand Down