20 #define ADD_ATT_FOR_NEW_VAR
25 CHARACTER(len=256) :: sfcdata_path
26 INTEGER :: cs_res, ncsmp, ncscp, i
27 INTEGER :: res_x, res_y
29 INTEGER*1,
ALLOCATABLE :: lakestatus(:)
30 INTEGER*2,
ALLOCATABLE :: lakedepth(:)
31 REAL,
ALLOCATABLE :: cs_grid(:,:)
32 REAL,
ALLOCATABLE :: cs_lakestatus(:), cs_lakedepth(:)
33 REAL,
ALLOCATABLE :: src_grid_lon(:), src_grid_lat(:)
35 INTEGER :: tile_req, tile_beg, tile_end, binary_lake
38 INTEGER,
PARAMETER :: nlat = 21600, nlon = 43200
39 REAL,
PARAMETER :: d2r = acos(-1.0) / 180.0
40 REAL,
PARAMETER :: r2d = 180.0 /acos(-1.0)
41 REAL,
PARAMETER :: pi = acos(-1.0)
42 REAL*8,
PARAMETER :: gppdeg = 120.0
43 REAL*8,
PARAMETER :: delta = 1.0 / 120.0
45 CHARACTER(len=32) :: arg, lakestatus_srce, lakedepth_srce
46 CHARACTER(len=256) :: lakedata_path
50 IF (iargc() /= 5 .AND. iargc() /= 7)
THEN
51 print*,
'Usage: ', trim(arg), &
52 ' [tile_num (0:all tiles, 7:regional)] [resolution (48,96, ...)] &
53 [lake data path] [lake status source] [lake depth source]'
54 print*,
'Or: ', trim(arg), &
55 ' [tile_num (0:all tiles, 7:regional)] [resolution (48,96, ...)] &
56 [lake data path] [lake status source] [lake depth source] [lake_cutoff] [binary_lake]'
60 READ(arg,*,iostat=stat) tile_req
62 READ(arg,*,iostat=stat) cs_res
63 CALL getarg(3, lakedata_path)
64 CALL getarg(4, lakestatus_srce)
65 CALL getarg(5, lakedepth_srce)
67 IF (iargc() == 5)
THEN
72 READ(arg,*,iostat=stat) lake_cutoff
74 READ(arg,*,iostat=stat) binary_lake
77 print*,
'lake status source:', trim(lakestatus_srce)
78 print*,
'lake depth source:', trim(lakedepth_srce)
79 print*,
'lake cutoff:', lake_cutoff
80 print*,
'binary lake:', binary_lake
82 IF (tile_req == 0)
THEN
83 tile_beg = 1; tile_end = 6
84 print*,
'Process tile 1 - 6 at resolution C',cs_res
85 ELSE IF (tile_req /= 7)
THEN
86 tile_beg = tile_req; tile_end = tile_req
87 print*,
'Process tile',tile_req,
'at resolution C',cs_res
89 tile_beg = 1; tile_end = 1
90 print*,
'Process regional tile (tile', tile_req,
') at resolution C',cs_res
95 ncsmp = (2*cs_res+1)*(2*cs_res+1)*6
96 IF (tile_req /= 7)
THEN
97 print*,
'Read in cubed sphere grid information ... ',ncsmp,
'pairs of lat/lons'
100 IF (tile_req /= 7)
THEN
101 ALLOCATE(cs_grid(ncsmp, 2))
108 ALLOCATE(src_grid_lon(nlon), src_grid_lat(nlat))
110 IF (lakestatus_srce ==
"GLDBV2" .OR. lakestatus_srce ==
"GLDBV3")
THEN
113 src_grid_lon(i) = -180.0 + delta * i
116 src_grid_lat(i) = 90.0 - delta * i
120 IF (lakestatus_srce ==
"MODISP" .OR. lakestatus_srce ==
"VIIRS")
THEN
123 src_grid_lon(i) = -180.0 + delta * (i-1)
126 src_grid_lat(i) = 90.0 - delta * (i-1)
131 lakedata_path = trim(lakedata_path) //
"/"
132 ALLOCATE(lakestatus(nlon*nlat),lakedepth(nlon*nlat))
133 print*,
'Read in lake data file ...'
134 CALL
read_lakedata(lakedata_path,lakestatus,lakedepth,nlat,nlon)
137 ncscp = cs_res*cs_res*6
138 ALLOCATE(cs_lakestatus(ncscp))
139 ALLOCATE(cs_lakedepth(ncscp))
141 print*,
'Calculate lake fraction and average depth for cubed-sphere cells ...'
146 IF (tile_req /= 7)
THEN
147 print*,
'Write lake fraction/depth on cubed sphere grid cells to netCDF files ...'
150 print*,
'Write lake fraction/depth on regional FV3 grid cells to a netCDF file ...'
154 DEALLOCATE(cs_lakestatus,cs_lakedepth)
156 DEALLOCATE(lakestatus,lakedepth)
157 DEALLOCATE(src_grid_lat, src_grid_lon)
171 INTEGER*1,
INTENT(IN) :: lakestat(:)
172 INTEGER*2,
INTENT(IN) :: lakedpth(:)
173 REAL,
INTENT(OUT) :: cs_lakestat(:), cs_lakedpth(:)
175 REAL*8 lolf(2), lort(2), uplf(2), uprt(2), sd_ltmn(4), sd_ltmx(4)
176 REAL*8 :: v(2,4), p(2)
177 REAL :: latmin1, latmax1
178 REAL :: latmin, latmax, lonmin, lonmax, lontmp, lat_sz_max, lon_sz_max
179 INTEGER :: tile_num, i, j, gp, row, col, cs_grid_idx, cs_data_idx
180 INTEGER :: sidex_res, sidey_res, sidex_sz, sidey_sz
181 INTEGER :: stride_lat, stride_lon
182 INTEGER :: src_grid_lat_beg,src_grid_lat_end,src_grid_lon_beg,src_grid_lon_end
183 INTEGER :: src_grid_lon_beg1,src_grid_lon_end1,src_grid_lon_beg2,src_grid_lon_end2
184 INTEGER :: grid_ct, lake_ct, co_gc, tmp
188 REAL*8 :: lake_dpth_sum, lake_avg_frac
191 IF (tile_req /= 7)
THEN
192 sidex_res = cs_res; sidey_res = cs_res
194 sidex_res = res_x; sidey_res = res_y
197 sidex_sz = 2*sidex_res+1; sidey_sz = 2*sidey_res+1
206 DO tile_num = tile_beg, tile_end
207 row = 2 + sidex_sz*(tile_num-1); col = 2
208 DO gp = 1, sidex_res*sidey_res
209 two_section = .false.
210 cs_grid_idx = (row-1)*sidex_sz+col
211 cs_data_idx = (tile_num-1)*sidex_res*sidey_res+gp
212 IF (abs(cs_grid(cs_grid_idx,1)) > 80.0 )
THEN
213 cs_lakestat(cs_data_idx) = 0.0
214 cs_lakedpth(cs_data_idx) = 0.0
217 IF (col > sidex_sz)
THEN
224 lolf(1) = cs_grid(cs_grid_idx-sidex_sz-1, 1)
225 lolf(2) = cs_grid(cs_grid_idx-sidex_sz-1, 2)
226 IF (lolf(2) > 180.0) lolf(2) = lolf(2) - 360.0
227 lort(1) = cs_grid(cs_grid_idx-sidex_sz+1, 1)
228 lort(2) = cs_grid(cs_grid_idx-sidex_sz+1, 2)
229 IF (lort(2) > 180.0) lort(2) = lort(2) - 360.0
230 uplf(1) = cs_grid(cs_grid_idx+sidex_sz-1,1)
231 uplf(2) = cs_grid(cs_grid_idx+sidex_sz-1,2)
232 IF (uplf(2) > 180.0) uplf(2) = uplf(2) - 360.0
233 uprt(1) = cs_grid(cs_grid_idx+sidex_sz+1,1)
234 uprt(2) = cs_grid(cs_grid_idx+sidex_sz+1,2)
236 v(1,1) = lolf(1); v(2,1) = lolf(2)
237 v(1,2) = lort(1); v(2,2) = lort(2)
238 v(1,3) = uprt(1); v(2,3) = uprt(2)
239 v(1,4) = uplf(1); v(2,4) = uplf(2)
240 v(:,:) = v(:,:) * d2r
242 IF (uprt(2) > 180.0) uprt(2) = uprt(2) - 360.0
245 CALL
find_limit(lolf, lort, sd_ltmn(1), sd_ltmx(1))
246 CALL
find_limit(lort, uprt, sd_ltmn(2), sd_ltmx(2))
247 CALL
find_limit(uprt, uplf, sd_ltmn(3), sd_ltmx(3))
248 CALL
find_limit(uplf, lolf, sd_ltmn(4), sd_ltmx(4))
249 latmin = min(sd_ltmn(1),min(sd_ltmn(2),min(sd_ltmn(3),sd_ltmn(4))))
250 latmax = max(sd_ltmx(1),max(sd_ltmx(2),max(sd_ltmx(3),sd_ltmx(4))))
252 latmin = min(lolf(1),min(lort(1),min(uplf(1),uprt(1))))
253 latmax = max(lolf(1),max(lort(1),max(uplf(1),uprt(1))))
254 lonmin = min(lolf(2),min(lort(2),min(uplf(2),uprt(2))))
255 lonmax = max(lolf(2),max(lort(2),max(uplf(2),uprt(2))))
259 src_grid_lat_beg = nint((90.0-latmax)*gppdeg+0.5)
260 src_grid_lat_end = nint((90.0-latmin)*gppdeg+0.5)
261 src_grid_lon_beg = nint((180.0+lonmin)*gppdeg+0.5)
262 src_grid_lon_end = nint((180.0+lonmax)*gppdeg+0.5)
264 IF (src_grid_lat_beg > src_grid_lat_end)
THEN
265 tmp = src_grid_lat_beg
266 src_grid_lat_beg = src_grid_lat_end
267 src_grid_lat_end = tmp
269 IF (src_grid_lon_beg > src_grid_lon_end)
THEN
270 tmp = src_grid_lon_beg
271 src_grid_lon_beg = src_grid_lon_end
272 src_grid_lon_end = tmp
274 IF ((src_grid_lon_end - src_grid_lon_beg) > nlon*0.75)
THEN
276 src_grid_lon_beg1 = src_grid_lon_end
277 src_grid_lon_end1 = nlon
278 src_grid_lon_beg2 = 1
279 src_grid_lon_end2 = src_grid_lon_beg
282 #ifdef DIAG_N_VERBOSE
283 print*,
'cell centre lat/lon =', &
284 gp, cs_res*cs_res, cs_grid(cs_grid_idx,1),cs_grid(cs_grid_idx,2)
285 print*,
'lat index range and stride', &
286 src_grid_lat_beg,src_grid_lat_end,stride_lat
287 print*,
'lat range ', &
288 src_grid_lat(src_grid_lat_beg),src_grid_lat(src_grid_lat_end)
290 lake_ct = 0; grid_ct = 0
293 DO j = src_grid_lat_beg, src_grid_lat_end, stride_lat
294 stride_lon = int(1.0/cos(src_grid_lat(j)*d2r)*
REAL(stride_lat))
295 #ifdef DIAG_N_VERBOSE
296 IF (j == src_grid_lat_beg)
THEN
297 print*,
'lon index range and stride', &
298 src_grid_lon_beg,src_grid_lon_end,stride_lon
299 print*,
'lon range ', &
300 src_grid_lon(src_grid_lon_beg),src_grid_lon(src_grid_lon_end)
301 IF (two_section .eqv. .true.)
THEN
302 print*,
'section1 index lon range and stride', &
303 src_grid_lon_beg1,src_grid_lon_end1,stride_lon
304 print*,
'section1 lon range ', &
305 src_grid_lon(src_grid_lon_beg1),src_grid_lon(src_grid_lon_end1)
306 print*,
'section2 index lon range and stride', &
307 src_grid_lon_beg2,src_grid_lon_end2,stride_lon
308 print*,
'section2 lon range ', &
309 src_grid_lon(src_grid_lon_beg2),src_grid_lon(src_grid_lon_end2)
313 IF (two_section .eqv. .false.)
THEN
314 DO i = src_grid_lon_beg, src_grid_lon_end, stride_lon
315 p(1) = src_grid_lat(j); p(2) = src_grid_lon(i)
319 lkst = lakestat((j-1)*nlon+i); lkdp = lakedpth((j-1)*nlon+i)
320 CALL
lake_cell_comp(lkst, lkdp, lake_ct, lake_avg_frac, lake_dpth_sum)
324 DO i = src_grid_lon_beg1, src_grid_lon_end1, stride_lon
325 p(1) = src_grid_lat(j); p(2) = src_grid_lon(i)
329 lkst = lakestat((j-1)*nlon+i); lkdp = lakedpth((j-1)*nlon+i)
330 CALL
lake_cell_comp(lkst, lkdp, lake_ct, lake_avg_frac, lake_dpth_sum)
333 DO i = src_grid_lon_beg2, src_grid_lon_end2, stride_lon
334 p(1) = src_grid_lat(j); p(2) = src_grid_lon(i)
338 lkst = lakestat((j-1)*nlon+i); lkdp = lakedpth((j-1)*nlon+i)
339 CALL
lake_cell_comp(lkst, lkdp, lake_ct, lake_avg_frac, lake_dpth_sum)
344 IF (lakestatus_srce ==
"GLDBV3" .OR. lakestatus_srce ==
"GLDBV2" .OR. &
345 lakestatus_srce ==
"VIIRS" )
THEN
346 cs_lakestat(cs_data_idx)=
REAL(lake_ct)/
REAL(grid_ct)
348 IF (lakestatus_srce ==
"MODISP")
THEN
349 cs_lakestat(cs_data_idx)=lake_avg_frac/
REAL(grid_ct)
351 IF (lake_ct /= 0)
THEN
352 cs_lakedpth(cs_data_idx)=lake_dpth_sum/
REAL(lake_ct)/10.0
354 cs_lakedpth(cs_data_idx)=0.0
356 #ifdef DIAG_N_VERBOSE
357 print*,
'tile_num, row, col:', tile_num, row, col
358 print*,
'grid_ct, lake_ct = ', grid_ct, lake_ct
359 print*,
'lake_frac= ', cs_lakestat(cs_data_idx)
360 print*,
'lake_depth (avg) = ', cs_lakedpth(cs_data_idx)
365 IF (col > sidex_sz)
THEN
385 INTEGER*1,
INTENT(IN) :: lkst
386 INTEGER*2,
INTENT(IN) :: lkdp
387 INTEGER,
INTENT(OUT) :: lake_ct
388 REAL*8,
INTENT(OUT) :: lake_avg_frac, lake_dpth_sum
392 IF (lakestatus_srce ==
"GLDBV3" .OR. lakestatus_srce ==
"GLDBV2")
THEN
395 lake_dpth_sum = lake_dpth_sum+30.0
397 lake_dpth_sum = lake_dpth_sum+100.0
400 lake_dpth_sum = lake_dpth_sum+
REAL(lkdp)
403 IF (lakestatus_srce ==
"MODISP" .OR. lakestatus_srce ==
"VIIRS")
THEN
404 lake_avg_frac = lake_avg_frac +
REAL(lkst) / 100.0
406 lake_dpth_sum = lake_dpth_sum+100.0
408 lake_dpth_sum = lake_dpth_sum+
REAL(lkdp)
424 INTEGER,
INTENT(IN) :: res
425 REAL,
INTENT(OUT) :: grid(:,:)
427 REAL*8,
ALLOCATABLE :: lat(:), lon(:)
428 INTEGER :: side_sz, tile_sz, ncid, varid
429 INTEGER :: i, tile_beg, tile_end, stat
430 CHARACTER(len=256) :: gridfile_path,gridfile
432 CHARACTER(len=4) res_ch
435 tile_sz = side_sz*side_sz
436 ALLOCATE(lat(tile_sz), lon(tile_sz))
438 IF (tile_req == 0)
THEN
439 tile_beg = 1; tile_end = 6
441 tile_beg = tile_req; tile_end = tile_req
443 WRITE(res_ch,
'(I4)') res
445 DO i = tile_beg, tile_end
447 gridfile = trim(gridfile_path)//
"C"//trim(adjustl(res_ch))//
"_grid.tile"//ich//
".nc"
449 print*,
'Open cubed sphere grid spec file ', trim(gridfile)
451 stat = nf90_open(trim(gridfile), nf90_nowrite, ncid)
454 stat = nf90_inq_varid(ncid,
'x',varid)
456 stat = nf90_get_var(ncid,varid,lon,start=(/1,1/),count=(/side_sz,side_sz/))
457 CALL
nc_opchk(stat,
'nf90_get_var_lon')
458 stat = nf90_inq_varid(ncid,
'y',varid)
460 stat = nf90_get_var(ncid,varid,lat,start=(/1,1/),count=(/side_sz,side_sz/))
461 CALL
nc_opchk(stat,
'nf90_get_var_lat')
462 stat = nf90_close(ncid)
465 tile_beg = (i-1)*tile_sz+1
467 grid(tile_beg:tile_end,1) = lat(1:tile_sz)
468 grid(tile_beg:tile_end,2) = lon(1:tile_sz)
484 INTEGER,
INTENT(IN) :: res, halo_depth
485 INTEGER,
INTENT(OUT) :: res_x, res_y
486 REAL,
ALLOCATABLE,
INTENT(OUT) :: grid(:,:)
488 REAL*8,
ALLOCATABLE :: lat(:), lon(:)
489 INTEGER :: sidex_sz, sidey_sz, tile_sz, ncid, varid, dimid
490 INTEGER :: x_start, y_start
491 INTEGER :: nxp, nyp, stat
492 CHARACTER(len=256) :: gridfile_path,gridfile
494 CHARACTER(len=4) res_ch
495 CHARACTER(len=8) dimname
497 WRITE(res_ch,
'(I4)') res
499 gridfile = trim(gridfile_path)//
"C"//trim(adjustl(res_ch))//
"_grid.tile7.nc"
501 print*,
'Open cubed sphere grid spec file ', trim(gridfile)
503 stat = nf90_open(trim(gridfile), nf90_nowrite, ncid)
506 stat = nf90_inq_dimid(ncid,
'nxp',dimid)
507 CALL
nc_opchk(stat,
'nf90_inq_dimid: nxp')
508 stat = nf90_inquire_dimension(ncid,dimid,dimname,len=nxp)
509 CALL
nc_opchk(stat,
'nf90_inquire_dimension: nxp')
511 stat = nf90_inq_dimid(ncid,
'nyp',dimid)
512 CALL
nc_opchk(stat,
'nf90_inq_dimid: nyp')
513 stat = nf90_inquire_dimension(ncid,dimid,dimname,len=nyp)
514 CALL
nc_opchk(stat,
'nf90_inquire_dimension: nyp')
518 tile_sz = sidex_sz*sidey_sz
519 ALLOCATE(lat(tile_sz), lon(tile_sz))
522 x_start = 1; y_start = 1
523 stat = nf90_inq_varid(ncid,
'x',varid)
525 stat = nf90_get_var(ncid,varid,lon,start=(/x_start,y_start/),count=(/sidex_sz,sidey_sz/))
526 CALL
nc_opchk(stat,
'nf90_get_var_lon')
527 stat = nf90_inq_varid(ncid,
'y',varid)
529 stat = nf90_get_var(ncid,varid,lat,start=(/x_start,y_start/),count=(/sidex_sz,sidey_sz/))
530 CALL
nc_opchk(stat,
'nf90_get_var_lat')
531 stat = nf90_close(ncid)
534 ALLOCATE(grid(tile_sz,2))
535 grid(1:tile_sz,1) = lat(1:tile_sz)
536 grid(1:tile_sz,2) = lon(1:tile_sz)
538 res_x = sidex_sz/2; res_y = sidey_sz/2
554 CHARACTER(len=256),
INTENT(IN) :: lakedata_path
555 INTEGER*1,
INTENT(OUT) :: lake_stat(:)
556 INTEGER*2,
INTENT(OUT) :: lake_dpth(:)
557 INTEGER,
INTENT(IN) :: nlat, nlon
559 CHARACTER(len=256) lakefile
560 INTEGER :: data_sz, i
566 lakefile = trim(lakedata_path) //
"GlobalLakeStatus_GLDBv3release.dat"
567 IF (lakestatus_srce ==
"GLDBV2")
THEN
568 lakefile = trim(lakedata_path) //
"GlobalLakeStatus.dat"
570 IF (lakestatus_srce ==
"GLDBV3")
THEN
571 lakefile = trim(lakedata_path) //
"GlobalLakeStatus_GLDBv3release.dat"
573 IF (lakestatus_srce ==
"MODISP")
THEN
575 lakefile = trim(lakedata_path) //
"GlobalLakeStatus_MODISp.dat"
577 IF (lakestatus_srce ==
"VIIRS")
THEN
578 lakefile = trim(lakedata_path) //
"GlobalLakeStatus_VIIRS.dat"
580 OPEN(10, file=lakefile, form=
'unformatted', access=
'direct', status=
'old', recl=data_sz*1)
581 READ(10,rec=1) lake_stat
585 lakefile = trim(lakedata_path) //
"GlobalLakeDepth_GLDBv3release.dat"
586 IF (lakedepth_srce ==
"GLDBV2")
THEN
587 lakefile = trim(lakedata_path) //
"GlobalLakeDepth.dat"
589 IF (lakedepth_srce ==
"GLDBV3")
THEN
590 lakefile = trim(lakedata_path) //
"GlobalLakeDepth_GLDBv3release.dat"
592 IF (lakedepth_srce ==
"GLOBATHY")
THEN
593 lakefile = trim(lakedata_path) //
"GlobalLakeDepth_GLOBathy.dat"
595 OPEN(10, file=lakefile, form=
'unformatted', access=
'direct', status=
'old', recl=data_sz*2)
596 READ(10,rec=1) lake_dpth
611 INTEGER,
INTENT(IN) :: cs_res
612 REAL,
INTENT(IN) :: cs_lakestat(:)
613 REAL,
INTENT(IN) :: cs_lakedpth(:)
615 INTEGER :: tile_sz, tile_num
616 INTEGER :: stat, ncid, x_dimid, y_dimid, varid, dimids(2)
617 INTEGER :: lake_frac_id, lake_depth_id
618 INTEGER :: land_frac_id, slmsk_id, inland_id, geolon_id, geolat_id
619 CHARACTER(len=256) :: filename,string,lakeinfo
620 CHARACTER(len=1) :: ich
621 CHARACTER(len=4) res_ch
622 REAL ::
lake_frac(cs_res*cs_res),lake_depth(cs_res*cs_res)
623 REAL :: geolon(cs_res*cs_res),geolat(cs_res*cs_res)
624 REAL :: land_frac(cs_res*cs_res),slmsk(cs_res*cs_res),inland(cs_res*cs_res)
625 real,
parameter :: epsil=1.e-6
626 real :: land_cutoff=1.e-4
630 tile_sz = cs_res*cs_res
632 WRITE(res_ch,
'(I4)') cs_res
633 DO tile_num = tile_beg, tile_end
634 WRITE(ich,
'(I1)') tile_num
637 filename =
"oro.C" // trim(adjustl(res_ch)) //
".tile" // ich //
".nc"
638 print *,
'Read, update, and write ',trim(filename)
639 stat = nf90_open(filename, nf90_write, ncid)
640 CALL
nc_opchk(stat,
"nf90_open oro_data.nc")
641 stat = nf90_inq_dimid(ncid,
"lon", x_dimid)
642 CALL
nc_opchk(stat,
"nf90_inq_dim: x")
643 stat = nf90_inq_dimid(ncid,
"lat", y_dimid)
644 CALL
nc_opchk(stat,
"nf90_inq_dim: y")
647 dimids = (/ x_dimid, y_dimid /)
648 stat = nf90_inq_varid(ncid,
"land_frac", land_frac_id)
649 CALL
nc_opchk(stat,
"nf90_inq_varid: land_frac")
650 stat = nf90_inq_varid(ncid,
"slmsk", slmsk_id)
651 CALL
nc_opchk(stat,
"nf90_inq_varid: slmsk")
653 stat = nf90_redef(ncid)
655 stat = nf90_def_var(ncid,
"lake_frac",nf90_float,dimids,lake_frac_id)
656 CALL
nc_opchk(stat,
"nf90_def_var: lake_frac")
657 #ifdef ADD_ATT_FOR_NEW_VAR
658 stat = nf90_put_att(ncid, lake_frac_id,
'coordinates',
'geolon geolat')
659 CALL
nc_opchk(stat,
"nf90_put_att: lake_frac:coordinates")
660 stat = nf90_put_att(ncid, lake_frac_id,
'long_name',
'lake fraction')
661 CALL
nc_opchk(stat,
"nf90_put_att: lake_frac:long_name")
662 stat = nf90_put_att(ncid, lake_frac_id,
'unit',
'fraction')
663 CALL
nc_opchk(stat,
"nf90_put_att: lake_frac:unit")
664 write(lakeinfo,
'(a,f4.2,a,i1)')
' lake_frac cutoff=',lake_cutoff,
'; binary_lake=',binary_lake
665 IF (lakestatus_srce ==
"GLDBV3")
THEN
666 write(string,
'(2a)')
'based on GLDBv3 (Choulga et al. 2019); missing lakes & added based on land_frac in this dataset;', &
668 ELSE IF (lakestatus_srce ==
"GLDBV2")
THEN
669 write(string,
'(2a)')
'based on GLDBv2 (Choulga et al. 2014); missing lakes & added based on land_frac in this dataset;', &
671 ELSE IF (lakestatus_srce ==
"MODISP")
THEN
672 write(string,
'(2a)')
'based on MODIS (2011-2015) product updated with two &
673 Landsat products: the JRC water product (2016-2020) and the GLC-FCS30 (2020); &
674 the source data set was created by Chengquan Huang of UMD;',trim(lakeinfo)
675 ELSE IF (lakestatus_srce ==
"VIIRS")
THEN
676 write(string,
'(2a)')
'based on multi-year VIIRS global surface type &
677 classification map (2012-2019); the source data set was created by &
678 Chengquan Huang of UMD and Michael Barlage of NOAA;',trim(lakeinfo)
680 stat = nf90_put_att(ncid, lake_frac_id,
'description',trim(string))
681 CALL
nc_opchk(stat,
"nf90_put_att: lake_frac:description")
683 stat = nf90_def_var(ncid,
"lake_depth",nf90_float,dimids,lake_depth_id)
684 CALL
nc_opchk(stat,
"nf90_def_var: lake_depth")
685 #ifdef ADD_ATT_FOR_NEW_VAR
686 stat = nf90_put_att(ncid, lake_depth_id,
'coordinates',
'geolon geolat')
687 CALL
nc_opchk(stat,
"nf90_put_att: lake_depth:coordinates")
688 stat = nf90_put_att(ncid, lake_depth_id,
'long_name',
'lake depth')
689 CALL
nc_opchk(stat,
"nf90_put_att: lake_depth:long_name")
690 stat = nf90_put_att(ncid, lake_depth_id,
'unit',
'meter')
691 CALL
nc_opchk(stat,
"nf90_put_att: lake_depth:long_name")
692 IF (lakedepth_srce ==
"GLDBV3")
THEN
693 stat = nf90_put_att(ncid, lake_depth_id,
'description', &
694 'based on GLDBv3 (Choulga et al. 2019); missing depth set to 10m &
695 (except to 211m in Caspian Sea); spurious large pos. depths are left unchanged.')
696 ELSE IF (lakedepth_srce ==
"GLDBV2")
THEN
697 stat = nf90_put_att(ncid, lake_depth_id,
'description', &
698 'based on GLDBv2 (Choulga et al. 2014); missing depth set to 10m &
699 (except to 211m in Caspian Sea); spurious large pos. depths are left unchanged.')
700 ELSE IF (lakedepth_srce ==
"GLOBATHY")
THEN
701 stat = nf90_put_att(ncid, lake_depth_id,
'description', &
702 'based on GLOBathy data resampled and projected to the MODIS domain.')
704 CALL
nc_opchk(stat,
"nf90_put_att: lake_depth:description")
707 write(string,
'(a,es8.1)')
'land_frac and lake_frac are adjusted such that &
708 their sum is 1 at points where inland=1; land_frac cutoff is',land_cutoff
709 stat = nf90_put_att(ncid, land_frac_id,
'description',trim(string))
710 CALL
nc_opchk(stat,
"nf90_put_att: land_frac:description")
712 write(string,
'(a)')
'slmsk = nint(land_frac)'
713 stat = nf90_put_att(ncid, slmsk_id,
'description',trim(string))
714 CALL
nc_opchk(stat,
"nf90_put_att: slmsk:description")
716 stat = nf90_enddef(ncid)
720 stat = nf90_inq_varid(ncid,
"geolon", geolon_id)
721 CALL
nc_opchk(stat,
"nf90_inq_varid: geolon")
722 stat = nf90_inq_varid(ncid,
"geolat", geolat_id)
723 CALL
nc_opchk(stat,
"nf90_inq_varid: geolat")
724 stat = nf90_inq_varid(ncid,
"land_frac", land_frac_id)
725 CALL
nc_opchk(stat,
"nf90_inq_varid: land_frac")
726 stat = nf90_inq_varid(ncid,
"slmsk", slmsk_id)
727 CALL
nc_opchk(stat,
"nf90_inq_varid: slmsk")
728 stat = nf90_inq_varid(ncid,
"inland", inland_id)
729 CALL
nc_opchk(stat,
"nf90_inq_varid: inland")
730 stat = nf90_get_var(ncid, geolon_id, geolon, &
731 start = (/ 1, 1 /), count = (/ cs_res, cs_res /) )
732 CALL
nc_opchk(stat,
"nf90_get_var: geolon")
733 stat = nf90_get_var(ncid, geolat_id, geolat, &
734 start = (/ 1, 1 /), count = (/ cs_res, cs_res /) )
735 CALL
nc_opchk(stat,
"nf90_get_var: geolat")
736 stat = nf90_get_var(ncid, land_frac_id, land_frac, &
737 start = (/ 1, 1 /), count = (/ cs_res, cs_res /) )
738 CALL
nc_opchk(stat,
"nf90_get_var: land_frac")
739 stat = nf90_get_var(ncid, slmsk_id, slmsk, &
740 start = (/ 1, 1 /), count = (/ cs_res, cs_res /) )
741 CALL
nc_opchk(stat,
"nf90_get_var: slmsk")
742 stat = nf90_get_var(ncid, inland_id, inland, &
743 start = (/ 1, 1 /), count = (/ cs_res, cs_res /) )
744 CALL
nc_opchk(stat,
"nf90_get_var: inland")
746 lake_frac(:) = cs_lakestat((tile_num-1)*tile_sz+1:tile_num*tile_sz)
747 lake_depth(:) = cs_lakedepth((tile_num-1)*tile_sz+1:tile_num*tile_sz)
755 IF (min(lake_cutoff,land_cutoff) < epsil)
then
756 print *,
'lake_cutoff/land_cutoff cannot be smaller than epsil, reset...'
757 lake_cutoff=max(epsil,lake_cutoff)
758 land_cutoff=max(epsil,land_cutoff)
763 if (land_frac(i)< land_cutoff) land_frac(i)=0.
764 if (land_frac(i)>1.-land_cutoff) land_frac(i)=1.
766 if (inland(i) /= 1.)
then
772 elseif (binary_lake == 1)
then
780 if (inland(i) == 1.)
then
786 elseif (
lake_frac(i) >= lake_cutoff .and. lake_depth(i)==0.)
then
789 slmsk(i) = nint(land_frac(i))
793 stat = nf90_put_var(ncid, lake_frac_id,
lake_frac, &
794 start = (/ 1, 1 /), count = (/ cs_res, cs_res /) )
795 CALL
nc_opchk(stat,
"nf90_put_var: lake_frac")
797 stat = nf90_put_var(ncid, lake_depth_id, lake_depth, &
798 start = (/ 1, 1 /), count = (/ cs_res, cs_res /) )
799 CALL
nc_opchk(stat,
"nf90_put_var: lake_depth")
802 stat = nf90_put_var(ncid, land_frac_id, land_frac, &
803 start = (/ 1, 1 /), count = (/ cs_res, cs_res /) )
804 CALL
nc_opchk(stat,
"nf90_put_var: land_frac")
806 stat = nf90_put_var(ncid, slmsk_id, slmsk, &
807 start = (/ 1, 1 /), count = (/ cs_res, cs_res /) )
808 CALL
nc_opchk(stat,
"nf90_put_var: slmsk")
810 stat = nf90_close(ncid)
828 INTEGER,
INTENT(IN) :: cs_res, tile_x_dim, tile_y_dim
829 REAL,
INTENT(IN) :: cs_lakestat(:)
830 REAL,
INTENT(IN) :: cs_lakedpth(:)
832 INTEGER :: tile_sz, tile_num
833 INTEGER :: stat, ncid, x_dimid, y_dimid, varid, dimids(2)
834 INTEGER :: lake_frac_id, lake_depth_id
835 INTEGER :: land_frac_id, slmsk_id, geolon_id, geolat_id, inland_id
836 CHARACTER(len=256) :: filename,string
837 CHARACTER(len=1) :: ich
838 CHARACTER(len=4) res_ch
840 REAL,
ALLOCATABLE ::
lake_frac(:), lake_depth(:), geolon(:), geolat(:)
841 REAL,
ALLOCATABLE :: land_frac(:), slmsk(:), inland(:)
843 real,
parameter :: epsil=1.e-6
844 real :: land_cutoff=1.e-6
846 INTEGER :: i, j, var_id
849 tile_sz = tile_x_dim*tile_y_dim
851 ALLOCATE(
lake_frac(tile_sz), lake_depth(tile_sz))
852 ALLOCATE(geolon(tile_sz), geolat(tile_sz))
853 ALLOCATE(land_frac(tile_sz), slmsk(tile_sz), inland(tile_sz))
855 WRITE(res_ch,
'(I4)') cs_res
857 WRITE(ich,
'(I1)') tile_num
859 filename =
"oro.C" // trim(adjustl(res_ch)) //
".tile" // ich //
".nc"
860 print*,
'Open and write regional orography data netCDF file ', trim(filename)
861 stat = nf90_open(filename, nf90_write, ncid)
862 CALL
nc_opchk(stat,
"nf90_open oro_data.nc")
863 stat = nf90_inq_dimid(ncid,
"lon", x_dimid)
864 CALL
nc_opchk(stat,
"nf90_inq_dim: x")
865 stat = nf90_inq_dimid(ncid,
"lat", y_dimid)
866 CALL
nc_opchk(stat,
"nf90_inq_dim: y")
867 dimids = (/ x_dimid, y_dimid /)
869 stat = nf90_inq_varid(ncid,
"land_frac", land_frac_id)
870 CALL
nc_opchk(stat,
"nf90_inq_varid: land_frac")
871 stat = nf90_inq_varid(ncid,
"slmsk", slmsk_id)
872 CALL
nc_opchk(stat,
"nf90_inq_varid: slmsk")
875 stat = nf90_redef(ncid)
878 IF (nf90_inq_varid(ncid,
"lake_frac",lake_frac_id) /= 0)
THEN
879 stat = nf90_def_var(ncid,
"lake_frac",nf90_float,dimids,lake_frac_id)
880 CALL
nc_opchk(stat,
"nf90_def_var: lake_frac")
881 #ifdef ADD_ATT_FOR_NEW_VAR
882 stat = nf90_put_att(ncid, lake_frac_id,
'coordinates',
'geolon geolat')
883 CALL
nc_opchk(stat,
"nf90_put_att: lake_frac:coordinates")
884 stat = nf90_put_att(ncid, lake_frac_id,
'long_name',
'lake fraction')
885 CALL
nc_opchk(stat,
"nf90_put_att: lake_frac:long_name")
886 stat = nf90_put_att(ncid, lake_frac_id,
'unit',
'fraction')
887 CALL
nc_opchk(stat,
"nf90_put_att: lake_frac:unit")
888 IF (lakestatus_srce ==
"GLDBV3")
THEN
889 write(string,
'(a,es8.1)')
'based on GLDBv3 (Choulga et al. 2019); missing lakes &
890 added based on land_frac in this dataset; lake_frac cutoff:',lake_cutoff
891 ELSE IF (lakestatus_srce ==
"GLDBV2")
THEN
892 write(string,
'(a,es8.1)')
'based on GLDBv2 (Choulga et al. 2014); missing lakes &
893 added based on land_frac in this dataset; lake_frac cutoff:',lake_cutoff
894 ELSE IF (lakestatus_srce ==
"MODISP")
THEN
895 write(string,
'(a,es8.1)')
'based on MODIS (2011-2015) product updated with two &
896 Landsat products: the JRC water product (2016-2020) and the GLC-FCS30 (2020); &
897 the source data set was created by Chengquan Huang of UMD; &
898 lake_frac cutoff:',lake_cutoff
899 ELSE IF (lakestatus_srce ==
"VIIRS")
THEN
900 write(string,
'(a,es8.1)')
'based on multi-year VIIRS global surface type &
901 classification map (2012-2019); the source data set was created by &
902 Chengquan Huang of UMD and Michael Barlage of NOAA; &
903 lake_frac cutoff:',lake_cutoff
905 stat = nf90_put_att(ncid, lake_frac_id,
'description',trim(string))
906 CALL
nc_opchk(stat,
"nf90_put_att: lake_frac:description")
909 IF (nf90_inq_varid(ncid,
"lake_depth",lake_depth_id) /= 0)
THEN
910 stat = nf90_def_var(ncid,
"lake_depth",nf90_float,dimids,lake_depth_id)
911 CALL
nc_opchk(stat,
"nf90_def_var: lake_depth")
912 #ifdef ADD_ATT_FOR_NEW_VAR
913 stat = nf90_put_att(ncid, lake_depth_id,
'coordinates',
'geolon geolat')
914 CALL
nc_opchk(stat,
"nf90_put_att: lake_depth:coordinates")
915 stat = nf90_put_att(ncid, lake_depth_id,
'long_name',
'lake depth')
916 CALL
nc_opchk(stat,
"nf90_put_att: lake_depth:long_name")
917 stat = nf90_put_att(ncid, lake_depth_id,
'unit',
'meter')
918 CALL
nc_opchk(stat,
"nf90_put_att: lake_depth:long_name")
919 IF (lakedepth_srce ==
"GLDBV3")
THEN
920 stat = nf90_put_att(ncid, lake_depth_id,
'description', &
921 'based on GLDBv3 (Choulga et al. 2019); missing depth set to 10m &
922 (except to 211m in Caspian Sea); spurious large pos. depths are left unchanged.')
923 ELSE IF (lakedepth_srce ==
"GLDBV2")
THEN
924 stat = nf90_put_att(ncid, lake_depth_id,
'description', &
925 'based on GLDBv2 (Choulga et al. 2014); missing depth set to 10m &
926 (except to 211m in Caspian Sea); spurious large pos. depths are left unchanged.')
927 ELSE IF (lakedepth_srce ==
"GLOBATHY")
THEN
928 stat = nf90_put_att(ncid, lake_depth_id,
'description', &
929 'based on GLOBathy data resampled and projected to the MODIS domain.')
931 CALL
nc_opchk(stat,
"nf90_put_att: lake_depth:description")
934 write(string,
'(a,es8.1)')
'land_frac and lake_frac are adjusted such that &
935 their sum is 1 at points where inland=1; land_frac cutoff is',land_cutoff
936 stat = nf90_put_att(ncid, land_frac_id,
'description',trim(string))
937 CALL
nc_opchk(stat,
"nf90_put_att: land_frac:description")
938 write(string,
'(a)')
'slmsk = nint(land_frac)'
939 stat = nf90_put_att(ncid, slmsk_id,
'description',trim(string))
940 CALL
nc_opchk(stat,
"nf90_put_att: slmsk:description")
942 stat = nf90_enddef(ncid)
946 stat = nf90_inq_varid(ncid,
"geolon", geolon_id)
947 CALL
nc_opchk(stat,
"nf90_inq_varid: geolon")
948 stat = nf90_inq_varid(ncid,
"geolat", geolat_id)
949 CALL
nc_opchk(stat,
"nf90_inq_varid: geolat")
950 stat = nf90_inq_varid(ncid,
"land_frac", land_frac_id)
951 CALL
nc_opchk(stat,
"nf90_inq_varid: land_frac")
952 stat = nf90_inq_varid(ncid,
"slmsk", slmsk_id)
953 CALL
nc_opchk(stat,
"nf90_inq_varid: slmsk")
954 stat = nf90_inq_varid(ncid,
"inland", inland_id)
955 CALL
nc_opchk(stat,
"nf90_inq_varid: inland")
957 stat = nf90_get_var(ncid, geolon_id, geolon, &
958 start = (/ 1, 1 /), count = (/ tile_x_dim, tile_y_dim /) )
959 CALL
nc_opchk(stat,
"nf90_get_var: geolon")
960 stat = nf90_get_var(ncid, geolat_id, geolat, &
961 start = (/ 1, 1 /), count = (/ tile_x_dim, tile_y_dim /) )
962 CALL
nc_opchk(stat,
"nf90_get_var: geolat")
963 stat = nf90_get_var(ncid, land_frac_id, land_frac, &
964 start = (/ 1, 1 /), count = (/ tile_x_dim, tile_y_dim /) )
965 CALL
nc_opchk(stat,
"nf90_get_var: land_frac")
966 stat = nf90_get_var(ncid, slmsk_id, slmsk, &
967 start = (/ 1, 1 /), count = (/ tile_x_dim, tile_y_dim /) )
968 CALL
nc_opchk(stat,
"nf90_get_var: slmsk")
969 stat = nf90_get_var(ncid, inland_id, inland, &
970 start = (/ 1, 1 /), count = (/ tile_x_dim, tile_y_dim /) )
971 CALL
nc_opchk(stat,
"nf90_get_var: inland")
974 lake_frac(:) = cs_lakestat((tile_num-1)*tile_sz+1:tile_num*tile_sz)
975 lake_depth(:) = cs_lakedepth((tile_num-1)*tile_sz+1:tile_num*tile_sz)
983 IF (min(lake_cutoff,land_cutoff) < epsil)
then
984 print *,
'lake_cutoff/land_cutoff cannot be smaller than epsil, reset...'
985 lake_cutoff=max(epsil,lake_cutoff)
986 land_cutoff=max(epsil,land_cutoff)
991 if (land_frac(i)< land_cutoff) land_frac(i)=0.
992 if (land_frac(i)>1.-land_cutoff) land_frac(i)=1.
994 if (inland(i) /= 1.)
then
1003 if (inland(i) == 1.)
then
1009 elseif (
lake_frac(i) >= lake_cutoff .and. lake_depth(i)==0.)
then
1012 slmsk(i) = nint(land_frac(i))
1016 stat = nf90_put_var(ncid, lake_frac_id,
lake_frac, &
1017 start = (/ 1, 1 /), count = (/ tile_x_dim, tile_y_dim /) )
1018 CALL
nc_opchk(stat,
"nf90_put_var: lake_frac")
1020 stat = nf90_put_var(ncid, lake_depth_id, lake_depth, &
1021 start = (/ 1, 1 /), count = (/ tile_x_dim, tile_y_dim /) )
1022 CALL
nc_opchk(stat,
"nf90_put_var: lake_depth")
1025 stat = nf90_put_var(ncid, land_frac_id, land_frac, &
1026 start = (/ 1, 1 /), count = (/ tile_x_dim, tile_y_dim /) )
1027 CALL
nc_opchk(stat,
"nf90_put_var: land_frac")
1029 stat = nf90_put_var(ncid, slmsk_id, slmsk, &
1030 start = (/ 1, 1 /), count = (/ tile_x_dim, tile_y_dim /) )
1031 CALL
nc_opchk(stat,
"nf90_put_var: slmsk")
1033 stat = nf90_close(ncid)
1049 REAL,
INTENT(INOUT) ::
lake_frac(cs_res*cs_res), lake_depth(cs_res*cs_res)
1050 REAL,
INTENT(IN) :: land_frac(cs_res*cs_res)
1051 REAL,
INTENT(IN) :: geolat(cs_res*cs_res), geolon(cs_res*cs_res)
1052 INTEGER,
INTENT(IN) :: tile_num
1056 tile_sz = cs_res*cs_res
1058 IF (tile_num == 2 .OR. tile_num == 3 .OR. tile_num == 7)
THEN
1059 IF (lakedepth_srce ==
"GLDBV3" .OR. lakedepth_srce ==
"GLDBV2")
THEN
1060 IF (lakestatus_srce ==
"GLDBV3" .OR. lakestatus_srce ==
"GLDBV2")
THEN
1062 IF (land_frac(i) < 0.9 .AND.
lake_frac(i) < 0.1)
THEN
1063 IF (geolat(i) > 35.0 .AND. geolat(i) <= 50.0 .AND. &
1064 geolon(i) > 45.0 .AND. geolon(i) <= 55.0)
THEN
1066 lake_depth(i) = 211.0
1068 IF (geolat(i) > 35.0 .AND. geolat(i) <= 50.0 .AND. &
1069 geolon(i) > 57.0 .AND. geolon(i) <= 63.0)
THEN
1071 lake_depth(i) = 10.0
1076 IF (lakestatus_srce ==
"MODISP" .OR. lakestatus_srce ==
"VIIRS")
THEN
1078 IF (land_frac(i) < 0.9)
THEN
1079 IF (geolat(i) > 35.0 .AND. geolat(i) <= 50.0 .AND. &
1080 geolon(i) > 45.0 .AND. geolon(i) <= 55.0)
THEN
1081 lake_depth(i) = 211.0
1083 IF (geolat(i) > 35.0 .AND. geolat(i) <= 50.0 .AND. &
1084 geolon(i) > 57.0 .AND. geolon(i) <= 63.0)
THEN
1085 lake_depth(i) = 10.0
1093 IF (tile_num == 6)
THEN
1094 IF (lakestatus_srce ==
"MODISP" .OR. lakestatus_srce ==
"VIIRS")
THEN
1096 IF (geolat(i) < -60.0)
THEN
1115 CHARACTER(len=*) opname
1118 IF (stat .NE.0)
THEN
1119 msg = trim(opname) //
' Error, status code and message:'
1120 print*,trim(msg), stat, nf90_strerror(stat)
subroutine read_lakedata(lakedata_path, lake_stat, lake_dpth, nlat, nlon)
Read a high-resolution lake depth dataset, and a corresponding lake status dataset which provides a s...
subroutine lake_cell_comp(lkst, lkdp, lake_ct, lake_avg_frac, lake_dpth_sum)
Compute cumulatively the lake fraction and lake depth for a cell.
subroutine write_lakedata_to_orodata(cs_res, cs_lakestat, cs_lakedpth)
Write lake depth and fraction to an existing model orography file.
subroutine read_cubed_sphere_reg_grid(res, grid, halo_depth, res_x, res_y)
Read the latitude and longitude for a regional grid from the 'grid' file.
subroutine read_cubed_sphere_grid(res, grid)
Read the latitude and longitude for a cubed-sphere grid from the 'grid' files.
subroutine find_limit(p1_in, p2_in, latmin, latmax)
Given two points on a cubed-sphere grid, compute the maximum and minimum latitudinal extent of the re...
program lake_frac
This program computes lake fraction and depth numbers for FV3 cubed sphere grid cells, from a high resolution lat/lon data set.
subroutine include_exclude_lakes(lake_frac, land_frac, lake_depth, geolat, geolon, tile_num)
Include Caspian Sea and Aral Sea if GLDB dataset is used, and exclude lakes in the coastal areas of A...
subroutine write_reg_lakedata_to_orodata(cs_res, tile_x_dim, tile_y_dim, cs_lakestat, cs_lakedpth)
Write lake depth and fraction to an existing model orography file.
LOGICAL function enclosure_cnvx(v, n, p, co_gc)
Test whether a given point 'p' is inside a convex spherical polygon defined with a series of 'n' vert...
subroutine cal_lake_frac_depth(lakestat, cs_lakestat, lakedpth, cs_lakedpth)
Calculate lake fraction and depth on the model grid from high-resolution data.
subroutine nc_opchk(stat, opname)
Check NetCDF return code and print error message.