cpld_gridgen  1.13.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_ang, find_angq, find_angchk
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, atmres, dirsrc, dirout, 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,n,i2,j2
52  integer :: ii
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 mosaic directory '//trim(fv3dir)
93  print '(a)',' MOM6 topography file '//trim(topofile)
94  print '(a)',' MOM6 edits file '//trim(editsfile)
95  print *,'editmask flag ',editmask
96  print *,'debug flag ',debug
97  print *,'do_postwgts flag ',do_postwgts
98  print *
99 
100  call allocate_all
101 
102  call esmf_logwrite("Starting gen_fixgrid", esmf_logmsg_info)
103  !---------------------------------------------------------------------
104  ! set up the arrays to retrieve the vertices
105  !---------------------------------------------------------------------
106 
107  ivertcu = ivertct + 1; jvertcu = jvertct + 0
108  ivertcv = ivertct + 0; jvertcv = jvertct + 1
109  ivertbu = ivertct + 1; jvertbu = jvertct + 1
110 
111  print '(a8,4i6)','iVertCt ',(ivertct(i),i=1,4)
112  print '(a8,4i6)','jVertCt ',(jvertct(i),i=1,4)
113  print *
114  print '(a8,4i6)','iVertCu ',(ivertcu(i),i=1,4)
115  print '(a8,4i6)','jVertCu ',(jvertcu(i),i=1,4)
116  print *
117  print '(a8,4i6)','iVertCv ',(ivertcv(i),i=1,4)
118  print '(a8,4i6)','jVertCv ',(jvertcv(i),i=1,4)
119  print *
120  print '(a8,4i6)','iVertBu ',(ivertbu(i),i=1,4)
121  print '(a8,4i6)','jVertBu ',(jvertbu(i),i=1,4)
122  print *
123 
124  latct_vert = -9999.0 ; lonct_vert = -9999.0
125  latcu_vert = -9999.0 ; loncu_vert = -9999.0
126  latcv_vert = -9999.0 ; loncv_vert = -9999.0
127  latbu_vert = -9999.0 ; lonbu_vert = -9999.0
128 
129  !---------------------------------------------------------------------
130  ! read the MOM6 land mask
131  !---------------------------------------------------------------------
132 
133  fsrc = trim(dirsrc)//'/'//trim(maskfile)
134 
135  rc = nf90_open(fsrc, nf90_nowrite, ncid)
136  print '(a)', 'reading ocean mask from '//trim(fsrc)
137  if(rc .ne. 0)print '(a)', 'nf90_open = '//trim(nf90_strerror(rc))
138 
139  wet4 = 0.0; wet8 = 0.0
140  rc = nf90_inq_varid(ncid, trim(maskname), id)
141  rc = nf90_inquire_variable(ncid, id, xtype=xtype)
142  if(xtype .eq. 5)rc = nf90_get_var(ncid, id, wet4)
143  if(xtype .eq. 6)rc = nf90_get_var(ncid, id, wet8)
144  rc = nf90_close(ncid)
145 
146  if(xtype.eq. 6)wet4 = real(wet8,4)
147 
148  !---------------------------------------------------------------------
149  ! read the MOM6 depth file
150  !---------------------------------------------------------------------
151 
152  fsrc = trim(dirsrc)//'/'//trim(topofile)
153 
154  rc = nf90_open(fsrc, nf90_nowrite, ncid)
155  print '(a)', 'reading ocean topography from '//trim(fsrc)
156  if(rc .ne. 0)print '(a)', 'nf90_open = '//trim(nf90_strerror(rc))
157 
158  dp4 = 0.0; dp8 = 0.0
159  rc = nf90_inq_varid(ncid, trim(toponame), id)
160  rc = nf90_inquire_variable(ncid, id, xtype=xtype)
161  if(xtype .eq. 5)rc = nf90_get_var(ncid, id, dp4)
162  if(xtype .eq. 6)rc = nf90_get_var(ncid, id, dp8)
163  rc = nf90_close(ncid)
164 
165  if(xtype.eq. 6)dp4 = real(dp8,4)
166 
167  if(editmask)then
168  !---------------------------------------------------------------------
169  ! apply topoedits run time mask changes if required for this config
170  ! this will create a modified topoedits file which accounts for any
171  ! land mask changes created at run time by MOM6
172  !---------------------------------------------------------------------
173 
174  if(trim(editsfile) == 'none')then
175  print '(a)', 'Need a valid editsfile to make mask edits '
176  call abort()
177  end if
178  inquire(file=trim(dirsrc)//'/'//trim(editsfile),exist=fexist)
179  if (.not. fexist) then
180  print '(a)', 'Required topoedits file '//trim(editsfile) &
181  //'for land mask changes is missing '
182  call abort()
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  ! If the the topoedits file has been modified to account for MOM6 run
195  ! time land mask changes (above), then the depth will be created using
196  ! this modified topoedits file
197  !---------------------------------------------------------------------
198 
199  fsrc = trim(dirsrc)//'/'//trim(editsfile)
200  if(editmask)fsrc = trim(dirout)//'/'//'ufs.'//trim(editsfile)
201 
202  if (trim(editsfile) /= 'none') then
203  inquire(file=trim(fsrc),exist=fexist)
204  if (.not. fexist) then
205  print '(a)', 'Required topoedits file '//trim(fsrc)//' is missing '
206  call abort()
207  end if
208  end if
209  call apply_topoedits(fsrc)
210 
211  !---------------------------------------------------------------------
212  ! read MOM6 supergrid file
213  !---------------------------------------------------------------------
214 
215  fsrc = trim(dirsrc)//'/'//'ocean_hgrid.nc'
216 
217  rc = nf90_open(fsrc, nf90_nowrite, ncid)
218  print '(a)', 'reading supergrid from '//trim(fsrc)
219  if(rc .ne. 0)print '(a)', 'nf90_open = '//trim(nf90_strerror(rc))
220 
221  rc = nf90_inq_varid(ncid, 'x', id) !lon
222  rc = nf90_get_var(ncid, id, x)
223 
224  rc = nf90_inq_varid(ncid, 'y', id) !lat
225  rc = nf90_get_var(ncid, id, y)
226 
227  rc = nf90_inq_varid(ncid, 'dx', id)
228  rc = nf90_get_var(ncid, id, dx)
229 
230  rc = nf90_inq_varid(ncid, 'dy', id)
231  rc = nf90_get_var(ncid, id, dy)
232 
233  rc = nf90_close(ncid)
234  sg_maxlat = maxval(y)
235  write(logmsg,'(a,f12.2)')'max lat in super grid ',maxval(y)
236  print '(a)',trim(logmsg)
237 
238  !---------------------------------------------------------------------
239  ! fill grid variables
240  !---------------------------------------------------------------------
241 
242  do j = 1,nj
243  do i = 1,ni
244  i2 = 2*i ; j2 = 2*j
245  !deg->rad
246  ulon(i,j) = x(i2,j2)*deg2rad
247  ulat(i,j) = y(i2,j2)*deg2rad
248  !m->cm
249  htn(i,j) = (dx(i2-1,j2) + dx(i2,j2))*100._dbl_kind
250  hte(i,j) = (dy(i2,j2-1) + dy(i2,j2))*100._dbl_kind
251  !deg
252  lonbu(i,j) = x(i2,j2)
253  latbu(i,j) = y(i2,j2)
254  !deg
255  lonct(i,j) = x(i2-1,j2-1)
256  loncu(i,j) = x(i2, j2-1)
257  loncv(i,j) = x(i2-1,j2 )
258  !deg
259  latct(i,j) = y(i2-1,j2-1)
260  latcu(i,j) = y(i2, j2-1)
261  latcv(i,j) = y(i2-1,j2 )
262  !m2
263  dxt = dx(i2-1,j2-1) + dx(i2,j2-1)
264  dyt = dy(i2-1,j2-1) + dy(i2-1,j2)
265  areact(i,j) = dxt*dyt
266  enddo
267  enddo
268 
269  !---------------------------------------------------------------------
270  ! locate the ith index of the two poles on j=nj
271  ! the corner points must lie on the pole
272  !---------------------------------------------------------------------
273 
274  ipole = -1
275  j = nj
276  do i = 1,ni/2
277  if(latbu(i,j) .eq. sg_maxlat)ipole(1) = i
278  enddo
279  do i = ni/2+1,ni
280  if(latbu(i,j) .eq. sg_maxlat)ipole(2) = i
281  enddo
282  write(logmsg,'(a,2i6,2f12.2)')'poles found at i = ',ipole, latbu(ipole(1),nj), &
283  latbu(ipole(2),nj)
284  print '(a)',trim(logmsg)
285 
286  !---------------------------------------------------------------------
287  ! find the angle on centers using the same procedure as MOM6
288  !---------------------------------------------------------------------
289 
290  call find_ang((/1,ni/),(/1,nj/),lonbu,latbu,lonct,anglet)
291  write(logmsg,'(a,2f12.2)')'ANGLET min,max: ',minval(anglet),maxval(anglet)
292  print '(a)',trim(logmsg)
293  write(logmsg,'(a,2f12.2)')'ANGLET edges i=1,i=ni: ',anglet(1,nj),anglet(ni,nj)
294  print '(a)',trim(logmsg)
295 
296  xangct(:) = 0.0
297  do i = 1,ni
298  i2 = ipole(2)+(ipole(1)-i)+1
299  xangct(i) = -anglet(i2,nj) ! angle changes sign across seam
300  end do
301 
302  !---------------------------------------------------------------------
303  ! find the angle on corners using the same procedure as CICE6
304  !---------------------------------------------------------------------
305 
306  call find_angq((/1,ni/),(/1,nj/),xangct,anglet,angle)
307  angle(ni,:) = -angle(1,:)
308  ! reverse angle for CICE
309  angle = -angle
310  write(logmsg,'(a,2f12.2)')'ANGLE min,max: ',minval(angle),maxval(angle)
311  print '(a)',trim(logmsg)
312  write(logmsg,'(a,2f12.2)')'ANGLE edges i=1,i=ni: ',angle(1,nj),angle(ni,nj)
313  print '(a)',trim(logmsg)
314 
315  !---------------------------------------------------------------------
316  ! check the Bu angle
317  !---------------------------------------------------------------------
318 
319  call find_angchk((/1,ni/),(/1,nj/),angle,angchk)
320  angchk(1,:) = -angchk(ni,:)
321  ! reverse angle for MOM6
322  angchk = -angchk
323  write(logmsg,'(a,2f12.2)')'ANGCHK min,max: ',minval(angchk),maxval(angchk)
324  print '(a)',trim(logmsg)
325  write(logmsg,'(a,2f12.2)')'ANGCHK edges i=1,i=ni: ',angchk(1,nj),angchk(ni,nj)
326  print '(a)',trim(logmsg)
327 
328  !---------------------------------------------------------------------
329  ! For the 1/4deg grid, hte at j=720 and j = 1440 is identically=0.0 for
330  ! j > 840 (64.0N). These are land points, but since CICE uses hte to
331  ! generate remaining variables, setting them to zero will cause problems
332  ! For 1deg grid, hte at ni/2 and ni are very small O~10-12, so test for
333  ! hte < 1.0
334  !---------------------------------------------------------------------
335 
336  write(logmsg,'(a,2e12.5)')'min vals of hte at folds ', minval(hte(ni/2,:)),minval(hte(ni,:))
337  print '(a)',trim(logmsg)
338  do j = 1,nj
339  ii = ni/2
340  if(hte(ii,j) .le. 1.0)hte(ii,j) = 0.5*(hte(ii-1,j) + hte(ii+1,j))
341  ii = ni
342  if(hte(ii,j) .le. 1.0)hte(ii,j) = 0.5*(hte(ii-1,j) + hte( 1,j))
343  enddo
344  write(logmsg,'(a,2e12.5)')'min vals of hte at folds ', minval(hte(ni/2,:)),minval(hte(ni,:))
345  print '(a)',trim(logmsg)
346 
347  !---------------------------------------------------------------------
348  ! find required extended values for setting all vertices
349  !---------------------------------------------------------------------
350 
351  if(debug)call checkseam
352 
353  do i = 1,ni
354  i2 = ipole(2)+(ipole(1)-i)+1
355  xlonct(i) = lonct(i2,nj)
356  xlatct(i) = latct(i2,nj)
357  enddo
358 
359  do i = 1,ni
360  i2 = ipole(2)+(ipole(1)-i)
361  if(i2 .lt. 1)i2 = ni
362  xloncu(i) = loncu(i2,nj)
363  xlatcu(i) = latcu(i2,nj)
364  enddo
365 
366  if(debug)call checkxlatlon
367 
368  !approx lat at grid bottom
369  do i = 1,ni
370  dlatbu(i) = latbu(i,1) + 2.0*(latcu(i,1) - latbu(i,1))
371  dlatcv(i) = latct(i,1) + 2.0*(latct(i,1) - latcv(i,1))
372  enddo
373 
374  !---------------------------------------------------------------------
375  ! fill grid vertices variables
376  !---------------------------------------------------------------------
377 
378  !Ct and Cu grids align in j
379  call fill_vertices(2,nj , ivertct,jvertct, latbu,lonbu, latct_vert,lonct_vert)
380  call fill_bottom(ivertct,jvertct, latbu,lonbu, latct_vert,lonct_vert,dlatbu)
381 
382  call fill_vertices(2,nj , ivertcu,jvertcu, latcv,loncv, latcu_vert,loncu_vert)
383  call fill_bottom(ivertcu,jvertcu, latcv,loncv, latcu_vert,loncu_vert,dlatcv)
384 
385  !Cv and Bu grids align in j
386  call fill_vertices(1,nj-1, ivertcv,jvertcv, latcu,loncu, latcv_vert,loncv_vert)
387  call fill_top(ivertcv,jvertcv, latcu,loncu, latcv_vert,loncv_vert, xlatcu, xloncu)
388 
389  call fill_vertices(1,nj-1, ivertbu,jvertbu, latct,lonct, latbu_vert,lonbu_vert)
390  call fill_top(ivertbu,jvertbu, latct,lonct, latbu_vert,lonbu_vert, xlatct, xlonct)
391 
392  if(debug)call checkpoint
393 
394  if(minval(latct_vert) .lt. -1.e3)stop
395  if(minval(lonct_vert) .lt. -1.e3)stop
396  if(minval(latcu_vert) .lt. -1.e3)stop
397  if(minval(loncu_vert) .lt. -1.e3)stop
398  if(minval(latcv_vert) .lt. -1.e3)stop
399  if(minval(loncv_vert) .lt. -1.e3)stop
400  if(minval(latbu_vert) .lt. -1.e3)stop
401  if(minval(lonbu_vert) .lt. -1.e3)stop
402  deallocate(xlonct, xlatct, xloncu, xlatcu, dlatbu, dlatcv)
403 
404  !---------------------------------------------------------------------
405  ! write out grid file files
406  !---------------------------------------------------------------------
407 
408  ! create a history attribute
409  call date_and_time(date=cdate)
410  history = 'created on '//trim(cdate)//' from '//trim(fsrc)
411 
412  ! write fix grid
413  fdst = trim(dirout)//'/'//'tripole.mx'//trim(res)//'.nc'
414  call write_tripolegrid(trim(fdst))
415 
416  ! write cice grid
417  fdst = trim(dirout)//'/'//'grid_cice_NEMS_mx'//trim(res)//'.nc'
418  call write_cicegrid(trim(fdst))
419  deallocate(ulon, ulat, htn, hte)
420 
421  ! write SCRIP files for generation of positional weights
422  do k = 1,nv
423  cstagger = trim(staggerlocs(k))
424  fdst = trim(dirout)//'/'//trim(cstagger)//'.mx'//trim(res)//'_SCRIP.nc'
425  logmsg = 'creating SCRIP file '//trim(fdst)
426  print '(a)',trim(logmsg)
427  call write_scripgrid(trim(fdst),trim(cstagger))
428  end do
429  deallocate(latcv_vert, loncv_vert)
430  deallocate(latcu_vert, loncu_vert)
431  deallocate(latbu_vert, lonbu_vert)
432 
433  ! write SCRIP file with land mask, used for mapped ocean mask
434  ! and mesh creation
435  cstagger = trim(staggerlocs(1))
436  fdst= trim(dirout)//'/'//trim(cstagger)//'.mx'//trim(res)//'_SCRIP_land.nc'
437  logmsg = 'creating SCRIP file '//trim(fdst)
438  print '(a)',trim(logmsg)
439  call write_scripgrid(trim(fdst),trim(cstagger),imask=int(wet4))
440  deallocate(latct_vert, lonct_vert)
441 
442  !---------------------------------------------------------------------
443  ! write lat,lon,depth and mask arrays required by ww3 in creating
444  ! mod_def file
445  !---------------------------------------------------------------------
446 
447  write(cnx,i4fmt)nx
448  write(form1,'(a)')'('//trim(cnx)//'f14.8)'
449  write(form2,'(a)')'('//trim(cnx)//'i2)'
450 
451  allocate(ww3mask(1:ni,1:nj)); ww3mask = int(wet4)
452  allocate(ww3dpth(1:ni,1:nj)); ww3dpth = dp4
453 
454  where(latct .ge. maximum_lat)ww3mask = 3
455  !close last row
456  ww3mask(:,nj) = 3
457 
458  open(unit=21,file=trim(dirout)//'/'//'ww3.mx'//trim(res)//'_x.inp',form='formatted')
459  open(unit=22,file=trim(dirout)//'/'//'ww3.mx'//trim(res)//'_y.inp',form='formatted')
460  open(unit=23,file=trim(dirout)//'/'//'ww3.mx'//trim(res)//'_bottom.inp',form='formatted')
461  open(unit=24,file=trim(dirout)//'/'//'ww3.mx'//trim(res)//'_mapsta.inp',form='formatted')
462  ! cice0 .ne. cicen requires obstruction map, should be initialized as zeros (w3grid,ln3032)
463  open(unit=25,file=trim(dirout)//'/'//'ww3.mx'//trim(res)//'_obstr.inp',form='formatted')
464 
465  do j = 1,nj
466  write( 21,trim(form1))lonct(:,j)
467  write( 22,trim(form1))latct(:,j)
468  end do
469  do j = 1,nj
470  write( 23,trim(form1))ww3dpth(:,j)
471  write( 24,trim(form2))ww3mask(:,j)
472  !'obsx' and 'obsy' arrays ???
473  write( 25,trim(form2))ww3mask(:,j)*0
474  write( 25,trim(form2))ww3mask(:,j)*0
475  end do
476 
477  close(21); close(22); close(23); close(24); close(25)
478  deallocate(ww3mask); deallocate(ww3dpth)
479  deallocate(wet4, wet8)
480 
481  !---------------------------------------------------------------------
482  ! use ESMF regridding to produce mapped ocean mask; first generate
483  ! conservative regrid weights from ocean to tiles; then generate the
484  ! tiled files containing the mapped ocean mask
485  !---------------------------------------------------------------------
486 
487  do n = 1,nar
488  npx = catm(n)
489  if (npx < 100) then
490  write(atmres,'(a,i2)')'C',npx
491  elseif (npx < 1000) then
492  write(atmres,'(a,i3)')'C',npx
493  else
494  write(atmres,'(a,i4)')'C',npx
495  end if
496 
497  method=esmf_regridmethod_conserve
498  fsrc = trim(dirout)//'/'//'Ct.mx'//trim(res)//'_SCRIP_land.nc'
499  fdst = trim(fv3dir)//'/'//trim(atmres)//'/'//trim(atmres)//'_mosaic.nc'
500  fwgt = trim(dirout)//'/'//'Ct.mx'//trim(res)//'.to.'//trim(atmres)//'.nc'
501  logmsg = 'creating weight file '//trim(fwgt)
502  print '(a)',trim(logmsg)
503 
504  call esmf_regridweightgen(srcfile=trim(fsrc),dstfile=trim(fdst), &
505  weightfile=trim(fwgt), regridmethod=method, &
506  unmappedaction=esmf_unmappedaction_ignore, ignoredegenerate=.true., &
507  netcdf4fileflag=.true., tilefilepath=trim(fv3dir)//'/'//trim(atmres)//'/', rc=rc)
508  if (esmf_logfounderror(rctocheck=rc, msg=esmf_logerr_passthru, &
509  line=__line__, file=__file__)) call esmf_finalize(endflag=esmf_end_abort)
510 
511  logmsg = 'creating mapped ocean mask for '//trim(atmres)
512  print '(a)',trim(logmsg)
513  call make_frac_land(trim(fsrc), trim(fwgt))
514  end do
515 
516  !---------------------------------------------------------------------
517  ! use ESMF to create positional weights for mapping a field from its
518  ! native stagger location (Cu,Cv,Bu) onto the center (Ct) grid location
519  !---------------------------------------------------------------------
520 
521  method=esmf_regridmethod_bilinear
522  fdst = trim(dirout)//'/'//'Ct.mx'//trim(res)//'_SCRIP.nc'
523  do k = 2,nv
524  cstagger = trim(staggerlocs(k))
525  fsrc = trim(dirout)//'/'//trim(cstagger)//'.mx'//trim(res)//'_SCRIP.nc'
526  fwgt = trim(dirout)//'/'//'tripole.mx'//trim(res)//'.'//trim(cstagger)//'.to.Ct.bilinear.nc'
527  logmsg = 'creating weight file '//trim(fwgt)
528  print '(a)',trim(logmsg)
529 
530  call esmf_regridweightgen(srcfile=trim(fsrc),dstfile=trim(fdst), &
531  weightfile=trim(fwgt), regridmethod=method, &
532  ignoredegenerate=.true., &
533  unmappedaction=esmf_unmappedaction_ignore, rc=rc)
534  if (esmf_logfounderror(rctocheck=rc, msg=esmf_logerr_passthru, &
535  line=__line__, file=__file__)) call esmf_finalize(endflag=esmf_end_abort)
536  end do
537 
538  !---------------------------------------------------------------------
539  ! use ESMF to create positional weights for mapping a field from the
540  ! center (Ct) grid location back to the native stagger location
541  ! (Cu,Cv,Bu). The destination is never mx025.
542  !---------------------------------------------------------------------
543 
544  if(trim(res) .ne. '025') then
545  method=esmf_regridmethod_bilinear
546  fsrc = trim(dirout)//'/'//'Ct.mx'//trim(res)//'_SCRIP.nc'
547  do k = 2,nv
548  cstagger = trim(staggerlocs(k))
549  fdst = trim(dirout)//'/'//trim(cstagger)//'.mx'//trim(res)//'_SCRIP.nc'
550  fwgt = trim(dirout)//'/'//'tripole.mx'//trim(res)//'.Ct.to.'//trim(cstagger)//'.bilinear.nc'
551  logmsg = 'creating weight file '//trim(fwgt)
552  print '(a)',trim(logmsg)
553 
554  call esmf_regridweightgen(srcfile=trim(fsrc),dstfile=trim(fdst), &
555  weightfile=trim(fwgt), regridmethod=method, &
556  ignoredegenerate=.true., &
557  unmappedaction=esmf_unmappedaction_ignore, rc=rc)
558  if (esmf_logfounderror(rctocheck=rc, msg=esmf_logerr_passthru, &
559  line=__line__, file=__file__)) call esmf_finalize(endflag=esmf_end_abort)
560  end do
561  end if
562 
563  !---------------------------------------------------------------------
564  !
565  !---------------------------------------------------------------------
566 
567  if(do_postwgts)call make_postwgts
568 
569  !---------------------------------------------------------------------
570  ! clean up
571  !---------------------------------------------------------------------
572 
573  deallocate(x, y, dx, dy)
574  deallocate(areact, anglet, angle, angchk)
575  deallocate(latct, lonct)
576  deallocate(latcv, loncv)
577  deallocate(latcu, loncu)
578  deallocate(latbu, lonbu)
579 
580 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