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

Added Princeton V3 met data reader #5

Merged
merged 3 commits into from
Dec 21, 2018
Merged
Show file tree
Hide file tree
Changes from all 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
8 changes: 6 additions & 2 deletions lis/metforcing/princeton/0Intro_princeton.txt
Original file line number Diff line number Diff line change
Expand Up @@ -8,7 +8,11 @@
!
!BOP
!\section{Princeton}
!This is the 50 year forcing data~\cite{princeton} produced at Princeton
!University. The data is global at 1 degree spatial resolution, at 3 hourly
!This is the 50 year forcing data (Sheffield et al., 2006) produced at Princeton
!University.
!Versions 2 and 2.2 of the data are global at 1 degree spatial resolution, at 3 hourly
!intervals.
!Version 3 is at 0.25 degrees and is psedo-global (southern boundary at 59.875 deg S), at
!3 hourly intervals. Unlike versions 2 and 2.2, version 3 does not have precipitation data
!at 3-hour intervals.
!EOP
38 changes: 28 additions & 10 deletions lis/metforcing/princeton/princeton_forcingMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -113,6 +113,7 @@ module princeton_forcingMod
! !REVISION HISTORY:
! 26Jan2007: Hiroko Kato; Initial Specification
! 15 May 2017: Bailing Li; Added changes for reading in version 2.2 data
! 22 Oct 2018: Daniel Sarmiento; Added changes to support version 3 data
!
! !INTERFACE:
subroutine init_PRINCETON(findex)
Expand Down Expand Up @@ -158,10 +159,16 @@ subroutine init_PRINCETON(findex)

LIS_rc%met_nf(findex) = 9

princeton_struc(:)%ncold = 360
princeton_struc(:)%nrold = 180

do n=1,LIS_rc%nnest
!Set dataset dimensions
if (princeton_struc(n)%version == "2" .OR. princeton_struc(n)%version == "2.2") then
princeton_struc(:)%ncold = 360
princeton_struc(:)%nrold = 180
elseif (princeton_struc(n)%version == "3") then
!Dimensions of driver data changed from versions 2.x to version 3
princeton_struc(:)%ncold = 1440
princeton_struc(:)%nrold = 600
endif

! Forecast mode:
if(LIS_rc%forecastMode.eq.1) then
Expand Down Expand Up @@ -207,13 +214,24 @@ subroutine init_PRINCETON(findex)
gridDesci(1) = 0
gridDesci(2) = princeton_struc(n)%ncold
gridDesci(3) = princeton_struc(n)%nrold
gridDesci(4) = -89.50
gridDesci(5) = -179.50
gridDesci(6) = 128
gridDesci(7) = 89.50
gridDesci(8) = 179.50
gridDesci(9) = 1.00
gridDesci(10) = 1.00
!Define driver data domains
if (princeton_struc(n)%version == "2" .OR. princeton_struc(n)%version == "2.2") then
gridDesci(4) = -89.50
gridDesci(5) = -179.50
gridDesci(6) = 128
gridDesci(7) = 89.50
gridDesci(8) = 179.50
gridDesci(9) = 1.00
gridDesci(10) = 1.00
elseif (princeton_struc(n)%version == "3") then
gridDesci(4) = -59.875
gridDesci(5) = -179.875
gridDesci(6) = 128
gridDesci(7) = 89.875
gridDesci(8) = 179.875
gridDesci(9) = 0.25
gridDesci(10) = 0.25
endif
gridDesci(20) = 0
princeton_struc(n)%mi = princeton_struc(n)%ncold*princeton_struc(n)%nrold

Expand Down
47 changes: 34 additions & 13 deletions lis/metforcing/princeton/read_princeton.F90
Original file line number Diff line number Diff line change
Expand Up @@ -19,6 +19,7 @@
! correct individual fields.
! 15 May 2017: Bailing Li; Added changes for reading in version 2.2 data
! that is in 3D array (4D in version 1 & 2).
! 22 Oct 2018: Daniel Sarmiento; Added changes to support version 3 data
!
! !INTERFACE:
subroutine read_princeton( order, n, findex, yr, mon, da, hr, ferror )
Expand Down Expand Up @@ -134,6 +135,7 @@ subroutine read_princeton( order, n, findex, yr, mon, da, hr, ferror )
real :: gridDesco(50) ! Input,output grid info arrays

