22 use gengrid_kinds,
only: cl, cs, dbl_kind, real_kind, int_kind
23 use angles,
only: find_angq, find_ang
24 use vertices,
only: fill_vertices, fill_bottom, fill_top
25 use mapped_mask,
only: make_frac_land
26 use postwgts,
only: make_postwgts
27 use tripolegrid,
only: write_tripolegrid
28 use cicegrid,
only: write_cicegrid
29 use scripgrid,
only: write_scripgrid
30 use topoedits,
only: add_topoedits, apply_topoedits
31 use charstrings,
only: logmsg, res, dirsrc, dirout, atmres, fv3dir, editsfile
32 use charstrings,
only: maskfile, maskname, topofile, toponame, editsfile, staggerlocs, cdate, history
33 use debugprint,
only: checkseam, checkxlatlon, checkpoint
39 real(dbl_kind) :: dxt, dyt
41 real(kind=dbl_kind),
parameter :: pi = 3.14159265358979323846_dbl_kind
42 real(kind=dbl_kind),
parameter :: deg2rad = pi/180.0_dbl_kind
44 real(real_kind),
allocatable,
dimension(:,:) :: ww3dpth
45 integer(int_kind),
allocatable,
dimension(:,:) :: ww3mask
47 character(len=CL) :: fsrc, fdst, fwgt
48 character(len= 2) :: cstagger
50 integer :: rc,ncid,id,xtype
51 integer :: i,j,k,i2,j2
53 integer :: localpet, npet
54 logical :: fexist = .false.
56 type(esmf_regridmethod_flag) :: method
60 character(len= 6) :: i4fmt =
'(i4.4)' 61 character(len=CS) :: form1
62 character(len=CS) :: form2
63 character(len= 6) :: cnx
69 call esmf_vmgetglobal(vm, rc=rc)
70 call esmf_initialize(vm=vm, logkindflag=esmf_logkind_multi, rc=rc)
71 call esmf_vmget(vm, localpet=localpet, pecount=npet, rc=rc)
72 if (esmf_logfounderror(rctocheck=rc, msg=esmf_logerr_passthru, &
73 line=__line__, file=__file__))
call esmf_finalize(endflag=esmf_end_abort)
75 if (localpet == 0) roottask =.true.
77 print *,npet,
' More than one task specified; Aborting ' 78 call esmf_finalize(endflag=esmf_end_abort)
85 call read_inputnml(
'grid.nml')
87 print
'(a,2i6)',
' output grid requested ',ni,nj
88 print
'(a,2i6)',
' supergrid size used ', nx,ny
89 print
'(a)',
' output grid tag '//trim(res)
90 print
'(a)',
' supergrid source directory '//trim(dirsrc)
91 print
'(a)',
' output grid directory '//trim(dirout)
92 print
'(a)',
' atm resolution '//trim(atmres)
93 print
'(a,i6)',
' fv3 tile grid size ',npx
94 print
'(a)',
' atm mosaic directory '//trim(fv3dir)
95 print
'(a)',
' MOM6 topography file '//trim(topofile)
96 print
'(a)',
' MOM6 edits file '//trim(editsfile)
97 print *,
'editmask flag ',editmask
98 print *,
'debug flag ',debug
99 print *,
'do_postwgts flag ',do_postwgts
104 call esmf_logwrite(
"Starting gen_fixgrid", esmf_logmsg_info)
109 ivertcu = ivertct + 1; jvertcu = jvertct + 0
110 ivertcv = ivertct + 0; jvertcv = jvertct + 1
111 ivertbu = ivertct + 1; jvertbu = jvertct + 1
113 print
'(a8,4i6)',
'iVertCt ',(ivertct(i),i=1,4)
114 print
'(a8,4i6)',
'jVertCt ',(jvertct(i),i=1,4)
116 print
'(a8,4i6)',
'iVertCu ',(ivertcu(i),i=1,4)
117 print
'(a8,4i6)',
'jVertCu ',(jvertcu(i),i=1,4)
119 print
'(a8,4i6)',
'iVertCv ',(ivertcv(i),i=1,4)
120 print
'(a8,4i6)',
'jVertCv ',(jvertcv(i),i=1,4)
122 print
'(a8,4i6)',
'iVertBu ',(ivertbu(i),i=1,4)
123 print
'(a8,4i6)',
'jVertBu ',(jvertbu(i),i=1,4)
126 latct_vert = -9999.0 ; lonct_vert = -9999.0
127 latcu_vert = -9999.0 ; loncu_vert = -9999.0
128 latcv_vert = -9999.0 ; loncv_vert = -9999.0
129 latbu_vert = -9999.0 ; lonbu_vert = -9999.0
135 fsrc = trim(dirsrc)//
'/'//trim(maskfile)
137 rc = nf90_open(fsrc, nf90_nowrite, ncid)
138 print
'(a)',
'reading ocean mask from '//trim(fsrc)
139 if(rc .ne. 0)print
'(a)',
'nf90_open = '//trim(nf90_strerror(rc))
141 wet4 = 0.0; wet8 = 0.0
142 rc = nf90_inq_varid(ncid, trim(maskname), id)
143 rc = nf90_inquire_variable(ncid, id, xtype=xtype)
144 if(xtype .eq. 5)rc = nf90_get_var(ncid, id, wet4)
145 if(xtype .eq. 6)rc = nf90_get_var(ncid, id, wet8)
146 rc = nf90_close(ncid)
148 if(xtype.eq. 6)wet4 =
real(wet8,4)
157 fsrc = trim(dirsrc)//
'/'//trim(topofile)
159 rc = nf90_open(fsrc, nf90_nowrite, ncid)
160 print
'(a)',
'reading ocean topography from '//trim(fsrc)
161 if(rc .ne. 0)print
'(a)',
'nf90_open = '//trim(nf90_strerror(rc))
164 rc = nf90_inq_varid(ncid, trim(toponame), id)
165 rc = nf90_inquire_variable(ncid, id, xtype=xtype)
166 if(xtype .eq. 5)rc = nf90_get_var(ncid, id, dp4)
167 if(xtype .eq. 6)rc = nf90_get_var(ncid, id, dp8)
168 rc = nf90_close(ncid)
170 if(xtype.eq. 6)dp4 =
real(dp8,4)
180 if(trim(editsfile) ==
'none')
then 181 print
'(a)',
'Need a valid editsfile to make mask edits ' 185 fsrc = trim(dirsrc)//
'/'//trim(editsfile)
186 fdst = trim(dirout)//
'/'//
'ufs.'//trim(editsfile)
187 call add_topoedits(fsrc,fdst)
196 fsrc = trim(dirsrc)//
'/'//trim(editsfile)
197 if(editmask)fsrc = trim(dirout)//
'/'//
'ufs.'//trim(editsfile)
198 call apply_topoedits(fsrc)
204 fsrc = trim(dirsrc)//
'/'//
'ocean_hgrid.nc' 206 rc = nf90_open(fsrc, nf90_nowrite, ncid)
207 print
'(a)',
'reading supergrid from '//trim(fsrc)
208 if(rc .ne. 0)print
'(a)',
'nf90_open = '//trim(nf90_strerror(rc))
210 rc = nf90_inq_varid(ncid,
'x', id)
211 rc = nf90_get_var(ncid, id, x)
213 rc = nf90_inq_varid(ncid,
'y', id)
214 rc = nf90_get_var(ncid, id, y)
216 rc = nf90_inq_varid(ncid,
'dx', id)
217 rc = nf90_get_var(ncid, id, dx)
219 rc = nf90_inq_varid(ncid,
'dy', id)
220 rc = nf90_get_var(ncid, id, dy)
222 rc = nf90_close(ncid)
225 sg_maxlat = maxval(y)
241 ulon(i,j) = x(i2,j2)*deg2rad
242 ulat(i,j) = y(i2,j2)*deg2rad
244 angle(i,j) = -angq(i2,j2)
246 htn(i,j) = (dx(i2-1,j2) + dx(i2,j2))*100._dbl_kind
247 hte(i,j) = (dy(i2,j2-1) + dy(i2,j2))*100._dbl_kind
249 lonbu(i,j) = x(i2,j2)
250 latbu(i,j) = y(i2,j2)
252 lonct(i,j) = x(i2-1,j2-1)
253 loncu(i,j) = x(i2, j2-1)
254 loncv(i,j) = x(i2-1,j2 )
256 latct(i,j) = y(i2-1,j2-1)
257 latcu(i,j) = y(i2, j2-1)
258 latcv(i,j) = y(i2-1,j2 )
260 dxt = dx(i2-1,j2-1) + dx(i2,j2-1)
261 dyt = dy(i2-1,j2-1) + dy(i2-1,j2)
262 areact(i,j) = dxt*dyt
271 print *,
'ANGLET ',minval(anglet),maxval(anglet)
272 print *,
'ANGLE ',minval(angle),maxval(angle)
282 write(logmsg,
'(a,2e12.5)')
'min vals of hte at folds ', &
283 minval(hte(ni/2,:)),minval(hte(ni,:))
284 print
'(a)',trim(logmsg)
287 if(hte(ii,j) .le. 1.0)hte(ii,j) = 0.5*(hte(ii-1,j) + hte(ii+1,j))
289 if(hte(ii,j) .le. 1.0)hte(ii,j) = 0.5*(hte(ii-1,j) + hte( 1,j))
291 write(logmsg,
'(a,2e12.5)')
'min vals of hte at folds ', &
292 minval(hte(ni/2,:)),minval(hte(ni,:))
293 print
'(a)',trim(logmsg)
299 where(lonct .lt. 0.0)lonct = lonct + 360._dbl_kind
300 where(loncu .lt. 0.0)loncu = loncu + 360._dbl_kind
301 where(loncv .lt. 0.0)loncv = loncv + 360._dbl_kind
302 where(lonbu .lt. 0.0)lonbu = lonbu + 360._dbl_kind
313 if(latbu(i,j) .eq. sg_maxlat)ipole(1) = i
316 if(latbu(i,j) .eq. sg_maxlat)ipole(2) = i
318 write(logmsg,
'(a,2i6,2f12.2)')
'poles found at i = ',ipole,latbu(ipole(1),nj), &
320 print
'(a)',trim(logmsg)
322 if(debug)
call checkseam
325 i2 = ipole(2)+(ipole(1)-i)+1
326 xlonct(i) = lonct(i2,nj)
327 xlatct(i) = latct(i2,nj)
331 i2 = ipole(2)+(ipole(1)-i)
333 xloncu(i) = loncu(i2,nj)
334 xlatcu(i) = latcu(i2,nj)
337 if(debug)
call checkxlatlon
341 dlatbu(i) = latbu(i,1) + 2.0*(latcu(i,1) - latbu(i,1))
342 dlatcv(i) = latct(i,1) + 2.0*(latct(i,1) - latcv(i,1))
350 call fill_vertices(2,nj , ivertct,jvertct, latbu,lonbu, latct_vert,lonct_vert)
351 call fill_bottom(ivertct,jvertct, latbu,lonbu, latct_vert,lonct_vert,dlatbu)
353 call fill_vertices(2,nj , ivertcu,jvertcu, latcv,loncv, latcu_vert,loncu_vert)
354 call fill_bottom(ivertcu,jvertcu, latcv,loncv, latcu_vert,loncu_vert,dlatcv)
357 call fill_vertices(1,nj-1, ivertcv,jvertcv, latcu,loncu, latcv_vert,loncv_vert)
358 call fill_top(ivertcv,jvertcv, latcu,loncu, latcv_vert,loncv_vert, xlatcu, xloncu)
360 call fill_vertices(1,nj-1, ivertbu,jvertbu, latct,lonct, latbu_vert,lonbu_vert)
361 call fill_top(ivertbu,jvertbu, latct,lonct, latbu_vert,lonbu_vert, xlatct, xlonct)
363 if(debug)
call checkpoint
365 if(minval(latct_vert) .lt. -1.e3)stop
366 if(minval(lonct_vert) .lt. -1.e3)stop
367 if(minval(latcu_vert) .lt. -1.e3)stop
368 if(minval(loncu_vert) .lt. -1.e3)stop
369 if(minval(latcv_vert) .lt. -1.e3)stop
370 if(minval(loncv_vert) .lt. -1.e3)stop
371 if(minval(latbu_vert) .lt. -1.e3)stop
372 if(minval(lonbu_vert) .lt. -1.e3)stop
373 deallocate(xlonct, xlatct, xloncu, xlatcu, dlatbu, dlatcv)
380 call date_and_time(date=cdate)
381 history =
'created on '//trim(cdate)//
' from '//trim(fsrc)
384 fdst = trim(dirout)//
'/'//
'tripole.mx'//trim(res)//
'.nc' 385 call write_tripolegrid(trim(fdst))
388 fdst = trim(dirout)//
'/'//
'grid_cice_NEMS_mx'//trim(res)//
'.nc' 389 call write_cicegrid(trim(fdst))
390 deallocate(ulon, ulat, htn, hte)
394 cstagger = trim(staggerlocs(k))
395 fdst = trim(dirout)//
'/'//trim(cstagger)//
'.mx'//trim(res)//
'_SCRIP.nc' 396 logmsg =
'creating SCRIP file '//trim(fdst)
397 print
'(a)',trim(logmsg)
398 call write_scripgrid(trim(fdst),trim(cstagger))
400 deallocate(latcv_vert, loncv_vert)
401 deallocate(latcu_vert, loncu_vert)
402 deallocate(latbu_vert, lonbu_vert)
406 cstagger = trim(staggerlocs(1))
407 fdst= trim(dirout)//
'/'//trim(cstagger)//
'.mx'//trim(res)//
'_SCRIP_land.nc' 408 logmsg =
'creating SCRIP file '//trim(fdst)
409 print
'(a)',trim(logmsg)
410 call write_scripgrid(trim(fdst),trim(cstagger),imask=int(wet4))
411 deallocate(latct_vert, lonct_vert)
419 write(form1,
'(a)')
'('//trim(cnx)//
'f14.8)' 420 write(form2,
'(a)')
'('//trim(cnx)//
'i2)' 422 allocate(ww3mask(1:ni,1:nj)); ww3mask = wet4
423 allocate(ww3dpth(1:ni,1:nj)); ww3dpth = dp4
425 where(latct .ge. maximum_lat)ww3mask = 3
429 open(unit=21,file=trim(dirout)//
'/'//
'ww3.mx'//trim(res)//
'_x.inp',form=
'formatted')
430 open(unit=22,file=trim(dirout)//
'/'//
'ww3.mx'//trim(res)//
'_y.inp',form=
'formatted')
431 open(unit=23,file=trim(dirout)//
'/'//
'ww3.mx'//trim(res)//
'_bottom.inp',form=
'formatted')
432 open(unit=24,file=trim(dirout)//
'/'//
'ww3.mx'//trim(res)//
'_mapsta.inp',form=
'formatted')
434 open(unit=25,file=trim(dirout)//
'/'//
'ww3.mx'//trim(res)//
'_obstr.inp',form=
'formatted')
437 write( 21,trim(form1))lonct(:,j)
438 write( 22,trim(form1))latct(:,j)
441 write( 23,trim(form1))ww3dpth(:,j)
442 write( 24,trim(form2))ww3mask(:,j)
444 write( 25,trim(form2))ww3mask(:,j)*0
445 write( 25,trim(form2))ww3mask(:,j)*0
448 close(21);
close(22);
close(23);
close(24);
close(25)
449 deallocate(ww3mask);
deallocate(ww3dpth)
450 deallocate(wet4, wet8)
458 method=esmf_regridmethod_conserve
459 fsrc = trim(dirout)//
'/'//
'Ct.mx'//trim(res)//
'_SCRIP_land.nc' 460 fdst = trim(fv3dir)//
'/'//trim(atmres)//
'/'//trim(atmres)//
'_mosaic.nc' 461 fwgt = trim(dirout)//
'/'//
'Ct.mx'//trim(res)//
'.to.'//trim(atmres)//
'.nc' 462 logmsg =
'creating weight file '//trim(fwgt)
463 print
'(a)',trim(logmsg)
465 call esmf_regridweightgen(srcfile=trim(fsrc),dstfile=trim(fdst), &
466 weightfile=trim(fwgt), regridmethod=method, &
467 unmappedaction=esmf_unmappedaction_ignore, ignoredegenerate=.true., &
468 netcdf4fileflag=.true., tilefilepath=trim(fv3dir)//
'/'//trim(atmres)//
'/', rc=rc)
469 if (esmf_logfounderror(rctocheck=rc, msg=esmf_logerr_passthru, &
470 line=__line__, file=__file__))
call esmf_finalize(endflag=esmf_end_abort)
472 logmsg =
'creating mapped ocean mask for '//trim(atmres)
473 print
'(a)',trim(logmsg)
474 call make_frac_land(trim(fsrc), trim(fwgt))
482 if(trim(res) .ne.
'025')
then 483 fsrc = trim(dirout)//
'/'//
'Ct.mx025_SCRIP.nc' 484 inquire(file=trim(fsrc), exist=fexist)
486 method=esmf_regridmethod_nearest_stod
487 fdst = trim(dirout)//
'/'//
'Ct.mx'//trim(res)//
'_SCRIP.nc' 488 fwgt = trim(dirout)//
'/'//
'tripole.mx025.Ct.to.mx'//trim(res)//
'.Ct.neareststod.nc' 489 logmsg =
'creating weight file '//trim(fwgt)
490 print
'(a)',trim(logmsg)
492 call esmf_regridweightgen(srcfile=trim(fsrc),dstfile=trim(fdst), &
493 weightfile=trim(fwgt), regridmethod=method, &
494 ignoredegenerate=.true., unmappedaction=esmf_unmappedaction_ignore, rc=rc)
495 if (esmf_logfounderror(rctocheck=rc, msg=esmf_logerr_passthru, &
496 line=__line__, file=__file__))
call esmf_finalize(endflag=esmf_end_abort)
498 logmsg =
'ERROR: '//trim(fsrc)//
' is required to generate tripole:triple weights' 499 print
'(a)',trim(logmsg)
508 if(do_postwgts)
call make_postwgts
514 deallocate(x,y, angq, dx, dy, xsgp1, ysgp1)
515 deallocate(areact, anglet, angle)
516 deallocate(latct, lonct)
517 deallocate(latcv, loncv)
518 deallocate(latcu, loncu)
519 deallocate(latbu, lonbu)
program gen_fixgrid
Generate fixed grid files required for coupled model using the MOM6 super grid file and ocean mask fi...