Skip to content

Commit

Permalink
Add missing mpp_land_bcast_real8() and read_rst_crt_nc8() subroutines
Browse files Browse the repository at this point in the history
  • Loading branch information
rcabell committed Sep 22, 2022
1 parent 2f91cb7 commit c6b8a46
Show file tree
Hide file tree
Showing 2 changed files with 56 additions and 7 deletions.
11 changes: 11 additions & 0 deletions src/MPP/mpp_land.F
Original file line number Diff line number Diff line change
Expand Up @@ -1096,6 +1096,17 @@ subroutine mpp_land_bcast_real(size1,inout)
return
end subroutine mpp_land_bcast_real

subroutine mpp_land_bcast_real8(size1,inout)
integer size1
! real inout(size1)
real*8 , dimension(:) :: inout
integer ierr, len
call mpi_bcast(inout,size1,MPI_REAL8, &
IO_id,HYDRO_COMM_WORLD,ierr)
call mpp_land_sync()
return
end subroutine mpp_land_bcast_real8

subroutine mpp_land_bcast_int2d(inout)
integer length1, k,length2
integer inout(:,:)
Expand Down
52 changes: 45 additions & 7 deletions src/Routing/module_HYDRO_io.F
Original file line number Diff line number Diff line change
Expand Up @@ -6700,7 +6700,7 @@ subroutine read_rst_gwbucket_real(ncid,outV,numbasns,&
character(len=*) :: vName
integer i, j,k
real, dimension(gnumbasns) :: buf
call read_rst_crt_nc(ncid,dble(buf),gnumbasns,vName)
call read_rst_crt_nc(ncid,buf,gnumbasns,vName)
do k = 1, numbasns
outV(k) = buf(basnsInd(k))
end do
Expand Down Expand Up @@ -6875,9 +6875,9 @@ subroutine RESTART_IN_NC(inFile,did)
endif
if(rt_domain(did)%NLAKES .gt. 0) then
call read_rst_crt_nc(ncid,rt_domain(did)%RESHT,rt_domain(did)%NLAKES,"resht")
call read_rst_crt_nc(ncid,rt_domain(did)%QLAKEO,rt_domain(did)%NLAKES,"qlakeo")
call read_rst_crt_nc(ncid,rt_domain(did)%QLAKEI,rt_domain(did)%NLAKES,"qlakei")
call read_rst_crt_nc8(ncid,rt_domain(did)%RESHT,rt_domain(did)%NLAKES,"resht")
call read_rst_crt_nc8(ncid,rt_domain(did)%QLAKEO,rt_domain(did)%NLAKES,"qlakeo")
call read_rst_crt_nc8(ncid,rt_domain(did)%QLAKEI,rt_domain(did)%NLAKES,"qlakei")
endif
if( nlst(did)%channel_only .eq. 0 .and. &
Expand Down Expand Up @@ -7152,7 +7152,7 @@ subroutine read_rst_crt_nc(ncid,var,n,varStr)
implicit none
integer :: ireg, ncid, varid, n, iret
!real,dimension(n) :: var
real*8,dimension(n) :: var
real,dimension(n) :: var
character(len=*) :: varStr
if( n .le. 0) return
Expand All @@ -7179,14 +7179,52 @@ subroutine read_rst_crt_nc(ncid,var,n,varStr)
#ifdef MPP_LAND
endif
if(n .gt. 0) then
call mpp_land_bcast(var)
!call mpp_land_bcast_real(n,var)
call mpp_land_bcast_real(n,var)
endif
#endif
return
end subroutine read_rst_crt_nc
subroutine read_rst_crt_nc8(ncid,var,n,varStr)
implicit none
integer :: ireg, ncid, varid, n, iret
!real,dimension(n) :: var
real*8,dimension(n) :: var
character(len=*) :: varStr
if( n .le. 0) return
#ifdef MPP_LAND
if(my_id .eq. IO_id) &
#endif
iret = nf90_inq_varid(ncid, trim(varStr), varid)
#ifdef MPP_LAND
call mpp_land_bcast_int1(iret)
#endif
if (iret /= 0) then
#ifdef HYDRO_D
print*, 'variable not found: name = "', trim(varStr)//'"'
#endif
return
endif
#ifdef HYDRO_D
print*, "read restart variable ", varStr
#endif
#ifdef MPP_LAND
if(my_id .eq. IO_id) then
#endif
iret = nf90_get_var(ncid, varid, var)
#ifdef MPP_LAND
endif
if(n .gt. 0) then
call mpp_land_bcast_real8(n,var)
endif
#endif
return
end subroutine read_rst_crt_nc8
subroutine read_rst_crt_stream_nc(ncid,var_out,n,varStr,gnlinks,map_l2g)
implicit none
integer :: ncid, varid, n, iret, gnlinks
Expand Down

0 comments on commit c6b8a46

Please sign in to comment.