integer :: kk ! forecast index
integer :: num_met ! number of fields, V3 will have 6, 2 and 2.2 have 7 (precip)

real,allocatable :: datain(:,:) ! input data (longitude,latitude)
real,allocatable :: temp2princeton(:,:,:)
Expand All @@ -152,6 +154,7 @@ subroutine read_princeton( order, n, findex, yr, mon, da, hr, ferror )
! If a problem, ferror is set to zero
ferror = 1
mo = LIS_rc%lnc(n)*LIS_rc%lnr(n)
num_met = 7

! Allocate memory
allocate(datain(princeton_struc(n)%ncold,princeton_struc(n)%nrold))
Expand Down Expand Up @@ -188,6 +191,13 @@ subroutine read_princeton( order, n, findex, yr, mon, da, hr, ferror )

!=== Open PRINCETON forcing files ===

! If version 3 is being used, set num_met to six because
! the 3-hr precip fields are not available. Turn off
! if future versions include 3-hr precip fields.
if(princeton_struc(n)%version == "3") then
num_met = 6
endif

! Loop over forecast index:
do kk= princeton_struc(n)%st_iterid, princeton_struc(n)%en_iterid

Expand All @@ -202,7 +212,7 @@ subroutine read_princeton( order, n, findex, yr, mon, da, hr, ferror )
endif

! Forcing variable loop:
do v = 1, 7 ! N_PF
do v = 1, num_met ! Number of met fields in Princeton data

! File name for data year/variable(v)_3hourly_year-year.nc
infile=trim(princeton_struc(n)%princetondir)//'/'//cyr//'/'//&
Expand All @@ -223,7 +233,7 @@ subroutine read_princeton( order, n, findex, yr, mon, da, hr, ferror )
if(LIS_masterproc) write(LIS_logunit,*)'[INFO] Opened file: ',infile
end if

if (princeton_struc(n)%version=="2.2")then
if (princeton_struc(n)%version=="2.2" .OR. princeton_struc(n)%version=="3")then
status = nf90_get_var(ncid, varid, datain, &
start=(/1,1,timestep/), &
count=(/princeton_struc(n)%ncold,princeton_struc(n)%nrold,1/))
Expand Down Expand Up @@ -389,8 +399,17 @@ subroutine read_princeton( order, n, findex, yr, mon, da, hr, ferror )
do i = 1, LIS_rc%lnc(n)
if(LIS_domain(n)%gindex(i,j).ne.-1) then
if ((tg(i,j) .lt. 0) .and. (tg(i,j) .ne. LIS_rc%udef)) then
write(LIS_logunit,*)'[ERR] No nearest neighbors, v, i, j',v,i,j,tg(i,j)
stop
!For version 3, replace the value with an undefined value instead of ending the run.
!For all other versions, keep old call to end the program.
if(princeton_struc(n)%version == "3") then
write(LIS_logunit,*)'[WARN] No nearest neighbors, v, i, j',v,i,j,tg(i,j)
write(LIS_logunit,*)'[WARN] Check output to make sure the missing neighbor'
write(LIS_logunit,*)'is isolated and does not distort the output.'
tg(i,j) = LIS_rc%udef ! New code to change data to undefined
else
write(LIS_logunit,*)'[WARN] No nearest neighbors, v, i, j',v,i,j,tg(i,j)
call LIS_endrun() ! Old code to call a fatal error
endif
endif
templdas(i,j,v) = tg(i,j)
endif
Expand Down Expand Up @@ -445,6 +464,7 @@ end subroutine read_princeton
!
! !INTERFACE:
subroutine princetongrid_2_lisgrid( nx, ny, grid_data )
use LIS_logMod, only : LIS_logunit, LIS_endrun

