chgres_cube  1.12.0
 All Data Structures Files Functions Variables
model_grid.F90
Go to the documentation of this file.
1 
4 
9  module model_grid
10 
11  use esmf
12  use esmf_logpublicmod
13 
15  implicit none
16 
17  private
18 
19  character(len=5), allocatable, public :: tiles_target_grid(:)
20 
21  character(len=50), public :: input_grid_type = "latlon"
22 
23 
24  ! Made lsoil_target non-parameter to allow for RAP land surface initiation
25  integer, public :: lsoil_target = 4 ! # soil layers
26 
27  integer, public :: i_input
28 
30  integer, public :: j_input
31 
33  integer, public :: ip1_input
34 
35  integer, public :: jp1_input
36 
37  integer, public :: i_target
38 
40  integer, public :: j_target
41 
43  integer, public :: ip1_target
44 
45  integer, public :: jp1_target
46 
47  integer, public :: num_tiles_input_grid
48 
49  integer, public :: num_tiles_target_grid
50 
51 
52  type(esmf_grid), public :: input_grid
53 
54  type(esmf_grid), public :: target_grid
55 
56 
57  type(esmf_field), public :: latitude_input_grid
58 
59  type(esmf_field), public :: longitude_input_grid
60 
61  type(esmf_field), public :: latitude_s_input_grid
62 
64  type(esmf_field), public :: longitude_s_input_grid
65 
67  type(esmf_field), public :: latitude_w_input_grid
68 
70  type(esmf_field), public :: longitude_w_input_grid
71 
73 
74  type(esmf_field), public :: landmask_target_grid
75 
77  type(esmf_field), public :: latitude_target_grid
78 
79  type(esmf_field), public :: latitude_s_target_grid
80 
82  type(esmf_field), public :: latitude_w_target_grid
83 
85  type(esmf_field), public :: longitude_target_grid
86 
87  type(esmf_field), public :: longitude_s_target_grid
88 
90  type(esmf_field), public :: longitude_w_target_grid
91 
93  type(esmf_field), public :: seamask_target_grid
94 
96  type(esmf_field), public :: terrain_target_grid
97 
98 
99  public :: define_target_grid
100  public :: define_input_grid
102 
103  contains
104 
116  subroutine define_input_grid(localpet, npets)
117 
118  use program_setup, only : input_type
119 
120  implicit none
121 
122  integer, intent(in) :: localpet, npets
123 
124  if (trim(input_type) == "gaussian_nemsio" .or. &
125  trim(input_type) == "gfs_gaussian_nemsio" .or. &
126  trim(input_type) == "gfs_sigio" .or. &
127  trim(input_type) == "gaussian_netcdf") then
128  call define_input_grid_gaussian(npets)
129  elseif (trim(input_type) == "grib2") then
130  call define_input_grid_grib2(npets)
131  else
132  call define_input_grid_mosaic(localpet, npets)
133  endif
134 
135  end subroutine define_input_grid
136 
147  subroutine define_input_grid_gaussian(npets)
148 
149  use nemsio_module
150 
151  use program_setup, only : data_dir_input_grid, &
152  atm_files_input_grid, &
153  sfc_files_input_grid, &
154  input_type, &
155  convert_atm, convert_sfc
156 
157  use sfcio_module
158  use sigio_module
159  use netcdf
160 
161  implicit none
162 
163  integer, intent(in) :: npets
164 
165  character(len=250) :: the_file
166 
167  integer :: i, j, rc, clb(2), cub(2), ncid, id_grid
168  integer(sfcio_intkind) :: rc2
169  integer(sigio_intkind) :: rc3
170 
171  real(esmf_kind_r8), allocatable :: latitude(:,:)
172  real(esmf_kind_r8), allocatable :: longitude(:,:)
173  real(esmf_kind_r8), pointer :: lat_src_ptr(:,:)
174  real(esmf_kind_r8), pointer :: lon_src_ptr(:,:)
175  real(esmf_kind_r8), pointer :: lat_corner_src_ptr(:,:)
176  real(esmf_kind_r8), pointer :: lon_corner_src_ptr(:,:)
177  real(esmf_kind_r8) :: deltalon
178  real(esmf_kind_r8), allocatable :: slat(:), wlat(:)
179 
180  type(nemsio_gfile) :: gfile
181  type(esmf_polekind_flag) :: polekindflag(2)
182  type(sfcio_head) :: sfchead
183  type(sigio_head) :: sighead
184 
185  print*,"- DEFINE INPUT GRID OBJECT FOR GAUSSIAN DATA."
186 
187  num_tiles_input_grid = 1
188 
189  if (convert_sfc) then
190  the_file=trim(data_dir_input_grid) // "/" // trim(sfc_files_input_grid(1))
191  elseif (convert_atm) then
192  the_file=trim(data_dir_input_grid) // "/" // trim(atm_files_input_grid(1))
193  endif
194 
195  if (trim(input_type) == "gfs_sigio") then ! sigio/sfcio format, used by
196  ! spectral gfs prior to 7/19/2017.
197 
198  if (convert_sfc) then ! sfcio format
199  print*,"- OPEN AND READ ", trim(the_file)
200  call sfcio_sropen(21, trim(the_file), rc2)
201  if (rc2 /= 0) call error_handler("OPENING FILE", rc2)
202  call sfcio_srhead(21, sfchead, rc2)
203  if (rc2 /= 0) call error_handler("READING FILE", rc2)
204  call sfcio_sclose(21, rc2)
205  i_input = sfchead%lonb
206  j_input = sfchead%latb
207  elseif (convert_atm) then ! sigio format
208  print*,"- OPEN AND READ ", trim(the_file)
209  call sigio_sropen(21, trim(the_file), rc3)
210  if (rc3 /= 0) call error_handler("OPENING FILE", rc3)
211  call sigio_srhead(21, sighead, rc3)
212  if (rc3 /= 0) call error_handler("READING FILE", rc3)
213  call sigio_sclose(21, rc3)
214  i_input = sighead%lonb
215  j_input = sighead%latb
216  endif
217 
218  elseif (trim(input_type) == "gaussian_netcdf") then
219 
220  print*,'- OPEN AND READ: ',trim(the_file)
221  rc=nf90_open(trim(the_file),nf90_nowrite,ncid)
222  call netcdf_err(rc, 'opening file')
223 
224  print*,"- READ grid_xt"
225  rc=nf90_inq_dimid(ncid, 'grid_xt', id_grid)
226  call netcdf_err(rc, 'reading grid_xt id')
227  rc=nf90_inquire_dimension(ncid,id_grid,len=i_input)
228  call netcdf_err(rc, 'reading grid_xt')
229 
230  print*,"- READ grid_yt"
231  rc=nf90_inq_dimid(ncid, 'grid_yt', id_grid)
232  call netcdf_err(rc, 'reading grid_yt id')
233  rc=nf90_inquire_dimension(ncid,id_grid,len=j_input)
234  call netcdf_err(rc, 'reading grid_yt')
235 
236  rc = nf90_close(ncid)
237 
238  else ! nemsio format
239 
240  call nemsio_init(iret=rc)
241 
242  print*,"- OPEN AND READ ", trim(the_file)
243  call nemsio_open(gfile, the_file, "read", iret=rc)
244  if (rc /= 0) call error_handler("OPENING FILE", rc)
245 
246  call nemsio_getfilehead(gfile, iret=rc, dimx=i_input, dimy=j_input)
247  if (rc /= 0) call error_handler("READING FILE", rc)
248 
249  call nemsio_close(gfile)
250 
251  endif
252 
253  ip1_input = i_input + 1
254  jp1_input = j_input + 1
255 
256  polekindflag(1:2) = esmf_polekind_monopole
257 
258  print*,"- CALL GridCreate1PeriDim FOR INPUT GRID."
259  input_grid = esmf_gridcreate1peridim(minindex=(/1,1/), &
260  maxindex=(/i_input,j_input/), &
261  polekindflag=polekindflag, &
262  periodicdim=1, &
263  poledim=2, &
264  coordsys=esmf_coordsys_sph_deg, &
265  regdecomp=(/1,npets/), &
266  indexflag=esmf_index_global, rc=rc)
267  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
268  call error_handler("IN GridCreate1PeriDim", rc)
269 
270  print*,"- CALL FieldCreate FOR INPUT GRID LATITUDE."
271  latitude_input_grid = esmf_fieldcreate(input_grid, &
272  typekind=esmf_typekind_r8, &
273  staggerloc=esmf_staggerloc_center, &
274  name="input_grid_latitude", rc=rc)
275  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
276  call error_handler("IN FieldCreate", rc)
277 
278  print*,"- CALL FieldCreate FOR INPUT GRID LONGITUDE."
279  longitude_input_grid = esmf_fieldcreate(input_grid, &
280  typekind=esmf_typekind_r8, &
281  staggerloc=esmf_staggerloc_center, &
282  name="input_grid_longitude", rc=rc)
283  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
284  call error_handler("IN FieldCreate", rc)
285 
286  allocate(longitude(i_input,j_input))
287  allocate(latitude(i_input,j_input))
288 
289  deltalon = 360.0_esmf_kind_r8 / real(i_input,kind=esmf_kind_r8)
290  do i = 1, i_input
291  longitude(i,:) = real((i-1),kind=esmf_kind_r8) * deltalon
292  enddo
293 
294  allocate(slat(j_input))
295  allocate(wlat(j_input))
296  call splat(4, j_input, slat, wlat)
297 
298  do i = 1, j_input
299  latitude(:,i) = 90.0_esmf_kind_r8 - (acos(slat(i))* 180.0_esmf_kind_r8 / &
300  (4.0_esmf_kind_r8*atan(1.0_esmf_kind_r8)))
301  enddo
302 
303  deallocate(slat, wlat)
304 
305  print*,"- CALL FieldScatter FOR INPUT GRID LONGITUDE."
306  call esmf_fieldscatter(longitude_input_grid, longitude, rootpet=0, rc=rc)
307  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
308  call error_handler("IN FieldScatter", rc)
309 
310  print*,"- CALL FieldScatter FOR INPUT GRID LATITUDE."
311  call esmf_fieldscatter(latitude_input_grid, latitude, rootpet=0, rc=rc)
312  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
313  call error_handler("IN FieldScatter", rc)
314 
315  print*,"- CALL GridAddCoord FOR INPUT GRID."
316  call esmf_gridaddcoord(input_grid, &
317  staggerloc=esmf_staggerloc_center, rc=rc)
318  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
319  call error_handler("IN GridAddCoord", rc)
320 
321  print*,"- CALL GridGetCoord FOR INPUT GRID X-COORD."
322  nullify(lon_src_ptr)
323  call esmf_gridgetcoord(input_grid, &
324  staggerloc=esmf_staggerloc_center, &
325  coorddim=1, &
326  farrayptr=lon_src_ptr, rc=rc)
327  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
328  call error_handler("IN GridGetCoord", rc)
329 
330  print*,"- CALL GridGetCoord FOR INPUT GRID Y-COORD."
331  nullify(lat_src_ptr)
332  call esmf_gridgetcoord(input_grid, &
333  staggerloc=esmf_staggerloc_center, &
334  coorddim=2, &
335  computationallbound=clb, &
336  computationalubound=cub, &
337  farrayptr=lat_src_ptr, rc=rc)
338  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
339  call error_handler("IN GridGetCoord", rc)
340 
341  do j = clb(2), cub(2)
342  do i = clb(1), cub(1)
343  lon_src_ptr(i,j) = longitude(i,j)
344  if (lon_src_ptr(i,j) > 360.0_esmf_kind_r8) lon_src_ptr(i,j) = lon_src_ptr(i,j) - 360.0_esmf_kind_r8
345  lat_src_ptr(i,j) = latitude(i,j)
346  enddo
347  enddo
348 
349  print*,"- CALL GridAddCoord FOR INPUT GRID."
350  call esmf_gridaddcoord(input_grid, &
351  staggerloc=esmf_staggerloc_corner, rc=rc)
352  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
353  call error_handler("IN GridAddCoord", rc)
354 
355  print*,"- CALL GridGetCoord FOR INPUT GRID X-COORD."
356  nullify(lon_corner_src_ptr)
357  call esmf_gridgetcoord(input_grid, &
358  staggerloc=esmf_staggerloc_corner, &
359  coorddim=1, &
360  farrayptr=lon_corner_src_ptr, rc=rc)
361  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
362  call error_handler("IN GridGetCoord", rc)
363 
364  print*,"- CALL GridGetCoord FOR INPUT GRID Y-COORD."
365  nullify(lat_corner_src_ptr)
366  call esmf_gridgetcoord(input_grid, &
367  staggerloc=esmf_staggerloc_corner, &
368  coorddim=2, &
369  computationallbound=clb, &
370  computationalubound=cub, &
371  farrayptr=lat_corner_src_ptr, rc=rc)
372  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
373  call error_handler("IN GridGetCoord", rc)
374 
375  do j = clb(2), cub(2)
376  do i = clb(1), cub(1)
377  lon_corner_src_ptr(i,j) = longitude(i,1) - (0.5_esmf_kind_r8*deltalon)
378  if (lon_corner_src_ptr(i,j) > 360.0_esmf_kind_r8) lon_corner_src_ptr(i,j) = lon_corner_src_ptr(i,j) - 360.0_esmf_kind_r8
379  if (j == 1) then
380  lat_corner_src_ptr(i,j) = 90.0_esmf_kind_r8
381  cycle
382  endif
383  if (j == jp1_input) then
384  lat_corner_src_ptr(i,j) = -90.0_esmf_kind_r8
385  cycle
386  endif
387  lat_corner_src_ptr(i,j) = 0.5_esmf_kind_r8 * (latitude(i,j-1)+ latitude(i,j))
388  enddo
389  enddo
390 
391  deallocate(latitude,longitude)
392 
393  end subroutine define_input_grid_gaussian
394 
401  subroutine define_input_grid_mosaic(localpet, npets)
402 
403  use netcdf
404  use program_setup, only : mosaic_file_input_grid, &
405  orog_dir_input_grid, &
406  orog_files_input_grid
407 
408  implicit none
409 
410  character(len=500) :: the_file
411 
412  integer, intent(in) :: localpet, npets
413 
414  integer :: id_tiles, id_dim, tile
415  integer :: extra, error, ncid
416  integer, allocatable :: decomptile(:,:)
417 
418  integer(esmf_kind_i8), allocatable :: landmask_one_tile(:,:)
419 
420  real(esmf_kind_r8), allocatable :: latitude_one_tile(:,:)
421  real(esmf_kind_r8), allocatable :: latitude_s_one_tile(:,:)
422  real(esmf_kind_r8), allocatable :: latitude_w_one_tile(:,:)
423  real(esmf_kind_r8), allocatable :: longitude_one_tile(:,:)
424  real(esmf_kind_r8), allocatable :: longitude_s_one_tile(:,:)
425  real(esmf_kind_r8), allocatable :: longitude_w_one_tile(:,:)
426 
427  print*,'- OPEN INPUT GRID MOSAIC FILE: ',trim(mosaic_file_input_grid)
428  error=nf90_open(trim(mosaic_file_input_grid),nf90_nowrite,ncid)
429  call netcdf_err(error, 'opening grid mosaic file')
430 
431  print*,"- READ NUMBER OF TILES"
432  error=nf90_inq_dimid(ncid, 'ntiles', id_tiles)
433  call netcdf_err(error, 'reading ntiles id')
434  error=nf90_inquire_dimension(ncid,id_tiles,len=num_tiles_input_grid)
435  call netcdf_err(error, 'reading ntiles')
436 
437  error = nf90_close(ncid)
438 
439  print*,'- NUMBER OF TILES, INPUT MODEL GRID IS ', num_tiles_input_grid
440 
441  if (mod(npets,num_tiles_input_grid) /= 0) then
442  call error_handler("MUST RUN WITH A TASK COUNT THAT IS A MULTIPLE OF 6.", 1)
443  endif
444 
445 !-----------------------------------------------------------------------
446 ! Create ESMF grid object for the model grid.
447 !-----------------------------------------------------------------------
448 
449  extra = npets / num_tiles_input_grid
450 
451  allocate(decomptile(2,num_tiles_input_grid))
452 
453  do tile = 1, num_tiles_input_grid
454  decomptile(:,tile)=(/1,extra/)
455  enddo
456 
457  print*,"- CALL GridCreateMosaic FOR INPUT MODEL GRID"
458  input_grid = esmf_gridcreatemosaic(filename=trim(mosaic_file_input_grid), &
459  regdecompptile=decomptile, &
460  staggerloclist=(/esmf_staggerloc_center, esmf_staggerloc_corner, &
461  esmf_staggerloc_edge1, esmf_staggerloc_edge2/), &
462  indexflag=esmf_index_global, &
463  tilefilepath=trim(orog_dir_input_grid), &
464  rc=error)
465  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
466  call error_handler("IN GridCreateMosaic", error)
467 
468 !-----------------------------------------------------------------------
469 ! Read the mask and lat/lons.
470 !-----------------------------------------------------------------------
471 
472  print*,"- CALL FieldCreate FOR INPUT GRID LATITUDE."
473  latitude_input_grid = esmf_fieldcreate(input_grid, &
474  typekind=esmf_typekind_r8, &
475  staggerloc=esmf_staggerloc_center, &
476  name="input_grid_latitude", &
477  rc=error)
478  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
479  call error_handler("IN FieldCreate", error)
480 
481  print*,"- CALL FieldCreate FOR INPUT GRID LONGITUDE."
482  longitude_input_grid = esmf_fieldcreate(input_grid, &
483  typekind=esmf_typekind_r8, &
484  staggerloc=esmf_staggerloc_center, &
485  name="input_grid_longitude", &
486  rc=error)
487  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
488  call error_handler("IN FieldCreate", error)
489 
490  print*,"- CALL FieldCreate FOR INPUT GRID LATITUDE_S."
491  latitude_s_input_grid = esmf_fieldcreate(input_grid, &
492  typekind=esmf_typekind_r8, &
493  staggerloc=esmf_staggerloc_edge2, &
494  name="input_grid_latitude_s", &
495  rc=error)
496  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
497  call error_handler("IN FieldCreate", error)
498 
499  print*,"- CALL FieldCreate FOR INPUT GRID LONGITUDE_S."
500  longitude_s_input_grid = esmf_fieldcreate(input_grid, &
501  typekind=esmf_typekind_r8, &
502  staggerloc=esmf_staggerloc_edge2, &
503  name="input_grid_longitude_s", &
504  rc=error)
505  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
506  call error_handler("IN FieldCreate", error)
507 
508  print*,"- CALL FieldCreate FOR INPUT GRID LATITUDE_W."
509  latitude_w_input_grid = esmf_fieldcreate(input_grid, &
510  typekind=esmf_typekind_r8, &
511  staggerloc=esmf_staggerloc_edge1, &
512  name="input_grid_latitude_w", &
513  rc=error)
514  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
515  call error_handler("IN FieldCreate", error)
516 
517  print*,"- CALL FieldCreate FOR INPUT GRID LONGITUDE_W."
518  longitude_w_input_grid = esmf_fieldcreate(input_grid, &
519  typekind=esmf_typekind_r8, &
520  staggerloc=esmf_staggerloc_edge1, &
521  name="input_grid_longitude_w", &
522  rc=error)
523  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
524  call error_handler("IN FieldCreate", error)
525 
526  the_file = trim(orog_dir_input_grid) // trim(orog_files_input_grid(1))
527 
528  print*,'- OPEN FIRST INPUT GRID OROGRAPHY FILE: ',trim(the_file)
529  error=nf90_open(trim(the_file),nf90_nowrite,ncid)
530  call netcdf_err(error, 'opening ororgraphy file')
531  print*,"- READ GRID DIMENSIONS"
532  error=nf90_inq_dimid(ncid, 'lon', id_dim)
533  call netcdf_err(error, 'reading lon id')
534  error=nf90_inquire_dimension(ncid,id_dim,len=i_input)
535  call netcdf_err(error, 'reading lon')
536  error=nf90_inq_dimid(ncid, 'lat', id_dim)
537  call netcdf_err(error, 'reading lat id')
538  error=nf90_inquire_dimension(ncid,id_dim,len=j_input)
539  call netcdf_err(error, 'reading lat')
540  error = nf90_close(ncid)
541 
542  print*,"- I/J DIMENSIONS OF THE INPUT GRID TILES ", i_input, j_input
543 
544  ip1_input = i_input + 1
545  jp1_input = j_input + 1
546 
547  if (localpet == 0) then
548  allocate(longitude_one_tile(i_input,j_input))
549  allocate(longitude_s_one_tile(i_input,jp1_input))
550  allocate(longitude_w_one_tile(ip1_input,j_input))
551  allocate(latitude_one_tile(i_input,j_input))
552  allocate(latitude_s_one_tile(i_input,jp1_input))
553  allocate(latitude_w_one_tile(ip1_input,j_input))
554  allocate(landmask_one_tile(i_input,j_input))
555  else
556  allocate(longitude_one_tile(0,0))
557  allocate(longitude_s_one_tile(0,0))
558  allocate(longitude_w_one_tile(0,0))
559  allocate(latitude_one_tile(0,0))
560  allocate(latitude_s_one_tile(0,0))
561  allocate(latitude_w_one_tile(0,0))
562  allocate(landmask_one_tile(0,0))
563  endif
564 
565  do tile = 1, num_tiles_input_grid
566  if (localpet == 0) then
567  call get_model_latlons(mosaic_file_input_grid, orog_dir_input_grid, num_tiles_input_grid, tile, &
568  i_input, j_input, ip1_input, jp1_input, latitude_one_tile, &
569  latitude_s_one_tile, latitude_w_one_tile, longitude_one_tile, &
570  longitude_s_one_tile, longitude_w_one_tile)
571  endif
572  print*,"- CALL FieldScatter FOR INPUT GRID LATITUDE. TILE IS: ", tile
573  call esmf_fieldscatter(latitude_input_grid, latitude_one_tile, rootpet=0, tile=tile, rc=error)
574  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
575  call error_handler("IN FieldScatter", error)
576  print*,"- CALL FieldScatter FOR INPUT GRID LONGITUDE. TILE IS: ", tile
577  call esmf_fieldscatter(longitude_input_grid, longitude_one_tile, rootpet=0, tile=tile, rc=error)
578  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
579  call error_handler("IN FieldScatter", error)
580  print*,"- CALL FieldScatter FOR INPUT GRID LATITUDE_S. TILE IS: ", tile
581  call esmf_fieldscatter(latitude_s_input_grid, latitude_s_one_tile, rootpet=0, tile=tile, rc=error)
582  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
583  call error_handler("IN FieldScatter", error)
584  print*,"- CALL FieldScatter FOR INPUT GRID LONGITUDE_S. TILE IS: ", tile
585  call esmf_fieldscatter(longitude_s_input_grid, longitude_s_one_tile, rootpet=0, tile=tile, rc=error)
586  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
587  call error_handler("IN FieldScatter", error)
588  print*,"- CALL FieldScatter FOR INPUT GRID LATITUDE_W. TILE IS: ", tile
589  call esmf_fieldscatter(latitude_w_input_grid, latitude_w_one_tile, rootpet=0, tile=tile, rc=error)
590  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
591  call error_handler("IN FieldScatter", error)
592  print*,"- CALL FieldScatter FOR INPUT GRID LONGITUDE_W. TILE IS: ", tile
593  call esmf_fieldscatter(longitude_w_input_grid, longitude_w_one_tile, rootpet=0, tile=tile, rc=error)
594  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
595  call error_handler("IN FieldScatter", error)
596  enddo
597 
598  deallocate(longitude_one_tile)
599  deallocate(longitude_s_one_tile)
600  deallocate(longitude_w_one_tile)
601  deallocate(latitude_one_tile)
602  deallocate(latitude_s_one_tile)
603  deallocate(latitude_w_one_tile)
604  deallocate(landmask_one_tile)
605 
606  end subroutine define_input_grid_mosaic
607 
614  subroutine define_input_grid_grib2(npets)
615 
616  use grib_mod
617  use gdswzd_mod
618  use program_setup, only : grib2_file_input_grid, data_dir_input_grid
619 
620  implicit none
621 
622  integer, intent(in) :: npets
623 
624  character(len=500) :: the_file
625 
626  integer :: i, j, k, jdisc, jgdtn, jpdtn, lugb, lugi
627  integer :: jids(200), jgdt(200), jpdt(200), rc
628  integer :: kgds(200), nret, clb(2), cub(2)
629 
630  logical :: unpack
631 
632  real :: res
633  real, allocatable :: rlon(:,:),rlat(:,:),xpts(:,:),ypts(:,:)
634  real, allocatable :: rlon_corner(:,:),rlat_corner(:,:)
635  real, allocatable :: rlon_diff(:,:),rlat_diff(:,:)
636  real, allocatable :: xpts_corner(:,:),ypts_corner(:,:)
637  real(esmf_kind_r8), allocatable :: latitude(:,:)
638  real(esmf_kind_r8), allocatable :: longitude(:,:)
639  real(esmf_kind_r8), allocatable :: latitude_corner(:,:)
640  real(esmf_kind_r8), allocatable :: longitude_corner(:,:)
641  real(esmf_kind_r8), pointer :: lat_src_ptr(:,:)
642  real(esmf_kind_r8), pointer :: lat_corner_src_ptr(:,:)
643  real(esmf_kind_r8), pointer :: lon_src_ptr(:,:)
644  real(esmf_kind_r8), pointer :: lon_corner_src_ptr(:,:)
645 
646  type(esmf_polekind_flag) :: polekindflag(2)
647 
648  type(gribfield) :: gfld
649 
650  the_file = trim(data_dir_input_grid) // "/" // grib2_file_input_grid
651 
652  lugb=12
653 
654  print*,"- OPEN AND READ INPUT DATA GRIB2 FILE: ", trim(the_file)
655  call baopenr(lugb,the_file,rc)
656  if (rc /= 0) call error_handler("OPENING FILE", rc)
657 
658 ! Read the first record and get the grid definition template.
659 
660  j = 0 ! Search at beginning of file
661  lugi = 0 ! No grib index file
662  jdisc = -1 ! Search for any discipline
663  jpdtn = -1 ! Search for any product definition template number
664  jgdtn = -1 ! Search for any grid definition template number
665  jids = -9999 ! Array of values in identification section, set to wildcard.
666  jgdt = -9999 ! Array of values in grid definition template, set to wildcard.
667  jpdt = -9999 ! Array of values in product definition template, set to wildcard.
668  unpack = .false. ! unpack data
669 
670  call getgb2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
671  unpack, k, gfld, rc)
672  if (rc /= 0) call error_handler("DEGRIBBING INPUT FILE.", rc)
673 
674  call baclose(lugb,rc)
675 
676  if (gfld%igdtnum == 0) then
677  print*,"- INPUT DATA ON LAT/LON GRID."
678  input_grid_type = 'latlon'
679  elseif (gfld%igdtnum == 30) then
680  print*,"- INPUT DATA ON LAMBERT CONFORMAL GRID."
681  input_grid_type = 'lambert'
682  elseif (gfld%igdtnum == 32769) then
683  print*,"- INPUT DATA ON ROTATED LAT/LON GRID."
684  input_grid_type = 'rotated_latlon'
685  else
686  call error_handler("INPUT GRID TEMPLATE NOT SUPPORTED.", 2)
687  endif
688 
689  kgds = 0
690  call gdt_to_gds(gfld%igdtnum, gfld%igdtmpl, gfld%igdtlen, kgds, i_input, j_input, res)
691 
692  ip1_input = i_input + 1
693  jp1_input = j_input + 1
694 
695  allocate(rlat(i_input,j_input))
696  allocate(rlon(i_input,j_input))
697  allocate(rlat_diff(i_input,j_input))
698  allocate(rlon_diff(i_input,j_input))
699  allocate(xpts(i_input,j_input))
700  allocate(ypts(i_input,j_input))
701  allocate(rlat_corner(ip1_input,jp1_input))
702  allocate(rlon_corner(ip1_input,jp1_input))
703  allocate(xpts_corner(ip1_input,jp1_input))
704  allocate(ypts_corner(ip1_input,jp1_input))
705 
706  do j = 1, j_input
707  do i = 1, i_input
708  xpts(i,j) = float(i)
709  ypts(i,j) = float(j)
710  enddo
711  enddo
712 
713  print*,"- COMPUTE GRID CELL CENTER COORDINATES."
714  call gdswzd(kgds,1,(i_input*j_input),-9999.,xpts,ypts,rlon,rlat,nret)
715 
716  if (nret /= (i_input*j_input)) then
717  call error_handler("GDSWZD RETURNED WRONG NUMBER OF POINTS.", 2)
718  endif
719 
720  deallocate(xpts, ypts)
721 
722  do j = 1, jp1_input
723  do i = 1, ip1_input
724  xpts_corner(i,j) = float(i) - 0.5
725  ypts_corner(i,j) = float(j) - 0.5
726  enddo
727  enddo
728 
729  print*,"- COMPUTE GRID CELL CORNER COORDINATES."
730  call gdswzd(kgds,1,(ip1_input*jp1_input),-9999.,xpts_corner,ypts_corner,rlon_corner,rlat_corner,nret)
731 
732  if (nret /= (ip1_input*jp1_input)) then
733  call error_handler("GDSWZD RETURNED WRONG NUMBER OF POINTS.", 2)
734  endif
735 
736  deallocate(xpts_corner, ypts_corner)
737 
738  if (gfld%igdtnum == 0) then ! gfs lat/lon data
739 
740  print*,"- CALL GridCreate1PeriDim FOR INPUT GRID."
741 
742  polekindflag(1:2) = esmf_polekind_monopole
743 
744  input_grid = esmf_gridcreate1peridim(minindex=(/1,1/), &
745  maxindex=(/i_input,j_input/), &
746  polekindflag=polekindflag, &
747  periodicdim=1, &
748  poledim=2, &
749  coordsys=esmf_coordsys_sph_deg, &
750  regdecomp=(/1,npets/), &
751  indexflag=esmf_index_global, rc=rc)
752 
753  else
754 
755  print*,"- CALL GridCreateNoPeriDim FOR INPUT GRID."
756 
757  input_grid = esmf_gridcreatenoperidim(maxindex=(/i_input,j_input/), &
758  indexflag=esmf_index_global, &
759  rc=rc)
760 
761  endif
762 
763  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
764  call error_handler("IN GridCreate1PeriDim", rc)
765 
766  print*,"- CALL FieldCreate FOR INPUT GRID LATITUDE."
767  latitude_input_grid = esmf_fieldcreate(input_grid, &
768  typekind=esmf_typekind_r8, &
769  staggerloc=esmf_staggerloc_center, &
770  name="input_grid_latitude", rc=rc)
771  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
772  call error_handler("IN FieldCreate", rc)
773 
774  print*,"- CALL FieldCreate FOR INPUT GRID LONGITUDE."
775  longitude_input_grid = esmf_fieldcreate(input_grid, &
776  typekind=esmf_typekind_r8, &
777  staggerloc=esmf_staggerloc_center, &
778  name="input_grid_longitude", rc=rc)
779  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
780  call error_handler("IN FieldCreate", rc)
781 
782  allocate(latitude(i_input,j_input))
783  allocate(longitude(i_input,j_input))
784 
785  latitude = rlat
786  longitude = rlon
787 
788  deallocate (rlat, rlon)
789 
790  print*,"- CALL FieldScatter FOR INPUT GRID LONGITUDE."
791  call esmf_fieldscatter(longitude_input_grid, longitude, rootpet=0, rc=rc)
792  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
793  call error_handler("IN FieldScatter", rc)
794 
795  print*,"- CALL FieldScatter FOR INPUT GRID LATITUDE."
796  call esmf_fieldscatter(latitude_input_grid, latitude, rootpet=0, rc=rc)
797  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
798  call error_handler("IN FieldScatter", rc)
799 
800  print*,"- CALL GridAddCoord FOR INPUT GRID."
801  call esmf_gridaddcoord(input_grid, &
802  staggerloc=esmf_staggerloc_center, rc=rc)
803  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
804  call error_handler("IN GridAddCoord", rc)
805 
806  print*,"- CALL GridGetCoord FOR INPUT GRID X-COORD."
807  nullify(lon_src_ptr)
808  call esmf_gridgetcoord(input_grid, &
809  staggerloc=esmf_staggerloc_center, &
810  coorddim=1, &
811  farrayptr=lon_src_ptr, rc=rc)
812  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
813  call error_handler("IN GridGetCoord", rc)
814 
815  print*,"- CALL GridGetCoord FOR INPUT GRID Y-COORD."
816  nullify(lat_src_ptr)
817  call esmf_gridgetcoord(input_grid, &
818  staggerloc=esmf_staggerloc_center, &
819  coorddim=2, &
820  computationallbound=clb, &
821  computationalubound=cub, &
822  farrayptr=lat_src_ptr, rc=rc)
823  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
824  call error_handler("IN GridGetCoord", rc)
825 
826  do j = clb(2), cub(2)
827  do i = clb(1), cub(1)
828  lon_src_ptr(i,j) = longitude(i,j)
829  if (lon_src_ptr(i,j) > 360.0_esmf_kind_r8) lon_src_ptr(i,j) = lon_src_ptr(i,j) - 360.0_esmf_kind_r8
830  lat_src_ptr(i,j) = latitude(i,j)
831  enddo
832  enddo
833 
834  deallocate(latitude, longitude)
835 
836  allocate(latitude_corner(ip1_input,jp1_input))
837  allocate(longitude_corner(ip1_input,jp1_input))
838 
839  latitude_corner = rlat_corner
840  longitude_corner = rlon_corner
841 
842  deallocate (rlat_corner, rlon_corner)
843 
844  print*,"- CALL GridAddCoord FOR INPUT GRID."
845  call esmf_gridaddcoord(input_grid, &
846  staggerloc=esmf_staggerloc_corner, rc=rc)
847  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
848  call error_handler("IN GridAddCoord", rc)
849 
850  print*,"- CALL GridGetCoord FOR INPUT GRID X-COORD."
851  nullify(lon_corner_src_ptr)
852  call esmf_gridgetcoord(input_grid, &
853  staggerloc=esmf_staggerloc_corner, &
854  coorddim=1, &
855  farrayptr=lon_corner_src_ptr, rc=rc)
856  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
857  call error_handler("IN GridGetCoord", rc)
858 
859  print*,"- CALL GridGetCoord FOR INPUT GRID Y-COORD."
860  nullify(lat_corner_src_ptr)
861  call esmf_gridgetcoord(input_grid, &
862  staggerloc=esmf_staggerloc_corner, &
863  coorddim=2, &
864  computationallbound=clb, &
865  computationalubound=cub, &
866  farrayptr=lat_corner_src_ptr, rc=rc)
867  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
868  call error_handler("IN GridGetCoord", rc)
869 
870  do j = clb(2), cub(2)
871  do i = clb(1), cub(1)
872  lon_corner_src_ptr(i,j) = longitude_corner(i,j)
873  if (lon_corner_src_ptr(i,j) > 360.0_esmf_kind_r8) lon_corner_src_ptr(i,j) = lon_corner_src_ptr(i,j) - 360.0_esmf_kind_r8
874  lat_corner_src_ptr(i,j) = latitude_corner(i,j)
875  enddo
876  enddo
877 
878  deallocate(latitude_corner, longitude_corner)
879 
880  end subroutine define_input_grid_grib2
881 
887  subroutine define_target_grid(localpet, npets)
888 
889  use netcdf
890  use program_setup, only : mosaic_file_target_grid, &
891  orog_dir_target_grid, &
892  orog_files_target_grid, &
893  nsoill_out
894 
895  implicit none
896 
897  integer, intent(in) :: localpet, npets
898 
899  character(len=500) :: the_file
900 
901  integer :: error, ncid, extra
902  integer :: id_tiles
903  integer :: id_dim, id_grid_tiles
904  integer :: tile
905  integer, allocatable :: decomptile(:,:)
906  integer(esmf_kind_i8), allocatable :: landmask_one_tile(:,:)
907  integer(esmf_kind_i8), allocatable :: seamask_one_tile(:,:)
908 
909  real(esmf_kind_r8), allocatable :: latitude_one_tile(:,:)
910  real(esmf_kind_r8), allocatable :: latitude_s_one_tile(:,:)
911  real(esmf_kind_r8), allocatable :: latitude_w_one_tile(:,:)
912  real(esmf_kind_r8), allocatable :: longitude_one_tile(:,:)
913  real(esmf_kind_r8), allocatable :: longitude_s_one_tile(:,:)
914  real(esmf_kind_r8), allocatable :: longitude_w_one_tile(:,:)
915  real(esmf_kind_r8), allocatable :: terrain_one_tile(:,:)
916 
917  lsoil_target = nsoill_out
918 
919  print*,'- OPEN TARGET GRID MOSAIC FILE: ',trim(mosaic_file_target_grid)
920  error=nf90_open(trim(mosaic_file_target_grid),nf90_nowrite,ncid)
921  call netcdf_err(error, 'opening grid mosaic file')
922 
923  print*,"- READ NUMBER OF TILES"
924  error=nf90_inq_dimid(ncid, 'ntiles', id_tiles)
925  call netcdf_err(error, 'reading ntile id')
926  error=nf90_inquire_dimension(ncid,id_tiles,len=num_tiles_target_grid)
927  call netcdf_err(error, 'reading ntiles')
928  error=nf90_inq_varid(ncid, 'gridtiles', id_grid_tiles)
929  call netcdf_err(error, 'reading gridtiles id')
930  allocate(tiles_target_grid(num_tiles_target_grid))
931  tiles_target_grid="NULL"
932  print*,"- READ TILE NAMES"
933  error=nf90_get_var(ncid, id_grid_tiles, tiles_target_grid)
934  call netcdf_err(error, 'reading gridtiles')
935 
936  error = nf90_close(ncid)
937 
938  print*,'- NUMBER OF TILES, TARGET MODEL GRID IS ', num_tiles_target_grid
939 
940  if (mod(npets,num_tiles_target_grid) /= 0) then
941  call error_handler("MUST RUN WITH TASK COUNT THAT IS A MULTIPLE OF # OF TILES.", 1)
942  endif
943 
944 !-----------------------------------------------------------------------
945 ! Get the model grid specs and land mask from the orography files.
946 !-----------------------------------------------------------------------
947 
948  the_file = trim(orog_dir_target_grid) // trim(orog_files_target_grid(1))
949 
950  print*,'- OPEN FIRST TARGET GRID OROGRAPHY FILE: ',trim(the_file)
951  error=nf90_open(trim(the_file),nf90_nowrite,ncid)
952  call netcdf_err(error, 'opening orography file')
953  print*,"- READ GRID DIMENSIONS"
954  error=nf90_inq_dimid(ncid, 'lon', id_dim)
955  call netcdf_err(error, 'reading lon id')
956  error=nf90_inquire_dimension(ncid,id_dim,len=i_target)
957  call netcdf_err(error, 'reading lon')
958  error=nf90_inq_dimid(ncid, 'lat', id_dim)
959  call netcdf_err(error, 'reading lat id')
960  error=nf90_inquire_dimension(ncid,id_dim,len=j_target)
961  call netcdf_err(error, 'reading lat')
962  error = nf90_close(ncid)
963 
964  print*,"- I/J DIMENSIONS OF THE TARGET GRID TILES ", i_target, j_target
965 
966  ip1_target = i_target + 1
967  jp1_target = j_target + 1
968 
969 !-----------------------------------------------------------------------
970 ! Create ESMF grid object for the model grid.
971 !-----------------------------------------------------------------------
972 
973  extra = npets / num_tiles_target_grid
974 
975  allocate(decomptile(2,num_tiles_target_grid))
976 
977  do tile = 1, num_tiles_target_grid
978  decomptile(:,tile)=(/1,extra/)
979  enddo
980 
981  print*,"- CALL GridCreateMosaic FOR TARGET GRID"
982  target_grid = esmf_gridcreatemosaic(filename=trim(mosaic_file_target_grid), &
983  regdecompptile=decomptile, &
984  staggerloclist=(/esmf_staggerloc_center, esmf_staggerloc_corner, &
985  esmf_staggerloc_edge1, esmf_staggerloc_edge2/), &
986  indexflag=esmf_index_global, &
987  tilefilepath=trim(orog_dir_target_grid), rc=error)
988  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
989  call error_handler("IN GridCreateMosaic", error)
990 
991 !-----------------------------------------------------------------------
992 ! Set target model landmask (1 - land, 0 - not land) and
993 ! seamask (1 - non-land, 0 -land). Read lat/lon on target grid.
994 !-----------------------------------------------------------------------
995 
996  print*,"- CALL FieldCreate FOR TARGET GRID LANDMASK."
997  landmask_target_grid = esmf_fieldcreate(target_grid, &
998  typekind=esmf_typekind_i8, &
999  staggerloc=esmf_staggerloc_center, &
1000  name="target_grid_landmask", rc=error)
1001  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1002  call error_handler("IN FieldCreate", error)
1003 
1004  print*,"- CALL FieldCreate FOR TARGET GRID SEAMASK."
1005  seamask_target_grid = esmf_fieldcreate(target_grid, &
1006  typekind=esmf_typekind_i8, &
1007  staggerloc=esmf_staggerloc_center, &
1008  name="target_grid_seamask", rc=error)
1009  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1010  call error_handler("IN FieldCreate", error)
1011 
1012  print*,"- CALL FieldCreate FOR TARGET GRID LATITUDE."
1013  latitude_target_grid = esmf_fieldcreate(target_grid, &
1014  typekind=esmf_typekind_r8, &
1015  staggerloc=esmf_staggerloc_center, &
1016  name="target_grid_latitude", rc=error)
1017  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1018  call error_handler("IN FieldCreate", error)
1019 
1020  print*,"- CALL FieldCreate FOR TARGET GRID LATITUDE_S."
1021  latitude_s_target_grid = esmf_fieldcreate(target_grid, &
1022  typekind=esmf_typekind_r8, &
1023  staggerloc=esmf_staggerloc_edge2, &
1024  name="target_grid_latitude_s", rc=error)
1025  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1026  call error_handler("IN FieldCreate", error)
1027 
1028  print*,"- CALL FieldCreate FOR TARGET GRID LATITUDE_W."
1029  latitude_w_target_grid = esmf_fieldcreate(target_grid, &
1030  typekind=esmf_typekind_r8, &
1031  staggerloc=esmf_staggerloc_edge1, &
1032  name="target_grid_latitude_w", &
1033  rc=error)
1034  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1035  call error_handler("IN FieldCreate", error)
1036 
1037  print*,"- CALL FieldCreate FOR TARGET GRID LONGITUDE."
1038  longitude_target_grid = esmf_fieldcreate(target_grid, &
1039  typekind=esmf_typekind_r8, &
1040  staggerloc=esmf_staggerloc_center, &
1041  name="target_grid_longitude", &
1042  rc=error)
1043  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1044  call error_handler("IN FieldCreate", error)
1045 
1046  print*,"- CALL FieldCreate FOR TARGET GRID LONGITUDE_S."
1047  longitude_s_target_grid = esmf_fieldcreate(target_grid, &
1048  typekind=esmf_typekind_r8, &
1049  staggerloc=esmf_staggerloc_edge2, &
1050  name="target_grid_longitude_s", &
1051  rc=error)
1052  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1053  call error_handler("IN FieldCreate", error)
1054 
1055  print*,"- CALL FieldCreate FOR TARGET GRID LONGITUDE_W."
1056  longitude_w_target_grid = esmf_fieldcreate(target_grid, &
1057  typekind=esmf_typekind_r8, &
1058  staggerloc=esmf_staggerloc_edge1, &
1059  name="target_grid_longitude_w", &
1060  rc=error)
1061  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1062  call error_handler("IN FieldCreate", error)
1063 
1064  print*,"- CALL FieldCreate FOR TARGET GRID TERRAIN."
1065  terrain_target_grid = esmf_fieldcreate(target_grid, &
1066  typekind=esmf_typekind_r8, &
1067  staggerloc=esmf_staggerloc_center, &
1068  name="target_grid_terrain", &
1069  rc=error)
1070  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1071  call error_handler("IN FieldCreate", error)
1072 
1073  if (localpet == 0) then
1074  allocate(landmask_one_tile(i_target,j_target))
1075  allocate(seamask_one_tile(i_target,j_target))
1076  allocate(latitude_one_tile(i_target,j_target))
1077  allocate(latitude_s_one_tile(i_target,jp1_target))
1078  allocate(latitude_w_one_tile(ip1_target,j_target))
1079  allocate(longitude_one_tile(i_target,j_target))
1080  allocate(longitude_s_one_tile(i_target,jp1_target))
1081  allocate(longitude_w_one_tile(ip1_target,j_target))
1082  allocate(terrain_one_tile(i_target,j_target))
1083  else
1084  allocate(landmask_one_tile(0,0))
1085  allocate(seamask_one_tile(0,0))
1086  allocate(longitude_one_tile(0,0))
1087  allocate(longitude_s_one_tile(0,0))
1088  allocate(longitude_w_one_tile(0,0))
1089  allocate(latitude_one_tile(0,0))
1090  allocate(latitude_s_one_tile(0,0))
1091  allocate(latitude_w_one_tile(0,0))
1092  allocate(terrain_one_tile(0,0))
1093  endif
1094 
1095  do tile = 1, num_tiles_target_grid
1096  if (localpet == 0) then
1097  the_file = trim(orog_dir_target_grid) // trim(orog_files_target_grid(tile))
1098  call get_model_mask_terrain(trim(the_file), i_target, j_target, landmask_one_tile, &
1099  terrain_one_tile)
1100  seamask_one_tile = 0
1101  where(landmask_one_tile == 0) seamask_one_tile = 1
1102  call get_model_latlons(mosaic_file_target_grid, orog_dir_target_grid, num_tiles_target_grid, tile, &
1103  i_target, j_target, ip1_target, jp1_target, latitude_one_tile, &
1104  latitude_s_one_tile, latitude_w_one_tile, longitude_one_tile, &
1105  longitude_s_one_tile, longitude_w_one_tile)
1106  endif
1107  print*,"- CALL FieldScatter FOR TARGET GRID LANDMASK. TILE IS: ", tile
1108  call esmf_fieldscatter(landmask_target_grid, landmask_one_tile, rootpet=0, tile=tile, rc=error)
1109  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1110  call error_handler("IN FieldScatter", error)
1111  print*,"- CALL FieldScatter FOR TARGET GRID SEAMASK. TILE IS: ", tile
1112  call esmf_fieldscatter(seamask_target_grid, seamask_one_tile, rootpet=0, tile=tile, rc=error)
1113  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1114  call error_handler("IN FieldScatter", error)
1115  print*,"- CALL FieldScatter FOR TARGET GRID LONGITUDE. TILE IS: ", tile
1116  call esmf_fieldscatter(longitude_target_grid, longitude_one_tile, rootpet=0, tile=tile, rc=error)
1117  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1118  call error_handler("IN FieldScatter", error)
1119  print*,"- CALL FieldScatter FOR TARGET GRID LONGITUDE_S. TILE IS: ", tile
1120  call esmf_fieldscatter(longitude_s_target_grid, longitude_s_one_tile, rootpet=0, tile=tile, rc=error)
1121  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1122  call error_handler("IN FieldScatter", error)
1123  print*,"- CALL FieldScatter FOR TARGET GRID LONGITUDE_W. TILE IS: ", tile
1124  call esmf_fieldscatter(longitude_w_target_grid, longitude_w_one_tile, rootpet=0, tile=tile, rc=error)
1125  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1126  call error_handler("IN FieldScatter", error)
1127  print*,"- CALL FieldScatter FOR TARGET GRID LATITUDE. TILE IS: ", tile
1128  call esmf_fieldscatter(latitude_target_grid, latitude_one_tile, rootpet=0, tile=tile, rc=error)
1129  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1130  call error_handler("IN FieldScatter", error)
1131  print*,"- CALL FieldScatter FOR TARGET GRID LATITUDE_S. TILE IS: ", tile
1132  call esmf_fieldscatter(latitude_s_target_grid, latitude_s_one_tile, rootpet=0, tile=tile, rc=error)
1133  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1134  call error_handler("IN FieldScatter", error)
1135  print*,"- CALL FieldScatter FOR TARGET GRID LATITUDE_W. TILE IS: ", tile
1136  call esmf_fieldscatter(latitude_w_target_grid, latitude_w_one_tile, rootpet=0, tile=tile, rc=error)
1137  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1138  call error_handler("IN FieldScatter", error)
1139  print*,"- CALL FieldScatter FOR TARGET GRID TERRAIN. TILE IS: ", tile
1140  call esmf_fieldscatter(terrain_target_grid, terrain_one_tile, rootpet=0, tile=tile, rc=error)
1141  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1142  call error_handler("IN FieldScatter", error)
1143  enddo
1144 
1145  deallocate(landmask_one_tile)
1146  deallocate(seamask_one_tile)
1147  deallocate(longitude_one_tile)
1148  deallocate(longitude_s_one_tile)
1149  deallocate(longitude_w_one_tile)
1150  deallocate(latitude_one_tile)
1151  deallocate(latitude_s_one_tile)
1152  deallocate(latitude_w_one_tile)
1153  deallocate(terrain_one_tile)
1154 
1155  end subroutine define_target_grid
1156 
1175  subroutine get_model_latlons(mosaic_file, orog_dir, num_tiles, tile, &
1176  i_tile, j_tile, ip1_tile, jp1_tile, &
1177  latitude, latitude_s, latitude_w, &
1178  longitude, longitude_s, longitude_w)
1179 
1180  use netcdf
1181 
1182  implicit none
1183 
1184  character(len=*), intent(in) :: mosaic_file, orog_dir
1185 
1186  integer, intent(in) :: num_tiles, tile
1187  integer, intent(in) :: i_tile, j_tile
1188  integer, intent(in) :: ip1_tile, jp1_tile
1189 
1190  real(esmf_kind_r8), intent(out) :: latitude(i_tile, j_tile)
1191  real(esmf_kind_r8), intent(out) :: latitude_s(i_tile, jp1_tile)
1192  real(esmf_kind_r8), intent(out) :: latitude_w(ip1_tile, j_tile)
1193  real(esmf_kind_r8), intent(out) :: longitude(i_tile, j_tile)
1194  real(esmf_kind_r8), intent(out) :: longitude_s(i_tile, jp1_tile)
1195  real(esmf_kind_r8), intent(out) :: longitude_w(ip1_tile, j_tile)
1196 
1197  character(len=50) :: grid_files(num_tiles)
1198  character(len=255) :: grid_file
1199 
1200  integer :: error, id_var, ncid
1201  integer :: id_dim, nxp, nyp, i, j, ii, jj
1202 
1203  real(esmf_kind_r8), allocatable :: tmpvar(:,:)
1204 
1205  print*,"- READ MODEL GRID FILE"
1206 
1207  print*,'- OPEN MOSAIC FILE: ', trim(mosaic_file)
1208  error=nf90_open(trim(mosaic_file), nf90_nowrite, ncid)
1209  call netcdf_err(error, 'opening mosaic file')
1210 
1211  print*,"- READ GRID FILE NAMES"
1212  error=nf90_inq_varid(ncid, 'gridfiles', id_var)
1213  call netcdf_err(error, 'reading gridfiles id')
1214  error=nf90_get_var(ncid, id_var, grid_files)
1215  call netcdf_err(error, 'reading gridfiles')
1216 
1217  error = nf90_close(ncid)
1218 
1219  grid_file = trim(orog_dir) // trim(grid_files(tile))
1220 
1221  print*,'- OPEN GRID FILE: ', trim(grid_file)
1222  error=nf90_open(trim(grid_file), nf90_nowrite, ncid)
1223  call netcdf_err(error, 'opening grid file')
1224 
1225  print*,'- READ NXP ID'
1226  error=nf90_inq_dimid(ncid, 'nxp', id_dim)
1227  call netcdf_err(error, 'reading nxp id')
1228 
1229  print*,'- READ NXP'
1230  error=nf90_inquire_dimension(ncid,id_dim,len=nxp)
1231  call netcdf_err(error, 'reading nxp')
1232 
1233  print*,'- READ NYP ID'
1234  error=nf90_inq_dimid(ncid, 'nyp', id_dim)
1235  call netcdf_err(error, 'reading nyp id')
1236 
1237  print*,'- READ NYP'
1238  error=nf90_inquire_dimension(ncid,id_dim,len=nyp)
1239  call netcdf_err(error, 'reading nyp')
1240 
1241  if ((nxp/2 /= i_tile) .or. (nyp/2 /= j_tile)) then
1242  call error_handler("DIMENSION MISMATCH IN GRID FILE.", 1)
1243  endif
1244 
1245  allocate(tmpvar(nxp,nyp))
1246 
1247  print*,'- READ LONGITUDE ID'
1248  error=nf90_inq_varid(ncid, 'x', id_var)
1249  call netcdf_err(error, 'reading longitude id')
1250 
1251  print*,'- READ LONGITUDE'
1252  error=nf90_get_var(ncid, id_var, tmpvar)
1253  call netcdf_err(error, 'reading longitude')
1254 
1255  do j = 1, j_tile
1256  do i = 1, i_tile
1257  ii = 2*i
1258  jj = 2*j
1259  longitude(i,j) = tmpvar(ii,jj)
1260  enddo
1261  enddo
1262 
1263  do j = 1, jp1_tile
1264  do i = 1, i_tile
1265  ii = 2*i
1266  jj = (2*j) - 1
1267  longitude_s(i,j) = tmpvar(ii,jj)
1268  enddo
1269  enddo
1270 
1271  do j = 1, j_tile
1272  do i = 1, ip1_tile
1273  ii = (2*i) - 1
1274  jj = 2*j
1275  longitude_w(i,j) = tmpvar(ii,jj)
1276  enddo
1277  enddo
1278 
1279  print*,'- READ LATITUDE ID'
1280  error=nf90_inq_varid(ncid, 'y', id_var)
1281  call netcdf_err(error, 'reading latitude id')
1282 
1283  print*,'- READ LATIITUDE'
1284  error=nf90_get_var(ncid, id_var, tmpvar)
1285  call netcdf_err(error, 'reading latitude')
1286 
1287  do j = 1, j_tile
1288  do i = 1, i_tile
1289  ii = 2*i
1290  jj = 2*j
1291  latitude(i,j) = tmpvar(ii,jj)
1292  enddo
1293  enddo
1294 
1295  do j = 1, jp1_tile
1296  do i = 1, i_tile
1297  ii = 2*i
1298  jj = (2*j) - 1
1299  latitude_s(i,j) = tmpvar(ii,jj)
1300  enddo
1301  enddo
1302 
1303  do j = 1, j_tile
1304  do i = 1, ip1_tile
1305  ii = (2*i) - 1
1306  jj = 2*j
1307  latitude_w(i,j) = tmpvar(ii,jj)
1308  enddo
1309  enddo
1310 
1311  deallocate(tmpvar)
1312 
1313  error = nf90_close(ncid)
1314 
1315  end subroutine get_model_latlons
1316 
1326  subroutine get_model_mask_terrain(orog_file, idim, jdim, mask, terrain)
1327 
1328  use netcdf
1329 
1330  implicit none
1331 
1332  character(len=*), intent(in) :: orog_file
1333 
1334  integer, intent(in) :: idim, jdim
1335  integer(esmf_kind_i8), intent(out) :: mask(idim,jdim)
1336 
1337  real(esmf_kind_i8), intent(out) :: terrain(idim,jdim)
1338 
1339  integer :: error, lat, lon
1340  integer :: ncid, id_dim, id_var
1341 
1342  real(kind=4), allocatable :: dummy(:,:)
1343 
1344  print*,"- READ MODEL LAND MASK FILE"
1345 
1346  print*,'- OPEN LAND MASK FILE: ', orog_file
1347  error=nf90_open(orog_file,nf90_nowrite,ncid)
1348  call netcdf_err(error, 'opening land mask file')
1349 
1350  print*,"- READ I-DIMENSION"
1351  error=nf90_inq_dimid(ncid, 'lon', id_dim)
1352  call netcdf_err(error, 'reading idim id')
1353  error=nf90_inquire_dimension(ncid,id_dim,len=lon)
1354  call netcdf_err(error, 'reading idim')
1355 
1356  print*,"- READ J-DIMENSION"
1357  error=nf90_inq_dimid(ncid, 'lat', id_dim)
1358  call netcdf_err(error, 'reading jdim id')
1359  error=nf90_inquire_dimension(ncid,id_dim,len=lat)
1360  call netcdf_err(error, 'reading jdim')
1361 
1362  print*,"- I/J DIMENSIONS: ", lon, lat
1363 
1364  if ((lon /= idim) .or. (lat /= jdim)) then
1365  call error_handler("MISMATCH IN DIMENSIONS.", 1)
1366  endif
1367 
1368  allocate(dummy(idim,jdim))
1369 
1370  print*,"- READ LAND MASK"
1371  error=nf90_inq_varid(ncid, 'slmsk', id_var)
1372  call netcdf_err(error, 'reading slmsk id')
1373  error=nf90_get_var(ncid, id_var, dummy)
1374  call netcdf_err(error, 'reading slmsk')
1375  mask = nint(dummy)
1376 
1377  print*,"- READ RAW OROGRAPHY."
1378  error=nf90_inq_varid(ncid, 'orog_raw', id_var)
1379  call netcdf_err(error, 'reading orog_raw id')
1380  error=nf90_get_var(ncid, id_var, dummy)
1381  call netcdf_err(error, 'reading orog_raw')
1382  terrain = dummy
1383 
1384  error = nf90_close(ncid)
1385 
1386  deallocate (dummy)
1387 
1388  end subroutine get_model_mask_terrain
1389 
1394 
1395  implicit none
1396 
1397  integer :: rc
1398 
1399  print*,"- DESTROY MODEL DATA."
1400 
1401  call esmf_fielddestroy(latitude_input_grid,rc=rc)
1402  call esmf_fielddestroy(longitude_input_grid,rc=rc)
1403  if (esmf_fieldiscreated(latitude_s_input_grid)) then
1404  call esmf_fielddestroy(latitude_s_input_grid, rc=rc)
1405  endif
1406  if (esmf_fieldiscreated(latitude_w_input_grid)) then
1407  call esmf_fielddestroy(latitude_w_input_grid, rc=rc)
1408  endif
1409  if (esmf_fieldiscreated(longitude_s_input_grid)) then
1410  call esmf_fielddestroy(longitude_s_input_grid, rc=rc)
1411  endif
1412  if (esmf_fieldiscreated(longitude_w_input_grid)) then
1413  call esmf_fielddestroy(longitude_w_input_grid, rc=rc)
1414  endif
1415  call esmf_fielddestroy(landmask_target_grid, rc=rc)
1416  call esmf_fielddestroy(latitude_target_grid, rc=rc)
1417  if (esmf_fieldiscreated(latitude_s_target_grid)) then
1418  call esmf_fielddestroy(latitude_s_target_grid, rc=rc)
1419  endif
1420  if (esmf_fieldiscreated(latitude_w_target_grid)) then
1421  call esmf_fielddestroy(latitude_w_target_grid, rc=rc)
1422  endif
1423  call esmf_fielddestroy(longitude_target_grid, rc=rc)
1424  if (esmf_fieldiscreated(longitude_s_target_grid)) then
1425  call esmf_fielddestroy(longitude_s_target_grid, rc=rc)
1426  endif
1427  if (esmf_fieldiscreated(longitude_w_target_grid)) then
1428  call esmf_fielddestroy(longitude_w_target_grid, rc=rc)
1429  endif
1430  call esmf_fielddestroy(seamask_target_grid, rc=rc)
1431  call esmf_fielddestroy(terrain_target_grid, rc=rc)
1432  call esmf_griddestroy(input_grid, rc=rc)
1433  call esmf_griddestroy(target_grid, rc=rc)
1434 
1435  end subroutine cleanup_input_target_grid_data
1436 
1448  subroutine gdt_to_gds(igdtnum, igdstmpl, igdtlen, kgds, ni, nj, res)
1449 
1450  implicit none
1451 
1452  integer, intent(in ) :: igdtnum, igdtlen, igdstmpl(igdtlen)
1453  integer, intent( out) :: kgds(200), ni, nj
1454  integer :: iscale
1455 
1456  real, intent( out) :: res
1457 
1458  kgds=0
1459 
1460  if (igdtnum.eq.32769) then ! rot lat/lon b grid
1461 
1462  iscale=igdstmpl(10)*igdstmpl(11)
1463  if (iscale == 0) iscale = 1e6
1464  kgds(1)=205 ! oct 6, rotated lat/lon for Non-E
1465  ! Stagger grid
1466  kgds(2)=igdstmpl(8) ! octs 7-8, Ni
1467  ni = kgds(2)
1468  kgds(3)=igdstmpl(9) ! octs 9-10, Nj
1469  nj = kgds(3)
1470  kgds(4)=nint(float(igdstmpl(12))/float(iscale)*1000.) ! octs 11-13, Lat of
1471  ! 1st grid point
1472  kgds(5)=nint(float(igdstmpl(13))/float(iscale)*1000.) ! octs 14-16, Lon of
1473  ! 1st grid point
1474 
1475  kgds(6)=0 ! oct 17, resolution and component flags
1476  if (igdstmpl(1)==2 ) kgds(6)=64
1477  if ( btest(igdstmpl(14),4).OR.btest(igdstmpl(14),5) ) kgds(6)=kgds(6)+128
1478  if ( btest(igdstmpl(14),3) ) kgds(6)=kgds(6)+8
1479 
1480  kgds(7)=nint(float(igdstmpl(15))/float(iscale)*1000.) ! octs 18-20,
1481  ! Lat of cent of rotation
1482  kgds(8)=nint(float(igdstmpl(16))/float(iscale)*1000.) ! octs 21-23,
1483  ! Lon of cent of rotation
1484  kgds(9)=nint(float(igdstmpl(17))/float(iscale)*1000.) ! octs 24-25,
1485  ! Di
1486  kgds(10)=nint(float(igdstmpl(18))/float(iscale)*1000.) ! octs 26-27,
1487  ! Dj
1488 
1489  kgds(11) = 0 ! oct 28, scan mode
1490  if (btest(igdstmpl(19),7)) kgds(11) = 128
1491  if (btest(igdstmpl(19),6)) kgds(11) = kgds(11) + 64
1492  if (btest(igdstmpl(19),5)) kgds(11) = kgds(11) + 32
1493 
1494  kgds(12)=nint(float(igdstmpl(20))/float(iscale)*1000.) ! octs 29-31, Lat of
1495  ! last grid point
1496  kgds(13)=nint(float(igdstmpl(21))/float(iscale)*1000.) ! octs 32-34, Lon of
1497  ! last grid point
1498 
1499  kgds(19)=0 ! oct 4, # vert coordinate parameters
1500  kgds(20)=255 ! oct 5, used for thinned grids, set to 255
1501 
1502  res = ((float(kgds(9)) / 1000.0) + (float(kgds(10)) / 1000.0)) &
1503  * 0.5 * 111.0
1504 
1505  elseif(igdtnum==30) then
1506 
1507  kgds(1)=3 ! oct 6, lambert conformal
1508  kgds(2)=igdstmpl(8) ! octs 7-8, Ni
1509  ni = kgds(2)
1510  kgds(3)=igdstmpl(9) ! octs 9-10, Nj
1511  nj = kgds(3)
1512 
1513  iscale = 1e6
1514  kgds(4) = nint(float(igdstmpl(10))/1000.0)
1515  kgds(5) = nint(float(igdstmpl(11))/1000.0)
1516 
1517  kgds(6)=0 ! oct 17, resolution and component flags
1518  if (igdstmpl(1)==2 ) kgds(6)=64
1519  if ( btest(igdstmpl(12),4).OR.btest(igdstmpl(12),5) ) kgds(6)=kgds(6)+128
1520  if ( btest(igdstmpl(12),3) ) kgds(6)=kgds(6)+8
1521 
1522  kgds(7) = nint(float(igdstmpl(14))/1000.0)
1523  kgds(8) = nint(float(igdstmpl(15))/1000.0)
1524  kgds(9) = nint(float(igdstmpl(16))/1000.0)
1525  kgds(10) = 0
1526 
1527  kgds(11) = 0 ! oct 28, scan mode
1528  if (btest(igdstmpl(18),7)) kgds(11) = 128
1529  if (btest(igdstmpl(18),6)) kgds(11) = kgds(11) + 64
1530  if (btest(igdstmpl(18),5)) kgds(11) = kgds(11) + 32
1531 
1532  kgds(12) = nint(float(igdstmpl(19))/1000.0)
1533  kgds(13) = nint(float(igdstmpl(20))/1000.0)
1534  kgds(14) = -90
1535  kgds(15) = 0
1536 
1537  elseif(igdtnum==0) then ! lat/lon grid
1538 
1539  iscale=igdstmpl(10)*igdstmpl(11)
1540  if (iscale == 0) iscale = 1e6
1541  kgds(1)=0 ! oct 6, data representation type.
1542  kgds(2)=igdstmpl(8) ! octs 7-8, Ni
1543  ni = kgds(2)
1544  kgds(3)=igdstmpl(9) ! octs 9-10, Nj
1545  nj = kgds(3)
1546  kgds(4)=nint(float(igdstmpl(12))/float(iscale)*1000.) ! octs 11-13, Lat of 1st grid point
1547  kgds(5)=nint(float(igdstmpl(13))/float(iscale)*1000.) ! octs 14-16, Lon of 1st grid point
1548 
1549  kgds(6)=0 ! oct 17, resolution and component flags
1550  if (igdstmpl(1)==2 ) kgds(6)=64
1551  if ( btest(igdstmpl(14),4).OR.btest(igdstmpl(14),5) ) kgds(6)=kgds(6)+128
1552  if ( btest(igdstmpl(14),3) ) kgds(6)=kgds(6)+8
1553 
1554  kgds(7)=nint(float(igdstmpl(15))/float(iscale)*1000.) ! octs 18-20, Lat of last grid point
1555  kgds(8)=nint(float(igdstmpl(16))/float(iscale)*1000.) ! octs 21-23, Lon of last grid point
1556  kgds(9)=nint(float(igdstmpl(17))/float(iscale)*1000.) ! octs 24-25, "i" resolution.
1557  kgds(10)=nint(float(igdstmpl(18))/float(iscale)*1000.) ! octs 26-27, "j" resolution.
1558 
1559  kgds(11) = 0 ! oct 28, scan mode
1560  if (btest(igdstmpl(19),7)) kgds(11) = 128
1561  if (btest(igdstmpl(19),6)) kgds(11) = kgds(11) + 64
1562  if (btest(igdstmpl(19),5)) kgds(11) = kgds(11) + 32
1563 
1564  kgds(12)=0 ! octs 29-32, reserved
1565  kgds(19)=0 ! oct 4, # vert coordinate parameters
1566  kgds(20)=255 ! oct 5, used for thinned grids, set to 255
1567 
1568  else
1569 
1570  call error_handler("UNRECOGNIZED INPUT GRID TYPE ", 1)
1571 
1572  endif
1573 
1574  end subroutine gdt_to_gds
1575 
1576  end module model_grid
subroutine get_model_mask_terrain(orog_file, idim, jdim, mask, terrain)
Read the model land mask and terrain for a single tile from the orography file.
subroutine, public cleanup_input_target_grid_data
Deallocate all esmf grid objects.
subroutine define_input_grid_grib2(npets)
Define input grid object for grib2 input data.
Definition: model_grid.F90:614
subroutine gdt_to_gds(igdtnum, igdstmpl, igdtlen, kgds, ni, nj, res)
Convert the GRIB2 grid description template to to the GRIB1 grid description section.
subroutine define_input_grid_gaussian(npets)
Define grid object for input data on global gaussian grids.
Definition: model_grid.F90:147
subroutine netcdf_err(err, string)
Error handler for netcdf.
Definition: utils.F90:34
Sets up the ESMF grid objects for the input data grid and target FV3 grid.
Definition: model_grid.F90:9
subroutine get_model_latlons(mosaic_file, orog_dir, num_tiles, tile, i_tile, j_tile, ip1_tile, jp1_tile, latitude, latitude_s, latitude_w, longitude, longitude_s, longitude_w)
Read model lat/lons for a single tile from the "grid" specificaton file.
subroutine error_handler(string, rc)
General error handler.
Definition: utils.F90:12
subroutine define_input_grid_mosaic(localpet, npets)
Define input grid for tiled data using the 'mosaic', 'grid' and orography files.
Definition: model_grid.F90:401
This module contains code to read the setup namelist file, handle the varmap file for GRIB2 data...
subroutine, public define_input_grid(localpet, npets)
Driver routine to setup the esmf grid object for the input grid.
Definition: model_grid.F90:116
subroutine, public define_target_grid(localpet, npets)
Setup the esmf grid object for the target grid.
Definition: model_grid.F90:887