chgres_cube  1.7.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 
14  implicit none
15 
16  private
17 
18  character(len=5), allocatable, public :: tiles_target_grid(:)
19 
20  character(len=50), public :: input_grid_type = "latlon"
21 
22 
23  ! Made lsoil_target non-parameter to allow for RAP land surface initiation
24  integer, public :: lsoil_target = 4 ! # soil layers
25 
26  integer, public :: i_input
27 
29  integer, public :: j_input
30 
32  integer, public :: ip1_input
33 
34  integer, public :: jp1_input
35 
36  integer, public :: i_target
37 
39  integer, public :: j_target
40 
42  integer, public :: ip1_target
43 
44  integer, public :: jp1_target
45 
46  integer, public :: num_tiles_input_grid
47 
48  integer, public :: num_tiles_target_grid
49 
50 
51  type(esmf_grid), public :: input_grid
52 
53  type(esmf_grid), public :: target_grid
54 
55 
56  type(esmf_field), public :: latitude_input_grid
57 
58  type(esmf_field), public :: longitude_input_grid
59 
60  type(esmf_field), public :: latitude_s_input_grid
61 
63  type(esmf_field), public :: longitude_s_input_grid
64 
66  type(esmf_field), public :: latitude_w_input_grid
67 
69  type(esmf_field), public :: longitude_w_input_grid
70 
72 
73  type(esmf_field), public :: landmask_target_grid
74 
76  type(esmf_field), public :: latitude_target_grid
77 
78  type(esmf_field), public :: latitude_s_target_grid
79 
81  type(esmf_field), public :: latitude_w_target_grid
82 
84  type(esmf_field), public :: longitude_target_grid
85 
86  type(esmf_field), public :: longitude_s_target_grid
87 
89  type(esmf_field), public :: longitude_w_target_grid
90 
92  type(esmf_field), public :: seamask_target_grid
93 
95  type(esmf_field), public :: terrain_target_grid
96 
97 
98  public :: define_target_grid
99  public :: define_input_grid
101 
102  contains
103 
115  subroutine define_input_grid(localpet, npets)
116 
117  use program_setup, only : input_type
118 
119  implicit none
120 
121  integer, intent(in) :: localpet, npets
122 
123  if (trim(input_type) == "gaussian_nemsio" .or. &
124  trim(input_type) == "gfs_gaussian_nemsio" .or. &
125  trim(input_type) == "gfs_sigio" .or. &
126  trim(input_type) == "gaussian_netcdf") then
127  call define_input_grid_gaussian(localpet, npets)
128  elseif (trim(input_type) == "grib2") then
129  call define_input_grid_grib2(localpet,npets)
130  else
131  call define_input_grid_mosaic(localpet, npets)
132  endif
133 
134  end subroutine define_input_grid
135 
147  subroutine define_input_grid_gaussian(localpet, 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) :: localpet, 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 
615  subroutine define_input_grid_grib2(localpet,npets)
616 
617  use grib_mod
618  use gdswzd_mod
619  use program_setup, only : grib2_file_input_grid, data_dir_input_grid
620 
621  implicit none
622 
623  integer, intent(in) :: localpet, npets
624 
625  character(len=500) :: the_file
626 
627  integer :: i, j, k, jdisc, jgdtn, jpdtn, lugb, lugi
628  integer :: jids(200), jgdt(200), jpdt(200), rc
629  integer :: kgds(200), nret, clb(2), cub(2)
630 
631  logical :: unpack
632 
633  real :: res
634  real, allocatable :: rlon(:,:),rlat(:,:),xpts(:,:),ypts(:,:)
635  real, allocatable :: rlon_corner(:,:),rlat_corner(:,:)
636  real, allocatable :: rlon_diff(:,:),rlat_diff(:,:)
637  real, allocatable :: xpts_corner(:,:),ypts_corner(:,:)
638  real(esmf_kind_r8), allocatable :: latitude(:,:)
639  real(esmf_kind_r8), allocatable :: longitude(:,:)
640  real(esmf_kind_r8), allocatable :: latitude_corner(:,:)
641  real(esmf_kind_r8), allocatable :: longitude_corner(:,:)
642  real(esmf_kind_r8), pointer :: lat_src_ptr(:,:)
643  real(esmf_kind_r8), pointer :: lat_corner_src_ptr(:,:)
644  real(esmf_kind_r8), pointer :: lon_src_ptr(:,:)
645  real(esmf_kind_r8), pointer :: lon_corner_src_ptr(:,:)
646 
647  type(esmf_polekind_flag) :: polekindflag(2)
648 
649  type(gribfield) :: gfld
650 
651  the_file = trim(data_dir_input_grid) // "/" // grib2_file_input_grid
652 
653  lugb=12
654 
655  print*,"- OPEN AND READ INPUT DATA GRIB2 FILE: ", trim(the_file)
656  call baopenr(lugb,the_file,rc)
657  if (rc /= 0) call error_handler("OPENING FILE", rc)
658 
659 ! Read the first record and get the grid definition template.
660 
661  j = 0 ! Search at beginning of file
662  lugi = 0 ! No grib index file
663  jdisc = -1 ! Search for any discipline
664  jpdtn = -1 ! Search for any product definition template number
665  jgdtn = -1 ! Search for any grid definition template number
666  jids = -9999 ! Array of values in identification section, set to wildcard.
667  jgdt = -9999 ! Array of values in grid definition template, set to wildcard.
668  jpdt = -9999 ! Array of values in product definition template, set to wildcard.
669  unpack = .false. ! unpack data
670 
671  call getgb2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
672  unpack, k, gfld, rc)
673  if (rc /= 0) call error_handler("DEGRIBBING INPUT FILE.", rc)
674 
675  call baclose(lugb,rc)
676 
677  if (gfld%igdtnum == 0) then
678  print*,"- INPUT DATA ON LAT/LON GRID."
679  input_grid_type = 'latlon'
680  elseif (gfld%igdtnum == 30) then
681  print*,"- INPUT DATA ON LAMBERT CONFORMAL GRID."
682  input_grid_type = 'lambert'
683  elseif (gfld%igdtnum == 32769) then
684  print*,"- INPUT DATA ON ROTATED LAT/LON GRID."
685  input_grid_type = 'rotated_latlon'
686  else
687  call error_handler("INPUT GRID TEMPLATE NOT SUPPORTED.", 2)
688  endif
689 
690  kgds = 0
691  call gdt_to_gds(gfld%igdtnum, gfld%igdtmpl, gfld%igdtlen, kgds, i_input, j_input, res)
692 
693  ip1_input = i_input + 1
694  jp1_input = j_input + 1
695 
696  allocate(rlat(i_input,j_input))
697  allocate(rlon(i_input,j_input))
698  allocate(rlat_diff(i_input,j_input))
699  allocate(rlon_diff(i_input,j_input))
700  allocate(xpts(i_input,j_input))
701  allocate(ypts(i_input,j_input))
702  allocate(rlat_corner(ip1_input,jp1_input))
703  allocate(rlon_corner(ip1_input,jp1_input))
704  allocate(xpts_corner(ip1_input,jp1_input))
705  allocate(ypts_corner(ip1_input,jp1_input))
706 
707  do j = 1, j_input
708  do i = 1, i_input
709  xpts(i,j) = float(i)
710  ypts(i,j) = float(j)
711  enddo
712  enddo
713 
714  print*,"- COMPUTE GRID CELL CENTER COORDINATES."
715  call gdswzd(kgds,1,(i_input*j_input),-9999.,xpts,ypts,rlon,rlat,nret)
716 
717  if (nret /= (i_input*j_input)) then
718  call error_handler("GDSWZD RETURNED WRONG NUMBER OF POINTS.", 2)
719  endif
720 
721  deallocate(xpts, ypts)
722 
723  do j = 1, jp1_input
724  do i = 1, ip1_input
725  xpts_corner(i,j) = float(i) - 0.5
726  ypts_corner(i,j) = float(j) - 0.5
727  enddo
728  enddo
729 
730  print*,"- COMPUTE GRID CELL CORNER COORDINATES."
731  call gdswzd(kgds,1,(ip1_input*jp1_input),-9999.,xpts_corner,ypts_corner,rlon_corner,rlat_corner,nret)
732 
733  if (nret /= (ip1_input*jp1_input)) then
734  call error_handler("GDSWZD RETURNED WRONG NUMBER OF POINTS.", 2)
735  endif
736 
737  deallocate(xpts_corner, ypts_corner)
738 
739  if (gfld%igdtnum == 0) then ! gfs lat/lon data
740 
741  print*,"- CALL GridCreate1PeriDim FOR INPUT GRID."
742 
743  polekindflag(1:2) = esmf_polekind_monopole
744 
745  input_grid = esmf_gridcreate1peridim(minindex=(/1,1/), &
746  maxindex=(/i_input,j_input/), &
747  polekindflag=polekindflag, &
748  periodicdim=1, &
749  poledim=2, &
750  coordsys=esmf_coordsys_sph_deg, &
751  regdecomp=(/1,npets/), &
752  indexflag=esmf_index_global, rc=rc)
753 
754  else
755 
756  print*,"- CALL GridCreateNoPeriDim FOR INPUT GRID."
757 
758  input_grid = esmf_gridcreatenoperidim(maxindex=(/i_input,j_input/), &
759  indexflag=esmf_index_global, &
760  rc=rc)
761 
762  endif
763 
764  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
765  call error_handler("IN GridCreate1PeriDim", rc)
766 
767  print*,"- CALL FieldCreate FOR INPUT GRID LATITUDE."
768  latitude_input_grid = esmf_fieldcreate(input_grid, &
769  typekind=esmf_typekind_r8, &
770  staggerloc=esmf_staggerloc_center, &
771  name="input_grid_latitude", rc=rc)
772  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
773  call error_handler("IN FieldCreate", rc)
774 
775  print*,"- CALL FieldCreate FOR INPUT GRID LONGITUDE."
776  longitude_input_grid = esmf_fieldcreate(input_grid, &
777  typekind=esmf_typekind_r8, &
778  staggerloc=esmf_staggerloc_center, &
779  name="input_grid_longitude", rc=rc)
780  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
781  call error_handler("IN FieldCreate", rc)
782 
783  allocate(latitude(i_input,j_input))
784  allocate(longitude(i_input,j_input))
785 
786  latitude = rlat
787  longitude = rlon
788 
789  deallocate (rlat, rlon)
790 
791  print*,"- CALL FieldScatter FOR INPUT GRID LONGITUDE."
792  call esmf_fieldscatter(longitude_input_grid, longitude, rootpet=0, rc=rc)
793  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
794  call error_handler("IN FieldScatter", rc)
795 
796  print*,"- CALL FieldScatter FOR INPUT GRID LATITUDE."
797  call esmf_fieldscatter(latitude_input_grid, latitude, rootpet=0, rc=rc)
798  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
799  call error_handler("IN FieldScatter", rc)
800 
801  print*,"- CALL GridAddCoord FOR INPUT GRID."
802  call esmf_gridaddcoord(input_grid, &
803  staggerloc=esmf_staggerloc_center, rc=rc)
804  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
805  call error_handler("IN GridAddCoord", rc)
806 
807  print*,"- CALL GridGetCoord FOR INPUT GRID X-COORD."
808  nullify(lon_src_ptr)
809  call esmf_gridgetcoord(input_grid, &
810  staggerloc=esmf_staggerloc_center, &
811  coorddim=1, &
812  farrayptr=lon_src_ptr, rc=rc)
813  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
814  call error_handler("IN GridGetCoord", rc)
815 
816  print*,"- CALL GridGetCoord FOR INPUT GRID Y-COORD."
817  nullify(lat_src_ptr)
818  call esmf_gridgetcoord(input_grid, &
819  staggerloc=esmf_staggerloc_center, &
820  coorddim=2, &
821  computationallbound=clb, &
822  computationalubound=cub, &
823  farrayptr=lat_src_ptr, rc=rc)
824  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
825  call error_handler("IN GridGetCoord", rc)
826 
827  do j = clb(2), cub(2)
828  do i = clb(1), cub(1)
829  lon_src_ptr(i,j) = longitude(i,j)
830  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
831  lat_src_ptr(i,j) = latitude(i,j)
832  enddo
833  enddo
834 
835  deallocate(latitude, longitude)
836 
837  allocate(latitude_corner(ip1_input,jp1_input))
838  allocate(longitude_corner(ip1_input,jp1_input))
839 
840  latitude_corner = rlat_corner
841  longitude_corner = rlon_corner
842 
843  deallocate (rlat_corner, rlon_corner)
844 
845  print*,"- CALL GridAddCoord FOR INPUT GRID."
846  call esmf_gridaddcoord(input_grid, &
847  staggerloc=esmf_staggerloc_corner, rc=rc)
848  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
849  call error_handler("IN GridAddCoord", rc)
850 
851  print*,"- CALL GridGetCoord FOR INPUT GRID X-COORD."
852  nullify(lon_corner_src_ptr)
853  call esmf_gridgetcoord(input_grid, &
854  staggerloc=esmf_staggerloc_corner, &
855  coorddim=1, &
856  farrayptr=lon_corner_src_ptr, rc=rc)
857  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
858  call error_handler("IN GridGetCoord", rc)
859 
860  print*,"- CALL GridGetCoord FOR INPUT GRID Y-COORD."
861  nullify(lat_corner_src_ptr)
862  call esmf_gridgetcoord(input_grid, &
863  staggerloc=esmf_staggerloc_corner, &
864  coorddim=2, &
865  computationallbound=clb, &
866  computationalubound=cub, &
867  farrayptr=lat_corner_src_ptr, rc=rc)
868  if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
869  call error_handler("IN GridGetCoord", rc)
870 
871  do j = clb(2), cub(2)
872  do i = clb(1), cub(1)
873  lon_corner_src_ptr(i,j) = longitude_corner(i,j)
874  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
875  lat_corner_src_ptr(i,j) = latitude_corner(i,j)
876  enddo
877  enddo
878 
879  deallocate(latitude_corner, longitude_corner)
880 
881  end subroutine define_input_grid_grib2
882 
888  subroutine define_target_grid(localpet, npets)
889 
890  use netcdf
891  use program_setup, only : mosaic_file_target_grid, &
892  orog_dir_target_grid, &
893  orog_files_target_grid, &
894  nsoill_out
895 
896  implicit none
897 
898  integer, intent(in) :: localpet, npets
899 
900  character(len=500) :: the_file
901 
902  integer :: error, ncid, extra
903  integer :: id_tiles
904  integer :: id_dim, id_grid_tiles
905  integer :: tile
906  integer, allocatable :: decomptile(:,:)
907  integer(esmf_kind_i8), allocatable :: landmask_one_tile(:,:)
908  integer(esmf_kind_i8), allocatable :: seamask_one_tile(:,:)
909 
910  real(esmf_kind_r8), allocatable :: latitude_one_tile(:,:)
911  real(esmf_kind_r8), allocatable :: latitude_s_one_tile(:,:)
912  real(esmf_kind_r8), allocatable :: latitude_w_one_tile(:,:)
913  real(esmf_kind_r8), allocatable :: longitude_one_tile(:,:)
914  real(esmf_kind_r8), allocatable :: longitude_s_one_tile(:,:)
915  real(esmf_kind_r8), allocatable :: longitude_w_one_tile(:,:)
916  real(esmf_kind_r8), allocatable :: terrain_one_tile(:,:)
917 
918  lsoil_target = nsoill_out
919 
920  print*,'- OPEN TARGET GRID MOSAIC FILE: ',trim(mosaic_file_target_grid)
921  error=nf90_open(trim(mosaic_file_target_grid),nf90_nowrite,ncid)
922  call netcdf_err(error, 'opening grid mosaic file')
923 
924  print*,"- READ NUMBER OF TILES"
925  error=nf90_inq_dimid(ncid, 'ntiles', id_tiles)
926  call netcdf_err(error, 'reading ntile id')
927  error=nf90_inquire_dimension(ncid,id_tiles,len=num_tiles_target_grid)
928  call netcdf_err(error, 'reading ntiles')
929  error=nf90_inq_varid(ncid, 'gridtiles', id_grid_tiles)
930  call netcdf_err(error, 'reading gridtiles id')
931  allocate(tiles_target_grid(num_tiles_target_grid))
932  tiles_target_grid="NULL"
933  print*,"- READ TILE NAMES"
934  error=nf90_get_var(ncid, id_grid_tiles, tiles_target_grid)
935  call netcdf_err(error, 'reading gridtiles')
936 
937  error = nf90_close(ncid)
938 
939  print*,'- NUMBER OF TILES, TARGET MODEL GRID IS ', num_tiles_target_grid
940 
941  if (mod(npets,num_tiles_target_grid) /= 0) then
942  call error_handler("MUST RUN WITH TASK COUNT THAT IS A MULTIPLE OF # OF TILES.", 1)
943  endif
944 
945 !-----------------------------------------------------------------------
946 ! Get the model grid specs and land mask from the orography files.
947 !-----------------------------------------------------------------------
948 
949  the_file = trim(orog_dir_target_grid) // trim(orog_files_target_grid(1))
950 
951  print*,'- OPEN FIRST TARGET GRID OROGRAPHY FILE: ',trim(the_file)
952  error=nf90_open(trim(the_file),nf90_nowrite,ncid)
953  call netcdf_err(error, 'opening orography file')
954  print*,"- READ GRID DIMENSIONS"
955  error=nf90_inq_dimid(ncid, 'lon', id_dim)
956  call netcdf_err(error, 'reading lon id')
957  error=nf90_inquire_dimension(ncid,id_dim,len=i_target)
958  call netcdf_err(error, 'reading lon')
959  error=nf90_inq_dimid(ncid, 'lat', id_dim)
960  call netcdf_err(error, 'reading lat id')
961  error=nf90_inquire_dimension(ncid,id_dim,len=j_target)
962  call netcdf_err(error, 'reading lat')
963  error = nf90_close(ncid)
964 
965  print*,"- I/J DIMENSIONS OF THE TARGET GRID TILES ", i_target, j_target
966 
967  ip1_target = i_target + 1
968  jp1_target = j_target + 1
969 
970 !-----------------------------------------------------------------------
971 ! Create ESMF grid object for the model grid.
972 !-----------------------------------------------------------------------
973 
974  extra = npets / num_tiles_target_grid
975 
976  allocate(decomptile(2,num_tiles_target_grid))
977 
978  do tile = 1, num_tiles_target_grid
979  decomptile(:,tile)=(/1,extra/)
980  enddo
981 
982  print*,"- CALL GridCreateMosaic FOR TARGET GRID"
983  target_grid = esmf_gridcreatemosaic(filename=trim(mosaic_file_target_grid), &
984  regdecompptile=decomptile, &
985  staggerloclist=(/esmf_staggerloc_center, esmf_staggerloc_corner, &
986  esmf_staggerloc_edge1, esmf_staggerloc_edge2/), &
987  indexflag=esmf_index_global, &
988  tilefilepath=trim(orog_dir_target_grid), rc=error)
989  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
990  call error_handler("IN GridCreateMosaic", error)
991 
992 !-----------------------------------------------------------------------
993 ! Set target model landmask (1 - land, 0 - not land) and
994 ! seamask (1 - non-land, 0 -land). Read lat/lon on target grid.
995 !-----------------------------------------------------------------------
996 
997  print*,"- CALL FieldCreate FOR TARGET GRID LANDMASK."
998  landmask_target_grid = esmf_fieldcreate(target_grid, &
999  typekind=esmf_typekind_i8, &
1000  staggerloc=esmf_staggerloc_center, &
1001  name="target_grid_landmask", rc=error)
1002  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1003  call error_handler("IN FieldCreate", error)
1004 
1005  print*,"- CALL FieldCreate FOR TARGET GRID SEAMASK."
1006  seamask_target_grid = esmf_fieldcreate(target_grid, &
1007  typekind=esmf_typekind_i8, &
1008  staggerloc=esmf_staggerloc_center, &
1009  name="target_grid_seamask", rc=error)
1010  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1011  call error_handler("IN FieldCreate", error)
1012 
1013  print*,"- CALL FieldCreate FOR TARGET GRID LATITUDE."
1014  latitude_target_grid = esmf_fieldcreate(target_grid, &
1015  typekind=esmf_typekind_r8, &
1016  staggerloc=esmf_staggerloc_center, &
1017  name="target_grid_latitude", rc=error)
1018  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1019  call error_handler("IN FieldCreate", error)
1020 
1021  print*,"- CALL FieldCreate FOR TARGET GRID LATITUDE_S."
1022  latitude_s_target_grid = esmf_fieldcreate(target_grid, &
1023  typekind=esmf_typekind_r8, &
1024  staggerloc=esmf_staggerloc_edge2, &
1025  name="target_grid_latitude_s", rc=error)
1026  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1027  call error_handler("IN FieldCreate", error)
1028 
1029  print*,"- CALL FieldCreate FOR TARGET GRID LATITUDE_W."
1030  latitude_w_target_grid = esmf_fieldcreate(target_grid, &
1031  typekind=esmf_typekind_r8, &
1032  staggerloc=esmf_staggerloc_edge1, &
1033  name="target_grid_latitude_w", &
1034  rc=error)
1035  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1036  call error_handler("IN FieldCreate", error)
1037 
1038  print*,"- CALL FieldCreate FOR TARGET GRID LONGITUDE."
1039  longitude_target_grid = esmf_fieldcreate(target_grid, &
1040  typekind=esmf_typekind_r8, &
1041  staggerloc=esmf_staggerloc_center, &
1042  name="target_grid_longitude", &
1043  rc=error)
1044  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1045  call error_handler("IN FieldCreate", error)
1046 
1047  print*,"- CALL FieldCreate FOR TARGET GRID LONGITUDE_S."
1048  longitude_s_target_grid = esmf_fieldcreate(target_grid, &
1049  typekind=esmf_typekind_r8, &
1050  staggerloc=esmf_staggerloc_edge2, &
1051  name="target_grid_longitude_s", &
1052  rc=error)
1053  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1054  call error_handler("IN FieldCreate", error)
1055 
1056  print*,"- CALL FieldCreate FOR TARGET GRID LONGITUDE_W."
1057  longitude_w_target_grid = esmf_fieldcreate(target_grid, &
1058  typekind=esmf_typekind_r8, &
1059  staggerloc=esmf_staggerloc_edge1, &
1060  name="target_grid_longitude_w", &
1061  rc=error)
1062  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1063  call error_handler("IN FieldCreate", error)
1064 
1065  print*,"- CALL FieldCreate FOR TARGET GRID TERRAIN."
1066  terrain_target_grid = esmf_fieldcreate(target_grid, &
1067  typekind=esmf_typekind_r8, &
1068  staggerloc=esmf_staggerloc_center, &
1069  name="target_grid_terrain", &
1070  rc=error)
1071  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1072  call error_handler("IN FieldCreate", error)
1073 
1074  if (localpet == 0) then
1075  allocate(landmask_one_tile(i_target,j_target))
1076  allocate(seamask_one_tile(i_target,j_target))
1077  allocate(latitude_one_tile(i_target,j_target))
1078  allocate(latitude_s_one_tile(i_target,jp1_target))
1079  allocate(latitude_w_one_tile(ip1_target,j_target))
1080  allocate(longitude_one_tile(i_target,j_target))
1081  allocate(longitude_s_one_tile(i_target,jp1_target))
1082  allocate(longitude_w_one_tile(ip1_target,j_target))
1083  allocate(terrain_one_tile(i_target,j_target))
1084  else
1085  allocate(landmask_one_tile(0,0))
1086  allocate(seamask_one_tile(0,0))
1087  allocate(longitude_one_tile(0,0))
1088  allocate(longitude_s_one_tile(0,0))
1089  allocate(longitude_w_one_tile(0,0))
1090  allocate(latitude_one_tile(0,0))
1091  allocate(latitude_s_one_tile(0,0))
1092  allocate(latitude_w_one_tile(0,0))
1093  allocate(terrain_one_tile(0,0))
1094  endif
1095 
1096  do tile = 1, num_tiles_target_grid
1097  if (localpet == 0) then
1098  the_file = trim(orog_dir_target_grid) // trim(orog_files_target_grid(tile))
1099  call get_model_mask_terrain(trim(the_file), i_target, j_target, landmask_one_tile, &
1100  terrain_one_tile)
1101  seamask_one_tile = 0
1102  where(landmask_one_tile == 0) seamask_one_tile = 1
1103  call get_model_latlons(mosaic_file_target_grid, orog_dir_target_grid, num_tiles_target_grid, tile, &
1104  i_target, j_target, ip1_target, jp1_target, latitude_one_tile, &
1105  latitude_s_one_tile, latitude_w_one_tile, longitude_one_tile, &
1106  longitude_s_one_tile, longitude_w_one_tile)
1107  endif
1108  print*,"- CALL FieldScatter FOR TARGET GRID LANDMASK. TILE IS: ", tile
1109  call esmf_fieldscatter(landmask_target_grid, landmask_one_tile, rootpet=0, tile=tile, rc=error)
1110  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1111  call error_handler("IN FieldScatter", error)
1112  print*,"- CALL FieldScatter FOR TARGET GRID SEAMASK. TILE IS: ", tile
1113  call esmf_fieldscatter(seamask_target_grid, seamask_one_tile, rootpet=0, tile=tile, rc=error)
1114  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1115  call error_handler("IN FieldScatter", error)
1116  print*,"- CALL FieldScatter FOR TARGET GRID LONGITUDE. TILE IS: ", tile
1117  call esmf_fieldscatter(longitude_target_grid, longitude_one_tile, rootpet=0, tile=tile, rc=error)
1118  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1119  call error_handler("IN FieldScatter", error)
1120  print*,"- CALL FieldScatter FOR TARGET GRID LONGITUDE_S. TILE IS: ", tile
1121  call esmf_fieldscatter(longitude_s_target_grid, longitude_s_one_tile, rootpet=0, tile=tile, rc=error)
1122  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1123  call error_handler("IN FieldScatter", error)
1124  print*,"- CALL FieldScatter FOR TARGET GRID LONGITUDE_W. TILE IS: ", tile
1125  call esmf_fieldscatter(longitude_w_target_grid, longitude_w_one_tile, rootpet=0, tile=tile, rc=error)
1126  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1127  call error_handler("IN FieldScatter", error)
1128  print*,"- CALL FieldScatter FOR TARGET GRID LATITUDE. TILE IS: ", tile
1129  call esmf_fieldscatter(latitude_target_grid, latitude_one_tile, rootpet=0, tile=tile, rc=error)
1130  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1131  call error_handler("IN FieldScatter", error)
1132  print*,"- CALL FieldScatter FOR TARGET GRID LATITUDE_S. TILE IS: ", tile
1133  call esmf_fieldscatter(latitude_s_target_grid, latitude_s_one_tile, rootpet=0, tile=tile, rc=error)
1134  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1135  call error_handler("IN FieldScatter", error)
1136  print*,"- CALL FieldScatter FOR TARGET GRID LATITUDE_W. TILE IS: ", tile
1137  call esmf_fieldscatter(latitude_w_target_grid, latitude_w_one_tile, rootpet=0, tile=tile, rc=error)
1138  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1139  call error_handler("IN FieldScatter", error)
1140  print*,"- CALL FieldScatter FOR TARGET GRID TERRAIN. TILE IS: ", tile
1141  call esmf_fieldscatter(terrain_target_grid, terrain_one_tile, rootpet=0, tile=tile, rc=error)
1142  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
1143  call error_handler("IN FieldScatter", error)
1144  enddo
1145 
1146  deallocate(landmask_one_tile)
1147  deallocate(seamask_one_tile)
1148  deallocate(longitude_one_tile)
1149  deallocate(longitude_s_one_tile)
1150  deallocate(longitude_w_one_tile)
1151  deallocate(latitude_one_tile)
1152  deallocate(latitude_s_one_tile)
1153  deallocate(latitude_w_one_tile)
1154  deallocate(terrain_one_tile)
1155 
1156  end subroutine define_target_grid
1157 
1176  subroutine get_model_latlons(mosaic_file, orog_dir, num_tiles, tile, &
1177  i_tile, j_tile, ip1_tile, jp1_tile, &
1178  latitude, latitude_s, latitude_w, &
1179  longitude, longitude_s, longitude_w)
1180 
1181  use netcdf
1182 
1183  implicit none
1184 
1185  character(len=*), intent(in) :: mosaic_file, orog_dir
1186 
1187  integer, intent(in) :: num_tiles, tile
1188  integer, intent(in) :: i_tile, j_tile
1189  integer, intent(in) :: ip1_tile, jp1_tile
1190 
1191  real(esmf_kind_r8), intent(out) :: latitude(i_tile, j_tile)
1192  real(esmf_kind_r8), intent(out) :: latitude_s(i_tile, jp1_tile)
1193  real(esmf_kind_r8), intent(out) :: latitude_w(ip1_tile, j_tile)
1194  real(esmf_kind_r8), intent(out) :: longitude(i_tile, j_tile)
1195  real(esmf_kind_r8), intent(out) :: longitude_s(i_tile, jp1_tile)
1196  real(esmf_kind_r8), intent(out) :: longitude_w(ip1_tile, j_tile)
1197 
1198  character(len=25) :: grid_files(num_tiles)
1199  character(len=255) :: grid_file
1200 
1201  integer :: error, id_var, ncid
1202  integer :: id_dim, nxp, nyp, i, j, ii, jj
1203 
1204  real(esmf_kind_r8), allocatable :: tmpvar(:,:)
1205 
1206  print*,"- READ MODEL GRID FILE"
1207 
1208  print*,'- OPEN MOSAIC FILE: ', trim(mosaic_file)
1209  error=nf90_open(trim(mosaic_file), nf90_nowrite, ncid)
1210  call netcdf_err(error, 'opening mosaic file')
1211 
1212  print*,"- READ GRID FILE NAMES"
1213  error=nf90_inq_varid(ncid, 'gridfiles', id_var)
1214  call netcdf_err(error, 'reading gridfiles id')
1215  error=nf90_get_var(ncid, id_var, grid_files)
1216  call netcdf_err(error, 'reading gridfiles')
1217 
1218  error = nf90_close(ncid)
1219 
1220  grid_file = trim(orog_dir) // trim(grid_files(tile))
1221 
1222  print*,'- OPEN GRID FILE: ', trim(grid_file)
1223  error=nf90_open(trim(grid_file), nf90_nowrite, ncid)
1224  call netcdf_err(error, 'opening grid file')
1225 
1226  print*,'- READ NXP ID'
1227  error=nf90_inq_dimid(ncid, 'nxp', id_dim)
1228  call netcdf_err(error, 'reading nxp id')
1229 
1230  print*,'- READ NXP'
1231  error=nf90_inquire_dimension(ncid,id_dim,len=nxp)
1232  call netcdf_err(error, 'reading nxp')
1233 
1234  print*,'- READ NYP ID'
1235  error=nf90_inq_dimid(ncid, 'nyp', id_dim)
1236  call netcdf_err(error, 'reading nyp id')
1237 
1238  print*,'- READ NYP'
1239  error=nf90_inquire_dimension(ncid,id_dim,len=nyp)
1240  call netcdf_err(error, 'reading nyp')
1241 
1242  if ((nxp/2 /= i_tile) .or. (nyp/2 /= j_tile)) then
1243  call error_handler("DIMENSION MISMATCH IN GRID FILE.", 1)
1244  endif
1245 
1246  allocate(tmpvar(nxp,nyp))
1247 
1248  print*,'- READ LONGITUDE ID'
1249  error=nf90_inq_varid(ncid, 'x', id_var)
1250  call netcdf_err(error, 'reading longitude id')
1251 
1252  print*,'- READ LONGITUDE'
1253  error=nf90_get_var(ncid, id_var, tmpvar)
1254  call netcdf_err(error, 'reading longitude')
1255 
1256  do j = 1, j_tile
1257  do i = 1, i_tile
1258  ii = 2*i
1259  jj = 2*j
1260  longitude(i,j) = tmpvar(ii,jj)
1261  enddo
1262  enddo
1263 
1264  do j = 1, jp1_tile
1265  do i = 1, i_tile
1266  ii = 2*i
1267  jj = (2*j) - 1
1268  longitude_s(i,j) = tmpvar(ii,jj)
1269  enddo
1270  enddo
1271 
1272  do j = 1, j_tile
1273  do i = 1, ip1_tile
1274  ii = (2*i) - 1
1275  jj = 2*j
1276  longitude_w(i,j) = tmpvar(ii,jj)
1277  enddo
1278  enddo
1279 
1280  print*,'- READ LATITUDE ID'
1281  error=nf90_inq_varid(ncid, 'y', id_var)
1282  call netcdf_err(error, 'reading latitude id')
1283 
1284  print*,'- READ LATIITUDE'
1285  error=nf90_get_var(ncid, id_var, tmpvar)
1286  call netcdf_err(error, 'reading latitude')
1287 
1288  do j = 1, j_tile
1289  do i = 1, i_tile
1290  ii = 2*i
1291  jj = 2*j
1292  latitude(i,j) = tmpvar(ii,jj)
1293  enddo
1294  enddo
1295 
1296  do j = 1, jp1_tile
1297  do i = 1, i_tile
1298  ii = 2*i
1299  jj = (2*j) - 1
1300  latitude_s(i,j) = tmpvar(ii,jj)
1301  enddo
1302  enddo
1303 
1304  do j = 1, j_tile
1305  do i = 1, ip1_tile
1306  ii = (2*i) - 1
1307  jj = 2*j
1308  latitude_w(i,j) = tmpvar(ii,jj)
1309  enddo
1310  enddo
1311 
1312  deallocate(tmpvar)
1313 
1314  error = nf90_close(ncid)
1315 
1316  end subroutine get_model_latlons
1317 
1327  subroutine get_model_mask_terrain(orog_file, idim, jdim, mask, terrain)
1328 
1329  use netcdf
1330 
1331  implicit none
1332 
1333  character(len=*), intent(in) :: orog_file
1334 
1335  integer, intent(in) :: idim, jdim
1336  integer(esmf_kind_i8), intent(out) :: mask(idim,jdim)
1337 
1338  real(esmf_kind_i8), intent(out) :: terrain(idim,jdim)
1339 
1340  integer :: error, lat, lon
1341  integer :: ncid, id_dim, id_var
1342 
1343  real(kind=4), allocatable :: dummy(:,:)
1344 
1345  print*,"- READ MODEL LAND MASK FILE"
1346 
1347  print*,'- OPEN LAND MASK FILE: ', orog_file
1348  error=nf90_open(orog_file,nf90_nowrite,ncid)
1349  call netcdf_err(error, 'opening land mask file')
1350 
1351  print*,"- READ I-DIMENSION"
1352  error=nf90_inq_dimid(ncid, 'lon', id_dim)
1353  call netcdf_err(error, 'reading idim id')
1354  error=nf90_inquire_dimension(ncid,id_dim,len=lon)
1355  call netcdf_err(error, 'reading idim')
1356 
1357  print*,"- READ J-DIMENSION"
1358  error=nf90_inq_dimid(ncid, 'lat', id_dim)
1359  call netcdf_err(error, 'reading jdim id')
1360  error=nf90_inquire_dimension(ncid,id_dim,len=lat)
1361  call netcdf_err(error, 'reading jdim')
1362 
1363  print*,"- I/J DIMENSIONS: ", lon, lat
1364 
1365  if ((lon /= idim) .or. (lat /= jdim)) then
1366  call error_handler("MISMATCH IN DIMENSIONS.", 1)
1367  endif
1368 
1369  allocate(dummy(idim,jdim))
1370 
1371  print*,"- READ LAND MASK"
1372  error=nf90_inq_varid(ncid, 'slmsk', id_var)
1373  call netcdf_err(error, 'reading slmsk id')
1374  error=nf90_get_var(ncid, id_var, dummy)
1375  call netcdf_err(error, 'reading slmsk')
1376  mask = nint(dummy)
1377 
1378  print*,"- READ RAW OROGRAPHY."
1379  error=nf90_inq_varid(ncid, 'orog_raw', id_var)
1380  call netcdf_err(error, 'reading orog_raw id')
1381  error=nf90_get_var(ncid, id_var, dummy)
1382  call netcdf_err(error, 'reading orog_raw')
1383  terrain = dummy
1384 
1385  error = nf90_close(ncid)
1386 
1387  deallocate (dummy)
1388 
1389  end subroutine get_model_mask_terrain
1390 
1395 
1396  implicit none
1397 
1398  integer :: rc
1399 
1400  print*,"- DESTROY MODEL DATA."
1401 
1402  call esmf_fielddestroy(latitude_input_grid,rc=rc)
1403  call esmf_fielddestroy(longitude_input_grid,rc=rc)
1404  if (esmf_fieldiscreated(latitude_s_input_grid)) then
1405  call esmf_fielddestroy(latitude_s_input_grid, rc=rc)
1406  endif
1407  if (esmf_fieldiscreated(latitude_w_input_grid)) then
1408  call esmf_fielddestroy(latitude_w_input_grid, rc=rc)
1409  endif
1410  if (esmf_fieldiscreated(longitude_s_input_grid)) then
1411  call esmf_fielddestroy(longitude_s_input_grid, rc=rc)
1412  endif
1413  if (esmf_fieldiscreated(longitude_w_input_grid)) then
1414  call esmf_fielddestroy(longitude_w_input_grid, rc=rc)
1415  endif
1416  call esmf_fielddestroy(landmask_target_grid, rc=rc)
1417  call esmf_fielddestroy(latitude_target_grid, rc=rc)
1418  if (esmf_fieldiscreated(latitude_s_target_grid)) then
1419  call esmf_fielddestroy(latitude_s_target_grid, rc=rc)
1420  endif
1421  if (esmf_fieldiscreated(latitude_w_target_grid)) then
1422  call esmf_fielddestroy(latitude_w_target_grid, rc=rc)
1423  endif
1424  call esmf_fielddestroy(longitude_target_grid, rc=rc)
1425  if (esmf_fieldiscreated(longitude_s_target_grid)) then
1426  call esmf_fielddestroy(longitude_s_target_grid, rc=rc)
1427  endif
1428  if (esmf_fieldiscreated(longitude_w_target_grid)) then
1429  call esmf_fielddestroy(longitude_w_target_grid, rc=rc)
1430  endif
1431  call esmf_fielddestroy(seamask_target_grid, rc=rc)
1432  call esmf_fielddestroy(terrain_target_grid, rc=rc)
1433  call esmf_griddestroy(input_grid, rc=rc)
1434  call esmf_griddestroy(target_grid, rc=rc)
1435 
1436  end subroutine cleanup_input_target_grid_data
1437 
1449  subroutine gdt_to_gds(igdtnum, igdstmpl, igdtlen, kgds, ni, nj, res)
1450 
1451  implicit none
1452 
1453  integer, intent(in ) :: igdtnum, igdtlen, igdstmpl(igdtlen)
1454  integer, intent( out) :: kgds(200), ni, nj
1455  integer :: iscale
1456 
1457  real, intent( out) :: res
1458 
1459  kgds=0
1460 
1461  if (igdtnum.eq.32769) then ! rot lat/lon b grid
1462 
1463  iscale=igdstmpl(10)*igdstmpl(11)
1464  if (iscale == 0) iscale = 1e6
1465  kgds(1)=205 ! oct 6, rotated lat/lon for Non-E
1466  ! Stagger grid
1467  kgds(2)=igdstmpl(8) ! octs 7-8, Ni
1468  ni = kgds(2)
1469  kgds(3)=igdstmpl(9) ! octs 9-10, Nj
1470  nj = kgds(3)
1471  kgds(4)=nint(float(igdstmpl(12))/float(iscale)*1000.) ! octs 11-13, Lat of
1472  ! 1st grid point
1473  kgds(5)=nint(float(igdstmpl(13))/float(iscale)*1000.) ! octs 14-16, Lon of
1474  ! 1st grid point
1475 
1476  kgds(6)=0 ! oct 17, resolution and component flags
1477  if (igdstmpl(1)==2 ) kgds(6)=64
1478  if ( btest(igdstmpl(14),4).OR.btest(igdstmpl(14),5) ) kgds(6)=kgds(6)+128
1479  if ( btest(igdstmpl(14),3) ) kgds(6)=kgds(6)+8
1480 
1481  kgds(7)=nint(float(igdstmpl(15))/float(iscale)*1000.) ! octs 18-20,
1482  ! Lat of cent of rotation
1483  kgds(8)=nint(float(igdstmpl(16))/float(iscale)*1000.) ! octs 21-23,
1484  ! Lon of cent of rotation
1485  kgds(9)=nint(float(igdstmpl(17))/float(iscale)*1000.) ! octs 24-25,
1486  ! Di
1487  kgds(10)=nint(float(igdstmpl(18))/float(iscale)*1000.) ! octs 26-27,
1488  ! Dj
1489 
1490  kgds(11) = 0 ! oct 28, scan mode
1491  if (btest(igdstmpl(19),7)) kgds(11) = 128
1492  if (btest(igdstmpl(19),6)) kgds(11) = kgds(11) + 64
1493  if (btest(igdstmpl(19),5)) kgds(11) = kgds(11) + 32
1494 
1495  kgds(12)=nint(float(igdstmpl(20))/float(iscale)*1000.) ! octs 29-31, Lat of
1496  ! last grid point
1497  kgds(13)=nint(float(igdstmpl(21))/float(iscale)*1000.) ! octs 32-34, Lon of
1498  ! last grid point
1499 
1500  kgds(19)=0 ! oct 4, # vert coordinate parameters
1501  kgds(20)=255 ! oct 5, used for thinned grids, set to 255
1502 
1503  res = ((float(kgds(9)) / 1000.0) + (float(kgds(10)) / 1000.0)) &
1504  * 0.5 * 111.0
1505 
1506  elseif(igdtnum==30) then
1507 
1508  kgds(1)=3 ! oct 6, lambert conformal
1509  kgds(2)=igdstmpl(8) ! octs 7-8, Ni
1510  ni = kgds(2)
1511  kgds(3)=igdstmpl(9) ! octs 9-10, Nj
1512  nj = kgds(3)
1513 
1514  iscale = 1e6
1515  kgds(4) = nint(float(igdstmpl(10))/1000.0)
1516  kgds(5) = nint(float(igdstmpl(11))/1000.0)
1517 
1518  kgds(6)=0 ! oct 17, resolution and component flags
1519  if (igdstmpl(1)==2 ) kgds(6)=64
1520  if ( btest(igdstmpl(12),4).OR.btest(igdstmpl(12),5) ) kgds(6)=kgds(6)+128
1521  if ( btest(igdstmpl(12),3) ) kgds(6)=kgds(6)+8
1522 
1523  kgds(7) = nint(float(igdstmpl(14))/1000.0)
1524  kgds(8) = nint(float(igdstmpl(15))/1000.0)
1525  kgds(9) = nint(float(igdstmpl(16))/1000.0)
1526  kgds(10) = 0
1527 
1528  kgds(11) = 0 ! oct 28, scan mode
1529  if (btest(igdstmpl(18),7)) kgds(11) = 128
1530  if (btest(igdstmpl(18),6)) kgds(11) = kgds(11) + 64
1531  if (btest(igdstmpl(18),5)) kgds(11) = kgds(11) + 32
1532 
1533  kgds(12) = nint(float(igdstmpl(19))/1000.0)
1534  kgds(13) = nint(float(igdstmpl(20))/1000.0)
1535  kgds(14) = -90
1536  kgds(15) = 0
1537 
1538  elseif(igdtnum==0) then ! lat/lon grid
1539 
1540  iscale=igdstmpl(10)*igdstmpl(11)
1541  if (iscale == 0) iscale = 1e6
1542  kgds(1)=0 ! oct 6, data representation type.
1543  kgds(2)=igdstmpl(8) ! octs 7-8, Ni
1544  ni = kgds(2)
1545  kgds(3)=igdstmpl(9) ! octs 9-10, Nj
1546  nj = kgds(3)
1547  kgds(4)=nint(float(igdstmpl(12))/float(iscale)*1000.) ! octs 11-13, Lat of 1st grid point
1548  kgds(5)=nint(float(igdstmpl(13))/float(iscale)*1000.) ! octs 14-16, Lon of 1st grid point
1549 
1550  kgds(6)=0 ! oct 17, resolution and component flags
1551  if (igdstmpl(1)==2 ) kgds(6)=64
1552  if ( btest(igdstmpl(14),4).OR.btest(igdstmpl(14),5) ) kgds(6)=kgds(6)+128
1553  if ( btest(igdstmpl(14),3) ) kgds(6)=kgds(6)+8
1554 
1555  kgds(7)=nint(float(igdstmpl(15))/float(iscale)*1000.) ! octs 18-20, Lat of last grid point
1556  kgds(8)=nint(float(igdstmpl(16))/float(iscale)*1000.) ! octs 21-23, Lon of last grid point
1557  kgds(9)=nint(float(igdstmpl(17))/float(iscale)*1000.) ! octs 24-25, "i" resolution.
1558  kgds(10)=nint(float(igdstmpl(18))/float(iscale)*1000.) ! octs 26-27, "j" resolution.
1559 
1560  kgds(11) = 0 ! oct 28, scan mode
1561  if (btest(igdstmpl(19),7)) kgds(11) = 128
1562  if (btest(igdstmpl(19),6)) kgds(11) = kgds(11) + 64
1563  if (btest(igdstmpl(19),5)) kgds(11) = kgds(11) + 32
1564 
1565  kgds(12)=0 ! octs 29-32, reserved
1566  kgds(19)=0 ! oct 4, # vert coordinate parameters
1567  kgds(20)=255 ! oct 5, used for thinned grids, set to 255
1568 
1569  else
1570 
1571  call error_handler("UNRECOGNIZED INPUT GRID TYPE ", 1)
1572 
1573  endif
1574 
1575  end subroutine gdt_to_gds
1576 
1577  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 gdt_to_gds(igdtnum, igdstmpl, igdtlen, kgds, ni, nj, res)
Convert the GRIB2 grid description template to to the GRIB1 grid description section.
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 define_input_grid_grib2(localpet, npets)
Define input grid object for grib2 input data.
Definition: model_grid.F90:615
subroutine netcdf_err(err, string)
Error handler for netcdf.
Definition: utils.F90:31
subroutine define_input_grid_gaussian(localpet, npets)
Define grid object for input data on global gaussian grids.
Definition: model_grid.F90:147
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
subroutine error_handler(string, rc)
General error handler.
Definition: utils.F90:9
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:115
subroutine, public define_target_grid(localpet, npets)
Setup the esmf grid object for the target grid.
Definition: model_grid.F90:888