FV3-JEDI
fv3jedi_io_geos_mod.f90
Go to the documentation of this file.
1 ! (C) Copyright 2017-2020 UCAR
2 !
3 ! This software is licensed under the terms of the Apache Licence Version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
5 
7 
8 ! iso
9 use iso_c_binding
10 
11 ! libs
12 use mpi
13 use netcdf
14 
15 ! fckit
16 use fckit_configuration_module, only: fckit_configuration
17 use fckit_mpi_module
18 
19 ! oops
20 use datetime_mod
21 
22 ! fv3-jedi
29 
30 implicit none
31 private
32 public fv3jedi_io_geos
33 
34 integer, parameter :: numfiles = 5 ! May change as more restarts are used
35 
37  ! File names and paths
38  integer :: ncid(numfiles)
39  integer, allocatable :: ncid_forfield(:)
40  logical :: ncid_isneeded(numfiles)
41  ! Whether to expect/use the tile dimension
42  logical :: tiledim(numfiles) = .false.
43  logical:: restart(numfiles) = .false.
44  ! Filename
45  character(len=maxstring) :: datapath = ''
46  character(len=maxstring) :: filenames(numfiles) = ''
47  character(len=maxstring) :: filenames_conf(numfiles) = ''
48  character(len=maxstring) :: filenames_default(numfiles) = ''
49  ! Write extra meta data needed for GEOS ingest
50  logical :: geosingestmeta = .false.
51  ! Comms
52  logical :: iam_io_proc
53  type(fckit_mpi_comm) :: ccomm
54  integer :: tcomm, ocomm !Communicator for each tile and for output
55  integer :: trank, tsize !Tile come info
56  integer :: crank, csize !Component comm info
57  integer :: orank, osize !Output comm info
58  ! IO ranges
59  integer :: is_r3_tile(5), ic_r3_tile(5)
60  integer :: is_r2_tile(4), ic_r2_tile(4)
61  integer :: is_r3_noti(4), ic_r3_noti(4)
62  integer :: is_r2_noti(3), ic_r2_noti(3)
63  integer :: vindex_tile = 4
64  integer :: vindex_noti = 3
65  ! Clobber option
66  logical :: clobber = .true.
67  integer :: x_dimid, y_dimid, n_dimid, z_dimid, e_dimid, t_dimid, f_dimid, c_dimid, o_dimid
68  ! Ps in file option
69  logical :: ps_in_file = .false.
70  contains
71  procedure :: setup_conf
72  procedure :: setup_date
73  procedure :: delete
74  procedure :: read_meta
75  procedure :: read_fields
76  procedure :: write
77  final :: dummy_final
78 end type fv3jedi_io_geos
79 
80 ! ------------------------------------------------------------------------------
81 
82 contains
83 
84 ! ------------------------------------------------------------------------------
85 
86 subroutine setup_conf(self, geom, f_conf)
87 
88 class(fv3jedi_io_geos), intent(inout) :: self
89 type(fv3jedi_geom), intent(in) :: geom
90 type(fckit_configuration), intent(in) :: f_conf
91 
92 integer :: ierr, n, var
93 
94 integer :: tileoffset, dt_in_name
95 character(len=4) :: yyyy
96 character(len=2) :: mm, dd, hh, min, ss
97 
98 ! Default file names
99 ! ------------------
100 n = 0
101 n = n+1; self%filenames_default(n) = 'bkg'
102 n = n+1; self%filenames_default(n) = 'crtmsrf'
103 n = n+1; self%filenames_default(n) = 'fvcore_internal_rst'
104 n = n+1; self%filenames_default(n) = 'moist_internal_rst'
105 n = n+1; self%filenames_default(n) = 'surf_import_rst'
106 
107 if (n .ne. numfiles) &
108  call abor1_ftn("fv3jedi_io_geos_mod.setup: number of potential restart files &
109  does not match numfiles")
110 
111 ! Set default tile dim and restart flags
112 ! --------------------------------------
113 self%tiledim(1) = .true. ! History defult tile dim
114 self%tiledim(2) = .true. ! History defult tile dim
115 self%tiledim(3) = .false. ! Restarts do not use tiledim
116 self%tiledim(4) = .false. ! Restarts do not use tiledim
117 self%tiledim(5) = .false. ! Restarts do not use tiledim
118 
119 self%restart(1) = .false. !Is not a restart file
120 self%restart(2) = .false. !Is not a restart file
121 self%restart(3) = .true. !Is a restart file
122 self%restart(4) = .true. !Is a restart file
123 self%restart(5) = .true. !Is a restart file
124 
125 ! Get configuration
126 ! -----------------
127 call get_conf(self, f_conf)
128 
129 ! Config filenames to filenames
130 ! -----------------------------
131 do n = 1,numfiles
132  self%filenames(n) = trim(self%filenames_conf(n))
133 enddo
134 
135 ! Component communicator / get the main communicator from fv3 geometry
136 ! --------------------------------------------------------------------
137 self%ccomm = geom%f_comm
138 self%csize = self%ccomm%size()
139 self%crank = self%ccomm%rank()
140 
141 ! Split comm and set flag for IO procs
142 ! ------------------------------------
143 self%iam_io_proc = .true.
144 
145 if (self%csize > 6) then
146 
147  ! Communicator for all procs on each tile. To communicate to tiles write proc
148  call mpi_comm_split(self%ccomm%communicator(), geom%ntile, self%ccomm%rank(), self%tcomm, ierr)
149  call mpi_comm_rank(self%tcomm, self%trank, ierr)
150  call mpi_comm_size(self%tcomm, self%tsize, ierr)
151 
152  ! Communicator for procs that will write, one per tile
153  call mpi_comm_split(self%ccomm%communicator(), self%trank, geom%ntile, self%ocomm, ierr)
154  call mpi_comm_rank(self%ocomm, self%orank, ierr)
155  call mpi_comm_size(self%ocomm, self%osize, ierr)
156 
157  if (self%trank .ne. 0) self%iam_io_proc = .false.
158 
159 else
160 
161  ! For 6 cores output comm is componenent comm
162  call mpi_comm_dup(self%ccomm%communicator(), self%ocomm, ierr)
163 
164 endif
165 
166 ! Create local to this proc start/count
167 ! -------------------------------------
168 if (self%iam_io_proc) then
169 
170  ! Starts/counts with tile dimension
171  self%is_r3_tile(1) = 1; self%ic_r3_tile(1) = geom%npx-1 !X
172  self%is_r3_tile(2) = 1; self%ic_r3_tile(2) = geom%npy-1 !Y
173  self%is_r3_tile(3) = geom%ntile; self%ic_r3_tile(3) = 1 !Tile
174  self%is_r3_tile(4) = 1; self%ic_r3_tile(4) = 1 !Lev
175  self%is_r3_tile(5) = 1; self%ic_r3_tile(5) = 1 !Time
176  self%is_r2_tile(1) = 1; self%ic_r2_tile(1) = geom%npx-1
177  self%is_r2_tile(2) = 1; self%ic_r2_tile(2) = geom%npy-1
178  self%is_r2_tile(3) = geom%ntile; self%ic_r2_tile(3) = 1
179  self%is_r2_tile(4) = 1; self%ic_r2_tile(4) = 1
180 
181  ! Starts/counts with no tile dimension
182  tileoffset = (geom%ntile-1)*(6*(geom%npy-1)/geom%ntiles)
183  self%is_r3_noti(1) = 1; self%ic_r3_noti(1) = geom%npx-1
184  self%is_r3_noti(2) = tileoffset+1; self%ic_r3_noti(2) = geom%npy-1
185  self%is_r3_noti(3) = 1; self%ic_r3_noti(3) = 1
186  self%is_r3_noti(4) = 1; self%ic_r3_noti(4) = 1
187  self%is_r2_noti(1) = 1; self%ic_r2_noti(1) = geom%npx-1
188  self%is_r2_noti(2) = tileoffset+1; self%ic_r2_noti(2) = geom%npy-1
189  self%is_r2_noti(3) = 1; self%ic_r2_noti(3) = 1
190 
191 endif
192 
193 end subroutine setup_conf
194 
195 ! ------------------------------------------------------------------------------
196 
197 subroutine setup_date(self, vdate)
198 
199 class(fv3jedi_io_geos), intent(inout) :: self
200 type(datetime), intent(in) :: vdate
201 
202 integer :: n
203 character(len=4) :: yyyy
204 character(len=2) :: mm, dd, hh, min, ss
205 
206 ! Datetime to strings
207 ! -------------------
208 call vdate_to_datestring(vdate, yyyy=yyyy, mm=mm, dd=dd, hh=hh, min=min, ss=ss)
209 
210 do n = 1, numfiles
211 
212  ! Config filenames to filenames
213  ! -----------------------------
214  self%filenames(n) = trim(self%filenames_conf(n))
215 
216  ! Swap out datetime templates if needed
217  if (index(self%filenames(n),"%yyyy") > 0) &
218  self%filenames(n) = replace_text(self%filenames(n),'%yyyy',yyyy)
219  if (index(self%filenames(n),"%mm" ) > 0) &
220  self%filenames(n) = replace_text(self%filenames(n),'%mm' ,mm )
221  if (index(self%filenames(n),"%dd" ) > 0) &
222  self%filenames(n) = replace_text(self%filenames(n),'%dd' ,dd )
223  if (index(self%filenames(n),"%hh" ) > 0) &
224  self%filenames(n) = replace_text(self%filenames(n),'%hh' ,hh )
225  if (index(self%filenames(n),"%MM" ) > 0) &
226  self%filenames(n) = replace_text(self%filenames(n),'%MM' ,min )
227  if (index(self%filenames(n),"%ss" ) > 0) &
228  self%filenames(n) = replace_text(self%filenames(n),'%ss' ,ss )
229 
230 enddo
231 
232 end subroutine setup_date
233 
234 ! ------------------------------------------------------------------------------
235 
236 subroutine get_conf(self,f_conf)
237 
238 implicit none
239 class(fv3jedi_io_geos), intent(inout) :: self
240 type(fckit_configuration), intent(in) :: f_conf
241 
242 integer :: n
243 
244 ! Path where files are read from or saved to
245 ! ------------------------------------------
246 call string_from_conf(f_conf,"datapath",self%datapath,'Data',memberswap=.true.)
247 
248 ! User can ask for extra meta data needed for GEOS ingest
249 ! -------------------------------------------------------
250 if (.not. f_conf%get('geosingestmeta',self%geosingestmeta)) self%geosingestmeta = .false.
251 
252 ! Whether to expect/use the tile dimenstion in the file
253 ! -----------------------------------------------------
254 if (.not. f_conf%get('clobber',self%clobber)) self%clobber = .true.
255 
256 ! Whether to expect/use the tile dimenstion in the file
257 ! -----------------------------------------------------
258 if (.not. f_conf%get('tiledim',self%tiledim(1))) self%tiledim(1) = .true.
259 if (.not. f_conf%get('tiledim',self%tiledim(2))) self%tiledim(2) = .true.
260 
261 ! User can optionally specify the file names
262 ! ------------------------------------------
263 n = 0
264 n = n+1; call string_from_conf(f_conf,"filename_bkgd",self%filenames_conf(1), &
265  self%filenames_default(1),memberswap=.true.)
266 n = n+1; call string_from_conf(f_conf,"filename_crtm",self%filenames_conf(2), &
267  self%filenames_default(2),memberswap=.true.)
268 n = n+1; call string_from_conf(f_conf,"filename_core",self%filenames_conf(3), &
269  self%filenames_default(3),memberswap=.true.)
270 n = n+1; call string_from_conf(f_conf,"filename_mois",self%filenames_conf(4), &
271  self%filenames_default(4),memberswap=.true.)
272 n = n+1; call string_from_conf(f_conf,"filename_surf",self%filenames_conf(5), &
273  self%filenames_default(5),memberswap=.true.)
274 
275 ! Sanity check
276 ! ------------
277 if (n .ne. numfiles) &
278  call abor1_ftn("fv3jedi_io_geos_mod.get_conf: number of potential restart files &
279  does not match numfiles")
280 
281 ! Option to allow for ps infile
282 if (f_conf%has("psinfile")) call f_conf%get_or_die("psinfile",self%ps_in_file)
283 
284 end subroutine get_conf
285 
286 ! ------------------------------------------------------------------------------
287 
288 subroutine delete(self)
289 
290 implicit none
291 class(fv3jedi_io_geos), intent(inout) :: self
292 
293 integer :: ierr, n
294 
295 ! Deallocate
296 ! ----------
297 if (allocated(self%ncid_forfield)) deallocate(self%ncid_forfield)
298 
299 ! Release split comms
300 ! -------------------
301 if (self%csize > 6) call mpi_comm_free(self%tcomm, ierr)
302 call mpi_comm_free(self%ocomm, ierr)
303 
304 end subroutine delete
305 
306 ! ------------------------------------------------------------------------------
307 
308 subroutine read_meta(self, geom, vdate, calendar_type, date_init, fields)
309 
310 implicit none
311 class(fv3jedi_io_geos), intent(inout) :: self
312 type(fv3jedi_geom), intent(inout) :: geom !< Geometry
313 type(datetime), intent(inout) :: vdate !< DateTime
314 integer, intent(inout) :: calendar_type !< Calendar type
315 integer, intent(inout) :: date_init(6) !< Date intialized
316 type(fv3jedi_field), intent(in) :: fields(:)
317 
318 integer :: varid, date(6), intdate, inttime, idate, isecs
319 character(len=8) :: cdate
320 character(len=6) :: ctime
321 integer(kind=c_int) :: cidate, cisecs, n, df
322 
323 ! Set filenames
324 ! -------------
325 call set_file_names(self, fields)
326 
327 ! Open files
328 ! ----------
329 call open_files(self)
330 
331 calendar_type = -1
332 date_init = 0
333 
334 idate = 0
335 isecs = 0
336 
337 if (self%iam_io_proc) then
338 
339  ! Get time attributes
340  do n = 1,numfiles
341  if (self%ncid_isneeded(n)) then
342  df = n
343  exit
344  endif
345  enddo
346 
347  call nccheck ( nf90_inq_varid(self%ncid(df), "time", varid), "nf90_inq_varid time" )
348  call nccheck ( nf90_get_att(self%ncid(df), varid, "begin_date", intdate), &
349  "nf90_get_att begin_date" )
350  call nccheck ( nf90_get_att(self%ncid(df), varid, "begin_time", inttime), &
351  "nf90_get_att begin_time" )
352 
353  ! Pad with leading zeros if need be
354  write(cdate,"(I0.8)") intdate
355  write(ctime,"(I0.6)") inttime
356 
357  ! Back to integer
358  read(cdate(1:4),*) date(1)
359  read(cdate(5:6),*) date(2)
360  read(cdate(7:8),*) date(3)
361  read(ctime(1:2),*) date(4)
362  read(ctime(3:4),*) date(5)
363  read(ctime(5:6),*) date(6)
364 
365  ! To idate/isecs for Jedi
366  idate = date(1)*10000 + date(2)*100 + date(3)
367  isecs = date(4)*3600 + date(5)*60 + date(6)
368 
369 endif
370 
371 ! Set the object date from the date of the file
372 call self%ccomm%broadcast(idate,0)
373 call self%ccomm%broadcast(isecs,0)
374 cidate = idate
375 cisecs = isecs
376 call datetime_from_ifs(vdate, cidate, cisecs)
377 
378 ! Close files
379 ! -----------
380 call close_files(self)
381 
382 end subroutine read_meta
383 
384 ! ------------------------------------------------------------------------------
385 
386 subroutine read_fields(self, geom, fields)
387 
388 implicit none
389 class(fv3jedi_io_geos), target, intent(inout) :: self
390 type(fv3jedi_geom), intent(inout) :: geom
391 type(fv3jedi_field), intent(inout) :: fields(:)
392 
393 integer :: varid, var, lev, ncid
394 logical :: tiledim
395 integer, pointer :: istart(:), icount(:)
396 integer, allocatable, target :: is_r3_tile(:), is_r3_noti(:)
397 real(kind=kind_real), allocatable :: arrayg(:,:), delp(:,:,:)
398 
399 ! Set filenames
400 ! -------------
401 call set_file_names(self, fields)
402 
403 ! Open files
404 ! ----------
405 call open_files(self)
406 
407 ! Local copy of starts for rank 3 in order to do one level at a time
408 ! ------------------------------------------------------------------
409 allocate(is_r3_tile(size(self%is_r3_tile)))
410 is_r3_tile = self%is_r3_tile
411 allocate(is_r3_noti(size(self%is_r3_noti)))
412 is_r3_noti = self%is_r3_noti
413 
414 ! Array for level of whole tile
415 ! -----------------------------
416 allocate(arrayg(1:geom%npx-1,1:geom%npy-1))
417 
418 
419 ! Loop over fields
420 ! ----------------
421 do var = 1,size(fields)
422 
423  ! ncid for this variable
424  ncid = self%ncid(self%ncid_forfield(var))
425 
426  !If ps then switch to delp
427  if (.not. self%ps_in_file) then
428  if (trim(fields(var)%fv3jedi_name) == 'ps') then
429  fields(var)%short_name = 'delp'
430  fields(var)%npz = geom%npz
431  allocate(delp(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
432  endif
433  endif
434 
435  ! Set pointers to the appropriate array ranges
436  ! --------------------------------------------
437  if (self%iam_io_proc) then
438 
439  tiledim = self%tiledim(self%ncid_forfield(var))
440 
441  if (associated(istart)) nullify(istart)
442  if (associated(icount)) nullify(icount)
443 
444  if (fields(var)%npz == 1) then
445  if (tiledim) then
446  istart => self%is_r2_tile
447  icount => self%ic_r2_tile
448  else
449  istart => self%is_r2_noti
450  icount => self%ic_r2_noti
451  endif
452  elseif (fields(var)%npz > 1) then
453  if (tiledim) then
454  istart => is_r3_tile;
455  icount => self%ic_r3_tile
456  else
457  istart => is_r3_noti
458  icount => self%ic_r3_noti
459  endif
460  endif
461  endif
462 
463  if (self%iam_io_proc) then
464  call nccheck ( nf90_inq_varid(ncid, trim(fields(var)%short_name), varid), &
465  "nf90_inq_varid "//trim(fields(var)%short_name) )
466  endif
467 
468  ! Loop over level and read the data
469  ! ---------------------------------
470  do lev = 1,fields(var)%npz
471 
472  arrayg = 0.0_kind_real
473 
474  if (self%iam_io_proc) then
475 
476  !Set start to current level
477  is_r3_tile(self%vindex_tile) = lev
478  is_r3_noti(self%vindex_noti) = lev
479 
480  ! Read the level
481  call nccheck ( nf90_get_var( ncid, varid, arrayg, istart, icount), &
482  "nf90_get_var "//trim(fields(var)%short_name) )
483  endif
484 
485 
486  ! Scatter the field to all processors on the tile
487  ! -----------------------------------------------
488  if (.not. self%ps_in_file .and. trim(fields(var)%fv3jedi_name) == 'ps') then
489  if (self%csize > 6) then
490  call scatter_tile(geom, self%tcomm, 1, arrayg, delp(geom%isc:geom%iec,geom%jsc:geom%jec,lev))
491  else
492  delp(geom%isc:geom%iec,geom%jsc:geom%jec,lev) = arrayg(geom%isc:geom%iec,geom%jsc:geom%jec)
493  endif
494  else
495  if (self%csize > 6) then
496  call scatter_tile(geom, self%tcomm, 1, arrayg, &
497  fields(var)%array(geom%isc:geom%iec,geom%jsc:geom%jec,lev))
498  else
499  fields(var)%array(geom%isc:geom%iec,geom%jsc:geom%jec,lev) = &
500  arrayg(geom%isc:geom%iec,geom%jsc:geom%jec)
501  endif
502  endif
503 
504  enddo
505 
506  if (.not. self%ps_in_file .and. trim(fields(var)%fv3jedi_name) == 'ps') then
507  fields(var)%short_name = 'ps'
508  fields(var)%npz = 1
509  fields(var)%array(:,:,1) = sum(delp,3)
510  deallocate(delp)
511  endif
512 
513 enddo
514 
515 
516 ! Deallocate locals
517 ! -----------------
518 deallocate(is_r3_tile,is_r3_noti)
519 deallocate(arrayg)
520 
521 ! Close files
522 ! -----------
523 call close_files(self)
524 
525 end subroutine read_fields
526 
527 ! ------------------------------------------------------------------------------
528 
529 subroutine write(self, geom, fields, vdate)
530 
531 implicit none
532 class(fv3jedi_io_geos), intent(inout) :: self
533 type(fv3jedi_geom), intent(inout) :: geom !< Geom
534 type(fv3jedi_field), intent(in) :: fields(:) !< Fields to be written
535 type(datetime), intent(in) :: vdate !< DateTime
536 
537 ! Set filenames
538 ! -------------
539 call set_file_names(self, fields)
540 
541 ! Open/create files
542 ! -----------------
543 call create_files(self)
544 
545 ! Write meta data
546 ! ---------------
547 if (self%clobber) call write_meta(self, geom, fields, vdate)
548 
549 ! Write fields
550 ! ------------
551 call write_fields(self, geom, fields, vdate)
552 
553 ! Close files
554 ! -----------
555 call close_files(self)
556 
557 end subroutine write
558 
559 ! ------------------------------------------------------------------------------
560 
561 subroutine write_meta(self, geom, fields, vdate)
562 
563 implicit none
564 class(fv3jedi_io_geos), target, intent(inout) :: self
565 type(fv3jedi_geom), intent(inout) :: geom !< Geom
566 type(fv3jedi_field), intent(in) :: fields(:) !< Fields to be written
567 type(datetime), intent(in) :: vdate !< DateTime
568 
569 integer :: var, ymult, k, n, vc
570 character(len=15) :: datefile
571 integer :: date(6), date8, time6
572 character(len=8) :: date8s, cubesize
573 character(len=6) :: time6s
574 character(len=4) :: XdimVar, YdimVar
575 integer :: varid(10000)
576 
577 integer, pointer :: istart(:), icount(:)
578 integer, allocatable :: dimidsg(:), tiles(:), levels(:)
579 real(kind=kind_real), allocatable :: latg(:,:), long(:,:), xdimydim(:)
580 
581 
582 ! Gathered lats/lons
583 ! ------------------
584 allocate(latg(1:geom%npx-1,1:geom%npy-1))
585 allocate(long(1:geom%npx-1,1:geom%npy-1))
586 
587 if (self%csize > 6) then
588  call gather_tile(geom, self%tcomm, 1, rad2deg*geom%grid_lat(geom%isc:geom%iec,geom%jsc:geom%jec), latg)
589  call gather_tile(geom, self%tcomm, 1, rad2deg*geom%grid_lon(geom%isc:geom%iec,geom%jsc:geom%jec), long)
590 else
591  latg = rad2deg*geom%grid_lat(geom%isc:geom%iec,geom%jsc:geom%jec)
592  long = rad2deg*geom%grid_lon(geom%isc:geom%iec,geom%jsc:geom%jec)
593 endif
594 
595 write(cubesize,'(I8)') geom%npx-1
596 
597 ! IO processors write the metadata
598 ! --------------------------------
599 if (self%iam_io_proc) then
600 
601  ! Get datetime information ready to write
602  ! ---------------------------------------
603  call vdate_to_datestring(vdate,datest=datefile,date=date)
604  write(date8s,'(I4,I0.2,I0.2)') date(1),date(2),date(3)
605  write(time6s,'(I0.2,I0.2,I0.2)') date(4),date(5),date(6)
606  read(date8s,*) date8
607  read(time6s,*) time6
608 
609 
610  ! Xdim, Ydim arrays
611  ! -----------------
612  allocate(xdimydim(geom%npx-1))
613  do k = 1,geom%npx-1
614  xdimydim(k) = real(k,kind_real)
615  enddo
616 
617 
618  ! Tile array
619  ! ----------
620  allocate(tiles(6))
621  do k = 1,6
622  tiles(k) = k
623  enddo
624 
625 
626  ! Level and edge arrays
627  ! ---------------------
628  allocate(levels(geom%npz+1))
629  do k = 1,geom%npz+1
630  levels(k) = k
631  enddo
632 
633 
634  ! Loop over all files to be created/written to
635  ! --------------------------------------------
636  do n = 1, numfiles
637 
638  if (self%ncid_isneeded(n)) then
639 
640  ! Multiplication factor when no tile dimension
641  ! --------------------------------------------
642  ymult = 1
643  if (.not. self%tiledim(n)) ymult = 6
644 
645  if ( self%tiledim(n) ) then
646  xdimvar = 'Xdim'
647  ydimvar = 'Ydim'
648  else
649  xdimvar = 'lon'
650  ydimvar = 'lat'
651  endif
652 
653  ! Main dimensions
654  call nccheck ( nf90_def_dim(self%ncid(n), trim(xdimvar), geom%npx-1, self%x_dimid), &
655  "nf90_def_dim "//trim(xdimvar) )
656  call nccheck ( nf90_def_dim(self%ncid(n), trim(ydimvar), ymult*(geom%npy-1), self%y_dimid), &
657  "nf90_def_dim "//trim(ydimvar) )
658  if (self%tiledim(n)) &
659  call nccheck ( nf90_def_dim(self%ncid(n), "n", geom%ntiles, self%n_dimid), "nf90_def_dim n" )
660  call nccheck ( nf90_def_dim(self%ncid(n), "lev", geom%npz, self%z_dimid), "nf90_def_dim lev" )
661  call nccheck ( nf90_def_dim(self%ncid(n), "edge", geom%npz+1, self%e_dimid), "nf90_def_dim edge" )
662  call nccheck ( nf90_def_dim(self%ncid(n), "time", 1, self%t_dimid), "nf90_def_dim time" )
663 
664  ! In case the four level surface fields need to be written
665  do var = 1,size(fields)
666  if (fields(var)%npz == 4) then
667  call nccheck ( nf90_def_dim(self%ncid(n), "lev4", 4, self%f_dimid), "nf90_def_dim lev" )
668  exit
669  endif
670  enddo
671 
672  !Needed by GEOS for ingesting cube sphere field
673  call nccheck ( nf90_def_dim(self%ncid(n), "ncontact", 4, self%c_dimid), &
674  "nf90_def_dim ncontact" )
675  call nccheck ( nf90_def_dim(self%ncid(n), "orientationStrLen", 5, self%o_dimid), &
676  "nf90_def_dim orientationStrLend" )
677 
678  ! Dimension ID array for lat/lon arrays
679  if (allocated(dimidsg)) deallocate(dimidsg)
680  if ( self%tiledim(n) ) then
681  allocate(dimidsg(3))
682  dimidsg(:) = (/ self%x_dimid, self%y_dimid, self%n_dimid /)
683  else
684  allocate(dimidsg(2))
685  dimidsg(:) = (/ self%x_dimid, self%y_dimid /)
686  endif
687 
688  ! Define fields to be written (geom)
689  vc=0;
690 
691  if (self%tiledim(n)) then
692  vc=vc+1;
693  call nccheck( nf90_def_var(self%ncid(n), "n", nf90_int, self%n_dimid, varid(vc)), "nf90_def_var n" )
694  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "long_name", "cubed-sphere face") )
695  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "axis", "e") )
696  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "grads_dim", "e") )
697  endif
698 
699  vc=vc+1;
700  call nccheck( nf90_def_var(self%ncid(n), trim(xdimvar), nf90_double, self%x_dimid, varid(vc)), &
701  "nf90_def_var "//trim(xdimvar) )
702  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "long_name", "Fake Longitude for GrADS Compatibility") )
703  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "units", "degrees_east") )
704 
705  vc=vc+1;
706  call nccheck( nf90_def_var(self%ncid(n), trim(ydimvar), nf90_double, self%y_dimid, varid(vc)), &
707  "nf90_def_var "//trim(ydimvar) )
708  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "long_name", "Fake Latitude for GrADS Compatibility") )
709  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "units", "degrees_north") )
710 
711  vc=vc+1;
712  call nccheck( nf90_def_var(self%ncid(n), "lons", nf90_double, dimidsg, varid(vc)), "nf90_def_var lons" )
713  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "long_name", "longitude") )
714  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "units", "degrees_east") )
715 
716  vc=vc+1;
717  call nccheck( nf90_def_var(self%ncid(n), "lats", nf90_double, dimidsg, varid(vc)), "nf90_def_var lats" )
718  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "long_name", "latitude") )
719  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "units", "degrees_north") )
720 
721  vc=vc+1;
722  call nccheck( nf90_def_var(self%ncid(n), "lev", nf90_double, self%z_dimid, varid(vc)), "nf90_def_var lev" )
723  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "long_name", "vertical level") )
724  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "units", "layer") )
725  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "positive", "down") )
726  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "coordinate", "eta") )
727  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "standard_name", "model_layers") )
728 
729  vc=vc+1;
730  call nccheck( nf90_def_var(self%ncid(n), "edge", nf90_double, self%e_dimid, varid(vc)), "nf90_def_var edge" )
731  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "long_name", "vertical level edges") )
732  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "units", "layer") )
733  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "positive", "down") )
734  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "coordinate", "eta") )
735  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "standard_name", "model_layers") )
736 
737  vc=vc+1;
738  call nccheck( nf90_def_var(self%ncid(n), "time", nf90_int, self%t_dimid, varid(vc)), "nf90_def_var time" )
739  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "long_name", "time"), "nf90_def_var time long_name" )
740  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "begin_date", date8), "nf90_def_var time begin_date" )
741  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "begin_time", time6), "nf90_def_var time begin_time" )
742 
743  vc=vc+1; !(Needed by GEOS to ingest cube sphere analysis)
744  call nccheck( nf90_def_var(self%ncid(n), "cubed_sphere", nf90_char, varid(vc)), &
745  "nf90_def_var cubed_sphere" )
746  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "grid_mapping_name", "gnomonic cubed-sphere"), &
747  "nf90_def_var time grid_mapping_name" )
748  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "file_format_version", "2.90"), &
749  "nf90_def_var time file_format_version" )
750  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "additional_vars", "contacts,orientation,anchor"), &
751  "nf90_def_var time additional_vars" )
752  call nccheck( nf90_put_att(self%ncid(n), varid(vc), "gridspec_file", "C"//trim(cubesize)//"_gridspec.nc4"), &
753  "nf90_def_var gridspec_file" )
754 
755  !vc=vc+1; !(Needed by GEOS to ingest cube sphere analysis)
756  !call nccheck( nf90_def_var(self%ncid(n), "ncontact", NF90_INT, varid(vc)), "nf90_def_var ncontact" )
757 
758 
759  ! End define mode
760  ! ---------------
761  call nccheck( nf90_enddef(self%ncid(n)), "nf90_enddef" )
762 
763 
764  ! Reset counter
765  ! -------------
766  vc=0
767 
768 
769  ! Write metadata arrays
770  ! ---------------------
771 
772  ! Tiles
773  if (self%tiledim(n)) then
774  vc=vc+1
775  call nccheck( nf90_put_var( self%ncid(n), varid(vc), tiles ), "nf90_put_var n" )
776  endif
777 
778  ! Xdim & Ydim arrays
779  vc=vc+1;call nccheck( nf90_put_var( self%ncid(n), varid(vc), xdimydim ), "nf90_put_var "//trim(xdimvar) )
780  vc=vc+1;call nccheck( nf90_put_var( self%ncid(n), varid(vc), xdimydim ), "nf90_put_var "//trim(ydimvar) )
781 
782  ! Start/counts
783  if (associated(istart)) nullify(istart)
784  if (associated(icount)) nullify(icount)
785  if (self%tiledim(n)) then
786  istart => self%is_r2_tile(1:3); icount => self%ic_r2_tile(1:3)
787  else
788  istart => self%is_r2_noti(1:2); icount => self%ic_r2_noti(1:2)
789  endif
790 
791  ! Lat/lon arrays
792  vc=vc+1;call nccheck( nf90_put_var( self%ncid(n), varid(vc), long, &
793  start = istart, &
794  count = icount ), &
795  "nf90_put_var lons" )
796 
797  vc=vc+1;call nccheck( nf90_put_var( self%ncid(n), varid(vc), latg, &
798  start = istart, &
799  count = icount ), &
800  "nf90_put_var lats" )
801 
802  ! Write model levels & time
803  vc=vc+1;call nccheck( nf90_put_var( self%ncid(n), varid(vc), levels(1:geom%npz) ), "nf90_put_var lev" )
804  vc=vc+1;call nccheck( nf90_put_var( self%ncid(n), varid(vc), levels ), "nf90_put_var edge" )
805 
806  ! Time
807  vc=vc+1;call nccheck( nf90_put_var( self%ncid(n), varid(vc), 0 ), "nf90_put_var time" )
808 
809  endif
810 
811  enddo
812 
813 endif
814 
815 end subroutine write_meta
816 
817 ! ------------------------------------------------------------------------------
818 
819 subroutine write_fields(self, geom, fields, vdate)
820 
821 implicit none
822 class(fv3jedi_io_geos), target, intent(inout) :: self
823 type(fv3jedi_geom), intent(inout) :: geom !< Geom
824 type(fv3jedi_field), intent(in) :: fields(:) !< Fields to be written
825 type(datetime), intent(in) :: vdate !< DateTime
826 
827 integer :: var, lev, n, ncid, filei, varid
828 integer, target :: dimids2_tile(4), dimids3_tile(5), dimidse_tile(5), dimids4_tile(5)
829 integer, target :: dimids2_noti(3), dimids3_noti(4), dimidse_noti(4), dimids4_noti(4)
830 integer, pointer :: dimids2(:), dimids3(:), dimidse(:), dimids4(:), dimids(:)
831 integer, allocatable, target :: is_r3_tile(:), is_r3_noti(:)
832 real(kind=kind_real), allocatable :: arrayg(:,:)
833 integer, pointer :: istart(:), icount(:)
834 
835 
836 ! Whole level of tile array
837 ! -------------------------
838 allocate(arrayg(1:geom%npx-1,1:geom%npy-1))
839 
840 
841 ! Local counts for 3 to allow changing start point
842 ! ------------------------------------------------
843 allocate(is_r3_tile(size(self%is_r3_tile)))
844 is_r3_tile = self%is_r3_tile
845 allocate(is_r3_noti(size(self%is_r3_noti)))
846 is_r3_noti = self%is_r3_noti
847 
848 
849 ! Dimension ID arrays for the various fields with and without tile dimension
850 ! --------------------------------------------------------------------------
851 dimids2_tile = (/ self%x_dimid, self%y_dimid, self%n_dimid, self%t_dimid /)
852 dimids3_tile = (/ self%x_dimid, self%y_dimid, self%n_dimid, self%z_dimid, self%t_dimid /)
853 dimidse_tile = (/ self%x_dimid, self%y_dimid, self%n_dimid, self%e_dimid, self%t_dimid /)
854 dimids4_tile = (/ self%x_dimid, self%y_dimid, self%n_dimid, self%f_dimid, self%t_dimid /)
855 
856 dimids2_noti = (/ self%x_dimid, self%y_dimid, self%t_dimid /)
857 dimids3_noti = (/ self%x_dimid, self%y_dimid, self%z_dimid, self%t_dimid /)
858 dimidse_noti = (/ self%x_dimid, self%y_dimid, self%e_dimid, self%t_dimid /)
859 dimids4_noti = (/ self%x_dimid, self%y_dimid, self%f_dimid, self%t_dimid /)
860 
861 
862 ! Loop over the fields
863 ! --------------------
864 do var = 1,size(fields)
865 
866  if (self%iam_io_proc) then
867 
868  ! ncid for this field
869  ! -------------------
870  filei = self%ncid_forfield(var)
871  ncid = self%ncid(filei)
872 
873  ! Redefine
874  ! --------
875  if (self%clobber) then
876 
877  call nccheck( nf90_redef(ncid), "nf90_enddef" )
878 
879  ! Dimension IDs for this field
880  ! ----------------------------
881  if (associated(dimids2)) nullify(dimids2)
882  if (associated(dimids3)) nullify(dimids3)
883  if (associated(dimidse)) nullify(dimidse)
884  if (associated(dimids4)) nullify(dimids4)
885 
886  if (self%tiledim(filei)) then
887  dimids2 => dimids2_tile
888  dimids3 => dimids3_tile
889  dimidse => dimidse_tile
890  dimids4 => dimids4_tile
891  else
892  dimids2 => dimids2_noti
893  dimids3 => dimids3_noti
894  dimidse => dimidse_noti
895  dimids4 => dimids4_noti
896  endif
897 
898  if (associated(dimids)) nullify (dimids)
899 
900  if (fields(var)%npz == 1) then
901  dimids => dimids2
902  elseif (fields(var)%npz == geom%npz) then
903  dimids => dimids3
904  elseif (fields(var)%npz == geom%npz+1) then
905  dimids => dimidse
906  elseif (fields(var)%npz == 4) then
907  dimids => dimids4
908  else
909  call abor1_ftn("write_geos: vertical dimension not supported")
910  endif
911 
912  ! Define field
913  call nccheck( nf90_def_var(ncid, trim(fields(var)%short_name), nf90_double, dimids, varid), &
914  "nf90_def_var"//trim(fields(var)%short_name))
915 
916  ! Write attributes if clobbering
917  if (self%clobber) then
918 
919  ! Long name and units
920  call nccheck( nf90_put_att(ncid, varid, "long_name" , trim(fields(var)%long_name) ), "nf90_put_att" )
921  call nccheck( nf90_put_att(ncid, varid, "units" , trim(fields(var)%units) ), "nf90_put_att" )
922 
923  ! Additional attributes for history and or plotting compatibility
924  if (.not.self%restart(filei)) then
925  call nccheck( nf90_put_att(ncid, varid, "standard_name", trim(fields(var)%long_name) ), "nf90_put_att" )
926  call nccheck( nf90_put_att(ncid, varid, "coordinates" , "lons lats" ), "nf90_put_att" )
927  call nccheck( nf90_put_att(ncid, varid, "grid_mapping" , "cubed_sphere" ), "nf90_put_att" )
928  endif
929 
930  endif
931 
932  ! End define mode
933  call nccheck( nf90_enddef(ncid), "nf90_enddef" )
934 
935  else
936 
937  ! Get existing variable id to write to
938  call nccheck ( nf90_inq_varid(ncid, trim(fields(var)%short_name), varid), &
939  "nf90_inq_varid "//trim(fields(var)%short_name) )
940 
941  endif
942 
943  ! Set starts and counts based on levels and tiledim flag
944  ! ------------------------------------------------------
945 
946  if (associated(istart)) nullify(istart)
947  if (associated(icount)) nullify(icount)
948 
949  if (fields(var)%npz == 1) then
950  if (self%tiledim(filei)) then
951  istart => self%is_r2_tile; icount => self%ic_r2_tile
952  else
953  istart => self%is_r2_noti; icount => self%ic_r2_noti
954  endif
955  elseif (fields(var)%npz > 1) then
956  if (self%tiledim(filei)) then
957  istart => is_r3_tile; icount => self%ic_r3_tile
958  else
959  istart => is_r3_noti; icount => self%ic_r3_noti
960  endif
961  endif
962 
963  endif
964 
965  do lev = 1,fields(var)%npz
966 
967  if (self%csize > 6) then
968  call gather_tile(geom, self%tcomm, 1, fields(var)%array(geom%isc:geom%iec,geom%jsc:geom%jec,lev), arrayg)
969  else
970  arrayg = fields(var)%array(geom%isc:geom%iec,geom%jsc:geom%jec,lev)
971  endif
972 
973  if (self%iam_io_proc) then
974 
975  is_r3_tile(self%vindex_tile) = lev
976  is_r3_noti(self%vindex_noti) = lev
977 
978  call nccheck( nf90_put_var( ncid, varid, arrayg, start = istart, count = icount ), &
979  "nf90_put_var "//trim(fields(var)%short_name) )
980 
981  endif
982 
983  enddo
984 
985 enddo
986 
987 ! Deallocate locals
988 ! -----------------
989 deallocate(is_r3_tile,is_r3_noti)
990 deallocate(arrayg)
991 
992 end subroutine write_fields
993 
994 ! ------------------------------------------------------------------------------
995 
996 subroutine gather_tile(geom, comm, nlev, array_l, array_g)
997 
998 implicit none
999 
1000 type(fv3jedi_geom), intent(in) :: geom
1001 integer, intent(in) :: comm
1002 integer, intent(in) :: nlev
1003 real(kind=kind_real), intent(in) :: array_l(geom%isc:geom%iec,geom%jsc:geom%jec,1:nlev) ! Local array
1004 real(kind=kind_real), intent(inout) :: array_g(1:geom%npx-1,1:geom%npy-1,1:nlev) ! Gathered array (only valid on root)
1005 
1006 real(kind=kind_real), allocatable :: vector_g(:), vector_l(:)
1007 integer :: comm_size, ierr
1008 integer :: ji, jj, jk, jc, n
1009 integer :: npx_g, npy_g, npx_l, npy_l
1010 integer, allocatable :: isc_l(:), iec_l(:), jsc_l(:), jec_l(:)
1011 integer, allocatable :: counts(:), displs(:), vectorcounts(:), vectordispls(:)
1012 
1013 !Get comm size
1014 call mpi_comm_size(comm, comm_size, ierr)
1015 
1016 !Array of counts and displacement
1017 allocate(counts(comm_size), displs(comm_size))
1018 do n = 1,comm_size
1019  displs(n) = n-1
1020  counts(n) = 1
1021 enddo
1022 
1023 !Horizontal size for global and local
1024 npx_g = geom%npx-1
1025 npy_g = geom%npy-1
1026 npx_l = geom%iec-geom%isc+1
1027 npy_l = geom%jec-geom%jsc+1
1028 
1029 !Gather local dimensions
1030 allocate(isc_l(comm_size), iec_l(comm_size), jsc_l(comm_size), jec_l(comm_size))
1031 call mpi_allgatherv(geom%isc, 1, mpi_int, isc_l, counts, displs, mpi_int, comm, ierr)
1032 call mpi_allgatherv(geom%iec, 1, mpi_int, iec_l, counts, displs, mpi_int, comm, ierr)
1033 call mpi_allgatherv(geom%jsc, 1, mpi_int, jsc_l, counts, displs, mpi_int, comm, ierr)
1034 call mpi_allgatherv(geom%jec, 1, mpi_int, jec_l, counts, displs, mpi_int, comm, ierr)
1035 deallocate(counts,displs)
1036 
1037 ! Pack whole tile array into vector
1038 allocate(vectorcounts(comm_size), vectordispls(comm_size))
1039 
1040 !Gather counts and displacement
1041 n = 0
1042 do jc = 1,comm_size
1043  vectordispls(jc) = n
1044  do jk = 1,nlev
1045  do jj = jsc_l(jc),jec_l(jc)
1046  do ji = isc_l(jc),iec_l(jc)
1047  n = n+1
1048  enddo
1049  enddo
1050  enddo
1051  vectorcounts(jc) = n - vectordispls(jc)
1052 enddo
1053 
1054 ! Pack local array into vector
1055 allocate(vector_l(npx_l*npy_l*nlev))
1056 n = 0
1057 do jk = 1,nlev
1058  do jj = geom%jsc,geom%jec
1059  do ji = geom%isc,geom%iec
1060  n = n+1
1061  vector_l(n) = array_l(ji,jj,jk)
1062  enddo
1063  enddo
1064 enddo
1065 
1066 ! Gather the full field
1067 allocate(vector_g(npx_g*npy_g*nlev))
1068 call mpi_gatherv( vector_l, npx_l*npy_l, mpi_double_precision, &
1069  vector_g, vectorcounts, vectordispls, mpi_double_precision, &
1070  0, comm, ierr)
1071 deallocate(vector_l,vectorcounts,vectordispls)
1072 
1073 !Unpack global vector into array
1074 n = 0
1075 do jc = 1,comm_size
1076  do jk = 1,nlev
1077  do jj = jsc_l(jc),jec_l(jc)
1078  do ji = isc_l(jc),iec_l(jc)
1079  n = n+1
1080  array_g(ji,jj,jk) = vector_g(n)
1081  enddo
1082  enddo
1083  enddo
1084 enddo
1085 deallocate(isc_l, iec_l, jsc_l, jec_l)
1086 
1087 deallocate(vector_g)
1088 
1089 end subroutine gather_tile
1090 
1091 ! ------------------------------------------------------------------------------
1092 
1093 subroutine scatter_tile(geom, comm, nlev, array_g, array_l)
1094 
1095 implicit none
1096 
1097 type(fv3jedi_geom), intent(in) :: geom
1098 integer, intent(in) :: comm
1099 integer, intent(in) :: nlev
1100 real(kind=kind_real), intent(in) :: array_g(1:geom%npx-1,1:geom%npy-1,nlev) ! Gathered array (only valid on root)
1101 real(kind=kind_real), intent(inout) :: array_l(geom%isc:geom%iec,geom%jsc:geom%jec,nlev) ! Local array
1102 
1103 real(kind=kind_real), allocatable :: vector_g(:), vector_l(:)
1104 integer :: comm_size, ierr
1105 integer :: ji, jj, jk, jc, n
1106 integer :: npx_g, npy_g, npx_l, npy_l
1107 integer, allocatable :: isc_l(:), iec_l(:), jsc_l(:), jec_l(:)
1108 integer, allocatable :: counts(:), displs(:), vectorcounts(:), vectordispls(:)
1109 
1110 !Get comm size
1111 call mpi_comm_size(comm, comm_size, ierr)
1112 
1113 !Array of counts and displacement
1114 allocate(counts(comm_size), displs(comm_size))
1115 do n = 1,comm_size
1116  displs(n) = n-1
1117  counts(n) = 1
1118 enddo
1119 
1120 !Horizontal size for global and local
1121 npx_g = geom%npx-1
1122 npy_g = geom%npy-1
1123 npx_l = geom%iec-geom%isc+1
1124 npy_l = geom%jec-geom%jsc+1
1125 
1126 !Gather local dimensions
1127 allocate(isc_l(comm_size), iec_l(comm_size), jsc_l(comm_size), jec_l(comm_size))
1128 call mpi_allgatherv(geom%isc, 1, mpi_int, isc_l, counts, displs, mpi_int, comm, ierr)
1129 call mpi_allgatherv(geom%iec, 1, mpi_int, iec_l, counts, displs, mpi_int, comm, ierr)
1130 call mpi_allgatherv(geom%jsc, 1, mpi_int, jsc_l, counts, displs, mpi_int, comm, ierr)
1131 call mpi_allgatherv(geom%jec, 1, mpi_int, jec_l, counts, displs, mpi_int, comm, ierr)
1132 deallocate(counts,displs)
1133 
1134 ! Pack whole tile array into vector
1135 allocate(vector_g(npx_g*npy_g*nlev))
1136 allocate(vectorcounts(comm_size), vectordispls(comm_size))
1137 
1138 n = 0
1139 do jc = 1,comm_size
1140  vectordispls(jc) = n
1141  do jk = 1,nlev
1142  do jj = jsc_l(jc),jec_l(jc)
1143  do ji = isc_l(jc),iec_l(jc)
1144  n = n+1
1145  vector_g(n) = array_g(ji,jj,jk)
1146  enddo
1147  enddo
1148  enddo
1149  vectorcounts(jc) = n - vectordispls(jc)
1150 enddo
1151 deallocate(isc_l, iec_l, jsc_l, jec_l)
1152 
1153 ! Scatter tile array to processors
1154 allocate(vector_l(npx_l*npy_l*nlev))
1155 
1156 call mpi_scatterv( vector_g, vectorcounts, vectordispls, mpi_double_precision, &
1157  vector_l, npx_l*npy_l, mpi_double_precision, &
1158  0, comm, ierr )
1159 
1160 deallocate(vector_g,vectorcounts,vectordispls)
1161 
1162 ! Unpack local vector into array
1163 n = 0
1164 do jk = 1,nlev
1165  do jj = geom%jsc,geom%jec
1166  do ji = geom%isc,geom%iec
1167  n = n+1
1168  array_l(ji,jj,jk) = vector_l(n)
1169  enddo
1170  enddo
1171 enddo
1172 deallocate(vector_l)
1173 
1174 end subroutine scatter_tile
1175 
1176 ! --------------------------------------------------------------------------------------------------
1177 
1178 ! Not really needed but prevents gnu compiler bug
1179 subroutine dummy_final(self)
1180 type(fv3jedi_io_geos), intent(inout) :: self
1181 end subroutine dummy_final
1182 
1183 ! --------------------------------------------------------------------------------------------------
1184 
1185 subroutine open_files(self)
1186 
1187 type(fv3jedi_io_geos), intent(inout) :: self
1188 
1189 integer :: n
1190 
1191 self%ncid = -1
1192 
1193 ! Open files for reading
1194 do n = 1,numfiles
1195  if (self%ncid_isneeded(n)) then
1196  call nccheck ( nf90_open( trim(self%datapath)//'/'//trim(self%filenames(n)), nf90_nowrite, &
1197  self%ncid(n) ), "nf90_open"//trim(self%filenames(n)) )
1198  endif
1199 enddo
1200 
1201 end subroutine open_files
1202 
1203 ! --------------------------------------------------------------------------------------------------
1204 
1205 subroutine create_files(self)
1206 
1207 type(fv3jedi_io_geos), intent(inout) :: self
1208 
1209 integer :: n, fileopts
1210 
1211 self%ncid = -1
1212 
1213 if (self%iam_io_proc) then
1214 
1215  if (self%clobber) then
1216 
1217  fileopts = ior(nf90_netcdf4, nf90_mpiio)
1218 
1219  ! Create/open files for writing
1220  do n = 1,numfiles
1221  if (self%ncid_isneeded(n)) then
1222  call nccheck( nf90_create( trim(self%datapath)//'/'//trim(self%filenames(n)), fileopts, &
1223  self%ncid(n), comm = self%ocomm, info = mpi_info_null), &
1224  "nf90_create"//trim(self%filenames(n)) )
1225  endif
1226  enddo
1227 
1228  else
1229 
1230  ! Open files for writing
1231  do n = 1,numfiles
1232 
1233  if (self%ncid_isneeded(n)) then
1234  call nccheck ( nf90_open( trim(self%datapath)//'/'//trim(self%filenames(n)), nf90_write, &
1235  self%ncid(n) ), "nf90_open"//trim(self%filenames(n)) )
1236  endif
1237  enddo
1238 
1239  endif
1240 
1241 endif
1242 
1243 end subroutine create_files
1244 
1245 ! --------------------------------------------------------------------------------------------------
1246 
1247 subroutine close_files(self)
1248 
1249 type(fv3jedi_io_geos), intent(inout) :: self
1250 
1251 integer :: n
1252 
1253 ! Close the files
1254 ! ---------------
1255 if (self%iam_io_proc) then
1256  do n = 1,numfiles
1257  if (self%ncid_isneeded(n)) then
1258  call nccheck ( nf90_close(self%ncid(n)), "nf90_close" )
1259  endif
1260  enddo
1261 endif
1262 
1263 end subroutine close_files
1264 
1265 ! --------------------------------------------------------------------------------------------------
1266 
1267 subroutine set_file_names(self, fields)
1268 
1269 type(fv3jedi_io_geos), intent(inout) :: self
1270 type(fv3jedi_field), intent(in) :: fields(:)
1271 
1272 integer :: var, grp, indfile
1273 
1274 if (allocated(self%ncid_forfield)) deallocate(self%ncid_forfield)
1275 allocate(self%ncid_forfield(size(fields)))
1276 self%ncid_forfield = -1
1277 self%ncid_isneeded = .false.
1278 
1279 ! Set files for the fields
1280 ! ------------------------
1281 do var = 1,size(fields)
1282 
1283  select case (trim(fields(var)%short_name))
1284 
1285  case("vtype","stype","vfrac") ! CRTM surface
1286  grp = 2
1287  case("U","V","W","PT","PKZ","PE","DZ") ! GEOS fv_core restart
1288  grp = 3
1289  case("Q","QILS","QICN","QLLS","QLCN","CLLS","CLCN") ! GEOS moist restart
1290  grp = 4
1291  case("PHIS") ! GEOS surface restart
1292  grp = 5
1293  case default ! Background file
1294  grp = 1
1295 
1296  endselect
1297 
1298  self%ncid_forfield(var) = grp
1299  self%ncid_isneeded(grp) = .true.
1300 
1301 enddo
1302 
1303 end subroutine set_file_names
1304 
1305 ! --------------------------------------------------------------------------------------------------
1306 
1307 end module fv3jedi_io_geos_mod
1308 
fv3jedi_io_utils_mod::replace_text
character(len(inputstr)+100) function replace_text(inputstr, search, replace)
Definition: fv3jedi_io_utils_mod.f90:93
fv3jedi_field_mod
Definition: fv3jedi_field_mod.f90:6
fv3jedi_io_geos_mod::set_file_names
subroutine set_file_names(self, fields)
Definition: fv3jedi_io_geos_mod.f90:1268
fv3jedi_netcdf_utils_mod
Definition: fv3jedi_netcdf_utils_mod.F90:6
fv3jedi_io_geos_mod::write_fields
subroutine write_fields(self, geom, fields, vdate)
Definition: fv3jedi_io_geos_mod.f90:820
fv3jedi_io_geos_mod::get_conf
subroutine get_conf(self, f_conf)
Definition: fv3jedi_io_geos_mod.f90:237
fv3jedi_io_geos_mod::delete
subroutine delete(self)
Definition: fv3jedi_io_geos_mod.f90:289
fv3jedi_io_geos_mod::close_files
subroutine close_files(self)
Definition: fv3jedi_io_geos_mod.f90:1248
fv3jedi_constants_mod::rad2deg
real(kind=kind_real), parameter, public rad2deg
Definition: fv3jedi_constants_mod.f90:13
fv3jedi_io_geos_mod::setup_date
subroutine setup_date(self, vdate)
Definition: fv3jedi_io_geos_mod.f90:198
fv3jedi_io_utils_mod
Definition: fv3jedi_io_utils_mod.f90:6
fv3jedi_geom_mod
Fortran module handling geometry for the FV3 model.
Definition: fv3jedi_geom_mod.f90:8
fv3jedi_io_geos_mod::dummy_final
subroutine dummy_final(self)
Definition: fv3jedi_io_geos_mod.f90:1180
fv3jedi_io_geos_mod::numfiles
integer, parameter numfiles
Definition: fv3jedi_io_geos_mod.f90:34
fv3jedi_geom_mod::fv3jedi_geom
Fortran derived type to hold geometry data for the FV3JEDI model.
Definition: fv3jedi_geom_mod.f90:46
fv3jedi_netcdf_utils_mod::nccheck
subroutine, public nccheck(status, iam)
Definition: fv3jedi_netcdf_utils_mod.F90:21
fv3jedi_constants_mod
Definition: fv3jedi_constants_mod.f90:6
fv3jedi_io_geos_mod::open_files
subroutine open_files(self)
Definition: fv3jedi_io_geos_mod.f90:1186
fv3jedi_io_geos_mod::gather_tile
subroutine gather_tile(geom, comm, nlev, array_l, array_g)
Definition: fv3jedi_io_geos_mod.f90:997
fv3jedi_io_geos_mod::scatter_tile
subroutine scatter_tile(geom, comm, nlev, array_g, array_l)
Definition: fv3jedi_io_geos_mod.f90:1094
fv3jedi_io_geos_mod::create_files
subroutine create_files(self)
Definition: fv3jedi_io_geos_mod.f90:1206
fv3jedi_kinds_mod::kind_real
integer, parameter, public kind_real
Definition: fv3jedi_kinds_mod.f90:14
fv3jedi_io_geos_mod
Definition: fv3jedi_io_geos_mod.f90:6
fv3jedi_io_utils_mod::vdate_to_datestring
subroutine vdate_to_datestring(vdate, datest, date, yyyy, mm, dd, hh, min, ss)
Definition: fv3jedi_io_utils_mod.f90:48
fv3jedi_io_geos_mod::write
subroutine write(self, geom, fields, vdate)
Definition: fv3jedi_io_geos_mod.f90:530
fv3jedi_field_mod::fv3jedi_field
Definition: fv3jedi_field_mod.f90:36
fv3jedi_io_geos_mod::read_meta
subroutine read_meta(self, geom, vdate, calendar_type, date_init, fields)
Definition: fv3jedi_io_geos_mod.f90:309
fv3jedi_io_geos_mod::fv3jedi_io_geos
Definition: fv3jedi_io_geos_mod.f90:36
fv3jedi_kinds_mod
Definition: fv3jedi_kinds_mod.f90:6
fv3jedi_io_utils_mod::string_from_conf
subroutine string_from_conf(f_conf, varstring, var, default, memberswap)
Definition: fv3jedi_io_utils_mod.f90:117
fv3jedi_io_geos_mod::write_meta
subroutine write_meta(self, geom, fields, vdate)
Definition: fv3jedi_io_geos_mod.f90:562
fv3jedi_io_geos_mod::setup_conf
subroutine setup_conf(self, geom, f_conf)
Definition: fv3jedi_io_geos_mod.f90:87
fv3jedi_io_geos_mod::read_fields
subroutine read_fields(self, geom, fields)
Definition: fv3jedi_io_geos_mod.f90:387