implicit none
! !ARGUMENTS:
Expand All @@ -471,17 +491,18 @@ subroutine princetongrid_2_lisgrid( nx, ny, grid_data )
! ------------------------------------------------------------------
! Some checks
! ------------------------------------------------------------------
if ((nx /= 360) .or. (ny /= 180)) then
write (*,*) '[ERR] Rprincetongrid_2_gldasgrid(): This routine has only been'
write (*,*) ' checked for nx=360 and ny=180. Make sure you know'
write (*,*) ' what you are doing. STOPPING.'
stop
if( ((nx /= 360) .or. (ny /= 180)) .and. ((nx /= 1440) .or. (ny /= 600)) ) then
write(LIS_logunit,*) '[ERR] Rprincetongrid_2_gldasgrid(): This routine has only been'
write(LIS_logunit,*) ' checked for nx=360 and ny=180 (version 2 and 2.2) and for '
write(LIS_logunit,*) ' nx=1440 and ny=600 (version 3). '
write(LIS_logunit,*) ' Make sure you know what you are doing. STOPPING.'
call LIS_endrun()
end if
if ((mod(nx,2) /= 0) .or. (mod(ny,2) /= 0)) then
write (*,*) '[ERR] Rprincetongrid_2_gldasgrid(): This routine can only work'
write (*,*) ' for even nx and ny. Make sure you know'
write (*,*) ' what you are doing. STOPPING.'
stop
write(LIS_logunit,*) '[ERR] Rprincetongrid_2_gldasgrid(): This routine can only work'
write(LIS_logunit,*) ' for even nx and ny. Make sure you know'
write(LIS_logunit,*) ' what you are doing. STOPPING.'
call LIS_endrun()
end if

!-------------------------------------------------------------------
Expand Down
22 changes: 21 additions & 1 deletion lis/metforcing/princeton/readcrd_princeton.F90
Original file line number Diff line number Diff line change
Expand Up @@ -14,13 +14,14 @@
! 26 Jan 2007; Hiroko Kato, Initial Code adopted from readcrd_princeton.F90
! 25 Jun 2007; Hiroko Kato, upgraded to LISv5.0
! 15 May 2017: Bailing Li; Added changes for reading in version 2.2 data
! 22 Oct 2018: Daniel Sarmiento; Added changes to support version 3 data
!
! !INTERFACE:
subroutine readcrd_princeton()
! !USES:
use ESMF
use LIS_coreMod, only : LIS_rc, LIS_config
use LIS_logMod, only : LIS_logunit
use LIS_logMod, only : LIS_logunit, LIS_endrun, LIS_verify
use princeton_forcingMod, only : princeton_struc
! !DESCRIPTION:
!
Expand All @@ -35,16 +36,35 @@ subroutine readcrd_princeton()
write(LIS_logunit,*)'[INFO] Using PRINCETON forcing'

call ESMF_ConfigFindLabel(LIS_config,"PRINCETON forcing directory:",rc=rc)
call LIS_verify(rc, 'PRINCETON forcing directory: not defined')

do n=1,LIS_rc%nnest
call ESMF_ConfigGetAttribute(LIS_config,princeton_struc(n)%princetondir,rc=rc)
enddo

call ESMF_ConfigFindLabel(LIS_config,"PRINCETON forcing version:",rc=rc)
call LIS_verify(rc, 'PRINCETON forcing version: not defined')

do n=1,LIS_rc%nnest
call ESMF_ConfigGetAttribute(LIS_config,princeton_struc(n)%version,rc=rc)
if (princeton_struc(n)%version .ne. "2" .AND. princeton_struc(n)%version .ne. "2.2" .AND. princeton_struc(n)%version .ne. "3") then
write(LIS_logunit,*) "[ERR] The current version of the Princeton LIS reader only"
write(LIS_logunit,*) "[ERR] supports versions 2 (1 deg), 2.2 (1 deg), or 3 (0.25 deg)."
write(LIS_logunit,*) "[ERR] Please input a valid version in the configuration file."
call LIS_endrun()
endif
if (princeton_struc(n)%version .eq. "3") then
write(LIS_logunit,*) "[WARN] The version 3 of the Princeton driver data is not global."
write(LIS_logunit,*) "[WARN] The southern latitude boundary is at -59.875 deg. LIS has"
write(LIS_logunit,*) "[WARN] not been tested when the domain defined in LDT extends beyond"
write(LIS_logunit,*) "[WARN] this boundary. Redefine the domain in LDT if you encounter"
write(LIS_logunit,*) "[WARN] any issues."
endif
enddo

do n=1,LIS_rc%nnest
write(LIS_logunit,*)'[INFO] PRINCETON forcing directory :',princeton_struc(n)%princetonDIR
write(LIS_logunit,*)'[INFO] PRINCETON forcing version :',princeton_struc(n)%version
princeton_struc(n)%princetontime1 = 3000.0
princeton_struc(n)%princetontime2 = 0.0
enddo
Expand Down