Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / hydro / Routing / module_reservoir_routing.F90
blobb2b20b459efb7669a970c2fa562976b45b543227
1 ! Intended purpose is to provide a module for all subroutines related to
2 ! reservoir routing, including active management, level pool, and integrating live 
3 ! data feeds. As of NWMv2.0, this module stub can read in a timeslice file
4 ! to incorporate data from external sources, should a data service become available.
6 ! NOTE: THIS CODE IS NOT ACTIVE AND CANNOT BE ACTIVATED THROUGH
7 !       NAMELIST OPTIONS AT THIS TIME - NWM Version 2.0.
9 ! Logan Karsten
10 ! National Center for Atmospheric Research
11 ! Research Applications Laboratory
12 ! karsten at ucar dot edu
14 module module_reservoir_routing
15 implicit none
17 contains
19 subroutine read_reservoir_obs(domainId)
20   use config_base, only: nlst
21    use netcdf
22    use module_hydro_stop, only: HYDRO_stop
23 #ifdef MPP_LAND
24    use module_mpp_land
25 #endif
26    implicit none
28    ! Pass in domain ID value from parent calling program.
29    integer, intent(in) :: domainId
31    ! Local variables
32    integer :: ftn ! NetCDF file handle.
33    integer :: nLakesNc ! Number of lake objects in the input file.
34    real, allocatable, dimension(:) :: minRelease ! Specified minimum reservoir release (cms)
35    real, allocatable, dimension(:) :: poolElev ! Value of time-varying water surface elevation
36    real, allocatable, dimension(:) :: resFlow ! Reservoir discharge in cms
37    real, allocatable, dimension(:) :: qType ! Reservoir discharge type
38    real*8, allocatable, dimension(:) :: ifd ! Initial fractional depth
39    real*8, allocatable, dimension(:) :: lkArea ! Gridded lake area (sq. km)
40    real*8, allocatable, dimension(:) :: lkMxE ! Maximum lake elevation (m ASL)
41    real*8, allocatable, dimension(:) :: orificeA ! Orifice cross-sectional area (sq. m)
42    real*8, allocatable, dimension(:) :: orificeC ! Orifice coefficient
43    real*8, allocatable, dimension(:) :: orificeE ! Orifice elevation (m ASL)
44    real*8, allocatable, dimension(:) :: weirC ! Weir coefficient
45    real*8, allocatable, dimension(:) :: weirE ! Weir height (m ASL)
46    real*8, allocatable, dimension(:) :: weirL ! Weir length (m)
47    integer, allocatable, dimension(:) :: ascIndex ! Index to sort lake objects by ascending ID
48    integer, allocatable, dimension(:) :: lakeID ! Lake index value
49    real, allocatable, dimension(:) :: lakeLat ! Lake latitude values.
50    real, allocatable, dimension(:) :: lakeLon ! Lake longitude values
51    character(len=1024) :: inFlnm ! NetCDF file with discharge data.
52    integer :: myId ! MPI processor ID
53    integer :: ierr ! MPI return status
54    integer :: diagFlag
55    logical :: file_exists
56    integer :: missingFlag
57    integer :: varTmpId
58    integer :: mppFlag
59    integer :: iret
60    integer :: lakeDimId
62 #ifdef HYDRO_D
63    diagFlag = 1
64 #else
65    diagFlag = 0
66 #endif
68    ! Sync up processes.
69    if(mppFlag .eq. 1) then
70 #ifdef MPP_LAND
71       call mpp_land_sync()
72 #endif
73    endif
75    ! Check to ensure the namelist option for reading in the reservoir discharge data
76    ! has been set to 1. If not, return back to the main calling program.
77    if(nlst(domainId)%reservoir_data_ingest .eq. 0) then
78        ! No reservoir realtime data requested.
79        return
80    endif
82    ! If we are running over MPI, determine which processor number we are on.
83    ! If not MPI, then default to 0, which is the I/O ID.
84    if(mppFlag .eq. 1) then
85 #ifdef MPP_LAND
86       call MPI_COMM_RANK( HYDRO_COMM_WORLD, myId, ierr )
87       call nwmCheck(diagFlag,ierr,'ERROR: Unable to determine MPI process ID.')
88 #endif
89    else
90       myId = 0
91    endif
93    ! Open up and read in the NetCDF file containing disharge data.
94    if(myId .eq. 0) then
95       ! Initialize our missing flag to 0. If at any point we don't find a file, 
96       ! the flag value will go to 1 to indicate no files were found.
97       missingFlag = 0
99       ! We are on the I/O processor.
100       write(inFlnm,'(A,"/LAKEFILE_DISCHARGE_",A12,".nc")') nlst(domainId)%reservoir_obs_dir,nlst(domainId)%olddate(1:4)//&
101             nlst(domainId)%olddate(6:7)//nlst(domainId)%olddate(9:10)//&
102             nlst(domainId)%olddate(12:13)//nlst(domainId)%olddate(15:16)
104       ! Check to see if the file exists.
105       INQUIRE(FILE=inFlnm,EXIST=file_exists)
107       if(file_exists) then
108          iret = nf90_open(trim(inFlnm),NF90_NOWRITE,ncid=ftn)
109          call nwmCheck(diagFlag,iret,"ERROR: Unable to open LAKE reservoir discharge file.")
110          iret = nf90_inq_dimid(ftn,'nlakes',lakeDimId)
111          call nwmCheck(diagFlag,iret,'ERROR: Unable to find nlakes dimension in LAKE reservoir discharge file.')
112          iret = nf90_inquire_dimension(ftn,lakeDimId,len=nLakesNc)
113          call nwmCheck(diagFlag,iret,'ERROR: Unable to find nlakes size in LAKE reservoir discharge file.')
115          ! Allocate the lake variables based on the nlakes dimension
116          allocate(minRelease(nLakesNc))
117          allocate(poolElev(nLakesNc))
118          allocate(resFlow(nLakesNc))
119          allocate(qType(nLakesNc))
120          allocate(ifd(nLakesNc))
121          allocate(lkArea(nLakesNc))
122          allocate(lkMxE(nLakesNc))
123          allocate(orificeA(nLakesNc))
124          allocate(orificeC(nLakesNc))
125          allocate(orificeE(nLakesNc))
126          allocate(weirC(nLakesNc))
127          allocate(weirE(nLakesNc))
128          allocate(weirL(nLakesNc))
129          allocate(ascIndex(nLakesNc))
130          allocate(lakeId(nLakesNc))
131          allocate(lakeLat(nLakesNc))
132          allocate(lakeLon(nLakesNc))
134          ! Read in data.
135          iret = nf90_inq_varid(ftn,'MIN_RELEASE',varTmpId)
136          call nwmCheck(diagFlag,iret,'ERROR: Unable to find MIN_RELEASE in LAKE reservoir discharge file.')
137          iret = nf90_get_var(ftn,varTmpId,minRelease)
138          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract MIN_RELEASE in LAKE reservoir discharge file.')
140          iret = nf90_inq_varid(ftn,'POOL_ELEV',varTmpId)
141          call nwmCheck(diagFlag,iret,'ERROR: Unable to find POOL_ELEV in LAKE reservoir discharge file.')
142          iret = nf90_get_var(ftn,varTmpId,poolElev)
143          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract POOL_ELEV in LAKE reservoir discharge file.')
145          iret = nf90_inq_varid(ftn,'RES_FLOW',varTmpId)
146          call nwmCheck(diagFlag,iret,'ERROR: Unable to find RES_FLOW in LAKE reservoir discharge file.')
147          iret = nf90_get_var(ftn,varTmpId,resFlow)
148          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract RES_FLOW in LAKE reservoir discharge file.')
150          iret = nf90_inq_varid(ftn,'Q_TYPE',varTmpId)
151          call nwmCheck(diagFlag,iret,'ERROR: Unable to find Q_TYPE in LAKE reservoir discharge file.')
152          iret = nf90_get_var(ftn,varTmpId,qType)
153          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract Q_TYPE in LAKE reservoir discharge file.')
155          iret = nf90_inq_varid(ftn,'ifd',varTmpId)
156          call nwmCheck(diagFlag,iret,'ERROR: Unable to find ifd in LAKE reservoir discharge file.')
157          iret = nf90_get_var(ftn,varTmpId,ifd)
158          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract ifd in LAKE reservoir discharge file.')
160          iret = nf90_inq_varid(ftn,'LkArea',varTmpId)
161          call nwmCheck(diagFlag,iret,'ERROR: Unable to find lkArea in LAKE reservoir discharge file.')
162          iret = nf90_get_var(ftn,varTmpId,lkArea)
163          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract lkArea in LAKE reservoir discharge file.')
165          iret = nf90_inq_varid(ftn,'LkMxE',varTmpId)
166          call nwmCheck(diagFlag,iret,'ERROR: Unable to find LkMxE in LAKE reservoir discharge file.')
167          iret = nf90_get_var(ftn,varTmpId,lkMxE)
168          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract lkMxE in LAKE reservoir discharge file.')
170          iret = nf90_inq_varid(ftn,'OrificeA',varTmpId)
171          call nwmCheck(diagFlag,iret,'ERROR: Unable to find OrificeA in LAKE reservoir discharge file.')
172          iret = nf90_get_var(ftn,varTmpId,orificeA)
173          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract OrificeA in LAKE reservoir discharge file.')
175          iret = nf90_inq_varid(ftn,'OrificeC',varTmpId)
176          call nwmCheck(diagFlag,iret,'ERROR: Unable to find OrificeC in LAKE reservoir discharge file.')
177          iret = nf90_get_var(ftn,varTmpId,orificeC)
178          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract OrificeC in LAKE reservoir discharge file.')
180          iret = nf90_inq_varid(ftn,'OrificeE',varTmpId)
181          call nwmCheck(diagFlag,iret,'ERROR: Unable to find OrificeE in LAKE reservoir discharge file.')
182          iret = nf90_get_var(ftn,varTmpId,orificeE)
183          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract OrificeE in LAKE reservoir discharge file.')
185          iret = nf90_inq_varid(ftn,'WeirC',varTmpId)
186          call nwmCheck(diagFlag,iret,'ERROR: Unable to find WeirC in LAKE reservoir discharge file.')
187          iret = nf90_get_var(ftn,varTmpId,weirC)
188          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract WeirC in LAKE reservoir discharge file.')
190          iret = nf90_inq_varid(ftn,'WeirE',varTmpId)
191          call nwmCheck(diagFlag,iret,'ERROR: Unable to find WeirE in LAKE reservoir discharge file.')
192          iret = nf90_get_var(ftn,varTmpId,weirE)
193          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract WeirE in LAKE reservoir discharge file.')
195          iret = nf90_inq_varid(ftn,'WeirL',varTmpId)
196          call nwmCheck(diagFlag,iret,'ERROR: Unable to find WeirL in LAKE reservoir discharge file.')
197          iret = nf90_get_var(ftn,varTmpId,weirL)
198          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract WeirL in LAKE reservoir discharge file.')
200          iret = nf90_inq_varid(ftn,'ascendingIndex',varTmpId)
201          call nwmCheck(diagFlag,iret,'ERROR: Unable to find ascendingIndex in LAKE reservoir discharge file.')
202          iret = nf90_get_var(ftn,varTmpId,ascIndex)
203          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract ascendingIndex in LAKE reservoir discharge file.')
205          iret = nf90_inq_varid(ftn,'lake_id',varTmpId)
206          call nwmCheck(diagFlag,iret,'ERROR: Unable to find lake_id in LAKE reservoir discharge file.')
207          iret = nf90_get_var(ftn,varTmpId,lakeId)
208          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract lake_id in LAKE reservoir discharge file.')
210          iret = nf90_inq_varid(ftn,'lat',varTmpId)
211          call nwmCheck(diagFlag,iret,'ERROR: Unable to find lat in LAKE reservoir discharge file.')
212          iret = nf90_get_var(ftn,varTmpId,lakeLat)
213          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract lat in LAKE reservoir discharge file.')
215          iret = nf90_inq_varid(ftn,'lon',varTmpId)
216          call nwmCheck(diagFlag,iret,'ERROR: Unable to find lon in LAKE reservoir discharge file.')
217          iret = nf90_get_var(ftn,varTmpId,lakeLon)
218          call nwmCheck(diagFlag,iret,'ERROR: Unable to extract lon in LAKE reservoir discharge file.')
220          ! Close the NetCDF file
221          iret = nf90_close(ftn)
222          call nwmCheck(diagFlag,iret,'ERROR: Unable to close LAKE reservoir discharge file.')
224          ! Deallocate memory appropriately
225          deallocate(minRelease)
226          deallocate(poolElev)
227          deallocate(resFlow)
228          deallocate(qType)
229          deallocate(ifd)
230          deallocate(lkArea)
231          deallocate(lkMxE)
232          deallocate(orificeA)
233          deallocate(orificeC)
234          deallocate(orificeE)
235          deallocate(weirC)
236          deallocate(weirE)
237          deallocate(weirL)
238          deallocate(ascIndex)
239          deallocate(lakeId)
240          deallocate(lakeLat)
241          deallocate(lakeLon)
243       else
244          missingFlag = 1
245          return
246       endif
247    endif
249 end subroutine read_reservoir_obs
251 subroutine postDiagMsg(diagFlag,diagMsg)
252    implicit none
254    ! Subroutine arguments.
255    integer, intent(in) :: diagFlag
256    character(len=*), intent(in) :: diagMsg
258    ! Only write out message if the diagnostic WRF_HYDRO_D flag was
259    ! set to 1
260    if (diagFlag .eq. 1) then
261       print*, trim(diagMsg)
262    end if
263 end subroutine postDiagMsg
265 subroutine nwmCheck(diagFlag,iret,msg)
266    implicit none
268    ! Subroutine arguments.
269    integer, intent(in) :: diagFlag,iret
270    character(len=*), intent(in) :: msg
272    ! Check status. If status of command is not 0, then post the error message
273    ! if WRF_HYDRO_D was set to be 1.
274    if (iret .ne. 0) then
275       call hydro_stop(trim(msg))
276    end if
278 end subroutine nwmCheck
281 end module module_reservoir_routing