FMS  2025.04
Flexible Modeling System
fms_io_utils.F90
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 !> @defgroup fms_io_utils_mod fms_io_utils_mod
19 !> @ingroup fms2_io
20 !> @brief Misc. utility routines for use in @ref fms2_io
21 
22 !> @addtogroup fms_io_utils_mod
23 !> @{
24 module fms_io_utils_mod
25 use, intrinsic :: iso_fortran_env, only: error_unit
26 !use mpp_mod, only : get_ascii_file_num_lines_and_length, read_ascii_file
27 #ifdef _OPENMP
28 use omp_lib
29 #endif
30 use mpp_mod
31 use mpp_domains_mod, only: domain2d, domainug, mpp_get_ntile_count, &
33  mpp_get_ug_domain_ntiles, mpp_get_ug_domain_tile_id
34 use platform_mod
35 use fms_string_utils_mod, only: string_copy
36 implicit none
37 private
38 
39 character(len=32), save :: filename_appendix = '' !< Appendix added to the restart filename
40 
41 public :: char_linked_list
42 public :: error
43 public :: file_exists
44 public :: openmp_thread_trap
45 public :: string_copy
46 public :: is_in_list
47 public :: append_to_list
48 public :: destroy_list
51 public :: allocate_array
52 public :: put_array_section
53 public :: get_array_section
54 public :: get_data_type_string
55 public :: open_check
56 public :: string_compare
58 public :: ascii_read
59 public :: parse_mask_table
60 public :: get_mosaic_tile_file
61 public :: get_filename_appendix
62 public :: set_filename_appendix
63 public :: get_instance_filename
65 
66 !> @}
67 
68 !> @brief A linked list of strings
69 !> @ingroup fms_io_utils_mod
71  character(len=128) :: string
72  type(char_linked_list), pointer :: head => null()
73 endtype char_linked_list
74 
75 !> @ingroup fms_io_utils_mod
77  module procedure parse_mask_table_2d
78  module procedure parse_mask_table_3d
79 end interface parse_mask_table
80 
81 !> @ingroup fms_io_utils_mod
83  module procedure get_mosaic_tile_file_sg
84  module procedure get_mosaic_tile_file_ug
85 end interface get_mosaic_tile_file
86 
87 !> @ingroup fms_io_utils_mod
88 interface allocate_array
89  module procedure allocate_array_i4_kind_1d
90  module procedure allocate_array_i4_kind_2d
91  module procedure allocate_array_i4_kind_3d
92  module procedure allocate_array_i4_kind_4d
93  module procedure allocate_array_i4_kind_5d
94  module procedure allocate_array_i8_kind_1d
95  module procedure allocate_array_i8_kind_2d
96  module procedure allocate_array_i8_kind_3d
97  module procedure allocate_array_i8_kind_4d
98  module procedure allocate_array_i8_kind_5d
99  module procedure allocate_array_r4_kind_1d
100  module procedure allocate_array_r4_kind_2d
101  module procedure allocate_array_r4_kind_3d
102  module procedure allocate_array_r4_kind_4d
103  module procedure allocate_array_r4_kind_5d
104  module procedure allocate_array_r8_kind_1d
105  module procedure allocate_array_r8_kind_2d
106  module procedure allocate_array_r8_kind_3d
107  module procedure allocate_array_r8_kind_4d
108  module procedure allocate_array_r8_kind_5d
109  module procedure allocate_array_char_1d
110  module procedure allocate_array_char_2d
111  module procedure allocate_array_char_3d
112  module procedure allocate_array_char_4d
113  module procedure allocate_array_char_5d
114  module procedure allocate_array_char_6d
115 end interface allocate_array
116 
117 
118 !> @ingroup fms_io_utils_mod
120  module procedure put_array_section_i4_kind_1d
121  module procedure put_array_section_i4_kind_2d
122  module procedure put_array_section_i4_kind_3d
123  module procedure put_array_section_i4_kind_4d
124  module procedure put_array_section_i4_kind_5d
125  module procedure put_array_section_i8_kind_1d
126  module procedure put_array_section_i8_kind_2d
127  module procedure put_array_section_i8_kind_3d
128  module procedure put_array_section_i8_kind_4d
129  module procedure put_array_section_i8_kind_5d
130  module procedure put_array_section_r4_kind_1d
131  module procedure put_array_section_r4_kind_2d
132  module procedure put_array_section_r4_kind_3d
133  module procedure put_array_section_r4_kind_4d
134  module procedure put_array_section_r4_kind_5d
135  module procedure put_array_section_r8_kind_1d
136  module procedure put_array_section_r8_kind_2d
137  module procedure put_array_section_r8_kind_3d
138  module procedure put_array_section_r8_kind_4d
139  module procedure put_array_section_r8_kind_5d
140 end interface put_array_section
141 
142 
143 !> @ingroup fms_io_utils_mod
145  module procedure get_array_section_i4_kind_1d
146  module procedure get_array_section_i4_kind_2d
147  module procedure get_array_section_i4_kind_3d
148  module procedure get_array_section_i4_kind_4d
149  module procedure get_array_section_i4_kind_5d
150  module procedure get_array_section_i8_kind_1d
151  module procedure get_array_section_i8_kind_2d
152  module procedure get_array_section_i8_kind_3d
153  module procedure get_array_section_i8_kind_4d
154  module procedure get_array_section_i8_kind_5d
155  module procedure get_array_section_r4_kind_1d
156  module procedure get_array_section_r4_kind_2d
157  module procedure get_array_section_r4_kind_3d
158  module procedure get_array_section_r4_kind_4d
159  module procedure get_array_section_r4_kind_5d
160  module procedure get_array_section_r8_kind_1d
161  module procedure get_array_section_r8_kind_2d
162  module procedure get_array_section_r8_kind_3d
163  module procedure get_array_section_r8_kind_4d
164  module procedure get_array_section_r8_kind_5d
165 end interface get_array_section
166 
167 
168 !> @ingroup fms_io_utils_mod
170  module procedure get_data_type_string_0d
171  module procedure get_data_type_string_1d
172  module procedure get_data_type_string_2d
173  module procedure get_data_type_string_3d
174  module procedure get_data_type_string_4d
175  module procedure get_data_type_string_5d
176 end interface get_data_type_string
177 
178 !> @addtogroup fms_io_utils_mod
179 !> @{
180 contains
181 
182 
183 !> @brief Print a message to stderr, then stop the program.
184 subroutine error(mesg)
185 
186  character(len=*), intent(in) :: mesg !< Message that will be printed to
187  !! stderr.
188 
189  call mpp_error(fatal, trim(mesg))
190 end subroutine error
191 
192 
193 !> @brief Determine if a file exists.
194 !! @return Flag telling if the file exists.
195 function file_exists(path) &
196  result(exists)
197 
198  character(len=*), intent(in) :: path !< Path to file.
199  logical :: exists
200 
201 !$omp critical (file_existence_inquiry)
202  inquire(file=trim(path), exist=exists)
203 !$omp end critical (file_existence_inquiry)
204 end function file_exists
205 
206 
207 !> @brief Catch OpenMP parallel machines.
208 subroutine openmp_thread_trap()
209 
210 #ifdef _OPENMP
211  if (omp_get_level() .gt. 0) then
212  call error("this routine is not thread-safe. Please do not" &
213  //" call it in an OpenMP threaded region.")
214  endif
215 #endif
216 end subroutine openmp_thread_trap
217 
218 !> @brief Compare strings.
219 !! @return Flag telling if the strings are the same.
220 function string_compare(string1, string2, ignore_case) &
221  result(same)
222 
223  character(len=*), intent(in) :: string1
224  character(len=*), intent(in) :: string2
225  logical,intent(in), optional :: ignore_case
226  logical :: same
227 
228  if (len_trim(string1) .ne. len_trim(string2)) then
229  same = .false.
230  return
231  endif
232  if (present(ignore_case)) then
233  if (ignore_case) then
234  same = trim(lowercase(string1)) .eq. trim(lowercase(string2))
235  return
236  endif
237  endif
238  same = trim(string1) .eq. trim(string2)
239 end function string_compare
240 
241 
242 !> @brief Determine if a string exists in a character linked list.
243 !! @return Flag telling if the string is found in the list.
244 function is_in_list(list, string, ignore_case) &
245  result(in_list)
246 
247  type(char_linked_list), pointer, intent(in) :: list !< Linked list.
248  character(len=*), intent(in) :: string !< Input string.
249  logical, intent(in), optional :: ignore_case !< Flag to ignore case.
250  logical :: in_list
251 
252  type(char_linked_list), pointer :: p
253 
254  in_list = .false.
255  p => list
256  do while(associated(p))
257  if (string_compare(string, p%string, ignore_case)) then
258  in_list = .true.
259  return
260  else
261  p => p%head
262  endif
263  enddo
264 end function is_in_list
265 
266 
267 !> @brief Add node to character linked list.
268 subroutine append_to_list(list, string)
269  type(char_linked_list), pointer, intent(inout) :: list !< Linked list.
270  character(len=*), intent(in) :: string !< Input string.
271 
272  type(char_linked_list), pointer :: node
273  type(char_linked_list), pointer :: p
274 
275  allocate(node)
276  call string_copy(node%string, string)
277  node%head => null()
278  if (associated(list)) then
279  p => list
280  do while (associated(p%head))
281  p => p%head
282  enddo
283  p%head => node
284  else
285  list => node
286  endif
287 end subroutine append_to_list
288 
289 
290 !> @brief Deallocate all nodes on a character linked list.
291 subroutine destroy_list(list)
292 
293  type(char_linked_list), pointer, intent(inout) :: list !< Linked list.
294 
295  type(char_linked_list), pointer :: p
296  type(char_linked_list), pointer :: p2
297  p => list
298  do while (associated(p))
299  p2 => p%head
300  deallocate(p)
301  p => p2
302  enddo
303  list => null()
304 end subroutine destroy_list
305 
306 
307 !> @brief Determine if the "domain tile string" (.tilex.) exists in the input filename.
308 !! @internal
309 function has_domain_tile_string(string) &
310  result(has_string)
311 
312  character(len=*), intent(in) :: string !< Input string.
313  logical :: has_string
314 
315  integer :: l
316  integer :: i, j
317 
318  has_string = .false.
319 ! Assigns i to the index where ".tile" starts
320  i = index(trim(string), ".tile", back=.true.)
321  if (i .ne. 0) then
322  l = len_trim(string)
323 ! Sets i to the index after .tile
324  i = i + 5
325  j = i
326  do while (i .le. l)
327 ! If the ith characters is a dot but i not equal to the index after .tile set has_string to true
328  if (verify(string(i:i), ".") .eq. 0 .and. j .ne. i) then
329  has_string = .true.
330  exit
331 ! If the ith characters is NOT a number exit function and has_string will stay as false
332  elseif (verify(string(i:i), "0123456789") .ne. 0) then
333  exit
334  endif
335  i = i + 1
336  enddo
337  endif
338 end function has_domain_tile_string
339 
340 
341 !> @brief Add the domain tile id to an input filepath.
342 !! @internal
343 subroutine domain_tile_filepath_mangle(dest, source, domain_tile_id)
344 
345  character(len=*), intent(inout) :: dest !< Output filepath.
346  character(len=*), intent(in) :: source !< Input filepath.
347  integer, intent(in) :: domain_tile_id !< Domain tile id.
348 
349  integer :: i
350 
351  if (has_domain_tile_string(source)) then
352  call error("The file "//trim(source)//" has a domain tile id (tileX) added. Check your open_file call")
353  endif
354  i = index(trim(source), ".nc", back=.true.)
355  if (i .eq. 0) then
356  call error("The file "//trim(source)//" does not contain .nc. Check your open_file call")
357  endif
358  write(dest, '(a,i0,a)') source(1:i-1)//".tile", &
359  domain_tile_id, source(i:len_trim(source))
360 end subroutine domain_tile_filepath_mangle
361 
362 
363 !> @brief Determine if the "I/O domain tile string" (.nc.xxxx) exists in the input filename.
364 !! @internal
365 function has_io_domain_tile_string(string) &
366  result(has_string)
367 
368  character(len=*), intent(in) :: string !< Input string.
369  logical :: has_string
370 
371  integer :: i
372  integer :: l
373 
374  has_string = .false.
375  i = index(trim(string), ".nc.", back=.true.)
376  if (i .ne. 0) then
377  l = len_trim(string)
378  i = i + 1
379  do while (i .le. l)
380  if (verify(string(i:i), "0123456789") .ne. 0) then
381  return
382  endif
383  i = i + 1
384  enddo
385  has_string = .true.
386  endif
387 end function has_io_domain_tile_string
388 
389 
390 !> @brief Add the I/O domain tile id to an input filepath.
391 subroutine io_domain_tile_filepath_mangle(dest, source, io_domain_tile_id)
392 
393  character(len=*), intent(inout) :: dest !< Output filepath.
394  character(len=*), intent(in) :: source !< Input filepath.
395  integer, intent(in) :: io_domain_tile_id !< I/O domain tile id.
396 
397  if (has_io_domain_tile_string(source)) then
398  call error("The file "//trim(source)// &
399  & " has already had a domain tile id (.nc.XXXX) added. Check your open_file call.")
400  endif
401  write(dest,'(a,i4.4)') trim(source)//".", io_domain_tile_id
402 end subroutine io_domain_tile_filepath_mangle
403 
404 
405 !> @brief Determine if the "restart string" (.res.) exists in the input filename.
406 !! @internal
407 function has_restart_string(string) &
408  result(has_string)
409 
410  character(len=*), intent(in) :: string !< Input string.
411  logical :: has_string
412 
413  has_string = index(trim(string), ".res.", back=.true.) .ne. 0
414 end function has_restart_string
415 
416 
417 !> @brief Add ".res" to an input file path.
418 subroutine restart_filepath_mangle(dest, source)
419 
420  character(len=*), intent(inout) :: dest
421  character(len=*), intent(in) :: source
422 
423  integer :: i
424 
425  if (has_restart_string(source)) then
426  call string_copy(dest, source)
427  return
428  endif
429  if (has_domain_tile_string(source)) then
430  i = index(trim(source), ".tile", back=.true.)
431  else
432  i = index(trim(source), ".nc", back=.true.)
433  if (i .eq. 0) then
434  call error("The file "//trim(source)//" does not contain .nc. Check your open_file call")
435  endif
436  endif
437  call string_copy(dest, source(1:i-1)//".res"//source(i:len_trim(source)))
438 end subroutine restart_filepath_mangle
439 
440 subroutine open_check(flag, fname)
441 
442  logical, intent(in) :: flag
443  character(len=*), intent(in), optional :: fname !< The file name
444 
445  if (.not. flag) then
446  if (present(fname)) then
447  call mpp_error(fatal, "Error occurred while opening file "//trim(fname))
448  else
449  call mpp_error(fatal, "Error occurred while opening file.")
450  endif
451  endif
452 end subroutine open_check
453 
454 !> @brief Read the ascii text from filename `ascii_filename`into string array
455 !! `ascii_var`
456 subroutine ascii_read(ascii_filename, ascii_var, num_lines, max_length)
457  character(len=*), intent(in) :: ascii_filename !< The file name to be read
458  character(len=:), dimension(:), allocatable, intent(out) :: ascii_var !< The
459  !! string
460  !! array
461  integer, optional, intent(out) :: num_lines !< Optional argument to return number of lines in file
462  integer, optional, intent(out) :: max_length !< Optional argument to return max_length of line in file
463  integer, dimension(2) :: lines_and_length !< lines = 1, length = 2
464  if(allocated(ascii_var)) deallocate(ascii_var)
465  lines_and_length = get_ascii_file_num_lines_and_length(ascii_filename)
466  allocate(character(len=lines_and_length(2))::ascii_var(lines_and_length(1)))
467  call read_ascii_file(ascii_filename, lines_and_length(2), ascii_var)
468  if(present(num_lines)) num_lines = lines_and_length(1)
469  if(present(max_length)) max_length = lines_and_length(2)
470 end subroutine ascii_read
471 
472 !> @brief Populate 2D maskmap from mask_table given a model
473 subroutine parse_mask_table_2d(mask_table, maskmap, modelname)
474 
475  character(len=*), intent(in) :: mask_table !< Mask table to be read in
476  logical, intent(out) :: maskmap(:,:) !< 2D Mask output
477  character(len=*), intent(in) :: modelname !< Model to which this applies
478 
479  integer :: nmask, layout(2)
480  integer, allocatable :: mask_list(:,:)
481  character(len=:), dimension(:), allocatable :: mask_table_contents
482  integer :: iocheck, n, stdoutunit, offset
483  character(len=128) :: record
484 
485  maskmap = .true.
486  nmask = 0
487  stdoutunit = stdout()
488  call ascii_read(mask_table, mask_table_contents)
489  if( mpp_pe() == mpp_root_pe() ) then
490  read(mask_table_contents(1), fmt=*, iostat=iocheck) nmask
491  if (iocheck > 0) then
492  call mpp_error(fatal, "fms2_io(parse_mask_table_2d): Error in reading nmask from file variable")
493  elseif (iocheck < 0) then
494  call mpp_error(fatal, "fms2_io(parse_mask_table_2d): Error: nmask not completely read from file variable")
495  endif
496  write(stdoutunit,*)"parse_mask_table: Number of domain regions masked in ", trim(modelname), " = ", nmask
497  if( nmask > 0 ) then
498  !--- read layout from mask_table and confirm it matches the shape of maskmap
499  read(mask_table_contents(2), fmt=*, iostat=iocheck) layout
500  if (iocheck > 0) then
501  call mpp_error(fatal, "fms2_io(parse_mask_table_2d): Error in reading layout from file variable")
502  elseif (iocheck < 0) then
503  call mpp_error(fatal, "fms2_io(parse_mask_table_2d): Error: layout not completely read from file variable")
504  endif
505  if( (layout(1) .NE. size(maskmap,1)) .OR. (layout(2) .NE. size(maskmap,2)) )then
506  write(stdoutunit,*)"layout=", layout, ", size(maskmap) = ", size(maskmap,1), size(maskmap,2)
507  call mpp_error(fatal, "fms2_io(parse_mask_table_2d): layout in file "//trim(mask_table)// &
508  "does not match size of maskmap for "//trim(modelname))
509  endif
510  !--- make sure mpp_npes() == layout(1)*layout(2) - nmask
511  if( mpp_npes() .NE. layout(1)*layout(2) - nmask ) call mpp_error(fatal, &
512  .NE."fms2_io(parse_mask_table_2d): mpp_npes() layout(1)*layout(2) - nmask for "//trim(modelname))
513  endif
514  endif
515 
516  call mpp_broadcast(nmask, mpp_root_pe())
517 
518  if(nmask==0) return
519 
520  allocate(mask_list(nmask,2))
521 
522  if( mpp_pe() == mpp_root_pe() ) then
523  n = 0
524  offset = 3
525  do while (offset + n < size(mask_table_contents)+1)
526  read(mask_table_contents(n+offset),'(a)',iostat=iocheck) record
527  if (iocheck > 0) then
528  call mpp_error(fatal, "fms2_io(parse_mask_table_2d): Error in reading record from file variable")
529  elseif (iocheck < 0) then
530  call mpp_error(fatal, "fms2_io(parse_mask_table_2d): Error: record not completely read from file variable")
531  endif
532  if (record(1:1) == '#') then
533  offset = offset + 1
534  cycle
535  elseif (record(1:10) == ' ') then
536  offset = offset + 1
537  cycle
538  endif
539  n = n + 1
540  if( n > nmask ) then
541  call mpp_error(fatal, "fms2_io(parse_mask_table_2d): number of mask_list entry "// &
542  "is greater than nmask in file "//trim(mask_table) )
543  endif
544  read(record,*,iostat=iocheck) mask_list(n,1), mask_list(n,2)
545  if (iocheck > 0) then
546  call mpp_error(fatal, "fms2_io(parse_mask_table_2d): Error in reading mask_list from record variable")
547  elseif (iocheck < 0) then
548  call mpp_error(fatal, &
549  & "fms2_io(parse_mask_table_2d): Error: mask_list not completely read from record variable")
550  endif
551  enddo
552 
553  !--- make sure the number of entry for mask_list is nmask
554  if( n .NE. nmask) call mpp_error(fatal, &
555  "fms2_io(parse_mask_table_2d): number of mask_list entry does not match nmask in file "//trim(mask_table))
556  endif
557 
558  call mpp_broadcast(mask_list, 2*nmask, mpp_root_pe())
559  do n = 1, nmask
560  maskmap(mask_list(n,1),mask_list(n,2)) = .false.
561  enddo
562 
563  deallocate(mask_list, mask_table_contents)
564 
565 end subroutine parse_mask_table_2d
566 
567 
568 !> @brief Populate 3D maskmap from mask_table given a model
569 subroutine parse_mask_table_3d(mask_table, maskmap, modelname)
570 
571  character(len=*), intent(in) :: mask_table !< Mask table to be read in
572  logical, intent(out) :: maskmap(:,:,:) !< 2D Mask output
573  character(len=*), intent(in) :: modelname !< Model to which this applies
574 
575  integer :: nmask, layout(2)
576  integer, allocatable :: mask_list(:,:)
577  character(len=:), dimension(:), allocatable :: mask_table_contents
578  integer :: iocheck, n, stdoutunit, ntiles, offset
579  character(len=128) :: record
580 
581  maskmap = .true.
582  nmask = 0
583  stdoutunit = stdout()
584  call ascii_read(mask_table, mask_table_contents)
585  if( mpp_pe() == mpp_root_pe() ) then
586  read(mask_table_contents(1), fmt=*, iostat=iocheck) nmask
587  if (iocheck > 0) then
588  call mpp_error(fatal, "fms2_io(parse_mask_table_3d): Error in reading nmask from file variable")
589  elseif (iocheck < 0) then
590  call mpp_error(fatal, "fms2_io(parse_mask_table_3d): Error: nmask not completely read from file variable")
591  endif
592  write(stdoutunit,*)"parse_mask_table: Number of domain regions masked in ", trim(modelname), " = ", nmask
593  if( nmask > 0 ) then
594  !--- read layout from mask_table and confirm it matches the shape of maskmap
595  read(mask_table_contents(2), fmt=*, iostat=iocheck) layout(1), layout(2), ntiles
596  if (iocheck > 0) then
597  call mpp_error(fatal, "fms2_io(parse_mask_table_3d): Error in reading layout from file variable")
598  elseif (iocheck < 0) then
599  call mpp_error(fatal, "fms2_io(parse_mask_table_3d): Error: layout not completely read from file variable")
600  endif
601  if( (layout(1) .NE. size(maskmap,1)) .OR. (layout(2) .NE. size(maskmap,2)) )then
602  write(stdoutunit,*)"layout=", layout, ", size(maskmap) = ", size(maskmap,1), size(maskmap,2)
603  call mpp_error(fatal, "fms2_io(parse_mask_table_3d): layout in file "//trim(mask_table)// &
604  "does not match size of maskmap for "//trim(modelname))
605  endif
606  if( ntiles .NE. size(maskmap,3) ) then
607  write(stdoutunit,*)"ntiles=", ntiles, ", size(maskmap,3) = ", size(maskmap,3)
608  call mpp_error(fatal, "fms2_io(parse_mask_table_3d): ntiles in file "//trim(mask_table)// &
609  "does not match size of maskmap for "//trim(modelname))
610  endif
611  !--- make sure mpp_npes() == layout(1)*layout(2) - nmask
612  if( mpp_npes() .NE. layout(1)*layout(2)*ntiles - nmask ) then
613  print*, "layout=", layout, nmask, mpp_npes()
614  call mpp_error(fatal, &
615  .NE."fms2_io(parse_mask_table_3d): mpp_npes() layout(1)*layout(2) - nmask for "//trim(modelname))
616  endif
617  endif
618  endif
619 
620  call mpp_broadcast(nmask, mpp_root_pe())
621 
622  if(nmask==0) return
623 
624  allocate(mask_list(nmask,3))
625 
626  if( mpp_pe() == mpp_root_pe() ) then
627  n = 0
628  offset = 3
629  do while (offset + n < size(mask_table_contents)+1)
630  read(mask_table_contents(n+offset),'(a)',iostat=iocheck) record
631  if (iocheck > 0) then
632  call mpp_error(fatal, "fms2_io(parse_mask_table_3d): Error in reading record from file variable")
633  elseif (iocheck < 0) then
634  call mpp_error(fatal, "fms2_io(parse_mask_table_3d): Error: record not completely read from file variable")
635  endif
636  if (record(1:1) == '#') then
637  offset = offset + 1
638  cycle
639  elseif (record(1:10) == ' ') then
640  offset = offset + 1
641  cycle
642  endif
643  n = n + 1
644  if( n > nmask ) then
645  call mpp_error(fatal, "fms2_io(parse_mask_table_3d): number of mask_list entry "// &
646  "is greater than nmask in file "//trim(mask_table) )
647  endif
648  read(record,*,iostat=iocheck) mask_list(n,1), mask_list(n,2), mask_list(n,3)
649  if (iocheck > 0) then
650  call mpp_error(fatal, "fms2_io(parse_mask_table_3d): Error in reading mask_list from record variable")
651  elseif (iocheck < 0) then
652  call mpp_error(fatal, &
653  & "fms2_io(parse_mask_table_3d): Error: mask_list not completely read from record variable")
654  endif
655  enddo
656 
657  !--- make sure the number of entry for mask_list is nmask
658  if( n .NE. nmask) call mpp_error(fatal, &
659  "fms2_io(parse_mask_table_3d): number of mask_list entry does not match nmask in file "//trim(mask_table))
660 ! call mpp_close(unit)
661  endif
662 
663  call mpp_broadcast(mask_list, 3*nmask, mpp_root_pe())
664  do n = 1, nmask
665  maskmap(mask_list(n,1),mask_list(n,2),mask_list(n,3)) = .false.
666  enddo
667  deallocate(mask_list, mask_table_contents)
668 end subroutine parse_mask_table_3d
669 
670 !> @brief Determine tile_file for structured grid based on filename and current
671 !! tile on mpp_domain (this is mostly used for ongrid data_overrides)
672 subroutine get_mosaic_tile_file_sg(file_in, file_out, is_no_domain, domain, tile_count)
673  character(len=*), intent(in) :: file_in !< name of 'base' file
674  character(len=*), intent(out) :: file_out !< name of tile_file
675  logical, intent(in) :: is_no_domain !< are we providing a
676  !! domain
677  type(domain2d), intent(in), optional, target :: domain !< domain provided
678  integer, intent(in), optional :: tile_count !< tile count
679 
680  character(len=FMS_FILE_LEN) :: basefile
681  character(len=6) :: tilename
682  character(len=2) :: my_tile_str
683  integer :: lens, ntiles, ntileMe, tile, my_tile_id
684  integer, dimension(:), allocatable :: tile_id
685  type(domain2d), pointer, save :: d_ptr =>null()
686  logical :: domain_exist
687 
688  if(index(file_in, '.nc', back=.true.)==0) then
689  basefile = trim(file_in)
690  else
691  lens = len_trim(file_in)
692  if(file_in(lens-2:lens) .NE. '.nc') call mpp_error(fatal, &
693  'get_mosaic_tile_file_sg: .nc should be at the end of file '//trim(file_in))
694  basefile = file_in(1:lens-3)
695  end if
696 
697  !--- get the tile name
698  ntiles = 1
699  my_tile_id = 1
700  domain_exist = .false.
701  if(PRESENT(domain))then
702  domain_exist = .true.
703  ntiles = mpp_get_ntile_count(domain)
704  d_ptr => domain
705  endif
706 
707  if(domain_exist) then
708  ntileme = mpp_get_current_ntile(d_ptr)
709  allocate(tile_id(ntileme))
710  tile_id = mpp_get_tile_id(d_ptr)
711  tile = 1
712  if(present(tile_count)) tile = tile_count
713  my_tile_id = tile_id(tile)
714  endif
715 
716  if(ntiles > 1 .or. my_tile_id > 1 )then
717  write(my_tile_str, '(I0)') my_tile_id
718  tilename = 'tile'//trim(my_tile_str)
719  if(index(basefile,'.'//trim(tilename),back=.true.) == 0)then
720  basefile = trim(basefile)//'.'//trim(tilename);
721  end if
722  end if
723  if(allocated(tile_id)) deallocate(tile_id)
724 
725  file_out = trim(basefile)//'.nc'
726 
727  d_ptr =>null()
728 
729 end subroutine get_mosaic_tile_file_sg
730 
731 !> @brief Determine tile_file for unstructured grid based on filename and current
732 !! tile on mpp_domain (this is mostly used for ongrid data_overrides)
733 subroutine get_mosaic_tile_file_ug(file_in, file_out, domain)
734  character(len=*), intent(in) :: file_in !< name of base file
735  character(len=*), intent(out) :: file_out !< name of tile file
736  type(domainug), intent(in), optional :: domain !< domain provided
737 
738  character(len=FMS_FILE_LEN) :: basefile
739  character(len=6) :: tilename
740  character(len=2) :: my_tile_str
741  integer :: lens, ntiles, my_tile_id
742 
743  if(index(file_in, '.nc', back=.true.)==0) then
744  basefile = trim(file_in)
745  else
746  lens = len_trim(file_in)
747  if(file_in(lens-2:lens) .NE. '.nc') call mpp_error(fatal, &
748  'fms_io_mod: .nc should be at the end of file '//trim(file_in))
749  basefile = file_in(1:lens-3)
750  end if
751 
752  !--- get the tile name
753  ntiles = 1
754  my_tile_id = 1
755  if(PRESENT(domain))then
756  ntiles = mpp_get_ug_domain_ntiles(domain)
757  my_tile_id = mpp_get_ug_domain_tile_id(domain)
758  endif
759 
760  if(ntiles > 1 .or. my_tile_id > 1 )then
761  write(my_tile_str, '(I0)') my_tile_id
762  tilename = 'tile'//trim(my_tile_str)
763  if(index(basefile,'.'//trim(tilename),back=.true.) == 0)then
764  basefile = trim(basefile)//'.'//trim(tilename);
765  end if
766  end if
767 
768  file_out = trim(basefile)//'.nc'
769 
770 end subroutine get_mosaic_tile_file_ug
771 
772 !> @brief Writes filename appendix to "string_out"
773 subroutine get_filename_appendix(string_out)
774  character(len=*) , intent(out) :: string_out !< String to write the filename_appendix to
775 
776  string_out = trim(filename_appendix)
777 
778 end subroutine get_filename_appendix
779 
780 !> @brief Clears the filename_appendix module variable
782 
783  filename_appendix = ''
784 
785 end subroutine nullify_filename_appendix
786 
787 !> @brief Save "string_in" as a module variable that will added to the filename
788 !! of the restart files
789 subroutine set_filename_appendix(string_in)
790  character(len=*) , intent(in) :: string_in !< String that will be saved as a module variable
791 
792  ! Check if string has already been added
793  if (len_trim(filename_appendix) > 0) then
794  call error("Set_filename_appendix: The filename appendix has already be set " &
795  //"call 'nullify_filename_appendix' first")
796  endif
797 
798  filename_appendix = trim(string_in)
799 
800 end subroutine set_filename_appendix
801 
802 !> @brief Adds the filename_appendix to name_in and sets it as name_out
803 subroutine get_instance_filename(name_in,name_out)
804  character(len=*) , intent(in) :: name_in !< Buffer to add the filename_appendix to
805  character(len=*), intent(inout) :: name_out !< name_in with the filename_appendix
806 
807  integer :: length !< Length of name_in
808  integer :: i !< no description
809 
810  length = len_trim(name_in)
811  name_out = name_in(1:length)
812 
813  if(len_trim(filename_appendix) > 0) then
814  !< If .tileXX is in the filename add the appendix before it
815  if (has_domain_tile_string(name_in)) then
816  i = index(trim(name_in), ".tile", back=.true.)
817  name_out = name_in(1:i-1) //'.'//trim(filename_appendix)//name_in(i:length)
818  return
819  endif
820 
821  !< If .nc is in the filename add the appendix before it
822  i = index(trim(name_in), ".nc", back=.true.)
823  if ( i .ne. 0 ) then
824  name_out = name_in(1:i-1)//'.'//trim(filename_appendix)//name_in(i:length)
825  else
826  i = index(trim(name_in), ".yaml", back=.true.)
827  if (i .ne. 0) then
828  !< If .yaml is in the filename add the appendix before it
829  name_out = name_in(1:i-1)//'.'//trim(filename_appendix)//name_in(i:length)
830  else
831  !< If .nc and .yaml are not in the name, add the appendix at the end of the file
832  name_out = name_in(1:length) //'.'//trim(filename_appendix)
833  endif
834  end if
835  end if
836 
837 end subroutine get_instance_filename
838 
839 include "array_utils.inc"
840 include "array_utils_char.inc"
841 include "get_data_type_string.inc"
842 
843 
844 end module fms_io_utils_mod
845 !> @}
846 ! close documentation grouping
subroutine allocate_array_r8_kind_1d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine get_data_type_string_4d(sdata, type_string)
Return a string describing the input data's type.
subroutine put_array_section_i8_kind_2d(section, array, s, c)
Put a section of an array into a larger array.
subroutine put_array_section_r8_kind_4d(section, array, s, c)
Put a section of an array into a larger array.
subroutine allocate_array_char_5d(buf, sizes, initialize)
Allocate character arrays using an input array of sizes.
subroutine put_array_section_i4_kind_2d(section, array, s, c)
Put a section of an array into a larger array.
Definition: array_utils.inc:76
subroutine put_array_section_i4_kind_1d(section, array, s, c)
Put a section of an array into a larger array.
Definition: array_utils.inc:39
subroutine, public get_instance_filename(name_in, name_out)
Adds the filename_appendix to name_in and sets it as name_out.
subroutine allocate_array_char_3d(buf, sizes, initialize)
Allocate character arrays using an input array of sizes.
subroutine get_array_section_i8_kind_4d(section, array, s, c)
Get a section of larger array.
logical function has_io_domain_tile_string(string)
Determine if the "I/O domain tile string" (.nc.xxxx) exists in the input filename.
subroutine get_array_section_r4_kind_3d(section, array, s, c)
Get a section of larger array.
subroutine, public ascii_read(ascii_filename, ascii_var, num_lines, max_length)
Read the ascii text from filename ascii_filenameinto string array ascii_var
subroutine get_mosaic_tile_file_ug(file_in, file_out, domain)
Determine tile_file for unstructured grid based on filename and current tile on mpp_domain (this is m...
subroutine get_array_section_r8_kind_1d(section, array, s, c)
Get a section of larger array.
subroutine get_array_section_i8_kind_3d(section, array, s, c)
Get a section of larger array.
subroutine allocate_array_i8_kind_2d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine put_array_section_i4_kind_3d(section, array, s, c)
Put a section of an array into a larger array.
subroutine, public restart_filepath_mangle(dest, source)
Add ".res" to an input file path.
subroutine put_array_section_r8_kind_2d(section, array, s, c)
Put a section of an array into a larger array.
subroutine put_array_section_r8_kind_3d(section, array, s, c)
Put a section of an array into a larger array.
subroutine get_array_section_i8_kind_5d(section, array, s, c)
Get a section of larger array.
subroutine get_data_type_string_5d(sdata, type_string)
Return a string describing the input data's type.
subroutine get_array_section_r4_kind_1d(section, array, s, c)
Get a section of larger array.
logical function, public file_exists(path)
Determine if a file exists.
subroutine put_array_section_i8_kind_1d(section, array, s, c)
Put a section of an array into a larger array.
subroutine put_array_section_i8_kind_5d(section, array, s, c)
Put a section of an array into a larger array.
subroutine allocate_array_r4_kind_3d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine allocate_array_i8_kind_1d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine parse_mask_table_2d(mask_table, maskmap, modelname)
Populate 2D maskmap from mask_table given a model.
subroutine allocate_array_r8_kind_4d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine allocate_array_r8_kind_3d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine, public open_check(flag, fname)
logical function has_restart_string(string)
Determine if the "restart string" (.res.) exists in the input filename.
subroutine, public nullify_filename_appendix()
Clears the filename_appendix module variable.
subroutine, public io_domain_tile_filepath_mangle(dest, source, io_domain_tile_id)
Add the I/O domain tile id to an input filepath.
subroutine allocate_array_i8_kind_5d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine allocate_array_i4_kind_2d(buf, sizes)
Allocate arrays using an input array of sizes.
Definition: array_utils.inc:63
subroutine, public get_filename_appendix(string_out)
Writes filename appendix to "string_out".
subroutine allocate_array_r4_kind_1d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine, public set_filename_appendix(string_in)
Save "string_in" as a module variable that will added to the filename of the restart files.
subroutine get_array_section_r8_kind_5d(section, array, s, c)
Get a section of larger array.
logical function, public string_compare(string1, string2, ignore_case)
Compare strings.
subroutine get_data_type_string_1d(sdata, type_string)
Return a string describing the input data's type.
subroutine put_array_section_r4_kind_5d(section, array, s, c)
Put a section of an array into a larger array.
subroutine allocate_array_char_4d(buf, sizes, initialize)
Allocate character arrays using an input array of sizes.
subroutine allocate_array_r8_kind_5d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine allocate_array_r4_kind_4d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine get_array_section_i4_kind_5d(section, array, s, c)
Get a section of larger array.
subroutine allocate_array_char_6d(buf, sizes, initialize)
Allocate character arrays using an input array of sizes.
subroutine get_array_section_i4_kind_3d(section, array, s, c)
Get a section of larger array.
subroutine allocate_array_r4_kind_2d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine allocate_array_i4_kind_4d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine put_array_section_i8_kind_4d(section, array, s, c)
Put a section of an array into a larger array.
subroutine allocate_array_i4_kind_3d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine put_array_section_r4_kind_3d(section, array, s, c)
Put a section of an array into a larger array.
subroutine parse_mask_table_3d(mask_table, maskmap, modelname)
Populate 3D maskmap from mask_table given a model.
subroutine get_mosaic_tile_file_sg(file_in, file_out, is_no_domain, domain, tile_count)
Determine tile_file for structured grid based on filename and current tile on mpp_domain (this is mos...
subroutine get_data_type_string_0d(sdata, type_string)
Return a string describing the input data's type.
subroutine get_array_section_i8_kind_1d(section, array, s, c)
Get a section of larger array.
subroutine get_array_section_r4_kind_4d(section, array, s, c)
Get a section of larger array.
subroutine put_array_section_r8_kind_1d(section, array, s, c)
Put a section of an array into a larger array.
subroutine get_array_section_i4_kind_2d(section, array, s, c)
Get a section of larger array.
Definition: array_utils.inc:88
subroutine get_data_type_string_3d(sdata, type_string)
Return a string describing the input data's type.
subroutine get_array_section_r8_kind_2d(section, array, s, c)
Get a section of larger array.
subroutine, public destroy_list(list)
Deallocate all nodes on a character linked list.
subroutine get_array_section_r8_kind_4d(section, array, s, c)
Get a section of larger array.
logical function has_domain_tile_string(string)
Determine if the "domain tile string" (.tilex.) exists in the input filename.
subroutine allocate_array_i8_kind_3d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine get_data_type_string_2d(sdata, type_string)
Return a string describing the input data's type.
subroutine put_array_section_i4_kind_5d(section, array, s, c)
Put a section of an array into a larger array.
subroutine get_array_section_i4_kind_4d(section, array, s, c)
Get a section of larger array.
subroutine, public append_to_list(list, string)
Add node to character linked list.
subroutine put_array_section_i4_kind_4d(section, array, s, c)
Put a section of an array into a larger array.
subroutine, public openmp_thread_trap()
Catch OpenMP parallel machines.
subroutine allocate_array_i4_kind_5d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine allocate_array_char_1d(buf, sizes, initialize)
Allocate character arrays using an input array of sizes.
subroutine allocate_array_r4_kind_5d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine allocate_array_i4_kind_1d(buf, sizes)
Allocate arrays using an input array of sizes.
Definition: array_utils.inc:26
subroutine put_array_section_r8_kind_5d(section, array, s, c)
Put a section of an array into a larger array.
subroutine get_array_section_i4_kind_1d(section, array, s, c)
Get a section of larger array.
Definition: array_utils.inc:51
subroutine put_array_section_i8_kind_3d(section, array, s, c)
Put a section of an array into a larger array.
subroutine get_array_section_i8_kind_2d(section, array, s, c)
Get a section of larger array.
subroutine put_array_section_r4_kind_1d(section, array, s, c)
Put a section of an array into a larger array.
subroutine allocate_array_i8_kind_4d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine put_array_section_r4_kind_4d(section, array, s, c)
Put a section of an array into a larger array.
subroutine allocate_array_r8_kind_2d(buf, sizes)
Allocate arrays using an input array of sizes.
subroutine allocate_array_char_2d(buf, sizes, initialize)
Allocate character arrays using an input array of sizes.
subroutine, public domain_tile_filepath_mangle(dest, source, domain_tile_id)
Add the domain tile id to an input filepath.
logical function, public is_in_list(list, string, ignore_case)
Determine if a string exists in a character linked list.
subroutine get_array_section_r4_kind_2d(section, array, s, c)
Get a section of larger array.
subroutine get_array_section_r8_kind_3d(section, array, s, c)
Get a section of larger array.
subroutine get_array_section_r4_kind_5d(section, array, s, c)
Get a section of larger array.
subroutine put_array_section_r4_kind_2d(section, array, s, c)
Put a section of an array into a larger array.
A linked list of strings.
subroutine, public string_copy(dest, source, check_for_null)
Safely copy a string from one buffer to another.
integer function, dimension(size(domain%tile_id(:))) mpp_get_tile_id(domain)
Returns the tile_id on current pe.
integer function mpp_get_current_ntile(domain)
Returns number of tile on current pe.
integer function mpp_get_ntile_count(domain)
Returns number of tiles in mosaic.
The domain2D type contains all the necessary information to define the global, compute and data domai...
Domain information for managing data on unstructured grids.
integer function stdout()
This function returns the current standard fortran unit numbers for output.
Definition: mpp_util.inc:42
subroutine read_ascii_file(FILENAME, LENGTH, Content, PELIST)
Reads any ascii file into a character array and broadcasts it to the non-root mpi-tasks....
Definition: mpp_util.inc:1446
integer function mpp_npes()
Returns processor count for current pelist.
Definition: mpp_util.inc:420
integer function, dimension(2) get_ascii_file_num_lines_and_length(FILENAME, PELIST)
Function to determine the maximum line length and number of lines from an ascii file.
Definition: mpp_util.inc:1353
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:406
Perform parallel broadcasts.
Definition: mpp.F90:1104
Error handler.
Definition: mpp.F90:381