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

Remove pftmask from history and PFTDATA_MASK from fsurdat files #1866

Merged
merged 1 commit into from
Oct 8, 2022
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
12 changes: 0 additions & 12 deletions src/main/histFileMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -3346,17 +3346,6 @@ subroutine htape_timeconst(t, mode)
long_name='land/ocean mask (0.=ocean and 1.=land)', ncid=nfid(t), &
imissing_value=ispval, ifill_value=ispval)
end if
if (ldomain%isgrid2d) then
call ncd_defvar(varname='pftmask' , xtype=ncd_int, &
dim1name='lon', dim2name='lat', &
long_name='pft real/fake mask (0.=fake and 1.=real)', ncid=nfid(t), &
imissing_value=ispval, ifill_value=ispval)
else
call ncd_defvar(varname='pftmask' , xtype=ncd_int, &
dim1name=grlnd, &
long_name='pft real/fake mask (0.=fake and 1.=real)', ncid=nfid(t), &
imissing_value=ispval, ifill_value=ispval)
end if
if (ldomain%isgrid2d) then
call ncd_defvar(varname='nbedrock' , xtype=ncd_int, &
dim1name='lon', dim2name='lat', &
Expand Down Expand Up @@ -3384,7 +3373,6 @@ subroutine htape_timeconst(t, mode)
call ncd_io(varname='area' , data=ldomain%area, dim1name=grlnd, ncid=nfid(t), flag='write')
call ncd_io(varname='landfrac', data=ldomain%frac, dim1name=grlnd, ncid=nfid(t), flag='write')
call ncd_io(varname='landmask', data=ldomain%mask, dim1name=grlnd, ncid=nfid(t), flag='write')
call ncd_io(varname='pftmask' , data=ldomain%pftm, dim1name=grlnd, ncid=nfid(t), flag='write')
call ncd_io(varname='nbedrock' , data=grc%nbedrock, dim1name=grlnd, ncid=nfid(t), flag='write')

end if ! (define/write mode
Expand Down
5 changes: 1 addition & 4 deletions tools/mksurfdata_esmf/src/mkfileMod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -604,10 +604,7 @@ subroutine mkfile_define_vars(pioid, dynlanduse)
end if

call mkpio_def_spatial_var(pioid=pioid, varname='LANDFRAC_PFT', xtype=xtype, &
long_name='land fraction from pft dataset', units='unitless')

call mkpio_def_spatial_var(pioid=pioid, varname='PFTDATA_MASK', xtype=PIO_INT, &
long_name='land mask from pft dataset, indicative of real/fake points', units='unitless')
long_name='land fraction from pft dataset (DIFF FROM landfrac USED IN SIMULATION, SHOWN IN HISTORY)', units='unitless')

if (.not. dynlanduse) then
call mkpio_def_spatial_var(pioid=pioid, varname='PCT_NATVEG', xtype=xtype, &
Expand Down
24 changes: 2 additions & 22 deletions tools/mksurfdata_esmf/src/mksurfdata.F90
Original file line number Diff line number Diff line change
Expand Up @@ -145,7 +145,6 @@ program mksurfdata
! pct vegetation data
real(r8), allocatable :: landfrac_pft(:) ! PFT data: % land per gridcell
real(r8), allocatable :: pctlnd_pft(:) ! PFT data: % of gridcell for PFTs
integer , allocatable :: pftdata_mask(:) ! mask indicating real or fake land type
type(pct_pft_type), allocatable :: pctnatpft(:) ! % of grid cell that is nat veg, and breakdown into PFTs
type(pct_pft_type), allocatable :: pctcft(:) ! % of grid cell that is crop, and breakdown into CFTs

Expand Down Expand Up @@ -376,13 +375,12 @@ program mksurfdata

! -----------------------------------
! Make PFTs [pctnatpft, pctcft] from dataset [fvegtyp]
! Make landfrac_pft and pftdata_mask
! Make landfrac_pft
! -----------------------------------
! Determine fractional land from pft dataset
allocate(pctlnd_pft(lsize_o)); pctlnd_pft(:) = spval
allocate(pctnatpft(lsize_o)) ;
allocate(pctcft(lsize_o)) ;
allocate(pftdata_mask(lsize_o)) ; pftdata_mask(:) = -999
allocate(landfrac_pft(lsize_o)) ; landfrac_pft(:) = spval
call mkpft( mksrf_fvegtyp_mesh, mksrf_fvegtyp, mesh_model, &
pctlnd_o=pctlnd_pft, pctnatpft_o=pctnatpft, pctcft_o=pctcft, rc=rc)
Expand All @@ -400,16 +398,8 @@ program mksurfdata
call pctcft(n)%set_pct_l2g(0._r8)
end if
landfrac_pft(n) = pctlnd_pft(n)/100._r8
if (pctlnd_pft(n) < 1.e-6_r8) then
pftdata_mask(n) = 0
else
pftdata_mask(n) = 1
end if
end do
if (fsurdat /= ' ') then
if (root_task) write(ndiag, '(a)') trim(subname)//" writing land mask from pft dataset"
call mkfile_output(pioid, mesh_model, 'PFTDATA_MASK', pftdata_mask, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) call shr_sys_abort('error in calling mkfile_output')
if (root_task) write(ndiag, '(a)') trim(subname)//" writing land fraction from pft dataset"
call mkfile_output(pioid, mesh_model, 'LANDFRAC_PFT', landfrac_pft, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) call shr_sys_abort('error in calling mkfile_output')
Expand Down Expand Up @@ -789,12 +779,6 @@ program mksurfdata
rcode = pio_inq_varid(pioid, 'cft', pio_varid)
rcode = pio_put_var(pioid, pio_varid, (/(n,n=cft_lb,cft_ub)/))

! Write out PFTDATA_MASK
! pftdata_mask was calculated ABOVE
if (root_task) write(ndiag, '(a)') trim(subname)//" writing land mask (calculated in furdata calc)"
call mkfile_output(pioid, mesh_model, 'PFTDATA_MASK', pftdata_mask, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) call shr_sys_abort('error in calling mkfile_output')

! Write out LANDFRAC_PFT
! landfrac_pft was calculated ABOVE
if (root_task) write(ndiag, '(a)') trim(subname)//" writing land fraction calculated in fsurdata calc)"
Expand Down Expand Up @@ -1097,10 +1081,8 @@ subroutine normalize_and_check_landuse(ns_o)

! Assume wetland, glacier and/or lake when dataset landmask implies ocean
! (assume medium soil color (15) and loamy texture).
! Also set pftdata_mask here

if (pctlnd_pft(n) < 1.e-6_r8) then
pftdata_mask(n) = 0
if (pctgla(n) < 1.e-6_r8) then
pctwet(n) = 100._r8 - pctlak(n)
pctgla(n) = 0._r8
Expand All @@ -1110,8 +1092,6 @@ subroutine normalize_and_check_landuse(ns_o)
pcturb(n) = 0._r8
call pctnatpft(n)%set_pct_l2g(0._r8)
call pctcft(n)%set_pct_l2g(0._r8)
else
pftdata_mask(n) = 1
end if

! Make sure sum of all land cover types except natural vegetation does
Expand Down Expand Up @@ -1265,7 +1245,7 @@ subroutine normalize_and_check_landuse(ns_o)

! Make sure that there is no vegetation outside the pft mask
do n = 1,ns_o
if (pftdata_mask(n) == 0 .and. (pctnatpft(n)%get_pct_l2g() > 0 .or. pctcft(n)%get_pct_l2g() > 0)) then
if (pctlnd_pft(n) < 1.e-6_r8 .and. (pctnatpft(n)%get_pct_l2g() > 0 .or. pctcft(n)%get_pct_l2g() > 0)) then
write (6,*)'vegetation found outside the pft mask at n=',n
write (6,*)'pctnatveg,pctcrop=', pctnatpft(n)%get_pct_l2g(), pctcft(n)%get_pct_l2g()
call shr_sys_abort()
Expand Down