cpld_gridgen  1.11.0
All Data Structures Files Functions Variables Pages
gen_fixgrid.F90
Go to the documentation of this file.
1 
5 
16 program gen_fixgrid
17 
18  use esmf
19 
20  use grdvars
21  use inputnml
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
34  use netcdf
35 
36  implicit none
37 
38  ! local variables
39  real(dbl_kind) :: dxt, dyt
40 
41  real(kind=dbl_kind), parameter :: pi = 3.14159265358979323846_dbl_kind
42  real(kind=dbl_kind), parameter :: deg2rad = pi/180.0_dbl_kind
43 
44  real(real_kind), allocatable, dimension(:,:) :: ww3dpth
45  integer(int_kind), allocatable, dimension(:,:) :: ww3mask
46 
47  character(len=CL) :: fsrc, fdst, fwgt
48  character(len= 2) :: cstagger
49 
50  integer :: rc,ncid,id,xtype
51  integer :: i,j,k,i2,j2
52  integer :: ii,jj
53  integer :: localpet, npet
54  logical :: fexist = .false.
55 
56  type(esmf_regridmethod_flag) :: method
57  type(esmf_vm) :: vm
58 
59  !WW3 file format for mod_def generation
60  character(len= 6) :: i4fmt = '(i4.4)'
61  character(len=CS) :: form1
62  character(len=CS) :: form2
63  character(len= 6) :: cnx
64 
65  !-------------------------------------------------------------------------
66  ! Initialize esmf environment.
67  !-------------------------------------------------------------------------
68 
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)
74  roottask = .false.
75  if (localpet == 0) roottask =.true.
76  if (npet /= 1) then
77  print *,npet,' More than one task specified; Aborting '
78  call esmf_finalize(endflag=esmf_end_abort)
79  end if
80 
81  !---------------------------------------------------------------------
82  !
83  !---------------------------------------------------------------------
84 
85  call read_inputnml('grid.nml')
86 
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
100  print *
101 
102  call allocate_all
103 
104  call esmf_logwrite("Starting gen_fixgrid", esmf_logmsg_info)
105  !---------------------------------------------------------------------
106  ! set up the arrays to retrieve the vertices
107  !---------------------------------------------------------------------
108 
109  ivertcu = ivertct + 1; jvertcu = jvertct + 0
110  ivertcv = ivertct + 0; jvertcv = jvertct + 1
111  ivertbu = ivertct + 1; jvertbu = jvertct + 1
112 
113  print '(a8,4i6)','iVertCt ',(ivertct(i),i=1,4)
114  print '(a8,4i6)','jVertCt ',(jvertct(i),i=1,4)
115  print *
116  print '(a8,4i6)','iVertCu ',(ivertcu(i),i=1,4)
117  print '(a8,4i6)','jVertCu ',(jvertcu(i),i=1,4)
118  print *
119  print '(a8,4i6)','iVertCv ',(ivertcv(i),i=1,4)
120  print '(a8,4i6)','jVertCv ',(jvertcv(i),i=1,4)
121  print *
122  print '(a8,4i6)','iVertBu ',(ivertbu(i),i=1,4)
123  print '(a8,4i6)','jVertBu ',(jvertbu(i),i=1,4)
124  print *
125 
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
130 
131  !---------------------------------------------------------------------
132  ! read the MOM6 land mask
133  !---------------------------------------------------------------------
134 
135  fsrc = trim(dirsrc)//'/'//trim(maskfile)
136 
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))
140 
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)
147 
148  if(xtype.eq. 6)wet4 = real(wet8,4)
149 
150  !print *,minval(wet8),maxval(wet8)
151  !print *,minval(wet4),maxval(wet4)
152 
153  !---------------------------------------------------------------------
154  ! read the MOM6 depth file
155  !---------------------------------------------------------------------
156 
157  fsrc = trim(dirsrc)//'/'//trim(topofile)
158 
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))
162 
163  dp4 = 0.0; dp8 = 0.0
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)
169 
170  if(xtype.eq. 6)dp4 = real(dp8,4)
171 
172  !print *,minval(dp8),maxval(dp8)
173  !print *,minval(dp4),maxval(dp4)
174 
175  if(editmask)then
176  !---------------------------------------------------------------------
177  ! apply topoedits run time mask changes if required for this config
178  !---------------------------------------------------------------------
179 
180  if(trim(editsfile) == 'none')then
181  print '(a)', 'Need a valid editsfile to make mask edits '
182  stop
183  end if
184 
185  fsrc = trim(dirsrc)//'/'//trim(editsfile)
186  fdst = trim(dirout)//'/'//'ufs.'//trim(editsfile)
187  call add_topoedits(fsrc,fdst)
188  endif
189 
190  !---------------------------------------------------------------------
191  ! MOM6 reads the depth file, applies the topo edits and then adjusts
192  ! depth using masking_depth and min/max depth. This call mimics
193  ! MOM6 routines apply_topography_edits_from_file and limit_topography
194  !---------------------------------------------------------------------
195 
196  fsrc = trim(dirsrc)//'/'//trim(editsfile)
197  if(editmask)fsrc = trim(dirout)//'/'//'ufs.'//trim(editsfile)
198  call apply_topoedits(fsrc)
199 
200  !---------------------------------------------------------------------
201  ! read MOM6 supergrid file
202  !---------------------------------------------------------------------
203 
204  fsrc = trim(dirsrc)//'/'//'ocean_hgrid.nc'
205 
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))
209 
210  rc = nf90_inq_varid(ncid, 'x', id) !lon
211  rc = nf90_get_var(ncid, id, x)
212 
213  rc = nf90_inq_varid(ncid, 'y', id) !lat
214  rc = nf90_get_var(ncid, id, y)
215 
216  rc = nf90_inq_varid(ncid, 'dx', id)
217  rc = nf90_get_var(ncid, id, dx)
218 
219  rc = nf90_inq_varid(ncid, 'dy', id)
220  rc = nf90_get_var(ncid, id, dy)
221 
222  rc = nf90_close(ncid)
223  !print *,'super grid size ',size(y,1),size(y,2)
224  !print *,'max lat in super grid ',maxval(y)
225  sg_maxlat = maxval(y)
226 
227  !---------------------------------------------------------------------
228  ! find the angle on corners---this requires the supergrid
229  !---------------------------------------------------------------------
230 
231  call find_angq
232 
233  !---------------------------------------------------------------------
234  ! fill grid variables
235  !---------------------------------------------------------------------
236 
237  do j = 1,nj
238  do i = 1,ni
239  i2 = 2*i ; j2 = 2*j
240  !deg->rad
241  ulon(i,j) = x(i2,j2)*deg2rad
242  ulat(i,j) = y(i2,j2)*deg2rad
243  !in rad already
244  angle(i,j) = -angq(i2,j2)
245  !m->cm
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
248  !deg
249  lonbu(i,j) = x(i2,j2)
250  latbu(i,j) = y(i2,j2)
251  !deg
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 )
255  !deg
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 )
259  !m2
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
263  enddo
264  enddo
265 
266  !---------------------------------------------------------------------
267  ! find the angle on centers---this does not requires the supergrid
268  !---------------------------------------------------------------------
269 
270  call find_ang
271  print *,'ANGLET ',minval(anglet),maxval(anglet)
272  print *,'ANGLE ',minval(angle),maxval(angle)
273 
274  !---------------------------------------------------------------------
275  ! For the 1/4deg grid, hte at j=720 and j = 1440 is identically=0.0 for
276  ! j > 840 (64.0N). These are land points, but since CICE uses hte to
277  ! generate remaining variables, setting them to zero will cause problems
278  ! For 1deg grid, hte at ni/2 and ni are very small O~10-12, so test for
279  ! hte < 1.0
280  !---------------------------------------------------------------------
281 
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)
285  do j = 1,nj
286  ii = ni/2
287  if(hte(ii,j) .le. 1.0)hte(ii,j) = 0.5*(hte(ii-1,j) + hte(ii+1,j))
288  ii = ni
289  if(hte(ii,j) .le. 1.0)hte(ii,j) = 0.5*(hte(ii-1,j) + hte( 1,j))
290  enddo
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)
294 
295  !---------------------------------------------------------------------
296  !
297  !---------------------------------------------------------------------
298 
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
303 
304  !---------------------------------------------------------------------
305  ! some basic error checking
306  ! find the i-th index of the poles at j= nj
307  ! the corner points must lie on the pole
308  !---------------------------------------------------------------------
309 
310  ipole = -1
311  j = nj
312  do i = 1,ni/2
313  if(latbu(i,j) .eq. sg_maxlat)ipole(1) = i
314  enddo
315  do i = ni/2+1,ni
316  if(latbu(i,j) .eq. sg_maxlat)ipole(2) = i
317  enddo
318  write(logmsg,'(a,2i6,2f12.2)')'poles found at i = ',ipole,latbu(ipole(1),nj), &
319  latbu(ipole(2),nj)
320  print '(a)',trim(logmsg)
321 
322  if(debug)call checkseam
323 
324  do i = 1,ni
325  i2 = ipole(2)+(ipole(1)-i)+1
326  xlonct(i) = lonct(i2,nj)
327  xlatct(i) = latct(i2,nj)
328  enddo
329 
330  do i = 1,ni
331  i2 = ipole(2)+(ipole(1)-i)
332  if(i2 .lt. 1)i2 = ni
333  xloncu(i) = loncu(i2,nj)
334  xlatcu(i) = latcu(i2,nj)
335  enddo
336 
337  if(debug)call checkxlatlon
338 
339  !approx lat at grid bottom
340  do i = 1,ni
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))
343  enddo
344 
345  !---------------------------------------------------------------------
346  ! fill grid vertices variables
347  !---------------------------------------------------------------------
348 
349  !Ct and Cu grids align in j
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)
352 
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)
355 
356  !Cv and Bu grids align in j
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)
359 
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)
362 
363  if(debug)call checkpoint
364 
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)
374 
375  !---------------------------------------------------------------------
376  ! write out grid file files
377  !---------------------------------------------------------------------
378 
379  ! create a history attribute
380  call date_and_time(date=cdate)
381  history = 'created on '//trim(cdate)//' from '//trim(fsrc)
382 
383  ! write fix grid
384  fdst = trim(dirout)//'/'//'tripole.mx'//trim(res)//'.nc'
385  call write_tripolegrid(trim(fdst))
386 
387  ! write cice grid
388  fdst = trim(dirout)//'/'//'grid_cice_NEMS_mx'//trim(res)//'.nc'
389  call write_cicegrid(trim(fdst))
390  deallocate(ulon, ulat, htn, hte)
391  ! write scrip grids; only the Ct is required, the remaining
392  ! staggers are used only in the postweights generation
393  do k = 1,nv
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))
399  end do
400  deallocate(latcv_vert, loncv_vert)
401  deallocate(latcu_vert, loncu_vert)
402  deallocate(latbu_vert, lonbu_vert)
403 
404  ! write SCRIP file with land mask, used for mapped ocean mask
405  ! and mesh creation
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)
412 
413  !---------------------------------------------------------------------
414  ! write lat,lon,depth and mask arrays required by ww3 in creating
415  ! mod_def file
416  !---------------------------------------------------------------------
417 
418  write(cnx,i4fmt)nx
419  write(form1,'(a)')'('//trim(cnx)//'f14.8)'
420  write(form2,'(a)')'('//trim(cnx)//'i2)'
421 
422  allocate(ww3mask(1:ni,1:nj)); ww3mask = wet4
423  allocate(ww3dpth(1:ni,1:nj)); ww3dpth = dp4
424 
425  where(latct .ge. maximum_lat)ww3mask = 3
426  !close last row
427  ww3mask(:,nj) = 3
428 
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')
433  ! cice0 .ne. cicen requires obstruction map, should be initialized as zeros (w3grid,ln3032)
434  open(unit=25,file=trim(dirout)//'/'//'ww3.mx'//trim(res)//'_obstr.inp',form='formatted')
435 
436  do j = 1,nj
437  write( 21,trim(form1))lonct(:,j)
438  write( 22,trim(form1))latct(:,j)
439  end do
440  do j = 1,nj
441  write( 23,trim(form1))ww3dpth(:,j)
442  write( 24,trim(form2))ww3mask(:,j)
443  !'obsx' and 'obsy' arrays ???
444  write( 25,trim(form2))ww3mask(:,j)*0
445  write( 25,trim(form2))ww3mask(:,j)*0
446  end do
447 
448  close(21); close(22); close(23); close(24); close(25)
449  deallocate(ww3mask); deallocate(ww3dpth)
450  deallocate(wet4, wet8)
451 
452  !---------------------------------------------------------------------
453  ! use ESMF regridding to produce mapped ocean mask; first generate
454  ! conservative regrid weights from ocean to tiles; then generate the
455  ! tiled files containing the mapped ocean mask
456  !---------------------------------------------------------------------
457 
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)
464 
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)
471 
472  logmsg = 'creating mapped ocean mask for '//trim(atmres)
473  print '(a)',trim(logmsg)
474  call make_frac_land(trim(fsrc), trim(fwgt))
475 
476  !---------------------------------------------------------------------
477  ! use ESMF to find the tripole:tripole weights for creation
478  ! of CICE ICs; the source grid is always mx025; don't create this
479  ! file if destination is also mx025
480  !---------------------------------------------------------------------
481 
482  if(trim(res) .ne. '025') then
483  fsrc = trim(dirout)//'/'//'Ct.mx025_SCRIP.nc'
484  inquire(file=trim(fsrc), exist=fexist)
485  if (fexist ) then
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)
491 
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)
497  else
498  logmsg = 'ERROR: '//trim(fsrc)//' is required to generate tripole:triple weights'
499  print '(a)',trim(logmsg)
500  stop
501  end if
502  end if
503 
504  !---------------------------------------------------------------------
505  !
506  !---------------------------------------------------------------------
507 
508  if(do_postwgts)call make_postwgts
509 
510  !---------------------------------------------------------------------
511  ! clean up
512  !---------------------------------------------------------------------
513 
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)
520 
521 end program gen_fixgrid
program gen_fixgrid
Generate fixed grid files required for coupled model using the MOM6 super grid file and ocean mask fi...
Definition: gen_fixgrid.F90:16