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